Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - kernel/none - level1.h (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 20291-5fbfea9) Lines: 546 707 77.2 %
Date: 2017-02-25 05:49:34 Functions: 195 259 75.3 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : #line 2 "../src/kernel/none/level1.h"
       2             : /* Copyright (C) 2000  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             : /* This file defines "level 1" kernel functions.
      16             :  * These functions can be inline; they are also defined externally in
      17             :  * mpinl.c, which includes this file and never needs to be changed */
      18             : 
      19             : INLINE long
      20 38366259679 : evallg(long x)
      21             : {
      22 38366259679 :   if (x & ~LGBITS) pari_err_OVERFLOW("lg()");
      23 38365849964 :   return _evallg(x);
      24             : }
      25             : INLINE long
      26    13657008 : evalvalp(long x)
      27             : {
      28    13657008 :   long v = _evalvalp(x);
      29    13657008 :   if (v & ~VALPBITS) pari_err_OVERFLOW("valp()");
      30    13657005 :   return v;
      31             : }
      32             : INLINE long
      33  5007248154 : evalexpo(long x)
      34             : {
      35  5007248154 :   long v = _evalexpo(x);
      36  5007248154 :   if (v & ~EXPOBITS) pari_err_OVERFLOW("expo()");
      37  5007248154 :   return v;
      38             : }
      39             : INLINE long
      40     1844364 : evalprecp(long x)
      41             : {
      42     1844364 :   long v = _evalprecp(x);
      43     1844364 :   if (x & ~((1UL<<(BITS_IN_LONG-VALPnumBITS))-1)) pari_err_OVERFLOW("precp()");
      44     1844364 :   return v;
      45             : }
      46             : 
      47             : INLINE int
      48   137918800 : varncmp(long x, long y)
      49             : {
      50   137918800 :   if (varpriority[x] < varpriority[y]) return  1;
      51   107796314 :   if (varpriority[x] > varpriority[y]) return -1;
      52    57022787 :   return 0;
      53             : }
      54             : INLINE long
      55           0 : varnmin(long x, long y)
      56           0 : { return (varpriority[x] <= varpriority[y])? x: y; }
      57             : INLINE long
      58         196 : varnmax(long x, long y)
      59         196 : { return (varpriority[x] >= varpriority[y])? x: y; }
      60             : 
      61             : /* Inhibit some area gerepile-wise: declare it to be a non recursive
      62             :  * type, of length l. Thus gerepile won't inspect the zone, just copy it.
      63             :  * For the following situation:
      64             :  *   z = cgetg(t,a); av = avma; garbage(); ltop = avma;
      65             :  *   for (i=1; i<HUGE; i++) gel(z,i) = blah();
      66             :  *   stackdummy(av,ltop);
      67             :  * loses (av-ltop) words but save a costly gerepile. */
      68             : INLINE void
      69  1562913108 : stackdummy(pari_sp av, pari_sp ltop) {
      70  1562913108 :   long l = ((GEN)av) - ((GEN)ltop);
      71  1562913108 :   if (l > 0) {
      72   565186115 :     GEN z = (GEN)ltop;
      73   565186115 :     z[0] = evaltyp(t_VECSMALL) | evallg(l);
      74             : #ifdef DEBUG
      75             :     { long i; for (i = 1; i < l; i++) z[i] = 0; }
      76             : #endif
      77             :   }
      78  1562312597 : }
      79             : INLINE void
      80    38458414 : fixlg(GEN x, long ly) {
      81    38458414 :   long lx = lg(x), l = lx - ly;
      82    38458414 :   if (l > 0)
      83             :   { /* stackdummy(x+lx, x+ly) */
      84    23598400 :     GEN z = x + ly;
      85    23598400 :     z[0] = evaltyp(t_VECSMALL) | evallg(l);
      86    23598395 :     setlg(x, ly);
      87             : #ifdef DEBUG
      88             :     { long i; for (i = 1; i < l; i++) z[i] = 0; }
      89             : #endif
      90             :   }
      91    38458397 : }
      92             : /* update lg(z) before affrr(y, z)  [ to cater for precision loss ]*/
      93             : INLINE void
      94    14862453 : affrr_fixlg(GEN y, GEN z) { fixlg(z, lg(y)); affrr(y, z); }
      95             : 
      96             : /*******************************************************************/
      97             : /*                                                                 */
      98             : /*                       ALLOCATE ON STACK                         */
      99             : /*                                                                 */
     100             : /*******************************************************************/
     101             : INLINE GEN
     102 35318264109 : new_chunk(size_t x) /* x is a number of longs */
     103             : {
     104 35318264109 :   GEN z = ((GEN) avma) - x;
     105             :   CHECK_CTRLC
     106 35318264109 :   if (x > (avma-pari_mainstack->bot) / sizeof(long))
     107           7 :     new_chunk_resize(x);
     108 35326150357 :   avma = (pari_sp)z;
     109             : #ifdef MEMSTEP
     110             :   if (DEBUGMEM && pari_mainstack->memused != DISABLE_MEMUSED) {
     111             :     long d = (long)pari_mainstack->memused - (long)z;
     112             :     if (labs(d) > 4*MEMSTEP)
     113             :     {
     114             :       pari_mainstack->memused = (pari_sp)z;
     115             :       err_printf("...%4.0lf Mbytes used\n",
     116             :                 (pari_mainstack->top-pari_mainstack->memused)/1048576.);
     117             :     }
     118             :   }
     119             : #endif
     120 35326150357 :   return z;
     121             : }
     122             : 
     123             : INLINE char *
     124     3684552 : stack_malloc(size_t N)
     125             : {
     126     3684552 :   long n = nchar2nlong(N);
     127     3684552 :   return (char*)new_chunk(n);
     128             : }
     129             : 
     130             : INLINE char *
     131      551054 : stack_malloc_align(size_t N, long k)
     132             : {
     133      551054 :   ulong d = ((ulong)avma) % k;
     134      551054 :   if (d) (void)new_chunk(d/sizeof(long));
     135      551054 :   return (char*) new_chunk(nchar2nlong(N));
     136             : }
     137             : 
     138             : INLINE char *
     139       23312 : stack_calloc(size_t N)
     140             : {
     141       23312 :   char *p = stack_malloc(N);
     142       23312 :   memset(p, 0, N); return p;
     143             : }
     144             : 
     145             : /* cgetg(lg(x), typ(x)), set *lx. Implicit unsetisclone() */
     146             : INLINE GEN
     147   341939428 : cgetg_copy(GEN x, long *plx) {
     148             :   GEN y;
     149   341939428 :   *plx = lg(x); y = new_chunk((size_t)*plx);
     150   341939050 :   y[0] = x[0] & (TYPBITS|LGBITS); return y;
     151             : }
     152             : INLINE GEN
     153       95411 : cgetg_block(long x, long y)
     154             : {
     155       95411 :   GEN z = newblock((size_t)x);
     156       95602 :   z[0] = evaltyp(y) | evallg(x);
     157       95624 :   return z;
     158             : }
     159             : INLINE GEN
     160  4085863003 : cgetg(long x, long y)
     161             : {
     162  4085863003 :   GEN z = new_chunk((size_t)x);
     163  4084106344 :   z[0] = evaltyp(y) | evallg(x);
     164  4081365718 :   return z;
     165             : }
     166             : INLINE GEN
     167 12446399239 : cgeti(long x)
     168             : {
     169 12446399239 :   GEN z = new_chunk((size_t)x);
     170 12433420179 :   z[0] = evaltyp(t_INT) | evallg(x);
     171 12423807064 :   return z;
     172             : }
     173             : INLINE GEN
     174  7464632994 : cgetipos(long x)
     175             : {
     176  7464632994 :   GEN z = cgeti(x);
     177  7463482533 :   z[1] = evalsigne(1) | evallgefint(x);
     178  7463482533 :   return z;
     179             : }
     180             : INLINE GEN
     181   211045390 : cgetineg(long x)
     182             : {
     183   211045390 :   GEN z = cgeti(x);
     184   211045443 :   z[1] = evalsigne(-1) | evallgefint(x);
     185   211045443 :   return z;
     186             : }
     187             : INLINE GEN
     188        2571 : cgetr_block(long x)
     189             : {
     190        2571 :   GEN z = newblock((size_t)x);
     191        2571 :   z[0] = evaltyp(t_REAL) | evallg(x);
     192        2571 :   return z;
     193             : }
     194             : INLINE GEN
     195  5185031380 : cgetr(long x)
     196             : {
     197  5185031380 :   GEN z = new_chunk((size_t)x);
     198  5185031380 :   z[0] = evaltyp(t_REAL) | evallg(x);
     199  5185031380 :   return z;
     200             : }
     201             : 
     202             : /*******************************************************************/
     203             : /*                                                                 */
     204             : /*                     COPY, NEGATION, ABSOLUTE VALUE              */
     205             : /*                                                                 */
     206             : /*******************************************************************/
     207             : /* cannot do memcpy because sometimes x and y overlap */
     208             : INLINE GEN
     209  1138421789 : leafcopy(GEN x)
     210             : {
     211  1138421789 :   register long lx = lg(x);
     212  1138421789 :   GEN y = new_chunk(lx); /* can't use cgetg_copy, in case x,y overlap */
     213  1138851449 :   while (--lx > 0) y[lx] = x[lx];
     214  1138851449 :   y[0] = x[0] & (TYPBITS|LGBITS); return y;
     215             : }
     216             : INLINE GEN
     217  4284535246 : icopy(GEN x)
     218             : {
     219  4284535246 :   long i = lgefint(x), lx = i;
     220  4284535246 :   GEN y = new_chunk(lx); /* can't use cgeti, in case x,y overlap */
     221  4287373882 :   while (--i > 0) y[i] = x[i];
     222  4287373882 :   y[0] = evaltyp(t_INT) | evallg(lx);
     223  4286097514 :   return y;
     224             : }
     225             : INLINE GEN
     226    27732268 : icopyspec(GEN x, long nx)
     227             : {
     228    27732268 :   long i = nx+2, lx = i;
     229    27732268 :   GEN y = new_chunk(lx); /* can't use cgeti, in case x,y overlap */
     230    27732235 :   x -= 2; while (--i >= 2) y[i] = x[i];
     231    27732235 :   y[1] = evalsigne(1) | evallgefint(lx);
     232    27732235 :   y[0] = evaltyp(t_INT) | evallg(lx);
     233    27732226 :   return y;
     234             : }
     235   180913990 : INLINE GEN rcopy(GEN x) { return leafcopy(x); }
     236     1505054 : INLINE GEN mpcopy(GEN x) { return leafcopy(x); }
     237             : 
     238             : INLINE GEN
     239   237386048 : mpabs(GEN x) { GEN y = leafcopy(x); setabssign(y); return y; }
     240             : INLINE GEN
     241       11508 : mpabs_shallow(GEN x) { return signe(x) < 0? mpabs(x): x; }
     242   198943894 : INLINE GEN absi(GEN x) { return mpabs(x); }
     243      103131 : INLINE GEN absi_shallow(GEN x) { return signe(x) < 0? negi(x): x; }
     244    24190815 : INLINE GEN absr(GEN x) { return mpabs(x); }
     245             : 
     246             : INLINE GEN
     247   322954689 : mpneg(GEN x) { GEN y = leafcopy(x); togglesign(y); return y; }
     248   226143884 : INLINE GEN negi(GEN x) { return mpneg(x); }
     249     1198652 : INLINE GEN negr(GEN x) { return mpneg(x); }
     250             : 
     251             : /* negate in place */
     252             : INLINE void
     253   770946638 : togglesign(GEN x) { if (x[1] & SIGNBITS) { x[1] ^= HIGHBIT; } }
     254             : INLINE void
     255   261288173 : setabssign(GEN x) { x[1] &= ~HIGHBIT; }
     256             : /* negate in place, except universal constants */
     257             : INLINE void
     258    54207942 : togglesign_safe(GEN *px)
     259             : {
     260    54207942 :   switch(*px - gen_1) /* gen_1, gen_2, gen_m1, gen_m2 */
     261             :   {
     262      838166 :     case 0: *px = gen_m1; break;
     263           0 :     case 3: *px = gen_m2;  break;
     264      195298 :     case 6: *px = gen_1; break;
     265           0 :     case 9: *px = gen_2;  break;
     266    53174478 :     default: togglesign(*px);
     267             :   }
     268    54207942 : }
     269             : /* setsigne(y, signe(x)) */
     270             : INLINE void
     271           0 : affectsign(GEN x, GEN y)
     272             : {
     273           0 :   y[1] = (x[1] & SIGNBITS) | (y[1] & ~SIGNBITS);
     274           0 : }
     275             : /* copies sign in place, except for universal constants */
     276             : INLINE void
     277     5217965 : affectsign_safe(GEN x, GEN *py)
     278             : {
     279     5217965 :   if (((*py)[1] ^ x[1]) & HIGHBIT) togglesign_safe(py);
     280     5217965 : }
     281             : /*******************************************************************/
     282             : /*                                                                 */
     283             : /*                     GEN -> LONG, LONG -> GEN                    */
     284             : /*                                                                 */
     285             : /*******************************************************************/
     286             : /* assume x != 0, return -x as a t_INT */
     287             : INLINE GEN
     288   210883603 : utoineg(ulong x) { GEN y = cgetineg(3); y[2] = x; return y; }
     289             : /* assume x != 0, return utoi(x) */
     290             : INLINE GEN
     291  6393491972 : utoipos(ulong x) { GEN y = cgetipos(3); y[2] = x; return y; }
     292             : INLINE GEN
     293  5300071409 : utoi(ulong x) { return x? utoipos(x): gen_0; }
     294             : INLINE GEN
     295   306100948 : stoi(long x)
     296             : {
     297   306100948 :   if (!x) return gen_0;
     298   191722142 :   return x > 0? utoipos((ulong)x): utoineg((ulong)-x);
     299             : }
     300             : 
     301             : /* x 2^BIL + y */
     302             : INLINE GEN
     303  3939287433 : uutoi(ulong x, ulong y)
     304             : {
     305             :   GEN z;
     306  3939287433 :   if (!x) return utoi(y);
     307   429101841 :   z = cgetipos(4);
     308   428637282 :   *int_W_lg(z, 1, 4) = x;
     309   428637282 :   *int_W_lg(z, 0, 4) = y; return z;
     310             : }
     311             : /* - (x 2^BIL + y) */
     312             : INLINE GEN
     313       60187 : uutoineg(ulong x, ulong y)
     314             : {
     315             :   GEN z;
     316       60187 :   if (!x) return y? utoineg(y): gen_0;
     317           0 :   z = cgetineg(4);
     318           0 :   *int_W_lg(z, 1, 4) = x;
     319           0 :   *int_W_lg(z, 0, 4) = y; return z;
     320             : }
     321             : 
     322             : INLINE long
     323   199109006 : itos(GEN x)
     324             : {
     325   199109006 :   long s = signe(x);
     326             :   long u;
     327             : 
     328   199109006 :   if (!s) return 0;
     329   188060377 :   u = x[2];
     330   188060377 :   if (lgefint(x) > 3 || u < 0)
     331          19 :     pari_err_OVERFLOW("t_INT-->long assignment");
     332   188060379 :   return (s>0) ? u : -u;
     333             : }
     334             : /* as itos, but return 0 if too large. Cf is_bigint */
     335             : INLINE long
     336     9115655 : itos_or_0(GEN x) {
     337             :   long n;
     338     9115655 :   if (lgefint(x) != 3 || (n = x[2]) & HIGHBIT) return 0;
     339     9096388 :   return signe(x) > 0? n: -n;
     340             : }
     341             : INLINE ulong
     342    25446263 : itou(GEN x)
     343             : {
     344    25446263 :   switch(lgefint(x)) {
     345     5399483 :     case 2: return 0;
     346    20046780 :     case 3: return x[2];
     347             :     default:
     348           0 :       pari_err_OVERFLOW("t_INT-->ulong assignment");
     349             :       return 0; /* LCOV_EXCL_LINE */
     350             :   }
     351             : }
     352             : 
     353             : /* as itou, but return 0 if too large. Cf is_bigint */
     354             : INLINE ulong
     355     4162688 : itou_or_0(GEN x) {
     356     4162688 :   if (lgefint(x) != 3) return 0;
     357     4152664 :   return (ulong)x[2];
     358             : }
     359             : 
     360             : INLINE GEN
     361    78108771 : real_0_bit(long bitprec) { GEN x=cgetr(2); x[1]=evalexpo(bitprec); return x; }
     362             : INLINE GEN
     363    25351306 : real_0(long prec) { return real_0_bit(-prec2nbits(prec)); }
     364             : INLINE GEN
     365      487063 : real_1_bit(long bit) { return real_1(nbits2prec(bit)); }
     366             : INLINE GEN
     367    27327790 : real_1(long prec) {
     368    27327790 :   GEN x = cgetr(prec);
     369             :   long i;
     370    27327790 :   x[1] = evalsigne(1) | _evalexpo(0);
     371    27327790 :   x[2] = (long)HIGHBIT; for (i=3; i<prec; i++) x[i] = 0;
     372    27327790 :   return x;
     373             : }
     374             : INLINE GEN
     375          91 : real_m1(long prec) {
     376          91 :   GEN x = cgetr(prec);
     377             :   long i;
     378          91 :   x[1] = evalsigne(-1) | _evalexpo(0);
     379          91 :   x[2] = (long)HIGHBIT; for (i=3; i<prec; i++) x[i] = 0;
     380          91 :   return x;
     381             : }
     382             : 
     383             : /* 2.^n */
     384             : INLINE GEN
     385       56041 : real2n(long n, long prec) { GEN z = real_1(prec); setexpo(z, n); return z; }
     386             : INLINE GEN
     387           0 : real_m2n(long n, long prec) { GEN z = real_m1(prec); setexpo(z, n); return z; }
     388             : INLINE GEN
     389   143115372 : stor(long s, long prec) { GEN z = cgetr(prec); affsr(s,z); return z; }
     390             : INLINE GEN
     391     5417984 : utor(ulong s, long prec){ GEN z = cgetr(prec); affur(s,z); return z; }
     392             : INLINE GEN
     393   140355667 : itor(GEN x, long prec) { GEN z = cgetr(prec); affir(x,z); return z; }
     394             : INLINE GEN
     395    72922202 : rtor(GEN x, long prec) { GEN z = cgetr(prec); affrr(x,z); return z; }
     396             : 
     397     5248760 : INLINE ulong int_bit(GEN x, long n)
     398             : {
     399     5248760 :   long r, q = dvmdsBIL(n, &r);
     400     5250248 :   return q < lgefint(x)-2?((ulong)*int_W(x,q) >> r) & 1UL:0;
     401             : }
     402             : 
     403             : /*******************************************************************/
     404             : /*                                                                 */
     405             : /*                           COMPARISON                            */
     406             : /*                                                                 */
     407             : /*******************************************************************/
     408             : INLINE int
     409      696466 : cmpir(GEN x, GEN y)
     410             : {
     411             :   pari_sp av;
     412             :   GEN z;
     413             : 
     414      696466 :   if (!signe(x)) return -signe(y);
     415      300564 :   if (!signe(y))
     416             :   {
     417        4594 :     if (expo(y) >= expi(x)) return 0;
     418        4566 :     return signe(x);
     419             :   }
     420      295970 :   av=avma; z = itor(x, realprec(y)); avma=av;
     421      295970 :   return cmprr(z,y); /* cmprr does no memory adjustment */
     422             : }
     423             : INLINE int
     424      437156 : cmpri(GEN x, GEN y) { return -cmpir(y,x); }
     425             : INLINE int
     426       61200 : cmpsr(long x, GEN y)
     427             : {
     428             :   pari_sp av;
     429             :   GEN z;
     430             : 
     431       61200 :   if (!x) return -signe(y);
     432       61200 :   av=avma; z = stor(x, LOWDEFAULTPREC); avma=av;
     433       61200 :   return cmprr(z,y);
     434             : }
     435             : INLINE int
     436       38668 : cmprs(GEN x, long y) { return -cmpsr(y,x); }
     437             : /* compare x and |y| */
     438             : INLINE int
     439    22583019 : abscmpui(ulong x, GEN y)
     440             : {
     441    22583019 :   long l = lgefint(y);
     442             :   ulong p;
     443             : 
     444    22583019 :   if (!x) return (l > 2)? -1: 0;
     445    22583005 :   if (l == 2) return 1;
     446    22392112 :   if (l > 3) return -1;
     447    22381972 :   p = y[2]; if (p == x) return 0;
     448    21992395 :   return p < x ? 1 : -1;
     449             : }
     450             : INLINE int
     451    22530475 : abscmpiu(GEN x, ulong y) { return -abscmpui(y,x); }
     452             : INLINE int
     453     4126562 : cmpsi(long x, GEN y)
     454             : {
     455             :   ulong p;
     456             : 
     457     4126562 :   if (!x) return -signe(y);
     458             : 
     459     4125946 :   if (x > 0)
     460             :   {
     461     4125610 :     if (signe(y)<=0) return 1;
     462     4125414 :     if (lgefint(y)>3) return -1;
     463     4112283 :     p = y[2]; if (p == (ulong)x) return 0;
     464     4078459 :     return p < (ulong)x ? 1 : -1;
     465             :   }
     466             : 
     467         336 :   if (signe(y)>=0) return -1;
     468          63 :   if (lgefint(y)>3) return 1;
     469          63 :   p = y[2]; if (p == (ulong)-x) return 0;
     470          14 :   return p < (ulong)(-x) ? -1 : 1;
     471             : }
     472             : INLINE int
     473     4124868 : cmpis(GEN x, long y) { return -cmpsi(y,x); }
     474             : INLINE int
     475      535615 : mpcmp(GEN x, GEN y)
     476             : {
     477      535615 :   if (typ(x)==t_INT)
     478        6337 :     return (typ(y)==t_INT) ? cmpii(x,y) : cmpir(x,y);
     479      529278 :   return (typ(y)==t_INT) ? -cmpir(y,x) : cmprr(x,y);
     480             : }
     481             : 
     482             : /* x == y ? */
     483             : INLINE int
     484      102707 : equalsi(long x, GEN y)
     485             : {
     486      102707 :   if (!x) return !signe(y);
     487      102707 :   if (x > 0)
     488             :   {
     489       33260 :     if (signe(y) <= 0 || lgefint(y) != 3) return 0;
     490       30299 :     return ((ulong)y[2] == (ulong)x);
     491             :   }
     492       69447 :   if (signe(y) >= 0 || lgefint(y) != 3) return 0;
     493       69447 :   return ((ulong)y[2] == (ulong)-x);
     494             : }
     495             : /* x == |y| ? */
     496             : INLINE int
     497    11598613 : absequalui(ulong x, GEN y)
     498             : {
     499    11598613 :   if (!x) return !signe(y);
     500    11598613 :   return (lgefint(y) == 3 && (ulong)y[2] == x);
     501             : }
     502             : INLINE int
     503    11495417 : absequaliu(GEN x, ulong y) { return absequalui(y,x); }
     504             : INLINE int
     505      102588 : equalis(GEN x, long y) { return equalsi(y,x); }
     506             : 
     507             : /* assume x != 0, is |x| == 2^n ? */
     508             : INLINE int
     509      760649 : absrnz_equal2n(GEN x) {
     510      760649 :   if ((ulong)x[2]==HIGHBIT)
     511             :   {
     512       41570 :     long i, lx = lg(x);
     513      623645 :     for (i = 3; i < lx; i++)
     514      587831 :       if (x[i]) return 0;
     515       35814 :     return 1;
     516             :   }
     517      719079 :   return 0;
     518             : }
     519             : /* assume x != 0, is |x| == 1 ? */
     520             : INLINE int
     521     1280125 : absrnz_equal1(GEN x) { return !expo(x) && absrnz_equal2n(x); }
     522             : 
     523             : INLINE long
     524  2507276842 : maxss(long x, long y) { return x>y?x:y; }
     525             : INLINE long
     526   827956993 : minss(long x, long y) { return x<y?x:y; }
     527             : INLINE long
     528     3288789 : minuu(ulong x, ulong y) { return x<y?x:y; }
     529             : INLINE long
     530     6389144 : maxuu(ulong x, ulong y) { return x>y?x:y; }
     531             : INLINE double
     532      963033 : maxdd(double x, double y) { return x>y?x:y; }
     533             : INLINE double
     534      377280 : mindd(double x, double y) { return x<y?x:y; }
     535             : 
     536             : /*******************************************************************/
     537             : /*                                                                 */
     538             : /*                             ADD / SUB                           */
     539             : /*                                                                 */
     540             : /*******************************************************************/
     541             : INLINE GEN
     542       25060 : subuu(ulong x, ulong y)
     543             : {
     544             :   ulong z;
     545             :   LOCAL_OVERFLOW;
     546       25060 :   z = subll(x, y);
     547       25060 :   return overflow? utoineg(-z): utoi(z);
     548             : }
     549             : INLINE GEN
     550  1520440294 : adduu(ulong x, ulong y) { ulong t = x+y; return uutoi((t < x), t); }
     551             : 
     552             : INLINE GEN
     553       25060 : addss(long x, long y)
     554             : {
     555       25060 :   if (!x) return stoi(y);
     556       25060 :   if (!y) return stoi(x);
     557       25060 :   if (x > 0) return y > 0? adduu(x,y): subuu(x, -y);
     558             : 
     559       25060 :   if (y > 0) return subuu(y, -x);
     560             :   else { /* - adduu(-x, -y) */
     561           0 :     ulong t = (-x)+(-y); return uutoineg((t < (ulong)(-x)), t);
     562             :   }
     563             : }
     564       25060 : INLINE GEN subss(long x, long y) { return addss(-y,x); }
     565             : 
     566             : INLINE GEN
     567  3653349792 : subii(GEN x, GEN y)
     568             : {
     569  3653349792 :   if (x==y) return gen_0; /* frequent with x = y = gen_0 */
     570  2654733897 :   return addii_sign(x, signe(x), y, -signe(y));
     571             : }
     572             : INLINE GEN
     573  4552060010 : addii(GEN x, GEN y) { return addii_sign(x, signe(x), y, signe(y)); }
     574             : INLINE GEN
     575   694111621 : addrr(GEN x, GEN y) { return addrr_sign(x, signe(x), y, signe(y)); }
     576             : INLINE GEN
     577   818252590 : subrr(GEN x, GEN y) { return addrr_sign(x, signe(x), y, -signe(y)); }
     578             : INLINE GEN
     579    97925097 : addir(GEN x, GEN y) { return addir_sign(x, signe(x), y, signe(y)); }
     580             : INLINE GEN
     581    75515353 : subir(GEN x, GEN y) { return addir_sign(x, signe(x), y, -signe(y)); }
     582             : INLINE GEN
     583     1359252 : subri(GEN x, GEN y) { return addir_sign(y, -signe(y), x, signe(x)); }
     584             : INLINE GEN
     585    49557537 : addsi(long x, GEN y) { return addsi_sign(x, y, signe(y)); }
     586             : INLINE GEN
     587    30199745 : addui(ulong x, GEN y) { return addui_sign(x, y, signe(y)); }
     588             : INLINE GEN
     589     3390292 : subsi(long x, GEN y) { return addsi_sign(x, y, -signe(y)); }
     590             : INLINE GEN
     591    22693120 : subui(ulong x, GEN y) { return addui_sign(x, y, -signe(y)); }
     592             : 
     593             : /*******************************************************************/
     594             : /*                                                                 */
     595             : /*                           MOD, REM, DIV                         */
     596             : /*                                                                 */
     597             : /*******************************************************************/
     598    76725605 : INLINE ulong mod2BIL(GEN x) { return *int_LSW(x); }
     599           0 : INLINE long mod64(GEN x) { return mod2BIL(x) & 63; }
     600         252 : INLINE long mod32(GEN x) { return mod2BIL(x) & 31; }
     601       57603 : INLINE long mod16(GEN x) { return mod2BIL(x) & 15; }
     602    39318649 : INLINE long mod8(GEN x)  { return mod2BIL(x) & 7; }
     603    10879345 : INLINE long mod4(GEN x)  { return mod2BIL(x) & 3; }
     604    15783762 : INLINE long mod2(GEN x)  { return mod2BIL(x) & 1; }
     605             : INLINE int
     606    22240630 : mpodd(GEN x) { return signe(x) && mod2(x); }
     607             : /* x mod 2^n, n < BITS_IN_LONG */
     608             : INLINE ulong
     609     5443893 : umodi2n(GEN x, long n)
     610             : {
     611     5443893 :   long s = signe(x);
     612     5443893 :   const ulong _2n = 1UL << n;
     613             :   ulong m;
     614     5443893 :   if (!s) return 0;
     615     5434058 :   m = *int_LSW(x) & (_2n - 1);
     616     5434058 :   if (s < 0 && m) m = _2n - m;
     617     5434058 :   return m;
     618             : }
     619           0 : INLINE ulong Mod64(GEN x){ return umodi2n(x,6); }
     620      166579 : INLINE ulong Mod32(GEN x){ return umodi2n(x,5); }
     621      167461 : INLINE ulong Mod16(GEN x){ return umodi2n(x,4); }
     622     1307992 : INLINE ulong Mod8(GEN x) { return umodi2n(x,3); }
     623     1323406 : INLINE ulong Mod4(GEN x) { return umodi2n(x,2); }
     624     2478385 : INLINE ulong Mod2(GEN x) { return umodi2n(x,1); }
     625             : 
     626             : INLINE GEN
     627    34818189 : truedivii(GEN a,GEN b) { return truedvmdii(a,b,NULL); }
     628             : INLINE GEN
     629         659 : truedivis(GEN a, long b) { return truedvmdis(a,b,NULL); }
     630             : INLINE GEN
     631     6084042 : truedivsi(long a, GEN b) { return truedvmdsi(a,b,NULL); }
     632             : 
     633             : INLINE GEN
     634     7291345 : divii(GEN a, GEN b) { return dvmdii(a,b,NULL); }
     635             : INLINE GEN
     636  1161747874 : remii(GEN a, GEN b) { return dvmdii(a,b,ONLY_REM); }
     637             : 
     638             : INLINE GEN
     639           0 : divss(long x, long y) { return stoi(x / y); }
     640             : INLINE GEN
     641      599053 : modss(long x, long y) { return stoi(smodss(x, y)); }
     642             : INLINE GEN
     643           0 : remss(long x, long y) { return stoi(x % y); }
     644             : INLINE long
     645     6306183 : smodss(long x, long y)
     646             : {
     647     6306183 :   long r = x%y;
     648     6306183 :   return (r >= 0)? r: labs(y) + r;
     649             : }
     650             : INLINE ulong
     651    11844359 : umodsu(long x, ulong y)
     652             : {
     653    11844359 :   return x>=0 ? x%y: Fl_neg((-x)%y, y);
     654             : }
     655             : 
     656             : INLINE long
     657           0 : sdivss_rem(long x, long y, long *r)
     658             : {
     659             :   long q;
     660             :   LOCAL_HIREMAINDER;
     661           0 :   if (!y) pari_err_INV("sdivss_rem",gen_0);
     662           0 :   hiremainder = 0; q = divll((ulong)labs(x),(ulong)labs(y));
     663           0 :   if (x < 0) { hiremainder = -((long)hiremainder); q = -q; }
     664           0 :   if (y < 0) q = -q;
     665           0 :   *r = hiremainder; return q;
     666             : }
     667             : INLINE GEN
     668           0 : divss_rem(long x, long y, long *r) { return stoi(sdivss_rem(x,y,r)); }
     669             : INLINE ulong
     670         483 : udivuu_rem(ulong x, ulong y, ulong *r)
     671             : {
     672         483 :   if (!y) pari_err_INV("udivuu_rem",gen_0);
     673         483 :   *r = x % y; return x / y;
     674             : }
     675             : INLINE ulong
     676       77672 : ceildivuu(ulong a, ulong b)
     677             : {
     678       77672 :   ulong c = a/b;
     679       77672 :   return (a%b)? c+1: c;
     680             : }
     681             : 
     682             : 
     683             : INLINE ulong
     684       10780 : udivui_rem(ulong x, GEN y, ulong *r)
     685             : {
     686       10780 :   long q, s = signe(y);
     687             :   LOCAL_HIREMAINDER;
     688             : 
     689       10780 :   if (!s) pari_err_INV("udivui_rem",gen_0);
     690       10780 :   if (!x || lgefint(y)>3) { *r = x; return 0; }
     691       10500 :   hiremainder=0; q = (long)divll(x, (ulong)y[2]);
     692       10500 :   if (s < 0) q = -q;
     693       10500 :   *r = hiremainder; return q;
     694             : }
     695             : 
     696             : /* assume d != 0 and |n| / d can be represented as an ulong.
     697             :  * Return |n|/d, set *r = |n| % d */
     698             : INLINE ulong
     699    12480475 : udiviu_rem(GEN n, ulong d, ulong *r)
     700             : {
     701    12480475 :   switch(lgefint(n))
     702             :   {
     703           0 :     case 2: *r = 0; return 0;
     704             :     case 3:
     705             :     {
     706    12480475 :       ulong nn = n[2];
     707    12480475 :       *r = nn % d; return nn / d;
     708             :     }
     709             :     default: /* 4 */
     710             :     {
     711             :       ulong n1, n0, q;
     712             :       LOCAL_HIREMAINDER;
     713           0 :       n0 = *int_W(n,0);
     714           0 :       n1 = *int_W(n,1);
     715           0 :       hiremainder = n1;
     716           0 :       q = divll(n0, d);
     717           0 :       *r = hiremainder; return q;
     718             :     }
     719             :   }
     720             : }
     721             : 
     722             : INLINE long
     723    13548508 : sdivsi_rem(long x, GEN y, long *r)
     724             : {
     725    13548508 :   long q, s = signe(y);
     726             :   LOCAL_HIREMAINDER;
     727             : 
     728    13548508 :   if (!s) pari_err_INV("sdivsi_rem",gen_0);
     729    13548507 :   if (!x || lgefint(y)>3 || ((long)y[2]) < 0) { *r = x; return 0; }
     730    11873292 :   hiremainder=0; q = (long)divll(labs(x), (ulong)y[2]);
     731    11873292 :   if (x < 0) { hiremainder = -((long)hiremainder); q = -q; }
     732    11873292 :   if (s < 0) q = -q;
     733    11873292 :   *r = hiremainder; return q;
     734             : }
     735             : INLINE GEN
     736           0 : divsi_rem(long s, GEN y, long *r) { return stoi(sdivsi_rem(s,y,r)); }
     737             : 
     738             : INLINE long
     739          98 : sdivsi(long x, GEN y)
     740             : {
     741          98 :   long q, s = signe(y);
     742             : 
     743          98 :   if (!s) pari_err_INV("sdivsi",gen_0);
     744          98 :   if (!x || lgefint(y)>3 || ((long)y[2]) < 0) return 0;
     745          70 :   q = labs(x) / y[2];
     746          70 :   if (x < 0) q = -q;
     747          70 :   if (s < 0) q = -q;
     748          70 :   return q;
     749             : }
     750             : 
     751             : INLINE GEN
     752           0 : dvmdss(long x, long y, GEN *z)
     753             : {
     754             :   long r;
     755           0 :   GEN q = divss_rem(x,y, &r);
     756           0 :   *z = stoi(r); return q;
     757             : }
     758             : INLINE long
     759  2496765520 : dvmdsBIL(long n, long *r) { *r = remsBIL(n); return divsBIL(n); }
     760             : INLINE ulong
     761   159152622 : dvmduBIL(ulong n, ulong *r) { *r = remsBIL(n); return divsBIL(n); }
     762             : INLINE GEN
     763           0 : dvmdsi(long x, GEN y, GEN *z)
     764             : {
     765             :   long r;
     766           0 :   GEN q = divsi_rem(x,y, &r);
     767           0 :   *z = stoi(r); return q;
     768             : }
     769             : INLINE GEN
     770           0 : dvmdis(GEN x, long y, GEN *z)
     771             : {
     772             :   long r;
     773           0 :   GEN q = divis_rem(x,y, &r);
     774           0 :   *z = stoi(r); return q;
     775             : }
     776             : 
     777             : INLINE long
     778      750097 : smodis(GEN x, long y)
     779             : {
     780      750097 :   pari_sp av = avma;
     781             :   long r;
     782      750097 :   (void)divis_rem(x,y, &r); avma = av; return (r >= 0) ? r: labs(y) + r;
     783             : }
     784             : INLINE GEN
     785       36892 : modis(GEN x, long y) { return stoi(smodis(x,y)); }
     786             : INLINE GEN
     787     7464220 : modsi(long x, GEN y) {
     788             :   long r;
     789     7464220 :   (void)sdivsi_rem(x, y, &r);
     790     7464220 :   return (r >= 0)? stoi(r): addsi_sign(r, y, 1);
     791             : }
     792             : 
     793             : INLINE ulong
     794     1194820 : umodui(ulong x, GEN y)
     795             : {
     796     1194820 :   if (!signe(y)) pari_err_INV("umodui",gen_0);
     797     1194820 :   if (!x || lgefint(y) > 3) return x;
     798      204147 :   return x % (ulong)y[2];
     799             : }
     800             : 
     801             : INLINE GEN
     802           0 : remsi(long x, GEN y)
     803           0 : { long r; (void)sdivsi_rem(x,y, &r); return stoi(r); }
     804             : INLINE GEN
     805           0 : remis(GEN x, long y)
     806             : {
     807           0 :   pari_sp av = avma;
     808             :   long r;
     809           0 :   (void)divis_rem(x,y, &r); avma = av; return stoi(r);
     810             : }
     811             : 
     812             : INLINE GEN
     813           0 : rdivis(GEN x, long y, long prec)
     814             : {
     815           0 :   GEN z = cgetr(prec);
     816           0 :   pari_sp av = avma;
     817           0 :   affrr(divrs(itor(x,prec), y),z);
     818           0 :   avma = av; return z;
     819             : }
     820             : INLINE GEN
     821           0 : rdivsi(long x, GEN y, long prec)
     822             : {
     823           0 :   GEN z = cgetr(prec);
     824           0 :   pari_sp av = avma;
     825           0 :   affrr(divsr(x, itor(y,prec)), z);
     826           0 :   avma = av; return z;
     827             : }
     828             : INLINE GEN
     829      839647 : rdivss(long x, long y, long prec)
     830             : {
     831      839647 :   GEN z = cgetr(prec);
     832      839647 :   pari_sp av = avma;
     833      839647 :   affrr(divrs(stor(x, prec), y), z);
     834      839647 :   avma = av; return z;
     835             : }
     836             : 
     837             : INLINE void
     838          14 : rdiviiz(GEN x, GEN y, GEN z)
     839             : {
     840          14 :   pari_sp av = avma;
     841          14 :   long prec = realprec(z);
     842          14 :   affir(x, z);
     843          14 :   if (!is_bigint(y)) {
     844           0 :     affrr(divrs(z, y[2]), z);
     845           0 :     if (signe(y) < 0) togglesign(z);
     846             :   }
     847             :   else
     848          14 :     affrr(divrr(z, itor(y,prec)), z);
     849          14 :   avma = av;
     850          14 : }
     851             : INLINE GEN
     852    11498669 : rdivii(GEN x, GEN y, long prec)
     853             : {
     854    11498669 :   GEN z = cgetr(prec);
     855    11498669 :   pari_sp av = avma;
     856    11498669 :   affir(x, z);
     857    11498669 :   if (lg(y) == 3) {
     858     7409869 :     affrr(divru(z, y[2]), z);
     859     7409869 :     if (signe(y) < 0) togglesign(z);
     860             :   }
     861             :   else
     862     4088800 :     affrr(divrr(z, itor(y,prec)), z);
     863    11498669 :   avma = av; return z;
     864             : }
     865             : INLINE GEN
     866    11445069 : fractor(GEN x, long prec) { return rdivii(gel(x,1), gel(x,2), prec); }
     867             : 
     868             : INLINE int
     869     1340087 : dvdii(GEN x, GEN y)
     870             : {
     871     1340087 :   pari_sp av=avma;
     872     1340087 :   GEN r = remii(x,y);
     873     1340087 :   avma = av; return r == gen_0;
     874             : }
     875             : INLINE int
     876          77 : dvdsi(long x, GEN y)
     877             : {
     878          77 :   if (!signe(y)) return x == 0;
     879          77 :   if (lgefint(y) != 3) return 0;
     880          77 :   return x % y[2] == 0;
     881             : }
     882             : INLINE int
     883         147 : dvdui(ulong x, GEN y)
     884             : {
     885         147 :   if (!signe(y)) return x == 0;
     886         147 :   if (lgefint(y) != 3) return 0;
     887         147 :   return x % y[2] == 0;
     888             : }
     889             : INLINE int
     890       31512 : dvdis(GEN x, long y)
     891       31512 : { return y? smodis(x, y) == 0: signe(x) == 0; }
     892             : INLINE int
     893        1561 : dvdiu(GEN x, ulong y)
     894        1561 : { return y? umodiu(x, y) == 0: signe(x) == 0; }
     895             : 
     896             : INLINE int
     897           0 : dvdisz(GEN x, long y, GEN z)
     898             : {
     899           0 :   const pari_sp av = avma;
     900             :   long r;
     901           0 :   GEN p1 = divis_rem(x,y, &r);
     902           0 :   avma = av; if (r) return 0;
     903           0 :   affii(p1,z); return 1;
     904             : }
     905             : INLINE int
     906           0 : dvdiuz(GEN x, ulong y, GEN z)
     907             : {
     908           0 :   const pari_sp av = avma;
     909             :   ulong r;
     910           0 :   GEN p1 = diviu_rem(x,y, &r);
     911           0 :   avma = av; if (r) return 0;
     912           0 :   affii(p1,z); return 1;
     913             : }
     914             : INLINE int
     915        3705 : dvdiiz(GEN x, GEN y, GEN z)
     916             : {
     917        3705 :   const pari_sp av=avma;
     918             :   GEN p2;
     919        3705 :   const GEN p1=dvmdii(x,y,&p2);
     920             : 
     921        3705 :   if (signe(p2)) { avma=av; return 0; }
     922        3233 :   affii(p1,z); avma=av; return 1;
     923             : }
     924             : 
     925             : INLINE ulong
     926     6254883 : remlll_pre(ulong u2, ulong u1, ulong u0, ulong n, ulong ninv)
     927             : {
     928     6254883 :   u1 = remll_pre(u2, u1, n, ninv);
     929     6254883 :   return remll_pre(u1, u0, n, ninv);
     930             : }
     931             : 
     932             : INLINE ulong
     933  1749634741 : Fl_sqr_pre(ulong a, ulong p, ulong pi)
     934             : {
     935             :   register ulong x;
     936             :   LOCAL_HIREMAINDER;
     937  1749634741 :   x = mulll(a,a);
     938  1749634741 :   return remll_pre(hiremainder, x, p, pi);
     939             : }
     940             : 
     941             : INLINE ulong
     942  2272683739 : Fl_mul_pre(ulong a, ulong b, ulong p, ulong pi)
     943             : {
     944             :   register ulong x;
     945             :   LOCAL_HIREMAINDER;
     946  2272683739 :   x = mulll(a,b);
     947  2272683739 :   return remll_pre(hiremainder, x, p, pi);
     948             : }
     949             : 
     950             : INLINE ulong
     951  4233753987 : Fl_addmul_pre(ulong x0, ulong x1, ulong y0, ulong p, ulong pi)
     952             : {
     953             :   ulong l0, h0;
     954             :   LOCAL_HIREMAINDER;
     955  4233753987 :   hiremainder = y0;
     956  4233753987 :   l0 = addmul(x0, x1); h0 = hiremainder;
     957  4233753987 :   return remll_pre(h0, l0, p, pi);
     958             : }
     959             : 
     960             : INLINE ulong
     961    45225169 : Fl_addmulmul_pre(ulong x0, ulong y0, ulong x1, ulong y1, ulong p, ulong pi)
     962             : {
     963             :   ulong l0, l1, h0, h1;
     964             :   LOCAL_OVERFLOW;
     965             :   LOCAL_HIREMAINDER;
     966    45225169 :   l0 = mulll(x0, y0); h0 = hiremainder;
     967    45225169 :   l1 = mulll(x1, y1); h1 = hiremainder;
     968    45225169 :   l0 = addll(l0, l1); h0 = addllx(h0, h1);
     969    45225169 :   return overflow ? remlll_pre(1, h0, l0, p, pi): remll_pre(h0, l0, p, pi);
     970             : }
     971             : 
     972             : INLINE ulong
     973     2705052 : Fl_ellj_pre(ulong a4, ulong a6, ulong p, ulong pi)
     974             : {
     975             :   /* a43 = 4 a4^3 */
     976     2705052 :   ulong a43 = Fl_double(Fl_double(
     977             :               Fl_mul_pre(a4, Fl_sqr_pre(a4, p, pi), p, pi), p), p);
     978             :   /* a62 = 27 a6^2 */
     979     2705055 :   ulong a62 = Fl_mul_pre(Fl_sqr_pre(a6, p, pi), 27 % p, p, pi);
     980     2705058 :   ulong z1 = Fl_mul_pre(a43, 1728 % p, p, pi);
     981     2705059 :   ulong z2 = Fl_add(a43, a62, p);
     982     2705059 :   return Fl_div(z1, z2, p);
     983             : }
     984             : 
     985             : /*******************************************************************/
     986             : /*                                                                 */
     987             : /*                        MP (INT OR REAL)                         */
     988             : /*                                                                 */
     989             : /*******************************************************************/
     990             : INLINE GEN
     991          21 : mptrunc(GEN x) { return typ(x)==t_INT? icopy(x): truncr(x); }
     992             : INLINE GEN
     993           0 : mpfloor(GEN x) { return typ(x)==t_INT? icopy(x): floorr(x); }
     994             : INLINE GEN
     995           0 : mpceil(GEN x) { return typ(x)==t_INT? icopy(x): ceilr(x); }
     996             : INLINE GEN
     997      357796 : mpround(GEN x) { return typ(x) == t_INT? icopy(x): roundr(x); }
     998             : 
     999             : INLINE long
    1000      203267 : mpexpo(GEN x) { return typ(x) == t_INT? expi(x): expo(x); }
    1001             : 
    1002             : INLINE GEN
    1003    31561674 : mpadd(GEN x, GEN y)
    1004             : {
    1005    31561674 :   if (typ(x)==t_INT)
    1006     4722597 :     return (typ(y)==t_INT) ? addii(x,y) : addir(x,y);
    1007    26839077 :   return (typ(y)==t_INT) ? addir(y,x) : addrr(x,y);
    1008             : }
    1009             : INLINE GEN
    1010    14183289 : mpsub(GEN x, GEN y)
    1011             : {
    1012    14183289 :   if (typ(x)==t_INT)
    1013      235563 :     return (typ(y)==t_INT) ? subii(x,y) : subir(x,y);
    1014    13947726 :   return (typ(y)==t_INT) ? subri(x,y) : subrr(x,y);
    1015             : }
    1016             : INLINE GEN
    1017    46195731 : mpmul(GEN x, GEN y)
    1018             : {
    1019    46195731 :   if (typ(x)==t_INT)
    1020    13543028 :     return (typ(y)==t_INT) ? mulii(x,y) : mulir(x,y);
    1021    32652703 :   return (typ(y)==t_INT) ? mulir(y,x) : mulrr(x,y);
    1022             : }
    1023             : INLINE GEN
    1024     5769937 : mpsqr(GEN x) { return (typ(x)==t_INT) ? sqri(x) : sqrr(x); }
    1025             : INLINE GEN
    1026      143641 : mpdiv(GEN x, GEN y)
    1027             : {
    1028      143641 :   if (typ(x)==t_INT)
    1029       71913 :     return (typ(y)==t_INT) ? divii(x,y) : divir(x,y);
    1030       71728 :   return (typ(y)==t_INT) ? divri(x,y) : divrr(x,y);
    1031             : }
    1032             : 
    1033             : /*******************************************************************/
    1034             : /*                                                                 */
    1035             : /*                          Z/nZ, n ULONG                          */
    1036             : /*                                                                 */
    1037             : /*******************************************************************/
    1038             : INLINE ulong
    1039   640077121 : Fl_double(ulong a, ulong p)
    1040             : {
    1041   640077121 :   ulong res = a << 1;
    1042   640077121 :   return (res >= p || res < a) ? res - p : res;
    1043             : }
    1044             : INLINE ulong
    1045    88720678 : Fl_triple(ulong a, ulong p)
    1046             : {
    1047    88720678 :   ulong res = a << 1;
    1048    88720678 :   if (res >= p || res < a) res -= p;
    1049    88720678 :   res += a;
    1050    88720678 :   return (res >= p || res < a)? res - p: res;
    1051             : }
    1052             : INLINE ulong
    1053    12784303 : Fl_halve(ulong a, ulong p)
    1054             : {
    1055             :   ulong ap, ap2;
    1056    12784303 :   if ((a&1UL)==0) return a>>1;
    1057     6457800 :   ap = a + p; ap2 = ap>>1;
    1058     6457800 :   return ap>=a ? ap2: (ap2|HIGHBIT);
    1059             : }
    1060             : 
    1061             : INLINE ulong
    1062  2625144026 : Fl_add(ulong a, ulong b, ulong p)
    1063             : {
    1064  2625144026 :   ulong res = a + b;
    1065  2625144026 :   return (res >= p || res < a) ? res - p : res;
    1066             : }
    1067             : INLINE ulong
    1068    66697921 : Fl_neg(ulong x, ulong p) { return x ? p - x: 0; }
    1069             : 
    1070             : INLINE ulong
    1071  3877511948 : Fl_sub(ulong a, ulong b, ulong p)
    1072             : {
    1073  3877511948 :   ulong res = a - b;
    1074  3877511948 :   return (res > a) ? res + p: res;
    1075             : }
    1076             : 
    1077             : /* centerlift(u mod p) */
    1078             : INLINE long
    1079    19683503 : Fl_center(ulong u, ulong p, ulong ps2) { return (long) (u > ps2)? u - p: u; }
    1080             : 
    1081             : INLINE ulong
    1082  1046116627 : Fl_mul(ulong a, ulong b, ulong p)
    1083             : {
    1084             :   register ulong x;
    1085             :   LOCAL_HIREMAINDER;
    1086  1046116627 :   x = mulll(a,b);
    1087  1046116627 :   if (!hiremainder) return x % p;
    1088   168325406 :   (void)divll(x,p); return hiremainder;
    1089             : }
    1090             : INLINE ulong
    1091    36238677 : Fl_sqr(ulong a, ulong p)
    1092             : {
    1093             :   register ulong x;
    1094             :   LOCAL_HIREMAINDER;
    1095    36238677 :   x = mulll(a,a);
    1096    36238677 :   if (!hiremainder) return x % p;
    1097    14007882 :   (void)divll(x,p); return hiremainder;
    1098             : }
    1099             : INLINE ulong
    1100    10990392 : Fl_div(ulong a, ulong b, ulong p) { return Fl_mul(a, Fl_inv(b, p), p); }
    1101             : 
    1102             : /*******************************************************************/
    1103             : /*                                                                 */
    1104             : /*        DEFINED FROM EXISTING ONE EXPLOITING COMMUTATIVITY       */
    1105             : /*                                                                 */
    1106             : /*******************************************************************/
    1107             : INLINE GEN
    1108      181590 : addri(GEN x, GEN y) { return addir(y,x); }
    1109             : INLINE GEN
    1110    48754076 : addis(GEN x, long s) { return addsi(s,x); }
    1111             : INLINE GEN
    1112    29577094 : addiu(GEN x, ulong s) { return addui(s,x); }
    1113             : INLINE GEN
    1114     2592963 : addrs(GEN x, long s) { return addsr(s,x); }
    1115             : 
    1116             : INLINE GEN
    1117    19300982 : subiu(GEN x, long y) { GEN z = subui(y, x); togglesign(z); return z; }
    1118             : INLINE GEN
    1119        4268 : subis(GEN x, long y) { return addsi(-y,x); }
    1120             : INLINE GEN
    1121     2801000 : subrs(GEN x, long y) { return addsr(-y,x); }
    1122             : 
    1123             : INLINE GEN
    1124   115448748 : mulis(GEN x, long s) { return mulsi(s,x); }
    1125             : INLINE GEN
    1126   554601405 : muliu(GEN x, ulong s) { return mului(s,x); }
    1127             : INLINE GEN
    1128     4245258 : mulru(GEN x, ulong s) { return mulur(s,x); }
    1129             : INLINE GEN
    1130    13883157 : mulri(GEN x, GEN s) { return mulir(s,x); }
    1131             : INLINE GEN
    1132    10015339 : mulrs(GEN x, long s) { return mulsr(s,x); }
    1133             : 
    1134             : /*******************************************************************/
    1135             : /*                                                                 */
    1136             : /*                  VALUATION, EXPONENT, SHIFTS                    */
    1137             : /*                                                                 */
    1138             : /*******************************************************************/
    1139             : INLINE long
    1140    50891141 : vali(GEN x)
    1141             : {
    1142             :   long i;
    1143             :   GEN xp;
    1144             : 
    1145    50891141 :   if (!signe(x)) return -1;
    1146    50891141 :   xp=int_LSW(x);
    1147    50891141 :   for (i=0; !*xp; i++) xp=int_nextW(xp);
    1148    50891141 :   return vals(*xp) + i * BITS_IN_LONG;
    1149             : }
    1150             : 
    1151             : 
    1152             : /* assume x > 0 */
    1153             : INLINE long
    1154   321933145 : expu(ulong x) { return (BITS_IN_LONG-1) - (long)bfffo(x); }
    1155             : 
    1156             : INLINE long
    1157   725489583 : expi(GEN x)
    1158             : {
    1159   725489583 :   const long lx=lgefint(x);
    1160   725489583 :   return lx==2? -(long)HIGHEXPOBIT: bit_accuracy(lx)-(long)bfffo(*int_MSW(x))-1;
    1161             : }
    1162             : 
    1163             : INLINE GEN
    1164    47787100 : shiftr(GEN x, long n)
    1165             : {
    1166    47787100 :   const long e = evalexpo(expo(x)+n);
    1167    47787100 :   const GEN y = rcopy(x);
    1168             : 
    1169    47787100 :   if (e & ~EXPOBITS) pari_err_OVERFLOW("expo()");
    1170    47787100 :   y[1] = (y[1]&~EXPOBITS) | e; return y;
    1171             : }
    1172             : INLINE GEN
    1173    10960567 : mpshift(GEN x,long s) { return (typ(x)==t_INT)?shifti(x,s):shiftr(x,s); }
    1174             : 
    1175             : /* FIXME: adapt/use mpn_[lr]shift instead */
    1176             : /* z2[imin..imax] := z1[imin..imax].f shifted left sh bits
    1177             :  * (feeding f from the right). Assume sh > 0 */
    1178             : INLINE void
    1179  2545982268 : shift_left(GEN z2, GEN z1, long imin, long imax, ulong f,  ulong sh)
    1180             : {
    1181  2545982268 :   GEN sb = z1 + imin, se = z1 + imax, te = z2 + imax;
    1182  2545982268 :   ulong l, m = BITS_IN_LONG - sh, k = f >> m;
    1183 12436143718 :   while (se > sb) {
    1184  7344179182 :     l     = *se--;
    1185  7344179182 :     *te-- = (l << sh) | k;
    1186  7344179182 :     k     = l >> m;
    1187             :   }
    1188  2545982268 :   *te = (((ulong)*se) << sh) | k;
    1189  2545982268 : }
    1190             : /* z2[imin..imax] := f.z1[imin..imax-1] shifted right sh bits
    1191             :  * (feeding f from the left). Assume sh > 0 */
    1192             : INLINE void
    1193  2333963362 : shift_right(GEN z2, GEN z1, long imin, long imax, ulong f, ulong sh)
    1194             : {
    1195  2333963362 :   GEN sb = z1 + imin, se = z1 + imax, tb = z2 + imin;
    1196  2333963362 :   ulong k, l = *sb++, m = BITS_IN_LONG - sh;
    1197  2333963362 :   *tb++ = (l >> sh) | (f << m);
    1198 11704751975 :   while (sb < se) {
    1199  7036825251 :     k     = l << m;
    1200  7036825251 :     l     = *sb++;
    1201  7036825251 :     *tb++ = (l >> sh) | k;
    1202             :   }
    1203  2333963362 : }
    1204             : 
    1205             : /* Backward compatibility. Inefficient && unused */
    1206             : extern ulong hiremainder;
    1207             : INLINE ulong
    1208           0 : shiftl(ulong x, ulong y)
    1209           0 : { hiremainder = x>>(BITS_IN_LONG-y); return (x<<y); }
    1210             : 
    1211             : INLINE ulong
    1212           0 : shiftlr(ulong x, ulong y)
    1213           0 : { hiremainder = x<<(BITS_IN_LONG-y); return (x>>y); }
    1214             : 
    1215             : INLINE void
    1216    78234244 : shiftr_inplace(GEN z, long d)
    1217             : {
    1218    78234244 :   setexpo(z, expo(z)+d);
    1219    78234244 : }
    1220             : 
    1221             : /*******************************************************************/
    1222             : /*                                                                 */
    1223             : /*                           ASSIGNMENT                            */
    1224             : /*                                                                 */
    1225             : /*******************************************************************/
    1226             : INLINE void
    1227    37993730 : affii(GEN x, GEN y)
    1228             : {
    1229    37993730 :   long lx = lgefint(x);
    1230    37993730 :   if (lg(y)<lx) pari_err_OVERFLOW("t_INT-->t_INT assignment");
    1231    37993721 :   while (--lx) y[lx] = x[lx];
    1232    37993721 : }
    1233             : INLINE void
    1234      291275 : affsi(long s, GEN x)
    1235             : {
    1236      291275 :   if (!s) x[1] = evalsigne(0) | evallgefint(2);
    1237             :   else
    1238             :   {
    1239      275947 :     if (s > 0) { x[1] = evalsigne( 1) | evallgefint(3); x[2] =  s; }
    1240       91969 :     else       { x[1] = evalsigne(-1) | evallgefint(3); x[2] = -s; }
    1241             :   }
    1242      291275 : }
    1243             : INLINE void
    1244     8404821 : affui(ulong u, GEN x)
    1245             : {
    1246     8404821 :   if (!u) x[1] = evalsigne(0) | evallgefint(2);
    1247     8403547 :   else  { x[1] = evalsigne(1) | evallgefint(3); x[2] = u; }
    1248     8404821 : }
    1249             : 
    1250             : INLINE void
    1251   143193072 : affsr(long x, GEN y)
    1252             : {
    1253   143193072 :   long sh, i, ly = lg(y);
    1254             : 
    1255   143193072 :   if (!x)
    1256             :   {
    1257     1366460 :     y[1] = evalexpo(-prec2nbits(ly));
    1258   144559532 :     return;
    1259             :   }
    1260   141826612 :   if (x < 0) {
    1261      169018 :     x = -x; sh = bfffo(x);
    1262      169018 :     y[1] = evalsigne(-1) | _evalexpo((BITS_IN_LONG-1)-sh);
    1263             :   }
    1264             :   else
    1265             :   {
    1266   141657594 :     sh = bfffo(x);
    1267   141657594 :     y[1] = evalsigne(1) | _evalexpo((BITS_IN_LONG-1)-sh);
    1268             :   }
    1269   141826612 :   y[2] = ((ulong)x)<<sh; for (i=3; i<ly; i++) y[i]=0;
    1270             : }
    1271             : 
    1272             : INLINE void
    1273     5487408 : affur(ulong x, GEN y)
    1274             : {
    1275     5487408 :   long sh, i, ly = lg(y);
    1276             : 
    1277     5487408 :   if (!x)
    1278             :   {
    1279         266 :     y[1] = evalexpo(-prec2nbits(ly));
    1280     5487674 :     return;
    1281             :   }
    1282     5487142 :   sh = bfffo(x);
    1283     5487142 :   y[1] = evalsigne(1) | _evalexpo((BITS_IN_LONG-1)-sh);
    1284     5487142 :   y[2] = x<<sh; for (i=3; i<ly; i++) y[i] = 0;
    1285             : }
    1286             : 
    1287             : INLINE void
    1288      106100 : affiz(GEN x, GEN y) { if (typ(y)==t_INT) affii(x,y); else affir(x,y); }
    1289             : INLINE void
    1290           0 : affsz(long x, GEN y) { if (typ(y)==t_INT) affsi(x,y); else affsr(x,y); }
    1291             : INLINE void
    1292      133373 : mpaff(GEN x, GEN y) { if (typ(x)==t_INT) affiz(x, y); else affrr(x,y); }
    1293             : 
    1294             : /*******************************************************************/
    1295             : /*                                                                 */
    1296             : /*                    OPERATION + ASSIGNMENT                       */
    1297             : /*                                                                 */
    1298             : /*******************************************************************/
    1299             : 
    1300           0 : INLINE void addiiz(GEN x, GEN y, GEN z)
    1301           0 : { pari_sp av = avma; affii(addii(x,y),z); avma = av; }
    1302           0 : INLINE void addirz(GEN x, GEN y, GEN z)
    1303           0 : { pari_sp av = avma; affrr(addir(x,y),z); avma = av; }
    1304           0 : INLINE void addriz(GEN x, GEN y, GEN z)
    1305           0 : { pari_sp av = avma; affrr(addri(x,y),z); avma = av; }
    1306     1411916 : INLINE void addrrz(GEN x, GEN y, GEN z)
    1307     1411916 : { pari_sp av = avma; affrr(addrr(x,y),z); avma = av; }
    1308           0 : INLINE void addsiz(long s, GEN y, GEN z)
    1309           0 : { pari_sp av = avma; affii(addsi(s,y),z); avma = av; }
    1310           0 : INLINE void addsrz(long s, GEN y, GEN z)
    1311           0 : { pari_sp av = avma; affrr(addsr(s,y),z); avma = av; }
    1312           0 : INLINE void addssz(long s, long y, GEN z)
    1313           0 : { pari_sp av = avma; affii(addss(s,y),z); avma = av; }
    1314             : 
    1315          21 : INLINE void diviiz(GEN x, GEN y, GEN z)
    1316          21 : { pari_sp av = avma; affii(divii(x,y),z); avma = av; }
    1317           0 : INLINE void divirz(GEN x, GEN y, GEN z)
    1318           0 : { pari_sp av = avma; mpaff(divir(x,y),z); avma = av; }
    1319           0 : INLINE void divisz(GEN x, long y, GEN z)
    1320           0 : { pari_sp av = avma; affii(divis(x,y),z); avma = av; }
    1321           0 : INLINE void divriz(GEN x, GEN y, GEN z)
    1322           0 : { pari_sp av = avma; affrr(divri(x,y),z); avma = av; }
    1323         354 : INLINE void divrrz(GEN x, GEN y, GEN z)
    1324         354 : { pari_sp av = avma; affrr(divrr(x,y),z); avma = av; }
    1325       90548 : INLINE void divrsz(GEN y, long s, GEN z)
    1326       90548 : { pari_sp av = avma; affrr(divrs(y,s),z); avma = av; }
    1327           0 : INLINE void divsiz(long x, GEN y, GEN z)
    1328           0 : { long junk; affsi(sdivsi_rem(x,y,&junk), z); }
    1329           0 : INLINE void divsrz(long s, GEN y, GEN z)
    1330           0 : { pari_sp av = avma; mpaff(divsr(s,y),z); avma = av; }
    1331           0 : INLINE void divssz(long x, long y, GEN z)
    1332           0 : { affsi(x/y, z); }
    1333             : 
    1334           0 : INLINE void modisz(GEN y, long s, GEN z)
    1335           0 : { pari_sp av = avma; affii(modis(y,s),z); avma = av; }
    1336           0 : INLINE void modsiz(long s, GEN y, GEN z)
    1337           0 : { pari_sp av = avma; affii(modsi(s,y),z); avma = av; }
    1338           0 : INLINE void modssz(long s, long y, GEN z)
    1339           0 : { pari_sp av = avma; affii(modss(s,y),z); avma = av; }
    1340             : 
    1341           0 : INLINE void mpaddz(GEN x, GEN y, GEN z)
    1342           0 : { pari_sp av = avma; mpaff(mpadd(x,y),z); avma = av; }
    1343           0 : INLINE void mpsubz(GEN x, GEN y, GEN z)
    1344           0 : { pari_sp av = avma; mpaff(mpsub(x,y),z); avma = av; }
    1345           0 : INLINE void mpmulz(GEN x, GEN y, GEN z)
    1346           0 : { pari_sp av = avma; mpaff(mpmul(x,y),z); avma = av; }
    1347             : 
    1348           0 : INLINE void muliiz(GEN x, GEN y, GEN z)
    1349           0 : { pari_sp av = avma; affii(mulii(x,y),z); avma = av; }
    1350           0 : INLINE void mulirz(GEN x, GEN y, GEN z)
    1351           0 : { pari_sp av = avma; mpaff(mulir(x,y),z); avma = av; }
    1352           0 : INLINE void mulriz(GEN x, GEN y, GEN z)
    1353           0 : { pari_sp av = avma; mpaff(mulri(x,y),z); avma = av; }
    1354       56940 : INLINE void mulrrz(GEN x, GEN y, GEN z)
    1355       56940 : { pari_sp av = avma; affrr(mulrr(x,y),z); avma = av; }
    1356           0 : INLINE void mulsiz(long s, GEN y, GEN z)
    1357           0 : { pari_sp av = avma; affii(mulsi(s,y),z); avma = av; }
    1358           0 : INLINE void mulsrz(long s, GEN y, GEN z)
    1359           0 : { pari_sp av = avma; mpaff(mulsr(s,y),z); avma = av; }
    1360           0 : INLINE void mulssz(long s, long y, GEN z)
    1361           0 : { pari_sp av = avma; affii(mulss(s,y),z); avma = av; }
    1362             : 
    1363           0 : INLINE void remiiz(GEN x, GEN y, GEN z)
    1364           0 : { pari_sp av = avma; affii(remii(x,y),z); avma = av; }
    1365           0 : INLINE void remisz(GEN y, long s, GEN z)
    1366           0 : { pari_sp av = avma; affii(remis(y,s),z); avma = av; }
    1367           0 : INLINE void remsiz(long s, GEN y, GEN z)
    1368           0 : { pari_sp av = avma; affii(remsi(s,y),z); avma = av; }
    1369           0 : INLINE void remssz(long s, long y, GEN z)
    1370           0 : { pari_sp av = avma; affii(remss(s,y),z); avma = av; }
    1371             : 
    1372           0 : INLINE void subiiz(GEN x, GEN y, GEN z)
    1373           0 : { pari_sp av = avma; affii(subii(x,y),z); avma = av; }
    1374           0 : INLINE void subirz(GEN x, GEN y, GEN z)
    1375           0 : { pari_sp av = avma; affrr(subir(x,y),z); avma = av; }
    1376           0 : INLINE void subisz(GEN y, long s, GEN z)
    1377           0 : { pari_sp av = avma; affii(addsi(-s,y),z); avma = av; }
    1378           0 : INLINE void subriz(GEN x, GEN y, GEN z)
    1379           0 : { pari_sp av = avma; affrr(subri(x,y),z); avma = av; }
    1380     1295642 : INLINE void subrrz(GEN x, GEN y, GEN z)
    1381     1295642 : { pari_sp av = avma; affrr(subrr(x,y),z); avma = av; }
    1382           0 : INLINE void subrsz(GEN y, long s, GEN z)
    1383           0 : { pari_sp av = avma; affrr(addsr(-s,y),z); avma = av; }
    1384           0 : INLINE void subsiz(long s, GEN y, GEN z)
    1385           0 : { pari_sp av = avma; affii(subsi(s,y),z); avma = av; }
    1386           0 : INLINE void subsrz(long s, GEN y, GEN z)
    1387           0 : { pari_sp av = avma; affrr(subsr(s,y),z); avma = av; }
    1388           0 : INLINE void subssz(long x, long y, GEN z) { addssz(x,-y,z); }
    1389             : 
    1390             : INLINE void
    1391           0 : dvmdssz(long x, long y, GEN z, GEN t) {
    1392           0 :   pari_sp av = avma;
    1393             :   long r;
    1394           0 :   affii(divss_rem(x,y, &r), z); avma = av; affsi(r,t);
    1395           0 : }
    1396             : INLINE void
    1397           0 : dvmdsiz(long x, GEN y, GEN z, GEN t) {
    1398           0 :   pari_sp av = avma;
    1399             :   long r;
    1400           0 :   affii(divsi_rem(x,y, &r), z); avma = av; affsi(r,t);
    1401           0 : }
    1402             : INLINE void
    1403           0 : dvmdisz(GEN x, long y, GEN z, GEN t) {
    1404           0 :   pari_sp av = avma;
    1405             :   long r;
    1406           0 :   affii(divis_rem(x,y, &r),z); avma = av; affsi(r,t);
    1407           0 : }
    1408             : INLINE void
    1409           0 : dvmdiiz(GEN x, GEN y, GEN z, GEN t) {
    1410           0 :   pari_sp av = avma;
    1411             :   GEN r;
    1412           0 :   affii(dvmdii(x,y,&r),z); affii(r,t); avma=av;
    1413           0 : }

Generated by: LCOV version 1.11