Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - kernel/gmp - mp.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.12.1 lcov report (development 24038-ebe36f6c4) Lines: 648 717 90.4 %
Date: 2019-07-23 05:53:17 Functions: 54 55 98.2 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : #line 2 "../src/kernel/gmp/mp.c"
       2             : /* Copyright (C) 2002-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             : /**                               GMP KERNEL                          **/
      18             : /** BA2002Sep24                                                       **/
      19             : /***********************************************************************/
      20             : /* GMP t_INT as just like normal t_INT, just the mantissa is the other way
      21             :  * round
      22             :  *
      23             :  *   `How would you like to live in Looking-glass House, Kitty?  I
      24             :  *   wonder if they'd give you milk in there?  Perhaps Looking-glass
      25             :  *   milk isn't good to drink--But oh, Kitty! now we come to the
      26             :  *   passage.  You can just see a little PEEP of the passage in
      27             :  *   Looking-glass House, if you leave the door of our drawing-room
      28             :  *   wide open:  and it's very like our passage as far as you can see,
      29             :  *   only you know it may be quite different on beyond.  Oh, Kitty!
      30             :  *   how nice it would be if we could only get through into Looking-
      31             :  *   glass House!  I'm sure it's got, oh! such beautiful things in it!
      32             :  *
      33             :  *  Through the Looking Glass,  Lewis Carrol
      34             :  *
      35             :  *  (pityful attempt to beat GN code/comments rate)
      36             :  *  */
      37             : 
      38             : #include <gmp.h>
      39             : #include "pari.h"
      40             : #include "paripriv.h"
      41             : #include "../src/kernel/none/tune-gen.h"
      42             : 
      43             : /*We need PARI invmod renamed to invmod_pari*/
      44             : #define INVMOD_PARI
      45             : 
      46           0 : static void *pari_gmp_realloc(void *ptr, size_t old_size, size_t new_size) {
      47           0 :   (void)old_size; return (void *) pari_realloc(ptr,new_size);
      48             : }
      49             : 
      50      436566 : static void pari_gmp_free(void *ptr, size_t old_size){
      51      436566 :   (void)old_size; pari_free(ptr);
      52      436566 : }
      53             : 
      54             : static void *(*old_gmp_malloc)(size_t new_size);
      55             : static void *(*old_gmp_realloc)(void *ptr, size_t old_size, size_t new_size);
      56             : static void (*old_gmp_free)(void *ptr, size_t old_size);
      57             : 
      58             : void
      59         912 : pari_kernel_init(void)
      60             : {
      61         912 :   mp_get_memory_functions (&old_gmp_malloc, &old_gmp_realloc, &old_gmp_free);
      62         912 :   mp_set_memory_functions((void *(*)(size_t)) pari_malloc, pari_gmp_realloc, pari_gmp_free);
      63         912 : }
      64             : 
      65             : void
      66         904 : pari_kernel_close(void)
      67             : {
      68             :   void *(*new_gmp_malloc)(size_t new_size);
      69             :   void *(*new_gmp_realloc)(void *ptr, size_t old_size, size_t new_size);
      70             :   void (*new_gmp_free)(void *ptr, size_t old_size);
      71         904 :   mp_get_memory_functions (&new_gmp_malloc, &new_gmp_realloc, &new_gmp_free);
      72         904 :   if (new_gmp_malloc==pari_malloc) new_gmp_malloc = old_gmp_malloc;
      73         904 :   if (new_gmp_realloc==pari_gmp_realloc) new_gmp_realloc = old_gmp_realloc;
      74         904 :   if (new_gmp_free==pari_gmp_free) new_gmp_free = old_gmp_free;
      75         904 :   mp_set_memory_functions(new_gmp_malloc, new_gmp_realloc, new_gmp_free);
      76         904 : }
      77             : 
      78             : #define LIMBS(x)  ((mp_limb_t *)((x)+2))
      79             : #define NLIMBS(x) (lgefint(x)-2)
      80             : /*This one is for t_REALs to emphasize they are not t_INTs*/
      81             : #define RLIMBS(x)  ((mp_limb_t *)((x)+2))
      82             : #define RNLIMBS(x) (lg(x)-2)
      83             : 
      84             : INLINE void
      85     1642043 : xmpn_copy(mp_limb_t *x, mp_limb_t *y, long n)
      86             : {
      87     1642043 :   while (--n >= 0) x[n]=y[n];
      88     1642043 : }
      89             : 
      90             : INLINE void
      91   189828352 : xmpn_mirror(mp_limb_t *x, long n)
      92             : {
      93             :   long i;
      94  1182578515 :   for(i=0;i<(n>>1);i++)
      95             :   {
      96   992750163 :     ulong m=x[i];
      97   992750163 :     x[i]=x[n-1-i];
      98   992750163 :     x[n-1-i]=m;
      99             :   }
     100   189828352 : }
     101             : 
     102             : INLINE void
     103   260198262 : xmpn_mirrorcopy(mp_limb_t *z, mp_limb_t *x, long n)
     104             : {
     105             :   long i;
     106  3156114776 :   for(i=0;i<n;i++)
     107  2895916514 :     z[i]=x[n-1-i];
     108   260198262 : }
     109             : 
     110             : INLINE void
     111    83958718 : xmpn_zero(mp_ptr x, mp_size_t n)
     112             : {
     113    83958718 :   while (--n >= 0) x[n]=0;
     114    83958718 : }
     115             : 
     116             : INLINE GEN
     117    17058777 : icopy_ef(GEN x, long l)
     118             : {
     119    17058777 :   register long lx = lgefint(x);
     120    17058777 :   const GEN y = cgeti(l);
     121             : 
     122    17058716 :   while (--lx > 0) y[lx]=x[lx];
     123    17058716 :   return y;
     124             : }
     125             : 
     126             : /* NOTE: arguments of "spec" routines (muliispec, addiispec, etc.) aren't
     127             :  * GENs but pairs (long *a, long na) representing a list of digits (in basis
     128             :  * BITS_IN_LONG) : a[0], ..., a[na-1]. [ In ordre to facilitate splitting: no
     129             :  * need to reintroduce codewords ]
     130             :  * Use speci(a,na) to visualize the corresponding GEN.
     131             :  */
     132             : 
     133             : /***********************************************************************/
     134             : /**                                                                   **/
     135             : /**                     ADDITION / SUBTRACTION                        **/
     136             : /**                                                                   **/
     137             : /***********************************************************************/
     138             : 
     139             : GEN
     140     2981043 : setloop(GEN a)
     141             : {
     142     2981043 :   pari_sp av = avma - 2 * sizeof(long);
     143     2981043 :   (void)cgetg(lgefint(a) + 3, t_VECSMALL);
     144     2981045 :   return icopy_avma(a, av); /* two cells of extra space after a */
     145             : }
     146             : 
     147             : /* we had a = setloop(?), then some incloops. Reset a to b */
     148             : GEN
     149      170696 : resetloop(GEN a, GEN b) {
     150      170696 :   a[0] = evaltyp(t_INT) | evallg(lgefint(b));
     151      170696 :   affii(b, a); return a;
     152             : }
     153             : 
     154             : /* assume a > 0, initialized by setloop. Do a++ */
     155             : static GEN
     156    31752505 : incpos(GEN a)
     157             : {
     158    31752505 :   long i, l = lgefint(a);
     159    31752510 :   for (i=2; i<l; i++)
     160    31752506 :     if (++uel(a,i)) return a;
     161           4 :   a[l] = 1; l++;
     162           4 :   a[0]=evaltyp(t_INT) | _evallg(l);
     163           4 :   a[1]=evalsigne(1) | evallgefint(l);
     164           4 :   return a;
     165             : }
     166             : 
     167             : /* assume a < 0, initialized by setloop. Do a++ */
     168             : static GEN
     169        8072 : incneg(GEN a)
     170             : {
     171        8072 :   long i, l = lgefint(a)-1;
     172        8072 :   if (uel(a,2)--)
     173             :   {
     174        8068 :     if (!a[l]) /* implies l = 2 */
     175             :     {
     176         580 :       a[0] = evaltyp(t_INT) | _evallg(2);
     177         580 :       a[1] = evalsigne(0) | evallgefint(2);
     178             :     }
     179        8068 :     return a;
     180             :   }
     181           5 :   for (i=3; i<=l; i++)
     182           5 :     if (uel(a,i)--) break;
     183           4 :   if (!a[l])
     184             :   {
     185           4 :     a[0] = evaltyp(t_INT) | _evallg(l);
     186           4 :     a[1] = evalsigne(-1) | evallgefint(l);
     187             :   }
     188           4 :   return a;
     189             : }
     190             : 
     191             : /* assume a initialized by setloop. Do a++ */
     192             : GEN
     193    32073379 : incloop(GEN a)
     194             : {
     195    32073379 :   switch(signe(a))
     196             :   {
     197             :     case 0:
     198      331708 :       a[0]=evaltyp(t_INT) | _evallg(3);
     199      331708 :       a[1]=evalsigne(1) | evallgefint(3);
     200      331708 :       a[2]=1; return a;
     201        8072 :     case -1: return incneg(a);
     202    31733599 :     default: return incpos(a);
     203             :   }
     204             : }
     205             : 
     206             : INLINE GEN
     207  1299260578 : adduispec(ulong s, GEN x, long nx)
     208             : {
     209             :   GEN  zd;
     210             :   long lz;
     211             : 
     212  1299260578 :   if (nx == 1) return adduu(uel(x,0), s);
     213   239534278 :   lz = nx+3; zd = cgeti(lz);
     214   238108893 :   if (mpn_add_1(LIMBS(zd),(mp_limb_t *)x,nx,s))
     215     5016469 :     zd[lz-1]=1;
     216             :   else
     217   233126732 :     lz--;
     218   238143201 :   zd[1] = evalsigne(1) | evallgefint(lz);
     219   238143201 :   return zd;
     220             : }
     221             : 
     222             : GEN
     223   273047769 : adduispec_offset(ulong s, GEN x, long offset, long nx)
     224             : {
     225   273047769 :   GEN xd=x+2+offset;
     226   273047769 :   while (nx && *(xd+nx-1)==0) nx--;
     227   273047769 :   if (!nx) return utoi(s);
     228   249549556 :   return adduispec(s,xd,nx);
     229             : }
     230             : 
     231             : INLINE GEN
     232  1775216668 : addiispec(GEN x, GEN y, long nx, long ny)
     233             : {
     234             :   GEN zd;
     235             :   long lz;
     236             : 
     237  1775216668 :   if (nx < ny) swapspec(x,y, nx,ny);
     238  1775216668 :   if (ny == 1) return adduispec(*y,x,nx);
     239   793268537 :   lz = nx+3; zd = cgeti(lz);
     240             : 
     241   789282878 :   if (mpn_add(LIMBS(zd),(mp_limb_t *)x,nx,(mp_limb_t *)y,ny))
     242    11466494 :     zd[lz-1]=1;
     243             :   else
     244   778627576 :     lz--;
     245             : 
     246   790094070 :   zd[1] = evalsigne(1) | evallgefint(lz);
     247   790094070 :   return zd;
     248             : }
     249             : 
     250             : /* assume x >= y */
     251             : INLINE GEN
     252   894011188 : subiuspec(GEN x, ulong s, long nx)
     253             : {
     254             :   GEN zd;
     255             :   long lz;
     256             : 
     257   894011188 :   if (nx == 1) return utoi(x[0] - s);
     258             : 
     259    85497633 :   lz = nx + 2; zd = cgeti(lz);
     260    85444096 :   mpn_sub_1 (LIMBS(zd), (mp_limb_t *)x, nx, s);
     261    85444860 :   if (! zd[lz - 1]) { --lz; }
     262             : 
     263    85444860 :   zd[1] = evalsigne(1) | evallgefint(lz);
     264    85444860 :   return zd;
     265             : }
     266             : 
     267             : /* assume x > y */
     268             : INLINE GEN
     269  1669617122 : subiispec(GEN x, GEN y, long nx, long ny)
     270             : {
     271             :   GEN zd;
     272             :   long lz;
     273  1669617122 :   if (ny==1) return subiuspec(x,*y,nx);
     274   785966631 :   lz = nx+2; zd = cgeti(lz);
     275             : 
     276   776511558 :   mpn_sub (LIMBS(zd), (mp_limb_t *)x, nx, (mp_limb_t *)y, ny);
     277   777428524 :   while (lz >= 3 && zd[lz - 1] == 0) { lz--; }
     278             : 
     279   777428524 :   zd[1] = evalsigne(1) | evallgefint(lz);
     280   777428524 :   return zd;
     281             : }
     282             : 
     283             : static void
     284   341315606 : roundr_up_ip(GEN x, long l)
     285             : {
     286   341315606 :   long i = l;
     287             :   for(;;)
     288             :   {
     289   345541226 :     if (++((ulong*)x)[--i]) break;
     290     2274274 :     if (i == 2) { x[2] = HIGHBIT; shiftr_inplace(x, 1); break; }
     291             :   }
     292   341315606 : }
     293             : 
     294             : void
     295   150372338 : affir(GEN x, GEN y)
     296             : {
     297   150372338 :   const long s = signe(x), ly = lg(y);
     298             :   long lx, sh, i;
     299             : 
     300   150372338 :   if (!s)
     301             :   {
     302     5127157 :     y[1] = evalexpo(-bit_accuracy(ly));
     303     5127157 :     return;
     304             :   }
     305   145245181 :   lx = lgefint(x); sh = bfffo(*int_MSW(x));
     306   145245181 :   y[1] = evalsigne(s) | evalexpo(bit_accuracy(lx)-sh-1);
     307   145245181 :   if (sh) {
     308   143486409 :     if (lx <= ly)
     309             :     {
     310    88944896 :       for (i=lx; i<ly; i++) y[i]=0;
     311    88944896 :       mpn_lshift(LIMBS(y),LIMBS(x),lx-2,sh);
     312    88944896 :       xmpn_mirror(LIMBS(y),lx-2);
     313    88944896 :       return;
     314             :     }
     315    54541513 :     mpn_lshift(LIMBS(y),LIMBS(x)+lx-ly,ly-2,sh);
     316    54541513 :     uel(y,2) |= uel(x,lx-ly+1) >> (BITS_IN_LONG-sh);
     317    54541513 :     xmpn_mirror(LIMBS(y),ly-2);
     318             :     /* lx > ly: round properly */
     319    54541513 :     if ((uel(x,lx-ly+1)<<sh) & HIGHBIT) roundr_up_ip(y, ly);
     320             :   }
     321             :   else {
     322     1758772 :     GEN xd=int_MSW(x);
     323     1758772 :     if (lx <= ly)
     324             :     {
     325      889975 :       for (i=2; i<lx; i++,xd=int_precW(xd)) y[i]=*xd;
     326      889975 :       for (   ; i<ly; i++) y[i]=0;
     327      889975 :       return;
     328             :     }
     329      868797 :     for (i=2; i<ly; i++,xd=int_precW(xd)) y[i]=*xd;
     330             :     /* lx > ly: round properly */
     331      868797 :     if (uel(x,lx-ly+1) & HIGHBIT) roundr_up_ip(y, ly);
     332             :   }
     333             : }
     334             : 
     335             : INLINE GEN
     336   427091274 : shiftispec(GEN x, long nx, long n)
     337             : {
     338             :   long ny,m;
     339             :   GEN yd, y;
     340             : 
     341   427091274 :   if (!n) return icopyspec(x, nx);
     342             : 
     343   415572926 :   if (n > 0)
     344             :   {
     345   263314219 :     long d = dvmdsBIL(n, &m);
     346             :     long i;
     347             : 
     348   263255705 :     ny = nx + d + (m!=0);
     349   263255705 :     y = cgeti(ny + 2); yd = y + 2;
     350   262665809 :     for (i=0; i<d; i++) yd[i] = 0;
     351             : 
     352   262665809 :     if (!m) xmpn_copy((mp_limb_t *) (yd + d), (mp_limb_t *) x, nx);
     353             :     else
     354             :     {
     355   262220527 :       ulong carryd = mpn_lshift((mp_limb_t *) (yd + d), (mp_limb_t *) x, nx, m);
     356   262576708 :       if (carryd) yd[ny - 1] = carryd;
     357   246764458 :       else ny--;
     358             :     }
     359             :   }
     360             :   else
     361             :   {
     362   152258707 :     long d = dvmdsBIL(-n, &m);
     363             : 
     364   152258625 :     ny = nx - d;
     365   152258625 :     if (ny < 1) return gen_0;
     366   151547493 :     y = cgeti(ny + 2); yd = y + 2;
     367             : 
     368   151544096 :     if (!m) xmpn_copy((mp_limb_t *) yd, (mp_limb_t *) (x + d), nx - d);
     369             :     else
     370             :     {
     371   150510260 :       mpn_rshift((mp_limb_t *) yd, (mp_limb_t *) (x + d), nx - d, m);
     372   150510217 :       if (yd[ny - 1] == 0)
     373             :       {
     374    16300391 :         if (ny == 1) { set_avma((pari_sp)(yd + 1)); return gen_0; }
     375    14835687 :         ny--;
     376             :       }
     377             :     }
     378             :   }
     379   412998263 :   y[1] = evalsigne(1)|evallgefint(ny + 2);
     380   412998263 :   return y;
     381             : }
     382             : 
     383             : GEN
     384    31881600 : mantissa2nr(GEN x, long n)
     385             : {
     386    31881600 :   long ly, i, m, s = signe(x), lx = lg(x);
     387             :   GEN y;
     388    31881600 :   if (!s) return gen_0;
     389    31880300 :   if (!n)
     390             :   {
     391    12033275 :     y = cgeti(lx);
     392    12033275 :     y[1] = evalsigne(s) | evallgefint(lx);
     393    12033275 :     xmpn_mirrorcopy(LIMBS(y),RLIMBS(x),lx-2);
     394    12033275 :     return y;
     395             :   }
     396    19847025 :   if (n > 0)
     397             :   {
     398      326286 :     GEN z = (GEN)avma;
     399      326286 :     long d = dvmdsBIL(n, &m);
     400             : 
     401      326286 :     ly = lx+d; y = new_chunk(ly);
     402      326286 :     for ( ; d; d--) *--z = 0;
     403      326286 :     if (!m) for (i=2; i<lx; i++) y[i]=x[i];
     404             :     else
     405             :     {
     406      325854 :       register const ulong sh = BITS_IN_LONG - m;
     407      325854 :       shift_left(y,x, 2,lx-1, 0,m);
     408      325854 :       i = uel(x,2) >> sh;
     409             :       /* Extend y on the left? */
     410      325854 :       if (i) { ly++; y = new_chunk(1); y[2] = i; }
     411             :     }
     412             :   }
     413             :   else
     414             :   {
     415    19520739 :     ly = lx - dvmdsBIL(-n, &m);
     416    19520739 :     if (ly<3) return gen_0;
     417    19520739 :     y = new_chunk(ly);
     418    19520739 :     if (m) {
     419    19471893 :       shift_right(y,x, 2,ly, 0,m);
     420    19471893 :       if (y[2] == 0)
     421             :       {
     422           0 :         if (ly==3) { set_avma((pari_sp)(y+3)); return gen_0; }
     423           0 :         ly--; set_avma((pari_sp)(++y));
     424             :       }
     425             :     } else {
     426       48846 :       for (i=2; i<ly; i++) y[i]=x[i];
     427             :     }
     428             :   }
     429    19847025 :   xmpn_mirror(LIMBS(y),ly-2);
     430    19847025 :   y[1] = evalsigne(s)|evallgefint(ly);
     431    19847025 :   y[0] = evaltyp(t_INT)|evallg(ly); return y;
     432             : }
     433             : 
     434             : GEN
     435     1982809 : truncr(GEN x)
     436             : {
     437             :   long s, e, d, m, i;
     438             :   GEN y;
     439     1982809 :   if ((s=signe(x)) == 0 || (e=expo(x)) < 0) return gen_0;
     440     1584465 :   d = nbits2lg(e+1); m = remsBIL(e);
     441     1584465 :   if (d > lg(x)) pari_err_PREC( "truncr (precision loss in truncation)");
     442             : 
     443     1584465 :   y=cgeti(d); y[1] = evalsigne(s) | evallgefint(d);
     444     1584465 :   if (++m == BITS_IN_LONG)
     445         244 :     for (i=2; i<d; i++) y[d-i+1]=x[i];
     446             :   else
     447             :   {
     448     1584221 :     GEN z=cgeti(d);
     449     1584221 :     for (i=2; i<d; i++) z[d-i+1]=x[i];
     450     1584221 :     mpn_rshift(LIMBS(y),LIMBS(z),d-2,BITS_IN_LONG-m);
     451     1584221 :     set_avma((pari_sp)y);
     452             :   }
     453     1584465 :   return y;
     454             : }
     455             : 
     456             : /* integral part */
     457             : GEN
     458      803406 : floorr(GEN x)
     459             : {
     460             :   long e, d, m, i, lx;
     461             :   GEN y;
     462      803406 :   if (signe(x) >= 0) return truncr(x);
     463      271627 :   if ((e=expo(x)) < 0) return gen_m1;
     464       97555 :   d = nbits2lg(e+1); m = remsBIL(e);
     465       97555 :   lx=lg(x); if (d>lx) pari_err_PREC( "floorr (precision loss in truncation)");
     466       97555 :   y = cgeti(d+1);
     467       97555 :   if (++m == BITS_IN_LONG)
     468             :   {
     469         126 :     for (i=2; i<d; i++) y[d-i+1]=x[i];
     470         126 :     i=d; while (i<lx && !x[i]) i++;
     471         126 :     if (i==lx) goto END;
     472             :   }
     473             :   else
     474             :   {
     475       97429 :     GEN z=cgeti(d);
     476       97429 :     for (i=2; i<d; i++) z[d-i+1]=x[i];
     477       97429 :     mpn_rshift(LIMBS(y),LIMBS(z),d-2,BITS_IN_LONG-m);
     478       97429 :     if (uel(x,d-1)<<m == 0)
     479             :     {
     480       25923 :       i=d; while (i<lx && !x[i]) i++;
     481       25923 :       if (i==lx) goto END;
     482             :     }
     483             :   }
     484       85181 :   if (mpn_add_1(LIMBS(y),LIMBS(y),d-2,1))
     485           0 :     y[d++]=1;
     486             : END:
     487       97555 :   y[1] = evalsigne(-1) | evallgefint(d);
     488       97555 :   return y;
     489             : }
     490             : 
     491             : INLINE int
     492  2074229662 : cmpiispec(GEN x, GEN y, long lx, long ly)
     493             : {
     494  2074229662 :   if (lx < ly) return -1;
     495  1979240681 :   if (lx > ly) return  1;
     496  1885581489 :   return mpn_cmp((mp_limb_t*)x,(mp_limb_t*)y, lx);
     497             : }
     498             : 
     499             : INLINE int
     500   169329402 : equaliispec(GEN x, GEN y, long lx, long ly)
     501             : {
     502   169329402 :   if (lx != ly) return 0;
     503   169312686 :   return !mpn_cmp((mp_limb_t*)x,(mp_limb_t*)y, lx);
     504             : }
     505             : 
     506             : /***********************************************************************/
     507             : /**                                                                   **/
     508             : /**                          MULTIPLICATION                           **/
     509             : /**                                                                   **/
     510             : /***********************************************************************/
     511             : /* assume ny > 0 */
     512             : INLINE GEN
     513  2591765556 : muluispec(ulong x, GEN y, long ny)
     514             : {
     515  2591765556 :   if (ny == 1)
     516  2008020564 :     return muluu(x, *y);
     517             :   else
     518             :   {
     519   583744992 :     long lz = ny+3;
     520   583744992 :     GEN z = cgeti(lz);
     521   580932903 :     ulong hi = mpn_mul_1 (LIMBS(z), (mp_limb_t *)y, ny, x);
     522   581106738 :     if (hi) { z[lz - 1] = hi; } else lz--;
     523   581106738 :     z[1] = evalsigne(1) | evallgefint(lz);
     524   581106738 :     return z;
     525             :   }
     526             : }
     527             : 
     528             : /* a + b*|y| */
     529             : GEN
     530      868513 : addumului(ulong a, ulong b, GEN y)
     531             : {
     532             :   GEN z;
     533             :   long i, lz;
     534             :   ulong hi;
     535      868513 :   if (!b || !signe(y)) return utoi(a);
     536      868430 :   lz = lgefint(y)+1;
     537      868430 :   z = cgeti(lz);
     538      868430 :   z[2]=a;
     539      868430 :   for(i=3;i<lz;i++) z[i]=0;
     540      868430 :   hi=mpn_addmul_1(LIMBS(z), LIMBS(y), NLIMBS(y), b);
     541      868430 :   if (hi) z[lz-1]=hi; else lz--;
     542      868430 :   z[1] = evalsigne(1) | evallgefint(lz);
     543      868430 :   set_avma((pari_sp)z); return z;
     544             : }
     545             : 
     546             : /***********************************************************************/
     547             : /**                                                                   **/
     548             : /**                          DIVISION                                 **/
     549             : /**                                                                   **/
     550             : /***********************************************************************/
     551             : 
     552             : ulong
     553  1088415095 : umodiu(GEN y, ulong x)
     554             : {
     555  1088415095 :   long sy=signe(y);
     556             :   ulong hi;
     557  1088415095 :   if (!x) pari_err_INV("umodiu",gen_0);
     558  1088626346 :   if (!sy) return 0;
     559   881947652 :   hi = mpn_mod_1(LIMBS(y),NLIMBS(y),x);
     560   882087379 :   if (!hi) return 0;
     561   871658110 :   return (sy > 0)? hi: x - hi;
     562             : }
     563             : 
     564             : /* return |y| \/ x */
     565             : GEN
     566   150936017 : absdiviu_rem(GEN y, ulong x, ulong *rem)
     567             : {
     568             :   long ly;
     569             :   GEN z;
     570             : 
     571   150936017 :   if (!x) pari_err_INV("absdiviu_rem",gen_0);
     572   150936512 :   if (!signe(y)) { *rem = 0; return gen_0; }
     573             : 
     574   150378195 :   ly = lgefint(y);
     575   150378195 :   if (ly == 3 && (ulong)x > uel(y,2)) { *rem = uel(y,2); return gen_0; }
     576             : 
     577   150301452 :   z = cgeti(ly);
     578   150239292 :   *rem = mpn_divrem_1(LIMBS(z), 0, LIMBS(y), NLIMBS(y), x);
     579   150274398 :   if (z [ly - 1] == 0) ly--;
     580   150274398 :   z[1] = evallgefint(ly) | evalsigne(1);
     581   150274398 :   return z;
     582             : }
     583             : 
     584             : GEN
     585    48302218 : divis_rem(GEN y, long x, long *rem)
     586             : {
     587    48302218 :   long sy=signe(y),ly,s;
     588             :   GEN z;
     589             : 
     590    48302218 :   if (!x) pari_err_INV("divis_rem",gen_0);
     591    48302251 :   if (!sy) { *rem = 0; return gen_0; }
     592    35557412 :   if (x<0) { s = -sy; x = -x; } else s = sy;
     593             : 
     594    35557412 :   ly = lgefint(y);
     595    35557412 :   if (ly == 3 && (ulong)x > uel(y,2)) { *rem = itos(y); return gen_0; }
     596             : 
     597    23022152 :   z = cgeti(ly);
     598    23022162 :   *rem = mpn_divrem_1(LIMBS(z), 0, LIMBS(y), NLIMBS(y), x);
     599    23022162 :   if (sy<0) *rem = -  *rem;
     600    23022162 :   if (z[ly - 1] == 0) ly--;
     601    23022162 :   z[1] = evallgefint(ly) | evalsigne(s);
     602    23022162 :   return z;
     603             : }
     604             : 
     605             : GEN
     606     1163046 : divis(GEN y, long x)
     607             : {
     608     1163046 :   long sy=signe(y),ly,s;
     609             :   GEN z;
     610             : 
     611     1163046 :   if (!x) pari_err_INV("divis",gen_0);
     612     1163046 :   if (!sy) return gen_0;
     613     1163046 :   if (x<0) { s = -sy; x = -x; } else s=sy;
     614             : 
     615     1163046 :   ly = lgefint(y);
     616     1163046 :   if (ly == 3 && (ulong)x > uel(y,2)) return gen_0;
     617             : 
     618     1163046 :   z = cgeti(ly);
     619     1163046 :   (void)mpn_divrem_1(LIMBS(z), 0, LIMBS(y), NLIMBS(y), x);
     620     1163046 :   if (z[ly - 1] == 0) ly--;
     621     1163046 :   z[1] = evallgefint(ly) | evalsigne(s);
     622     1163046 :   return z;
     623             : }
     624             : 
     625             : /* We keep llx bits of x and lly bits of y*/
     626             : static GEN
     627    48944924 : divrr_with_gmp(GEN x, GEN y)
     628             : {
     629    48944924 :   long lx=RNLIMBS(x),ly=RNLIMBS(y);
     630    48944924 :   long lw=minss(lx,ly);
     631    48944924 :   long lly=minss(lw+1,ly);
     632    48944924 :   GEN  w=cgetr(lw+2);
     633    48944924 :   long lu=lw+lly;
     634    48944924 :   long llx=minss(lu,lx);
     635    48944924 :   mp_limb_t *u=(mp_limb_t *)new_chunk(lu);
     636    48944924 :   mp_limb_t *z=(mp_limb_t *)new_chunk(lly);
     637             :   mp_limb_t *q,*r;
     638    48944924 :   long e=expo(x)-expo(y);
     639    48944924 :   long sx=signe(x),sy=signe(y);
     640    48944924 :   xmpn_mirrorcopy(z,RLIMBS(y),lly);
     641    48944924 :   xmpn_mirrorcopy(u+lu-llx,RLIMBS(x),llx);
     642    48944924 :   xmpn_zero(u,lu-llx);
     643    48944924 :   q = (mp_limb_t *)new_chunk(lw+1);
     644    48944924 :   r = (mp_limb_t *)new_chunk(lly);
     645             : 
     646    48944924 :   mpn_tdiv_qr(q,r,0,u,lu,z,lly);
     647             : 
     648             :   /*Round up: This is not exactly correct we should test 2*r>z*/
     649    48944924 :   if (uel(r,lly-1) > (uel(z,lly-1)>>1))
     650    24260078 :     mpn_add_1(q,q,lw+1,1);
     651             : 
     652    48944924 :   xmpn_mirrorcopy(RLIMBS(w),q,lw);
     653             : 
     654    48944924 :   if (q[lw] == 0) e--;
     655    19584928 :   else if (q[lw] == 1) { shift_right(w,w, 2,lw+2, 1,1); }
     656          18 :   else { w[2] = HIGHBIT; e++; }
     657    48944924 :   if (sy < 0) sx = -sx;
     658    48944924 :   w[1] = evalsigne(sx) | evalexpo(e);
     659    48944924 :   set_avma((pari_sp)w); return w;
     660             : }
     661             : 
     662             : /* We keep llx bits of x and lly bits of y*/
     663             : static GEN
     664     6312629 : divri_with_gmp(GEN x, GEN y)
     665             : {
     666     6312629 :   long llx=RNLIMBS(x),ly=NLIMBS(y);
     667     6312629 :   long lly=minss(llx+1,ly);
     668     6312629 :   GEN  w=cgetr(llx+2);
     669     6312629 :   long lu=llx+lly, ld=ly-lly;
     670     6312629 :   mp_limb_t *u=(mp_limb_t *)new_chunk(lu);
     671     6312629 :   mp_limb_t *z=(mp_limb_t *)new_chunk(lly);
     672             :   mp_limb_t *q,*r;
     673     6312629 :   long sh=bfffo(y[ly+1]);
     674     6312629 :   long e=expo(x)-expi(y);
     675     6312629 :   long sx=signe(x),sy=signe(y);
     676     6312629 :   if (sh) mpn_lshift(z,LIMBS(y)+ld,lly,sh);
     677      162922 :   else xmpn_copy(z,LIMBS(y)+ld,lly);
     678     6312629 :   xmpn_mirrorcopy(u+lu-llx,RLIMBS(x),llx);
     679     6312629 :   xmpn_zero(u,lu-llx);
     680     6312629 :   q = (mp_limb_t *)new_chunk(llx+1);
     681     6312629 :   r = (mp_limb_t *)new_chunk(lly);
     682             : 
     683     6312629 :   mpn_tdiv_qr(q,r,0,u,lu,z,lly);
     684             : 
     685             :   /*Round up: This is not exactly correct we should test 2*r>z*/
     686     6312629 :   if (uel(r,lly-1) > (uel(z,lly-1)>>1))
     687     2727353 :     mpn_add_1(q,q,llx+1,1);
     688             : 
     689     6312629 :   xmpn_mirrorcopy(RLIMBS(w),q,llx);
     690             : 
     691     6312629 :   if (q[llx] == 0) e--;
     692     3354833 :   else if (q[llx] == 1) { shift_right(w,w, 2,llx+2, 1,1); }
     693           0 :   else { w[2] = HIGHBIT; e++; }
     694     6312629 :   if (sy < 0) sx = -sx;
     695     6312629 :   w[1] = evalsigne(sx) | evalexpo(e);
     696     6312629 :   set_avma((pari_sp)w); return w;
     697             : }
     698             : 
     699             : GEN
     700    31591174 : divri(GEN x, GEN y)
     701             : {
     702    31591174 :   long  s = signe(y);
     703             : 
     704    31591174 :   if (!s) pari_err_INV("divri",gen_0);
     705    31591174 :   if (!signe(x)) return real_0_bit(expo(x) - expi(y));
     706    31448843 :   if (!is_bigint(y)) {
     707    25136214 :     GEN z = divru(x, y[2]);
     708    25136214 :     if (s < 0) togglesign(z);
     709    25136214 :     return z;
     710             :   }
     711     6312629 :   return divri_with_gmp(x,y);
     712             : }
     713             : 
     714             : GEN
     715   111086206 : divrr(GEN x, GEN y)
     716             : {
     717   111086206 :   long sx=signe(x), sy=signe(y), lx,ly,lr,e,i,j;
     718             :   ulong y0,y1;
     719             :   GEN r, r1;
     720             : 
     721   111086206 :   if (!sy) pari_err_INV("divrr",y);
     722   111086206 :   e = expo(x) - expo(y);
     723   111086206 :   if (!sx) return real_0_bit(e);
     724   109434835 :   if (sy<0) sx = -sx;
     725             : 
     726   109434835 :   lx=lg(x); ly=lg(y);
     727   109434835 :   if (ly==3)
     728             :   {
     729    60489911 :     ulong k = x[2], l = (lx>3)? x[3]: 0;
     730             :     LOCAL_HIREMAINDER;
     731    60489911 :     if (k < uel(y,2)) e--;
     732             :     else
     733             :     {
     734    27852871 :       l >>= 1; if (k&1) l |= HIGHBIT;
     735    27852871 :       k >>= 1;
     736             :     }
     737    60489911 :     hiremainder = k; k = divll(l,y[2]);
     738    60489911 :     if (hiremainder & HIGHBIT)
     739             :     {
     740    16794200 :       k++;
     741    16794200 :       if (!k) { k = HIGHBIT; e++; }
     742             :     }
     743    60489911 :     r = cgetr(3);
     744    60489911 :     r[1] = evalsigne(sx) | evalexpo(e);
     745    60489911 :     r[2] = k; return r;
     746             :   }
     747             : 
     748    48944924 :   if (ly>=DIVRR_GMP_LIMIT)
     749    48944924 :     return divrr_with_gmp(x,y);
     750             : 
     751           0 :   lr = minss(lx,ly); r = new_chunk(lr);
     752           0 :   r1 = r-1;
     753           0 :   r1[1] = 0; for (i=2; i<lr; i++) r1[i]=x[i];
     754           0 :   r1[lr] = (lx>ly)? x[lr]: 0;
     755           0 :   y0 = y[2]; y1 = y[3];
     756           0 :   for (i=0; i<lr-1; i++)
     757             :   { /* r1 = r + (i-1), OK up to r1[2] (accesses at most r[lr]) */
     758             :     ulong k, qp;
     759             :     LOCAL_HIREMAINDER;
     760             :     LOCAL_OVERFLOW;
     761             : 
     762           0 :     if (uel(r1,1) == y0)
     763             :     {
     764           0 :       qp = ULONG_MAX; k = addll(y0,r1[2]);
     765             :     }
     766             :     else
     767             :     {
     768           0 :       if (uel(r1,1) > y0) /* can't happen if i=0 */
     769             :       {
     770           0 :         GEN y1 = y+1;
     771           0 :         j = lr-i; r1[j] = subll(r1[j],y1[j]);
     772           0 :         for (j--; j>0; j--) r1[j] = subllx(r1[j],y1[j]);
     773           0 :         j=i; do uel(r,--j)++; while (j && !r[j]);
     774             :       }
     775           0 :       hiremainder = r1[1]; overflow = 0;
     776           0 :       qp = divll(r1[2],y0); k = hiremainder;
     777             :     }
     778           0 :     j = lr-i+1;
     779           0 :     if (!overflow)
     780             :     {
     781             :       long k3, k4;
     782           0 :       k3 = mulll(qp,y1);
     783           0 :       if (j == 3) /* i = lr - 2 maximal, r1[3] undefined -> 0 */
     784           0 :         k4 = subll(hiremainder,k);
     785             :       else
     786             :       {
     787           0 :         k3 = subll(k3, r1[3]);
     788           0 :         k4 = subllx(hiremainder,k);
     789             :       }
     790           0 :       while (!overflow && k4) { qp--; k3=subll(k3,y1); k4=subllx(k4,y0); }
     791             :     }
     792           0 :     if (j<ly) (void)mulll(qp,y[j]); else { hiremainder = 0 ; j = ly; }
     793           0 :     for (j--; j>1; j--)
     794             :     {
     795           0 :       r1[j] = subll(r1[j], addmul(qp,y[j]));
     796           0 :       hiremainder += overflow;
     797             :     }
     798           0 :     if (uel(r1,1) != hiremainder)
     799             :     {
     800           0 :       if (uel(r1,1) < hiremainder)
     801             :       {
     802           0 :         qp--;
     803           0 :         j = lr-i-(lr-i>=ly); r1[j] = addll(r1[j], y[j]);
     804           0 :         for (j--; j>1; j--) r1[j] = addllx(r1[j], y[j]);
     805             :       }
     806             :       else
     807             :       {
     808           0 :         uel(r1,1) -= hiremainder;
     809           0 :         while (r1[1])
     810             :         {
     811           0 :           qp++; if (!qp) { j=i; do uel(r,--j)++; while (j && !r[j]); }
     812           0 :           j = lr-i-(lr-i>=ly); r1[j] = subll(r1[j],y[j]);
     813           0 :           for (j--; j>1; j--) r1[j] = subllx(r1[j],y[j]);
     814           0 :           uel(r1,1) -= overflow;
     815             :         }
     816             :       }
     817             :     }
     818           0 :     *++r1 = qp;
     819             :   }
     820             :   /* i = lr-1 */
     821             :   /* round correctly */
     822           0 :   if (uel(r1,1) > (y0>>1))
     823             :   {
     824           0 :     j=i; do uel(r,--j)++; while (j && !r[j]);
     825             :   }
     826           0 :   r1 = r-1; for (j=i; j>=2; j--) r[j]=r1[j];
     827           0 :   if (r[0] == 0) e--;
     828           0 :   else if (r[0] == 1) { shift_right(r,r, 2,lr, 1,1); }
     829             :   else { /* possible only when rounding up to 0x2 0x0 ... */
     830           0 :     r[2] = (long)HIGHBIT; e++;
     831             :   }
     832           0 :   r[0] = evaltyp(t_REAL)|evallg(lr);
     833           0 :   r[1] = evalsigne(sx) | evalexpo(e);
     834           0 :   return r;
     835             : }
     836             : 
     837             : /* Integer division x / y: such that sign(r) = sign(x)
     838             :  *   if z = ONLY_REM return remainder, otherwise return quotient
     839             :  *   if z != NULL set *z to remainder
     840             :  *   *z is the last object on stack (and thus can be disposed of with cgiv
     841             :  *   instead of gerepile)
     842             :  * If *z is zero, we put gen_0 here and no copy.
     843             :  * space needed: lx + ly
     844             :  */
     845             : GEN
     846   914925006 : dvmdii(GEN x, GEN y, GEN *z)
     847             : {
     848   914925006 :   long sx=signe(x),sy=signe(y);
     849             :   long lx, ly, lq;
     850             :   pari_sp av;
     851             :   GEN r,q;
     852             : 
     853   914925006 :   if (!sy) pari_err_INV("dvmdii",y);
     854   919018762 :   if (!sx)
     855             :   {
     856    24036231 :     if (!z || z == ONLY_REM) return gen_0;
     857     9407074 :     *z=gen_0; return gen_0;
     858             :   }
     859   894982531 :   lx=lgefint(x);
     860   894982531 :   ly=lgefint(y); lq=lx-ly;
     861   894982531 :   if (lq <= 0)
     862             :   {
     863   554616914 :     if (lq == 0)
     864             :     {
     865   531198710 :       long s=mpn_cmp(LIMBS(x),LIMBS(y),NLIMBS(x));
     866   533443722 :       if (s>0) goto DIVIDE;
     867   202791044 :       if (s==0)
     868             :       {
     869    10135286 :         if (z == ONLY_REM) return gen_0;
     870     4498574 :         if (z) *z = gen_0;
     871     4498574 :         if (sx < 0) sy = -sy;
     872     4498574 :         return stoi(sy);
     873             :       }
     874             :     }
     875   216073962 :     if (z == ONLY_REM) return icopy(x);
     876     8653022 :     if (z) *z = icopy(x);
     877     8653019 :     return gen_0;
     878             :   }
     879             : DIVIDE: /* quotient is non-zero */
     880   671018295 :   av=avma; if (sx<0) sy = -sy;
     881   671018295 :   if (ly==3)
     882             :   {
     883   286138948 :     ulong lq = lx;
     884             :     ulong si;
     885   286138948 :     q  = cgeti(lq);
     886   285881613 :     si = mpn_divrem_1(LIMBS(q), 0, LIMBS(x), NLIMBS(x), y[2]);
     887   285959592 :     if (q[lq - 1] == 0) lq--;
     888   285959592 :     if (z == ONLY_REM)
     889             :     {
     890   239712822 :       set_avma(av); if (!si) return gen_0;
     891   211379548 :       r=cgeti(3);
     892   211235372 :       r[1] = evalsigne(sx) | evallgefint(3);
     893   211235372 :       r[2] = si; return r;
     894             :     }
     895    46246770 :     q[1] = evalsigne(sy) | evallgefint(lq);
     896    46246770 :     if (!z) return q;
     897    44512934 :     if (!si) { *z=gen_0; return q; }
     898    19307076 :     r=cgeti(3);
     899    19307049 :     r[1] = evalsigne(sx) | evallgefint(3);
     900    19307049 :     r[2] = si; *z=r; return q;
     901             :   }
     902   384879347 :   if (z == ONLY_REM)
     903             :   {
     904   377335622 :     ulong lr = lgefint(y);
     905   377335622 :     ulong lq = lgefint(x)-lgefint(y)+3;
     906   377335622 :     GEN r = cgeti(lr);
     907   368155924 :     GEN q = cgeti(lq);
     908   364403278 :     mpn_tdiv_qr(LIMBS(q), LIMBS(r),0, LIMBS(x), NLIMBS(x), LIMBS(y), NLIMBS(y));
     909   385968297 :     if (!r[lr - 1])
     910             :     {
     911    89486752 :       while(lr>2 && !r[lr - 1]) lr--;
     912    89486752 :       if (lr == 2) {set_avma(av); return gen_0;} /* exact division */
     913             :     }
     914   380509824 :     r[1] = evalsigne(sx) | evallgefint(lr);
     915   380509824 :     set_avma((pari_sp)r); return r;
     916             :   }
     917             :   else
     918             :   {
     919     7543725 :     ulong lq = lgefint(x)-lgefint(y)+3;
     920     7543725 :     ulong lr = lgefint(y);
     921     7543725 :     GEN q = cgeti(lq);
     922     7543725 :     GEN r = cgeti(lr);
     923     7543725 :     mpn_tdiv_qr(LIMBS(q), LIMBS(r),0, LIMBS(x), NLIMBS(x), LIMBS(y), NLIMBS(y));
     924     7543725 :     if (q[lq - 1] == 0) lq--;
     925     7543725 :     q[1] = evalsigne(sy) | evallgefint(lq);
     926     7543725 :     if (!z) { set_avma((pari_sp)q); return q; }
     927     7396114 :     if (!r[lr - 1])
     928             :     {
     929     2683986 :       while(lr>2 && !r[lr - 1]) lr--;
     930     2683986 :       if (lr == 2) { set_avma((pari_sp)q); *z=gen_0; return q; } /* exact division */
     931             :     }
     932     5249327 :     r[1] = evalsigne(sx) | evallgefint(lr);
     933     5249327 :     set_avma((pari_sp)r); *z = r; return q;
     934             :   }
     935             : }
     936             : 
     937             : /* Montgomery reduction.
     938             :  * N has k words, assume T >= 0 has less than 2k.
     939             :  * Return res := T / B^k mod N, where B = 2^BIL
     940             :  * such that 0 <= res < T/B^k + N  and  res has less than k words */
     941             : GEN
     942    13240861 : red_montgomery(GEN T, GEN N, ulong inv)
     943             : {
     944             :   pari_sp av;
     945             :   GEN Te, Td, Ne, Nd, scratch;
     946    13240861 :   ulong i, j, m, t, d, k = NLIMBS(N);
     947             :   int carry;
     948             :   LOCAL_HIREMAINDER;
     949             :   LOCAL_OVERFLOW;
     950             : 
     951    13240861 :   if (k == 0) return gen_0;
     952    13240861 :   d = NLIMBS(T); /* <= 2*k */
     953    13240861 :   if (d == 0) return gen_0;
     954             : #ifdef DEBUG
     955             :   if (d > 2*k) pari_err_BUG("red_montgomery");
     956             : #endif
     957    13240852 :   if (k == 1)
     958             :   { /* as below, special cased for efficiency */
     959         492 :     ulong n = uel(N,2);
     960         492 :     if (d == 1) {
     961         492 :       hiremainder = uel(T,2);
     962         492 :       m = hiremainder * inv;
     963         492 :       (void)addmul(m, n); /* t + m*n = 0 */
     964         492 :       return utoi(hiremainder);
     965             :     } else { /* d = 2 */
     966           0 :       hiremainder = uel(T,2);
     967           0 :       m = hiremainder * inv;
     968           0 :       (void)addmul(m, n); /* t + m*n = 0 */
     969           0 :       t = addll(hiremainder, uel(T,3));
     970           0 :       if (overflow) t -= n; /* t > n doesn't fit in 1 word */
     971           0 :       return utoi(t);
     972             :     }
     973             :   }
     974             :   /* assume k >= 2 */
     975    13240360 :   av = avma; scratch = new_chunk(k<<1); /* >= k + 2: result fits */
     976             : 
     977             :   /* copy T to scratch space (pad with zeroes to 2k words) */
     978    12962086 :   Td = scratch;
     979    12962086 :   Te = T + 2;
     980    12962086 :   for (i=0; i < d     ; i++) *Td++ = *Te++;
     981    12962086 :   for (   ; i < (k<<1); i++) *Td++ = 0;
     982             : 
     983    12962086 :   Te = scratch - 1; /* 1 beyond end of current T mantissa (in scratch) */
     984    12962086 :   Ne = N + 1;       /* 1 beyond end of N mantissa */
     985             : 
     986    12962086 :   carry = 0;
     987    95189404 :   for (i=0; i<k; i++) /* set T := T/B nod N, k times */
     988             :   {
     989    82227318 :     Td = Te; /* one beyond end of (new) T mantissa */
     990    82227318 :     Nd = Ne;
     991    82227318 :     hiremainder = *++Td;
     992    82227318 :     m = hiremainder * inv; /* solve T + m N = O(B) */
     993             : 
     994             :     /* set T := (T + mN) / B */
     995    82227318 :     Te = Td;
     996    82227318 :     (void)addmul(m, *++Nd); /* = 0 */
     997   760209202 :     for (j=1; j<k; j++)
     998             :     {
     999   677981884 :       t = addll(addmul(m, *++Nd), *++Td);
    1000   677981884 :       *Td = t;
    1001   677981884 :       hiremainder += overflow;
    1002             :     }
    1003    82227318 :     t = addll(hiremainder, *++Td); *Td = t + carry;
    1004    82227318 :     carry = (overflow || (carry && *Td == 0));
    1005             :   }
    1006    12962086 :   if (carry)
    1007             :   { /* Td > N overflows (k+1 words), set Td := Td - N */
    1008       35904 :     GEN NE = N + k+1;
    1009       35904 :     Td = Te;
    1010       35904 :     Nd = Ne;
    1011       35904 :     t = subll(*++Td, *++Nd); *Td = t;
    1012       35904 :     while (Nd < NE) { t = subllx(*++Td, *++Nd); *Td = t; }
    1013             :   }
    1014             : 
    1015             :   /* copy result */
    1016    12962086 :   Td = (GEN)av - 1; /* *Td = high word of final result */
    1017    12962086 :   while (*Td == 0 && Te < Td) Td--; /* strip leading 0s */
    1018    12962086 :   k = Td - Te; if (!k) { set_avma(av); return gen_0; }
    1019    12962086 :   Td = (GEN)av - k; /* will write mantissa there */
    1020    12962086 :   (void)memmove(Td, Te+1, k*sizeof(long));
    1021    12962086 :   Td -= 2;
    1022    12962086 :   Td[0] = evaltyp(t_INT) | evallg(k+2);
    1023    12948349 :   Td[1] = evalsigne(1) | evallgefint(k+2);
    1024             : #ifdef DEBUG
    1025             : {
    1026             :   long l = NLIMBS(N), s = BITS_IN_LONG*l;
    1027             :   GEN R = int2n(s);
    1028             :   GEN res = remii(mulii(T, Fp_inv(R, N)), N);
    1029             :   if (k > lgefint(N)
    1030             :     || !equalii(remii(Td,N),res)
    1031             :     || cmpii(Td, addii(shifti(T, -s), N)) >= 0) pari_err_BUG("red_montgomery");
    1032             : }
    1033             : #endif
    1034    12948349 :   set_avma((pari_sp)Td); return Td;
    1035             : }
    1036             : 
    1037             : /* EXACT INTEGER DIVISION */
    1038             : 
    1039             : /* use undocumented GMP interface */
    1040             : static void
    1041    39580966 : GEN2mpz(mpz_t X, GEN x)
    1042             : {
    1043    39580966 :   long l = lgefint(x)-2;
    1044    39580966 :   X->_mp_alloc = l;
    1045    39580966 :   X->_mp_size = signe(x) > 0? l: -l;
    1046    39580966 :   X->_mp_d = LIMBS(x);
    1047    39580966 : }
    1048             : static void
    1049    19790491 : mpz2GEN(GEN z, mpz_t Z)
    1050             : {
    1051    19790491 :   long l = Z->_mp_size;
    1052    19790491 :   z[1] = evalsigne(l > 0? 1: -1) | evallgefint(labs(l)+2);
    1053    19790491 : }
    1054             : 
    1055             : #ifdef mpn_divexact_1
    1056             : static GEN
    1057   278699973 : diviuexact_i(GEN x, ulong y)
    1058             : {
    1059   278699973 :   long l = lgefint(x);
    1060   278699973 :   GEN z = cgeti(l);
    1061   278679272 :   mpn_divexact_1(LIMBS(z), LIMBS(x), NLIMBS(x), y);
    1062   278680075 :   if (z[l-1] == 0) l--;
    1063   278680075 :   z[1] = evallgefint(l) | evalsigne(signe(x));
    1064   278680075 :   return z;
    1065             : }
    1066             : #elif 1 && !defined(_WIN64) /* mpz_divexact_ui is not LLP64 friendly   */
    1067             :                             /* assume y != 0 and the division is exact */
    1068             : static GEN
    1069             : diviuexact_i(GEN x, ulong y)
    1070             : {
    1071             :   long l = lgefint(x);
    1072             :   GEN z = cgeti(l);
    1073             :   mpz_t X, Z;
    1074             :   GEN2mpz(X, x);
    1075             :   Z->_mp_alloc = l-2;
    1076             :   Z->_mp_size  = l-2;
    1077             :   Z->_mp_d = LIMBS(z);
    1078             :   mpz_divexact_ui(Z, X, y);
    1079             :   mpz2GEN(z, Z); return z;
    1080             : }
    1081             : #else
    1082             : /* assume y != 0 and the division is exact */
    1083             : static GEN
    1084             : diviuexact_i(GEN x, ulong y)
    1085             : {
    1086             :   /*TODO: implement true exact division.*/
    1087             :   return divii(x,utoi(y));
    1088             : }
    1089             : #endif
    1090             : 
    1091             : GEN
    1092    37514552 : diviuexact(GEN x, ulong y)
    1093             : {
    1094             :   GEN z;
    1095    37514552 :   if (!signe(x)) return gen_0;
    1096    37242548 :   z = diviuexact_i(x, y);
    1097    37242737 :   if (lgefint(z) == 2) pari_err_OP("exact division", x, utoi(y));
    1098    37242694 :   return z;
    1099             : }
    1100             : 
    1101             : /* Find z such that x=y*z, knowing that y | x (unchecked) */
    1102             : GEN
    1103   317809209 : diviiexact(GEN x, GEN y)
    1104             : {
    1105             :   GEN z;
    1106   317809209 :   if (!signe(y)) pari_err_INV("diviiexact",y);
    1107   317825586 :   if (!signe(x)) return gen_0;
    1108   261239212 :   if (lgefint(y) == 3)
    1109             :   {
    1110   241448724 :     z = diviuexact_i(x, y[2]);
    1111   241438208 :     if (signe(y) < 0) togglesign(z);
    1112             :   }
    1113             :   else
    1114             :   {
    1115    19790488 :     long l = lgefint(x);
    1116             :     mpz_t X, Y, Z;
    1117    19790488 :     z = cgeti(l);
    1118    19790486 :     GEN2mpz(X, x);
    1119    19790485 :     GEN2mpz(Y, y);
    1120    19790482 :     Z->_mp_alloc = l-2;
    1121    19790482 :     Z->_mp_size  = l-2;
    1122    19790482 :     Z->_mp_d = LIMBS(z);
    1123    19790482 :     mpz_divexact(Z, X, Y);
    1124    19790491 :     mpz2GEN(z, Z);
    1125             :   }
    1126   261229051 :   if (lgefint(z) == 2) pari_err_OP("exact division", x, y);
    1127   261228836 :   return z;
    1128             : }
    1129             : 
    1130             : /* assume yz != and yz | x */
    1131             : GEN
    1132      133800 : diviuuexact(GEN x, ulong y, ulong z)
    1133             : {
    1134             :   long tmp[4];
    1135             :   ulong t;
    1136             :   LOCAL_HIREMAINDER;
    1137      133800 :   t = mulll(y, z);
    1138      133800 :   if (!hiremainder) return diviuexact(x, t);
    1139           0 :   tmp[0] = evaltyp(t_INT)|_evallg(4);
    1140           0 :   tmp[1] = evalsigne(1)|evallgefint(4);
    1141           0 :   tmp[2] = t;
    1142           0 :   tmp[3] = hiremainder;
    1143           0 :   return diviiexact(x, tmp);
    1144             : }
    1145             : 
    1146             : 
    1147             : /********************************************************************/
    1148             : /**                                                                **/
    1149             : /**               INTEGER MULTIPLICATION                           **/
    1150             : /**                                                                **/
    1151             : /********************************************************************/
    1152             : 
    1153             : /* nx >= ny = num. of digits of x, y (not GEN, see mulii) */
    1154             : GEN
    1155  2418622325 : muliispec(GEN x, GEN y, long nx, long ny)
    1156             : {
    1157             :   GEN zd;
    1158             :   long lz;
    1159             :   ulong hi;
    1160             : 
    1161  2418622325 :   if (nx < ny) swapspec(x,y, nx,ny);
    1162  2418622325 :   if (!ny) return gen_0;
    1163  2418622325 :   if (ny == 1) return muluispec((ulong)*y, x, nx);
    1164             : 
    1165   580262876 :   lz = nx+ny+2;
    1166   580262876 :   zd = cgeti(lz);
    1167   576186011 :   hi = mpn_mul(LIMBS(zd), (mp_limb_t *)x, nx, (mp_limb_t *)y, ny);
    1168   583772458 :   if (!hi) lz--;
    1169             :   /*else zd[lz-1]=hi; GH tell me it is not necessary.*/
    1170             : 
    1171   583772458 :   zd[1] = evalsigne(1) | evallgefint(lz);
    1172   583772458 :   return zd;
    1173             : }
    1174             : GEN
    1175      156024 : muluui(ulong x, ulong y, GEN z)
    1176             : {
    1177      156024 :   long t, s = signe(z);
    1178             :   GEN r;
    1179             :   LOCAL_HIREMAINDER;
    1180             : 
    1181      156024 :   if (!x || !y || !signe(z)) return gen_0;
    1182      155628 :   t = mulll(x,y);
    1183      155628 :   if (!hiremainder)
    1184      155628 :     r = muluispec(t, z+2, lgefint(z)-2);
    1185             :   else
    1186             :   {
    1187             :     long tmp[2];
    1188           0 :     tmp[1] = hiremainder;
    1189           0 :     tmp[0] = t;
    1190           0 :     r = muliispec(z+2,tmp, lgefint(z)-2, 2);
    1191             :   }
    1192      155628 :   setsigne(r,s); return r;
    1193             : }
    1194             : 
    1195             : GEN
    1196   769775175 : sqrispec(GEN x, long nx)
    1197             : {
    1198             :   GEN zd;
    1199             :   long lz;
    1200             : 
    1201   769775175 :   if (!nx) return gen_0;
    1202   229981356 :   if (nx==1) return sqru(*x);
    1203             : 
    1204   148362449 :   lz = (nx<<1)+2;
    1205   148362449 :   zd = cgeti(lz);
    1206             : #ifdef mpn_sqr
    1207   143992495 :   mpn_sqr(LIMBS(zd), (mp_limb_t *)x, nx);
    1208             : #else
    1209             :   mpn_mul_n(LIMBS(zd), (mp_limb_t *)x, (mp_limb_t *)x, nx);
    1210             : #endif
    1211   146311373 :   if (zd[lz-1]==0) lz--;
    1212             : 
    1213   146311373 :   zd[1] = evalsigne(1) | evallgefint(lz);
    1214   146311373 :   return zd;
    1215             : }
    1216             : 
    1217             : INLINE GEN
    1218     7032386 : sqrispec_mirror(GEN x, long nx)
    1219             : {
    1220     7032386 :   GEN cx=new_chunk(nx);
    1221             :   GEN z;
    1222     7032386 :   xmpn_mirrorcopy((mp_limb_t *)cx,(mp_limb_t *)x,nx);
    1223     7032386 :   z=sqrispec(cx, nx);
    1224     7032386 :   xmpn_mirror(LIMBS(z), NLIMBS(z));
    1225     7032386 :   return z;
    1226             : }
    1227             : 
    1228             : /* leaves garbage on the stack. */
    1229             : INLINE GEN
    1230    19462532 : muliispec_mirror(GEN x, GEN y, long nx, long ny)
    1231             : {
    1232             :   GEN cx, cy, z;
    1233    19462532 :   long s = 0;
    1234    19462532 :   while (nx && x[nx-1]==0) { nx--; s++; }
    1235    19462532 :   while (ny && y[ny-1]==0) { ny--; s++; }
    1236    19462532 :   cx=new_chunk(nx); cy=new_chunk(ny);
    1237    19462532 :   xmpn_mirrorcopy((mp_limb_t *)cx,(mp_limb_t *)x,nx);
    1238    19462532 :   xmpn_mirrorcopy((mp_limb_t *)cy,(mp_limb_t *)y,ny);
    1239    19462532 :   z =  nx>=ny ? muliispec(cx, cy, nx, ny): muliispec(cy, cx, ny, nx);
    1240    19462532 :   if (s)
    1241             :   {
    1242     2007314 :     long i, lz = lgefint(z) + s;
    1243     2007314 :     (void)new_chunk(s);
    1244     2007314 :     z -= s;
    1245     2007314 :     for (i=0; i<s; i++) z[2+i]=0;
    1246     2007314 :     z[1] = evalsigne(1) | evallgefint(lz);
    1247     2007314 :     z[0] = evaltyp(t_INT) | evallg(lz);
    1248             :   }
    1249    19462532 :   xmpn_mirror(LIMBS(z), NLIMBS(z));
    1250    19462532 :   return z;
    1251             : }
    1252             : 
    1253             : /* x % (2^n), assuming n >= 0 */
    1254             : GEN
    1255    17493746 : remi2n(GEN x, long n)
    1256             : {
    1257             :   ulong hi;
    1258    17493746 :   long l, k, lx, ly, sx = signe(x);
    1259             :   GEN z, xd, zd;
    1260             : 
    1261    17493746 :   if (!sx || !n) return gen_0;
    1262             : 
    1263    17339756 :   k = dvmdsBIL(n, &l);
    1264    17974844 :   lx = lgefint(x);
    1265    17974844 :   if (lx < k+3) return icopy(x);
    1266             : 
    1267    17662683 :   xd = x + (2 + k);
    1268             :   /* x = |k|...|1|#|... : copy the last l bits of # and the first k words
    1269             :    *              ^--- initial xd  */
    1270    17662683 :   hi = ((ulong)*xd) & ((1UL<<l)-1); /* last l bits of # = top bits of result */
    1271    17662683 :   if (!hi)
    1272             :   { /* strip leading zeroes from result */
    1273     1782062 :     xd--; while (k && !*xd) { k--; xd--; }
    1274     1782062 :     if (!k) return gen_0;
    1275     1022579 :     ly = k+2;
    1276             :   }
    1277             :   else
    1278    15880621 :     ly = k+3;
    1279             : 
    1280    16903200 :   zd = z = cgeti(ly);
    1281    16426069 :   *++zd = evalsigne(sx) | evallgefint(ly);
    1282    16426069 :   xd = x+1; for ( ;k; k--) *++zd = *++xd;
    1283    16426069 :   if (hi) *++zd = hi;
    1284    16426069 :   return z;
    1285             : }
    1286             : 
    1287             : /********************************************************************/
    1288             : /**                                                                **/
    1289             : /**                      INTEGER SQUARE ROOT                       **/
    1290             : /**                                                                **/
    1291             : /********************************************************************/
    1292             : 
    1293             : /* Return S (and set R) s.t S^2 + R = N, 0 <= R <= 2S.
    1294             :  * As for dvmdii, R is last on stack and guaranteed to be gen_0 in case the
    1295             :  * remainder is 0. R = NULL is allowed. */
    1296             : GEN
    1297     3697650 : sqrtremi(GEN a, GEN *r)
    1298             : {
    1299     3697650 :   long l, na = NLIMBS(a);
    1300             :   mp_size_t nr;
    1301             :   GEN S;
    1302     3697650 :   if (!na) {
    1303          24 :     if (r) *r = gen_0;
    1304          24 :     return gen_0;
    1305             :   }
    1306     3697626 :   l = (na + 5) >> 1; /* 2 + ceil(na/2) */
    1307     3697626 :   S = cgetipos(l);
    1308     3697586 :   if (r) {
    1309      349132 :     GEN R = cgeti(2 + na);
    1310      349132 :     nr = mpn_sqrtrem(LIMBS(S), LIMBS(R), LIMBS(a), na);
    1311      349132 :     if (nr) R[1] = evalsigne(1) | evallgefint(nr+2);
    1312       10927 :     else    { set_avma((pari_sp)S); R = gen_0; }
    1313      349132 :     *r = R;
    1314             :   }
    1315             :   else
    1316     3348454 :     (void)mpn_sqrtrem(LIMBS(S), NULL, LIMBS(a), na);
    1317     3697586 :   return S;
    1318             : }
    1319             : 
    1320             : /* compute sqrt(|a|), assuming a != 0 */
    1321             : GEN
    1322    28701165 : sqrtr_abs(GEN a)
    1323             : {
    1324             :   GEN res;
    1325             :   mp_limb_t *b, *c;
    1326    28701165 :   long l = RNLIMBS(a), e = expo(a), er = e>>1;
    1327             :   long n;
    1328    28701165 :   res = cgetr(2 + l);
    1329    28701165 :   res[1] = evalsigne(1) | evalexpo(er);
    1330    28701165 :   if (e&1)
    1331             :   {
    1332    14046342 :     b = (mp_limb_t *) new_chunk(l<<1);
    1333    14046342 :     xmpn_zero(b,l);
    1334    14046342 :     xmpn_mirrorcopy(b+l, RLIMBS(a), l);
    1335    14046342 :     c = (mp_limb_t *) new_chunk(l);
    1336    14046342 :     n = mpn_sqrtrem(c,b,b,l<<1); /* c <- sqrt; b <- rem */
    1337    14046342 :     if (n>l || (n==l && mpn_cmp(b,c,l) > 0)) mpn_add_1(c,c,l,1);
    1338             :   }
    1339             :   else
    1340             :   {
    1341             :     ulong u;
    1342    14654823 :     b = (mp_limb_t *) mantissa2nr(a,-1);
    1343    14654823 :     b[1] = uel(a,l+1)<<(BITS_IN_LONG-1);
    1344    14654823 :     b = (mp_limb_t *) new_chunk(l);
    1345    14654823 :     xmpn_zero(b,l+1); /* overwrites the former b[0] */
    1346    14654823 :     c = (mp_limb_t *) new_chunk(l + 1);
    1347    14654823 :     n = mpn_sqrtrem(c,b,b,(l<<1)+2); /* c <- sqrt; b <- rem */
    1348    14654823 :     u = (ulong)*c++;
    1349    14654823 :     if ( u&HIGHBIT || (u == ~HIGHBIT &&
    1350           0 :              (n>l || (n==l && mpn_cmp(b,c,l)>0))))
    1351     7186144 :       mpn_add_1(c,c,l,1);
    1352             :   }
    1353    28701165 :   xmpn_mirrorcopy(RLIMBS(res),c,l);
    1354    28701165 :   set_avma((pari_sp)res); return res;
    1355             : }
    1356             : 
    1357             : /* Normalize a non-negative integer */
    1358             : GEN
    1359    95627246 : int_normalize(GEN x, long known_zero_words)
    1360             : {
    1361    95627246 :   long i =  lgefint(x) - 1 - known_zero_words;
    1362  1229732980 :   for ( ; i > 1; i--)
    1363  1196698069 :     if (x[i]) { setlgefint(x, i+1); return x; }
    1364    33034911 :   x[1] = evalsigne(0) | evallgefint(2); return x;
    1365             : }
    1366             : 
    1367             : /********************************************************************
    1368             :  **                                                                **
    1369             :  **                           Base Conversion                      **
    1370             :  **                                                                **
    1371             :  ********************************************************************/
    1372             : 
    1373             : ulong *
    1374      328141 : convi(GEN x, long *l)
    1375             : {
    1376      328141 :   long n = nchar2nlong(2 + (long)(NLIMBS(x) * (BITS_IN_LONG * LOG10_2)));
    1377      328141 :   GEN str = cgetg(n+1, t_VECSMALL);
    1378      328141 :   unsigned char *res = (unsigned char*) GSTR(str);
    1379      328141 :   long llz = mpn_get_str(res, 10, LIMBS(icopy(x)), NLIMBS(x));
    1380             :   long lz;
    1381             :   ulong *z;
    1382             :   long i, j;
    1383             :   unsigned char *t;
    1384      328141 :   while (!*res) {res++; llz--;} /*Strip leading zeros*/
    1385      328141 :   lz  = (8+llz)/9;
    1386      328141 :   z = (ulong*)new_chunk(1+lz);
    1387      328141 :   t=res+llz+9;
    1388      706223 :   for(i=0;i<llz-8;i+=9)
    1389             :   {
    1390             :     ulong s;
    1391      378082 :     t-=18;
    1392      378082 :     s=*t++;
    1393     3402738 :     for (j=1; j<9;j++)
    1394     3024656 :       s=10*s+*t++;
    1395      378082 :     *z++=s;
    1396             :   }
    1397      328141 :   if (i<llz)
    1398             :   {
    1399      324650 :     unsigned char *t = res;
    1400      324650 :     ulong s=*t++;
    1401      845211 :     for (j=i+1; j<llz;j++)
    1402      520561 :       s=10*s+*t++;
    1403      324650 :     *z++=s;
    1404             :   }
    1405      328141 :   *l = lz;
    1406      328141 :   return z;
    1407             : }

Generated by: LCOV version 1.13