Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - kernel/none - mp.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 19825-b77c7f8) Lines: 1106 1148 96.3 %
Date: 2016-12-05 05:49:04 Functions: 69 69 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : #line 2 "../src/kernel/none/mp.c"
       2             : /* Copyright (C) 2000-2003 The PARI group.
       3             : 
       4             : This file is part of the PARI/GP package.
       5             : 
       6             : PARI/GP is free software; you can redistribute it and/or modify it under the
       7             : terms of the GNU General Public License as published by the Free Software
       8             : Foundation. It is distributed in the hope that it will be useful, but WITHOUT
       9             : ANY WARRANTY WHATSOEVER.
      10             : 
      11             : Check the License for details. You should have received a copy of it, along
      12             : with the package; see the file 'COPYING'. If not, write to the Free Software
      13             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      14             : 
      15             : /***********************************************************************/
      16             : /**                                                                   **/
      17             : /**                       MULTIPRECISION KERNEL                       **/
      18             : /**                                                                   **/
      19             : /***********************************************************************/
      20             : #include "pari.h"
      21             : #include "paripriv.h"
      22             : #include "../src/kernel/none/tune-gen.h"
      23             : 
      24         574 : void pari_kernel_init(void) { }
      25         573 : void pari_kernel_close(void) { }
      26             : 
      27             : /* NOTE: arguments of "spec" routines (muliispec, addiispec, etc.) aren't
      28             :  * GENs but pairs (long *a, long na) representing a list of digits (in basis
      29             :  * BITS_IN_LONG) : a[0], ..., a[na-1]. [ In ordre to facilitate splitting: no
      30             :  * need to reintroduce codewords ] */
      31             : 
      32             : #define LIMBS(x)  ((x)+2)
      33             : #define NLIMBS(x) (lgefint(x)-2)
      34             : 
      35             : /* Normalize a non-negative integer */
      36             : GEN
      37   189174705 : int_normalize(GEN x, long known_zero_words)
      38             : {
      39   189174705 :   long i, lx = lgefint(x);
      40             :   GEN x0;
      41   189174705 :   if (lx == 2) { x[1] = evalsigne(0) | evallgefint(2); return x; }
      42   189174705 :   if (!known_zero_words && x[2]) return x;
      43   774279222 :   for (i = 2+known_zero_words; i < lx; i++)
      44   726813606 :     if (x[i]) break;
      45   117292977 :   x0 = x; i -= 2; x += i;
      46   117292977 :   if (x0 == (GEN)avma) avma = (pari_sp)x;
      47    48215238 :   else stackdummy((pari_sp)(x0+i), (pari_sp)x0);
      48   117292977 :   lx -= i;
      49   117292977 :   x[0] = evaltyp(t_INT) | evallg(lx);
      50   117292977 :   if (lx == 2) x[1] = evalsigne(0) | evallgefint(lx);
      51    69827361 :   else         x[1] = evalsigne(1) | evallgefint(lx);
      52   117292977 :   return x;
      53             : }
      54             : 
      55             : /***********************************************************************/
      56             : /**                                                                   **/
      57             : /**                      ADDITION / SUBTRACTION                       **/
      58             : /**                                                                   **/
      59             : /***********************************************************************/
      60             : 
      61             : GEN
      62     1809246 : setloop(GEN a)
      63             : {
      64     1809246 :   pari_sp av = avma;
      65     1809246 :   (void)cgetg(lgefint(a) + 3, t_VECSMALL);
      66     1809246 :   return icopy_avma(a, av); /* two cells of extra space before a */
      67             : }
      68             : 
      69             : /* we had a = setloop(?), then some incloops. Reset a to b */
      70             : GEN
      71      117531 : resetloop(GEN a, GEN b) {
      72      117531 :   long lb = lgefint(b);
      73      117531 :   a += lgefint(a) - lb;
      74      117531 :   a[0] = evaltyp(t_INT) | evallg(lb);
      75      117531 :   affii(b, a); return a;
      76             : }
      77             : 
      78             : /* assume a > 0, initialized by setloop. Do a++ */
      79             : static GEN
      80    19084320 : incpos(GEN a)
      81             : {
      82    19084320 :   long i, l = lgefint(a);
      83    19084323 :   for (i=l-1; i>1; i--)
      84    19084320 :     if (++a[i]) return a;
      85           3 :   l++; a--; /* use extra cell */
      86           3 :   a[0]=evaltyp(t_INT) | _evallg(l);
      87           3 :   a[1]=evalsigne(1) | evallgefint(l);
      88           3 :   a[2]=1; return a;
      89             : }
      90             : 
      91             : /* assume a < 0, initialized by setloop. Do a++ */
      92             : static GEN
      93        4638 : incneg(GEN a)
      94             : {
      95        4638 :   long i, l = lgefint(a)-1;
      96        4638 :   if (uel(a,l)--)
      97             :   {
      98        4635 :     if (l == 2 && !a[2])
      99             :     {
     100         399 :       a++; /* save one cell */
     101         399 :       a[0] = evaltyp(t_INT) | _evallg(2);
     102         399 :       a[1] = evalsigne(0) | evallgefint(2);
     103             :     }
     104        4635 :     return a;
     105             :   }
     106           3 :   for (i = l-1;; i--) /* finishes since a[2] != 0 */
     107           3 :     if (uel(a,i)--) break;
     108           0 :   if (!a[2])
     109             :   {
     110           3 :     a++; /* save one cell */
     111           3 :     a[0] = evaltyp(t_INT) | _evallg(l);
     112           3 :     a[1] = evalsigne(-1) | evallgefint(l);
     113             :   }
     114           3 :   return a;
     115             : }
     116             : 
     117             : /* assume a initialized by setloop. Do a++ */
     118             : GEN
     119    19104489 : incloop(GEN a)
     120             : {
     121    19104489 :   switch(signe(a))
     122             :   {
     123       15531 :     case 0: a--; /* use extra cell */
     124       15531 :       a[0]=evaltyp(t_INT) | _evallg(3);
     125       15531 :       a[1]=evalsigne(1) | evallgefint(3);
     126       15531 :       a[2]=1; return a;
     127        4638 :     case -1: return incneg(a);
     128    19084320 :     default: return incpos(a);
     129             :   }
     130             : }
     131             : 
     132             : INLINE GEN
     133   834130366 : adduispec(ulong s, GEN x, long nx)
     134             : {
     135   834130366 :   GEN xd, zd = (GEN)avma;
     136             :   long lz;
     137             : 
     138   834130366 :   if (nx == 1) return adduu(s, uel(x,0));
     139   184801350 :   lz = nx+3; (void)new_chunk(lz);
     140   184801350 :   xd = x + nx;
     141   184801350 :   *--zd = (ulong)*--xd + s;
     142   184801350 :   if ((ulong)*zd < s)
     143             :     for(;;)
     144             :     {
     145    15941037 :       if (xd == x) { *--zd = 1; break; } /* enlarge z */
     146    12347511 :       *--zd = ((ulong)*--xd) + 1;
     147    12347511 :       if (*zd) { lz--; break; }
     148     5678862 :     }
     149   174539175 :   else lz--;
     150   184801350 :   while (xd > x) *--zd = *--xd;
     151   184801350 :   *--zd = evalsigne(1) | evallgefint(lz);
     152   184801350 :   *--zd = evaltyp(t_INT) | evallg(lz);
     153   184801350 :   avma=(pari_sp)zd; return zd;
     154             : }
     155             : 
     156             : GEN
     157   211825614 : adduispec_offset(ulong s, GEN x, long offset, long nx)
     158             : {
     159   211825614 :   GEN xd = x+lgefint(x)-nx-offset;
     160   211825614 :   while (nx && *xd==0) {xd++; nx--;}
     161   211825614 :   if (!nx) return utoi(s);
     162   198913578 :   return adduispec(s,xd,nx);
     163             : }
     164             : 
     165             : static GEN
     166  1683976260 : addiispec(GEN x, GEN y, long nx, long ny)
     167             : {
     168             :   GEN xd, yd, zd;
     169  1683976260 :   long lz, i = -2;
     170             :   LOCAL_OVERFLOW;
     171             : 
     172  1683976260 :   if (nx < ny) swapspec(x,y, nx,ny);
     173  1683976260 :   if (ny == 1) return adduispec(*y,x,nx);
     174  1075872126 :   zd = (GEN)avma;
     175  1075872126 :   lz = nx+3; (void)new_chunk(lz);
     176  1075872126 :   xd = x + nx;
     177  1075872126 :   yd = y + ny;
     178  1075872126 :   zd[-1] = addll(xd[-1], yd[-1]);
     179             : #ifdef addllx8
     180   623718894 :   for (  ; i-8 > -ny; i-=8)
     181   265094852 :     addllx8(xd+i, yd+i, zd+i, overflow);
     182             : #endif
     183  1075872126 :   for (  ; i >= -ny; i--) zd[i] = addllx(xd[i], yd[i]);
     184  1075872126 :   if (overflow)
     185             :     for(;;)
     186             :     {
     187   104123625 :       if (i < -nx) { zd[i] = 1; i--; break; } /* enlarge z */
     188    77236098 :       zd[i] = uel(xd,i) + 1;
     189    77236098 :       if (zd[i]) { i--; lz--; break; }
     190    26061699 :       i--;
     191    26061699 :     }
     192   997810200 :   else lz--;
     193  1075872126 :   for (; i >= -nx; i--) zd[i] = xd[i];
     194  1075872126 :   zd += i+1;
     195  1075872126 :   *--zd = evalsigne(1) | evallgefint(lz);
     196  1075872126 :   *--zd = evaltyp(t_INT) | evallg(lz);
     197  1075872126 :   avma=(pari_sp)zd; return zd;
     198             : }
     199             : 
     200             : /* assume x >= s */
     201             : INLINE GEN
     202   672601302 : subiuspec(GEN x, ulong s, long nx)
     203             : {
     204   672601302 :   GEN xd, zd = (GEN)avma;
     205             :   long lz;
     206             :   LOCAL_OVERFLOW;
     207             : 
     208   672601302 :   if (nx == 1) return utoi(x[0] - s);
     209             : 
     210    74248852 :   lz = nx+2; (void)new_chunk(lz);
     211    74248852 :   xd = x + nx;
     212    74248852 :   *--zd = subll(*--xd, s);
     213    74248852 :   if (overflow)
     214             :     for(;;)
     215             :     {
     216    34005402 :       *--zd = ((ulong)*--xd) - 1;
     217    34005402 :       if (*xd) break;
     218     1475043 :     }
     219    74248852 :   if (xd == x)
     220    32119182 :     while (*zd == 0) { zd++; lz--; } /* shorten z */
     221             :   else
     222   298556546 :     do  *--zd = *--xd; while (xd > x);
     223    74248852 :   *--zd = evalsigne(1) | evallgefint(lz);
     224    74248852 :   *--zd = evaltyp(t_INT) | evallg(lz);
     225    74248852 :   avma=(pari_sp)zd; return zd;
     226             : }
     227             : 
     228             : /* assume x > y */
     229             : static GEN
     230  1358473487 : subiispec(GEN x, GEN y, long nx, long ny)
     231             : {
     232             :   GEN xd,yd,zd;
     233  1358473487 :   long lz, i = -2;
     234             :   LOCAL_OVERFLOW;
     235             : 
     236  1358473487 :   if (ny==1) return subiuspec(x,*y,nx);
     237   701450300 :   zd = (GEN)avma;
     238   701450300 :   lz = nx+2; (void)new_chunk(lz);
     239   701450300 :   xd = x + nx;
     240   701450300 :   yd = y + ny;
     241   701450300 :   zd[-1] = subll(xd[-1], yd[-1]);
     242             : #ifdef subllx8
     243   420389393 :   for (  ; i-8 > -ny; i-=8)
     244   186572626 :     subllx8(xd+i, yd+i, zd+i, overflow);
     245             : #endif
     246   701450300 :   for (  ; i >= -ny; i--) zd[i] = subllx(xd[i], yd[i]);
     247   701450300 :   if (overflow)
     248             :     for(;;)
     249             :     {
     250    90150638 :       zd[i] = uel(xd,i) - 1;
     251    90150638 :       if (xd[i--]) break;
     252    46955250 :     }
     253   701450300 :   if (i>=-nx)
     254    63515924 :     for (; i >= -nx; i--) zd[i] = xd[i];
     255             :   else
     256   637934376 :     while (zd[i+1] == 0) { i++; lz--; } /* shorten z */
     257   701450300 :   zd += i+1;
     258   701450300 :   *--zd = evalsigne(1) | evallgefint(lz);
     259   701450300 :   *--zd = evaltyp(t_INT) | evallg(lz);
     260   701450300 :   avma=(pari_sp)zd; return zd;
     261             : }
     262             : 
     263             : static void
     264   206502668 : roundr_up_ip(GEN x, long l)
     265             : {
     266   206502668 :   long i = l;
     267             :   for(;;)
     268             :   {
     269   206628422 :     if (++uel(x,--i)) break;
     270      180810 :     if (i == 2) { x[2] = (long)HIGHBIT; shiftr_inplace(x, 1); break; }
     271      125754 :   }
     272   206502668 : }
     273             : 
     274             : void
     275    76183617 : affir(GEN x, GEN y)
     276             : {
     277    76183617 :   const long s = signe(x), ly = lg(y);
     278             :   long lx, sh, i;
     279             : 
     280    76183617 :   if (!s)
     281             :   {
     282      854928 :     y[1] = evalexpo(-prec2nbits(ly));
     283      854928 :     return;
     284             :   }
     285             : 
     286    75328689 :   lx = lgefint(x); sh = bfffo(x[2]);
     287    75328689 :   y[1] = evalsigne(s) | evalexpo(bit_accuracy(lx)-sh-1);
     288    75328689 :   if (sh) {
     289    74768325 :     if (lx <= ly)
     290             :     {
     291    52378269 :       for (i=lx; i<ly; i++) y[i]=0;
     292    52378269 :       shift_left(y,x,2,lx-1, 0,sh);
     293    52378269 :       return;
     294             :     }
     295    22390056 :     shift_left(y,x,2,ly-1, x[ly],sh);
     296             :     /* lx > ly: round properly */
     297    22390056 :     if ((uel(x,ly)<<sh) & HIGHBIT) roundr_up_ip(y, ly);
     298             :   }
     299             :   else {
     300      560364 :     if (lx <= ly)
     301             :     {
     302      358662 :       for (i=2; i<lx; i++) y[i]=x[i];
     303      358662 :       for (   ; i<ly; i++) y[i]=0;
     304      358662 :       return;
     305             :     }
     306      201702 :     for (i=2; i<ly; i++) y[i]=x[i];
     307             :     /* lx > ly: round properly */
     308      201702 :     if (uel(x,ly) & HIGHBIT) roundr_up_ip(y, ly);
     309             :   }
     310             : }
     311             : 
     312             : INLINE GEN
     313   329684811 : shiftispec(GEN x, long nx, long n)
     314             : {
     315             :   long ny, i, m;
     316             :   GEN y, yd;
     317   329684811 :   if (!n)  return icopyspec(x, nx);
     318             : 
     319   310407228 :   if (n > 0)
     320             :   {
     321   204242733 :     GEN z = (GEN)avma;
     322   204242733 :     long d = dvmdsBIL(n, &m);
     323             : 
     324   204242733 :     ny = nx+d; y = new_chunk(ny + 2); yd = y + 2;
     325   204242733 :     for ( ; d; d--) *--z = 0;
     326   204242733 :     if (!m) for (i=0; i<nx; i++) yd[i]=x[i];
     327             :     else
     328             :     {
     329   203163681 :       register const ulong sh = BITS_IN_LONG - m;
     330   203163681 :       shift_left(yd,x, 0,nx-1, 0,m);
     331   203163681 :       i = uel(x,0) >> sh;
     332             :       /* Extend y on the left? */
     333   203163681 :       if (i) { ny++; y = new_chunk(1); y[2] = i; }
     334             :     }
     335             :   }
     336             :   else
     337             :   {
     338   106164495 :     ny = nx - dvmdsBIL(-n, &m);
     339   106164495 :     if (ny<1) return gen_0;
     340   105962907 :     y = new_chunk(ny + 2); yd = y + 2;
     341   105962907 :     if (m) {
     342    97201311 :       shift_right(yd,x, 0,ny, 0,m);
     343    97201311 :       if (yd[0] == 0)
     344             :       {
     345     6900963 :         if (ny==1) { avma = (pari_sp)(y+3); return gen_0; }
     346     6165441 :         ny--; avma = (pari_sp)(++y);
     347             :       }
     348             :     } else {
     349     8761596 :       for (i=0; i<ny; i++) yd[i]=x[i];
     350             :     }
     351             :   }
     352   309470118 :   y[1] = evalsigne(1)|evallgefint(ny + 2);
     353   309470118 :   y[0] = evaltyp(t_INT)|evallg(ny + 2); return y;
     354             : }
     355             : 
     356             : GEN
     357     8409753 : mantissa2nr(GEN x, long n)
     358             : { /*This is a kludge since x is not an integer*/
     359     8409753 :   long s = signe(x);
     360             :   GEN y;
     361             : 
     362     8409753 :   if(s == 0) return gen_0;
     363     8408883 :   y = shiftispec(x + 2, lg(x) - 2, n);
     364     8408883 :   if (signe(y)) setsigne(y, s);
     365     8408883 :   return y;
     366             : }
     367             : 
     368             : GEN
     369     1102605 : truncr(GEN x)
     370             : {
     371             :   long d,e,i,s,m;
     372             :   GEN y;
     373             : 
     374     1102605 :   if ((s=signe(x)) == 0 || (e=expo(x)) < 0) return gen_0;
     375      968736 :   d = nbits2prec(e+1); m = remsBIL(e);
     376      968736 :   if (d > lg(x)) pari_err_PREC( "truncr (precision loss in truncation)");
     377             : 
     378      968736 :   y=cgeti(d); y[1] = evalsigne(s) | evallgefint(d);
     379      968736 :   if (++m == BITS_IN_LONG)
     380           3 :     for (i=2; i<d; i++) y[i]=x[i];
     381             :   else
     382      968733 :     shift_right(y,x, 2,d,0, BITS_IN_LONG - m);
     383      968736 :   return y;
     384             : }
     385             : 
     386             : /* integral part */
     387             : GEN
     388      325704 : floorr(GEN x)
     389             : {
     390             :   long d,e,i,lx,m;
     391             :   GEN y;
     392             : 
     393      325704 :   if (signe(x) >= 0) return truncr(x);
     394      134067 :   if ((e=expo(x)) < 0) return gen_m1;
     395       56115 :   d = nbits2prec(e+1); m = remsBIL(e);
     396       56115 :   lx=lg(x); if (d>lx) pari_err_PREC( "floorr (precision loss in truncation)");
     397       56115 :   y = new_chunk(d);
     398       56115 :   if (++m == BITS_IN_LONG)
     399             :   {
     400           0 :     for (i=2; i<d; i++) y[i]=x[i];
     401           0 :     i=d; while (i<lx && !x[i]) i++;
     402           0 :     if (i==lx) goto END;
     403             :   }
     404             :   else
     405             :   {
     406       56115 :     shift_right(y,x, 2,d,0, BITS_IN_LONG - m);
     407       56115 :     if (uel(x,d-1)<<m == 0)
     408             :     {
     409       14376 :       i=d; while (i<lx && !x[i]) i++;
     410       14376 :       if (i==lx) goto END;
     411             :     }
     412             :   }
     413             :   /* set y:=y+1 */
     414       50985 :   for (i=d-1; i>=2; i--) { uel(y,i)++; if (y[i]) goto END; }
     415           0 :   y=new_chunk(1); y[2]=1; d++;
     416             : END:
     417       56115 :   y[1] = evalsigne(-1) | evallgefint(d);
     418       56115 :   y[0] = evaltyp(t_INT) | evallg(d); return y;
     419             : }
     420             : 
     421             : INLINE int
     422  1873781559 : cmpiispec(GEN x, GEN y, long lx, long ly)
     423             : {
     424             :   long i;
     425  1873781559 :   if (lx < ly) return -1;
     426  1798659000 :   if (lx > ly) return  1;
     427  1694050054 :   i = 0; while (i<lx && x[i]==y[i]) i++;
     428  1694050054 :   if (i==lx) return 0;
     429  1612586056 :   return (uel(x,i) > uel(y,i))? 1: -1;
     430             : }
     431             : 
     432             : INLINE int
     433    41103726 : equaliispec(GEN x, GEN y, long lx, long ly)
     434             : {
     435             :   long i;
     436    41103726 :   if (lx != ly) return 0;
     437    41103693 :   i = ly-1; while (i>=0 && x[i]==y[i]) i--;
     438    41103693 :   return i < 0;
     439             : }
     440             : 
     441             : /***********************************************************************/
     442             : /**                                                                   **/
     443             : /**                          MULTIPLICATION                           **/
     444             : /**                                                                   **/
     445             : /***********************************************************************/
     446             : /* assume ny > 0 */
     447             : INLINE GEN
     448  1439865888 : muluispec(ulong x, GEN y, long ny)
     449             : {
     450  1439865888 :   GEN yd, z = (GEN)avma;
     451  1439865888 :   long lz = ny+3;
     452             :   LOCAL_HIREMAINDER;
     453             : 
     454  1439865888 :   (void)new_chunk(lz);
     455  1439865888 :   yd = y + ny; *--z = mulll(x, *--yd);
     456  1439865888 :   while (yd > y) *--z = addmul(x,*--yd);
     457  1439865888 :   if (hiremainder) *--z = hiremainder; else lz--;
     458  1439865888 :   *--z = evalsigne(1) | evallgefint(lz);
     459  1439865888 :   *--z = evaltyp(t_INT) | evallg(lz);
     460  1439865888 :   avma=(pari_sp)z; return z;
     461             : }
     462             : 
     463             : /* a + b*|Y| */
     464             : GEN
     465      428220 : addumului(ulong a, ulong b, GEN Y)
     466             : {
     467             :   GEN yd,y,z;
     468             :   long ny,lz;
     469             :   LOCAL_HIREMAINDER;
     470             :   LOCAL_OVERFLOW;
     471             : 
     472      428220 :   if (!b || !signe(Y)) return utoi(a);
     473             : 
     474      428217 :   y = LIMBS(Y); z = (GEN)avma;
     475      428217 :   ny = NLIMBS(Y);
     476      428217 :   lz = ny+3;
     477             : 
     478      428217 :   (void)new_chunk(lz);
     479      428217 :   yd = y + ny; *--z = addll(a, mulll(b, *--yd));
     480      428217 :   if (overflow) hiremainder++; /* can't overflow */
     481      428217 :   while (yd > y) *--z = addmul(b,*--yd);
     482      428217 :   if (hiremainder) *--z = hiremainder; else lz--;
     483      428217 :   *--z = evalsigne(1) | evallgefint(lz);
     484      428217 :   *--z = evaltyp(t_INT) | evallg(lz);
     485      428217 :   avma=(pari_sp)z; return z;
     486             : }
     487             : 
     488             : /***********************************************************************/
     489             : /**                                                                   **/
     490             : /**                          DIVISION                                 **/
     491             : /**                                                                   **/
     492             : /***********************************************************************/
     493             : 
     494             : ulong
     495   906841743 : umodiu(GEN y, ulong x)
     496             : {
     497   906841743 :   long sy=signe(y),ly,i;
     498             :   ulong xi;
     499             :   LOCAL_HIREMAINDER;
     500             : 
     501   906841743 :   if (!x) pari_err_INV("umodiu",gen_0);
     502   906841743 :   if (!sy) return 0;
     503   760939389 :   ly = lgefint(y);
     504   760939389 :   if (x <= uel(y,2))
     505             :   {
     506   312622060 :     hiremainder=0;
     507   312622060 :     if (ly==3)
     508             :     {
     509   302549701 :       hiremainder=uel(y,2)%x;
     510   302549701 :       if (!hiremainder) return 0;
     511   278272066 :       return (sy > 0)? hiremainder: x - hiremainder;
     512             :     }
     513             :   }
     514             :   else
     515             :   {
     516   448317329 :     if (ly==3) return (sy > 0)? uel(y,2): x - uel(y,2);
     517   102715785 :     hiremainder=uel(y,2); ly--; y++;
     518             :   }
     519   112788144 :   xi = get_Fl_red(x);
     520   112788144 :   for (i=2; i<ly; i++) (void)divll_pre(y[i],x,xi);
     521   112788144 :   if (!hiremainder) return 0;
     522   110529024 :   return (sy > 0)? hiremainder: x - hiremainder;
     523             : }
     524             : 
     525             : /* return |y| \/ x */
     526             : GEN
     527   104245122 : diviu_rem(GEN y, ulong x, ulong *rem)
     528             : {
     529             :   long ly,i;
     530             :   GEN z;
     531             :   ulong xi;
     532             :   LOCAL_HIREMAINDER;
     533             : 
     534   104245122 :   if (!x) pari_err_INV("diviu_rem",gen_0);
     535   104245122 :   if (!signe(y)) { *rem = 0; return gen_0; }
     536             : 
     537   103841493 :   ly = lgefint(y);
     538   103841493 :   if (x <= uel(y,2))
     539             :   {
     540    88447446 :     hiremainder=0;
     541    88447446 :     if (ly==3)
     542             :     {
     543    28859088 :       z = cgetipos(3);
     544    28859088 :       z[2] = divll(uel(y,2),x);
     545    28859088 :       *rem = hiremainder; return z;
     546             :     }
     547             :   }
     548             :   else
     549             :   {
     550    15394047 :     if (ly==3) { *rem = uel(y,2); return gen_0; }
     551    12080460 :     hiremainder = uel(y,2); ly--; y++;
     552             :   }
     553    71668818 :   xi = get_Fl_red(x);
     554    71668818 :   z = cgetipos(ly);
     555    71668818 :   for (i=2; i<ly; i++) z[i]=divll_pre(y[i],x,xi);
     556    71668818 :   *rem = hiremainder; return z;
     557             : }
     558             : 
     559             : GEN
     560    22099419 : divis_rem(GEN y, long x, long *rem)
     561             : {
     562    22099419 :   long sy=signe(y),ly,s,i;
     563             :   GEN z;
     564             :   ulong xi;
     565             :   LOCAL_HIREMAINDER;
     566             : 
     567    22099419 :   if (!x) pari_err_INV("divis_rem",gen_0);
     568    22099419 :   if (!sy) { *rem=0; return gen_0; }
     569    17831751 :   if (x<0) { s = -sy; x = -x; } else s = sy;
     570             : 
     571    17831751 :   ly = lgefint(y);
     572    17831751 :   if ((ulong)x <= uel(y,2))
     573             :   {
     574    14149051 :     hiremainder=0;
     575    14149051 :     if (ly==3)
     576             :     {
     577    13771282 :       z = cgeti(3); z[1] = evallgefint(3) | evalsigne(s);
     578    13771282 :       z[2] = divll(uel(y,2),x);
     579    13771282 :       if (sy<0) hiremainder = - ((long)hiremainder);
     580    13771282 :       *rem = (long)hiremainder; return z;
     581             :     }
     582             :   }
     583             :   else
     584             :   {
     585     3682700 :     if (ly==3) { *rem = itos(y); return gen_0; }
     586       75789 :     hiremainder = uel(y,2); ly--; y++;
     587             :   }
     588      453558 :   xi = get_Fl_red(x);
     589      453558 :   z = cgeti(ly); z[1] = evallgefint(ly) | evalsigne(s);
     590      453558 :   for (i=2; i<ly; i++) z[i]=divll_pre(y[i],x,xi);
     591      453558 :   if (sy<0) hiremainder = - ((long)hiremainder);
     592      453558 :   *rem = (long)hiremainder; return z;
     593             : }
     594             : 
     595             : GEN
     596      350229 : divis(GEN y, long x)
     597             : {
     598      350229 :   long sy=signe(y),ly,s,i;
     599             :   ulong xi;
     600             :   GEN z;
     601             :   LOCAL_HIREMAINDER;
     602             : 
     603      350229 :   if (!x) pari_err_INV("divis",gen_0);
     604      350229 :   if (!sy) return gen_0;
     605      350229 :   if (x<0) { s = -sy; x = -x; } else s = sy;
     606             : 
     607      350229 :   ly = lgefint(y);
     608      350229 :   if ((ulong)x <= uel(y,2))
     609             :   {
     610      329178 :     hiremainder=0;
     611      329178 :     if (ly==3)
     612             :     {
     613      132552 :       z = cgeti(3); z[1] = evallgefint(3) | evalsigne(s);
     614      132552 :       z[2] = divll(y[2],x);
     615      132552 :       return z;
     616             :     }
     617             :   }
     618             :   else
     619             :   {
     620       21051 :     if (ly==3) return gen_0;
     621       21051 :     hiremainder=y[2]; ly--; y++;
     622             :   }
     623      217677 :   xi = get_Fl_red(x);
     624      217677 :   z = cgeti(ly); z[1] = evallgefint(ly) | evalsigne(s);
     625      217677 :   for (i=2; i<ly; i++) z[i]=divll_pre(y[i],x, xi);
     626      217677 :   return z;
     627             : }
     628             : 
     629             : GEN
     630    63684375 : divrr(GEN x, GEN y)
     631             : {
     632    63684375 :   long sx=signe(x), sy=signe(y), lx,ly,lr,e,i,j;
     633             :   ulong y0,y1;
     634             :   GEN r, r1;
     635             : 
     636    63684375 :   if (!x) pari_err_INV("divrr",y);
     637    63684375 :   e = expo(x) - expo(y);
     638    63684375 :   if (!sx) return real_0_bit(e);
     639    62615022 :   if (sy<0) sx = -sx;
     640             : 
     641    62615022 :   lx=lg(x); ly=lg(y);
     642    62615022 :   if (ly==3)
     643             :   {
     644    40206696 :     ulong k = x[2], l = (lx>3)? x[3]: 0;
     645             :     LOCAL_HIREMAINDER;
     646    40206696 :     if (k < uel(y,2)) e--;
     647             :     else
     648             :     {
     649    18356247 :       l >>= 1; if (k&1) l |= HIGHBIT;
     650    18356247 :       k >>= 1;
     651             :     }
     652    40206696 :     hiremainder = k; k = divll(l,y[2]);
     653    40206696 :     if (hiremainder & HIGHBIT)
     654             :     {
     655    11159046 :       k++;
     656    11159046 :       if (!k) { k = HIGHBIT; e++; }
     657             :     }
     658    40206696 :     r = cgetr(3);
     659    40206696 :     r[1] = evalsigne(sx) | evalexpo(e);
     660    40206696 :     r[2] = k; return r;
     661             :   }
     662             : 
     663    22408326 :   lr = minss(lx,ly); r = new_chunk(lr);
     664    22408326 :   r1 = r-1;
     665    22408326 :   r1[1] = 0; for (i=2; i<lr; i++) r1[i]=x[i];
     666    22408326 :   r1[lr] = (lx>ly)? x[lr]: 0;
     667    22408326 :   y0 = y[2]; y1 = y[3];
     668   173040888 :   for (i=0; i<lr-1; i++)
     669             :   { /* r1 = r + (i-1), OK up to r1[2] (accesses at most r[lr]) */
     670             :     ulong k, qp;
     671             :     LOCAL_HIREMAINDER;
     672             :     LOCAL_OVERFLOW;
     673             : 
     674   150632562 :     if (uel(r1,1) == y0)
     675             :     {
     676      360324 :       qp = ULONG_MAX; k = addll(y0,r1[2]);
     677             :     }
     678             :     else
     679             :     {
     680   150272238 :       if (uel(r1,1) > y0) /* can't happen if i=0 */
     681             :       {
     682           0 :         GEN y1 = y+1;
     683           0 :         j = lr-i; r1[j] = subll(r1[j],y1[j]);
     684           0 :         for (j--; j>0; j--) r1[j] = subllx(r1[j],y1[j]);
     685           0 :         j=i; do uel(r,--j)++; while (j && !uel(r,j));
     686             :       }
     687   150272238 :       hiremainder = r1[1]; overflow = 0;
     688   150272238 :       qp = divll(r1[2],y0); k = hiremainder;
     689             :     }
     690   150632562 :     j = lr-i+1;
     691   150632562 :     if (!overflow)
     692             :     {
     693             :       long k3, k4;
     694   150353271 :       k3 = mulll(qp,y1);
     695   150353271 :       if (j == 3) /* i = lr - 2 maximal, r1[3] undefined -> 0 */
     696    22398939 :         k4 = subll(hiremainder,k);
     697             :       else
     698             :       {
     699   127954332 :         k3 = subll(k3, r1[3]);
     700   127954332 :         k4 = subllx(hiremainder,k);
     701             :       }
     702   150353271 :       while (!overflow && k4) { qp--; k3 = subll(k3,y1); k4 = subllx(k4,y0); }
     703             :     }
     704   150632562 :     if (j<ly) (void)mulll(qp,y[j]); else { hiremainder = 0 ; j = ly; }
     705  1452847419 :     for (j--; j>1; j--)
     706             :     {
     707  1302214857 :       r1[j] = subll(r1[j], addmul(qp,y[j]));
     708  1302214857 :       hiremainder += overflow;
     709             :     }
     710   150632562 :     if (uel(r1,1) != hiremainder)
     711             :     {
     712      136539 :       if (uel(r1,1) < hiremainder)
     713             :       {
     714      136539 :         qp--;
     715      136539 :         j = lr-i-(lr-i>=ly); r1[j] = addll(r1[j], y[j]);
     716      136539 :         for (j--; j>1; j--) r1[j] = addllx(r1[j], y[j]);
     717             :       }
     718             :       else
     719             :       {
     720           0 :         r1[1] -= hiremainder;
     721           0 :         while (r1[1])
     722             :         {
     723           0 :           qp++; if (!qp) { j=i; do uel(r,--j)++; while (j && !r[j]); }
     724           0 :           j = lr-i-(lr-i>=ly); r1[j] = subll(r1[j],y[j]);
     725           0 :           for (j--; j>1; j--) r1[j] = subllx(r1[j],y[j]);
     726           0 :           r1[1] -= overflow;
     727             :         }
     728             :       }
     729             :     }
     730   150632562 :     *++r1 = qp;
     731             :   }
     732             :   /* i = lr-1 */
     733             :   /* round correctly */
     734    22408326 :   if (uel(r1,1) > (y0>>1))
     735             :   {
     736    11007950 :     j=i; do uel(r,--j)++; while (j && !r[j]);
     737             :   }
     738    22408326 :   r1 = r-1; for (j=i; j>=2; j--) r[j]=r1[j];
     739    22408326 :   if (r[0] == 0) e--;
     740     9240423 :   else if (r[0] == 1) { shift_right(r,r, 2,lr, 1,1); }
     741             :   else { /* possible only when rounding up to 0x2 0x0 ... */
     742          24 :     r[2] = (long)HIGHBIT; e++;
     743             :   }
     744    22408326 :   r[0] = evaltyp(t_REAL)|evallg(lr);
     745    22408326 :   r[1] = evalsigne(sx) | evalexpo(e);
     746    22408326 :   return r;
     747             : }
     748             : 
     749             : GEN
     750    13334976 : divri(GEN x, GEN y)
     751             : {
     752    13334976 :   long lx, s = signe(y);
     753             :   pari_sp av;
     754             :   GEN z;
     755             : 
     756    13334976 :   if (!s) pari_err_INV("divri",y);
     757    13334976 :   if (!signe(x)) return real_0_bit(expo(x) - expi(y));
     758    13309770 :   if (!is_bigint(y)) {
     759    11227674 :     GEN z = divru(x, y[2]);
     760    11227674 :     if (s < 0) togglesign(z);
     761    11227674 :     return z;
     762             :   }
     763     2082096 :   lx = lg(x); z = cgetr(lx); av = avma;
     764     2082096 :   affrr(divrr(x, itor(y, lx+1)), z);
     765     2082096 :   avma = av; return z;
     766             : }
     767             : 
     768             : /* Integer division x / y: such that sign(r) = sign(x)
     769             :  *   if z = ONLY_REM return remainder, otherwise return quotient
     770             :  *   if z != NULL set *z to remainder
     771             :  *   *z is the last object on stack (and thus can be disposed of with cgiv
     772             :  *   instead of gerepile)
     773             :  * If *z is zero, we put gen_0 here and no copy.
     774             :  * space needed: lx + ly */
     775             : GEN
     776   540692997 : dvmdii(GEN x, GEN y, GEN *z)
     777             : {
     778   540692997 :   long sx = signe(x), sy = signe(y);
     779   540692997 :   long lx, ly = lgefint(y), lz, i, j, sh, lq, lr;
     780             :   pari_sp av;
     781             :   ulong y0,y0i,y1, *xd,*rd,*qd;
     782             :   GEN q, r, r1;
     783             : 
     784   540692997 :   if (!sx)
     785             :   {
     786    12334740 :     if (ly < 3) pari_err_INV("dvmdii",gen_0);
     787    12334737 :     if (!z || z == ONLY_REM) return gen_0;
     788     6301899 :     *z=gen_0; return gen_0;
     789             :   }
     790   528358257 :   if (ly <= 3)
     791             :   {
     792             :     ulong rem;
     793   345072423 :     if (ly < 3) pari_err_INV("dvmdii",gen_0);
     794   345072423 :     if (z == ONLY_REM)
     795             :     {
     796   308002923 :       rem = umodiu(x,uel(y,2));
     797   308002923 :       if (!rem) return gen_0;
     798   287468418 :       return (sx < 0)? utoineg(uel(y,2) - rem): utoipos(rem);
     799             :     }
     800    37069500 :     q = diviu_rem(x, uel(y,2), &rem);
     801    37069500 :     if (sx != sy) togglesign(q);
     802    37069500 :     if (!z) return q;
     803    33989244 :     if (!rem) *z = gen_0;
     804    16045962 :     else *z = sx < 0? utoineg(rem): utoipos(rem);
     805    33989244 :     return q;
     806             :   }
     807   183285834 :   lx=lgefint(x);
     808   183285834 :   lz=lx-ly;
     809   183285834 :   if (lz <= 0)
     810             :   {
     811    75033804 :     if (lz == 0)
     812             :     {
     813    66868824 :       for (i=2; i<lx; i++)
     814    66723987 :         if (x[i] != y[i])
     815             :         {
     816    64118643 :           if (uel(x,i) > uel(y,i)) goto DIVIDE;
     817    20915142 :           goto TRIVIAL;
     818             :         }
     819      144837 :       if (z == ONLY_REM) return gen_0;
     820        6948 :       if (z) *z = gen_0;
     821        6948 :       if (sx < 0) sy = -sy;
     822        6948 :       return stoi(sy);
     823             :     }
     824             : TRIVIAL:
     825    31685466 :     if (z == ONLY_REM) return icopy(x);
     826     2388369 :     if (z) *z = icopy(x);
     827     2388369 :     return gen_0;
     828             :   }
     829             : DIVIDE: /* quotient is non-zero */
     830   151455531 :   av=avma; if (sx<0) sy = -sy;
     831   151455531 :   r1 = new_chunk(lx); sh = bfffo(y[2]);
     832   151455531 :   if (sh)
     833             :   { /* normalize so that highbit(y) = 1 (shift left x and y by sh bits)*/
     834   143265492 :     register const ulong m = BITS_IN_LONG - sh;
     835   143265492 :     r = new_chunk(ly);
     836   143265492 :     shift_left(r, y,2,ly-1, 0,sh); y = r;
     837   143265492 :     shift_left(r1,x,2,lx-1, 0,sh);
     838   143265492 :     r1[1] = uel(x,2) >> m;
     839             :   }
     840             :   else
     841             :   {
     842     8190039 :     r1[1] = 0; for (j=2; j<lx; j++) r1[j] = x[j];
     843             :   }
     844   151455531 :   x = r1;
     845   151455531 :   y0 = y[2]; y0i = get_Fl_red(y0);
     846   151455531 :   y1 = y[3];
     847   651945693 :   for (i=0; i<=lz; i++)
     848             :   { /* r1 = x + i */
     849             :     ulong k, qp;
     850             :     LOCAL_HIREMAINDER;
     851             :     LOCAL_OVERFLOW;
     852             : 
     853   500490162 :     if (uel(r1,1) == y0)
     854             :     {
     855        5978 :       qp = ULONG_MAX; k = addll(y0,r1[2]);
     856             :     }
     857             :     else
     858             :     {
     859   500484184 :       hiremainder = r1[1]; overflow = 0;
     860   500484184 :       qp = divll_pre(r1[2],y0,y0i); k = hiremainder;
     861             :     }
     862   500490162 :     if (!overflow)
     863             :     {
     864   500484838 :       long k3 = subll(mulll(qp,y1), r1[3]);
     865   500484838 :       long k4 = subllx(hiremainder,k);
     866   500484838 :       while (!overflow && k4) { qp--; k3 = subll(k3,y1); k4 = subllx(k4,y0); }
     867             :     }
     868   500490162 :     hiremainder = 0; j = ly;
     869  5849771970 :     for (j--; j>1; j--)
     870             :     {
     871  5349281808 :       r1[j] = subll(r1[j], addmul(qp,y[j]));
     872  5349281808 :       hiremainder += overflow;
     873             :     }
     874   500490162 :     if (uel(r1,1) < hiremainder)
     875             :     {
     876      231624 :       qp--;
     877      231624 :       j = ly-1; r1[j] = addll(r1[j],y[j]);
     878      231624 :       for (j--; j>1; j--) r1[j] = addllx(r1[j],y[j]);
     879             :     }
     880   500490162 :     *++r1 = qp;
     881             :   }
     882             : 
     883   151455531 :   lq = lz+2;
     884   151455531 :   if (!z)
     885             :   {
     886      451794 :     qd = (ulong*)av;
     887      451794 :     xd = (ulong*)(x + lq);
     888      451794 :     if (x[1]) { lz++; lq++; }
     889      451794 :     while (lz--) *--qd = *--xd;
     890      451794 :     *--qd = evalsigne(sy) | evallgefint(lq);
     891      451794 :     *--qd = evaltyp(t_INT) | evallg(lq);
     892      451794 :     avma = (pari_sp)qd; return (GEN)qd;
     893             :   }
     894             : 
     895   151003737 :   j=lq; while (j<lx && !x[j]) j++;
     896   151003737 :   lz = lx-j;
     897   151003737 :   if (z == ONLY_REM)
     898             :   {
     899   113768064 :     if (lz==0) { avma = av; return gen_0; }
     900   111841608 :     rd = (ulong*)av; lr = lz+2;
     901   111841608 :     xd = (ulong*)(x + lx);
     902   111841608 :     if (!sh) while (lz--) *--rd = *--xd;
     903             :     else
     904             :     { /* shift remainder right by sh bits */
     905   103719744 :       const ulong shl = BITS_IN_LONG - sh;
     906             :       ulong l;
     907   103719744 :       xd--;
     908   469642866 :       while (--lz) /* fill r[3..] */
     909             :       {
     910   262203378 :         l = *xd >> sh;
     911   262203378 :         *--rd = l | (*--xd << shl);
     912             :       }
     913   103719744 :       l = *xd >> sh;
     914   103719744 :       if (l) *--rd = l; else lr--;
     915             :     }
     916   111841608 :     *--rd = evalsigne(sx) | evallgefint(lr);
     917   111841608 :     *--rd = evaltyp(t_INT) | evallg(lr);
     918   111841608 :     avma = (pari_sp)rd; return (GEN)rd;
     919             :   }
     920             : 
     921    37235673 :   lr = lz+2;
     922    37235673 :   rd = NULL; /* gcc -Wall */
     923    37235673 :   if (lz)
     924             :   { /* non zero remainder: initialize rd */
     925    35378229 :     xd = (ulong*)(x + lx);
     926    35378229 :     if (!sh)
     927             :     {
     928       21786 :       rd = (ulong*)avma; (void)new_chunk(lr);
     929       21786 :       while (lz--) *--rd = *--xd;
     930             :     }
     931             :     else
     932             :     { /* shift remainder right by sh bits */
     933    35356443 :       const ulong shl = BITS_IN_LONG - sh;
     934             :       ulong l;
     935    35356443 :       rd = (ulong*)x; /* overwrite shifted y */
     936    35356443 :       xd--;
     937   164113371 :       while (--lz)
     938             :       {
     939    93400485 :         l = *xd >> sh;
     940    93400485 :         *--rd = l | (*--xd << shl);
     941             :       }
     942    35356443 :       l = *xd >> sh;
     943    35356443 :       if (l) *--rd = l; else lr--;
     944             :     }
     945    35378229 :     *--rd = evalsigne(sx) | evallgefint(lr);
     946    35378229 :     *--rd = evaltyp(t_INT) | evallg(lr);
     947    35378229 :     rd += lr;
     948             :   }
     949    37235673 :   qd = (ulong*)av;
     950    37235673 :   xd = (ulong*)(x + lq);
     951    37235673 :   if (x[1]) lq++;
     952    37235673 :   j = lq-2; while (j--) *--qd = *--xd;
     953    37235673 :   *--qd = evalsigne(sy) | evallgefint(lq);
     954    37235673 :   *--qd = evaltyp(t_INT) | evallg(lq);
     955    37235673 :   q = (GEN)qd;
     956    37235673 :   if (lr==2) *z = gen_0;
     957             :   else
     958             :   { /* rd has been properly initialized: we had lz > 0 */
     959    35378229 :     while (lr--) *--qd = *--rd;
     960    35378229 :     *z = (GEN)qd;
     961             :   }
     962    37235673 :   avma = (pari_sp)qd; return q;
     963             : }
     964             : 
     965             : /* Montgomery reduction.
     966             :  * N has k words, assume T >= 0 has less than 2k.
     967             :  * Return res := T / B^k mod N, where B = 2^BIL
     968             :  * such that 0 <= res < T/B^k + N  and  res has less than k words */
     969             : GEN
     970     6833880 : red_montgomery(GEN T, GEN N, ulong inv)
     971             : {
     972             :   pari_sp av;
     973             :   GEN Te, Td, Ne, Nd, scratch;
     974     6833880 :   ulong i, j, m, t, d, k = NLIMBS(N);
     975             :   int carry;
     976             :   LOCAL_HIREMAINDER;
     977             :   LOCAL_OVERFLOW;
     978             : 
     979     6833880 :   if (k == 0) return gen_0;
     980     6833880 :   d = NLIMBS(T); /* <= 2*k */
     981     6833880 :   if (d == 0) return gen_0;
     982             : #ifdef DEBUG
     983             :   if (d > 2*k) pari_err_BUG("red_montgomery");
     984             : #endif
     985     6833871 :   if (k == 1)
     986             :   { /* as below, special cased for efficiency */
     987           0 :     ulong n = uel(N,2);
     988           0 :     if (d == 1) {
     989           0 :       hiremainder = uel(T,2);
     990           0 :       m = hiremainder * inv;
     991           0 :       (void)addmul(m, n); /* t + m*n = 0 */
     992           0 :       return utoi(hiremainder);
     993             :     } else { /* d = 2 */
     994           0 :       hiremainder = uel(T,3);
     995           0 :       m = hiremainder * inv;
     996           0 :       (void)addmul(m, n); /* t + m*n = 0 */
     997           0 :       t = addll(hiremainder, uel(T,2));
     998           0 :       if (overflow) t -= n; /* t > n doesn't fit in 1 word */
     999           0 :       return utoi(t);
    1000             :     }
    1001             :   }
    1002             :   /* assume k >= 2 */
    1003     6833871 :   av = avma; scratch = new_chunk(k<<1); /* >= k + 2: result fits */
    1004             : 
    1005             :   /* copy T to scratch space (pad with zeroes to 2k words) */
    1006     6833871 :   Td = (GEN)av;
    1007     6833871 :   Te = T + (d+2);
    1008     6833871 :   for (i=0; i < d     ; i++) *--Td = *--Te;
    1009     6833871 :   for (   ; i < (k<<1); i++) *--Td = 0;
    1010             : 
    1011     6833871 :   Te = (GEN)av; /* 1 beyond end of current T mantissa (in scratch) */
    1012     6833871 :   Ne = N + k+2; /* 1 beyond end of N mantissa */
    1013             : 
    1014     6833871 :   carry = 0;
    1015    35568384 :   for (i=0; i<k; i++) /* set T := T/B nod N, k times */
    1016             :   {
    1017    28734513 :     Td = Te; /* one beyond end of (new) T mantissa */
    1018    28734513 :     Nd = Ne;
    1019    28734513 :     hiremainder = *--Td;
    1020    28734513 :     m = hiremainder * inv; /* solve T + m N = O(B) */
    1021             : 
    1022             :     /* set T := (T + mN) / B */
    1023    28734513 :     Te = Td;
    1024    28734513 :     (void)addmul(m, *--Nd); /* = 0 */
    1025   274285335 :     for (j=1; j<k; j++)
    1026             :     {
    1027   245550822 :       t = addll(addmul(m, *--Nd), *--Td);
    1028   245550822 :       *Td = t;
    1029   245550822 :       hiremainder += overflow;
    1030             :     }
    1031    28734513 :     t = addll(hiremainder, *--Td); *Td = t + carry;
    1032    28734513 :     carry = (overflow || (carry && *Td == 0));
    1033             :   }
    1034     6833871 :   if (carry)
    1035             :   { /* Td > N overflows (k+1 words), set Td := Td - N */
    1036       89868 :     Td = Te;
    1037       89868 :     Nd = Ne;
    1038       89868 :     t = subll(*--Td, *--Nd); *Td = t;
    1039       89868 :     while (Td > scratch) { t = subllx(*--Td, *--Nd); *Td = t; }
    1040             :   }
    1041             : 
    1042             :   /* copy result */
    1043     6833871 :   Td = (GEN)av;
    1044     6833871 :   while (*scratch == 0 && Te > scratch) scratch++; /* strip leading 0s */
    1045     6833871 :   while (Te > scratch) *--Td = *--Te;
    1046     6833871 :   k = (GEN)av - Td; if (!k) { avma = av; return gen_0; }
    1047     6833871 :   k += 2;
    1048     6833871 :   *--Td = evalsigne(1) | evallgefint(k);
    1049     6833871 :   *--Td = evaltyp(t_INT) | evallg(k);
    1050             : #ifdef DEBUG
    1051             : {
    1052             :   long l = NLIMBS(N), s = BITS_IN_LONG*l;
    1053             :   GEN R = int2n(s);
    1054             :   GEN res = remii(mulii(T, Fp_inv(R, N)), N);
    1055             :   if (k > lgefint(N)
    1056             :     || !equalii(remii(Td,N),res)
    1057             :     || cmpii(Td, addii(shifti(T, -s), N)) >= 0) pari_err_BUG("red_montgomery");
    1058             : }
    1059             : #endif
    1060     6833871 :   avma = (pari_sp)Td; return Td;
    1061             : }
    1062             : 
    1063             : /* EXACT INTEGER DIVISION */
    1064             : 
    1065             : /* assume xy>0, the division is exact and y is odd. Destroy x */
    1066             : static GEN
    1067    15786150 : diviuexact_i(GEN x, ulong y)
    1068             : {
    1069             :   long i, lz, lx;
    1070             :   ulong q, yinv;
    1071             :   GEN z, z0, x0, x0min;
    1072             : 
    1073    15786150 :   if (y == 1) return icopy(x);
    1074    13225113 :   lx = lgefint(x);
    1075    13225113 :   if (lx == 3) return utoipos(uel(x,2) / y);
    1076    12845184 :   yinv = invmod2BIL(y);
    1077    12845184 :   lz = (y <= uel(x,2)) ? lx : lx-1;
    1078    12845184 :   z = new_chunk(lz);
    1079    12845184 :   z0 = z + lz;
    1080    12845184 :   x0 = x + lx; x0min = x + lx-lz+2;
    1081             : 
    1082    59551248 :   while (x0 > x0min)
    1083             :   {
    1084    33860880 :     *--z0 = q = yinv*uel(--x0,0); /* i-th quotient */
    1085    33860880 :     if (!q) continue;
    1086             :     /* x := x - q * y */
    1087             :     { /* update neither lowest word (could set it to 0) nor highest ones */
    1088    30694602 :       register GEN x1 = x0 - 1;
    1089             :       LOCAL_HIREMAINDER;
    1090    30694602 :       (void)mulll(q,y);
    1091    30694602 :       if (hiremainder)
    1092             :       {
    1093    26623782 :         if (uel(x1,0) < hiremainder)
    1094             :         {
    1095     1189542 :           uel(x1,0) -= hiremainder;
    1096     1191276 :           do uel(--x1,0)--; while (uel(x1,0) == ULONG_MAX);
    1097             :         }
    1098             :         else
    1099    25434240 :           uel(x1,0) -= hiremainder;
    1100             :       }
    1101             :     }
    1102             :   }
    1103    12845184 :   i=2; while(!z[i]) i++;
    1104    12845184 :   z += i-2; lz -= i-2;
    1105    12845184 :   z[0] = evaltyp(t_INT)|evallg(lz);
    1106    12845184 :   z[1] = evalsigne(1)|evallg(lz);
    1107    12845184 :   avma = (pari_sp)z; return z;
    1108             : }
    1109             : 
    1110             : /* assume y != 0 and the division is exact */
    1111             : GEN
    1112    19664298 : diviuexact(GEN x, ulong y)
    1113             : {
    1114             :   pari_sp av;
    1115    19664298 :   long lx, vy, s = signe(x);
    1116             :   GEN z;
    1117             : 
    1118    19664298 :   if (!s) return gen_0;
    1119    19182012 :   if (y == 1) return icopy(x);
    1120    19177176 :   lx = lgefint(x);
    1121    19177176 :   if (lx == 3) {
    1122    18600756 :     ulong q = uel(x,2) / y;
    1123    18600756 :     return (s > 0)? utoipos(q): utoineg(q);
    1124             :   }
    1125      576420 :   av = avma; (void)new_chunk(lx); vy = vals(y);
    1126      576420 :   if (vy) {
    1127      443313 :     y >>= vy;
    1128      443313 :     if (y == 1) { avma = av; return shifti(x, -vy); }
    1129      421167 :     x = shifti(x, -vy);
    1130      421167 :     if (lx == 3) {
    1131           0 :       ulong q = uel(x,2) / y;
    1132           0 :       avma = av;
    1133           0 :       return (s > 0)? utoipos(q): utoineg(q);
    1134             :     }
    1135      133107 :   } else x = icopy(x);
    1136      554274 :   avma = av;
    1137      554274 :   z = diviuexact_i(x, y);
    1138      554274 :   setsigne(z, s); return z;
    1139             : }
    1140             : 
    1141             : /* Find z such that x=y*z, knowing that y | x (unchecked)
    1142             :  * Method: y0 z0 = x0 mod B = 2^BITS_IN_LONG ==> z0 = 1/y0 mod B.
    1143             :  *    Set x := (x - z0 y) / B, updating only relevant words, and repeat */
    1144             : GEN
    1145   171202419 : diviiexact(GEN x, GEN y)
    1146             : {
    1147   171202419 :   long lx, ly, lz, vy, i, ii, sx = signe(x), sy = signe(y);
    1148             :   pari_sp av;
    1149             :   ulong y0inv,q;
    1150             :   GEN z;
    1151             : 
    1152   171202419 :   if (!sy) pari_err_INV("diviiexact",gen_0);
    1153   171202416 :   if (!sx) return gen_0;
    1154   137869134 :   lx = lgefint(x);
    1155   137869134 :   if (lx == 3) {
    1156   117563142 :     q = uel(x,2) / uel(y,2);
    1157   117563142 :     return (sx+sy) ? utoipos(q): utoineg(q);
    1158             :   }
    1159    20305992 :   vy = vali(y); av = avma;
    1160    20305992 :   (void)new_chunk(lx); /* enough room for z */
    1161    20305992 :   if (vy)
    1162             :   { /* make y odd */
    1163     9101820 :     y = shifti(y,-vy);
    1164     9101820 :     x = shifti(x,-vy); lx = lgefint(x);
    1165             :   }
    1166    11204172 :   else x = icopy(x); /* necessary because we destroy x */
    1167    20305992 :   avma = av; /* will erase our x,y when exiting */
    1168             :   /* now y is odd */
    1169    20305992 :   ly = lgefint(y);
    1170    20305992 :   if (ly == 3)
    1171             :   {
    1172    15231876 :     x = diviuexact_i(x,uel(y,2)); /* x != 0 */
    1173    15231876 :     setsigne(x, (sx+sy)? 1: -1); return x;
    1174             :   }
    1175     5074116 :   y0inv = invmod2BIL(y[ly-1]);
    1176     5074116 :   i=2; while (i<ly && y[i]==x[i]) i++;
    1177     5074116 :   lz = (i==ly || uel(y,i) < uel(x,i)) ? lx-ly+3 : lx-ly+2;
    1178     5074116 :   z = new_chunk(lz);
    1179             : 
    1180     5074116 :   y += ly - 1; /* now y[-i] = i-th word of y */
    1181    23725452 :   for (ii=lx-1,i=lz-1; i>=2; i--,ii--)
    1182             :   {
    1183             :     long limj;
    1184             :     LOCAL_HIREMAINDER;
    1185             :     LOCAL_OVERFLOW;
    1186             : 
    1187    18651336 :     z[i] = q = y0inv*uel(x,ii); /* i-th quotient */
    1188    18651336 :     if (!q) continue;
    1189             : 
    1190             :     /* x := x - q * y */
    1191    16552635 :     (void)mulll(q,y[0]); limj = maxss(lx - lz, ii+3-ly);
    1192             :     { /* update neither lowest word (could set it to 0) nor highest ones */
    1193    16552635 :       register GEN x0 = x + (ii - 1), y0 = y - 1, xlim = x + limj;
    1194    70738776 :       for (; x0 >= xlim; x0--, y0--)
    1195             :       {
    1196    54186141 :         *x0 = subll(*x0, addmul(q,*y0));
    1197    54186141 :         hiremainder += overflow;
    1198             :       }
    1199    16552635 :       if (hiremainder && limj != lx - lz)
    1200             :       {
    1201     9039558 :         if ((ulong)*x0 < hiremainder)
    1202             :         {
    1203       78009 :           *x0 -= hiremainder;
    1204       78009 :           do (*--x0)--; while ((ulong)*x0 == ULONG_MAX);
    1205             :         }
    1206             :         else
    1207     8961549 :           *x0 -= hiremainder;
    1208             :       }
    1209             :     }
    1210             :   }
    1211     5074116 :   i=2; while(!z[i]) i++;
    1212     5074116 :   z += i-2; lz -= (i-2);
    1213     5074116 :   z[0] = evaltyp(t_INT)|evallg(lz);
    1214     5074116 :   z[1] = evalsigne((sx+sy)? 1: -1) | evallg(lz);
    1215     5074116 :   avma = (pari_sp)z; return z;
    1216             : }
    1217             : 
    1218             : /* assume yz != and yz | x */
    1219             : GEN
    1220      395091 : diviuuexact(GEN x, ulong y, ulong z)
    1221             : {
    1222             :   long tmp[4];
    1223             :   ulong t;
    1224             :   LOCAL_HIREMAINDER;
    1225      395091 :   t = mulll(y, z);
    1226      395091 :   if (!hiremainder) return diviuexact(x, t);
    1227           0 :   tmp[0] = evaltyp(t_INT)|_evallg(4);
    1228           0 :   tmp[1] = evalsigne(1)|evallgefint(4);
    1229           0 :   tmp[2] = hiremainder;
    1230           0 :   tmp[3] = t;
    1231           0 :   return diviiexact(x, tmp);
    1232             : }
    1233             : 
    1234             : /********************************************************************/
    1235             : /**                                                                **/
    1236             : /**               INTEGER MULTIPLICATION (BASECASE)                **/
    1237             : /**                                                                **/
    1238             : /********************************************************************/
    1239             : /* nx >= ny = num. of digits of x, y (not GEN, see mulii) */
    1240             : INLINE GEN
    1241  1494757493 : muliispec_basecase(GEN x, GEN y, long nx, long ny)
    1242             : {
    1243             :   GEN z2e,z2d,yd,xd,ye,zd;
    1244             :   long p1,lz;
    1245             :   LOCAL_HIREMAINDER;
    1246             : 
    1247  1494757493 :   if (ny == 1) return muluispec((ulong)*y, x, nx);
    1248   479865222 :   if (ny == 0) return gen_0;
    1249   479241270 :   zd = (GEN)avma; lz = nx+ny+2;
    1250   479241270 :   (void)new_chunk(lz);
    1251   479241270 :   xd = x + nx;
    1252   479241270 :   yd = y + ny;
    1253   479241270 :   ye = yd; p1 = *--xd;
    1254             : 
    1255   479241270 :   *--zd = mulll(p1, *--yd); z2e = zd;
    1256   479241270 :   while (yd > y) *--zd = addmul(p1, *--yd);
    1257   479241270 :   *--zd = hiremainder;
    1258             : 
    1259  3628706997 :   while (xd > x)
    1260             :   {
    1261             :     LOCAL_OVERFLOW;
    1262  2670224457 :     yd = ye; p1 = *--xd;
    1263             : 
    1264  2670224457 :     z2d = --z2e;
    1265  2670224457 :     *z2d = addll(mulll(p1, *--yd), *z2d); z2d--;
    1266 34393221204 :     while (yd > y)
    1267             :     {
    1268 29052772290 :       hiremainder += overflow;
    1269 29052772290 :       *z2d = addll(addmul(p1, *--yd), *z2d); z2d--;
    1270             :     }
    1271  2670224457 :     *--zd = hiremainder + overflow;
    1272             :   }
    1273   479241270 :   if (*zd == 0) { zd++; lz--; } /* normalize */
    1274   479241270 :   *--zd = evalsigne(1) | evallgefint(lz);
    1275   479241270 :   *--zd = evaltyp(t_INT) | evallg(lz);
    1276   479241270 :   avma=(pari_sp)zd; return zd;
    1277             : }
    1278             : 
    1279             : INLINE GEN
    1280   512423988 : sqrispec_basecase(GEN x, long nx)
    1281             : {
    1282             :   GEN z2e,z2d,yd,xd,zd,x0,z0;
    1283             :   long p1,lz;
    1284             :   LOCAL_HIREMAINDER;
    1285             :   LOCAL_OVERFLOW;
    1286             : 
    1287   512423988 :   if (nx == 1) return sqru((ulong)*x);
    1288   443293413 :   if (nx == 0) return gen_0;
    1289    81603270 :   zd = (GEN)avma; lz = (nx+1) << 1;
    1290    81603270 :   z0 = new_chunk(lz);
    1291    81603270 :   if (nx == 1)
    1292             :   {
    1293           0 :     *--zd = mulll(*x, *x);
    1294           0 :     *--zd = hiremainder; goto END;
    1295             :   }
    1296    81603270 :   xd = x + nx;
    1297             : 
    1298             :   /* compute double products --> zd */
    1299    81603270 :   p1 = *--xd; yd = xd; --zd;
    1300    81603270 :   *--zd = mulll(p1, *--yd); z2e = zd;
    1301    81603270 :   while (yd > x) *--zd = addmul(p1, *--yd);
    1302    81603270 :   *--zd = hiremainder;
    1303             : 
    1304    81603270 :   x0 = x+1;
    1305   703862775 :   while (xd > x0)
    1306             :   {
    1307             :     LOCAL_OVERFLOW;
    1308   540656235 :     p1 = *--xd; yd = xd;
    1309             : 
    1310   540656235 :     z2e -= 2; z2d = z2e;
    1311   540656235 :     *z2d = addll(mulll(p1, *--yd), *z2d); z2d--;
    1312  5115791880 :     while (yd > x)
    1313             :     {
    1314  4034479410 :       hiremainder += overflow;
    1315  4034479410 :       *z2d = addll(addmul(p1, *--yd), *z2d); z2d--;
    1316             :     }
    1317   540656235 :     *--zd = hiremainder + overflow;
    1318             :   }
    1319             :   /* multiply zd by 2 (put result in zd - 1) */
    1320    81603270 :   zd[-1] = ((*zd & HIGHBIT) != 0);
    1321    81603270 :   shift_left(zd, zd, 0, (nx<<1)-3, 0, 1);
    1322             : 
    1323             :   /* add the squares */
    1324    81603270 :   xd = x + nx; zd = z0 + lz;
    1325    81603270 :   p1 = *--xd;
    1326    81603270 :   zd--; *zd = mulll(p1,p1);
    1327    81603270 :   zd--; *zd = addll(hiremainder, *zd);
    1328   785466045 :   while (xd > x)
    1329             :   {
    1330   622259505 :     p1 = *--xd;
    1331   622259505 :     zd--; *zd = addll(mulll(p1,p1)+ overflow, *zd);
    1332   622259505 :     zd--; *zd = addll(hiremainder + overflow, *zd);
    1333             :   }
    1334             : 
    1335             : END:
    1336    81603270 :   if (*zd == 0) { zd++; lz--; } /* normalize */
    1337    81603270 :   *--zd = evalsigne(1) | evallgefint(lz);
    1338    81603270 :   *--zd = evaltyp(t_INT) | evallg(lz);
    1339    81603270 :   avma=(pari_sp)zd; return zd;
    1340             : }
    1341             : 
    1342             : /********************************************************************/
    1343             : /**                                                                **/
    1344             : /**               INTEGER MULTIPLICATION (FFT)                     **/
    1345             : /**                                                                **/
    1346             : /********************************************************************/
    1347             : 
    1348             : /*
    1349             :  Compute parameters for FFT:
    1350             :    len: result length
    1351             :    k: FFT depth.
    1352             :    n: number of blocks (2^k)
    1353             :    bs: block size
    1354             :    mod: Modulus is M=2^(BIL*mod)+1
    1355             :    ord: order of 2 in Z/MZ.
    1356             :  We must have:
    1357             :    bs*n >= l
    1358             :    2^(BIL*mod) > nb*2^(2*BIL*bs)
    1359             :    2^k | 2*BIL*mod
    1360             : */
    1361             : static void
    1362        5298 : mulliifft_params(long len, long *k, long *mod, long *bs, long *n, ulong *ord)
    1363             : {
    1364             :   long r;
    1365        5298 :   *k = expu((3*len)>>2)-3;
    1366             :   do {
    1367        5298 :     (*k)--;
    1368        5298 :     r  = *k-(TWOPOTBITS_IN_LONG+2);
    1369        5298 :     *n = 1L<<*k;
    1370        5298 :     *bs = (len+*n-1)>>*k;
    1371        5298 :     *mod= 2**bs+1;
    1372        5298 :     if (r>0)
    1373         240 :       *mod=((*mod+(1L<<r)-1)>>r)<<r;
    1374        5298 :   } while(*mod>=3**bs);
    1375        5298 :   *ord= 4**mod*BITS_IN_LONG;
    1376        5298 : }
    1377             : 
    1378             : /* Zf_: arithmetic in ring Z/MZ where M= 2^(BITS_IN_LONG*mod)+1
    1379             :  * for some mod.
    1380             :  * Do not garbage collect.
    1381             :  */
    1382             : 
    1383             : static GEN
    1384    10051008 : Zf_add(GEN a, GEN b, GEN M)
    1385             : {
    1386    10051008 :   GEN y, z = addii(a,b);
    1387    10051008 :   long mod = lgefint(M)-3;
    1388    10051008 :   long l = NLIMBS(z);
    1389    10051008 :   if (l<=mod) return z;
    1390     3366957 :   y = subis(z, 1);
    1391     3366957 :   if (NLIMBS(y)<=mod) return z;
    1392     3366957 :   return int_normalize(y,1);
    1393             : }
    1394             : 
    1395             : static GEN
    1396    10209636 : Zf_sub(GEN a, GEN b, GEN M)
    1397             : {
    1398    10209636 :   GEN z = subii(a,b);
    1399    10209636 :   return signe(z)>=0? z: addii(M,z);
    1400             : }
    1401             : 
    1402             : /* destroy z */
    1403             : static GEN
    1404    19798509 : Zf_red_destroy(GEN z, GEN M)
    1405             : {
    1406    19798509 :   long mod = lgefint(M)-3;
    1407    19798509 :   long l = NLIMBS(z);
    1408             :   GEN y;
    1409    19798509 :   if (l<=mod) return z;
    1410     8593407 :   y = shifti(z, -mod*BITS_IN_LONG);
    1411     8593407 :   z = int_normalize(z, NLIMBS(y));
    1412     8593407 :   y = Zf_red_destroy(y, M);
    1413     8593407 :   z = subii(z, y);
    1414     8593407 :   if (signe(z)<0) z = addii(z, M);
    1415     8593407 :   return z;
    1416             : }
    1417             : 
    1418             : INLINE GEN
    1419    10278894 : Zf_shift(GEN a, ulong s, GEN M) { return Zf_red_destroy(shifti(a, s), M); }
    1420             : 
    1421             : /*
    1422             :  Multiply by sqrt(2)^s
    1423             :  We use the formula sqrt(2)=z_8*(1-z_4)) && z_8=2^(ord/16) [2^(ord/4)+1]
    1424             : */
    1425             : 
    1426             : static GEN
    1427    10051008 : Zf_mulsqrt2(GEN a, ulong s, ulong ord, GEN M)
    1428             : {
    1429    10051008 :   ulong hord = ord>>1;
    1430    10051008 :   if (!signe(a)) return gen_0;
    1431     9035430 :   if (odd(s)) /* Multiply by 2^(s/2) */
    1432             :   {
    1433      158628 :     GEN az8  = Zf_shift(a,   ord>>4, M);
    1434      158628 :     GEN az83 = Zf_shift(az8, ord>>3, M);
    1435      158628 :     a = Zf_sub(az8, az83, M);
    1436      158628 :     s--;
    1437             :   }
    1438     9035430 :   if (s < hord)
    1439     6648285 :     return Zf_shift(a, s>>1, M);
    1440             :   else
    1441     2387145 :     return subii(M,Zf_shift(a, (s-hord)>>1, M));
    1442             : }
    1443             : 
    1444             : INLINE GEN
    1445      157056 : Zf_sqr(GEN a, GEN M) { return Zf_red_destroy(sqri(a), M); }
    1446             : 
    1447             : INLINE GEN
    1448      769152 : Zf_mul(GEN a, GEN b, GEN M) { return Zf_red_destroy(mulii(a,b), M); }
    1449             : 
    1450             : /* In place, bit reversing FFT */
    1451             : static void
    1452     1685748 : muliifft_dit(ulong o, ulong ord, GEN M, GEN FFT, long d, long step)
    1453             : {
    1454     1685748 :   pari_sp av = avma;
    1455             :   long i;
    1456     1685748 :   ulong j, no = (o<<1)%ord;
    1457     1685748 :   long hstep=step>>1;
    1458     8190132 :   for (i = d+1, j = 0; i <= d+hstep; ++i, j =(j+o)%ord)
    1459             :   {
    1460     6504384 :     GEN a = Zf_add(gel(FFT,i), gel(FFT,i+hstep), M);
    1461     6504384 :     GEN b = Zf_mulsqrt2(Zf_sub(gel(FFT,i), gel(FFT,i+hstep), M), j, ord, M);
    1462     6504384 :     affii(a,gel(FFT,i));
    1463     6504384 :     affii(b,gel(FFT,i+hstep));
    1464     6504384 :     avma = av;
    1465             :   }
    1466     1685748 :   if (hstep>1)
    1467             :   {
    1468      838068 :     muliifft_dit(no, ord, M, FFT, d, hstep);
    1469      838068 :     muliifft_dit(no, ord, M, FFT, d+hstep, hstep);
    1470             :   }
    1471     1685748 : }
    1472             : 
    1473             : /* In place, bit reversed FFT, inverse of muliifft_dit */
    1474             : static void
    1475      920910 : muliifft_dis(ulong o, ulong ord, GEN M, GEN FFT, long d, long step)
    1476             : {
    1477      920910 :   pari_sp av = avma;
    1478             :   long i;
    1479      920910 :   ulong j, no = (o<<1)%ord;
    1480      920910 :   long hstep=step>>1;
    1481      920910 :   if (hstep>1)
    1482             :   {
    1483      457806 :     muliifft_dis(no, ord, M, FFT, d, hstep);
    1484      457806 :     muliifft_dis(no, ord, M, FFT, d+hstep, hstep);
    1485             :   }
    1486     4467534 :   for (i = d+1, j = 0; i <= d+hstep; ++i, j =(j+o)%ord)
    1487             :   {
    1488     3546624 :     GEN z = Zf_mulsqrt2(gel(FFT,i+hstep), j, ord, M);
    1489     3546624 :     GEN a = Zf_add(gel(FFT,i), z, M);
    1490     3546624 :     GEN b = Zf_sub(gel(FFT,i), z, M);
    1491     3546624 :     affii(a,gel(FFT,i));
    1492     3546624 :     affii(b,gel(FFT,i+hstep));
    1493     3546624 :     avma = av;
    1494             :   }
    1495      920910 : }
    1496             : 
    1497             : static GEN
    1498        9612 : muliifft_spliti(GEN a, long na, long bs, long n, long mod)
    1499             : {
    1500        9612 :   GEN ap = a+na-1;
    1501        9612 :   GEN c  = cgetg(n+1, t_VEC);
    1502             :   long i,j;
    1503     1704972 :   for(i=1;i<=n;i++)
    1504             :   {
    1505     1695360 :     GEN z = cgeti(mod+3);
    1506     1695360 :     if (na)
    1507             :     {
    1508      835131 :       long m = minss(bs, na), v=0;
    1509      835131 :       GEN zp, aa=ap-m+1;
    1510      835131 :       while (!*aa && v<m) {aa++; v++;}
    1511      835131 :       zp = z+m-v+1;
    1512    19260318 :       for (j=v; j < m; j++)
    1513    18425187 :         *zp-- = *ap--;
    1514      835131 :       ap -= v; na -= m;
    1515      835131 :       z[1] = evalsigne(m!=v) | evallgefint(m-v+2);
    1516             :     } else
    1517      860229 :       z[1] = evalsigne(0) | evallgefint(2);
    1518     1695360 :     gel(c, i) = z;
    1519             :   }
    1520        9612 :   return c;
    1521             : }
    1522             : 
    1523             : static GEN
    1524        5298 : muliifft_unspliti(GEN V, long bs, long len)
    1525             : {
    1526        5298 :   long s, i, j, l = lg(V);
    1527        5298 :   GEN a = cgeti(len);
    1528        5298 :   a[1] = evalsigne(1)|evallgefint(len);
    1529    29656521 :   for(i=2;i<len;i++)
    1530    29651223 :     a[i] = 0;
    1531      931506 :   for(i=1, s=0; i<l; i++, s+=bs)
    1532             :   {
    1533      926208 :     GEN  u = gel(V,i);
    1534      926208 :     if (signe(u))
    1535             :     {
    1536      793014 :       GEN ap = int_W(a,s);
    1537      793014 :       GEN up = int_LSW(u);
    1538      793014 :       long lu = NLIMBS(u);
    1539             :       LOCAL_OVERFLOW;
    1540      793014 :       *ap = addll(*ap, *up--); ap--;
    1541    48552153 :       for (j=1; j<lu; j++)
    1542    47759139 :        { *ap = addllx(*ap, *up--); ap--; }
    1543     1586190 :       while (overflow)
    1544         162 :        { *ap = addll(*ap, 1); ap--; }
    1545             :     }
    1546             :   }
    1547        5298 :   return int_normalize(a,0);
    1548             : }
    1549             : 
    1550             : static GEN
    1551         984 : sqrispec_fft(GEN a, long na)
    1552             : {
    1553         984 :   pari_sp av, ltop = avma;
    1554         984 :   long len = 2*na;
    1555             :   long k, mod, bs, n;
    1556             :   GEN  FFT, M;
    1557             :   long i;
    1558             :   ulong o, ord;
    1559             : 
    1560         984 :   mulliifft_params(len,&k,&mod,&bs,&n,&ord);
    1561         984 :   o = ord>>k;
    1562         984 :   M = int2n(mod*BITS_IN_LONG);
    1563         984 :   M[2+mod] = 1;
    1564         984 :   FFT = muliifft_spliti(a, na, bs, n, mod);
    1565         984 :   muliifft_dit(o, ord, M, FFT, 0, n);
    1566         984 :   av = avma;
    1567      158040 :   for(i=1; i<=n; i++)
    1568             :   {
    1569      157056 :     affii(Zf_sqr(gel(FFT,i), M), gel(FFT,i));
    1570      157056 :     avma=av;
    1571             :   }
    1572         984 :   muliifft_dis(ord-o, ord, M, FFT, 0, n);
    1573      158040 :   for(i=1; i<=n; i++)
    1574             :   {
    1575      157056 :     affii(Zf_shift(gel(FFT,i), (ord>>1)-k, M), gel(FFT,i));
    1576      157056 :     avma=av;
    1577             :   }
    1578         984 :   return gerepileuptoint(ltop, muliifft_unspliti(FFT,bs,2+len));
    1579             : }
    1580             : 
    1581             : static GEN
    1582        4314 : muliispec_fft(GEN a, GEN b, long na, long nb)
    1583             : {
    1584        4314 :   pari_sp av, av2, ltop = avma;
    1585        4314 :   long len = na+nb;
    1586             :   long k, mod, bs, n;
    1587             :   GEN FFT, FFTb, M;
    1588             :   long i;
    1589             :   ulong o, ord;
    1590             : 
    1591        4314 :   mulliifft_params(len,&k,&mod,&bs,&n,&ord);
    1592        4314 :   o = ord>>k;
    1593        4314 :   M = int2n(mod*BITS_IN_LONG);
    1594        4314 :   M[2+mod] = 1;
    1595        4314 :   FFT = muliifft_spliti(a, na, bs, n, mod);
    1596        4314 :   av=avma;
    1597        4314 :   muliifft_dit(o, ord, M, FFT, 0, n);
    1598        4314 :   FFTb = muliifft_spliti(b, nb, bs, n, mod);
    1599        4314 :   av2 = avma;
    1600        4314 :   muliifft_dit(o, ord, M, FFTb, 0, n);
    1601      773466 :   for(i=1; i<=n; i++)
    1602             :   {
    1603      769152 :     affii(Zf_mul(gel(FFT,i), gel(FFTb,i), M), gel(FFT,i));
    1604      769152 :     avma=av2;
    1605             :   }
    1606        4314 :   avma=av;
    1607        4314 :   muliifft_dis(ord-o, ord, M, FFT, 0, n);
    1608      773466 :   for(i=1; i<=n; i++)
    1609             :   {
    1610      769152 :     affii(Zf_shift(gel(FFT,i),(ord>>1)-k,M), gel(FFT,i));
    1611      769152 :     avma=av;
    1612             :   }
    1613        4314 :   return gerepileuptoint(ltop, muliifft_unspliti(FFT,bs,2+len));
    1614             : }
    1615             : 
    1616             : /********************************************************************/
    1617             : /**                                                                **/
    1618             : /**               INTEGER MULTIPLICATION (KARATSUBA)               **/
    1619             : /**                                                                **/
    1620             : /********************************************************************/
    1621             : 
    1622             : /* return (x shifted left d words) + y. Assume d > 0, x > 0 and y >= 0 */
    1623             : static GEN
    1624   102783729 : addshiftw(GEN x, GEN y, long d)
    1625             : {
    1626   102783729 :   GEN z,z0,y0,yd, zd = (GEN)avma;
    1627   102783729 :   long a,lz,ly = lgefint(y);
    1628             : 
    1629   102783729 :   z0 = new_chunk(d);
    1630   102783729 :   a = ly-2; yd = y+ly;
    1631   102783729 :   if (a >= d)
    1632             :   {
    1633   100675794 :     y0 = yd-d; while (yd > y0) *--zd = *--yd; /* copy last d words of y */
    1634   100675794 :     a -= d;
    1635   100675794 :     if (a)
    1636    68042949 :       z = addiispec(LIMBS(x), LIMBS(y), NLIMBS(x), a);
    1637             :     else
    1638    32632845 :       z = icopy(x);
    1639             :   }
    1640             :   else
    1641             :   {
    1642     2107935 :     y0 = yd-a; while (yd > y0) *--zd = *--yd; /* copy last a words of y */
    1643     2107935 :     while (zd > z0) *--zd = 0;    /* complete with 0s */
    1644     2107935 :     z = icopy(x);
    1645             :   }
    1646   102783729 :   lz = lgefint(z)+d;
    1647   102783729 :   z[1] = evalsigne(1) | evallgefint(lz);
    1648   102783729 :   z[0] = evaltyp(t_INT) | evallg(lz); return z;
    1649             : }
    1650             : 
    1651             : /* Fast product (Karatsuba) of integers. a and b are "special" GENs
    1652             :  * c,c0,c1,c2 are genuine GENs.
    1653             :  */
    1654             : GEN
    1655  1526683349 : muliispec(GEN a, GEN b, long na, long nb)
    1656             : {
    1657             :   GEN a0,c,c0;
    1658             :   long n0, n0a, i;
    1659             :   pari_sp av;
    1660             : 
    1661  1526683349 :   if (na < nb) swapspec(a,b, na,nb);
    1662  1526683349 :   if (nb < MULII_KARATSUBA_LIMIT) return muliispec_basecase(a,b,na,nb);
    1663    31925856 :   if (nb >= MULII_FFT_LIMIT)      return muliispec_fft(a,b,na,nb);
    1664    31921542 :   i=(na>>1); n0=na-i; na=i;
    1665    31921542 :   av=avma; a0=a+na; n0a=n0;
    1666    31921542 :   while (n0a && !*a0) { a0++; n0a--; }
    1667             : 
    1668    31921542 :   if (n0a && nb > n0)
    1669    30694299 :   { /* nb <= na <= n0 */
    1670             :     GEN b0,c1,c2;
    1671             :     long n0b;
    1672             : 
    1673    30694299 :     nb -= n0;
    1674    30694299 :     c = muliispec(a,b,na,nb);
    1675    30694299 :     b0 = b+nb; n0b = n0;
    1676    30694299 :     while (n0b && !*b0) { b0++; n0b--; }
    1677    30694299 :     if (n0b)
    1678             :     {
    1679    30059400 :       c0 = muliispec(a0,b0, n0a,n0b);
    1680             : 
    1681    30059400 :       c2 = addiispec(a0,a, n0a,na);
    1682    30059400 :       c1 = addiispec(b0,b, n0b,nb);
    1683    30059400 :       c1 = muliispec(LIMBS(c1),LIMBS(c2), NLIMBS(c1),NLIMBS(c2));
    1684    30059400 :       c2 = addiispec(LIMBS(c0),LIMBS(c),  NLIMBS(c0),NLIMBS(c));
    1685             : 
    1686    30059400 :       c1 = subiispec(LIMBS(c1),LIMBS(c2), NLIMBS(c1),NLIMBS(c2));
    1687             :     }
    1688             :     else
    1689             :     {
    1690      634899 :       c0 = gen_0;
    1691      634899 :       c1 = muliispec(a0,b, n0a,nb);
    1692             :     }
    1693    30694299 :     c = addshiftw(c,c1, n0);
    1694             :   }
    1695             :   else
    1696             :   {
    1697     1227243 :     c = muliispec(a,b,na,nb);
    1698     1227243 :     c0 = muliispec(a0,b,n0a,nb);
    1699             :   }
    1700    31921542 :   return gerepileuptoint(av, addshiftw(c,c0, n0));
    1701             : }
    1702             : GEN
    1703      411747 : muluui(ulong x, ulong y, GEN z)
    1704             : {
    1705      411747 :   long t, s = signe(z);
    1706             :   GEN r;
    1707             :   LOCAL_HIREMAINDER;
    1708             : 
    1709      411747 :   if (!x || !y || !signe(z)) return gen_0;
    1710      411450 :   t = mulll(x,y);
    1711      411450 :   if (!hiremainder)
    1712      411450 :     r = muluispec(t, z+2, lgefint(z)-2);
    1713             :   else
    1714             :   {
    1715             :     long tmp[2];
    1716           0 :     tmp[0] = hiremainder;
    1717           0 :     tmp[1] = t;
    1718           0 :     r = muliispec(z+2,tmp,lgefint(z)-2,2);
    1719             :   }
    1720      411450 :   setsigne(r,s); return r;
    1721             : }
    1722             : 
    1723             : #define sqrispec_mirror sqrispec
    1724             : #define muliispec_mirror muliispec
    1725             : 
    1726             : /* x % (2^n), assuming n >= 0 */
    1727             : GEN
    1728     5993247 : remi2n(GEN x, long n)
    1729             : {
    1730     5993247 :   long hi,l,k,lx,ly, sx = signe(x);
    1731             :   GEN z, xd, zd;
    1732             : 
    1733     5993247 :   if (!sx || !n) return gen_0;
    1734             : 
    1735     5975097 :   k = dvmdsBIL(n, &l);
    1736     5975097 :   lx = lgefint(x);
    1737     5975097 :   if (lx < k+3) return icopy(x);
    1738             : 
    1739     5941521 :   xd = x + (lx-k-1);
    1740             :   /* x = |_|...|#|1|...|k| : copy the last l bits of # and the last k words
    1741             :    *            ^--- initial xd  */
    1742     5941521 :   hi = ((ulong)*xd) & ((1UL<<l)-1); /* last l bits of # = top bits of result */
    1743     5941521 :   if (!hi)
    1744             :   { /* strip leading zeroes from result */
    1745      252699 :     xd++; while (k && !*xd) { k--; xd++; }
    1746      252699 :     if (!k) return gen_0;
    1747       53946 :     ly = k+2; xd--;
    1748             :   }
    1749             :   else
    1750     5688822 :     ly = k+3;
    1751             : 
    1752     5742768 :   zd = z = cgeti(ly);
    1753     5742768 :   *++zd = evalsigne(sx) | evallgefint(ly);
    1754     5742768 :   if (hi) *++zd = hi;
    1755     5742768 :   for ( ;k; k--) *++zd = *++xd;
    1756     5742768 :   return z;
    1757             : }
    1758             : 
    1759             : GEN
    1760   515965566 : sqrispec(GEN a, long na)
    1761             : {
    1762             :   GEN a0,c;
    1763             :   long n0, n0a, i;
    1764             :   pari_sp av;
    1765             : 
    1766   515965566 :   if (na < SQRI_KARATSUBA_LIMIT) return sqrispec_basecase(a,na);
    1767     3541578 :   if (na >= SQRI_FFT_LIMIT) return sqrispec_fft(a,na);
    1768     3540594 :   i=(na>>1); n0=na-i; na=i;
    1769     3540594 :   av=avma; a0=a+na; n0a=n0;
    1770     3540594 :   while (n0a && !*a0) { a0++; n0a--; }
    1771     3540594 :   c = sqrispec(a,na);
    1772     3540594 :   if (n0a)
    1773             :   {
    1774     3525411 :     GEN t, c1, c0 = sqrispec(a0,n0a);
    1775             : #if 0
    1776             :     c1 = shifti(muliispec(a0,a, n0a,na),1);
    1777             : #else /* faster */
    1778     3525411 :     t = addiispec(a0,a,n0a,na);
    1779     3525411 :     t = sqrispec(LIMBS(t),NLIMBS(t));
    1780     3525411 :     c1= addiispec(LIMBS(c0),LIMBS(c), NLIMBS(c0), NLIMBS(c));
    1781     3525411 :     c1= subiispec(LIMBS(t),LIMBS(c1), NLIMBS(t), NLIMBS(c1));
    1782             : #endif
    1783     3525411 :     c = addshiftw(c,c1, n0);
    1784     3525411 :     c = addshiftw(c,c0, n0);
    1785             :   }
    1786             :   else
    1787       15183 :     c = addshiftw(c,gen_0,n0<<1);
    1788     3540594 :   return gerepileuptoint(av, c);
    1789             : }
    1790             : 
    1791             : /********************************************************************/
    1792             : /**                                                                **/
    1793             : /**                    KARATSUBA SQUARE ROOT                       **/
    1794             : /**      adapted from Paul Zimmermann's implementation of          **/
    1795             : /**      his algorithm in GMP (mpn_sqrtrem)                        **/
    1796             : /**                                                                **/
    1797             : /********************************************************************/
    1798             : 
    1799             : /* Square roots table */
    1800             : static const unsigned char approx_tab[192] = {
    1801             :   128,128,129,130,131,132,133,134,135,136,137,138,139,140,141,142,
    1802             :   143,144,144,145,146,147,148,149,150,150,151,152,153,154,155,155,
    1803             :   156,157,158,159,160,160,161,162,163,163,164,165,166,167,167,168,
    1804             :   169,170,170,171,172,173,173,174,175,176,176,177,178,178,179,180,
    1805             :   181,181,182,183,183,184,185,185,186,187,187,188,189,189,190,191,
    1806             :   192,192,193,193,194,195,195,196,197,197,198,199,199,200,201,201,
    1807             :   202,203,203,204,204,205,206,206,207,208,208,209,209,210,211,211,
    1808             :   212,212,213,214,214,215,215,216,217,217,218,218,219,219,220,221,
    1809             :   221,222,222,223,224,224,225,225,226,226,227,227,228,229,229,230,
    1810             :   230,231,231,232,232,233,234,234,235,235,236,236,237,237,238,238,
    1811             :   239,240,240,241,241,242,242,243,243,244,244,245,245,246,246,247,
    1812             :   247,248,248,249,249,250,250,251,251,252,252,253,253,254,254,255
    1813             : };
    1814             : 
    1815             : /* N[0], assume N[0] >= 2^(BIL-2).
    1816             :  * Return r,s such that s^2 + r = N, 0 <= r <= 2s */
    1817             : static void
    1818    17022767 : p_sqrtu1(ulong *N, ulong *ps, ulong *pr)
    1819             : {
    1820    17022767 :   ulong prec, r, s, q, u, n0 = N[0];
    1821             : 
    1822    17022767 :   q = n0 >> (BITS_IN_LONG - 8);
    1823             :   /* 2^6 = 64 <= q < 256 = 2^8 */
    1824    17022767 :   s = approx_tab[q - 64];                                /* 128 <= s < 255 */
    1825    17022767 :   r = (n0 >> (BITS_IN_LONG - 16)) - s * s;                /* r <= 2*s */
    1826    17022767 :   if (r > (s << 1)) { r -= (s << 1) | 1; s++; }
    1827             : 
    1828             :   /* 8-bit approximation from the high 8-bits of N[0] */
    1829    17022767 :   prec = 8;
    1830    17022767 :   n0 <<= 2 * prec;
    1831    68091068 :   while (2 * prec < BITS_IN_LONG)
    1832             :   { /* invariant: s has prec bits, and r <= 2*s */
    1833    34045534 :     r = (r << prec) + (n0 >> (BITS_IN_LONG - prec));
    1834    34045534 :     n0 <<= prec;
    1835    34045534 :     u = 2 * s;
    1836    34045534 :     q = r / u; u = r - q * u;
    1837    34045534 :     s = (s << prec) + q;
    1838    34045534 :     u = (u << prec) + (n0 >> (BITS_IN_LONG - prec));
    1839    34045534 :     q = q * q;
    1840    34045534 :     r = u - q;
    1841    34045534 :     if (u < q) { s--; r += (s << 1) | 1; }
    1842    34045534 :     n0 <<= prec;
    1843    34045534 :     prec = 2 * prec;
    1844             :   }
    1845    17022767 :   *ps = s;
    1846    17022767 :   *pr = r;
    1847    17022767 : }
    1848             : 
    1849             : /* N[0..1], assume N[0] >= 2^(BIL-2).
    1850             :  * Return 1 if remainder overflows, 0 otherwise */
    1851             : static int
    1852    15123123 : p_sqrtu2(ulong *N, ulong *ps, ulong *pr)
    1853             : {
    1854    15123123 :   ulong cc, qhl, r, s, q, u, n1 = N[1];
    1855             :   LOCAL_OVERFLOW;
    1856             : 
    1857    15123123 :   p_sqrtu1(N, &s, &r); /* r <= 2s */
    1858    15123123 :   qhl = 0; while (r >= s) { qhl++; r -= s; }
    1859             :   /* now r < s < 2^(BIL/2) */
    1860    15123123 :   r = (r << BITS_IN_HALFULONG) | (n1 >> BITS_IN_HALFULONG);
    1861    15123123 :   u = s << 1;
    1862    15123123 :   q = r / u; u = r - q * u;
    1863    15123123 :   q += (qhl & 1) << (BITS_IN_HALFULONG - 1);
    1864    15123123 :   qhl >>= 1;
    1865             :   /* (initial r)<<(BIL/2) + n1>>(BIL/2) = (qhl<<(BIL/2) + q) * 2s + u */
    1866    15123123 :   s = ((s + qhl) << BITS_IN_HALFULONG) + q;
    1867    15123123 :   cc = u >> BITS_IN_HALFULONG;
    1868    15123123 :   r = (u << BITS_IN_HALFULONG) | (n1 & LOWMASK);
    1869    15123123 :   r = subll(r, q * q);
    1870    15123123 :   cc -= overflow + qhl;
    1871             :   /* now subtract 2*q*2^(BIL/2) + 2^BIL if qhl is set */
    1872    15123123 :   if ((long)cc < 0)
    1873             :   {
    1874     3884685 :     if (s) {
    1875     3879948 :       r = addll(r, s);
    1876     3879948 :       cc += overflow;
    1877     3879948 :       s--;
    1878             :     } else {
    1879        4737 :       cc++;
    1880        4737 :       s = ~0UL;
    1881             :     }
    1882     3884685 :     r = addll(r, s);
    1883     3884685 :     cc += overflow;
    1884             :   }
    1885    15123123 :   *ps = s;
    1886    15123123 :   *pr = r; return cc;
    1887             : }
    1888             : 
    1889             : static void
    1890    15037893 : xmpn_zero(GEN x, long n)
    1891             : {
    1892    15037893 :   while (--n >= 0) x[n]=0;
    1893    15037893 : }
    1894             : static void
    1895   154811331 : xmpn_copy(GEN z, GEN x, long n)
    1896             : {
    1897   154811331 :   long k = n;
    1898   154811331 :   while (--k >= 0) z[k] = x[k];
    1899   154811331 : }
    1900             : /* a[0..la-1] * 2^(lb BIL) | b[0..lb-1] */
    1901             : static GEN
    1902    66203766 : catii(GEN a, long la, GEN b, long lb)
    1903             : {
    1904    66203766 :   long l = la + lb + 2;
    1905    66203766 :   GEN z = cgetipos(l);
    1906    66203766 :   xmpn_copy(LIMBS(z), a, la);
    1907    66203766 :   xmpn_copy(LIMBS(z) + la, b, lb);
    1908    66203766 :   return int_normalize(z, 0);
    1909             : }
    1910             : 
    1911             : /* sqrt n[0..1], assume n normalized */
    1912             : static GEN
    1913    15078948 : sqrtispec2(GEN n, GEN *pr)
    1914             : {
    1915             :   ulong s, r;
    1916    15078948 :   int hi = p_sqrtu2((ulong*)n, &s, &r);
    1917    15078948 :   GEN S = utoi(s);
    1918    15078948 :   *pr = hi? uutoi(1,r): utoi(r);
    1919    15078948 :   return S;
    1920             : }
    1921             : 
    1922             : /* sqrt n[0], _dont_ assume n normalized */
    1923             : static GEN
    1924     1899644 : sqrtispec1_sh(GEN n, GEN *pr)
    1925             : {
    1926             :   GEN S;
    1927     1899644 :   ulong r, s, u0 = uel(n,0);
    1928     1899644 :   int sh = bfffo(u0) & ~1UL;
    1929     1899644 :   if (sh) u0 <<= sh;
    1930     1899644 :   p_sqrtu1(&u0, &s, &r);
    1931             :   /* s^2 + r = u0, s < 2^(BIL/2). Rescale back:
    1932             :    * 2^(2k) n = S^2 + R
    1933             :    * so 2^(2k) n = (S - s0)^2 + (2*S*s0 - s0^2 + R), s0 = S mod 2^k. */
    1934     1899644 :   if (sh) {
    1935     1895204 :     int k = sh >> 1;
    1936     1895204 :     ulong s0 = s & ((1L<<k) - 1);
    1937     1895204 :     r += s * (s0<<1);
    1938     1895204 :     s >>= k;
    1939     1895204 :     r >>= sh;
    1940             :   }
    1941     1899644 :   S = utoi(s);
    1942     1899644 :   if (pr) *pr = utoi(r);
    1943     1899644 :   return S;
    1944             : }
    1945             : 
    1946             : /* sqrt n[0..1], _dont_ assume n normalized */
    1947             : static GEN
    1948       44175 : sqrtispec2_sh(GEN n, GEN *pr)
    1949             : {
    1950             :   GEN S;
    1951       44175 :   ulong U[2], r, s, u0 = uel(n,0), u1 = uel(n,1);
    1952       44175 :   int hi, sh = bfffo(u0) & ~1UL;
    1953       44175 :   if (sh) {
    1954       43896 :     u0 = (u0 << sh) | (u1 >> (BITS_IN_LONG-sh));
    1955       43896 :     u1 <<= sh;
    1956             :   }
    1957       44175 :   U[0] = u0;
    1958       44175 :   U[1] = u1; hi = p_sqrtu2(U, &s, &r);
    1959             :   /* s^2 + R = u0|u1. Rescale back:
    1960             :    * 2^(2k) n = S^2 + R
    1961             :    * so 2^(2k) n = (S - s0)^2 + (2*S*s0 - s0^2 + R), s0 = S mod 2^k. */
    1962       44175 :   if (sh) {
    1963       43896 :     int k = sh >> 1;
    1964       43896 :     ulong s0 = s & ((1L<<k) - 1);
    1965             :     LOCAL_HIREMAINDER;
    1966             :     LOCAL_OVERFLOW;
    1967       43896 :     r = addll(r, mulll(s, (s0<<1)));
    1968       43896 :     if (overflow) hiremainder++;
    1969       43896 :     hiremainder += hi; /* + 0 or 1 */
    1970       43896 :     s >>= k;
    1971       43896 :     r = (r>>sh) | (hiremainder << (BITS_IN_LONG-sh));
    1972       43896 :     hi = (hiremainder & (1L<<sh));
    1973             :   }
    1974       44175 :   S = utoi(s);
    1975       44175 :   if (pr) *pr = hi? uutoi(1,r): utoi(r);
    1976       44175 :   return S;
    1977             : }
    1978             : 
    1979             : /* Let N = N[0..2n-1]. Return S (and set R) s.t S^2 + R = N, 0 <= R <= 2S
    1980             :  * Assume N normalized */
    1981             : static GEN
    1982    48180831 : sqrtispec(GEN N, long n, GEN *r)
    1983             : {
    1984             :   GEN S, R, q, z, u;
    1985             :   long l, h;
    1986             : 
    1987    48180831 :   if (n == 1) return sqrtispec2(N, r);
    1988    33101883 :   l = n >> 1;
    1989    33101883 :   h = n - l; /* N = a3(h) | a2(h) | a1(l) | a0(l words) */
    1990    33101883 :   S = sqrtispec(N, h, &R); /* S^2 + R = a3|a2 */
    1991             : 
    1992    33101883 :   z = catii(LIMBS(R), NLIMBS(R), N + 2*h, l); /* = R | a1(l) */
    1993    33101883 :   q = dvmdii(z, shifti(S,1), &u);
    1994    33101883 :   z = catii(LIMBS(u), NLIMBS(u), N + n + h, l); /* = u | a0(l) */
    1995             : 
    1996    33101883 :   S = addshiftw(S, q, l);
    1997    33101883 :   R = subii(z, sqri(q));
    1998    33101883 :   if (signe(R) < 0)
    1999             :   {
    2000     5287533 :     GEN S2 = shifti(S,1);
    2001     5287533 :     R = addis(subiispec(LIMBS(S2),LIMBS(R), NLIMBS(S2),NLIMBS(R)), -1);
    2002     5287533 :     S = addis(S, -1);
    2003             :   }
    2004    33101883 :   *r = R; return S;
    2005             : }
    2006             : 
    2007             : /* Return S (and set R) s.t S^2 + R = N, 0 <= R <= 2S.
    2008             :  * As for dvmdii, R is last on stack and guaranteed to be gen_0 in case the
    2009             :  * remainder is 0. R = NULL is allowed. */
    2010             : GEN
    2011     1984874 : sqrtremi(GEN N, GEN *r)
    2012             : {
    2013             :   pari_sp av;
    2014     1984874 :   GEN S, R, n = N+2;
    2015     1984874 :   long k, l2, ln = NLIMBS(N);
    2016             :   int sh;
    2017             : 
    2018     1984874 :   if (ln <= 2)
    2019             :   {
    2020     1943819 :     if (ln == 2) return sqrtispec2_sh(n, r);
    2021     1899644 :     if (ln == 1) return sqrtispec1_sh(n, r);
    2022           0 :     if (r) *r = gen_0;
    2023           0 :     return gen_0;
    2024             :   }
    2025       41055 :   av = avma;
    2026       41055 :   sh = bfffo(n[0]) >> 1;
    2027       41055 :   l2 = (ln + 1) >> 1;
    2028       81774 :   if (sh || (ln & 1)) { /* normalize n, so that n[0] >= 2^BIL / 4 */
    2029       40719 :     GEN s0, t = new_chunk(ln + 1);
    2030       40719 :     t[ln] = 0;
    2031       40719 :     if (sh)
    2032       40284 :       shift_left(t, n, 0,ln-1, 0, sh << 1);
    2033             :     else
    2034         435 :       xmpn_copy(t, n, ln);
    2035       40719 :     S = sqrtispec(t, l2, &R); /* t normalized, 2 * l2 words */
    2036             :     /* Rescale back:
    2037             :      * 2^(2k) n = S^2 + R, k = sh + (ln & 1)*BIL/2
    2038             :      * so 2^(2k) n = (S - s0)^2 + (2*S*s0 - s0^2 + R), s0 = S mod 2^k. */
    2039       40719 :     k = sh + (ln & 1) * (BITS_IN_LONG/2);
    2040       40719 :     s0 = remi2n(S, k);
    2041       40719 :     R = addii(shifti(R,-1), mulii(s0, S));
    2042       40719 :     R = shifti(R, 1 - (k<<1));
    2043       40719 :     S = shifti(S, -k);
    2044             :   }
    2045             :   else
    2046         336 :     S = sqrtispec(n, l2, &R);
    2047             : 
    2048       41055 :   if (!r) { avma = (pari_sp)S; return gerepileuptoint(av, S); }
    2049        2058 :   gerepileall(av, 2, &S, &R); *r = R; return S;
    2050             : }
    2051             : 
    2052             : /* compute sqrt(|a|), assuming a != 0 */
    2053             : 
    2054             : #if 1
    2055             : GEN
    2056    15037893 : sqrtr_abs(GEN x)
    2057             : {
    2058    15037893 :   long l = realprec(x) - 2, e = expo(x), er = e>>1;
    2059    15037893 :   GEN b, c, res = cgetr(2 + l);
    2060    15037893 :   res[1] = evalsigne(1) | evalexpo(er);
    2061    15037893 :   if (e&1) {
    2062     7365471 :     b = new_chunk(l << 1);
    2063     7365471 :     xmpn_copy(b, x+2, l);
    2064     7365471 :     xmpn_zero(b + l,l);
    2065     7365471 :     b = sqrtispec(b, l, &c);
    2066     7365471 :     xmpn_copy(res+2, b+2, l);
    2067     7365471 :     if (cmpii(c, b) > 0) roundr_up_ip(res, l+2);
    2068             :   } else {
    2069             :     ulong u;
    2070     7672422 :     b = new_chunk(2 + (l << 1));
    2071     7672422 :     shift_left(b+1, x+2, 0,l-1, 0, BITS_IN_LONG-1);
    2072     7672422 :     b[0] = uel(x,2)>>1;
    2073     7672422 :     xmpn_zero(b + l+1,l+1);
    2074     7672422 :     b = sqrtispec(b, l+1, &c);
    2075     7672422 :     xmpn_copy(res+2, b+2, l);
    2076     7672422 :     u = uel(b,l+2);
    2077     7672422 :     if ( u&HIGHBIT || (u == ~HIGHBIT && cmpii(c,b) > 0))
    2078     3741833 :       roundr_up_ip(res, l+2);
    2079             :   }
    2080    15037893 :   avma = (pari_sp)res; return res;
    2081             : }
    2082             : 
    2083             : #else /* use t_REAL: currently much slower (quadratic division) */
    2084             : 
    2085             : #ifdef LONG_IS_64BIT
    2086             : /* 64 bits of b = sqrt(a[0] * 2^64 + a[1])  [ up to 1ulp ] */
    2087             : static ulong
    2088             : sqrtu2(ulong *a)
    2089             : {
    2090             :   ulong c, b = dblmantissa( sqrt((double)a[0]) );
    2091             :   LOCAL_HIREMAINDER;
    2092             :   LOCAL_OVERFLOW;
    2093             : 
    2094             :   /* > 32 correct bits, 1 Newton iteration to reach 64 */
    2095             :   if (b <= a[0]) return HIGHBIT | (a[0] >> 1);
    2096             :   hiremainder = a[0]; c = divll(a[1], b);
    2097             :   return (addll(c, b) >> 1) | HIGHBIT;
    2098             : }
    2099             : /* 64 bits of sqrt(a[0] * 2^63) */
    2100             : static ulong
    2101             : sqrtu2_1(ulong *a)
    2102             : {
    2103             :   ulong t[2];
    2104             :   t[0] = (a[0] >> 1);
    2105             :   t[1] = (a[0] << (BITS_IN_LONG-1)) | (a[1] >> 1);
    2106             :   return sqrtu2(t);
    2107             : }
    2108             : #else
    2109             : /* 32 bits of sqrt(a[0] * 2^32) */
    2110             : static ulong
    2111             : sqrtu2(ulong *a)   { return dblmantissa( sqrt((double)a[0]) ); }
    2112             : /* 32 bits of sqrt(a[0] * 2^31) */
    2113             : static ulong
    2114             : sqrtu2_1(ulong *a) { return dblmantissa( sqrt(2. * a[0]) ); }
    2115             : #endif
    2116             : 
    2117             : GEN
    2118             : sqrtr_abs(GEN x)
    2119             : {
    2120             :   long l1, i, l = lg(x), ex = expo(x);
    2121             :   GEN a, t, y = cgetr(l);
    2122             :   pari_sp av, av0 = avma;
    2123             : 
    2124             :   a = rtor(x, l+1);
    2125             :   t = cgetr(l+1);
    2126             :   if (ex & 1) { /* odd exponent */
    2127             :     a[1] = evalsigne(1) | _evalexpo(1);
    2128             :     t[2] = (long)sqrtu2((ulong*)a + 2);
    2129             :   } else { /* even exponent */
    2130             :     a[1] = evalsigne(1) | _evalexpo(0);
    2131             :     t[2] = (long)sqrtu2_1((ulong*)a + 2);
    2132             :   }
    2133             :   t[1] = evalsigne(1) | _evalexpo(0);
    2134             :   for (i = 3; i <= l; i++) t[i] = 0;
    2135             : 
    2136             :   /* |x| = 2^(ex/2) a, t ~ sqrt(a) */
    2137             :   l--; l1 = 1; av = avma;
    2138             :   while (l1 < l) { /* let t := (t + a/t)/2 */
    2139             :     l1 <<= 1; if (l1 > l) l1 = l;
    2140             :     setlg(a, l1 + 2);
    2141             :     setlg(t, l1 + 2);
    2142             :     affrr(addrr(t, divrr(a,t)), t); shiftr_inplace(t, -1);
    2143             :     avma = av;
    2144             :   }
    2145             :   affrr(t,y); shiftr_inplace(y, (ex>>1));
    2146             :   avma = av0; return y;
    2147             : }
    2148             : 
    2149             : #endif
    2150             : 
    2151             : /*******************************************************************
    2152             :  *                                                                 *
    2153             :  *                           Base Conversion                       *
    2154             :  *                                                                 *
    2155             :  *******************************************************************/
    2156             : 
    2157             : static void
    2158     1071579 : convi_dac(GEN x, ulong l, ulong *res)
    2159             : {
    2160     1071579 :   pari_sp ltop=avma;
    2161             :   ulong m;
    2162             :   GEN x1,x2;
    2163     2143158 :   if (l==1) { *res=itou(x); return; }
    2164      490056 :   m=l>>1;
    2165      490056 :   x1=dvmdii(x,powuu(1000000000UL,m),&x2);
    2166      490056 :   convi_dac(x1,l-m,res+m);
    2167      490056 :   convi_dac(x2,m,res);
    2168      490056 :   avma=ltop;
    2169             : }
    2170             : 
    2171             : /* convert integer --> base 10^9 [not memory clean] */
    2172             : ulong *
    2173      239904 : convi(GEN x, long *l)
    2174             : {
    2175      239904 :   long lz, lx = lgefint(x);
    2176             :   ulong *z;
    2177      239904 :   if (lx == 3 && uel(x,2) < 1000000000UL) {
    2178      148437 :     z = (ulong*)new_chunk(1);
    2179      148437 :     *z = x[2];
    2180      148437 :     *l = 1; return z+1;
    2181             :   }
    2182       91467 :   lz = 1 + (long)bit_accuracy_mul(lx, LOG10_2/9);
    2183       91467 :   z = (ulong*)new_chunk(lz);
    2184       91467 :   convi_dac(x,(ulong)lz,z);
    2185       91467 :   while (z[lz-1]==0) lz--;
    2186       91467 :   *l=lz; return z+lz;
    2187             : }
    2188             : 

Generated by: LCOV version 1.11