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.18.1 lcov report (development 30213-f9b05da6c0) Lines: 683 719 95.0 %
Date: 2025-04-26 09:18:30 Functions: 54 56 96.4 %
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; either version 2 of the License, or (at your option) any later
       9             : version. It is distributed in the hope that it will be useful, but WITHOUT
      10             : ANY WARRANTY WHATSOEVER.
      11             : 
      12             : Check the License for details. You should have received a copy of it, along
      13             : with the package; see the file 'COPYING'. If not, write to the Free Software
      14             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      15             : 
      16             : /***********************************************************************/
      17             : /**                                                                   **/
      18             : /**                               GMP KERNEL                          **/
      19             : /** BA2002Sep24                                                       **/
      20             : /***********************************************************************/
      21             : /* GMP t_INT as just like normal t_INT, just the mantissa is the other way
      22             :  * round
      23             :  *
      24             :  *   `How would you like to live in Looking-glass House, Kitty?  I
      25             :  *   wonder if they'd give you milk in there?  Perhaps Looking-glass
      26             :  *   milk isn't good to drink--But oh, Kitty! now we come to the
      27             :  *   passage.  You can just see a little PEEP of the passage in
      28             :  *   Looking-glass House, if you leave the door of our drawing-room
      29             :  *   wide open:  and it's very like our passage as far as you can see,
      30             :  *   only you know it may be quite different on beyond.  Oh, Kitty!
      31             :  *   how nice it would be if we could only get through into Looking-
      32             :  *   glass House!  I'm sure it's got, oh! such beautiful things in it!
      33             :  *
      34             :  *  Through the Looking Glass,  Lewis Carrol
      35             :  *
      36             :  *  (pityful attempt to beat GN code/comments rate)
      37             :  *  */
      38             : 
      39             : #include <gmp.h>
      40             : #include "pari.h"
      41             : #include "paripriv.h"
      42             : #include "../src/kernel/none/tune-gen.h"
      43             : 
      44             : /*We need PARI invmod renamed to invmod_pari*/
      45             : #define INVMOD_PARI
      46             : 
      47           0 : static void *pari_gmp_realloc(void *ptr, size_t old_size, size_t new_size) {
      48           0 :   (void)old_size; return (void *) pari_realloc(ptr,new_size);
      49             : }
      50             : 
      51     1759236 : static void pari_gmp_free(void *ptr, size_t old_size){
      52     1759236 :   (void)old_size; pari_free(ptr);
      53     1759236 : }
      54             : 
      55             : static void *(*old_gmp_malloc)(size_t new_size);
      56             : static void *(*old_gmp_realloc)(void *ptr, size_t old_size, size_t new_size);
      57             : static void (*old_gmp_free)(void *ptr, size_t old_size);
      58             : 
      59             : void
      60        1112 : pari_kernel_init(void)
      61             : {
      62        1112 :   mp_get_memory_functions (&old_gmp_malloc, &old_gmp_realloc, &old_gmp_free);
      63        1112 :   mp_set_memory_functions((void *(*)(size_t)) pari_malloc, pari_gmp_realloc, pari_gmp_free);
      64        1112 : }
      65             : 
      66             : const char *
      67           4 : pari_kernel_version(void)
      68             : {
      69             : #ifdef gmp_version
      70           4 :   return gmp_version;
      71             : #else
      72             :   return "";
      73             : #endif
      74             : }
      75             : 
      76             : void
      77        1104 : pari_kernel_close(void)
      78             : {
      79             :   void *(*new_gmp_malloc)(size_t new_size);
      80             :   void *(*new_gmp_realloc)(void *ptr, size_t old_size, size_t new_size);
      81             :   void (*new_gmp_free)(void *ptr, size_t old_size);
      82        1104 :   mp_get_memory_functions (&new_gmp_malloc, &new_gmp_realloc, &new_gmp_free);
      83        1104 :   if (new_gmp_malloc==pari_malloc) new_gmp_malloc = old_gmp_malloc;
      84        1104 :   if (new_gmp_realloc==pari_gmp_realloc) new_gmp_realloc = old_gmp_realloc;
      85        1104 :   if (new_gmp_free==pari_gmp_free) new_gmp_free = old_gmp_free;
      86        1104 :   mp_set_memory_functions(new_gmp_malloc, new_gmp_realloc, new_gmp_free);
      87        1104 : }
      88             : 
      89             : #define LIMBS(x)  ((mp_limb_t *)((x)+2))
      90             : #define NLIMBS(x) (lgefint(x)-2)
      91             : /*This one is for t_REALs to emphasize they are not t_INTs*/
      92             : #define RLIMBS(x)  ((mp_limb_t *)((x)+2))
      93             : #define RNLIMBS(x) (lg(x)-2)
      94             : 
      95             : INLINE void
      96     6705743 : xmpn_copy(mp_limb_t *x, mp_limb_t *y, long n)
      97             : {
      98    55248937 :   while (--n >= 0) x[n]=y[n];
      99     6705743 : }
     100             : 
     101             : INLINE void
     102   583052100 : xmpn_mirror(mp_limb_t *x, long n)
     103             : {
     104             :   long i;
     105  3875413176 :   for(i=0;i<(n>>1);i++)
     106             :   {
     107  3292361076 :     ulong m=x[i];
     108  3292361076 :     x[i]=x[n-1-i];
     109  3292361076 :     x[n-1-i]=m;
     110             :   }
     111   583052100 : }
     112             : 
     113             : INLINE void
     114   712714232 : xmpn_mirrorcopy(mp_limb_t *z, mp_limb_t *x, long n)
     115             : {
     116             :   long i;
     117  9818627191 :   for(i=0;i<n;i++)
     118  9105912959 :     z[i]=x[n-1-i];
     119   712714232 : }
     120             : 
     121             : INLINE void
     122   235821622 : xmpn_zero(mp_ptr x, mp_size_t n)
     123             : {
     124  2056367613 :   while (--n >= 0) x[n]=0;
     125   235821622 : }
     126             : 
     127             : INLINE GEN
     128    41190219 : icopy_ef(GEN x, long l)
     129             : {
     130    41190219 :   long lx = lgefint(x);
     131    41190219 :   const GEN y = cgeti(l);
     132             : 
     133   293251406 :   while (--lx > 0) y[lx]=x[lx];
     134    41188650 :   return y;
     135             : }
     136             : 
     137             : /* NOTE: arguments of "spec" routines (muliispec, addiispec, etc.) aren't
     138             :  * GENs but pairs (long *a, long na) representing a list of digits (in basis
     139             :  * BITS_IN_LONG) : a[0], ..., a[na-1]. In order to facilitate splitting: no
     140             :  * need to reintroduce codewords. */
     141             : 
     142             : /***********************************************************************/
     143             : /**                                                                   **/
     144             : /**                     ADDITION / SUBTRACTION                        **/
     145             : /**                                                                   **/
     146             : /***********************************************************************/
     147             : 
     148             : GEN
     149     2997600 : setloop(GEN a)
     150             : {
     151     2997600 :   pari_sp av = avma - 2 * sizeof(long);
     152     2997600 :   (void)cgetg(lgefint(a) + 3, t_VECSMALL);
     153     2997600 :   return icopy_avma(a, av); /* two cells of extra space after a */
     154             : }
     155             : 
     156             : /* we had a = setloop(?), then some incloops. Reset a to b */
     157             : GEN
     158      174328 : resetloop(GEN a, GEN b) {
     159      174328 :   a[0] = evaltyp(t_INT) | evallg(lgefint(b));
     160      174328 :   affii(b, a); return a;
     161             : }
     162             : 
     163             : /* assume a > 0, initialized by setloop. Do a++ */
     164             : static GEN
     165   102877718 : incpos(GEN a)
     166             : {
     167   102877718 :   long i, l = lgefint(a);
     168   102877723 :   for (i=2; i<l; i++)
     169   102883521 :     if (++uel(a,i)) return a;
     170           3 :   a[l] = 1; l++;
     171           3 :   a[0]=evaltyp(t_INT) | _evallg(l);
     172           3 :   a[1]=evalsigne(1) | evallgefint(l);
     173           3 :   return a;
     174             : }
     175             : 
     176             : /* assume a < 0, initialized by setloop. Do a++ */
     177             : static GEN
     178       66684 : incneg(GEN a)
     179             : {
     180       66684 :   long i, l = lgefint(a)-1;
     181       66684 :   if (uel(a,2)--)
     182             :   {
     183       66680 :     if (!a[l]) /* implies l = 2 */
     184             :     {
     185        1980 :       a[0] = evaltyp(t_INT) | _evallg(2);
     186        1980 :       a[1] = evalsigne(0) | evallgefint(2);
     187             :     }
     188       66680 :     return a;
     189             :   }
     190           5 :   for (i=3; i<=l; i++)
     191           5 :     if (uel(a,i)--) break;
     192           4 :   if (!a[l])
     193             :   {
     194           4 :     a[0] = evaltyp(t_INT) | _evallg(l);
     195           4 :     a[1] = evalsigne(-1) | evallgefint(l);
     196             :   }
     197           4 :   return a;
     198             : }
     199             : 
     200             : /* assume a initialized by setloop. Do a++ */
     201             : GEN
     202   103281905 : incloop(GEN a)
     203             : {
     204   103281905 :   switch(signe(a))
     205             :   {
     206      336632 :     case 0:
     207      336632 :       a[0]=evaltyp(t_INT) | _evallg(3);
     208      336632 :       a[1]=evalsigne(1) | evallgefint(3);
     209      336632 :       a[2]=1; return a;
     210       66684 :     case -1: return incneg(a);
     211   102878589 :     default: return incpos(a);
     212             :   }
     213             : }
     214             : 
     215             : INLINE GEN
     216  2774025495 : adduispec(ulong s, GEN x, long nx)
     217             : {
     218             :   GEN  zd;
     219             :   long lz;
     220             : 
     221  2774025495 :   if (nx == 1) return adduu(uel(x,0), s);
     222   900394539 :   lz = nx+3; zd = cgeti(lz);
     223   903331599 :   if (mpn_add_1(LIMBS(zd),(mp_limb_t *)x,nx,s))
     224     5005329 :     zd[lz-1]=1;
     225             :   else
     226   898334682 :     lz--;
     227   903340011 :   zd[1] = evalsigne(1) | evallgefint(lz);
     228   903340011 :   return zd;
     229             : }
     230             : 
     231             : GEN
     232   765559386 : adduispec_offset(ulong s, GEN x, long offset, long nx)
     233             : {
     234   765559386 :   GEN xd=x+2+offset;
     235   997186515 :   while (nx && *(xd+nx-1)==0) nx--;
     236   765559386 :   if (!nx) return utoi(s);
     237   720852800 :   return adduispec(s,xd,nx);
     238             : }
     239             : 
     240             : INLINE GEN
     241  3313320403 : addiispec(GEN x, GEN y, long nx, long ny)
     242             : {
     243             :   GEN zd;
     244             :   long lz;
     245             : 
     246  3313320403 :   if (nx < ny) swapspec(x,y, nx,ny);
     247  3313320403 :   if (ny == 1) return adduispec(*y,x,nx);
     248  1347994594 :   lz = nx+3; zd = cgeti(lz);
     249             : 
     250  1348286744 :   if (mpn_add(LIMBS(zd),(mp_limb_t *)x,nx,(mp_limb_t *)y,ny))
     251    30339077 :     zd[lz-1]=1;
     252             :   else
     253  1318256056 :     lz--;
     254             : 
     255  1348595133 :   zd[1] = evalsigne(1) | evallgefint(lz);
     256  1348595133 :   return zd;
     257             : }
     258             : 
     259             : /* assume x >= y */
     260             : INLINE GEN
     261  1767934759 : subiuspec(GEN x, ulong s, long nx)
     262             : {
     263             :   GEN zd;
     264             :   long lz;
     265             : 
     266  1767934759 :   if (nx == 1) return utoi(x[0] - s);
     267             : 
     268   251033788 :   lz = nx + 2; zd = cgeti(lz);
     269   251542984 :   mpn_sub_1 (LIMBS(zd), (mp_limb_t *)x, nx, s);
     270   251543677 :   if (! zd[lz - 1]) { --lz; }
     271             : 
     272   251543677 :   zd[1] = evalsigne(1) | evallgefint(lz);
     273   251543677 :   return zd;
     274             : }
     275             : 
     276             : /* assume x > y */
     277             : INLINE GEN
     278  3002703912 : subiispec(GEN x, GEN y, long nx, long ny)
     279             : {
     280             :   GEN zd;
     281             :   long lz;
     282  3002703912 :   if (ny==1) return subiuspec(x,*y,nx);
     283  1264971450 :   lz = nx+2; zd = cgeti(lz);
     284             : 
     285  1262524942 :   mpn_sub (LIMBS(zd), (mp_limb_t *)x, nx, (mp_limb_t *)y, ny);
     286  1730181043 :   while (lz >= 3 && zd[lz - 1] == 0) { lz--; }
     287             : 
     288  1262808761 :   zd[1] = evalsigne(1) | evallgefint(lz);
     289  1262808761 :   return zd;
     290             : }
     291             : 
     292             : static void
     293   521939860 : roundr_up_ip(GEN x, long l)
     294             : {
     295   521939860 :   long i = l;
     296             :   for(;;)
     297             :   {
     298   523733118 :     if (++((ulong*)x)[--i]) break;
     299     2179383 :     if (i == 2) { x[2] = HIGHBIT; shiftr_inplace(x, 1); break; }
     300             :   }
     301   521991487 : }
     302             : 
     303             : void
     304   400728415 : affir(GEN x, GEN y)
     305             : {
     306   400728415 :   const long s = signe(x), ly = lg(y);
     307             :   long lx, sh, i;
     308             : 
     309   400728415 :   if (!s)
     310             :   {
     311    41781811 :     y[1] = evalexpo(-bit_accuracy(ly));
     312    41778671 :     return;
     313             :   }
     314   358946604 :   lx = lgefint(x); sh = bfffo(*int_MSW(x));
     315   358946604 :   y[1] = evalsigne(s) | evalexpo(bit_accuracy(lx)-sh-1);
     316   359137848 :   if (sh) {
     317   351629455 :     if (lx <= ly)
     318             :     {
     319   953138198 :       for (i=lx; i<ly; i++) y[i]=0;
     320   247288781 :       mpn_lshift(LIMBS(y),LIMBS(x),lx-2,sh);
     321   247371084 :       xmpn_mirror(LIMBS(y),lx-2);
     322   247636027 :       return;
     323             :     }
     324   104340674 :     mpn_lshift(LIMBS(y),LIMBS(x)+lx-ly,ly-2,sh);
     325   104342004 :     uel(y,2) |= uel(x,lx-ly+1) >> (BITS_IN_LONG-sh);
     326   104342004 :     xmpn_mirror(LIMBS(y),ly-2);
     327             :     /* lx > ly: round properly */
     328   104339217 :     if ((uel(x,lx-ly+1)<<sh) & HIGHBIT) roundr_up_ip(y, ly);
     329             :   }
     330             :   else {
     331     7508393 :     GEN xd=int_MSW(x);
     332     7508393 :     if (lx <= ly)
     333             :     {
     334     9421082 :       for (i=2; i<lx; i++,xd=int_precW(xd)) y[i]=*xd;
     335     5135276 :       for (   ; i<ly; i++) y[i]=0;
     336     2353535 :       return;
     337             :     }
     338    14692168 :     for (i=2; i<ly; i++,xd=int_precW(xd)) y[i]=*xd;
     339             :     /* lx > ly: round properly */
     340     5154858 :     if (uel(x,lx-ly+1) & HIGHBIT) roundr_up_ip(y, ly);
     341             :   }
     342             : }
     343             : 
     344             : INLINE GEN
     345   703217686 : shiftispec(GEN x, long nx, long n)
     346             : {
     347             :   long ny,m;
     348             :   GEN yd, y;
     349             : 
     350   703217686 :   if (!n) return icopyspec(x, nx);
     351             : 
     352   674097647 :   if (n > 0)
     353             :   {
     354   395126891 :     long d = dvmdsBIL(n, &m);
     355             :     long i;
     356             : 
     357   395025834 :     ny = nx + d + (m!=0);
     358   395025834 :     y = cgeti(ny + 2); yd = y + 2;
     359   548443186 :     for (i=0; i<d; i++) yd[i] = 0;
     360             : 
     361   394097578 :     if (!m) xmpn_copy((mp_limb_t *) (yd + d), (mp_limb_t *) x, nx);
     362             :     else
     363             :     {
     364   392566661 :       ulong carryd = mpn_lshift((mp_limb_t *) (yd + d), (mp_limb_t *) x, nx, m);
     365   392624349 :       if (carryd) yd[ny - 1] = carryd;
     366   369525198 :       else ny--;
     367             :     }
     368             :   }
     369             :   else
     370             :   {
     371   278970756 :     long d = dvmdsBIL(-n, &m);
     372             : 
     373   280965724 :     ny = nx - d;
     374   280965724 :     if (ny < 1) return gen_0;
     375   278122206 :     y = cgeti(ny + 2); yd = y + 2;
     376             : 
     377   277826246 :     if (!m) xmpn_copy((mp_limb_t *) yd, (mp_limb_t *) (x + d), nx - d);
     378             :     else
     379             :     {
     380   273335151 :       mpn_rshift((mp_limb_t *) yd, (mp_limb_t *) (x + d), nx - d, m);
     381   273349025 :       if (yd[ny - 1] == 0)
     382             :       {
     383    58773355 :         if (ny == 1) return gc_const((pari_sp)(yd + 1), gen_0);
     384    49281062 :         ny--;
     385             :       }
     386             :     }
     387             :   }
     388   660598500 :   y[1] = evalsigne(1)|evallgefint(ny + 2);
     389   660598500 :   return y;
     390             : }
     391             : 
     392             : GEN
     393   137224679 : mantissa2nr(GEN x, long n)
     394             : {
     395   137224679 :   long ly, i, m, s = signe(x), lx = lg(x);
     396             :   GEN y;
     397   137224679 :   if (!s) return gen_0;
     398   137223330 :   if (!n)
     399             :   {
     400    30203774 :     y = cgeti(lx);
     401    30200632 :     y[1] = evalsigne(s) | evallgefint(lx);
     402    30200632 :     xmpn_mirrorcopy(LIMBS(y),RLIMBS(x),lx-2);
     403    30200390 :     return y;
     404             :   }
     405   107019556 :   if (n > 0)
     406             :   {
     407      218284 :     GEN z = (GEN)avma;
     408      218284 :     long d = dvmdsBIL(n, &m);
     409             : 
     410      218284 :     ly = lx+d; y = new_chunk(ly);
     411      556800 :     for ( ; d; d--) *--z = 0;
     412      222719 :     if (!m) for (i=2; i<lx; i++) y[i]=x[i];
     413             :     else
     414             :     {
     415      216842 :       const ulong sh = BITS_IN_LONG - m;
     416      216842 :       shift_left(y,x, 2,lx-1, 0,m);
     417      216842 :       i = uel(x,2) >> sh;
     418             :       /* Extend y on the left? */
     419      216842 :       if (i) { ly++; y = new_chunk(1); y[2] = i; }
     420             :     }
     421             :   }
     422             :   else
     423             :   {
     424   106801272 :     ly = lx - dvmdsBIL(-n, &m);
     425   106800192 :     if (ly<3) return gen_0;
     426   106800192 :     y = new_chunk(ly);
     427   106761158 :     if (m) {
     428   106502875 :       shift_right(y,x, 2,ly, 0,m);
     429   106540612 :       if (y[2] == 0)
     430             :       {
     431           0 :         if (ly==3) return gc_const((pari_sp)(y+3), gen_0);
     432           0 :         ly--; set_avma((pari_sp)(++y));
     433             :       }
     434             :     } else {
     435      701242 :       for (i=2; i<ly; i++) y[i]=x[i];
     436             :     }
     437             :   }
     438   107017179 :   xmpn_mirror(LIMBS(y),ly-2);
     439   107063283 :   y[1] = evalsigne(s)|evallgefint(ly);
     440   107063283 :   y[0] = evaltyp(t_INT)|evallg(ly); return y;
     441             : }
     442             : 
     443             : GEN
     444     3546570 : truncr(GEN x)
     445             : {
     446             :   long s, e, d, m, i;
     447             :   GEN y;
     448     3546570 :   if ((s=signe(x)) == 0 || (e=expo(x)) < 0) return gen_0;
     449     1513824 :   d = nbits2lg(e+1); m = remsBIL(e);
     450     1513813 :   if (d > lg(x)) pari_err_PREC( "truncr (precision loss in truncation)");
     451             : 
     452     1513809 :   y=cgeti(d); y[1] = evalsigne(s) | evallgefint(d);
     453     1513746 :   if (++m == BITS_IN_LONG)
     454       95961 :     for (i=2; i<d; i++) y[d-i+1]=x[i];
     455             :   else
     456             :   {
     457     1466024 :     GEN z=cgeti(d);
     458     2998619 :     for (i=2; i<d; i++) z[d-i+1]=x[i];
     459     1465910 :     mpn_rshift(LIMBS(y),LIMBS(z),d-2,BITS_IN_LONG-m);
     460     1465910 :     set_avma((pari_sp)y);
     461             :   }
     462     1513612 :   return y;
     463             : }
     464             : 
     465             : /* integral part */
     466             : GEN
     467     7001065 : floorr(GEN x)
     468             : {
     469             :   long e, d, m, i, lx;
     470             :   GEN y;
     471     7001065 :   if (signe(x) >= 0) return truncr(x);
     472     4206854 :   if ((e=expo(x)) < 0) return gen_m1;
     473     3520211 :   d = nbits2lg(e+1); m = remsBIL(e);
     474     3518668 :   lx=lg(x); if (d>lx) pari_err_PREC( "floorr (precision loss in truncation)");
     475     3518668 :   y = cgeti(d+1);
     476     3507400 :   if (++m == BITS_IN_LONG)
     477             :   {
     478        3053 :     for (i=2; i<d; i++) y[d-i+1]=x[i];
     479        1450 :     i=d; while (i<lx && !x[i]) i++;
     480         906 :     if (i==lx) goto END;
     481             :   }
     482             :   else
     483             :   {
     484     3506494 :     GEN z=cgeti(d);
     485     7841933 :     for (i=2; i<d; i++) z[d-i+1]=x[i];
     486     3492746 :     mpn_rshift(LIMBS(y),LIMBS(z),d-2,BITS_IN_LONG-m);
     487     3493524 :     if (uel(x,d-1)<<m == 0)
     488             :     {
     489      513576 :       i=d; while (i<lx && !x[i]) i++;
     490      118094 :       if (i==lx) goto END;
     491             :     }
     492             :   }
     493     3418520 :   if (mpn_add_1(LIMBS(y),LIMBS(y),d-2,1))
     494           0 :     y[d++]=1;
     495     3418369 : END:
     496     3494279 :   y[1] = evalsigne(-1) | evallgefint(d);
     497     3494279 :   return y;
     498             : }
     499             : 
     500             : INLINE int
     501  3941502325 : cmpiispec(GEN x, GEN y, long lx, long ly)
     502             : {
     503  3941502325 :   if (lx < ly) return -1;
     504  3622308085 :   if (lx > ly) return  1;
     505  3473099965 :   return mpn_cmp((mp_limb_t*)x,(mp_limb_t*)y, lx);
     506             : }
     507             : 
     508             : INLINE int
     509   267251605 : equaliispec(GEN x, GEN y, long lx, long ly)
     510             : {
     511   267251605 :   if (lx != ly) return 0;
     512   267110596 :   return !mpn_cmp((mp_limb_t*)x,(mp_limb_t*)y, lx);
     513             : }
     514             : 
     515             : /***********************************************************************/
     516             : /**                                                                   **/
     517             : /**                          MULTIPLICATION                           **/
     518             : /**                                                                   **/
     519             : /***********************************************************************/
     520             : /* assume ny > 0 */
     521             : INLINE GEN
     522  5544490369 : muluispec(ulong x, GEN y, long ny)
     523             : {
     524  5544490369 :   if (ny == 1)
     525  4453609182 :     return muluu(x, *y);
     526             :   else
     527             :   {
     528  1090881187 :     long lz = ny+3;
     529  1090881187 :     GEN z = cgeti(lz);
     530  1100217052 :     ulong hi = mpn_mul_1 (LIMBS(z), (mp_limb_t *)y, ny, x);
     531  1101614977 :     if (hi) { z[lz - 1] = hi; } else lz--;
     532  1101614977 :     z[1] = evalsigne(1) | evallgefint(lz);
     533  1101614977 :     return z;
     534             :   }
     535             : }
     536             : 
     537             : /* a + b*|y| */
     538             : GEN
     539           0 : addumului(ulong a, ulong b, GEN y)
     540             : {
     541             :   GEN z;
     542             :   long i, lz;
     543             :   ulong hi;
     544           0 :   if (!b || !signe(y)) return utoi(a);
     545           0 :   lz = lgefint(y)+1;
     546           0 :   z = cgeti(lz);
     547           0 :   z[2]=a;
     548           0 :   for(i=3;i<lz;i++) z[i]=0;
     549           0 :   hi=mpn_addmul_1(LIMBS(z), LIMBS(y), NLIMBS(y), b);
     550           0 :   if (hi) z[lz-1]=hi; else lz--;
     551           0 :   z[1] = evalsigne(1) | evallgefint(lz);
     552           0 :   return gc_const((pari_sp)z, z);
     553             : }
     554             : 
     555             : /***********************************************************************/
     556             : /**                                                                   **/
     557             : /**                          DIVISION                                 **/
     558             : /**                                                                   **/
     559             : /***********************************************************************/
     560             : 
     561             : ulong
     562  1303173910 : umodiu(GEN y, ulong x)
     563             : {
     564  1303173910 :   long sy=signe(y);
     565             :   ulong hi;
     566  1303173910 :   if (!x) pari_err_INV("umodiu",gen_0);
     567  1304105443 :   if (!sy) return 0;
     568   895899903 :   hi = mpn_mod_1(LIMBS(y),NLIMBS(y),x);
     569   895899903 :   if (!hi) return 0;
     570   879702095 :   return (sy > 0)? hi: x - hi;
     571             : }
     572             : 
     573             : /* return |y| \/ x */
     574             : GEN
     575   109856686 : absdiviu_rem(GEN y, ulong x, ulong *rem)
     576             : {
     577             :   long ly;
     578             :   GEN z;
     579             : 
     580   109856686 :   if (!x) pari_err_INV("absdiviu_rem",gen_0);
     581   109862349 :   if (!signe(y)) { *rem = 0; return gen_0; }
     582             : 
     583    83122180 :   ly = lgefint(y);
     584    83122180 :   if (ly == 3 && (ulong)x > uel(y,2)) { *rem = uel(y,2); return gen_0; }
     585             : 
     586    68657570 :   z = cgeti(ly);
     587    68656718 :   *rem = mpn_divrem_1(LIMBS(z), 0, LIMBS(y), NLIMBS(y), x);
     588    68657043 :   if (z [ly - 1] == 0) ly--;
     589    68657043 :   z[1] = evallgefint(ly) | evalsigne(1);
     590    68657043 :   return z;
     591             : }
     592             : 
     593             : GEN
     594    85237945 : divis_rem(GEN y, long x, long *rem)
     595             : {
     596    85237945 :   long sy=signe(y),ly,s;
     597             :   GEN z;
     598             : 
     599    85237945 :   if (!x) pari_err_INV("divis_rem",gen_0);
     600    85247744 :   if (!sy) { *rem = 0; return gen_0; }
     601    59906489 :   if (x<0) { s = -sy; x = -x; } else s = sy;
     602             : 
     603    59906489 :   ly = lgefint(y);
     604    59906489 :   if (ly == 3 && (ulong)x > uel(y,2)) { *rem = itos(y); return gen_0; }
     605             : 
     606    41464346 :   z = cgeti(ly);
     607    41458479 :   *rem = mpn_divrem_1(LIMBS(z), 0, LIMBS(y), NLIMBS(y), x);
     608    41458605 :   if (sy<0) *rem = -  *rem;
     609    41458605 :   if (z[ly - 1] == 0) ly--;
     610    41458605 :   z[1] = evallgefint(ly) | evalsigne(s);
     611    41458605 :   return z;
     612             : }
     613             : 
     614             : GEN
     615      966290 : divis(GEN y, long x)
     616             : {
     617      966290 :   long sy=signe(y),ly,s;
     618             :   GEN z;
     619             : 
     620      966290 :   if (!x) pari_err_INV("divis",gen_0);
     621      966290 :   if (!sy) return gen_0;
     622      966242 :   if (x<0) { s = -sy; x = -x; } else s=sy;
     623             : 
     624      966242 :   ly = lgefint(y);
     625      966242 :   if (ly == 3 && (ulong)x > uel(y,2)) return gen_0;
     626             : 
     627      965930 :   z = cgeti(ly);
     628      965928 :   (void)mpn_divrem_1(LIMBS(z), 0, LIMBS(y), NLIMBS(y), x);
     629      965928 :   if (z[ly - 1] == 0) ly--;
     630      965928 :   z[1] = evallgefint(ly) | evalsigne(s);
     631      965928 :   return z;
     632             : }
     633             : 
     634             : /* We keep llx bits of x and lly bits of y*/
     635             : static GEN
     636    75709108 : divrr_with_gmp(GEN x, GEN y)
     637             : {
     638    75709108 :   long lx=RNLIMBS(x),ly=RNLIMBS(y);
     639    75709108 :   long lw=minss(lx,ly);
     640    75709930 :   long lly=minss(lw+1,ly);
     641    75710301 :   GEN  w = cgetg(lw+2, t_REAL);
     642    75692575 :   long lu=lw+lly;
     643    75692575 :   long llx=minss(lu,lx);
     644    75689683 :   mp_limb_t *u=(mp_limb_t *)new_chunk(lu);
     645    75669000 :   mp_limb_t *z=(mp_limb_t *)new_chunk(lly);
     646             :   mp_limb_t *q,*r;
     647    75636337 :   long e=expo(x)-expo(y);
     648    75636337 :   long sx=signe(x),sy=signe(y);
     649    75636337 :   xmpn_mirrorcopy(z,RLIMBS(y),lly);
     650    75650356 :   xmpn_mirrorcopy(u+lu-llx,RLIMBS(x),llx);
     651    75697773 :   xmpn_zero(u,lu-llx);
     652    75725560 :   q = (mp_limb_t *)new_chunk(lw+1);
     653    75709505 :   r = (mp_limb_t *)new_chunk(lly);
     654             : 
     655    75685507 :   mpn_tdiv_qr(q,r,0,u,lu,z,lly);
     656             : 
     657             :   /*Round up: This is not exactly correct we should test 2*r>z*/
     658    75735612 :   if (uel(r,lly-1) > (uel(z,lly-1)>>1))
     659    37522458 :     mpn_add_1(q,q,lw+1,1);
     660             : 
     661    75735639 :   xmpn_mirrorcopy(RLIMBS(w),q,lw);
     662             : 
     663    75734103 :   if (q[lw] == 0) e--;
     664    41920573 :   else if (q[lw] == 1) { shift_right(w,w, 2,lw+2, 1,1); }
     665           0 :   else { w[2] = HIGHBIT; e++; }
     666    75733083 :   if (sy < 0) sx = -sx;
     667    75733083 :   w[1] = evalsigne(sx) | evalexpo(e);
     668    75726421 :   return gc_const((pari_sp)w, w);
     669             : }
     670             : 
     671             : /* We keep llx bits of x and lly bits of y*/
     672             : static GEN
     673    35195532 : divri_with_gmp(GEN x, GEN y)
     674             : {
     675    35195532 :   long llx=RNLIMBS(x),ly=NLIMBS(y);
     676    35195532 :   long lly=minss(llx+1,ly);
     677    35195652 :   GEN  w = cgetg(llx+2, t_REAL);
     678    35191350 :   long lu=llx+lly, ld=ly-lly;
     679    35191350 :   mp_limb_t *u=(mp_limb_t *)new_chunk(lu);
     680    35187664 :   mp_limb_t *z=(mp_limb_t *)new_chunk(lly);
     681             :   mp_limb_t *q,*r;
     682    35181884 :   long sh=bfffo(y[ly+1]);
     683    35181884 :   long e=expo(x)-expi(y);
     684    35182733 :   long sx=signe(x),sy=signe(y);
     685    35182733 :   if (sh) mpn_lshift(z,LIMBS(y)+ld,lly,sh);
     686      683755 :   else xmpn_copy(z,LIMBS(y)+ld,lly);
     687    35184805 :   xmpn_mirrorcopy(u+lu-llx,RLIMBS(x),llx);
     688    35192279 :   xmpn_zero(u,lu-llx);
     689    35197915 :   q = (mp_limb_t *)new_chunk(llx+1);
     690    35195440 :   r = (mp_limb_t *)new_chunk(lly);
     691             : 
     692    35190700 :   mpn_tdiv_qr(q,r,0,u,lu,z,lly);
     693             : 
     694             :   /*Round up: This is not exactly correct we should test 2*r>z*/
     695    35200536 :   if (uel(r,lly-1) > (uel(z,lly-1)>>1))
     696    16017684 :     mpn_add_1(q,q,llx+1,1);
     697             : 
     698    35200535 :   xmpn_mirrorcopy(RLIMBS(w),q,llx);
     699             : 
     700    35200056 :   if (q[llx] == 0) e--;
     701    10600284 :   else if (q[llx] == 1) { shift_right(w,w, 2,llx+2, 1,1); }
     702           0 :   else { w[2] = HIGHBIT; e++; }
     703    35199962 :   if (sy < 0) sx = -sx;
     704    35199962 :   w[1] = evalsigne(sx) | evalexpo(e);
     705    35198995 :   return gc_const((pari_sp)w, w);
     706             : }
     707             : 
     708             : GEN
     709   150466168 : divri(GEN x, GEN y)
     710             : {
     711   150466168 :   long  s = signe(y);
     712             : 
     713   150466168 :   if (!s) pari_err_INV("divri",gen_0);
     714   150467004 :   if (!signe(x)) return real_0_bit(expo(x) - expi(y));
     715   150235430 :   if (!is_bigint(y)) {
     716   115045055 :     GEN z = divru(x, y[2]);
     717   115042727 :     if (s < 0) togglesign(z);
     718   115042737 :     return z;
     719             :   }
     720    35189063 :   return divri_with_gmp(x,y);
     721             : }
     722             : 
     723             : GEN
     724   141744246 : divrr(GEN x, GEN y)
     725             : {
     726   141744246 :   long sx=signe(x), sy=signe(y), lx,ly,lr,e,i,j;
     727             :   ulong y0,y1;
     728             :   GEN r, r1;
     729             : 
     730   141744246 :   if (!sy) pari_err_INV("divrr",y);
     731   141752054 :   e = expo(x) - expo(y);
     732   141752054 :   if (!sx) return real_0_bit(e);
     733   141247454 :   if (sy<0) sx = -sx;
     734             : 
     735   141247454 :   lx=lg(x); ly=lg(y);
     736   141247454 :   if (ly==3)
     737             :   {
     738    26291669 :     ulong k = x[2], l = (lx>3)? x[3]: 0;
     739             :     LOCAL_HIREMAINDER;
     740    26291669 :     if (k < uel(y,2)) e--;
     741             :     else
     742             :     {
     743     8278611 :       l >>= 1; if (k&1) l |= HIGHBIT;
     744     8278611 :       k >>= 1;
     745             :     }
     746    26291669 :     hiremainder = k; k = divll(l,y[2]);
     747    26291669 :     if (hiremainder > (uel(y,2) >> 1) && !++k) { k = HIGHBIT; e++; }
     748    26291669 :     r = cgetg(3, t_REAL);
     749    26291077 :     r[1] = evalsigne(sx) | evalexpo(e);
     750    26289036 :     r[2] = k; return r;
     751             :   }
     752             : 
     753   114955785 :   if (ly >= prec2lg(DIVRR_GMP_LIMIT))
     754    75709617 :     return divrr_with_gmp(x,y);
     755             : 
     756    39294483 :   lr = minss(lx,ly); r = new_chunk(lr);
     757    39295167 :   r1 = r-1;
     758   133461782 :   r1[1] = 0; for (i=2; i<lr; i++) r1[i]=x[i];
     759    39295167 :   r1[lr] = (lx>ly)? x[lr]: 0;
     760    39295167 :   y0 = y[2]; y1 = y[3];
     761   172713017 :   for (i=0; i<lr-1; i++)
     762             :   { /* r1 = r + (i-1), OK up to r1[2] (accesses at most r[lr]) */
     763             :     ulong k, qp;
     764             :     LOCAL_HIREMAINDER;
     765             :     LOCAL_OVERFLOW;
     766             : 
     767   133417850 :     if (uel(r1,1) == y0) { qp = ULONG_MAX; k = addll(y0,r1[2]); }
     768             :     else
     769             :     {
     770   133171523 :       if (uel(r1,1) > y0) /* can't happen if i=0 */
     771             :       {
     772           0 :         GEN y1 = y+1;
     773           0 :         j = lr-i; r1[j] = subll(r1[j],y1[j]);
     774           0 :         for (j--; j>0; j--) r1[j] = subllx(r1[j],y1[j]);
     775           0 :         j=i; do uel(r,--j)++; while (j && !r[j]);
     776             :       }
     777   133171523 :       hiremainder = r1[1]; overflow = 0;
     778   133171523 :       qp = divll(r1[2],y0); k = hiremainder;
     779             :     }
     780   133417850 :     j = lr-i+1;
     781   133417850 :     if (!overflow)
     782             :     {
     783             :       long k3, k4;
     784   133263535 :       k3 = mulll(qp,y1);
     785   133263535 :       if (j == 3) /* i = lr - 2 maximal, r1[3] undefined -> 0 */
     786    39350012 :         k4 = subll(hiremainder,k);
     787             :       else
     788             :       {
     789    93913523 :         k3 = subll(k3, r1[3]);
     790    93913523 :         k4 = subllx(hiremainder,k);
     791             :       }
     792   170960823 :       while (!overflow && k4) { qp--; k3=subll(k3,y1); k4=subllx(k4,y0); }
     793             :     }
     794   133417850 :     if (j<ly) (void)mulll(qp,y[j]); else { hiremainder = 0 ; j = ly; }
     795   395846622 :     for (j--; j>1; j--)
     796             :     {
     797   262428772 :       r1[j] = subll(r1[j], addmul(qp,y[j]));
     798   262428772 :       hiremainder += overflow;
     799             :     }
     800   133417850 :     if (uel(r1,1) != hiremainder)
     801             :     {
     802      182776 :       if (uel(r1,1) < hiremainder)
     803             :       {
     804      182776 :         qp--;
     805      182776 :         j = lr-i-(lr-i>=ly); r1[j] = addll(r1[j], y[j]);
     806      520872 :         for (j--; j>1; j--) r1[j] = addllx(r1[j], y[j]);
     807             :       }
     808             :       else
     809             :       {
     810           0 :         uel(r1,1) -= hiremainder;
     811           0 :         while (r1[1])
     812             :         {
     813           0 :           qp++; if (!qp) { j=i; do uel(r,--j)++; while (j && !r[j]); }
     814           0 :           j = lr-i-(lr-i>=ly); r1[j] = subll(r1[j],y[j]);
     815           0 :           for (j--; j>1; j--) r1[j] = subllx(r1[j],y[j]);
     816           0 :           uel(r1,1) -= overflow;
     817             :         }
     818             :       }
     819             :     }
     820   133417850 :     *++r1 = qp;
     821             :   }
     822             :   /* i = lr-1 */
     823             :   /* round correctly */
     824    39295167 :   if (uel(r1,1) > (y0>>1))
     825             :   {
     826    19808951 :     j=i; do uel(r,--j)++; while (j && !r[j]);
     827             :   }
     828   133640837 :   r1 = r-1; for (j=i; j>=2; j--) r[j]=r1[j];
     829    39295167 :   if (r[0] == 0) e--;
     830    14179625 :   else if (r[0] == 1) { shift_right(r,r, 2,lr, 1,1); }
     831             :   else { /* possible only when rounding up to 0x2 0x0 ... */
     832          18 :     r[2] = (long)HIGHBIT; e++;
     833             :   }
     834    39294386 :   r[0] = evaltyp(t_REAL)|evallg(lr);
     835    39351427 :   r[1] = evalsigne(sx) | evalexpo(e);
     836    39341393 :   return r;
     837             : }
     838             : 
     839             : /* Integer division x / y: such that sign(r) = sign(x)
     840             :  *   if z = ONLY_REM return remainder, otherwise return quotient
     841             :  *   if z != NULL set *z to remainder
     842             :  *   *z is the last object on stack and can be disposed of with cgiv
     843             :  * If *z is zero, we put gen_0 here and no copy.
     844             :  * space needed: lx + ly
     845             :  */
     846             : GEN
     847  2177334581 : dvmdii(GEN x, GEN y, GEN *z)
     848             : {
     849  2177334581 :   long sx=signe(x),sy=signe(y);
     850             :   long lx, ly, lq;
     851             :   pari_sp av;
     852             :   GEN r,q;
     853             : 
     854  2177334581 :   if (!sy) pari_err_INV("dvmdii",y);
     855  2178449162 :   if (!sx)
     856             :   {
     857    68105849 :     if (!z || z == ONLY_REM) return gen_0;
     858    40578178 :     *z=gen_0; return gen_0;
     859             :   }
     860  2110343313 :   lx=lgefint(x);
     861  2110343313 :   ly=lgefint(y); lq=lx-ly;
     862  2110343313 :   if (lq <= 0)
     863             :   {
     864  1249942345 :     if (lq == 0)
     865             :     {
     866  1061993285 :       long s=mpn_cmp(LIMBS(x),LIMBS(y),NLIMBS(x));
     867  1061993285 :       if (s>0) goto DIVIDE;
     868   363509408 :       if (s==0)
     869             :       {
     870    32675595 :         if (z == ONLY_REM) return gen_0;
     871    22136349 :         if (z) *z = gen_0;
     872    22136349 :         if (sx < 0) sy = -sy;
     873    22136349 :         return stoi(sy);
     874             :       }
     875             :     }
     876   518782873 :     if (z == ONLY_REM) return icopy(x);
     877    14028855 :     if (z) *z = icopy(x);
     878    14028855 :     return gen_0;
     879             :   }
     880   860400968 : DIVIDE: /* quotient is nonzero */
     881  1558884845 :   av=avma; if (sx<0) sy = -sy;
     882  1558884845 :   if (ly==3)
     883             :   {
     884   570693947 :     ulong lq = lx;
     885             :     ulong si;
     886   570693947 :     q  = cgeti(lq);
     887   570233302 :     si = mpn_divrem_1(LIMBS(q), 0, LIMBS(x), NLIMBS(x), y[2]);
     888   570776684 :     if (q[lq - 1] == 0) lq--;
     889   570776684 :     if (z == ONLY_REM)
     890             :     {
     891   331599145 :       if (!si) return gc_const(av, gen_0);
     892   286109642 :       set_avma(av); r = cgeti(3);
     893   285590980 :       r[1] = evalsigne(sx) | evallgefint(3);
     894   285590980 :       r[2] = si; return r;
     895             :     }
     896   239177539 :     q[1] = evalsigne(sy) | evallgefint(lq);
     897   239177539 :     if (!z) return q;
     898   235074111 :     if (!si) { *z=gen_0; return q; }
     899    63682072 :     r=cgeti(3);
     900    63708578 :     r[1] = evalsigne(sx) | evallgefint(3);
     901    63708578 :     r[2] = si; *z=r; return q;
     902             :   }
     903   988190898 :   if (z == ONLY_REM)
     904             :   {
     905   965157238 :     ulong lr = lgefint(y);
     906   965157238 :     ulong lq = lgefint(x)-lgefint(y)+3;
     907   965157238 :     GEN r = cgeti(lr);
     908   958954141 :     GEN q = cgeti(lq);
     909   952402287 :     mpn_tdiv_qr(LIMBS(q), LIMBS(r),0, LIMBS(x), NLIMBS(x), LIMBS(y), NLIMBS(y));
     910   965893092 :     if (!r[lr - 1])
     911             :     {
     912  1180357198 :       while(lr>2 && !r[lr - 1]) lr--;
     913   544823558 :       if (lr == 2) return gc_const(av, gen_0); /* exact division */
     914             :     }
     915   951004178 :     r[1] = evalsigne(sx) | evallgefint(lr);
     916   951004178 :     return gc_const((pari_sp)r, r);
     917             :   }
     918             :   else
     919             :   {
     920    23033660 :     ulong lq = lgefint(x)-lgefint(y)+3;
     921    23033660 :     ulong lr = lgefint(y);
     922    23033660 :     GEN q = cgeti(lq);
     923    28465442 :     GEN r = cgeti(lr);
     924    28451187 :     mpn_tdiv_qr(LIMBS(q), LIMBS(r),0, LIMBS(x), NLIMBS(x), LIMBS(y), NLIMBS(y));
     925    28477613 :     if (q[lq - 1] == 0) lq--;
     926    28477613 :     q[1] = evalsigne(sy) | evallgefint(lq);
     927    28477613 :     if (!z) return gc_const((pari_sp)q, q);
     928    24784997 :     if (!r[lr - 1])
     929             :     {
     930    38517656 :       while(lr>2 && !r[lr - 1]) lr--;
     931     6343535 :       if (lr == 2) { *z = gen_0; return gc_const((pari_sp)q, q); } /* exact */
     932             :     }
     933    19996014 :     r[1] = evalsigne(sx) | evallgefint(lr);
     934    19996014 :     *z = gc_const((pari_sp)r, r); return q;
     935             :   }
     936             : }
     937             : 
     938             : /* Montgomery reduction.
     939             :  * N has k words, assume T >= 0 has less than 2k.
     940             :  * Return res := T / B^k mod N, where B = 2^BIL
     941             :  * such that 0 <= res < T/B^k + N  and  res has less than k words */
     942             : GEN
     943    33681769 : red_montgomery(GEN T, GEN N, ulong inv)
     944             : {
     945             :   pari_sp av;
     946             :   GEN Te, Td, Ne, Nd, scratch;
     947    33681769 :   ulong i, j, m, t, d, k = NLIMBS(N);
     948             :   int carry;
     949             :   LOCAL_HIREMAINDER;
     950             :   LOCAL_OVERFLOW;
     951             : 
     952    33681769 :   if (k == 0) return gen_0;
     953    33681769 :   d = NLIMBS(T); /* <= 2*k */
     954    33681769 :   if (d == 0) return gen_0;
     955             : #ifdef DEBUG
     956             :   if (d > 2*k) pari_err_BUG("red_montgomery");
     957             : #endif
     958    33681760 :   if (k == 1)
     959             :   { /* as below, special cased for efficiency */
     960      163341 :     ulong n = uel(N,2);
     961      163341 :     if (d == 1) {
     962      163194 :       hiremainder = uel(T,2);
     963      163194 :       m = hiremainder * inv;
     964      163194 :       (void)addmul(m, n); /* t + m*n = 0 */
     965      163194 :       return utoi(hiremainder);
     966             :     } else { /* d = 2 */
     967         147 :       hiremainder = uel(T,2);
     968         147 :       m = hiremainder * inv;
     969         147 :       (void)addmul(m, n); /* t + m*n = 0 */
     970         147 :       t = addll(hiremainder, uel(T,3));
     971         147 :       if (overflow) t -= n; /* t > n doesn't fit in 1 word */
     972         147 :       return utoi(t);
     973             :     }
     974             :   }
     975             :   /* assume k >= 2 */
     976    33518419 :   av = avma; scratch = new_chunk(k<<1); /* >= k + 2: result fits */
     977             : 
     978             :   /* copy T to scratch space (pad with zeroes to 2k words) */
     979    33183186 :   Td = scratch;
     980    33183186 :   Te = T + 2;
     981   454486638 :   for (i=0; i < d     ; i++) *Td++ = *Te++;
     982    59926460 :   for (   ; i < (k<<1); i++) *Td++ = 0;
     983             : 
     984    33183186 :   Te = scratch - 1; /* 1 beyond end of current T mantissa (in scratch) */
     985    33183186 :   Ne = N + 1;       /* 1 beyond end of N mantissa */
     986             : 
     987    33183186 :   carry = 0;
     988   244250688 :   for (i=0; i<k; i++) /* set T := T/B nod N, k times */
     989             :   {
     990   211067502 :     Td = Te; /* one beyond end of (new) T mantissa */
     991   211067502 :     Nd = Ne;
     992   211067502 :     hiremainder = *++Td;
     993   211067502 :     m = hiremainder * inv; /* solve T + m N = O(B) */
     994             : 
     995             :     /* set T := (T + mN) / B */
     996   211067502 :     Te = Td;
     997   211067502 :     (void)addmul(m, *++Nd); /* = 0 */
     998  1803835408 :     for (j=1; j<k; j++)
     999             :     {
    1000  1592767906 :       t = addll(addmul(m, *++Nd), *++Td);
    1001  1592767906 :       *Td = t;
    1002  1592767906 :       hiremainder += overflow;
    1003             :     }
    1004   211067502 :     t = addll(hiremainder, *++Td); *Td = t + carry;
    1005   211067502 :     carry = (overflow || (carry && *Td == 0));
    1006             :   }
    1007    33183186 :   if (carry)
    1008             :   { /* Td > N overflows (k+1 words), set Td := Td - N */
    1009       62483 :     GEN NE = N + k+1;
    1010       62483 :     Td = Te;
    1011       62483 :     Nd = Ne;
    1012       62483 :     t = subll(*++Td, *++Nd); *Td = t;
    1013      607978 :     while (Nd < NE) { t = subllx(*++Td, *++Nd); *Td = t; }
    1014             :   }
    1015             : 
    1016             :   /* copy result */
    1017    33183186 :   Td = (GEN)av - 1; /* *Td = high word of final result */
    1018    36653041 :   while (*Td == 0 && Te < Td) Td--; /* strip leading 0s */
    1019    33183186 :   k = Td - Te; if (!k) return gc_const(av, gen_0);
    1020    33183186 :   Td = (GEN)av - k; /* will write mantissa there */
    1021    33183186 :   (void)memmove(Td, Te+1, k*sizeof(long));
    1022    33183186 :   Td -= 2;
    1023    33183186 :   Td[0] = evaltyp(t_INT) | evallg(k+2);
    1024    33441942 :   Td[1] = evalsigne(1) | evallgefint(k+2);
    1025             : #ifdef DEBUG
    1026             : {
    1027             :   long l = NLIMBS(N), s = BITS_IN_LONG*l;
    1028             :   GEN R = int2n(s);
    1029             :   GEN res = remii(mulii(T, Fp_inv(R, N)), N);
    1030             :   if (k > lgefint(N)
    1031             :     || !equalii(remii(Td,N),res)
    1032             :     || cmpii(Td, addii(shifti(T, -s), N)) >= 0) pari_err_BUG("red_montgomery");
    1033             : }
    1034             : #endif
    1035    33441942 :   return gc_const((pari_sp)Td, Td);
    1036             : }
    1037             : 
    1038             : /* EXACT INTEGER DIVISION */
    1039             : 
    1040             : /* use undocumented GMP interface */
    1041             : static void
    1042   115116054 : GEN2mpz(mpz_t X, GEN x)
    1043             : {
    1044   115116054 :   long l = lgefint(x)-2;
    1045   115116054 :   X->_mp_alloc = l;
    1046   115116054 :   X->_mp_size = signe(x) > 0? l: -l;
    1047   115116054 :   X->_mp_d = LIMBS(x);
    1048   115116054 : }
    1049             : static void
    1050    57559286 : mpz2GEN(GEN z, mpz_t Z)
    1051             : {
    1052    57559286 :   long l = Z->_mp_size;
    1053    57559286 :   z[1] = evalsigne(l > 0? 1: -1) | evallgefint(labs(l)+2);
    1054    57559286 : }
    1055             : 
    1056             : #ifdef mpn_divexact_1
    1057             : static GEN
    1058   413532738 : diviuexact_i(GEN x, ulong y)
    1059             : {
    1060   413532738 :   long l = lgefint(x);
    1061   413532738 :   GEN z = cgeti(l);
    1062   413060190 :   mpn_divexact_1(LIMBS(z), LIMBS(x), NLIMBS(x), y);
    1063   413121328 :   if (z[l-1] == 0) l--;
    1064   413121328 :   z[1] = evallgefint(l) | evalsigne(signe(x));
    1065   413121328 :   return z;
    1066             : }
    1067             : #elif 1 && !defined(_WIN64) /* mpz_divexact_ui is not LLP64 friendly   */
    1068             :                             /* assume y != 0 and the division is exact */
    1069             : static GEN
    1070             : diviuexact_i(GEN x, ulong y)
    1071             : {
    1072             :   long l = lgefint(x);
    1073             :   GEN z = cgeti(l);
    1074             :   mpz_t X, Z;
    1075             :   GEN2mpz(X, x);
    1076             :   Z->_mp_alloc = l-2;
    1077             :   Z->_mp_size  = l-2;
    1078             :   Z->_mp_d = LIMBS(z);
    1079             :   mpz_divexact_ui(Z, X, y);
    1080             :   mpz2GEN(z, Z); return z;
    1081             : }
    1082             : #else
    1083             : /* assume y != 0 and the division is exact */
    1084             : static GEN
    1085             : diviuexact_i(GEN x, ulong y)
    1086             : {
    1087             :   /*TODO: implement true exact division.*/
    1088             :   return divii(x,utoi(y));
    1089             : }
    1090             : #endif
    1091             : 
    1092             : GEN
    1093    39176416 : diviuexact(GEN x, ulong y)
    1094             : {
    1095             :   GEN z;
    1096    39176416 :   if (!signe(x)) return gen_0;
    1097    38047840 :   z = diviuexact_i(x, y);
    1098    38043767 :   if (lgefint(z) == 2) pari_err_OP("exact division", x, utoi(y));
    1099    38043804 :   return z;
    1100             : }
    1101             : 
    1102             : /* Find z such that x=y*z, knowing that y | x (unchecked) */
    1103             : GEN
    1104   525449894 : diviiexact(GEN x, GEN y)
    1105             : {
    1106             :   GEN z;
    1107   525449894 :   if (!signe(y)) pari_err_INV("diviiexact",y);
    1108   525625006 :   if (!signe(x)) return gen_0;
    1109   433247290 :   if (lgefint(y) == 3)
    1110             :   {
    1111   375709154 :     z = diviuexact_i(x, y[2]);
    1112   375286379 :     if (signe(y) < 0) togglesign(z);
    1113             :   }
    1114             :   else
    1115             :   {
    1116    57538136 :     long l = lgefint(x);
    1117             :     mpz_t X, Y, Z;
    1118    57538136 :     z = cgeti(l);
    1119    57558163 :     GEN2mpz(X, x);
    1120    57558095 :     GEN2mpz(Y, y);
    1121    57557982 :     Z->_mp_alloc = l-2;
    1122    57557982 :     Z->_mp_size  = l-2;
    1123    57557982 :     Z->_mp_d = LIMBS(z);
    1124    57557982 :     mpz_divexact(Z, X, Y);
    1125    57559301 :     mpz2GEN(z, Z);
    1126             :   }
    1127   432845650 :   if (lgefint(z) == 2) pari_err_OP("exact division", x, y);
    1128   432824013 :   return z;
    1129             : }
    1130             : 
    1131             : /* assume yz != and yz | x */
    1132             : GEN
    1133      200516 : diviuuexact(GEN x, ulong y, ulong z)
    1134             : {
    1135             :   long tmp[4];
    1136             :   ulong t;
    1137             :   LOCAL_HIREMAINDER;
    1138      200516 :   t = mulll(y, z);
    1139      200516 :   if (!hiremainder) return diviuexact(x, t);
    1140           0 :   tmp[0] = evaltyp(t_INT)|_evallg(4);
    1141           0 :   tmp[1] = evalsigne(1)|evallgefint(4);
    1142           0 :   tmp[2] = t;
    1143           0 :   tmp[3] = hiremainder;
    1144           0 :   return diviiexact(x, tmp);
    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  5877490403 : muliispec(GEN x, GEN y, long nx, long ny)
    1156             : {
    1157             :   GEN zd;
    1158             :   long lz;
    1159             :   ulong hi;
    1160             : 
    1161  5877490403 :   if (nx < ny) swapspec(x,y, nx,ny);
    1162  5877490403 :   if (!ny) return gen_0;
    1163  5877490403 :   if (ny == 1) return muluispec((ulong)*y, x, nx);
    1164             : 
    1165  1045050661 :   lz = nx+ny+2;
    1166  1045050661 :   zd = cgeti(lz);
    1167  1047964023 :   hi = mpn_mul(LIMBS(zd), (mp_limb_t *)x, nx, (mp_limb_t *)y, ny);
    1168  1054828186 :   if (!hi) lz--;
    1169             :   /*else zd[lz-1]=hi; GH tell me it is not necessary.*/
    1170             : 
    1171  1054828186 :   zd[1] = evalsigne(1) | evallgefint(lz);
    1172  1054828186 :   return zd;
    1173             : }
    1174             : GEN
    1175      222737 : muluui(ulong x, ulong y, GEN z)
    1176             : {
    1177      222737 :   long t, s = signe(z);
    1178             :   GEN r;
    1179             :   LOCAL_HIREMAINDER;
    1180             : 
    1181      222737 :   if (!x || !y || !signe(z)) return gen_0;
    1182      222355 :   t = mulll(x,y);
    1183      222355 :   if (!hiremainder)
    1184      222370 :     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      222350 :   setsigne(r,s); return r;
    1193             : }
    1194             : 
    1195             : GEN
    1196  1019258231 : sqrispec(GEN x, long nx)
    1197             : {
    1198             :   GEN zd;
    1199             :   long lz;
    1200             : 
    1201  1019258231 :   if (!nx) return gen_0;
    1202   482713614 :   if (nx==1) return sqru(*x);
    1203             : 
    1204   273562415 :   lz = (nx<<1)+2;
    1205   273562415 :   zd = cgeti(lz);
    1206             : #ifdef mpn_sqr
    1207   270369026 :   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   273825507 :   if (zd[lz-1]==0) lz--;
    1212             : 
    1213   273825507 :   zd[1] = evalsigne(1) | evallgefint(lz);
    1214   273825507 :   return zd;
    1215             : }
    1216             : 
    1217             : INLINE GEN
    1218    41388474 : sqrispec_mirror(GEN x, long nx)
    1219             : {
    1220    41388474 :   GEN cx=new_chunk(nx);
    1221             :   GEN z;
    1222    41339286 :   xmpn_mirrorcopy((mp_limb_t *)cx,(mp_limb_t *)x,nx);
    1223    41452140 :   z=sqrispec(cx, nx);
    1224    41523357 :   xmpn_mirror(LIMBS(z), NLIMBS(z));
    1225    41523189 :   return z;
    1226             : }
    1227             : 
    1228             : /* leaves garbage on the stack. */
    1229             : INLINE GEN
    1230    84061691 : muliispec_mirror(GEN x, GEN y, long nx, long ny)
    1231             : {
    1232             :   GEN cx, cy, z;
    1233    84061691 :   long s = 0;
    1234   115145867 :   while (nx && x[nx-1]==0) { nx--; s++; }
    1235   121338563 :   while (ny && y[ny-1]==0) { ny--; s++; }
    1236    84061691 :   cx=new_chunk(nx); cy=new_chunk(ny);
    1237    83579329 :   xmpn_mirrorcopy((mp_limb_t *)cx,(mp_limb_t *)x,nx);
    1238    84482787 :   xmpn_mirrorcopy((mp_limb_t *)cy,(mp_limb_t *)y,ny);
    1239    85045382 :   z =  nx>=ny ? muliispec(cx, cy, nx, ny): muliispec(cy, cx, ny, nx);
    1240    85031458 :   if (s)
    1241             :   {
    1242     7731311 :     long i, lz = lgefint(z) + s;
    1243     7731311 :     (void)new_chunk(s);
    1244     7731308 :     z -= s;
    1245    76092352 :     for (i=0; i<s; i++) z[2+i]=0;
    1246     7731308 :     z[1] = evalsigne(1) | evallgefint(lz);
    1247     7731308 :     z[0] = evaltyp(t_INT) | evallg(lz);
    1248             :   }
    1249    85031455 :   xmpn_mirror(LIMBS(z), NLIMBS(z));
    1250    85729291 :   return z;
    1251             : }
    1252             : 
    1253             : /* x % (2^n), assuming n >= 0 */
    1254             : GEN
    1255    37315169 : remi2n(GEN x, long n)
    1256             : {
    1257             :   ulong hi;
    1258    37315169 :   long l, k, lx, ly, sx = signe(x);
    1259             :   GEN z, xd, zd;
    1260             : 
    1261    37315169 :   if (!sx || !n) return gen_0;
    1262             : 
    1263    36990842 :   k = dvmdsBIL(n, &l);
    1264    37003721 :   lx = lgefint(x);
    1265    37003721 :   if (lx < k+3) return icopy(x);
    1266             : 
    1267    36170621 :   xd = x + (2 + k);
    1268             :   /* x = |k|...|1|#|... : copy the last l bits of # and the first k words
    1269             :    *              ^--- initial xd  */
    1270    36170621 :   hi = ((ulong)*xd) & ((1UL<<l)-1); /* last l bits of # = top bits of result */
    1271    36170621 :   if (!hi)
    1272             :   { /* strip leading zeroes from result */
    1273     4105924 :     xd--; while (k && !*xd) { k--; xd--; }
    1274     3926784 :     if (!k) return gen_0;
    1275     2059702 :     ly = k+2;
    1276             :   }
    1277             :   else
    1278    32243837 :     ly = k+3;
    1279             : 
    1280    34303539 :   zd = z = cgeti(ly);
    1281    34267303 :   *++zd = evalsigne(sx) | evallgefint(ly);
    1282   481726017 :   xd = x+1; for ( ;k; k--) *++zd = *++xd;
    1283    34267303 :   if (hi) *++zd = hi;
    1284    34267303 :   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     5114590 : sqrtremi(GEN a, GEN *r)
    1298             : {
    1299     5114590 :   long l, na = NLIMBS(a);
    1300             :   mp_size_t nr;
    1301             :   GEN S;
    1302     5114590 :   if (!na) {
    1303         724 :     if (r) *r = gen_0;
    1304         724 :     return gen_0;
    1305             :   }
    1306     5113866 :   l = (na + 5) >> 1; /* 2 + ceil(na/2) */
    1307     5113866 :   S = cgetipos(l);
    1308     5113836 :   if (r) {
    1309     1308768 :     GEN R = cgeti(2 + na);
    1310     1308768 :     nr = mpn_sqrtrem(LIMBS(S), LIMBS(R), LIMBS(a), na);
    1311     1308768 :     if (nr) R[1] = evalsigne(1) | evallgefint(nr+2);
    1312       25116 :     else    { set_avma((pari_sp)S); R = gen_0; }
    1313     1308768 :     *r = R;
    1314             :   }
    1315             :   else
    1316     3805068 :     (void)mpn_sqrtrem(LIMBS(S), NULL, LIMBS(a), na);
    1317     5113838 :   return S;
    1318             : }
    1319             : 
    1320             : /* compute sqrt(|a|), assuming a != 0 */
    1321             : GEN
    1322   125198760 : sqrtr_abs(GEN a)
    1323             : {
    1324             :   GEN res;
    1325             :   mp_limb_t *b, *c;
    1326   125198760 :   long l = RNLIMBS(a), e = expo(a), er = e>>1;
    1327             :   long n;
    1328   125198760 :   res = cgetg(2 + l, t_REAL);
    1329   125114551 :   res[1] = evalsigne(1) | evalexpo(er);
    1330   125190638 :   if (e&1)
    1331             :   {
    1332    52607617 :     b = (mp_limb_t *) new_chunk(l<<1);
    1333    52592466 :     xmpn_zero(b,l);
    1334    52593183 :     xmpn_mirrorcopy(b+l, RLIMBS(a), l);
    1335    52606914 :     c = (mp_limb_t *) new_chunk(l);
    1336    52601502 :     n = mpn_sqrtrem(c,b,b,l<<1); /* c <- sqrt; b <- rem */
    1337    52635879 :     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    72583021 :     b = (mp_limb_t *) mantissa2nr(a,-1);
    1343    72616936 :     b[1] = uel(a,l+1)<<(BITS_IN_LONG-1);
    1344    72616936 :     b = (mp_limb_t *) new_chunk(l);
    1345    72593327 :     xmpn_zero(b,l+1); /* overwrites the former b[0] */
    1346    72596547 :     c = (mp_limb_t *) new_chunk(l + 1);
    1347    72558310 :     n = mpn_sqrtrem(c,b,b,(l<<1)+2); /* c <- sqrt; b <- rem */
    1348    72658247 :     u = (ulong)*c++;
    1349    72658247 :     if ( u&HIGHBIT || (u == ~HIGHBIT &&
    1350           0 :              (n>l || (n==l && mpn_cmp(b,c,l)>0))))
    1351    35804753 :       mpn_add_1(c,c,l,1);
    1352             :   }
    1353   125303652 :   xmpn_mirrorcopy(RLIMBS(res),c,l);
    1354   125278129 :   return gc_const((pari_sp)res, res);
    1355             : }
    1356             : 
    1357             : /* Normalize a nonnegative integer */
    1358             : GEN
    1359   306571927 : int_normalize(GEN x, long known_zero_words)
    1360             : {
    1361   306571927 :   long i =  lgefint(x) - 1 - known_zero_words;
    1362  2772241185 :   for ( ; i > 1; i--)
    1363  2720164975 :     if (x[i]) { setlgefint(x, i+1); return x; }
    1364    52076210 :   x[1] = evalsigne(0) | evallgefint(2); return x;
    1365             : }
    1366             : 
    1367             : /********************************************************************
    1368             :  **                                                                **
    1369             :  **                           Base Conversion                      **
    1370             :  **                                                                **
    1371             :  ********************************************************************/
    1372             : 
    1373             : ulong *
    1374      438548 : convi(GEN x, long *l)
    1375             : {
    1376      438548 :   long n = nchar2nlong(2 + (long)(NLIMBS(x) * (BITS_IN_LONG * LOG10_2)));
    1377      438548 :   GEN str = cgetg(n+1, t_VECSMALL);
    1378      438548 :   unsigned char *res = (unsigned char*) GSTR(str);
    1379      438548 :   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      438548 :   while (!*res) {res++; llz--;} /*Strip leading zeros*/
    1385      438548 :   lz  = (8+llz)/9;
    1386      438548 :   z = (ulong*)new_chunk(1+lz);
    1387      438548 :   t=res+llz+9;
    1388      869678 :   for(i=0;i<llz-8;i+=9)
    1389             :   {
    1390             :     ulong s;
    1391      431130 :     t-=18;
    1392      431130 :     s=*t++;
    1393     3880170 :     for (j=1; j<9;j++)
    1394     3449040 :       s=10*s+*t++;
    1395      431130 :     *z++=s;
    1396             :   }
    1397      438548 :   if (i<llz)
    1398             :   {
    1399      434597 :     unsigned char *t = res;
    1400      434597 :     ulong s=*t++;
    1401     1228637 :     for (j=i+1; j<llz;j++)
    1402      794040 :       s=10*s+*t++;
    1403      434597 :     *z++=s;
    1404             :   }
    1405      438548 :   *l = lz;
    1406      438548 :   return z;
    1407             : }

Generated by: LCOV version 1.16