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.8.0 lcov report (development 19222-cc5b867) Lines: 524 684 76.6 %
Date: 2016-07-28 07:10:28 Functions: 186 248 75.0 %
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 41846918062 : evallg(long x)
      21             : {
      22 41846918062 :   if (x & ~LGBITS) pari_err_OVERFLOW("lg()");
      23 41851517753 :   return _evallg(x);
      24             : }
      25             : INLINE long
      26    11750903 : evalvalp(long x)
      27             : {
      28    11750903 :   long v = _evalvalp(x);
      29    11750903 :   if (v & ~VALPBITS) pari_err_OVERFLOW("valp()");
      30    11750901 :   return v;
      31             : }
      32             : INLINE long
      33  6475641648 : evalexpo(long x)
      34             : {
      35  6475641648 :   long v = _evalexpo(x);
      36  6475641648 :   if (v & ~EXPOBITS) pari_err_OVERFLOW("expo()");
      37  6475641648 :   return v;
      38             : }
      39             : INLINE long
      40     1278948 : evalprecp(long x)
      41             : {
      42     1278948 :   long v = _evalprecp(x);
      43     1278948 :   if (x & ~((1UL<<(BITS_IN_LONG-VALPnumBITS))-1)) pari_err_OVERFLOW("precp()");
      44     1278948 :   return v;
      45             : }
      46             : 
      47             : INLINE int
      48   134680885 : varncmp(long x, long y)
      49             : {
      50   134680885 :   if (varpriority[x] < varpriority[y]) return  1;
      51   105148832 :   if (varpriority[x] > varpriority[y]) return -1;
      52    55116986 :   return 0;
      53             : }
      54             : 
      55             : /* Inhibit some area gerepile-wise: declare it to be a non recursive
      56             :  * type, of length l. Thus gerepile won't inspect the zone, just copy it.
      57             :  * For the following situation:
      58             :  *   z = cgetg(t,a); av = avma; garbage(); ltop = avma;
      59             :  *   for (i=1; i<HUGE; i++) gel(z,i) = blah();
      60             :  *   stackdummy(av,ltop);
      61             :  * loses (av-ltop) words but save a costly gerepile. */
      62             : INLINE void
      63  1699079747 : stackdummy(pari_sp av, pari_sp ltop) {
      64  1699079747 :   long l = ((GEN)av) - ((GEN)ltop);
      65  1699079747 :   if (l > 0) {
      66   628419492 :     GEN z = (GEN)ltop;
      67   628419492 :     z[0] = evaltyp(t_VECSMALL) | evallg(l);
      68             : #ifdef DEBUG
      69             :     { long i; for (i = 1; i < l; i++) z[i] = 0; }
      70             : #endif
      71             :   }
      72  1699450429 : }
      73             : INLINE void
      74    14282555 : fixlg(GEN x, long ly) {
      75    14282555 :   long lx = lg(x), l = lx - ly;
      76    14282555 :   if (l > 0)
      77             :   { /* stackdummy(x+lx, x+ly) */
      78        7383 :     GEN z = x + ly;
      79        7383 :     z[0] = evaltyp(t_VECSMALL) | evallg(l);
      80        7383 :     setlg(x, ly);
      81             : #ifdef DEBUG
      82             :     { long i; for (i = 1; i < l; i++) z[i] = 0; }
      83             : #endif
      84             :   }
      85    14282555 : }
      86             : /* update lg(z) before affrr(y, z)  [ to cater for precision loss ]*/
      87             : INLINE void
      88    14277266 : affrr_fixlg(GEN y, GEN z) { fixlg(z, lg(y)); affrr(y, z); }
      89             : 
      90             : /*******************************************************************/
      91             : /*                                                                 */
      92             : /*                       ALLOCATE ON STACK                         */
      93             : /*                                                                 */
      94             : /*******************************************************************/
      95             : INLINE GEN
      96 39331792944 : new_chunk(size_t x) /* x is a number of longs */
      97             : {
      98 39331792944 :   GEN z = ((GEN) avma) - x;
      99             :   CHECK_CTRLC
     100 39331792944 :   if (x > (avma-pari_mainstack->bot) / sizeof(long))
     101           0 :     new_chunk_resize(x);
     102 39343238115 :   avma = (pari_sp)z;
     103             : #ifdef MEMSTEP
     104             :   if (DEBUGMEM && pari_mainstack->memused != DISABLE_MEMUSED) {
     105             :     long d = (long)pari_mainstack->memused - (long)z;
     106             :     if (labs(d) > 4*MEMSTEP)
     107             :     {
     108             :       pari_mainstack->memused = (pari_sp)z;
     109             :       err_printf("...%4.0lf Mbytes used\n",
     110             :                 (pari_mainstack->top-pari_mainstack->memused)/1048576.);
     111             :     }
     112             :   }
     113             : #endif
     114 39343238115 :   return z;
     115             : }
     116             : 
     117             : INLINE char *
     118     4096253 : stack_malloc(size_t N)
     119             : {
     120     4096253 :   long n = nchar2nlong(N);
     121     4096253 :   return (char*)new_chunk(n);
     122             : }
     123             : 
     124             : INLINE char *
     125       21646 : stack_calloc(size_t N)
     126             : {
     127       21646 :   char *p = stack_malloc(N);
     128       21646 :   memset(p, 0, N); return p;
     129             : }
     130             : 
     131             : /* cgetg(lg(x), typ(x)), set *lx. Implicit unsetisclone() */
     132             : INLINE GEN
     133   346804390 : cgetg_copy(GEN x, long *plx) {
     134             :   GEN y;
     135   346804390 :   *plx = lg(x); y = new_chunk((size_t)*plx);
     136   346804075 :   y[0] = x[0] & (TYPBITS|LGBITS); return y;
     137             : }
     138             : INLINE GEN
     139       31500 : cgetg_block(long x, long y)
     140             : {
     141       31500 :   GEN z = newblock((size_t)x);
     142       31500 :   z[0] = evaltyp(y) | evallg(x);
     143       31500 :   return z;
     144             : }
     145             : INLINE GEN
     146  4219972186 : cgetg(long x, long y)
     147             : {
     148  4219972186 :   GEN z = new_chunk((size_t)x);
     149  4220389971 :   z[0] = evaltyp(y) | evallg(x);
     150  4218310010 :   return z;
     151             : }
     152             : INLINE GEN
     153 13165271819 : cgeti(long x)
     154             : {
     155 13165271819 :   GEN z = new_chunk((size_t)x);
     156 13151750627 :   z[0] = evaltyp(t_INT) | evallg(x);
     157 13141771374 :   return z;
     158             : }
     159             : INLINE GEN
     160  8097905464 : cgetipos(long x)
     161             : {
     162  8097905464 :   GEN z = cgeti(x);
     163  8096684179 :   z[1] = evalsigne(1) | evallgefint(x);
     164  8096684179 :   return z;
     165             : }
     166             : INLINE GEN
     167   199386087 : cgetineg(long x)
     168             : {
     169   199386087 :   GEN z = cgeti(x);
     170   199386088 :   z[1] = evalsigne(-1) | evallgefint(x);
     171   199386088 :   return z;
     172             : }
     173             : INLINE GEN
     174        2702 : cgetr_block(long x)
     175             : {
     176        2702 :   GEN z = newblock((size_t)x);
     177        2702 :   z[0] = evaltyp(t_REAL) | evallg(x);
     178        2702 :   return z;
     179             : }
     180             : INLINE GEN
     181  6651157077 : cgetr(long x)
     182             : {
     183  6651157077 :   GEN z = new_chunk((size_t)x);
     184  6651157077 :   z[0] = evaltyp(t_REAL) | evallg(x);
     185  6651157077 :   return z;
     186             : }
     187             : 
     188             : /*******************************************************************/
     189             : /*                                                                 */
     190             : /*                     COPY, NEGATION, ABSOLUTE VALUE              */
     191             : /*                                                                 */
     192             : /*******************************************************************/
     193             : /* cannot do memcpy because sometimes x and y overlap */
     194             : INLINE GEN
     195  1086360730 : leafcopy(GEN x)
     196             : {
     197  1086360730 :   register long lx = lg(x);
     198  1086360730 :   GEN y = new_chunk(lx); /* can't use cgetg_copy, in case x,y overlap */
     199  1086365858 :   while (--lx > 0) y[lx] = x[lx];
     200  1086365858 :   y[0] = x[0] & (TYPBITS|LGBITS); return y;
     201             : }
     202             : INLINE GEN
     203  4199777099 : icopy(GEN x)
     204             : {
     205  4199777099 :   long i = lgefint(x), lx = i;
     206  4199777099 :   GEN y = new_chunk(lx); /* can't use cgeti, in case x,y overlap */
     207  4202597680 :   while (--i > 0) y[i] = x[i];
     208  4202597680 :   y[0] = evaltyp(t_INT) | evallg(lx);
     209  4200922646 :   return y;
     210             : }
     211             : INLINE GEN
     212    39672495 : icopyspec(GEN x, long nx)
     213             : {
     214    39672495 :   long i = nx+2, lx = i;
     215    39672495 :   GEN y = new_chunk(lx); /* can't use cgeti, in case x,y overlap */
     216    39672474 :   x -= 2; while (--i >= 2) y[i] = x[i];
     217    39672474 :   y[1] = evalsigne(1) | evallgefint(lx);
     218    39672474 :   y[0] = evaltyp(t_INT) | evallg(lx);
     219    39672460 :   return y;
     220             : }
     221   170380169 : INLINE GEN rcopy(GEN x) { return leafcopy(x); }
     222     1150649 : INLINE GEN mpcopy(GEN x) { return leafcopy(x); }
     223             : 
     224             : INLINE GEN
     225   264465618 : mpabs(GEN x) { GEN y = leafcopy(x); setabssign(y); return y; }
     226             : INLINE GEN
     227         294 : mpabs_shallow(GEN x) { return signe(x) < 0? mpabs(x): x; }
     228   221192625 : INLINE GEN absi(GEN x) { return mpabs(x); }
     229      108749 : INLINE GEN absi_shallow(GEN x) { return signe(x) < 0? negi(x): x; }
     230    29078552 : INLINE GEN absr(GEN x) { return mpabs(x); }
     231             : 
     232             : INLINE GEN
     233   313737754 : mpneg(GEN x) { GEN y = leafcopy(x); togglesign(y); return y; }
     234   218152716 : INLINE GEN negi(GEN x) { return mpneg(x); }
     235     1171476 : INLINE GEN negr(GEN x) { return mpneg(x); }
     236             : 
     237             : /* negate in place */
     238             : INLINE void
     239   774036417 : togglesign(GEN x) { if (x[1] & SIGNBITS) { x[1] ^= HIGHBIT; } }
     240             : INLINE void
     241   296055831 : setabssign(GEN x) { x[1] &= ~HIGHBIT; }
     242             : /* negate in place, except universal constants */
     243             : INLINE void
     244    42490852 : togglesign_safe(GEN *px)
     245             : {
     246    42490852 :   switch(*px - gen_1) /* gen_1, gen_2, gen_m1, gen_m2 */
     247             :   {
     248      825690 :     case 0: *px = gen_m1; break;
     249           0 :     case 3: *px = gen_m2;  break;
     250      191816 :     case 6: *px = gen_1; break;
     251           0 :     case 9: *px = gen_2;  break;
     252    41473346 :     default: togglesign(*px);
     253             :   }
     254    42490852 : }
     255             : /* setsigne(y, signe(x)) */
     256             : INLINE void
     257           0 : affectsign(GEN x, GEN y)
     258             : {
     259           0 :   y[1] = (x[1] & SIGNBITS) | (y[1] & ~SIGNBITS);
     260           0 : }
     261             : /* copies sign in place, except for universal constants */
     262             : INLINE void
     263     5172701 : affectsign_safe(GEN x, GEN *py)
     264             : {
     265     5172701 :   if (((*py)[1] ^ x[1]) & HIGHBIT) togglesign_safe(py);
     266     5172701 : }
     267             : /*******************************************************************/
     268             : /*                                                                 */
     269             : /*                     GEN -> LONG, LONG -> GEN                    */
     270             : /*                                                                 */
     271             : /*******************************************************************/
     272             : /* assume x != 0, return -x as a t_INT */
     273             : INLINE GEN
     274   199251589 : utoineg(ulong x) { GEN y = cgetineg(3); y[2] = x; return y; }
     275             : /* assume x != 0, return utoi(x) */
     276             : INLINE GEN
     277  7081698836 : utoipos(ulong x) { GEN y = cgetipos(3); y[2] = x; return y; }
     278             : INLINE GEN
     279  6093172993 : utoi(ulong x) { return x? utoipos(x): gen_0; }
     280             : INLINE GEN
     281   282578645 : stoi(long x)
     282             : {
     283   282578645 :   if (!x) return gen_0;
     284   186606179 :   return x > 0? utoipos((ulong)x): utoineg((ulong)-x);
     285             : }
     286             : 
     287             : /* x 2^BIL + y */
     288             : INLINE GEN
     289  4139594736 : uutoi(ulong x, ulong y)
     290             : {
     291             :   GEN z;
     292  4139594736 :   if (!x) return utoi(y);
     293   378885394 :   z = cgetipos(4);
     294   378436135 :   *int_W_lg(z, 1, 4) = x;
     295   378436135 :   *int_W_lg(z, 0, 4) = y; return z;
     296             : }
     297             : /* - (x 2^BIL + y) */
     298             : INLINE GEN
     299       71324 : uutoineg(ulong x, ulong y)
     300             : {
     301             :   GEN z;
     302       71324 :   if (!x) return y? utoineg(y): gen_0;
     303           0 :   z = cgetineg(4);
     304           0 :   *int_W_lg(z, 1, 4) = x;
     305           0 :   *int_W_lg(z, 0, 4) = y; return z;
     306             : }
     307             : 
     308             : INLINE long
     309   200360687 : itos(GEN x)
     310             : {
     311   200360687 :   long s = signe(x);
     312             :   long u;
     313             : 
     314   200360687 :   if (!s) return 0;
     315   189013714 :   u = x[2];
     316   189013714 :   if (lgefint(x) > 3 || u < 0)
     317          21 :     pari_err_OVERFLOW("t_INT-->long assignment");
     318   189013701 :   return (s>0) ? u : -u;
     319             : }
     320             : /* as itos, but return 0 if too large. Cf is_bigint */
     321             : INLINE long
     322     7200131 : itos_or_0(GEN x) {
     323             :   long n;
     324     7200131 :   if (lgefint(x) != 3 || (n = x[2]) & HIGHBIT) return 0;
     325     7181315 :   return signe(x) > 0? n: -n;
     326             : }
     327             : INLINE ulong
     328    25605864 : itou(GEN x)
     329             : {
     330    25605864 :   switch(lgefint(x)) {
     331     5179553 :     case 2: return 0;
     332    20426311 :     case 3: return x[2];
     333             :     default:
     334           0 :       pari_err_OVERFLOW("t_INT-->ulong assignment");
     335           0 :       return 0; /* not reached */
     336             :   }
     337             : }
     338             : 
     339             : /* as itou, but return 0 if too large. Cf is_bigint */
     340             : INLINE ulong
     341     4166017 : itou_or_0(GEN x) {
     342     4166017 :   if (lgefint(x) != 3) return 0;
     343     4156003 :   return (ulong)x[2];
     344             : }
     345             : 
     346             : INLINE GEN
     347    91159866 : real_0_bit(long bitprec) { GEN x=cgetr(2); x[1]=evalexpo(bitprec); return x; }
     348             : INLINE GEN
     349    30216733 : real_0(long prec) { return real_0_bit(-prec2nbits(prec)); }
     350             : INLINE GEN
     351      476897 : real_1_bit(long bit) { return real_1(nbits2prec(bit)); }
     352             : INLINE GEN
     353    26529492 : real_1(long prec) {
     354    26529492 :   GEN x = cgetr(prec);
     355             :   long i;
     356    26529492 :   x[1] = evalsigne(1) | _evalexpo(0);
     357    26529492 :   x[2] = (long)HIGHBIT; for (i=3; i<prec; i++) x[i] = 0;
     358    26529492 :   return x;
     359             : }
     360             : INLINE GEN
     361          84 : real_m1(long prec) {
     362          84 :   GEN x = cgetr(prec);
     363             :   long i;
     364          84 :   x[1] = evalsigne(-1) | _evalexpo(0);
     365          84 :   x[2] = (long)HIGHBIT; for (i=3; i<prec; i++) x[i] = 0;
     366          84 :   return x;
     367             : }
     368             : 
     369             : /* 2.^n */
     370             : INLINE GEN
     371       48733 : real2n(long n, long prec) { GEN z = real_1(prec); setexpo(z, n); return z; }
     372             : INLINE GEN
     373           0 : real_m2n(long n, long prec) { GEN z = real_m1(prec); setexpo(z, n); return z; }
     374             : INLINE GEN
     375   144380931 : stor(long s, long prec) { GEN z = cgetr(prec); affsr(s,z); return z; }
     376             : INLINE GEN
     377     1735380 : utor(ulong s, long prec){ GEN z = cgetr(prec); affur(s,z); return z; }
     378             : INLINE GEN
     379   176803465 : itor(GEN x, long prec) { GEN z = cgetr(prec); affir(x,z); return z; }
     380             : INLINE GEN
     381    71329795 : rtor(GEN x, long prec) { GEN z = cgetr(prec); affrr(x,z); return z; }
     382             : 
     383     5250688 : INLINE ulong int_bit(GEN x, long n)
     384             : {
     385     5250688 :   long r, q = dvmdsBIL(n, &r);
     386     5251793 :   return q < lgefint(x)-2?((ulong)*int_W(x,q) >> r) & 1UL:0;
     387             : }
     388             : 
     389             : /*******************************************************************/
     390             : /*                                                                 */
     391             : /*                           COMPARISON                            */
     392             : /*                                                                 */
     393             : /*******************************************************************/
     394             : INLINE int
     395      705940 : cmpir(GEN x, GEN y)
     396             : {
     397             :   pari_sp av;
     398             :   GEN z;
     399             : 
     400      705940 :   if (!signe(x)) return -signe(y);
     401      336768 :   if (!signe(y))
     402             :   {
     403        3915 :     if (expo(y) >= expi(x)) return 0;
     404        3887 :     return signe(x);
     405             :   }
     406      332853 :   av=avma; z = itor(x, realprec(y)); avma=av;
     407      332853 :   return cmprr(z,y); /* cmprr does no memory adjustment */
     408             : }
     409             : INLINE int
     410      474708 : cmpri(GEN x, GEN y) { return -cmpir(y,x); }
     411             : INLINE int
     412       61690 : cmpsr(long x, GEN y)
     413             : {
     414             :   pari_sp av;
     415             :   GEN z;
     416             : 
     417       61690 :   if (!x) return -signe(y);
     418       61690 :   av=avma; z = stor(x, LOWDEFAULTPREC); avma=av;
     419       61690 :   return cmprr(z,y);
     420             : }
     421             : INLINE int
     422       38640 : cmprs(GEN x, long y) { return -cmpsr(y,x); }
     423             : /* compare x and |y| */
     424             : INLINE int
     425    26400404 : cmpui(ulong x, GEN y)
     426             : {
     427    26400404 :   long l = lgefint(y);
     428             :   ulong p;
     429             : 
     430    26400404 :   if (!x) return (l > 2)? -1: 0;
     431    26400390 :   if (l == 2) return 1;
     432    26210527 :   if (l > 3) return -1;
     433    26200447 :   p = y[2]; if (p == x) return 0;
     434    25809815 :   return p < x ? 1 : -1;
     435             : }
     436             : INLINE int
     437    26381904 : cmpiu(GEN x, ulong y) { return -cmpui(y,x); }
     438             : INLINE int
     439     4152945 : cmpsi(long x, GEN y)
     440             : {
     441             :   ulong p;
     442             : 
     443     4152945 :   if (!x) return -signe(y);
     444             : 
     445     4152343 :   if (x > 0)
     446             :   {
     447     4152042 :     if (signe(y)<=0) return 1;
     448     4151923 :     if (lgefint(y)>3) return -1;
     449     4138743 :     p = y[2]; if (p == (ulong)x) return 0;
     450     4104982 :     return p < (ulong)x ? 1 : -1;
     451             :   }
     452             : 
     453         301 :   if (signe(y)>=0) return -1;
     454          42 :   if (lgefint(y)>3) return 1;
     455          42 :   p = y[2]; if (p == (ulong)-x) return 0;
     456          14 :   return p < (ulong)(-x) ? -1 : 1;
     457             : }
     458             : INLINE int
     459     4151377 : cmpis(GEN x, long y) { return -cmpsi(y,x); }
     460             : INLINE int
     461      506622 : mpcmp(GEN x, GEN y)
     462             : {
     463      506622 :   if (typ(x)==t_INT)
     464        5696 :     return (typ(y)==t_INT) ? cmpii(x,y) : cmpir(x,y);
     465      500926 :   return (typ(y)==t_INT) ? -cmpir(y,x) : cmprr(x,y);
     466             : }
     467             : 
     468             : /* x == y ? */
     469             : INLINE int
     470       79576 : equalsi(long x, GEN y)
     471             : {
     472       79576 :   if (!x) return !signe(y);
     473       79576 :   if (x > 0)
     474             :   {
     475       10129 :     if (signe(y) <= 0 || lgefint(y) != 3) return 0;
     476        8792 :     return ((ulong)y[2] == (ulong)x);
     477             :   }
     478       69447 :   if (signe(y) >= 0 || lgefint(y) != 3) return 0;
     479       69447 :   return ((ulong)y[2] == (ulong)-x);
     480             : }
     481             : /* x == |y| ? */
     482             : INLINE int
     483    11104289 : equalui(ulong x, GEN y)
     484             : {
     485    11104289 :   if (!x) return !signe(y);
     486    11104289 :   return (lgefint(y) == 3 && (ulong)y[2] == x);
     487             : }
     488             : INLINE int
     489    11000949 : equaliu(GEN x, ulong y) { return equalui(y,x); }
     490             : INLINE int
     491       79457 : equalis(GEN x, long y) { return equalsi(y,x); }
     492             : 
     493             : /* assume x != 0, is |x| == 2^n ? */
     494             : INLINE int
     495      738180 : absrnz_equal2n(GEN x) {
     496      738180 :   if ((ulong)x[2]==HIGHBIT)
     497             :   {
     498       29165 :     long i, lx = lg(x);
     499      385379 :     for (i = 3; i < lx; i++)
     500      361786 :       if (x[i]) return 0;
     501       23593 :     return 1;
     502             :   }
     503      709015 :   return 0;
     504             : }
     505             : /* assume x != 0, is |x| == 1 ? */
     506             : INLINE int
     507     1170004 : absrnz_equal1(GEN x) { return !expo(x) && absrnz_equal2n(x); }
     508             : 
     509             : INLINE long
     510  2877707306 : maxss(long x, long y) { return x>y?x:y; }
     511             : INLINE long
     512   829236154 : minss(long x, long y) { return x<y?x:y; }
     513             : INLINE long
     514     3275971 : minuu(ulong x, ulong y) { return x<y?x:y; }
     515             : INLINE long
     516     6280996 : maxuu(ulong x, ulong y) { return x>y?x:y; }
     517             : INLINE double
     518      932588 : maxdd(double x, double y) { return x>y?x:y; }
     519             : INLINE double
     520      357992 : mindd(double x, double y) { return x<y?x:y; }
     521             : 
     522             : /*******************************************************************/
     523             : /*                                                                 */
     524             : /*                             ADD / SUB                           */
     525             : /*                                                                 */
     526             : /*******************************************************************/
     527             : INLINE GEN
     528       24990 : subuu(ulong x, ulong y)
     529             : {
     530             :   ulong z;
     531             :   LOCAL_OVERFLOW;
     532       24990 :   z = subll(x, y);
     533       24990 :   return overflow? utoineg(-z): utoi(z);
     534             : }
     535             : INLINE GEN
     536  1838418645 : adduu(ulong x, ulong y) { ulong t = x+y; return uutoi((t < x), t); }
     537             : 
     538             : INLINE GEN
     539       24990 : addss(long x, long y)
     540             : {
     541       24990 :   if (!x) return stoi(y);
     542       24990 :   if (!y) return stoi(x);
     543       24990 :   if (x > 0) return y > 0? adduu(x,y): subuu(x, -y);
     544             : 
     545       24990 :   if (y > 0) return subuu(y, -x);
     546             :   else { /* - adduu(-x, -y) */
     547           0 :     ulong t = (-x)+(-y); return uutoineg((t < (ulong)(-x)), t);
     548             :   }
     549             : }
     550       24990 : INLINE GEN subss(long x, long y) { return addss(-y,x); }
     551             : 
     552             : INLINE GEN
     553  4500529765 : subii(GEN x, GEN y)
     554             : {
     555  4500529765 :   if (x==y) return gen_0; /* frequent with x = y = gen_0 */
     556  3296032269 :   return addii_sign(x, signe(x), y, -signe(y));
     557             : }
     558             : INLINE GEN
     559  4808295741 : addii(GEN x, GEN y) { return addii_sign(x, signe(x), y, signe(y)); }
     560             : INLINE GEN
     561   747315833 : addrr(GEN x, GEN y) { return addrr_sign(x, signe(x), y, signe(y)); }
     562             : INLINE GEN
     563  1492603222 : subrr(GEN x, GEN y) { return addrr_sign(x, signe(x), y, -signe(y)); }
     564             : INLINE GEN
     565    88774033 : addir(GEN x, GEN y) { return addir_sign(x, signe(x), y, signe(y)); }
     566             : INLINE GEN
     567   115241026 : subir(GEN x, GEN y) { return addir_sign(x, signe(x), y, -signe(y)); }
     568             : INLINE GEN
     569     1004455 : subri(GEN x, GEN y) { return addir_sign(y, -signe(y), x, signe(x)); }
     570             : INLINE GEN
     571    98116392 : addsi(long x, GEN y) { return addsi_sign(x, y, signe(y)); }
     572             : INLINE GEN
     573    11765347 : addui(ulong x, GEN y) { return addui_sign(x, y, signe(y)); }
     574             : INLINE GEN
     575     3397806 : subsi(long x, GEN y) { return addsi_sign(x, y, -signe(y)); }
     576             : INLINE GEN
     577    17497417 : subui(ulong x, GEN y) { return addui_sign(x, y, -signe(y)); }
     578             : 
     579             : /*******************************************************************/
     580             : /*                                                                 */
     581             : /*                           MOD, REM, DIV                         */
     582             : /*                                                                 */
     583             : /*******************************************************************/
     584    83032916 : INLINE ulong mod2BIL(GEN x) { return *int_LSW(x); }
     585           0 : INLINE long mod64(GEN x) { return mod2BIL(x) & 63; }
     586     1492659 : INLINE long mod32(GEN x) { return mod2BIL(x) & 31; }
     587     3981089 : INLINE long mod16(GEN x) { return mod2BIL(x) & 15; }
     588    39418784 : INLINE long mod8(GEN x)  { return mod2BIL(x) & 7; }
     589     9381320 : INLINE long mod4(GEN x)  { return mod2BIL(x) & 3; }
     590    18096644 : INLINE long mod2(GEN x)  { return mod2BIL(x) & 1; }
     591             : INLINE int
     592    27643202 : mpodd(GEN x) { return signe(x) && mod2(x); }
     593             : 
     594             : INLINE GEN
     595    36663301 : truedivii(GEN a,GEN b) { return truedvmdii(a,b,NULL); }
     596             : INLINE GEN
     597         652 : truedivis(GEN a, long b) { return truedvmdis(a,b,NULL); }
     598             : INLINE GEN
     599     6083720 : truedivsi(long a, GEN b) { return truedvmdsi(a,b,NULL); }
     600             : 
     601             : INLINE GEN
     602     7975380 : divii(GEN a, GEN b) { return dvmdii(a,b,NULL); }
     603             : INLINE GEN
     604  1206686913 : remii(GEN a, GEN b) { return dvmdii(a,b,ONLY_REM); }
     605             : 
     606             : INLINE GEN
     607           0 : divss(long x, long y) { return stoi(x / y); }
     608             : INLINE GEN
     609      600364 : modss(long x, long y) { return stoi(smodss(x, y)); }
     610             : INLINE GEN
     611           0 : remss(long x, long y) { return stoi(x % y); }
     612             : INLINE long
     613     6371406 : smodss(long x, long y)
     614             : {
     615     6371406 :   long r = x%y;
     616     6371406 :   return (r >= 0)? r: labs(y) + r;
     617             : }
     618             : INLINE ulong
     619    11846095 : umodsu(long x, ulong y)
     620             : {
     621    11846095 :   return x>=0 ? x%y: Fl_neg((-x)%y, y);
     622             : }
     623             : 
     624             : INLINE long
     625           0 : sdivss_rem(long x, long y, long *r)
     626             : {
     627             :   long q;
     628             :   LOCAL_HIREMAINDER;
     629           0 :   if (!y) pari_err_INV("sdivss_rem",gen_0);
     630           0 :   hiremainder = 0; q = divll((ulong)labs(x),(ulong)labs(y));
     631           0 :   if (x < 0) { hiremainder = -((long)hiremainder); q = -q; }
     632           0 :   if (y < 0) q = -q;
     633           0 :   *r = hiremainder; return q;
     634             : }
     635             : INLINE GEN
     636           0 : divss_rem(long x, long y, long *r) { return stoi(sdivss_rem(x,y,r)); }
     637             : INLINE ulong
     638         427 : udivuu_rem(ulong x, ulong y, ulong *r)
     639             : {
     640         427 :   if (!y) pari_err_INV("udivuu_rem",gen_0);
     641         427 :   *r = x % y; return x / y;
     642             : }
     643             : 
     644             : INLINE ulong
     645       10780 : udivui_rem(ulong x, GEN y, ulong *r)
     646             : {
     647       10780 :   long q, s = signe(y);
     648             :   LOCAL_HIREMAINDER;
     649             : 
     650       10780 :   if (!s) pari_err_INV("udivui_rem",gen_0);
     651       10780 :   if (!x || lgefint(y)>3) { *r = x; return 0; }
     652       10500 :   hiremainder=0; q = (long)divll(x, (ulong)y[2]);
     653       10500 :   if (s < 0) q = -q;
     654       10500 :   *r = hiremainder; return q;
     655             : }
     656             : 
     657             : /* assume d != 0 and |n| / d can be represented as an ulong.
     658             :  * Return |n|/d, set *r = |n| % d */
     659             : INLINE ulong
     660    12480475 : udiviu_rem(GEN n, ulong d, ulong *r)
     661             : {
     662    12480475 :   switch(lgefint(n))
     663             :   {
     664           0 :     case 2: *r = 0; return 0;
     665             :     case 3:
     666             :     {
     667    12480475 :       ulong nn = n[2];
     668    12480475 :       *r = nn % d; return nn / d;
     669             :     }
     670             :     default: /* 4 */
     671             :     {
     672             :       ulong n1, n0, q;
     673             :       LOCAL_HIREMAINDER;
     674           0 :       n0 = *int_W(n,0);
     675           0 :       n1 = *int_W(n,1);
     676           0 :       hiremainder = n1;
     677           0 :       q = divll(n0, d);
     678           0 :       *r = hiremainder; return q;
     679             :     }
     680             :   }
     681             : }
     682             : 
     683             : INLINE long
     684    13557273 : sdivsi_rem(long x, GEN y, long *r)
     685             : {
     686    13557273 :   long q, s = signe(y);
     687             :   LOCAL_HIREMAINDER;
     688             : 
     689    13557273 :   if (!s) pari_err_INV("sdivsi_rem",gen_0);
     690    13557279 :   if (!x || lgefint(y)>3 || ((long)y[2]) < 0) { *r = x; return 0; }
     691    11882239 :   hiremainder=0; q = (long)divll(labs(x), (ulong)y[2]);
     692    11882239 :   if (x < 0) { hiremainder = -((long)hiremainder); q = -q; }
     693    11882239 :   if (s < 0) q = -q;
     694    11882239 :   *r = hiremainder; return q;
     695             : }
     696             : INLINE GEN
     697           0 : divsi_rem(long s, GEN y, long *r) { return stoi(sdivsi_rem(s,y,r)); }
     698             : 
     699             : INLINE long
     700          98 : sdivsi(long x, GEN y)
     701             : {
     702          98 :   long q, s = signe(y);
     703             : 
     704          98 :   if (!s) pari_err_INV("sdivsi",gen_0);
     705          98 :   if (!x || lgefint(y)>3 || ((long)y[2]) < 0) return 0;
     706          70 :   q = labs(x) / y[2];
     707          70 :   if (x < 0) q = -q;
     708          70 :   if (s < 0) q = -q;
     709          70 :   return q;
     710             : }
     711             : 
     712             : INLINE GEN
     713           0 : dvmdss(long x, long y, GEN *z)
     714             : {
     715             :   long r;
     716           0 :   GEN q = divss_rem(x,y, &r);
     717           0 :   *z = stoi(r); return q;
     718             : }
     719             : INLINE long
     720  3253613335 : dvmdsBIL(long n, long *r) { *r = remsBIL(n); return divsBIL(n); }
     721             : INLINE ulong
     722    99999664 : dvmduBIL(ulong n, ulong *r) { *r = remsBIL(n); return divsBIL(n); }
     723             : INLINE GEN
     724           0 : dvmdsi(long x, GEN y, GEN *z)
     725             : {
     726             :   long r;
     727           0 :   GEN q = divsi_rem(x,y, &r);
     728           0 :   *z = stoi(r); return q;
     729             : }
     730             : INLINE GEN
     731           0 : dvmdis(GEN x, long y, GEN *z)
     732             : {
     733             :   long r;
     734           0 :   GEN q = divis_rem(x,y, &r);
     735           0 :   *z = stoi(r); return q;
     736             : }
     737             : 
     738             : INLINE long
     739     1212457 : smodis(GEN x, long y)
     740             : {
     741     1212457 :   pari_sp av = avma;
     742             :   long r;
     743     1212457 :   (void)divis_rem(x,y, &r); avma = av; return (r >= 0) ? r: labs(y) + r;
     744             : }
     745             : INLINE GEN
     746       36439 : modis(GEN x, long y) { return stoi(smodis(x,y)); }
     747             : INLINE GEN
     748     7473309 : modsi(long x, GEN y) {
     749             :   long r;
     750     7473309 :   (void)sdivsi_rem(x, y, &r);
     751     7473312 :   return (r >= 0)? stoi(r): addsi_sign(r, y, 1);
     752             : }
     753             : 
     754             : INLINE ulong
     755     1192901 : umodui(ulong x, GEN y)
     756             : {
     757     1192901 :   if (!signe(y)) pari_err_INV("umodui",gen_0);
     758     1192901 :   if (!x || lgefint(y) > 3) return x;
     759      203752 :   return x % (ulong)y[2];
     760             : }
     761             : 
     762             : INLINE GEN
     763           0 : remsi(long x, GEN y)
     764           0 : { long r; (void)sdivsi_rem(x,y, &r); return stoi(r); }
     765             : INLINE GEN
     766           0 : remis(GEN x, long y)
     767             : {
     768           0 :   pari_sp av = avma;
     769             :   long r;
     770           0 :   (void)divis_rem(x,y, &r); avma = av; return stoi(r);
     771             : }
     772             : 
     773             : INLINE GEN
     774           0 : rdivis(GEN x, long y, long prec)
     775             : {
     776           0 :   GEN z = cgetr(prec);
     777           0 :   pari_sp av = avma;
     778           0 :   affrr(divrs(itor(x,prec), y),z);
     779           0 :   avma = av; return z;
     780             : }
     781             : INLINE GEN
     782           0 : rdivsi(long x, GEN y, long prec)
     783             : {
     784           0 :   GEN z = cgetr(prec);
     785           0 :   pari_sp av = avma;
     786           0 :   affrr(divsr(x, itor(y,prec)), z);
     787           0 :   avma = av; return z;
     788             : }
     789             : INLINE GEN
     790      839647 : rdivss(long x, long y, long prec)
     791             : {
     792      839647 :   GEN z = cgetr(prec);
     793      839647 :   pari_sp av = avma;
     794      839647 :   affrr(divrs(stor(x, prec), y), z);
     795      839647 :   avma = av; return z;
     796             : }
     797             : 
     798             : INLINE void
     799          14 : rdiviiz(GEN x, GEN y, GEN z)
     800             : {
     801          14 :   pari_sp av = avma;
     802          14 :   long prec = realprec(z);
     803          14 :   affir(x, z);
     804          14 :   if (!is_bigint(y)) {
     805           0 :     affrr(divrs(z, y[2]), z);
     806           0 :     if (signe(y) < 0) togglesign(z);
     807             :   }
     808             :   else
     809          14 :     affrr(divrr(z, itor(y,prec)), z);
     810          14 :   avma = av;
     811          14 : }
     812             : INLINE GEN
     813    11631526 : rdivii(GEN x, GEN y, long prec)
     814             : {
     815    11631526 :   GEN z = cgetr(prec);
     816    11631526 :   pari_sp av = avma;
     817    11631526 :   affir(x, z);
     818    11631526 :   if (lg(y) == 3) {
     819     7548726 :     affrr(divru(z, y[2]), z);
     820     7548726 :     if (signe(y) < 0) togglesign(z);
     821             :   }
     822             :   else
     823     4082800 :     affrr(divrr(z, itor(y,prec)), z);
     824    11631526 :   avma = av; return z;
     825             : }
     826             : INLINE GEN
     827    11578176 : fractor(GEN x, long prec) { return rdivii(gel(x,1), gel(x,2), prec); }
     828             : 
     829             : INLINE int
     830     1514761 : dvdii(GEN x, GEN y)
     831             : {
     832     1514761 :   pari_sp av=avma;
     833     1514761 :   GEN r = remii(x,y);
     834     1514761 :   avma = av; return r == gen_0;
     835             : }
     836             : INLINE int
     837           7 : dvdsi(long x, GEN y)
     838             : {
     839           7 :   if (!signe(y)) return x == 0;
     840           7 :   if (lgefint(y) != 3) return 0;
     841           7 :   return x % y[2] == 0;
     842             : }
     843             : INLINE int
     844         126 : dvdui(ulong x, GEN y)
     845             : {
     846         126 :   if (!signe(y)) return x == 0;
     847         126 :   if (lgefint(y) != 3) return 0;
     848         126 :   return x % y[2] == 0;
     849             : }
     850             : INLINE int
     851       31435 : dvdis(GEN x, long y)
     852       31435 : { return y? smodis(x, y) == 0: signe(x) == 0; }
     853             : INLINE int
     854        1540 : dvdiu(GEN x, ulong y)
     855        1540 : { return y? umodiu(x, y) == 0: signe(x) == 0; }
     856             : 
     857             : INLINE int
     858           0 : dvdisz(GEN x, long y, GEN z)
     859             : {
     860           0 :   const pari_sp av = avma;
     861             :   long r;
     862           0 :   GEN p1 = divis_rem(x,y, &r);
     863           0 :   avma = av; if (r) return 0;
     864           0 :   affii(p1,z); return 1;
     865             : }
     866             : INLINE int
     867           0 : dvdiuz(GEN x, ulong y, GEN z)
     868             : {
     869           0 :   const pari_sp av = avma;
     870             :   ulong r;
     871           0 :   GEN p1 = diviu_rem(x,y, &r);
     872           0 :   avma = av; if (r) return 0;
     873           0 :   affii(p1,z); return 1;
     874             : }
     875             : INLINE int
     876        3782 : dvdiiz(GEN x, GEN y, GEN z)
     877             : {
     878        3782 :   const pari_sp av=avma;
     879             :   GEN p2;
     880        3782 :   const GEN p1=dvmdii(x,y,&p2);
     881             : 
     882        3782 :   if (signe(p2)) { avma=av; return 0; }
     883        3282 :   affii(p1,z); avma=av; return 1;
     884             : }
     885             : 
     886             : INLINE ulong
     887     5680554 : remlll_pre(ulong u2, ulong u1, ulong u0, ulong n, ulong ninv)
     888             : {
     889     5680554 :   u1 = remll_pre(u2, u1, n, ninv);
     890     5680554 :   return remll_pre(u1, u0, n, ninv);
     891             : }
     892             : 
     893             : INLINE ulong
     894  4252009531 : Fl_sqr_pre(ulong a, ulong p, ulong pi)
     895             : {
     896             :   register ulong x;
     897             :   LOCAL_HIREMAINDER;
     898  4252009531 :   x = mulll(a,a);
     899  4252009531 :   return remll_pre(hiremainder, x, p, pi);
     900             : }
     901             : 
     902             : INLINE ulong
     903  2806447997 : Fl_mul_pre(ulong a, ulong b, ulong p, ulong pi)
     904             : {
     905             :   register ulong x;
     906             :   LOCAL_HIREMAINDER;
     907  2806447997 :   x = mulll(a,b);
     908  2806447997 :   return remll_pre(hiremainder, x, p, pi);
     909             : }
     910             : 
     911             : INLINE ulong
     912  4626936361 : Fl_addmul_pre(ulong x0, ulong x1, ulong y0, ulong p, ulong pi)
     913             : {
     914             :   ulong l0, h0;
     915             :   LOCAL_HIREMAINDER;
     916  4626936361 :   hiremainder = y0;
     917  4626936361 :   l0 = addmul(x0, x1); h0 = hiremainder;
     918  4626936361 :   return remll_pre(h0, l0, p, pi);
     919             : }
     920             : 
     921             : INLINE ulong
     922    45163262 : Fl_addmulmul_pre(ulong x0, ulong y0, ulong x1, ulong y1, ulong p, ulong pi)
     923             : {
     924             :   ulong l0, l1, h0, h1;
     925             :   LOCAL_OVERFLOW;
     926             :   LOCAL_HIREMAINDER;
     927    45163262 :   l0 = mulll(x0, y0); h0 = hiremainder;
     928    45163262 :   l1 = mulll(x1, y1); h1 = hiremainder;
     929    45163262 :   l0 = addll(l0, l1); h0 = addllx(h0, h1);
     930    45163262 :   return overflow ? remlll_pre(1, h0, l0, p, pi): remll_pre(h0, l0, p, pi);
     931             : }
     932             : 
     933             : INLINE ulong
     934     2706448 : Fl_ellj_pre(ulong a4, ulong a6, ulong p, ulong pi)
     935             : {
     936             :   /* a43 = 4 a4^3 */
     937     2706448 :   ulong a43 = Fl_double(Fl_double(
     938             :               Fl_mul_pre(a4, Fl_sqr_pre(a4, p, pi), p, pi), p), p);
     939             :   /* a62 = 27 a6^2 */
     940     2706449 :   ulong a62 = Fl_mul_pre(Fl_sqr_pre(a6, p, pi), 27 % p, p, pi);
     941     2706450 :   ulong z1 = Fl_mul_pre(a43, 1728 % p, p, pi);
     942     2706448 :   ulong z2 = Fl_add(a43, a62, p);
     943     2706447 :   return Fl_div(z1, z2, p);
     944             : }
     945             : 
     946             : /*******************************************************************/
     947             : /*                                                                 */
     948             : /*                        MP (INT OR REAL)                         */
     949             : /*                                                                 */
     950             : /*******************************************************************/
     951             : INLINE GEN
     952          21 : mptrunc(GEN x) { return typ(x)==t_INT? icopy(x): truncr(x); }
     953             : INLINE GEN
     954           0 : mpfloor(GEN x) { return typ(x)==t_INT? icopy(x): floorr(x); }
     955             : INLINE GEN
     956           0 : mpceil(GEN x) { return typ(x)==t_INT? icopy(x): ceilr(x); }
     957             : INLINE GEN
     958      357868 : mpround(GEN x) { return typ(x) == t_INT? icopy(x): roundr(x); }
     959             : 
     960             : INLINE long
     961      203279 : mpexpo(GEN x) { return typ(x) == t_INT? expi(x): expo(x); }
     962             : 
     963             : INLINE GEN
     964    19996359 : mpadd(GEN x, GEN y)
     965             : {
     966    19996359 :   if (typ(x)==t_INT)
     967     4097236 :     return (typ(y)==t_INT) ? addii(x,y) : addir(x,y);
     968    15899123 :   return (typ(y)==t_INT) ? addir(y,x) : addrr(x,y);
     969             : }
     970             : INLINE GEN
     971    10465484 : mpsub(GEN x, GEN y)
     972             : {
     973    10465484 :   if (typ(x)==t_INT)
     974      259051 :     return (typ(y)==t_INT) ? subii(x,y) : subir(x,y);
     975    10206433 :   return (typ(y)==t_INT) ? subri(x,y) : subrr(x,y);
     976             : }
     977             : INLINE GEN
     978    28020299 : mpmul(GEN x, GEN y)
     979             : {
     980    28020299 :   if (typ(x)==t_INT)
     981     4876851 :     return (typ(y)==t_INT) ? mulii(x,y) : mulir(x,y);
     982    23143448 :   return (typ(y)==t_INT) ? mulir(y,x) : mulrr(x,y);
     983             : }
     984             : INLINE GEN
     985     3523111 : mpsqr(GEN x) { return (typ(x)==t_INT) ? sqri(x) : sqrr(x); }
     986             : INLINE GEN
     987      195972 : mpdiv(GEN x, GEN y)
     988             : {
     989      195972 :   if (typ(x)==t_INT)
     990      101416 :     return (typ(y)==t_INT) ? divii(x,y) : divir(x,y);
     991       94556 :   return (typ(y)==t_INT) ? divri(x,y) : divrr(x,y);
     992             : }
     993             : 
     994             : /*******************************************************************/
     995             : /*                                                                 */
     996             : /*                          Z/nZ, n ULONG                          */
     997             : /*                                                                 */
     998             : /*******************************************************************/
     999             : INLINE ulong
    1000   643832164 : Fl_double(ulong a, ulong p)
    1001             : {
    1002   643832164 :   ulong res = a << 1;
    1003   643832164 :   return (res >= p || res < a) ? res - p : res;
    1004             : }
    1005             : INLINE ulong
    1006   102856226 : Fl_triple(ulong a, ulong p)
    1007             : {
    1008   102856226 :   ulong res = a << 1;
    1009   102856226 :   if (res >= p || res < a) res -= p;
    1010   102856226 :   res += a;
    1011   102856226 :   return (res >= p || res < a)? res - p: res;
    1012             : }
    1013             : INLINE ulong
    1014    12750284 : Fl_halve(ulong a, ulong p)
    1015             : {
    1016             :   ulong ap, ap2;
    1017    12750284 :   if ((a&1UL)==0) return a>>1;
    1018     6437732 :   ap = a + p; ap2 = ap>>1;
    1019     6437732 :   return ap>=a ? ap2: (ap2|HIGHBIT);
    1020             : }
    1021             : 
    1022             : INLINE ulong
    1023  3872550318 : Fl_add(ulong a, ulong b, ulong p)
    1024             : {
    1025  3872550318 :   ulong res = a + b;
    1026  3872550318 :   return (res >= p || res < a) ? res - p : res;
    1027             : }
    1028             : INLINE ulong
    1029    70638500 : Fl_neg(ulong x, ulong p) { return x ? p - x: 0; }
    1030             : 
    1031             : INLINE ulong
    1032  4335009929 : Fl_sub(ulong a, ulong b, ulong p)
    1033             : {
    1034  4335009929 :   ulong res = a - b;
    1035  4335009929 :   return (res > a) ? res + p: res;
    1036             : }
    1037             : 
    1038             : /* centerlift(u mod p) */
    1039             : INLINE long
    1040    18540345 : Fl_center(ulong u, ulong p, ulong ps2) { return (long) (u > ps2)? u - p: u; }
    1041             : 
    1042             : INLINE ulong
    1043  1977400263 : Fl_mul(ulong a, ulong b, ulong p)
    1044             : {
    1045             :   register ulong x;
    1046             :   LOCAL_HIREMAINDER;
    1047  1977400263 :   x = mulll(a,b);
    1048  1977400263 :   if (!hiremainder) return x % p;
    1049  1001195958 :   (void)divll(x,p); return hiremainder;
    1050             : }
    1051             : INLINE ulong
    1052   124157060 : Fl_sqr(ulong a, ulong p)
    1053             : {
    1054             :   register ulong x;
    1055             :   LOCAL_HIREMAINDER;
    1056   124157060 :   x = mulll(a,a);
    1057   124157060 :   if (!hiremainder) return x % p;
    1058    47203705 :   (void)divll(x,p); return hiremainder;
    1059             : }
    1060             : INLINE ulong
    1061    15330020 : Fl_div(ulong a, ulong b, ulong p) { return Fl_mul(a, Fl_inv(b, p), p); }
    1062             : 
    1063             : /*******************************************************************/
    1064             : /*                                                                 */
    1065             : /*        DEFINED FROM EXISTING ONE EXPLOITING COMMUTATIVITY       */
    1066             : /*                                                                 */
    1067             : /*******************************************************************/
    1068             : INLINE GEN
    1069      262135 : addri(GEN x, GEN y) { return addir(y,x); }
    1070             : INLINE GEN
    1071    65983642 : addis(GEN x, long s) { return addsi(s,x); }
    1072             : INLINE GEN
    1073    11674979 : addiu(GEN x, ulong s) { return addui(s,x); }
    1074             : INLINE GEN
    1075     2046912 : addrs(GEN x, long s) { return addsr(s,x); }
    1076             : 
    1077             : INLINE GEN
    1078    14114262 : subiu(GEN x, long y) { GEN z = subui(y, x); togglesign(z); return z; }
    1079             : INLINE GEN
    1080    16528225 : subis(GEN x, long y) { return addsi(-y,x); }
    1081             : INLINE GEN
    1082     2254940 : subrs(GEN x, long y) { return addsr(-y,x); }
    1083             : 
    1084             : INLINE GEN
    1085   114501299 : mulis(GEN x, long s) { return mulsi(s,x); }
    1086             : INLINE GEN
    1087   555327115 : muliu(GEN x, ulong s) { return mului(s,x); }
    1088             : INLINE GEN
    1089     3937318 : mulru(GEN x, ulong s) { return mulur(s,x); }
    1090             : INLINE GEN
    1091    14004329 : mulri(GEN x, GEN s) { return mulir(s,x); }
    1092             : INLINE GEN
    1093    10273157 : mulrs(GEN x, long s) { return mulsr(s,x); }
    1094             : 
    1095             : /*******************************************************************/
    1096             : /*                                                                 */
    1097             : /*                  VALUATION, EXPONENT, SHIFTS                    */
    1098             : /*                                                                 */
    1099             : /*******************************************************************/
    1100             : INLINE long
    1101    66630566 : vali(GEN x)
    1102             : {
    1103             :   long i;
    1104             :   GEN xp;
    1105             : 
    1106    66630566 :   if (!signe(x)) return -1;
    1107    66630566 :   xp=int_LSW(x);
    1108    66630566 :   for (i=0; !*xp; i++) xp=int_nextW(xp);
    1109    66630566 :   return vals(*xp) + i * BITS_IN_LONG;
    1110             : }
    1111             : 
    1112             : 
    1113             : /* assume x > 0 */
    1114             : INLINE long
    1115   302801300 : expu(ulong x) { return (BITS_IN_LONG-1) - (long)bfffo(x); }
    1116             : 
    1117             : INLINE long
    1118   776088705 : expi(GEN x)
    1119             : {
    1120   776088705 :   const long lx=lgefint(x);
    1121   776088705 :   return lx==2? -(long)HIGHEXPOBIT: bit_accuracy(lx)-(long)bfffo(*int_MSW(x))-1;
    1122             : }
    1123             : 
    1124             : INLINE GEN
    1125    46578599 : shiftr(GEN x, long n)
    1126             : {
    1127    46578599 :   const long e = evalexpo(expo(x)+n);
    1128    46578599 :   const GEN y = rcopy(x);
    1129             : 
    1130    46578599 :   if (e & ~EXPOBITS) pari_err_OVERFLOW("expo()");
    1131    46578599 :   y[1] = (y[1]&~EXPOBITS) | e; return y;
    1132             : }
    1133             : INLINE GEN
    1134    10856999 : mpshift(GEN x,long s) { return (typ(x)==t_INT)?shifti(x,s):shiftr(x,s); }
    1135             : 
    1136             : /* FIXME: adapt/use mpn_[lr]shift instead */
    1137             : /* z2[imin..imax] := z1[imin..imax].f shifted left sh bits
    1138             :  * (feeding f from the right). Assume sh > 0 */
    1139             : INLINE void
    1140  2808934254 : shift_left(GEN z2, GEN z1, long imin, long imax, ulong f,  ulong sh)
    1141             : {
    1142  2808934254 :   GEN sb = z1 + imin, se = z1 + imax, te = z2 + imax;
    1143  2808934254 :   ulong l, m = BITS_IN_LONG - sh, k = f >> m;
    1144 13422891401 :   while (se > sb) {
    1145  7805022893 :     l     = *se--;
    1146  7805022893 :     *te-- = (l << sh) | k;
    1147  7805022893 :     k     = l >> m;
    1148             :   }
    1149  2808934254 :   *te = (((ulong)*se) << sh) | k;
    1150  2808934254 : }
    1151             : /* z2[imin..imax] := f.z1[imin..imax-1] shifted right sh bits
    1152             :  * (feeding f from the left). Assume sh > 0 */
    1153             : INLINE void
    1154  3315829160 : shift_right(GEN z2, GEN z1, long imin, long imax, ulong f, ulong sh)
    1155             : {
    1156  3315829160 :   GEN sb = z1 + imin, se = z1 + imax, tb = z2 + imin;
    1157  3315829160 :   ulong k, l = *sb++, m = BITS_IN_LONG - sh;
    1158  3315829160 :   *tb++ = (l >> sh) | (f << m);
    1159 14675171007 :   while (sb < se) {
    1160  8043512687 :     k     = l << m;
    1161  8043512687 :     l     = *sb++;
    1162  8043512687 :     *tb++ = (l >> sh) | k;
    1163             :   }
    1164  3315829160 : }
    1165             : 
    1166             : /* Backward compatibility. Inefficient && unused */
    1167             : extern ulong hiremainder;
    1168             : INLINE ulong
    1169           0 : shiftl(ulong x, ulong y)
    1170           0 : { hiremainder = x>>(BITS_IN_LONG-y); return (x<<y); }
    1171             : 
    1172             : INLINE ulong
    1173           0 : shiftlr(ulong x, ulong y)
    1174           0 : { hiremainder = x<<(BITS_IN_LONG-y); return (x>>y); }
    1175             : 
    1176             : INLINE void
    1177    76668669 : shiftr_inplace(GEN z, long d)
    1178             : {
    1179    76668669 :   setexpo(z, expo(z)+d);
    1180    76668669 : }
    1181             : 
    1182             : /*******************************************************************/
    1183             : /*                                                                 */
    1184             : /*                           ASSIGNMENT                            */
    1185             : /*                                                                 */
    1186             : /*******************************************************************/
    1187             : INLINE void
    1188    91280662 : affii(GEN x, GEN y)
    1189             : {
    1190    91280662 :   long lx = lgefint(x);
    1191    91280662 :   if (lg(y)<lx) pari_err_OVERFLOW("t_INT-->t_INT assignment");
    1192    91280645 :   while (--lx) y[lx] = x[lx];
    1193    91280645 : }
    1194             : INLINE void
    1195      264927 : affsi(long s, GEN x)
    1196             : {
    1197      264927 :   if (!s) x[1] = evalsigne(0) | evallgefint(2);
    1198             :   else
    1199             :   {
    1200      250660 :     if (s > 0) { x[1] = evalsigne( 1) | evallgefint(3); x[2] =  s; }
    1201       90158 :     else       { x[1] = evalsigne(-1) | evallgefint(3); x[2] = -s; }
    1202             :   }
    1203      264927 : }
    1204             : INLINE void
    1205     8398042 : affui(ulong u, GEN x)
    1206             : {
    1207     8398042 :   if (!u) x[1] = evalsigne(0) | evallgefint(2);
    1208     8396838 :   else  { x[1] = evalsigne(1) | evallgefint(3); x[2] = u; }
    1209     8398042 : }
    1210             : 
    1211             : INLINE void
    1212   144459282 : affsr(long x, GEN y)
    1213             : {
    1214   144459282 :   long sh, i, ly = lg(y);
    1215             : 
    1216   144459282 :   if (!x)
    1217             :   {
    1218     1366432 :     y[1] = evalexpo(-prec2nbits(ly));
    1219   145825714 :     return;
    1220             :   }
    1221   143092850 :   if (x < 0) {
    1222      162332 :     x = -x; sh = bfffo(x);
    1223      162332 :     y[1] = evalsigne(-1) | _evalexpo((BITS_IN_LONG-1)-sh);
    1224             :   }
    1225             :   else
    1226             :   {
    1227   142930518 :     sh = bfffo(x);
    1228   142930518 :     y[1] = evalsigne(1) | _evalexpo((BITS_IN_LONG-1)-sh);
    1229             :   }
    1230   143092850 :   y[2] = ((ulong)x)<<sh; for (i=3; i<ly; i++) y[i]=0;
    1231             : }
    1232             : 
    1233             : INLINE void
    1234     1816834 : affur(ulong x, GEN y)
    1235             : {
    1236     1816834 :   long sh, i, ly = lg(y);
    1237             : 
    1238     1816834 :   if (!x)
    1239             :   {
    1240         210 :     y[1] = evalexpo(-prec2nbits(ly));
    1241     1817044 :     return;
    1242             :   }
    1243     1816624 :   sh = bfffo(x);
    1244     1816624 :   y[1] = evalsigne(1) | _evalexpo((BITS_IN_LONG-1)-sh);
    1245     1816624 :   y[2] = x<<sh; for (i=3; i<ly; i++) y[i] = 0;
    1246             : }
    1247             : 
    1248             : INLINE void
    1249       95113 : affiz(GEN x, GEN y) { if (typ(y)==t_INT) affii(x,y); else affir(x,y); }
    1250             : INLINE void
    1251           0 : affsz(long x, GEN y) { if (typ(y)==t_INT) affsi(x,y); else affsr(x,y); }
    1252             : INLINE void
    1253       99322 : mpaff(GEN x, GEN y) { if (typ(x)==t_INT) affiz(x, y); else affrr(x,y); }
    1254             : 
    1255             : /*******************************************************************/
    1256             : /*                                                                 */
    1257             : /*                    OPERATION + ASSIGNMENT                       */
    1258             : /*                                                                 */
    1259             : /*******************************************************************/
    1260             : 
    1261           0 : INLINE void addiiz(GEN x, GEN y, GEN z)
    1262           0 : { pari_sp av = avma; affii(addii(x,y),z); avma = av; }
    1263           0 : INLINE void addirz(GEN x, GEN y, GEN z)
    1264           0 : { pari_sp av = avma; affrr(addir(x,y),z); avma = av; }
    1265           0 : INLINE void addriz(GEN x, GEN y, GEN z)
    1266           0 : { pari_sp av = avma; affrr(addri(x,y),z); avma = av; }
    1267     1436864 : INLINE void addrrz(GEN x, GEN y, GEN z)
    1268     1436864 : { pari_sp av = avma; affrr(addrr(x,y),z); avma = av; }
    1269           0 : INLINE void addsiz(long s, GEN y, GEN z)
    1270           0 : { pari_sp av = avma; affii(addsi(s,y),z); avma = av; }
    1271           0 : INLINE void addsrz(long s, GEN y, GEN z)
    1272           0 : { pari_sp av = avma; affrr(addsr(s,y),z); avma = av; }
    1273           0 : INLINE void addssz(long s, long y, GEN z)
    1274           0 : { pari_sp av = avma; affii(addss(s,y),z); avma = av; }
    1275             : 
    1276          35 : INLINE void diviiz(GEN x, GEN y, GEN z)
    1277          35 : { pari_sp av = avma; affii(divii(x,y),z); avma = av; }
    1278           0 : INLINE void divirz(GEN x, GEN y, GEN z)
    1279           0 : { pari_sp av = avma; mpaff(divir(x,y),z); avma = av; }
    1280           0 : INLINE void divisz(GEN x, long y, GEN z)
    1281           0 : { pari_sp av = avma; affii(divis(x,y),z); avma = av; }
    1282           0 : INLINE void divriz(GEN x, GEN y, GEN z)
    1283           0 : { pari_sp av = avma; affrr(divri(x,y),z); avma = av; }
    1284         460 : INLINE void divrrz(GEN x, GEN y, GEN z)
    1285         460 : { pari_sp av = avma; affrr(divrr(x,y),z); avma = av; }
    1286       90548 : INLINE void divrsz(GEN y, long s, GEN z)
    1287       90548 : { pari_sp av = avma; affrr(divrs(y,s),z); avma = av; }
    1288           0 : INLINE void divsiz(long x, GEN y, GEN z)
    1289           0 : { long junk; affsi(sdivsi_rem(x,y,&junk), z); }
    1290           0 : INLINE void divsrz(long s, GEN y, GEN z)
    1291           0 : { pari_sp av = avma; mpaff(divsr(s,y),z); avma = av; }
    1292           0 : INLINE void divssz(long x, long y, GEN z)
    1293           0 : { affsi(x/y, z); }
    1294             : 
    1295           0 : INLINE void modisz(GEN y, long s, GEN z)
    1296           0 : { pari_sp av = avma; affii(modis(y,s),z); avma = av; }
    1297           0 : INLINE void modsiz(long s, GEN y, GEN z)
    1298           0 : { pari_sp av = avma; affii(modsi(s,y),z); avma = av; }
    1299           0 : INLINE void modssz(long s, long y, GEN z)
    1300           0 : { pari_sp av = avma; affii(modss(s,y),z); avma = av; }
    1301             : 
    1302           0 : INLINE void mpaddz(GEN x, GEN y, GEN z)
    1303           0 : { pari_sp av = avma; mpaff(mpadd(x,y),z); avma = av; }
    1304           0 : INLINE void mpsubz(GEN x, GEN y, GEN z)
    1305           0 : { pari_sp av = avma; mpaff(mpsub(x,y),z); avma = av; }
    1306           0 : INLINE void mpmulz(GEN x, GEN y, GEN z)
    1307           0 : { pari_sp av = avma; mpaff(mpmul(x,y),z); avma = av; }
    1308             : 
    1309           0 : INLINE void muliiz(GEN x, GEN y, GEN z)
    1310           0 : { pari_sp av = avma; affii(mulii(x,y),z); avma = av; }
    1311           0 : INLINE void mulirz(GEN x, GEN y, GEN z)
    1312           0 : { pari_sp av = avma; mpaff(mulir(x,y),z); avma = av; }
    1313           0 : INLINE void mulriz(GEN x, GEN y, GEN z)
    1314           0 : { pari_sp av = avma; mpaff(mulri(x,y),z); avma = av; }
    1315       56940 : INLINE void mulrrz(GEN x, GEN y, GEN z)
    1316       56940 : { pari_sp av = avma; affrr(mulrr(x,y),z); avma = av; }
    1317           0 : INLINE void mulsiz(long s, GEN y, GEN z)
    1318           0 : { pari_sp av = avma; affii(mulsi(s,y),z); avma = av; }
    1319           0 : INLINE void mulsrz(long s, GEN y, GEN z)
    1320           0 : { pari_sp av = avma; mpaff(mulsr(s,y),z); avma = av; }
    1321           0 : INLINE void mulssz(long s, long y, GEN z)
    1322           0 : { pari_sp av = avma; affii(mulss(s,y),z); avma = av; }
    1323             : 
    1324           0 : INLINE void remiiz(GEN x, GEN y, GEN z)
    1325           0 : { pari_sp av = avma; affii(remii(x,y),z); avma = av; }
    1326           0 : INLINE void remisz(GEN y, long s, GEN z)
    1327           0 : { pari_sp av = avma; affii(remis(y,s),z); avma = av; }
    1328           0 : INLINE void remsiz(long s, GEN y, GEN z)
    1329           0 : { pari_sp av = avma; affii(remsi(s,y),z); avma = av; }
    1330           0 : INLINE void remssz(long s, long y, GEN z)
    1331           0 : { pari_sp av = avma; affii(remss(s,y),z); avma = av; }
    1332             : 
    1333           0 : INLINE void subiiz(GEN x, GEN y, GEN z)
    1334           0 : { pari_sp av = avma; affii(subii(x,y),z); avma = av; }
    1335           0 : INLINE void subirz(GEN x, GEN y, GEN z)
    1336           0 : { pari_sp av = avma; affrr(subir(x,y),z); avma = av; }
    1337           0 : INLINE void subisz(GEN y, long s, GEN z)
    1338           0 : { pari_sp av = avma; affii(addsi(-s,y),z); avma = av; }
    1339           0 : INLINE void subriz(GEN x, GEN y, GEN z)
    1340           0 : { pari_sp av = avma; affrr(subri(x,y),z); avma = av; }
    1341     1295642 : INLINE void subrrz(GEN x, GEN y, GEN z)
    1342     1295642 : { pari_sp av = avma; affrr(subrr(x,y),z); avma = av; }
    1343           0 : INLINE void subrsz(GEN y, long s, GEN z)
    1344           0 : { pari_sp av = avma; affrr(addsr(-s,y),z); avma = av; }
    1345           0 : INLINE void subsiz(long s, GEN y, GEN z)
    1346           0 : { pari_sp av = avma; affii(subsi(s,y),z); avma = av; }
    1347           0 : INLINE void subsrz(long s, GEN y, GEN z)
    1348           0 : { pari_sp av = avma; affrr(subsr(s,y),z); avma = av; }
    1349           0 : INLINE void subssz(long x, long y, GEN z) { addssz(x,-y,z); }
    1350             : 
    1351             : INLINE void
    1352           0 : dvmdssz(long x, long y, GEN z, GEN t) {
    1353           0 :   pari_sp av = avma;
    1354             :   long r;
    1355           0 :   affii(divss_rem(x,y, &r), z); avma = av; affsi(r,t);
    1356           0 : }
    1357             : INLINE void
    1358           0 : dvmdsiz(long x, GEN y, GEN z, GEN t) {
    1359           0 :   pari_sp av = avma;
    1360             :   long r;
    1361           0 :   affii(divsi_rem(x,y, &r), z); avma = av; affsi(r,t);
    1362           0 : }
    1363             : INLINE void
    1364           0 : dvmdisz(GEN x, long y, GEN z, GEN t) {
    1365           0 :   pari_sp av = avma;
    1366             :   long r;
    1367           0 :   affii(divis_rem(x,y, &r),z); avma = av; affsi(r,t);
    1368           0 : }
    1369             : INLINE void
    1370           0 : dvmdiiz(GEN x, GEN y, GEN z, GEN t) {
    1371           0 :   pari_sp av = avma;
    1372             :   GEN r;
    1373           0 :   affii(dvmdii(x,y,&r),z); affii(r,t); avma=av;
    1374           0 : }

Generated by: LCOV version 1.11