Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - language - init.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.18.0 lcov report (development 29712-7c8a932571) Lines: 1148 1493 76.9 %
Date: 2024-11-15 09:08:45 Functions: 135 162 83.3 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000-2003  The PARI group.
       2             : 
       3             : This file is part of the PARI/GP package.
       4             : 
       5             : PARI/GP is free software; you can redistribute it and/or modify it under the
       6             : terms of the GNU General Public License as published by the Free Software
       7             : Foundation; either version 2 of the License, or (at your option) any later
       8             : version. It is distributed in the hope that it will be useful, but WITHOUT
       9             : ANY WARRANTY WHATSOEVER.
      10             : 
      11             : Check the License for details. You should have received a copy of it, along
      12             : with the package; see the file 'COPYING'. If not, write to the Free Software
      13             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      14             : 
      15             : /*******************************************************************/
      16             : /*                                                                 */
      17             : /*        INITIALIZING THE SYSTEM, ERRORS, STACK MANAGEMENT        */
      18             : /*                                                                 */
      19             : /*******************************************************************/
      20             : /* _GNU_SOURCE is needed before first include to get RUSAGE_THREAD */
      21             : #undef _GNU_SOURCE /* avoid warning */
      22             : #define _GNU_SOURCE
      23             : #include <string.h>
      24             : #if defined(_WIN32) || defined(__CYGWIN32__)
      25             : #  include "../systems/mingw/mingw.h"
      26             : #  include <process.h>
      27             : #endif
      28             : #include "paricfg.h"
      29             : #if defined(STACK_CHECK) && !defined(__EMX__) && !defined(_WIN32)
      30             : #  include <sys/types.h>
      31             : #  include <sys/time.h>
      32             : #  include <sys/resource.h>
      33             : #endif
      34             : #if defined(HAS_WAITPID) && defined(HAS_SETSID)
      35             : #  include <sys/wait.h>
      36             : #endif
      37             : #ifdef HAS_MMAP
      38             : #  include <sys/mman.h>
      39             : #endif
      40             : #if defined(USE_GETTIMEOFDAY) || defined(USE_GETRUSAGE) || defined(USE_TIMES)
      41             : #  include <sys/time.h>
      42             : #endif
      43             : #if defined(USE_GETRUSAGE)
      44             : #  include <sys/resource.h>
      45             : #endif
      46             : #if defined(USE_FTIME) || defined(USE_FTIMEFORWALLTIME)
      47             : #  include <sys/timeb.h>
      48             : #endif
      49             : #if defined(USE_CLOCK_GETTIME) || defined(USE_TIMES)
      50             : #  include <time.h>
      51             : #endif
      52             : #if defined(USE_TIMES)
      53             : #  include <sys/times.h>
      54             : #endif
      55             : #define PARI_INIT
      56             : #include "pari.h"
      57             : #include "paripriv.h"
      58             : #include "anal.h"
      59             : 
      60             : const double LOG10_2 = 0.3010299956639812; /* log_10(2) */
      61             : const double LOG2_10 = 3.321928094887362;  /* log_2(10) */
      62             : 
      63             : GEN gnil, gen_0, gen_1, gen_m1, gen_2, gen_m2, ghalf, err_e_STACK;
      64             : 
      65             : static const ulong readonly_constants[] = {
      66             :   evaltyp(t_INT) | _evallg(2),  /* gen_0 */
      67             :   evallgefint(2),
      68             :   evaltyp(t_INT) | _evallg(2),  /* gnil */
      69             :   evallgefint(2),
      70             :   evaltyp(t_INT) | _evallg(3),  /* gen_1 */
      71             :   evalsigne(1) | evallgefint(3),
      72             :   1,
      73             :   evaltyp(t_INT) | _evallg(3),  /* gen_2 */
      74             :   evalsigne(1) | evallgefint(3),
      75             :   2,
      76             :   evaltyp(t_INT) | _evallg(3),  /* gen_m1 */
      77             :   evalsigne(-1) | evallgefint(3),
      78             :   1,
      79             :   evaltyp(t_INT) | _evallg(3),  /* gen_m2 */
      80             :   evalsigne(-1) | evallgefint(3),
      81             :   2,
      82             :   evaltyp(t_ERROR) | _evallg(2), /* err_e_STACK */
      83             :   e_STACK,
      84             :   evaltyp(t_FRAC) | _evallg(3), /* ghalf */
      85             :   (ulong)(readonly_constants+4),
      86             :   (ulong)(readonly_constants+7)
      87             : };
      88             : THREAD GEN zetazone, bernzone, eulerzone, primetab;
      89             : pari_prime *pari_PRIMES;
      90             : FILE    *pari_outfile, *pari_errfile, *pari_logfile, *pari_infile;
      91             : char    *current_logfile, *current_psfile, *pari_datadir;
      92             : long    gp_colors[c_LAST];
      93             : int     disable_color;
      94             : ulong   DEBUGLEVEL, DEBUGMEM;
      95             : THREAD  long    DEBUGVAR;
      96             : ulong   pari_mt_nbthreads;
      97             : long    precreal;
      98             : ulong   precdl, pari_logstyle;
      99             : gp_data *GP_DATA;
     100             : 
     101             : entree  **varentries;
     102             : THREAD long *varpriority;
     103             : 
     104             : THREAD pari_sp avma;
     105             : THREAD struct pari_mainstack *pari_mainstack;
     106             : 
     107             : static void ** MODULES;
     108             : static pari_stack s_MODULES;
     109             : const long functions_tblsz = 135; /* size of functions_hash */
     110             : entree **functions_hash, **defaults_hash;
     111             : 
     112             : void (*cb_pari_display_hist)(long n);
     113             : char *(*cb_pari_fgets_interactive)(char *s, int n, FILE *f);
     114             : int (*cb_pari_get_line_interactive)(const char*, const char*, filtre_t *F);
     115             : void (*cb_pari_quit)(long);
     116             : void (*cb_pari_init_histfile)(void);
     117             : void (*cb_pari_ask_confirm)(const char *);
     118             : int  (*cb_pari_handle_exception)(long);
     119             : int  (*cb_pari_err_handle)(GEN);
     120             : int  (*cb_pari_whatnow)(PariOUT *out, const char *, int);
     121             : void (*cb_pari_sigint)(void);
     122             : void (*cb_pari_pre_recover)(long);
     123             : void (*cb_pari_err_recover)(long);
     124             : int (*cb_pari_break_loop)(int);
     125             : int (*cb_pari_is_interactive)(void);
     126             : void (*cb_pari_start_output)(void);
     127             : void (*cb_pari_long_help)(const char *s, long num);
     128             : 
     129             : const char * pari_library_path = NULL;
     130             : 
     131             : static THREAD GEN global_err_data;
     132             : THREAD jmp_buf *iferr_env;
     133             : const long CATCH_ALL = -1;
     134             : 
     135             : static void pari_init_timer(void);
     136             : 
     137             : /*********************************************************************/
     138             : /*                                                                   */
     139             : /*                       BLOCKS & CLONES                             */
     140             : /*                                                                   */
     141             : /*********************************************************************/
     142             : /*#define DEBUG*/
     143             : static THREAD long next_block;
     144             : static THREAD GEN cur_block; /* current block in block list */
     145             : static THREAD GEN root_block; /* current block in block list */
     146             : 
     147             : static void
     148      326354 : pari_init_blocks(void)
     149             : {
     150      326354 :   next_block = 0; cur_block = NULL; root_block = NULL;
     151      326354 : }
     152             : 
     153             : static void
     154      320311 : pari_close_blocks(void)
     155             : {
     156     2177380 :   while (cur_block) killblock(cur_block);
     157      324378 : }
     158             : 
     159             : static long
     160 11430702463 : blockheight(GEN bl) { return bl? bl_height(bl): 0; }
     161             : 
     162             : static long
     163  2743616930 : blockbalance(GEN bl)
     164  2743616930 : { return bl ? blockheight(bl_left(bl)) - blockheight(bl_right(bl)): 0; }
     165             : 
     166             : static void
     167  2971793366 : fix_height(GEN bl)
     168  2971793366 : { bl_height(bl) = maxss(blockheight(bl_left(bl)), blockheight(bl_right(bl)))+1; }
     169             : 
     170             : static GEN
     171    60657891 : bl_rotright(GEN y)
     172             : {
     173    60657891 :   GEN x = bl_left(y), t = bl_right(x);
     174    60657891 :   bl_right(x) = y;
     175    60657891 :   bl_left(y)  = t;
     176    60657891 :   fix_height(y);
     177    60657881 :   fix_height(x);
     178    60657886 :   return x;
     179             : }
     180             : 
     181             : static GEN
     182    63730020 : bl_rotleft(GEN x)
     183             : {
     184    63730020 :   GEN y = bl_right(x), t = bl_left(y);
     185    63730020 :   bl_left(y)  = x;
     186    63730020 :   bl_right(x) = t;
     187    63730020 :   fix_height(x);
     188    63730237 :   fix_height(y);
     189    63730285 :   return y;
     190             : }
     191             : 
     192             : static GEN
     193  1664463399 : blockinsert(GEN x, GEN bl, long *d)
     194             : {
     195             :   long b, c;
     196  1664463399 :   if (!bl)
     197             :   {
     198   232337922 :     bl_left(x)=NULL; bl_right(x)=NULL;
     199   232337922 :     bl_height(x)=1; return x;
     200             :   }
     201  1432125477 :   c = cmpuu((ulong)x, (ulong)bl);
     202  1432128114 :   if (c < 0)
     203   579126362 :     bl_left(bl) = blockinsert(x, bl_left(bl), d);
     204   853001752 :   else if (c > 0)
     205   853001752 :     bl_right(bl) = blockinsert(x, bl_right(bl), d);
     206           0 :   else return bl; /* ??? Already exist in the tree ? */
     207  1432117989 :   fix_height(bl);
     208  1432095618 :   b = blockbalance(bl);
     209  1432105474 :   if (b > 1)
     210             :   {
     211    32659485 :     if (*d > 0) bl_left(bl) = bl_rotleft(bl_left(bl));
     212    32659516 :     return bl_rotright(bl);
     213             :   }
     214  1399445989 :   if (b < -1)
     215             :   {
     216    30996877 :     if (*d < 0) bl_right(bl) = bl_rotright(bl_right(bl));
     217    30996880 :     return bl_rotleft(bl);
     218             :   }
     219  1368449112 :   *d = c; return bl;
     220             : }
     221             : 
     222             : static GEN
     223  1523404897 : blockdelete(GEN x, GEN bl)
     224             : {
     225             :   long b;
     226  1523404897 :   if (!bl) return NULL; /* ??? Do not exist in the tree */
     227  1523404897 :   if (x < bl)
     228   533393174 :     bl_left(bl) = blockdelete(x, bl_left(bl));
     229   990011723 :   else if (x > bl)
     230   705554262 :     bl_right(bl) = blockdelete(x, bl_right(bl));
     231             :   else
     232             :   {
     233   284457461 :     if (!bl_left(bl) && !bl_right(bl)) return NULL;
     234    95091145 :     else if (!bl_left(bl)) return bl_right(bl);
     235    73774214 :     else if (!bl_right(bl)) return bl_left(bl);
     236             :     else
     237             :     {
     238    52120606 :       GEN r = bl_right(bl);
     239    75120058 :       while (bl_left(r)) r = bl_left(r);
     240    52120606 :       bl_right(r) = blockdelete(r, bl_right(bl));
     241    52125194 :       bl_left(r) = bl_left(bl);
     242    52125194 :       bl = r;
     243             :     }
     244             :   }
     245  1291049445 :   fix_height(bl);
     246  1291027572 :   b = blockbalance(bl);
     247  1291025484 :   if (b > 1)
     248             :   {
     249    14583146 :     if (blockbalance(bl_left(bl)) >= 0) return bl_rotright(bl);
     250             :     else
     251     2102062 :     { bl_left(bl) = bl_rotleft(bl_left(bl)); return bl_rotright(bl); }
     252             :   }
     253  1276442338 :   if (b < -1)
     254             :   {
     255     5991345 :     if (blockbalance(bl_right(bl)) <= 0) return bl_rotleft(bl);
     256             :     else
     257     1946010 :     { bl_right(bl) = bl_rotright(bl_right(bl)); return bl_rotleft(bl); }
     258             :   }
     259  1270450993 :   return bl;
     260             : }
     261             : 
     262             : static GEN
     263   753181047 : blocksearch(GEN x, GEN bl)
     264             : {
     265   753181047 :   if (isclone(x)) return x;
     266   578711982 :   if (isonstack(x) || is_universal_constant(x)) return NULL;
     267   895202375 :   while (bl)
     268             :   {
     269   892286026 :     if (x >= bl  && x < bl + bl_size(bl))
     270   229683674 :       return bl;
     271   662602352 :     bl = x < bl ? bl_left(bl): bl_right(bl);
     272             :   }
     273     2916349 :   return NULL; /* Unknown address */
     274             : }
     275             : 
     276             : static int
     277     1748362 : check_clone(GEN x)
     278             : {
     279     1748362 :   GEN bl = root_block;
     280     1748362 :   if (isonstack(x) || is_universal_constant(x)) return 1;
     281     2368007 :   while (bl)
     282             :   {
     283     2368007 :     if (x >= bl  && x < bl + bl_size(bl))
     284      276310 :       return 1;
     285     2091697 :     bl = x < bl ? bl_left(bl): bl_right(bl);
     286             :   }
     287           0 :   return 0; /* Unknown address */
     288             : }
     289             : 
     290             : void
     291   376784206 : clone_lock(GEN x)
     292             : {
     293   376784206 :   GEN y = blocksearch(x, root_block);
     294   376461002 :   if (y && isclone(y))
     295             :   {
     296   201853051 :     if (DEBUGMEM > 2)
     297           0 :       err_printf("locking block no %ld: %08lx from %08lx\n", bl_num(y), y, x);
     298   201853051 :     ++bl_refc(y);
     299             :   }
     300   376461002 : }
     301             : 
     302             : void
     303   317485352 : clone_unlock(GEN x)
     304             : {
     305   317485352 :   GEN y = blocksearch(x, root_block);
     306   317418075 :   if (y && isclone(y))
     307             :   {
     308   150172206 :     if (DEBUGMEM > 2)
     309           0 :       err_printf("unlocking block no %ld: %08lx from %08lx\n", bl_num(y), y, x);
     310   150172206 :     gunclone(y);
     311             :   }
     312   317418075 : }
     313             : 
     314             : void
     315    59316004 : clone_unlock_deep(GEN x)
     316             : {
     317    59316004 :   GEN y = blocksearch(x, root_block);
     318    59316004 :   if (y && isclone(y))
     319             :   {
     320    52125730 :     if (DEBUGMEM > 2)
     321           0 :       err_printf("unlocking deep block no %ld: %08lx from %08lx\n", bl_num(y), y, x);
     322    52125730 :     gunclone_deep(y);
     323             :   }
     324    59316004 : }
     325             : 
     326             : /* Return x, where:
     327             :  * x[-8]: AVL height
     328             :  * x[-7]: adress of left child or NULL
     329             :  * x[-6]: adress of right child or NULL
     330             :  * x[-5]: size
     331             :  * x[-4]: reference count
     332             :  * x[-3]: adress of next block
     333             :  * x[-2]: adress of preceding block.
     334             :  * x[-1]: number of allocated blocs.
     335             :  * x[0..n-1]: malloc-ed memory. */
     336             : GEN
     337   232336655 : newblock(size_t n)
     338             : {
     339   232336655 :   long d = 0;
     340   232336655 :   long *x = (long *) pari_malloc((n + BL_HEAD)*sizeof(long)) + BL_HEAD;
     341             : 
     342   232341362 :   bl_size(x) = n;
     343   232341362 :   bl_refc(x) = 1;
     344   232341362 :   bl_next(x) = NULL;
     345   232341362 :   bl_prev(x) = cur_block;
     346   232341362 :   bl_num(x)  = next_block++;
     347   232341362 :   if (cur_block) bl_next(cur_block) = x;
     348   232341362 :   root_block = blockinsert(x, root_block, &d);
     349   232338320 :   if (DEBUGMEM > 2)
     350           0 :     err_printf("new block, size %6lu (no %ld): %08lx\n", n, next_block-1, x);
     351   232338338 :   return cur_block = x;
     352             : }
     353             : 
     354             : GEN
     355       37887 : gcloneref(GEN x)
     356             : {
     357       37887 :   if (isclone(x)) { ++bl_refc(x); return x; }
     358       37369 :   else return gclone(x);
     359             : }
     360             : 
     361             : void
     362           0 : gclone_refc(GEN x) { ++bl_refc(x); }
     363             : 
     364             : void
     365   382499665 : gunclone(GEN x)
     366             : {
     367   382499665 :   if (--bl_refc(x) > 0) return;
     368   232327361 :   BLOCK_SIGINT_START;
     369   232336163 :   root_block = blockdelete(x, root_block);
     370   232316098 :   if (bl_next(x)) bl_prev(bl_next(x)) = bl_prev(x);
     371             :   else
     372             :   {
     373    38315905 :     cur_block = bl_prev(x);
     374    38315905 :     next_block = bl_num(x);
     375             :   }
     376   232316098 :   if (bl_prev(x)) bl_next(bl_prev(x)) = bl_next(x);
     377   232316098 :   if (DEBUGMEM > 2)
     378           0 :     err_printf("killing block (no %ld): %08lx\n", bl_num(x), x);
     379   232316098 :   free((void*)bl_base(x)); /* pari_free not needed: we already block */
     380   232316098 :   BLOCK_SIGINT_END;
     381             : }
     382             : 
     383             : static void
     384   128067491 : vec_gunclone_deep(GEN x)
     385             : {
     386   128067491 :   long i, l = lg(x);
     387  3124763112 :   for (i = 1; i < l; i++) gunclone_deep(gel(x,i));
     388   128067506 : }
     389             : /* Recursively look for clones in the container and kill them. Then kill
     390             :  * container if clone. */
     391             : void
     392  3250796015 : gunclone_deep(GEN x)
     393             : {
     394             :   GEN v;
     395  3250796015 :   if (isclone(x) && bl_refc(x) > 1) { --bl_refc(x); return; }
     396  3198669235 :   BLOCK_SIGINT_START;
     397  3198669235 :   switch(typ(x))
     398             :   {
     399   128066113 :     case t_VEC: case t_COL: case t_MAT:
     400   128066113 :       vec_gunclone_deep(x);
     401   128066110 :       break;
     402        6175 :     case t_LIST:
     403        6175 :       if ((v = list_data(x))) { vec_gunclone_deep(v); gunclone(v); }
     404        6175 :       break;
     405             :   }
     406  3198669232 :   if (isclone(x)) gunclone(x);
     407  3198669109 :   BLOCK_SIGINT_END;
     408             : }
     409             : 
     410             : int
     411      357980 : pop_entree_block(entree *ep, long loc)
     412             : {
     413      357980 :   GEN x = (GEN)ep->value;
     414      357980 :   if (bl_num(x) < loc) return 0; /* older */
     415         455 :   if (DEBUGMEM>2)
     416           0 :     err_printf("popping %s (block no %ld)\n", ep->name, bl_num(x));
     417         455 :   gunclone_deep(x); return 1;
     418             : }
     419             : 
     420             : /***************************************************************************
     421             :  **                                                                       **
     422             :  **                           Export                                      **
     423             :  **                                                                       **
     424             :  ***************************************************************************/
     425             : 
     426             : static hashtable *export_hash;
     427             : static void
     428        1858 : pari_init_export(void)
     429             : {
     430        1858 :   export_hash = hash_create_str(1,0);
     431        1858 : }
     432             : static void
     433        1848 : pari_close_export(void)
     434             : {
     435        1848 :   hash_destroy(export_hash);
     436        1848 : }
     437             : 
     438             : /* Exported values are blocks, but do not have the clone bit set so that they
     439             :  * are not affected by clone_lock and ensure_nb, etc. */
     440             : 
     441             : void
     442          59 : export_add(const char *str, GEN val)
     443             : {
     444             :   hashentry *h;
     445          59 :   val = gclone(val); unsetisclone(val);
     446          59 :   h = hash_search(export_hash, (void*) str);
     447          59 :   if (h)
     448             :   {
     449          21 :     GEN v = (GEN)h->val;
     450          21 :     h->val = val;
     451          21 :     setisclone(v); gunclone(v);
     452             :   }
     453             :   else
     454          38 :     hash_insert(export_hash,(void*)str, (void*) val);
     455          59 : }
     456             : 
     457             : void
     458          24 : export_del(const char *str)
     459             : {
     460          24 :   hashentry *h = hash_remove(export_hash,(void*)str);
     461          24 :   if (h)
     462             :   {
     463          24 :     GEN v = (GEN)h->val;
     464          24 :     setisclone(v); gunclone(v);
     465          24 :     pari_free(h);
     466             :   }
     467          24 : }
     468             : 
     469             : GEN
     470        1500 : export_get(const char *str)
     471             : {
     472        1500 :   return hash_haskey_GEN(export_hash,(void*)str);
     473             : }
     474             : 
     475             : void
     476           6 : unexportall(void)
     477             : {
     478           6 :   pari_sp av = avma;
     479           6 :   GEN keys = hash_keys(export_hash);
     480           6 :   long i, l = lg(keys);
     481          24 :   for (i = 1; i < l; i++) mt_export_del((const char *)keys[i]);
     482           6 :   set_avma(av);
     483           6 : }
     484             : 
     485             : void
     486           6 : exportall(void)
     487             : {
     488             :   long i;
     489         816 :   for (i = 0; i < functions_tblsz; i++)
     490             :   {
     491             :     entree *ep;
     492        9024 :     for (ep = functions_hash[i]; ep; ep = ep->next)
     493        8214 :       if (EpVALENCE(ep)==EpVAR) mt_export_add(ep->name, (GEN)ep->value);
     494             :   }
     495           6 : }
     496             : 
     497             : /*********************************************************************/
     498             : /*                                                                   */
     499             : /*                       C STACK SIZE CONTROL                        */
     500             : /*                                                                   */
     501             : /*********************************************************************/
     502             : /* Avoid core dump on deep recursion. Adapted Perl code by Dominic Dunlop */
     503             : THREAD void *PARI_stack_limit = NULL;
     504             : 
     505             : #ifdef STACK_CHECK
     506             : 
     507             : #  ifdef __EMX__                                /* Emulate */
     508             : void
     509             : pari_stackcheck_init(void *pari_stack_base)
     510             : {
     511             :   if (!pari_stack_base) { PARI_stack_limit = NULL; return; }
     512             :   PARI_stack_limit = get_stack(1./16, 32*1024);
     513             : }
     514             : #  elif _WIN32
     515             : void
     516             : pari_stackcheck_init(void *pari_stack_base)
     517             : {
     518             :   ulong size = 1UL << 21;
     519             :   if (!pari_stack_base) { PARI_stack_limit = NULL; return; }
     520             :   if (size > (ulong)pari_stack_base)
     521             :     PARI_stack_limit = (void*)(((ulong)pari_stack_base) / 16);
     522             :   else
     523             :     PARI_stack_limit = (void*)((ulong)pari_stack_base - (size/16)*15);
     524             : }
     525             : #  else /* !__EMX__ && !_WIN32 */
     526             : /* Set PARI_stack_limit to (a little above) the lowest safe address that can be
     527             :  * used on the stack. Leave PARI_stack_limit at its initial value (NULL) to
     528             :  * show no check should be made [init failed]. Assume stack grows downward. */
     529             : void
     530      328213 : pari_stackcheck_init(void *pari_stack_base)
     531             : {
     532             :   struct rlimit rip;
     533             :   ulong size;
     534      328213 :   if (!pari_stack_base) { PARI_stack_limit = NULL; return; }
     535      328213 :   if (getrlimit(RLIMIT_STACK, &rip)) return;
     536      328220 :   size = rip.rlim_cur;
     537      328220 :   if (size == (ulong)RLIM_INFINITY || size > (ulong)pari_stack_base)
     538           0 :     PARI_stack_limit = (void*)(((ulong)pari_stack_base) / 16);
     539             :   else
     540      328220 :     PARI_stack_limit = (void*)((ulong)pari_stack_base - (size/16)*15);
     541             : }
     542             : #  endif /* !__EMX__ */
     543             : 
     544             : #else
     545             : void
     546             : pari_stackcheck_init(void *pari_stack_base)
     547             : {
     548             :   (void) pari_stack_base; PARI_stack_limit = NULL;
     549             : }
     550             : #endif /* STACK_CHECK */
     551             : 
     552             : /*******************************************************************/
     553             : /*                         HEAP TRAVERSAL                          */
     554             : /*******************************************************************/
     555             : struct getheap_t { long n, l; };
     556             : /* x is a block, not necessarily a clone [x[0] may not be set] */
     557             : static void
     558        6664 : f_getheap(GEN x, void *D)
     559             : {
     560        6664 :   struct getheap_t *T = (struct getheap_t*)D;
     561        6664 :   T->n++;
     562        6664 :   T->l += bl_size(x) + BL_HEAD;
     563        6664 : }
     564             : GEN
     565          84 : getheap(void)
     566             : {
     567          84 :   struct getheap_t T = { 0, 0 };
     568          84 :   traverseheap(&f_getheap, &T); return mkvec2s(T.n, T.l);
     569             : }
     570             : 
     571             : static void
     572       13412 : traverseheap_r(GEN bl, void(*f)(GEN, void *), void *data)
     573             : {
     574       13412 :   if (!bl) return;
     575        6664 :   traverseheap_r(bl_left(bl), f, data);
     576        6664 :   traverseheap_r(bl_right(bl), f, data);
     577        6664 :   f(bl, data);
     578             : }
     579             : 
     580             : void
     581          84 : traverseheap( void(*f)(GEN, void *), void *data)
     582             : {
     583          84 :   traverseheap_r(root_block,f, data);
     584          84 : }
     585             : 
     586             : /*********************************************************************/
     587             : /*                          DAEMON / FORK                            */
     588             : /*********************************************************************/
     589             : #if defined(HAS_WAITPID) && defined(HAS_SETSID)
     590             : /* Properly fork a process, detaching from main process group without creating
     591             :  * zombies on exit. Parent returns 1, son returns 0 */
     592             : int
     593          76 : pari_daemon(void)
     594             : {
     595          76 :   pid_t pid = fork();
     596          76 :   switch(pid) {
     597           0 :       case -1: return 1; /* father, fork failed */
     598           0 :       case 0:
     599           0 :         (void)setsid(); /* son becomes process group leader */
     600           0 :         if (fork()) _exit(0); /* now son exits, also when fork fails */
     601           0 :         break; /* grandson: its father is the son, which exited,
     602             :                 * hence father becomes 'init', that'll take care of it */
     603          76 :       default: /* father, fork succeeded */
     604          76 :         (void)waitpid(pid,NULL,0); /* wait for son to exit, immediate */
     605          76 :         return 1;
     606             :   }
     607             :   /* grandson. The silly '!' avoids a gcc-8 warning (unused value) */
     608           0 :   (void)!freopen("/dev/null","r",stdin);
     609           0 :   return 0;
     610             : }
     611             : #else
     612             : int
     613             : pari_daemon(void)
     614             : {
     615             :   pari_err_IMPL("pari_daemon without waitpid & setsid");
     616             :   return 0;
     617             : }
     618             : #endif
     619             : 
     620             : /*********************************************************************/
     621             : /*                                                                   */
     622             : /*                       SYSTEM INITIALIZATION                       */
     623             : /*                                                                   */
     624             : /*********************************************************************/
     625             : static int try_to_recover = 0;
     626             : THREAD VOLATILE int PARI_SIGINT_block = 0, PARI_SIGINT_pending = 0;
     627             : 
     628             : /*********************************************************************/
     629             : /*                         SIGNAL HANDLERS                           */
     630             : /*********************************************************************/
     631             : static void
     632           0 : dflt_sigint_fun(void) { pari_err(e_MISC, "user interrupt"); }
     633             : 
     634             : #if defined(_WIN32) || defined(__CYGWIN32__)
     635             : int win32ctrlc = 0, win32alrm = 0;
     636             : void
     637             : dowin32ctrlc(void)
     638             : {
     639             :   win32ctrlc = 0;
     640             :   cb_pari_sigint();
     641             : }
     642             : #endif
     643             : 
     644             : static void
     645           0 : pari_handle_SIGINT(void)
     646             : {
     647             : #ifdef _WIN32
     648             :   if (++win32ctrlc >= 5) _exit(3);
     649             : #else
     650           0 :   cb_pari_sigint();
     651             : #endif
     652           0 : }
     653             : 
     654             : typedef void (*pari_sighandler_t)(int);
     655             : 
     656             : pari_sighandler_t
     657       20378 : os_signal(int sig, pari_sighandler_t f)
     658             : {
     659             : #ifdef HAS_SIGACTION
     660             :   struct sigaction sa, oldsa;
     661             : 
     662       20378 :   sa.sa_handler = f;
     663       20378 :   sigemptyset(&sa.sa_mask);
     664       20378 :   sa.sa_flags = SA_NODEFER;
     665             : 
     666       20378 :   if (sigaction(sig, &sa, &oldsa)) return NULL;
     667       20378 :   return oldsa.sa_handler;
     668             : #else
     669             :   return signal(sig,f);
     670             : #endif
     671             : }
     672             : 
     673             : void
     674           0 : pari_sighandler(int sig)
     675             : {
     676             :   const char *msg;
     677             : #ifndef HAS_SIGACTION
     678             :   /*SYSV reset the signal handler in the handler*/
     679             :   (void)os_signal(sig,pari_sighandler);
     680             : #endif
     681           0 :   switch(sig)
     682             :   {
     683             : #ifdef SIGBREAK
     684             :     case SIGBREAK:
     685             :       if (PARI_SIGINT_block==1)
     686             :       {
     687             :         PARI_SIGINT_pending=SIGBREAK;
     688             :         mt_sigint();
     689             :       }
     690             :       else pari_handle_SIGINT();
     691             :       return;
     692             : #endif
     693             : 
     694             : #ifdef SIGINT
     695           0 :     case SIGINT:
     696           0 :       if (PARI_SIGINT_block==1)
     697             :       {
     698           0 :         PARI_SIGINT_pending=SIGINT;
     699           0 :         mt_sigint();
     700             :       }
     701           0 :       else pari_handle_SIGINT();
     702           0 :       return;
     703             : #endif
     704             : 
     705             : #ifdef SIGSEGV
     706           0 :     case SIGSEGV:
     707           0 :       msg="PARI/GP (Segmentation Fault)"; break;
     708             : #endif
     709             : #ifdef SIGBUS
     710           0 :     case SIGBUS:
     711           0 :       msg="PARI/GP (Bus Error)"; break;
     712             : #endif
     713             : #ifdef SIGFPE
     714           0 :     case SIGFPE:
     715           0 :       msg="PARI/GP (Floating Point Exception)"; break;
     716             : #endif
     717             : 
     718             : #ifdef SIGPIPE
     719           0 :     case SIGPIPE:
     720             :     {
     721           0 :       pariFILE *f = GP_DATA->pp->file;
     722           0 :       if (f && pari_outfile == f->file)
     723             :       {
     724           0 :         GP_DATA->pp->file = NULL; /* to avoid oo recursion on error */
     725           0 :         pari_outfile = stdout; pari_fclose(f);
     726           0 :         pari_err(e_MISC, "Broken Pipe, resetting file stack...");
     727             :       }
     728             :       return; /* LCOV_EXCL_LINE */
     729             :     }
     730             : #endif
     731             : 
     732           0 :     default: msg="signal handling"; break;
     733             :   }
     734           0 :   pari_err_BUG(msg);
     735             : }
     736             : 
     737             : void
     738        3706 : pari_sig_init(void (*f)(int))
     739             : {
     740             : #ifdef SIGBUS
     741        3706 :   (void)os_signal(SIGBUS,f);
     742             : #endif
     743             : #ifdef SIGFPE
     744        3706 :   (void)os_signal(SIGFPE,f);
     745             : #endif
     746             : #ifdef SIGINT
     747        3706 :   (void)os_signal(SIGINT,f);
     748             : #endif
     749             : #ifdef SIGBREAK
     750             :   (void)os_signal(SIGBREAK,f);
     751             : #endif
     752             : #ifdef SIGPIPE
     753        3706 :   (void)os_signal(SIGPIPE,f);
     754             : #endif
     755             : #ifdef SIGSEGV
     756        3706 :   (void)os_signal(SIGSEGV,f);
     757             : #endif
     758        3706 : }
     759             : 
     760             : /*********************************************************************/
     761             : /*                      STACK AND UNIVERSAL CONSTANTS                */
     762             : /*********************************************************************/
     763             : static void
     764        1858 : init_universal_constants(void)
     765             : {
     766        1858 :   gen_0  = (GEN)readonly_constants;
     767        1858 :   gnil   = (GEN)readonly_constants+2;
     768        1858 :   gen_1  = (GEN)readonly_constants+4;
     769        1858 :   gen_2  = (GEN)readonly_constants+7;
     770        1858 :   gen_m1 = (GEN)readonly_constants+10;
     771        1858 :   gen_m2 = (GEN)readonly_constants+13;
     772        1858 :   err_e_STACK = (GEN)readonly_constants+16;
     773        1858 :   ghalf  = (GEN)readonly_constants+18;
     774        1858 : }
     775             : 
     776             : static void
     777      326783 : pari_init_errcatch(void)
     778             : {
     779      326783 :   iferr_env = NULL;
     780      326783 :   global_err_data = NULL;
     781      326783 : }
     782             : 
     783             : void
     784        1886 : setalldebug(long n)
     785             : {
     786        1886 :   long i, l = numberof(pari_DEBUGLEVEL_ptr);
     787      115046 :   for (i = 0; i < l; i++) *pari_DEBUGLEVEL_ptr[i] = n;
     788        1886 : }
     789             : 
     790             : /*********************************************************************/
     791             : /*                           INIT DEFAULTS                           */
     792             : /*********************************************************************/
     793             : void
     794        1858 : pari_init_defaults(void)
     795             : {
     796             :   long i;
     797        1858 :   initout(1);
     798             : 
     799        1858 :   precreal = 128;
     800        1858 :   precdl = 16;
     801        1858 :   DEBUGLEVEL = 0;
     802        1858 :   setalldebug(0);
     803        1858 :   DEBUGMEM = 1;
     804        1858 :   disable_color = 1;
     805        1858 :   pari_logstyle = logstyle_none;
     806             : 
     807        1858 :   current_psfile = pari_strdup("pari.ps");
     808        1858 :   current_logfile= pari_strdup("pari.log");
     809        1858 :   pari_logfile = NULL;
     810             : 
     811        1858 :   pari_datadir = os_getenv("GP_DATA_DIR");
     812        1858 :   if (!pari_datadir)
     813             :   {
     814             : #if defined(_WIN32) || defined(__CYGWIN32__)
     815             :     if (paricfg_datadir[0]=='@' && paricfg_datadir[1]==0)
     816             :       pari_datadir = win32_datadir();
     817             :     else
     818             : #endif
     819        1858 :       pari_datadir = pari_strdup(paricfg_datadir);
     820             :   }
     821           0 :   else pari_datadir= pari_strdup(pari_datadir);
     822       14864 :   for (i=0; i<c_LAST; i++) gp_colors[i] = c_NONE;
     823        1858 : }
     824             : 
     825             : /*********************************************************************/
     826             : /*                   FUNCTION HASHTABLES, MODULES                    */
     827             : /*********************************************************************/
     828             : extern entree functions_basic[], functions_default[];
     829             : static void
     830        1858 : pari_init_functions(void)
     831             : {
     832        1858 :   pari_stack_init(&s_MODULES, sizeof(*MODULES),(void**)&MODULES);
     833        1858 :   pari_stack_pushp(&s_MODULES,functions_basic);
     834        1858 :   functions_hash = (entree**) pari_calloc(sizeof(entree*)*functions_tblsz);
     835        1858 :   pari_fill_hashtable(functions_hash, functions_basic);
     836        1858 :   defaults_hash = (entree**) pari_calloc(sizeof(entree*)*functions_tblsz);
     837        1858 :   pari_add_defaults_module(functions_default);
     838        1858 : }
     839             : 
     840             : void
     841        1848 : pari_add_module(entree *ep)
     842             : {
     843        1848 :   pari_fill_hashtable(functions_hash, ep);
     844        1848 :   pari_stack_pushp(&s_MODULES, ep);
     845        1848 : }
     846             : 
     847             : void
     848        1858 : pari_add_defaults_module(entree *ep)
     849        1858 : { pari_fill_hashtable(defaults_hash, ep); }
     850             : 
     851             : /*********************************************************************/
     852             : /*                       PARI MAIN STACK                             */
     853             : /*********************************************************************/
     854             : 
     855             : #ifdef HAS_MMAP
     856             : #define PARI_STACK_ALIGN (sysconf(_SC_PAGE_SIZE))
     857             : #ifndef MAP_ANONYMOUS
     858             : #define MAP_ANONYMOUS MAP_ANON
     859             : #endif
     860             : #ifndef MAP_NORESERVE
     861             : #define MAP_NORESERVE 0
     862             : #endif
     863             : static void *
     864      326797 : pari_mainstack_malloc(size_t size)
     865             : {
     866             :   void *b;
     867             :   /* Check that the system allows reserving "size" bytes. This is just
     868             :    * a check, we immediately free the memory. */
     869      326797 :   BLOCK_SIGINT_START;
     870      326797 :   b = mmap(NULL, size, PROT_READ|PROT_WRITE, MAP_PRIVATE|MAP_ANONYMOUS, -1, 0);
     871      326797 :   BLOCK_SIGINT_END;
     872      326797 :   if (b == MAP_FAILED) return NULL;
     873      326797 :   BLOCK_SIGINT_START;
     874      326797 :   munmap(b, size);
     875             : 
     876             :   /* Map again, this time with MAP_NORESERVE. On some operating systems
     877             :    * like Cygwin, this is needed because remapping with PROT_NONE and
     878             :    * MAP_NORESERVE does not work as expected. */
     879      326797 :   b = mmap(NULL, size, PROT_READ|PROT_WRITE,
     880             :                        MAP_PRIVATE|MAP_ANONYMOUS|MAP_NORESERVE, -1, 0);
     881      326797 :   BLOCK_SIGINT_END;
     882      326797 :   if (b == MAP_FAILED) return NULL;
     883      326797 :   return b;
     884             : }
     885             : 
     886             : static void
     887      326787 : pari_mainstack_mfree(void *s, size_t size)
     888             : {
     889      326787 :   BLOCK_SIGINT_START;
     890      326787 :   munmap(s, size);
     891      326787 :   BLOCK_SIGINT_END;
     892      326787 : }
     893             : 
     894             : /* Completely discard the memory mapped between the addresses "from"
     895             :  * and "to" (which must be page-aligned).
     896             :  *
     897             :  * We use mmap() with PROT_NONE, which means that the underlying memory
     898             :  * is freed and that the kernel should not commit memory for it. We
     899             :  * still keep the mapping such that we can change the flags to
     900             :  * PROT_READ|PROT_WRITE later.
     901             :  *
     902             :  * NOTE: remapping with MAP_FIXED and PROT_NONE is not the same as
     903             :  * calling mprotect(..., PROT_NONE) because the latter will keep the
     904             :  * memory committed (this is in particular relevant on Linux with
     905             :  * vm.overcommit = 2). This remains true even when calling
     906             :  * madvise(..., MADV_DONTNEED). */
     907             : static void
     908      443574 : pari_mainstack_mreset(pari_sp from, pari_sp to)
     909             : {
     910      443574 :   size_t s = to - from;
     911             :   void *addr, *res;
     912      443574 :   if (!s) return;
     913             : 
     914          21 :   addr = (void*)from;
     915          21 :   BLOCK_SIGINT_START;
     916          21 :   res = mmap(addr, s, PROT_NONE,
     917             :              MAP_FIXED|MAP_PRIVATE|MAP_ANONYMOUS|MAP_NORESERVE, -1, 0);
     918          21 :   BLOCK_SIGINT_END;
     919          21 :   if (res != addr) pari_err(e_MEM);
     920             : }
     921             : 
     922             : /* Commit (make available) the virtual memory mapped between the
     923             :  * addresses "from" and "to" (which must be page-aligned).
     924             :  * Return 0 if successful, -1 if failed. */
     925             : static int
     926      443574 : pari_mainstack_mextend(pari_sp from, pari_sp to)
     927             : {
     928      443574 :   size_t s = to - from;
     929             :   int ret;
     930      443574 :   BLOCK_SIGINT_START;
     931      443574 :   ret = mprotect((void*)from, s, PROT_READ|PROT_WRITE);
     932      443574 :   BLOCK_SIGINT_END;
     933      443574 :   return ret;
     934             : }
     935             : 
     936             : /* Set actual stack size to the given size. This sets st->size and
     937             :  * st->bot. If not enough system memory is available, this can fail.
     938             :  * Return 1 if successful, 0 if failed (in that case, st->size is not
     939             :  * changed) */
     940             : static int
     941      443574 : pari_mainstack_setsize(struct pari_mainstack *st, size_t size)
     942             : {
     943      443574 :   pari_sp newbot = st->top - size;
     944             :   /* Align newbot to pagesize */
     945      443574 :   pari_sp alignbot = newbot & ~(pari_sp)(PARI_STACK_ALIGN - 1);
     946      443574 :   if (pari_mainstack_mextend(alignbot, st->top))
     947             :   {
     948             :     /* Making the memory available did not work: limit vsize to the
     949             :      * current actual stack size. */
     950           0 :     st->vsize = st->size;
     951           0 :     pari_warn(warnstack, st->vsize);
     952           0 :     return 0;
     953             :   }
     954      443574 :   pari_mainstack_mreset(st->vbot, alignbot);
     955      443574 :   st->bot = newbot;
     956      443574 :   st->size = size;
     957      443574 :   return 1;
     958             : }
     959             : 
     960             : #else
     961             : #define PARI_STACK_ALIGN (0x40UL)
     962             : static void *
     963             : pari_mainstack_malloc(size_t s)
     964             : {
     965             :   char * tmp;
     966             :   BLOCK_SIGINT_START;
     967             :   tmp = malloc(s); /* NOT pari_malloc, e_MEM would be deadly */
     968             :   BLOCK_SIGINT_END;
     969             :   return tmp;
     970             : }
     971             : 
     972             : static void
     973             : pari_mainstack_mfree(void *s, size_t size) { (void) size; pari_free(s); }
     974             : 
     975             : static int
     976             : pari_mainstack_setsize(struct pari_mainstack *st, size_t size)
     977             : {
     978             :   st->bot = st->top - size;
     979             :   st->size = size;
     980             :   return 1;
     981             : }
     982             : 
     983             : #endif
     984             : 
     985             : static const size_t MIN_STACK = 500032UL;
     986             : static size_t
     987      653584 : fix_size(size_t a)
     988             : {
     989      653584 :   size_t ps = PARI_STACK_ALIGN;
     990      653584 :   size_t b = a & ~(ps - 1); /* Align */
     991      653584 :   if (b < a && b < ~(ps - 1)) b += ps;
     992      653584 :   if (b < MIN_STACK) b = MIN_STACK;
     993      653584 :   return b;
     994             : }
     995             : 
     996             : static void
     997      326797 : pari_mainstack_alloc(int numerr, struct pari_mainstack *st, size_t rsize, size_t vsize)
     998             : {
     999      326797 :   size_t sizemax = vsize ? vsize: rsize, s = fix_size(sizemax);
    1000             :   for (;;)
    1001             :   {
    1002      326797 :     st->vbot = (pari_sp)pari_mainstack_malloc(s);
    1003      326797 :     if (st->vbot) break;
    1004           0 :     if (s == MIN_STACK) pari_err(e_MEM); /* no way out. Die */
    1005           0 :     s = fix_size(s >> 1);
    1006           0 :     pari_warn(numerr, s);
    1007             :   }
    1008      326797 :   st->vsize = vsize ? s: 0;
    1009      326797 :   st->rsize = minuu(rsize, s);
    1010      326797 :   st->top = st->vbot+s;
    1011      326797 :   if (!pari_mainstack_setsize(st, st->rsize))
    1012             :   {
    1013             :     /* This should never happen since we only decrease the allocated space */
    1014           0 :     pari_err(e_MEM);
    1015             :   }
    1016      326797 :   st->memused = 0;
    1017      326797 : }
    1018             : 
    1019             : static void
    1020      326787 : pari_mainstack_free(struct pari_mainstack *st)
    1021             : {
    1022      326787 :   pari_mainstack_mfree((void*)st->vbot, st->vsize ? st->vsize : fix_size(st->rsize));
    1023      326787 :   st->top = st->bot = st->vbot = 0;
    1024      326787 :   st->size = st->vsize = 0;
    1025      326787 : }
    1026             : 
    1027             : static void
    1028         431 : pari_mainstack_resize(struct pari_mainstack *st, size_t rsize, size_t vsize)
    1029             : {
    1030         431 :   BLOCK_SIGINT_START;
    1031         431 :   pari_mainstack_free(st);
    1032         431 :   pari_mainstack_alloc(warnstack, st, rsize, vsize);
    1033         431 :   BLOCK_SIGINT_END;
    1034         431 : }
    1035             : 
    1036             : static void
    1037      326356 : pari_mainstack_use(struct pari_mainstack *st)
    1038             : {
    1039      326356 :   pari_mainstack = st;
    1040      326356 :   avma = st->top; /* don't use set_avma */
    1041      326356 : }
    1042             : 
    1043             : static void
    1044        1858 : paristack_alloc(size_t rsize, size_t vsize)
    1045             : {
    1046        1858 :   pari_mainstack_alloc(warnstack, pari_mainstack, rsize, vsize);
    1047        1858 :   pari_mainstack_use(pari_mainstack);
    1048        1858 : }
    1049             : 
    1050             : void
    1051           0 : paristack_setsize(size_t rsize, size_t vsize)
    1052             : {
    1053           0 :   pari_mainstack_resize(pari_mainstack, rsize, vsize);
    1054           0 :   pari_mainstack_use(pari_mainstack);
    1055           0 : }
    1056             : 
    1057             : void
    1058           0 : parivstack_resize(ulong newsize)
    1059             : {
    1060             :   size_t s;
    1061           0 :   if (newsize && newsize < pari_mainstack->rsize)
    1062           0 :     pari_err_DIM("stack sizes [parisizemax < parisize]");
    1063           0 :   if (newsize == pari_mainstack->vsize) return;
    1064           0 :   evalstate_reset();
    1065           0 :   paristack_setsize(pari_mainstack->rsize, newsize);
    1066           0 :   s = pari_mainstack->vsize ? pari_mainstack->vsize : pari_mainstack->rsize;
    1067           0 :   if (DEBUGMEM)
    1068           0 :     pari_warn(warner,"new maximum stack size = %lu (%.3f Mbytes)",
    1069             :               s, s/1048576.);
    1070           0 :   pari_init_errcatch();
    1071           0 :   cb_pari_err_recover(-1);
    1072             : }
    1073             : 
    1074             : void
    1075         438 : paristack_newrsize(ulong newsize)
    1076             : {
    1077         438 :   size_t s, vsize = pari_mainstack->vsize;
    1078         438 :   if (!newsize) newsize = pari_mainstack->rsize << 1;
    1079         438 :   if (newsize != pari_mainstack->rsize)
    1080         431 :     pari_mainstack_resize(pari_mainstack, newsize, vsize);
    1081         438 :   evalstate_reset();
    1082         438 :   s = pari_mainstack->rsize;
    1083         438 :   if (DEBUGMEM)
    1084         438 :     pari_warn(warner,"new stack size = %lu (%.3f Mbytes)", s, s/1048576.);
    1085         438 :   pari_init_errcatch();
    1086         438 :   cb_pari_err_recover(-1);
    1087           0 : }
    1088             : 
    1089             : void
    1090           0 : paristack_resize(ulong newsize)
    1091             : {
    1092           0 :   long size = pari_mainstack->size;
    1093           0 :   if (!newsize)
    1094           0 :     newsize = 2 * size;
    1095           0 :   newsize = minuu(newsize, pari_mainstack->vsize);
    1096           0 :   if (newsize <= pari_mainstack->size) return;
    1097           0 :   if (pari_mainstack_setsize(pari_mainstack, newsize))
    1098             :   {
    1099           0 :     if (DEBUGMEM)
    1100           0 :       pari_warn(warner, "increasing stack size to %lu", pari_mainstack->size);
    1101             :   }
    1102             :   else
    1103             :   {
    1104           0 :     pari_mainstack_setsize(pari_mainstack, size);
    1105           0 :     pari_err(e_STACK);
    1106             :   }
    1107             : }
    1108             : 
    1109             : void
    1110      116777 : parivstack_reset(void)
    1111             : {
    1112      116777 :   pari_mainstack_setsize(pari_mainstack, pari_mainstack->rsize);
    1113      116777 :   if (avma < pari_mainstack->bot)
    1114           0 :     pari_err_BUG("parivstack_reset [avma < bot]");
    1115      116777 : }
    1116             : 
    1117             : /* Enlarge the stack if needed such that the unused portion of the stack
    1118             :  * (between bot and avma) is large enough to contain x longs. */
    1119             : void
    1120          14 : new_chunk_resize(size_t x)
    1121             : {
    1122          14 :   if (pari_mainstack->vsize==0
    1123          14 :     || x > (avma-pari_mainstack->vbot) / sizeof(long)) pari_err(e_STACK);
    1124           0 :   while (x > (avma-pari_mainstack->bot) / sizeof(long))
    1125           0 :     paristack_resize(0);
    1126           0 : }
    1127             : 
    1128             : /*********************************************************************/
    1129             : /*                       PARI THREAD                                 */
    1130             : /*********************************************************************/
    1131             : 
    1132             : /* Initial PARI thread structure t with a stack of size s and
    1133             :  * argument arg */
    1134             : 
    1135             : static void
    1136      324249 : pari_thread_set_global(struct pari_global_state *gs)
    1137             : {
    1138      324249 :   setdebugvar(gs->debugvar);
    1139      324238 :   push_localbitprec(gs->bitprec);
    1140      324404 :   pari_set_primetab(gs->primetab);
    1141      324307 :   pari_set_seadata(gs->seadata);
    1142      324296 :   pari_set_varstate(gs->varpriority, &gs->varstate);
    1143      322610 : }
    1144             : 
    1145             : static void
    1146      324508 : pari_thread_get_global(struct pari_global_state *gs)
    1147             : {
    1148      324508 :   gs->debugvar = getdebugvar();
    1149      324508 :   gs->bitprec = get_localbitprec();
    1150      324508 :   gs->primetab = primetab;
    1151      324508 :   gs->seadata = pari_get_seadata();
    1152      324508 :   varstate_save(&gs->varstate);
    1153      324508 :   gs->varpriority = varpriority;
    1154      324508 : }
    1155             : 
    1156             : void
    1157      324508 : pari_thread_alloc(struct pari_thread *t, size_t s, GEN arg)
    1158             : {
    1159      324508 :   pari_mainstack_alloc(warnstackthread, &t->st,s,0);
    1160      324508 :   pari_thread_get_global(&t->gs);
    1161      324508 :   t->data = arg;
    1162      324508 : }
    1163             : 
    1164             : /* Initial PARI thread structure t with a stack of size s and virtual size v
    1165             :  * and argument arg */
    1166             : 
    1167             : void
    1168           0 : pari_thread_valloc(struct pari_thread *t, size_t s, size_t v, GEN arg)
    1169             : {
    1170           0 :   pari_mainstack_alloc(warnstackthread, &t->st,s,v);
    1171           0 :   pari_thread_get_global(&t->gs);
    1172           0 :   t->data = arg;
    1173           0 : }
    1174             : 
    1175             : void
    1176      324508 : pari_thread_free(struct pari_thread *t)
    1177             : {
    1178      324508 :   pari_mainstack_free(&t->st);
    1179      324508 : }
    1180             : 
    1181             : void
    1182      326356 : pari_thread_init(void)
    1183             : {
    1184             :   long var;
    1185      326356 :   pari_stackcheck_init((void*)&var);
    1186      326357 :   pari_init_blocks();
    1187      326351 :   pari_init_errcatch();
    1188      326338 :   pari_init_rand();
    1189      326303 :   pari_init_floats();
    1190      326302 :   pari_init_hgm();
    1191      326306 :   pari_init_parser();
    1192      326255 :   pari_init_compiler();
    1193      326174 :   pari_init_evaluator();
    1194      326109 :   pari_init_files();
    1195      326107 :   pari_init_ellcondfile();
    1196      326133 : }
    1197             : 
    1198             : void
    1199      325142 : pari_thread_close(void)
    1200             : {
    1201      325142 :   pari_thread_close_files();
    1202      321713 :   pari_close_evaluator();
    1203      325018 :   pari_close_compiler();
    1204      320665 :   pari_close_parser();
    1205      324648 :   pari_close_floats();
    1206      321179 :   pari_close_hgm();
    1207      320148 :   pari_close_blocks();
    1208      324112 : }
    1209             : 
    1210             : GEN
    1211      324504 : pari_thread_start(struct pari_thread *t)
    1212             : {
    1213      324504 :   pari_mainstack_use(&t->st);
    1214      324499 :   pari_thread_init();
    1215      324261 :   pari_thread_set_global(&t->gs);
    1216      322702 :   mt_thread_init();
    1217      322691 :   return t->data;
    1218             : }
    1219             : 
    1220             : /*********************************************************************/
    1221             : /*                       LIBPARI INIT / CLOSE                        */
    1222             : /*********************************************************************/
    1223             : 
    1224             : static void
    1225           0 : pari_exit(void)
    1226             : {
    1227           0 :   err_printf("  ***   Error in the PARI system. End of program.\n");
    1228           0 :   exit(1);
    1229             : }
    1230             : 
    1231             : static void
    1232           0 : dflt_err_recover(long errnum) { (void) errnum; pari_exit(); }
    1233             : 
    1234             : static void
    1235           0 : dflt_pari_quit(long err) { (void)err; /*do nothing*/; }
    1236             : 
    1237             : /* initialize PARI data. Initialize [new|old]fun to NULL for default set. */
    1238             : void
    1239        1858 : pari_init_opts(size_t parisize, ulong maxprime, ulong init_opts)
    1240             : {
    1241             :   ulong u;
    1242             : 
    1243        1858 :   pari_mt_nbthreads = 0;
    1244        1858 :   cb_pari_quit = dflt_pari_quit;
    1245        1858 :   cb_pari_init_histfile = NULL;
    1246        1858 :   cb_pari_get_line_interactive = NULL;
    1247        1858 :   cb_pari_fgets_interactive = NULL;
    1248        1858 :   cb_pari_whatnow = NULL;
    1249        1858 :   cb_pari_handle_exception = NULL;
    1250        1858 :   cb_pari_err_handle = pari_err_display;
    1251        1858 :   cb_pari_pre_recover = NULL;
    1252        1858 :   cb_pari_break_loop = NULL;
    1253        1858 :   cb_pari_is_interactive = NULL;
    1254        1858 :   cb_pari_start_output = NULL;
    1255        1858 :   cb_pari_sigint = dflt_sigint_fun;
    1256        1858 :   cb_pari_long_help = NULL;
    1257        1858 :   if (init_opts&INIT_JMPm) cb_pari_err_recover = dflt_err_recover;
    1258             : 
    1259        1858 :   pari_stackcheck_init(&u);
    1260        1858 :   pari_init_homedir();
    1261        1858 :   if (init_opts&INIT_DFTm) {
    1262           0 :     pari_init_defaults();
    1263           0 :     GP_DATA = default_gp_data();
    1264           0 :     pari_init_paths();
    1265             :   }
    1266             : 
    1267        1858 :   pari_mainstack = (struct pari_mainstack *) malloc(sizeof(*pari_mainstack));
    1268        1858 :   paristack_alloc(parisize, 0);
    1269        1858 :   init_universal_constants();
    1270        1858 :   pari_PRIMES = NULL;
    1271        1858 :   if (!(init_opts&INIT_noPRIMEm))
    1272             :   {
    1273           0 :     GP_DATA->primelimit = maxprime;
    1274           0 :     pari_init_primes(GP_DATA->primelimit);
    1275             :   }
    1276        1858 :   if (!(init_opts&INIT_noINTGMPm)) pari_kernel_init();
    1277        1858 :   pari_init_graphics();
    1278        1858 :   pari_thread_init();
    1279        1858 :   pari_set_primetab(NULL);
    1280        1858 :   pari_set_seadata(NULL);
    1281        1858 :   pari_init_functions();
    1282        1858 :   pari_init_export();
    1283        1858 :   pari_var_init();
    1284        1858 :   pari_init_timer();
    1285        1858 :   pari_init_buffers();
    1286        1858 :   (void)getabstime();
    1287        1858 :   try_to_recover = 1;
    1288        1858 :   if (!(init_opts&INIT_noIMTm)) pari_mt_init();
    1289        1858 :   if ((init_opts&INIT_SIGm)) pari_sig_init(pari_sighandler);
    1290        1858 : }
    1291             : 
    1292             : void
    1293           0 : pari_init(size_t parisize, ulong maxprime)
    1294           0 : { pari_init_opts(parisize, maxprime, INIT_JMPm | INIT_SIGm | INIT_DFTm); }
    1295             : 
    1296             : void
    1297        1848 : pari_close_opts(ulong init_opts)
    1298             : {
    1299             :   long i;
    1300             : 
    1301        1848 :   BLOCK_SIGINT_START;
    1302        1848 :   if ((init_opts&INIT_SIGm)) pari_sig_init(SIG_DFL);
    1303        1848 :   if (!(init_opts&INIT_noIMTm)) pari_mt_close();
    1304             : 
    1305        1848 :   pari_var_close(); /* must come before destruction of functions_hash */
    1306      251328 :   for (i = 0; i < functions_tblsz; i++)
    1307             :   {
    1308      249480 :     entree *ep = functions_hash[i];
    1309     2807072 :     while (ep) {
    1310     2557592 :       entree *EP = ep->next;
    1311     2557592 :       if (!EpSTATIC(ep)) { freeep(ep); free(ep); }
    1312     2557592 :       ep = EP;
    1313             :     }
    1314             :   }
    1315        1848 :   pari_close_mf();
    1316        1848 :   pari_thread_close();
    1317        1848 :   pari_close_export();
    1318        1848 :   pari_close_files();
    1319        1848 :   pari_close_homedir();
    1320        1848 :   if (!(init_opts&INIT_noINTGMPm)) pari_kernel_close();
    1321             : 
    1322        1848 :   free((void*)functions_hash);
    1323        1848 :   free((void*)defaults_hash);
    1324        1848 :   if (pari_PRIMES) pari_close_primes();
    1325        1848 :   free(current_logfile);
    1326        1848 :   free(current_psfile);
    1327        1848 :   pari_mainstack_free(pari_mainstack);
    1328        1848 :   free((void*)pari_mainstack);
    1329        1848 :   pari_stack_delete(&s_MODULES);
    1330        1848 :   if (pari_datadir) free(pari_datadir);
    1331        1848 :   if (init_opts&INIT_DFTm)
    1332             :   { /* delete GP_DATA */
    1333        1848 :     pari_close_paths();
    1334        1848 :     if (GP_DATA->hist->v) free((void*)GP_DATA->hist->v);
    1335        1848 :     if (GP_DATA->pp->cmd) free((void*)GP_DATA->pp->cmd);
    1336        1848 :     if (GP_DATA->help) free((void*)GP_DATA->help);
    1337        1848 :     if (GP_DATA->plothsizes) free((void*)GP_DATA->plothsizes);
    1338        1848 :     if (GP_DATA->colormap) pari_free(GP_DATA->colormap);
    1339        1848 :     if (GP_DATA->graphcolors) pari_free(GP_DATA->graphcolors);
    1340        1848 :     free((void*)GP_DATA->prompt);
    1341        1848 :     free((void*)GP_DATA->prompt_cont);
    1342        1848 :     free((void*)GP_DATA->histfile);
    1343             :   }
    1344        1848 :   BLOCK_SIGINT_END;
    1345        1848 : }
    1346             : 
    1347             : void
    1348        1848 : pari_close(void)
    1349        1848 : { pari_close_opts(INIT_JMPm | INIT_SIGm | INIT_DFTm); }
    1350             : 
    1351             : /*******************************************************************/
    1352             : /*                                                                 */
    1353             : /*                         ERROR RECOVERY                          */
    1354             : /*                                                                 */
    1355             : /*******************************************************************/
    1356             : void
    1357      140456 : gp_context_save(struct gp_context* rec)
    1358             : {
    1359      140456 :   rec->prettyp = GP_DATA->fmt->prettyp;
    1360      140456 :   rec->listloc = next_block;
    1361      140456 :   rec->iferr_env = iferr_env;
    1362      140456 :   rec->err_data  = global_err_data;
    1363      140456 :   varstate_save(&rec->var);
    1364      140456 :   evalstate_save(&rec->eval);
    1365      140456 :   parsestate_save(&rec->parse);
    1366      140456 :   filestate_save(&rec->file);
    1367      140456 : }
    1368             : 
    1369             : void
    1370       12739 : gp_context_restore(struct gp_context* rec)
    1371             : {
    1372             :   long i;
    1373             : 
    1374       12739 :   if (!try_to_recover) return;
    1375             :   /* disable gp_context_restore() and SIGINT */
    1376       12739 :   try_to_recover = 0;
    1377       12739 :   BLOCK_SIGINT_START
    1378       12739 :   if (DEBUGMEM>2) err_printf("entering recover(), loc = %ld\n", rec->listloc);
    1379       12739 :   evalstate_restore(&rec->eval);
    1380       12739 :   parsestate_restore(&rec->parse);
    1381       12739 :   filestate_restore(&rec->file);
    1382       12739 :   global_err_data = rec->err_data;
    1383       12739 :   iferr_env = rec->iferr_env;
    1384       12739 :   GP_DATA->fmt->prettyp = rec->prettyp;
    1385             : 
    1386     1732504 :   for (i = 0; i < functions_tblsz; i++)
    1387             :   {
    1388     1719765 :     entree *ep = functions_hash[i];
    1389    20194232 :     while (ep)
    1390             :     {
    1391    18474467 :       entree *EP = ep->next;
    1392    18474467 :       switch(EpVALENCE(ep))
    1393             :       {
    1394      390952 :         case EpVAR:
    1395      391407 :           while (pop_val_if_newer(ep,rec->listloc)) /* empty */;
    1396      390952 :           break;
    1397      694768 :         case EpNEW: break;
    1398             :       }
    1399    18474467 :       ep = EP;
    1400             :     }
    1401             :   }
    1402       12739 :   varstate_restore(&rec->var);
    1403       12739 :   if (DEBUGMEM>2) err_printf("leaving recover()\n");
    1404       12739 :   BLOCK_SIGINT_END
    1405       12739 :   try_to_recover = 1;
    1406             : }
    1407             : 
    1408             : static void
    1409       12664 : err_recover(long numerr)
    1410             : {
    1411       12664 :   if (cb_pari_pre_recover)
    1412       12664 :     cb_pari_pre_recover(numerr);
    1413           0 :   evalstate_reset();
    1414           0 :   killallfiles();
    1415           0 :   pari_init_errcatch();
    1416           0 :   cb_pari_err_recover(numerr);
    1417           0 : }
    1418             : 
    1419             : static void
    1420       13438 : err_init(void)
    1421             : {
    1422             :   /* make sure pari_err msg starts at the beginning of line */
    1423       13438 :   if (!pari_last_was_newline()) pari_putc('\n');
    1424       13438 :   pariOut->flush();
    1425       13438 :   pariErr->flush();
    1426       13438 :   out_term_color(pariErr, c_ERR);
    1427       13438 : }
    1428             : 
    1429             : static void
    1430       13298 : err_init_msg(int user)
    1431             : {
    1432             :   const char *gp_function_name;
    1433       13298 :   out_puts(pariErr, "  *** ");
    1434       13298 :   if (!user && (gp_function_name = closure_func_err()))
    1435        9627 :     out_printf(pariErr, "%s: ", gp_function_name);
    1436             :   else
    1437        3671 :     out_puts(pariErr, "  ");
    1438       13298 : }
    1439             : 
    1440             : void
    1441         753 : pari_warn(int numerr, ...)
    1442             : {
    1443             :   char *ch1;
    1444             :   va_list ap;
    1445             : 
    1446         753 :   va_start(ap,numerr);
    1447             : 
    1448         753 :   err_init();
    1449         753 :   err_init_msg(numerr==warnuser || numerr==warnstack);
    1450         753 :   switch (numerr)
    1451             :   {
    1452           7 :     case warnuser:
    1453           7 :       out_puts(pariErr, "user warning: ");
    1454           7 :       out_print1(pariErr, NULL, va_arg(ap, GEN), f_RAW);
    1455           7 :       break;
    1456             : 
    1457           0 :     case warnmem:
    1458           0 :       out_puts(pariErr, "collecting garbage in "); ch1=va_arg(ap, char*);
    1459           0 :       out_vprintf(pariErr, ch1,ap); out_putc(pariErr, '.');
    1460           0 :       break;
    1461             : 
    1462         746 :     case warner:
    1463         746 :       out_puts(pariErr, "Warning: "); ch1=va_arg(ap, char*);
    1464         746 :       out_vprintf(pariErr, ch1,ap); out_putc(pariErr, '.');
    1465         746 :       break;
    1466             : 
    1467           0 :     case warnprec:
    1468           0 :       out_vprintf(pariErr, "Warning: increasing prec in %s; new prec = %ld",
    1469             :                       ap);
    1470           0 :       break;
    1471             : 
    1472           0 :     case warnfile:
    1473           0 :       out_puts(pariErr, "Warning: failed to "),
    1474           0 :       ch1 = va_arg(ap, char*);
    1475           0 :       out_printf(pariErr, "%s: %s", ch1, va_arg(ap, char*));
    1476           0 :       break;
    1477             : 
    1478           0 :     case warnstack:
    1479             :     case warnstackthread:
    1480             :     {
    1481           0 :       ulong  s = va_arg(ap, ulong);
    1482             :       char buf[128];
    1483           0 :       const char * stk = numerr == warnstackthread
    1484           0 :                          || mt_is_thread() ? "thread": "PARI";
    1485           0 :       sprintf(buf,"Warning: not enough memory, new %s stack %lu", stk, s);
    1486           0 :       out_puts(pariErr,buf);
    1487           0 :       break;
    1488             :     }
    1489             :   }
    1490         753 :   va_end(ap);
    1491         753 :   out_term_color(pariErr, c_NONE);
    1492         753 :   out_putc(pariErr, '\n');
    1493         753 :   pariErr->flush();
    1494         753 : }
    1495             : void
    1496           0 : pari_sigint(const char *time_s)
    1497             : {
    1498           0 :   int recover=0;
    1499           0 :   BLOCK_SIGALRM_START
    1500           0 :   err_init();
    1501           0 :   mt_break_recover();
    1502           0 :   closure_err(0);
    1503           0 :   err_init_msg(0);
    1504           0 :   out_puts(pariErr, "user interrupt after ");
    1505           0 :   out_puts(pariErr, time_s);
    1506           0 :   out_term_color(pariErr, c_NONE);
    1507           0 :   pariErr->flush();
    1508           0 :   if (cb_pari_handle_exception)
    1509           0 :     recover = cb_pari_handle_exception(-1);
    1510           0 :   if (!recover && !block)
    1511           0 :     PARI_SIGINT_pending = 0;
    1512           0 :   BLOCK_SIGINT_END
    1513           0 :   if (!recover) err_recover(e_MISC);
    1514           0 : }
    1515             : 
    1516             : #define retmkerr2(x,y)\
    1517             :   do { GEN _v = cgetg(3, t_ERROR);\
    1518             :        _v[1] = (x);\
    1519             :        gel(_v,2) = (y); return _v; } while(0)
    1520             : #define retmkerr3(x,y,z)\
    1521             :   do { GEN _v = cgetg(4, t_ERROR);\
    1522             :        _v[1] = (x);\
    1523             :        gel(_v,2) = (y);\
    1524             :        gel(_v,3) = (z); return _v; } while(0)
    1525             : #define retmkerr4(x,y,z,t)\
    1526             :   do { GEN _v = cgetg(5, t_ERROR);\
    1527             :        _v[1] = (x);\
    1528             :        gel(_v,2) = (y);\
    1529             :        gel(_v,3) = (z);\
    1530             :        gel(_v,4) = (t); return _v; } while(0)
    1531             : #define retmkerr5(x,y,z,t,u)\
    1532             :   do { GEN _v = cgetg(6, t_ERROR);\
    1533             :        _v[1] = (x);\
    1534             :        gel(_v,2) = (y);\
    1535             :        gel(_v,3) = (z);\
    1536             :        gel(_v,4) = (t);\
    1537             :        gel(_v,5) = (u); return _v; } while(0)
    1538             : #define retmkerr6(x,y,z,t,u,v)\
    1539             :   do { GEN _v = cgetg(7, t_ERROR);\
    1540             :        _v[1] = (x);\
    1541             :        gel(_v,2) = (y);\
    1542             :        gel(_v,3) = (z);\
    1543             :        gel(_v,4) = (t);\
    1544             :        gel(_v,5) = (u);\
    1545             :        gel(_v,6) = (v); return _v; } while(0)
    1546             : 
    1547             : static GEN
    1548       49396 : pari_err2GEN(long numerr, va_list ap)
    1549             : {
    1550       49396 :   switch ((enum err_list) numerr)
    1551             :   {
    1552         140 :   case e_SYNTAX:
    1553             :     {
    1554         140 :       const char *msg = va_arg(ap, char*);
    1555         140 :       const char *s = va_arg(ap,char *);
    1556         140 :       const char *entry = va_arg(ap,char *);
    1557         140 :       retmkerr3(numerr,strtoGENstr(msg), mkvecsmall2((long)s,(long)entry));
    1558             :     }
    1559         355 :   case e_MISC: case e_ALARM:
    1560             :     {
    1561         355 :       const char *ch1 = va_arg(ap, char*);
    1562         355 :       retmkerr2(numerr, gvsprintf(ch1,ap));
    1563             :     }
    1564        2723 :   case e_NOTFUNC:
    1565             :   case e_USER:
    1566        2723 :     retmkerr2(numerr,va_arg(ap, GEN));
    1567           0 :   case e_FILE:
    1568             :     {
    1569           0 :       const char *f = va_arg(ap, const char*);
    1570           0 :       retmkerr3(numerr, strtoGENstr(f), strtoGENstr(va_arg(ap, char*)));
    1571             :     }
    1572          36 :   case e_FILEDESC:
    1573             :     {
    1574          36 :       const char *f = va_arg(ap, const char*);
    1575          36 :       retmkerr3(numerr, strtoGENstr(f), stoi(va_arg(ap, long)));
    1576             :     }
    1577        1819 :   case e_OVERFLOW:
    1578             :   case e_IMPL:
    1579             :   case e_DIM:
    1580             :   case e_CONSTPOL:
    1581             :   case e_ROOTS0:
    1582             :   case e_FLAG:
    1583             :   case e_PREC:
    1584             :   case e_BUG:
    1585             :   case e_ARCH:
    1586             :   case e_PACKAGE:
    1587        1819 :     retmkerr2(numerr, strtoGENstr(va_arg(ap, char*)));
    1588        1687 :   case e_MODULUS:
    1589             :   case e_VAR:
    1590             :     {
    1591        1687 :       const char *f = va_arg(ap, const char*);
    1592        1687 :       GEN x = va_arg(ap, GEN);
    1593        1687 :       GEN y = va_arg(ap, GEN);
    1594        1687 :       retmkerr4(numerr, strtoGENstr(f), x,y);
    1595             :     }
    1596       35430 :   case e_INV:
    1597             :   case e_IRREDPOL:
    1598             :   case e_PRIME:
    1599             :   case e_SQRTN:
    1600             :   case e_TYPE:
    1601             :     {
    1602       35430 :       const char *f = va_arg(ap, const char*);
    1603       35430 :       GEN x = va_arg(ap, GEN);
    1604       35430 :       retmkerr3(numerr, strtoGENstr(f), x);
    1605             :     }
    1606        3990 :   case e_COPRIME: case e_OP: case e_TYPE2:
    1607             :     {
    1608        3990 :       const char *f = va_arg(ap, const char*);
    1609        3990 :       GEN x = va_arg(ap, GEN);
    1610        3990 :       GEN y = va_arg(ap, GEN);
    1611        3990 :       retmkerr4(numerr,strtoGENstr(f),x,y);
    1612             :     }
    1613         214 :   case e_COMPONENT:
    1614             :     {
    1615         214 :       const char *f= va_arg(ap, const char *);
    1616         214 :       const char *op = va_arg(ap, const char *);
    1617         214 :       GEN l = va_arg(ap, GEN);
    1618         214 :       GEN x = va_arg(ap, GEN);
    1619         214 :       retmkerr5(numerr,strtoGENstr(f),strtoGENstr(op),l,x);
    1620             :     }
    1621        2743 :   case e_DOMAIN:
    1622             :     {
    1623        2743 :       const char *f = va_arg(ap, const char*);
    1624        2743 :       const char *v = va_arg(ap, const char *);
    1625        2743 :       const char *op = va_arg(ap, const char *);
    1626        2743 :       GEN l = va_arg(ap, GEN);
    1627        2743 :       GEN x = va_arg(ap, GEN);
    1628        2743 :       retmkerr6(numerr,strtoGENstr(f),strtoGENstr(v),strtoGENstr(op),l,x);
    1629             :     }
    1630         245 :   case e_PRIORITY:
    1631             :     {
    1632         245 :       const char *f = va_arg(ap, const char*);
    1633         245 :       GEN x = va_arg(ap, GEN);
    1634         245 :       const char *op = va_arg(ap, const char *);
    1635         245 :       long v = va_arg(ap, long);
    1636         245 :       retmkerr5(numerr,strtoGENstr(f),x,strtoGENstr(op),stoi(v));
    1637             :     }
    1638           0 :   case e_MAXPRIME:
    1639           0 :     retmkerr2(numerr, utoi(va_arg(ap, ulong)));
    1640          14 :   case e_STACK:
    1641          14 :     return err_e_STACK;
    1642           0 :   case e_STACKTHREAD:
    1643           0 :     retmkerr3(numerr, utoi(va_arg(ap, ulong)), utoi(va_arg(ap, ulong)));
    1644           0 :   default:
    1645           0 :     return mkerr(numerr);
    1646             :   }
    1647             : }
    1648             : 
    1649             : static char *
    1650        7322 : type_dim(GEN x)
    1651             : {
    1652        7322 :   char *v = stack_malloc(64);
    1653        7322 :   switch(typ(x))
    1654             :   {
    1655         133 :     case t_MAT:
    1656             :     {
    1657         133 :       long l = lg(x), r = (l == 1)? 1: lgcols(x);
    1658         133 :       sprintf(v, "t_MAT (%ld x %ld)", r-1,l-1);
    1659         133 :       break;
    1660             :     }
    1661         133 :     case t_COL:
    1662         133 :       sprintf(v, "t_COL (%ld elts)", lg(x)-1);
    1663         133 :       break;
    1664         259 :     case t_VEC:
    1665         259 :       sprintf(v, "t_VEC (%ld elts)", lg(x)-1);
    1666         259 :       break;
    1667        6797 :     default:
    1668        6797 :       v = (char*)type_name(typ(x));
    1669             :   }
    1670        7322 :   return v;
    1671             : }
    1672             : 
    1673             : static char *
    1674        3696 : gdisplay(GEN x)
    1675             : {
    1676        3696 :   char *s = GENtostr_raw(x);
    1677        3696 :   if (strlen(s) < 1600) return s;
    1678          35 :   if (! GP_DATA->breakloop) return (char*)"(...)";
    1679           0 :   return stack_sprintf("\n  ***  (...) Huge %s omitted; you can access it via dbg_err()", type_name(typ(x)));
    1680             : }
    1681             : 
    1682             : char *
    1683       21463 : pari_err2str(GEN e)
    1684             : {
    1685       21463 :   long numerr = err_get_num(e);
    1686       21463 :   switch ((enum err_list) numerr)
    1687             :   {
    1688           0 :   case e_ALARM:
    1689           0 :     return pari_sprintf("alarm interrupt after %Ps.",gel(e,2));
    1690         343 :   case e_MISC:
    1691         343 :     return pari_sprintf("%Ps.",gel(e,2));
    1692             : 
    1693           0 :   case e_ARCH:
    1694           0 :     return pari_sprintf("sorry, '%Ps' not available on this system.",gel(e,2));
    1695          14 :   case e_BUG:
    1696          14 :     return pari_sprintf("bug in %Ps, please report.",gel(e,2));
    1697          21 :   case e_CONSTPOL:
    1698          21 :     return pari_sprintf("constant polynomial in %Ps.", gel(e,2));
    1699          84 :   case e_COPRIME:
    1700         252 :     return pari_sprintf("elements not coprime in %Ps:\n    %s\n    %s",
    1701          84 :                         gel(e,2), gdisplay(gel(e,3)), gdisplay(gel(e,4)));
    1702         711 :   case e_DIM:
    1703         711 :     return pari_sprintf("inconsistent dimensions in %Ps.", gel(e,2));
    1704           0 :   case e_FILE:
    1705           0 :     return pari_sprintf("error opening %Ps: `%Ps'.", gel(e,2), gel(e,3));
    1706          36 :   case e_FILEDESC:
    1707          36 :     return pari_sprintf("invalid file descriptor in %Ps [%Ps]", gel(e,2), gel(e,3));
    1708          91 :   case e_FLAG:
    1709          91 :     return pari_sprintf("invalid flag in %Ps.", gel(e,2));
    1710         462 :   case e_IMPL:
    1711         462 :     return pari_sprintf("sorry, %Ps is not yet implemented.", gel(e,2));
    1712           0 :   case e_PACKAGE:
    1713           0 :     return pari_sprintf("package %Ps is required, please install it.", gel(e,2));
    1714         651 :   case e_INV:
    1715         651 :     return pari_sprintf("impossible inverse in %Ps: %s.", gel(e,2),
    1716         651 :                         gdisplay(gel(e,3)));
    1717          63 :   case e_IRREDPOL:
    1718         126 :     return pari_sprintf("not an irreducible polynomial in %Ps: %s.",
    1719          63 :                         gel(e,2), gdisplay(gel(e,3)));
    1720           0 :   case e_MAXPRIME:
    1721             :     {
    1722           0 :       const char * msg = "not enough precomputed primes";
    1723           0 :       ulong c = itou(gel(e,2));
    1724           0 :       if (c) return pari_sprintf("%s, need primelimit ~ %lu.",msg, c);
    1725           0 :       else   return pari_strdup(msg);
    1726             :     }
    1727           0 :   case e_MEM:
    1728           0 :     return pari_strdup("not enough memory");
    1729        1316 :   case e_MODULUS:
    1730             :     {
    1731        1316 :       GEN x = gel(e,3), y = gel(e,4);
    1732        1316 :       return pari_sprintf("inconsistent moduli in %Ps: %s != %s",
    1733        1316 :                           gel(e,2), gdisplay(x), gdisplay(y));
    1734             :     }
    1735           0 :   case e_NONE: return NULL;
    1736        2709 :   case e_NOTFUNC:
    1737        2709 :     return pari_strdup("not a function in function call");
    1738        3661 :   case e_OP: case e_TYPE2:
    1739             :     {
    1740        3661 :       pari_sp av = avma;
    1741             :       char *v;
    1742        3661 :       const char *f, *op = GSTR(gel(e,2));
    1743        3661 :       const char *what = numerr == e_OP? "inconsistent": "forbidden";
    1744        3661 :       GEN x = gel(e,3);
    1745        3661 :       GEN y = gel(e,4);
    1746        3661 :       switch(*op)
    1747             :       {
    1748          21 :       case '+': f = "addition"; break;
    1749         175 :       case '*': f = "multiplication"; break;
    1750        2744 :       case '/': case '%': case '\\': f = "division"; break;
    1751           0 :       case '=': op = "-->"; f = "assignment"; break;
    1752         721 :       default:  f = op; op = ","; break;
    1753             :       }
    1754        3661 :       v = pari_sprintf("%s %s %s %s %s.", what,f,type_dim(x),op,type_dim(y));
    1755        3661 :       set_avma(av); return v;
    1756             :     }
    1757         214 :   case e_COMPONENT:
    1758             :     {
    1759         214 :       const char *f= GSTR(gel(e,2));
    1760         214 :       const char *op= GSTR(gel(e,3));
    1761         214 :       GEN l = gel(e,4);
    1762         214 :       if (!*f)
    1763         154 :         return pari_sprintf("nonexistent component: index %s %Ps",op,l);
    1764          60 :       return pari_sprintf("nonexistent component in %s: index %s %Ps",f,op,l);
    1765             :     }
    1766        2633 :   case e_DOMAIN:
    1767             :     {
    1768        2633 :       const char *f = GSTR(gel(e,2));
    1769        2633 :       const char *v = GSTR(gel(e,3));
    1770        2633 :       const char *op= GSTR(gel(e,4));
    1771        2633 :       GEN l = gel(e,5);
    1772        2633 :       if (!*op)
    1773          42 :         return pari_sprintf("domain error in %s: %s out of range",f,v);
    1774        2591 :       return pari_sprintf("domain error in %s: %s %s %Ps",f,v,op,l);
    1775             :     }
    1776         196 :   case e_PRIORITY:
    1777             :     {
    1778         196 :       const char *f = GSTR(gel(e,2));
    1779         196 :       long vx = gvar(gel(e,3));
    1780         196 :       const char *op= GSTR(gel(e,4));
    1781         196 :       long v = itos(gel(e,5));
    1782         196 :       return pari_sprintf("incorrect priority in %s: variable %Ps %s %Ps",f,
    1783             :              pol_x(vx), op, pol_x(v));
    1784             :     }
    1785         161 :   case e_OVERFLOW:
    1786         161 :     return pari_sprintf("overflow in %Ps.", gel(e,2));
    1787         231 :   case e_PREC:
    1788         231 :     return pari_sprintf("precision too low in %Ps.", gel(e,2));
    1789          98 :   case e_PRIME:
    1790         196 :     return pari_sprintf("not a prime number in %Ps: %s.",
    1791          98 :                         gel(e,2), gdisplay(gel(e,3)));
    1792          63 :   case e_ROOTS0:
    1793          63 :     return pari_sprintf("zero polynomial in %Ps.", gel(e,2));
    1794          84 :   case e_SQRTN:
    1795         168 :     return pari_sprintf("not an n-th power residue in %Ps: %s.",
    1796          84 :                         gel(e,2), gdisplay(gel(e,3)));
    1797          14 :   case e_STACK:
    1798             :   case e_STACKTHREAD:
    1799             :   {
    1800          14 :     const char *what = numerr == e_STACK? "PARI": "thread";
    1801          14 :     const char *var = numerr == e_STACK? "parisizemax": "threadsizemax";
    1802          14 :     size_t s = numerr == e_STACK? pari_mainstack->vsize: GP_DATA->threadsizemax;
    1803          14 :     char *hint = (char *) pari_malloc(512*sizeof(char));
    1804          14 :     char *buf = (char *) pari_malloc(512*sizeof(char));
    1805          14 :     if (s)
    1806           0 :       sprintf(hint,"you can increase '%s' using default()", var);
    1807             :     else
    1808             :     {
    1809          14 :       s = (numerr != e_STACK || !GP_DATA->threadsize)? pari_mainstack->rsize
    1810          28 :                                                      : GP_DATA->threadsize;
    1811          14 :       sprintf(hint,"set '%s' to a nonzero value in your GPRC", var);
    1812             :     }
    1813          14 :     sprintf(buf, "the %s stack overflows !\n"
    1814             :                  "  current stack size: %lu (%.3f Mbytes)\n  [hint] %s\n",
    1815          14 :                  what, (ulong)s, (double)s/1048576., hint);
    1816          14 :     return buf;
    1817             :   }
    1818           0 :   case e_SYNTAX:
    1819           0 :     return pari_strdup(GSTR(gel(e,2)));
    1820        7222 :   case e_TYPE:
    1821       14444 :     return pari_sprintf("incorrect type in %Ps (%s).",
    1822        7222 :                         gel(e,2), type_name(typ(gel(e,3))));
    1823          14 :   case e_USER:
    1824          14 :     return pari_sprint0("user error: ", gel(e,2), f_RAW);
    1825         371 :   case e_VAR:
    1826             :     {
    1827         371 :       GEN x = gel(e,3), y = gel(e,4);
    1828        1113 :       return pari_sprintf("inconsistent variables in %Ps, %Ps != %Ps.",
    1829         371 :                           gel(e,2), pol_x(varn(x)), pol_x(varn(y)));
    1830             :     }
    1831             :   }
    1832             :   return NULL; /*LCOV_EXCL_LINE*/
    1833             : }
    1834             : 
    1835             : int
    1836       12685 : pari_err_display(GEN err)
    1837             : {
    1838       12685 :   long numerr=err_get_num(err);
    1839       12685 :   err_init();
    1840       12685 :   if (numerr==e_SYNTAX)
    1841             :   {
    1842         140 :     const char *msg = GSTR(gel(err,2));
    1843         140 :     const char *s     = (const char *) gmael(err,3,1);
    1844         140 :     const char *entry = (const char *) gmael(err,3,2);
    1845         140 :     print_errcontext(pariErr, msg, s, entry);
    1846             :   }
    1847             :   else
    1848             :   {
    1849             :     char *s;
    1850       12545 :     closure_err(0);
    1851       12545 :     err_init_msg(numerr==e_USER);
    1852       12545 :     s = pari_err2str(err); pariErr->puts(s); pari_free(s);
    1853       12545 :     if (numerr==e_NOTFUNC)
    1854             :     {
    1855        2709 :       GEN fun = gel(err,2);
    1856        2709 :       if (gequalX(fun))
    1857             :       {
    1858        2709 :         entree *ep = varentries[varn(fun)];
    1859        2709 :         const char *t = ep->name;
    1860        2709 :         if (cb_pari_whatnow) cb_pari_whatnow(pariErr,t,1);
    1861             :       }
    1862             :     }
    1863             :   }
    1864       12671 :   out_term_color(pariErr, c_NONE);
    1865       12671 :   pariErr->flush(); return 0;
    1866             : }
    1867             : 
    1868             : void
    1869       49408 : pari_err(int numerr, ...)
    1870             : {
    1871             :   va_list ap;
    1872             :   GEN E;
    1873             : 
    1874       49408 :   va_start(ap,numerr);
    1875             : 
    1876       49408 :   if (numerr)
    1877       49396 :     E = pari_err2GEN(numerr,ap);
    1878             :   else
    1879             :   {
    1880          12 :     E = va_arg(ap,GEN);
    1881          12 :     numerr = err_get_num(E);
    1882             :   }
    1883       49407 :   global_err_data = E;
    1884       49407 :   if (*iferr_env) longjmp(*iferr_env, numerr);
    1885       12695 :   mt_err_recover(numerr);
    1886       12685 :   va_end(ap);
    1887       25356 :   if (cb_pari_err_handle &&
    1888       12685 :       cb_pari_err_handle(E)) return;
    1889       25333 :   if (cb_pari_handle_exception &&
    1890       12669 :       cb_pari_handle_exception(numerr)) return;
    1891       12664 :   err_recover(numerr);
    1892             : }
    1893             : 
    1894             : GEN
    1895       73440 : pari_err_last(void) { return global_err_data; }
    1896             : 
    1897             : const char *
    1898       27109 : numerr_name(long numerr)
    1899             : {
    1900       27109 :   switch ((enum err_list) numerr)
    1901             :   {
    1902           0 :   case e_ALARM:    return "e_ALARM";
    1903           0 :   case e_ARCH:     return "e_ARCH";
    1904           0 :   case e_BUG:      return "e_BUG";
    1905           0 :   case e_COMPONENT: return "e_COMPONENT";
    1906           0 :   case e_CONSTPOL: return "e_CONSTPOL";
    1907           0 :   case e_COPRIME:  return "e_COPRIME";
    1908           0 :   case e_DIM:      return "e_DIM";
    1909          56 :   case e_DOMAIN:   return "e_DOMAIN";
    1910           0 :   case e_FILE:     return "e_FILE";
    1911           0 :   case e_FILEDESC: return "e_FILEDESC";
    1912           7 :   case e_FLAG:     return "e_FLAG";
    1913          49 :   case e_IMPL:     return "e_IMPL";
    1914       19101 :   case e_INV:      return "e_INV";
    1915           0 :   case e_IRREDPOL: return "e_IRREDPOL";
    1916           0 :   case e_MAXPRIME: return "e_MAXPRIME";
    1917           0 :   case e_MEM:      return "e_MEM";
    1918           0 :   case e_MISC:     return "e_MISC";
    1919           0 :   case e_MODULUS:  return "e_MODULUS";
    1920           0 :   case e_NONE:     return "e_NONE";
    1921           0 :   case e_NOTFUNC:  return "e_NOTFUNC";
    1922           0 :   case e_OP:       return "e_OP";
    1923           0 :   case e_OVERFLOW: return "e_OVERFLOW";
    1924           0 :   case e_PACKAGE:  return "e_PACKAGE";
    1925           0 :   case e_PREC:     return "e_PREC";
    1926           0 :   case e_PRIME:    return "e_PRIME";
    1927          49 :   case e_PRIORITY: return "e_PRIORITY";
    1928           0 :   case e_ROOTS0:   return "e_ROOTS0";
    1929           0 :   case e_SQRTN:    return "e_SQRTN";
    1930           0 :   case e_STACK:    return "e_STACK";
    1931           0 :   case e_SYNTAX:   return "e_SYNTAX";
    1932           0 :   case e_STACKTHREAD:   return "e_STACKTHREAD";
    1933           0 :   case e_TYPE2:    return "e_TYPE2";
    1934        7847 :   case e_TYPE:     return "e_TYPE";
    1935           0 :   case e_USER:     return "e_USER";
    1936           0 :   case e_VAR:      return "e_VAR";
    1937             :   }
    1938           0 :   return "invalid error number";
    1939             : }
    1940             : 
    1941             : long
    1942           0 : name_numerr(const char *s)
    1943             : {
    1944           0 :   if (!strcmp(s,"e_ALARM"))    return e_ALARM;
    1945           0 :   if (!strcmp(s,"e_ARCH"))     return e_ARCH;
    1946           0 :   if (!strcmp(s,"e_BUG"))      return e_BUG;
    1947           0 :   if (!strcmp(s,"e_COMPONENT")) return e_COMPONENT;
    1948           0 :   if (!strcmp(s,"e_CONSTPOL")) return e_CONSTPOL;
    1949           0 :   if (!strcmp(s,"e_COPRIME"))  return e_COPRIME;
    1950           0 :   if (!strcmp(s,"e_DIM"))      return e_DIM;
    1951           0 :   if (!strcmp(s,"e_DOMAIN"))   return e_DOMAIN;
    1952           0 :   if (!strcmp(s,"e_FILE"))     return e_FILE;
    1953           0 :   if (!strcmp(s,"e_FILEDESC")) return e_FILEDESC;
    1954           0 :   if (!strcmp(s,"e_FLAG"))     return e_FLAG;
    1955           0 :   if (!strcmp(s,"e_IMPL"))     return e_IMPL;
    1956           0 :   if (!strcmp(s,"e_INV"))      return e_INV;
    1957           0 :   if (!strcmp(s,"e_IRREDPOL")) return e_IRREDPOL;
    1958           0 :   if (!strcmp(s,"e_MAXPRIME")) return e_MAXPRIME;
    1959           0 :   if (!strcmp(s,"e_MEM"))      return e_MEM;
    1960           0 :   if (!strcmp(s,"e_MISC"))     return e_MISC;
    1961           0 :   if (!strcmp(s,"e_MODULUS"))  return e_MODULUS;
    1962           0 :   if (!strcmp(s,"e_NONE"))     return e_NONE;
    1963           0 :   if (!strcmp(s,"e_NOTFUNC"))  return e_NOTFUNC;
    1964           0 :   if (!strcmp(s,"e_OP"))       return e_OP;
    1965           0 :   if (!strcmp(s,"e_OVERFLOW")) return e_OVERFLOW;
    1966           0 :   if (!strcmp(s,"e_PACKAGE"))  return e_PACKAGE;
    1967           0 :   if (!strcmp(s,"e_PREC"))     return e_PREC;
    1968           0 :   if (!strcmp(s,"e_PRIME"))    return e_PRIME;
    1969           0 :   if (!strcmp(s,"e_PRIORITY")) return e_PRIORITY;
    1970           0 :   if (!strcmp(s,"e_ROOTS0"))   return e_ROOTS0;
    1971           0 :   if (!strcmp(s,"e_SQRTN"))    return e_SQRTN;
    1972           0 :   if (!strcmp(s,"e_STACK"))    return e_STACK;
    1973           0 :   if (!strcmp(s,"e_SYNTAX"))   return e_SYNTAX;
    1974           0 :   if (!strcmp(s,"e_TYPE"))     return e_TYPE;
    1975           0 :   if (!strcmp(s,"e_TYPE2"))    return e_TYPE2;
    1976           0 :   if (!strcmp(s,"e_USER"))     return e_USER;
    1977           0 :   if (!strcmp(s,"e_VAR"))      return e_VAR;
    1978           0 :   pari_err(e_MISC,"unknown error name");
    1979             :   return -1; /* LCOV_EXCL_LINE */
    1980             : }
    1981             : 
    1982             : GEN
    1983       27109 : errname(GEN err)
    1984             : {
    1985       27109 :   if (typ(err)!=t_ERROR) pari_err_TYPE("errname",err);
    1986       27109 :   return strtoGENstr(numerr_name(err_get_num(err)));
    1987             : }
    1988             : 
    1989             : /*******************************************************************/
    1990             : /*                                                                */
    1991             : /*                       CLONING & COPY                            */
    1992             : /*                  Replicate an existing GEN                      */
    1993             : /*                                                                 */
    1994             : /*******************************************************************/
    1995             : /* lontyp[tx] = 0 (non recursive type) or number of codewords for type tx */
    1996             : const  long lontyp[] = { 0,0,0,1,1,2,1,2,1,1, 2,2,0,1,1,1,1,1,1,1, 2,0,0,2,2,1 };
    1997             : 
    1998             : static GEN
    1999        7019 : list_internal_copy(GEN z, long nmax)
    2000             : {
    2001             :   long i, l;
    2002             :   GEN a;
    2003        7019 :   if (!z) return NULL;
    2004        1336 :   l = lg(z);
    2005        1336 :   a = newblock(nmax+1);
    2006       50018 :   for (i = 1; i < l; i++) gel(a,i) = gel(z,i)? gclone(gel(z,i)): gen_0;
    2007        1336 :   a[0] = z[0]; setisclone(a); return a;
    2008             : }
    2009             : 
    2010             : static void
    2011        7019 : listassign(GEN x, GEN y)
    2012             : {
    2013        7019 :   long nmax = list_nmax(x);
    2014        7019 :   GEN L = list_data(x);
    2015        7019 :   if (!nmax && L) nmax = lg(L) + 32; /* not malloc'ed yet */
    2016        7019 :   y[1] = evaltyp(list_typ(x))|evallg(nmax);
    2017        7019 :   list_data(y) = list_internal_copy(L, nmax);
    2018        7019 : }
    2019             : 
    2020             : /* transform a non-malloced list (e.g. from gtolist or gtomap) to a malloced
    2021             :  * list suitable for listput */
    2022             : GEN
    2023           0 : listinit(GEN x)
    2024             : {
    2025           0 :   GEN y = cgetg(3, t_LIST);
    2026           0 :   listassign(x, y); return y;
    2027             : }
    2028             : 
    2029             : /* copy list on the PARI stack */
    2030             : GEN
    2031         730 : listcopy(GEN x)
    2032             : {
    2033         730 :   GEN y = mklist(), L = list_data(x);
    2034         730 :   if (L) list_data(y) = gcopy(L);
    2035         730 :   y[1] = evaltyp(list_typ(x));
    2036         730 :   return y;
    2037             : }
    2038             : 
    2039             : GEN
    2040  4510962226 : gcopy(GEN x)
    2041             : {
    2042  4510962226 :   long tx = typ(x), lx, i;
    2043             :   GEN y;
    2044  4510962226 :   switch(tx)
    2045             :   { /* non recursive types */
    2046  3770334958 :     case t_INT: return signe(x)? icopy(x): gen_0;
    2047   482721482 :     case t_REAL:
    2048             :     case t_STR:
    2049   482721482 :     case t_VECSMALL: return leafcopy(x);
    2050             :     /* one more special case */
    2051         730 :     case t_LIST: return listcopy(x);
    2052             :   }
    2053   257905056 :   y = cgetg_copy(x, &lx);
    2054   258253478 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    2055  1089721067 :   for (; i<lx; i++) gel(y,i) = gcopy(gel(x,i));
    2056   258253356 :   return y;
    2057             : }
    2058             : 
    2059             : /* as gcopy, but truncate to the first lx components if recursive type
    2060             :  * [ leaves use their own lg ]. No checks. */
    2061             : GEN
    2062         742 : gcopy_lg(GEN x, long lx)
    2063             : {
    2064         742 :   long tx = typ(x), i;
    2065             :   GEN y;
    2066         742 :   switch(tx)
    2067             :   { /* non recursive types */
    2068           0 :     case t_INT: return signe(x)? icopy(x): gen_0;
    2069           0 :     case t_REAL:
    2070             :     case t_STR:
    2071           0 :     case t_VECSMALL: return leafcopy(x);
    2072             :     /* one more special case */
    2073           0 :     case t_LIST: return listcopy(x);
    2074             :   }
    2075         742 :   y = cgetg(lx, tx);
    2076         742 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    2077        2051 :   for (; i<lx; i++) gel(y,i) = gcopy(gel(x,i));
    2078         742 :   return y;
    2079             : }
    2080             : 
    2081             : /* cf cgetg_copy: "allocate" (by updating first codeword only) for subsequent
    2082             :  * copy of x, as if avma = *AVMA */
    2083             : INLINE GEN
    2084  1000438333 : cgetg_copy_avma(GEN x, long *plx, pari_sp *AVMA) {
    2085             :   GEN z;
    2086  1000438333 :   *plx = lg(x);
    2087  1000438333 :   z = ((GEN)*AVMA) - *plx;
    2088  1000438333 :   z[0] = x[0] & (TYPBITS|LGBITS);
    2089  1000438333 :   *AVMA = (pari_sp)z; return z;
    2090             : }
    2091             : INLINE GEN
    2092         620 : cgetlist_avma(pari_sp *AVMA)
    2093             : {
    2094         620 :   GEN y = ((GEN)*AVMA) - 3;
    2095         620 :   y[0] = _evallg(3) | evaltyp(t_LIST);
    2096         620 :   *AVMA = (pari_sp)y; return y;
    2097             : }
    2098             : 
    2099             : /* copy x as if avma = *AVMA, update *AVMA */
    2100             : GEN
    2101  3277163919 : gcopy_avma(GEN x, pari_sp *AVMA)
    2102             : {
    2103  3277163919 :   long i, lx, tx = typ(x);
    2104             :   GEN y;
    2105             : 
    2106  3277163919 :   switch(typ(x))
    2107             :   { /* non recursive types */
    2108  3070887088 :     case t_INT:
    2109  3070887088 :       if (lgefint(x) == 2) return gen_0;
    2110  2554853245 :       *AVMA = (pari_sp)icopy_avma(x, *AVMA);
    2111  2554853570 :       return (GEN)*AVMA;
    2112    62383936 :     case t_REAL: case t_STR: case t_VECSMALL:
    2113    62383936 :       *AVMA = (pari_sp)leafcopy_avma(x, *AVMA);
    2114    62222015 :       return (GEN)*AVMA;
    2115             : 
    2116             :     /* one more special case */
    2117         616 :     case t_LIST:
    2118         616 :       y = cgetlist_avma(AVMA);
    2119         616 :       listassign(x, y); return y;
    2120             : 
    2121             :   }
    2122   143892279 :   y = cgetg_copy_avma(x, &lx, AVMA);
    2123   143939217 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    2124   693327246 :   for (; i<lx; i++) gel(y,i) = gcopy_avma(gel(x,i), AVMA);
    2125   143942287 :   return y;
    2126             : }
    2127             : 
    2128             : /* [copy_bin/bin_copy:] same as gcopy_avma but use NULL to code an exact 0, and
    2129             :  * make shallow copies of finalized t_LISTs */
    2130             : static GEN
    2131  3860215076 : gcopy_av0(GEN x, pari_sp *AVMA)
    2132             : {
    2133  3860215076 :   long i, lx, tx = typ(x);
    2134             :   GEN y;
    2135             : 
    2136  3860215076 :   switch(tx)
    2137             :   { /* non recursive types */
    2138  2534886096 :     case t_INT:
    2139  2534886096 :       if (!signe(x)) return NULL; /* special marker */
    2140  1660892785 :       *AVMA = (pari_sp)icopy_avma(x, *AVMA);
    2141  1661098896 :       return (GEN)*AVMA;
    2142          49 :     case t_LIST:
    2143          49 :       if (list_data(x) && !list_nmax(x)) break; /* not finalized, need copy */
    2144             :       /* else finalized: shallow copy */
    2145             :     case t_REAL: case t_STR: case t_VECSMALL:
    2146   469568457 :       *AVMA = (pari_sp)leafcopy_avma(x, *AVMA);
    2147   469561948 :       return (GEN)*AVMA;
    2148             :   }
    2149   855760523 :   y = cgetg_copy_avma(x, &lx, AVMA);
    2150   856507336 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    2151  4259547010 :   for (; i<lx; i++) gel(y,i) = gcopy_av0(gel(x,i), AVMA);
    2152   856447061 :   return y;
    2153             : }
    2154             : 
    2155             : INLINE GEN
    2156          12 : icopy_avma_canon(GEN x, pari_sp AVMA)
    2157             : {
    2158          12 :   long i, lx = lgefint(x);
    2159          12 :   GEN y = ((GEN)AVMA) - lx;
    2160          12 :   y[0] = evaltyp(t_INT)|_evallg(lx); /* kills isclone */
    2161          12 :   y[1] = x[1]; x = int_MSW(x);
    2162          24 :   for (i=2; i<lx; i++, x = int_precW(x)) y[i] = *x;
    2163          12 :   return y;
    2164             : }
    2165             : 
    2166             : /* [copy_bin_canon:] same as gcopy_av0, but copy integers in
    2167             :  * canonical (native kernel) form and make a full copy of t_LISTs */
    2168             : static GEN
    2169          32 : gcopy_av0_canon(GEN x, pari_sp *AVMA)
    2170             : {
    2171          32 :   long i, lx, tx = typ(x);
    2172             :   GEN y;
    2173             : 
    2174          32 :   switch(tx)
    2175             :   { /* non recursive types */
    2176          20 :     case t_INT:
    2177          20 :       if (!signe(x)) return NULL; /* special marker */
    2178          12 :       *AVMA = (pari_sp)icopy_avma_canon(x, *AVMA);
    2179          12 :       return (GEN)*AVMA;
    2180           0 :     case t_REAL: case t_STR: case t_VECSMALL:
    2181           0 :       *AVMA = (pari_sp)leafcopy_avma(x, *AVMA);
    2182           0 :       return (GEN)*AVMA;
    2183             : 
    2184             :     /* one more special case */
    2185           4 :     case t_LIST:
    2186             :     {
    2187           4 :       long t = list_typ(x);
    2188           4 :       GEN y = cgetlist_avma(AVMA), z = list_data(x);
    2189           4 :       if (z) {
    2190           0 :         list_data(y) = gcopy_av0_canon(z, AVMA);
    2191           0 :         y[1] = evaltyp(t)|_evallg(lg(z)-1);
    2192             :       } else {
    2193           4 :         list_data(y) = NULL;
    2194           4 :         y[1] = evaltyp(t);
    2195             :       }
    2196           4 :       return y;
    2197             :     }
    2198             :   }
    2199           8 :   y = cgetg_copy_avma(x, &lx, AVMA);
    2200           8 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    2201          24 :   for (; i<lx; i++) gel(y,i) = gcopy_av0_canon(gel(x,i), AVMA);
    2202           8 :   return y;
    2203             : }
    2204             : 
    2205             : /* [copy_bin/bin_copy:] size (number of words) required for
    2206             :  * gcopy_av0_canon(x) */
    2207             : static long
    2208          32 : taille0_canon(GEN x)
    2209             : {
    2210          32 :   long i,n,lx, tx = typ(x);
    2211          32 :   switch(tx)
    2212             :   { /* non recursive types */
    2213          20 :     case t_INT: return signe(x)? lgefint(x): 0;
    2214           0 :     case t_REAL:
    2215             :     case t_STR:
    2216           0 :     case t_VECSMALL: return lg(x);
    2217             : 
    2218             :     /* one more special case */
    2219           4 :     case t_LIST:
    2220             :     {
    2221           4 :       GEN L = list_data(x);
    2222           4 :       return L? 3 + taille0_canon(L): 3;
    2223             :     }
    2224             :   }
    2225           8 :   n = lx = lg(x);
    2226          24 :   for (i=lontyp[tx]; i<lx; i++) n += taille0_canon(gel(x,i));
    2227           8 :   return n;
    2228             : }
    2229             : 
    2230             : /* [copy_bin/bin_copy:] size (number of words) required for gcopy_av0(x) */
    2231             : static long
    2232  3860436343 : taille0(GEN x)
    2233             : {
    2234  3860436343 :   long i,n,lx, tx = typ(x);
    2235  3860436343 :   switch(tx)
    2236             :   { /* non recursive types */
    2237  2535132606 :     case t_INT:
    2238  2535132606 :       lx = lgefint(x);
    2239  2535132606 :       return lx == 2? 0: lx;
    2240          49 :     case t_LIST:
    2241             :     {
    2242          49 :       GEN L = list_data(x);
    2243          49 :       if (L && !list_nmax(x)) break; /* not finalized, deep copy */
    2244             :     }
    2245             :     /* else finalized: shallow */
    2246             :     case t_REAL:
    2247             :     case t_STR:
    2248             :     case t_VECSMALL:
    2249   469565692 :       return lg(x);
    2250             :   }
    2251   855738045 :   n = lx = lg(x);
    2252  4259172094 :   for (i=lontyp[tx]; i<lx; i++) n += taille0(gel(x,i));
    2253   856062012 :   return n;
    2254             : }
    2255             : 
    2256             : static long
    2257  3363837079 : gsizeclone_i(GEN x)
    2258             : {
    2259  3363837079 :   long i,n,lx, tx = typ(x);
    2260  3363837079 :   switch(tx)
    2261             :   { /* non recursive types */
    2262  3070866306 :     case t_INT: lx = lgefint(x); return lx == 2? 0: lx;;
    2263    74632913 :     case t_REAL:
    2264             :     case t_STR:
    2265    74632913 :     case t_VECSMALL: return lg(x);
    2266             : 
    2267        7019 :     case t_LIST: return 3;
    2268   218330841 :     default:
    2269   218330841 :       n = lx = lg(x);
    2270  3495103113 :       for (i=lontyp[tx]; i<lx; i++) n += gsizeclone_i(gel(x,i));
    2271   218383469 :       return n;
    2272             :   }
    2273             : }
    2274             : 
    2275             : /* #words needed to clone x; t_LIST is a special case since list_data() is
    2276             :  * malloc'ed later, in list_internal_copy() */
    2277             : static long
    2278   231592676 : gsizeclone(GEN x) { return (typ(x) == t_INT)? lgefint(x): gsizeclone_i(x); }
    2279             : 
    2280             : long
    2281     2240315 : gsizeword(GEN x)
    2282             : {
    2283     2240315 :   long i, n, lx, tx = typ(x);
    2284     2240315 :   switch(tx)
    2285             :   { /* non recursive types */
    2286     1700433 :     case t_INT:
    2287             :     case t_REAL:
    2288             :     case t_STR:
    2289     1700433 :     case t_VECSMALL: return lg(x);
    2290             : 
    2291           7 :     case t_LIST:
    2292           7 :       x = list_data(x);
    2293           7 :       return x? 3 + gsizeword(x): 3;
    2294             : 
    2295      539875 :     default:
    2296      539875 :       n = lx = lg(x);
    2297     2778167 :       for (i=lontyp[tx]; i<lx; i++) n += gsizeword(gel(x,i));
    2298      539875 :       return n;
    2299             :   }
    2300             : }
    2301             : long
    2302         168 : gsizebyte(GEN x) { return gsizeword(x) * sizeof(long); }
    2303             : 
    2304             : /* return a clone of x structured as a gcopy */
    2305             : GENbin*
    2306   457242452 : copy_bin(GEN x)
    2307             : {
    2308   457242452 :   long t = taille0(x);
    2309   457343803 :   GENbin *p = (GENbin*)pari_malloc(sizeof(GENbin) + t*sizeof(long));
    2310   457392082 :   pari_sp AVMA = (pari_sp)(GENbinbase(p) + t);
    2311   457373065 :   p->rebase = &shiftaddress;
    2312   457373065 :   p->len = t;
    2313   457373065 :   p->x   = gcopy_av0(x, &AVMA);
    2314   457335893 :   p->base= (GEN)AVMA; return p;
    2315             : }
    2316             : 
    2317             : /* same, writing t_INT in canonical native form */
    2318             : GENbin*
    2319          16 : copy_bin_canon(GEN x)
    2320             : {
    2321          16 :   long t = taille0_canon(x);
    2322          16 :   GENbin *p = (GENbin*)pari_malloc(sizeof(GENbin) + t*sizeof(long));
    2323          16 :   pari_sp AVMA = (pari_sp)(GENbinbase(p) + t);
    2324          16 :   p->rebase = &shiftaddress_canon;
    2325          16 :   p->len = t;
    2326          16 :   p->x   = gcopy_av0_canon(x, &AVMA);
    2327          16 :   p->base= (GEN)AVMA; return p;
    2328             : }
    2329             : 
    2330             : GEN
    2331   231592861 : gclone(GEN x)
    2332             : {
    2333   231592861 :   long i,lx,tx = typ(x), t = gsizeclone(x);
    2334   231591819 :   GEN y = newblock(t);
    2335   231594970 :   switch(tx)
    2336             :   { /* non recursive types */
    2337   144551800 :     case t_INT:
    2338   144551800 :       lx = lgefint(x);
    2339   144551800 :       y[0] = evaltyp(t_INT)|_evallg(lx);
    2340   876421898 :       for (i=1; i<lx; i++) y[i] = x[i];
    2341   144551800 :       break;
    2342    12466796 :     case t_REAL:
    2343             :     case t_STR:
    2344             :     case t_VECSMALL:
    2345    12466796 :       lx = lg(x);
    2346   132184531 :       for (i=0; i<lx; i++) y[i] = x[i];
    2347    12466796 :       break;
    2348             : 
    2349             :     /* one more special case */
    2350        6403 :     case t_LIST:
    2351        6403 :       y[0] = evaltyp(t_LIST)|_evallg(3);
    2352        6403 :       listassign(x, y);
    2353        6403 :       break;
    2354    74569971 :     default: {
    2355    74569971 :       pari_sp AVMA = (pari_sp)(y + t);
    2356    74569971 :       lx = lg(x);
    2357    74569971 :       y[0] = x[0];
    2358    74569971 :       if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    2359  2802381868 :       for (; i<lx; i++) gel(y,i) = gcopy_avma(gel(x,i), &AVMA);
    2360             :     }
    2361             :   }
    2362   231594212 :   setisclone(y); return y;
    2363             : }
    2364             : 
    2365             : void
    2366  2987461449 : shiftaddress(GEN x, long dec)
    2367             : {
    2368  2987461449 :   long i, lx, tx = typ(x);
    2369  2987461449 :   if (is_recursive_t(tx))
    2370             :   {
    2371   856717404 :     if (tx == t_LIST)
    2372             :     {
    2373          49 :       if (!list_data(x) || list_nmax(x)) return; /* empty or finalized */
    2374             :       /* not finalized, update pointers  */
    2375             :     }
    2376   856717362 :     lx = lg(x);
    2377  4260888159 :     for (i=lontyp[tx]; i<lx; i++) {
    2378  3404304645 :       if (!x[i]) gel(x,i) = gen_0;
    2379             :       else
    2380             :       {
    2381  2556977821 :         x[i] += dec;
    2382  2556977821 :         shiftaddress(gel(x,i), dec);
    2383             :       }
    2384             :     }
    2385             :   }
    2386             : }
    2387             : 
    2388             : void
    2389          24 : shiftaddress_canon(GEN x, long dec)
    2390             : {
    2391          24 :   long i, lx, tx = typ(x);
    2392          24 :   switch(tx)
    2393             :   { /* non recursive types */
    2394          12 :     case t_INT: {
    2395             :       GEN y;
    2396          12 :       lx = lgefint(x); if (lx <= 3) return;
    2397           0 :       y = x + 2;
    2398           0 :       x = int_MSW(x);  if (x == y) return;
    2399           0 :       while (x > y) { lswap(*x, *y); x = int_precW(x); y++; }
    2400           0 :       break;
    2401             :     }
    2402           0 :     case t_REAL:
    2403             :     case t_STR:
    2404             :     case t_VECSMALL:
    2405           0 :       break;
    2406             : 
    2407             :     /* one more special case */
    2408           4 :     case t_LIST: {
    2409           4 :       GEN Lx = list_data(x);
    2410           4 :       if (Lx) {
    2411           0 :         GEN L = (GEN)((long)Lx+dec);
    2412           0 :         shiftaddress_canon(L, dec);
    2413           0 :         list_data(x) = gcopy(L);
    2414             :       }
    2415           4 :       break;
    2416             :     }
    2417           8 :     default: /* Fall through */
    2418           8 :       lx = lg(x);
    2419          24 :       for (i=lontyp[tx]; i<lx; i++) {
    2420          16 :         if (!x[i]) gel(x,i) = gen_0;
    2421             :         else
    2422             :         {
    2423           8 :           x[i] += dec;
    2424           8 :           shiftaddress_canon(gel(x,i), dec);
    2425             :         }
    2426             :       }
    2427             :   }
    2428             : }
    2429             : 
    2430             : /********************************************************************/
    2431             : /**                                                                **/
    2432             : /**                INSERT DYNAMIC OBJECT IN STRUCTURE              **/
    2433             : /**                                                                **/
    2434             : /********************************************************************/
    2435             : GEN
    2436          42 : obj_reinit(GEN S)
    2437             : {
    2438          42 :   GEN s, T = leafcopy(S);
    2439          42 :   long a = lg(T)-1;
    2440          42 :   s = gel(T,a);
    2441          42 :   gel(T,a) = zerovec(lg(s)-1);
    2442          42 :   return T;
    2443             : }
    2444             : 
    2445             : GEN
    2446     1622451 : obj_init(long d, long n)
    2447             : {
    2448     1622451 :   GEN S = cgetg(d+2, t_VEC);
    2449     1622451 :   gel(S, d+1) = zerovec(n);
    2450     1622451 :   return S;
    2451             : }
    2452             : 
    2453             : /* as obj_insert. WITHOUT cloning (for libpari, when creating a *new* obj
    2454             :  * from an existing one) */
    2455             : GEN
    2456     1748362 : obj_insert(GEN S, long K, GEN O)
    2457             : {
    2458     1748362 :   GEN o, v = veclast(S);
    2459     1748362 :   if (typ(v) != t_VEC) pari_err_TYPE("obj_insert", S);
    2460     1748362 :   if (!check_clone(v))
    2461             :   {
    2462           0 :     if (DEBUGLEVEL) pari_warn(warner,"trying to update parent object");
    2463           0 :     return gclone(O);
    2464             :   }
    2465     1748362 :   o = gel(v,K);
    2466     1748362 :   gel(v,K) = gclone(O); /*SIGINT: before unclone(o)*/
    2467     1748362 :   if (isclone(o)) gunclone(o); return gel(v,K);
    2468             : }
    2469             : 
    2470             : GEN
    2471      190085 : obj_insert_shallow(GEN S, long K, GEN O)
    2472             : {
    2473      190085 :   GEN v = veclast(S);
    2474      190085 :   if (typ(v) != t_VEC) pari_err_TYPE("obj_insert", S);
    2475      190085 :   gel(v,K) = O;
    2476      190085 :   return gel(v,K);
    2477             : }
    2478             : 
    2479             : /* Does S [last position] contain data at position K ? Return it, or NULL */
    2480             : GEN
    2481     4569903 : obj_check(GEN S, long K)
    2482             : {
    2483     4569903 :   GEN O, v = veclast(S);
    2484     4569900 :   if (typ(v) != t_VEC || K >= lg(v)) pari_err_TYPE("obj_check", S);
    2485     4569901 :   O = gel(v,K); return isintzero(O)? NULL: O;
    2486             : }
    2487             : 
    2488             : GEN
    2489     1032164 : obj_checkbuild(GEN S, long tag, GEN (*build)(GEN))
    2490             : {
    2491     1032164 :   GEN O = obj_check(S, tag);
    2492     1032164 :   if (!O)
    2493      896469 :   { pari_sp av = avma; O = obj_insert(S, tag, build(S)); set_avma(av); }
    2494     1032157 :   return O;
    2495             : }
    2496             : 
    2497             : GEN
    2498      530202 : obj_checkbuild_prec(GEN S, long tag, GEN (*build)(GEN,long),
    2499             :   long (*pr)(GEN), long prec)
    2500             : {
    2501      530202 :   pari_sp av = avma;
    2502      530202 :   GEN w = obj_check(S, tag);
    2503      530202 :   if (!w || pr(w) < prec) w = obj_insert(S, tag, build(S, prec));
    2504      530202 :   set_avma(av); return gcopy(w);
    2505             : }
    2506             : GEN
    2507      398512 : obj_checkbuild_realprec(GEN S, long tag, GEN (*build)(GEN,long), long prec)
    2508      398512 : { return obj_checkbuild_prec(S,tag,build,gprecision,prec); }
    2509             : GEN
    2510         546 : obj_checkbuild_padicprec(GEN S, long tag, GEN (*build)(GEN,long), long prec)
    2511         546 : { return obj_checkbuild_prec(S,tag,build,padicprec_relative,prec); }
    2512             : 
    2513             : /* Reset S [last position], freeing all clones */
    2514             : void
    2515      186354 : obj_free(GEN S)
    2516             : {
    2517      186354 :   GEN v = veclast(S);
    2518             :   long i;
    2519      186354 :   if (typ(v) != t_VEC) pari_err_TYPE("obj_free", S);
    2520     1206310 :   for (i = 1; i < lg(v); i++)
    2521             :   {
    2522     1019956 :     GEN o = gel(v,i);
    2523     1019956 :     gel(v,i) = gen_0;
    2524     1019956 :     gunclone_deep(o);
    2525             :   }
    2526      186354 : }
    2527             : 
    2528             : /*******************************************************************/
    2529             : /*                                                                 */
    2530             : /*                         STACK MANAGEMENT                        */
    2531             : /*                                                                 */
    2532             : /*******************************************************************/
    2533             : INLINE void
    2534  5039926728 : dec_gerepile(pari_sp *x, pari_sp av0, pari_sp av, pari_sp tetpil, size_t dec)
    2535             : {
    2536  5039926728 :   if (*x < av && *x >= av0)
    2537             :   { /* update address if in stack */
    2538  4420821033 :     if (*x < tetpil) *x += dec;
    2539        3126 :     else pari_err_BUG("gerepile, significant pointers lost");
    2540             :   }
    2541  5039923602 : }
    2542             : 
    2543             : void
    2544     6326237 : gerepileallsp(pari_sp av, pari_sp tetpil, int n, ...)
    2545             : {
    2546     6326237 :   const pari_sp av0 = avma;
    2547     6326237 :   const size_t dec = av-tetpil;
    2548             :   int i;
    2549     6326237 :   va_list a; va_start(a, n);
    2550     6326237 :   (void)gerepile(av,tetpil,NULL);
    2551    29597796 :   for (i=0; i<n; i++) dec_gerepile((pari_sp*)va_arg(a,GEN*), av0,av,tetpil,dec);
    2552     6326935 :   va_end(a);
    2553     6326935 : }
    2554             : 
    2555             : /* Takes an array of pointers to GENs, of length n.
    2556             :  * Cleans up the stack between av and tetpil, updating those GENs. */
    2557             : void
    2558    21423505 : gerepilemanysp(pari_sp av, pari_sp tetpil, GEN* gptr[], int n)
    2559             : {
    2560    21423505 :   const pari_sp av0 = avma;
    2561    21423505 :   const size_t dec = av-tetpil;
    2562             :   int i;
    2563    21423505 :   (void)gerepile(av,tetpil,NULL);
    2564    64314458 :   for (i=0; i<n; i++) dec_gerepile((pari_sp*)gptr[i], av0, av, tetpil, dec);
    2565    21427082 : }
    2566             : 
    2567             : /* Takes an array of GENs (cast to longs), of length n.
    2568             :  * Cleans up the stack between av and tetpil, updating those GENs. */
    2569             : void
    2570   307250777 : gerepilecoeffssp(pari_sp av, pari_sp tetpil, long *g, int n)
    2571             : {
    2572   307250777 :   const pari_sp av0 = avma;
    2573   307250777 :   const size_t dec = av-tetpil;
    2574             :   int i;
    2575   307250777 :   (void)gerepile(av,tetpil,NULL);
    2576   922298896 :   for (i=0; i<n; i++,g++) dec_gerepile((pari_sp*)g, av0, av, tetpil, dec);
    2577   307444760 : }
    2578             : 
    2579             : static int
    2580           0 : dochk_gerepileupto(GEN av, GEN x)
    2581             : {
    2582             :   long i,lx,tx;
    2583           0 :   if (!isonstack(x)) return 1;
    2584           0 :   if (x > av)
    2585             :   {
    2586           0 :     pari_warn(warner,"bad object %Ps",x);
    2587           0 :     return 0;
    2588             :   }
    2589           0 :   tx = typ(x);
    2590           0 :   if (! is_recursive_t(tx)) return 1;
    2591             : 
    2592           0 :   lx = lg(x);
    2593           0 :   for (i=lontyp[tx]; i<lx; i++)
    2594           0 :     if (!dochk_gerepileupto(av, gel(x,i)))
    2595             :     {
    2596           0 :       pari_warn(warner,"bad component %ld in object %Ps",i,x);
    2597           0 :       return 0;
    2598             :     }
    2599           0 :   return 1;
    2600             : }
    2601             : /* check that x and all its components are out of stack, or have been
    2602             :  * created after av */
    2603             : int
    2604           0 : chk_gerepileupto(GEN x) { return dochk_gerepileupto(x, x); }
    2605             : 
    2606             : /* print stack between avma & av */
    2607             : void
    2608           0 : dbg_gerepile(pari_sp av)
    2609             : {
    2610           0 :   GEN x = (GEN)avma;
    2611           0 :   while (x < (GEN)av)
    2612             :   {
    2613           0 :     const long tx = typ(x), lx = lg(x);
    2614             :     GEN *a;
    2615             : 
    2616           0 :     pari_printf(" [%ld] %Ps:", x - (GEN)avma, x);
    2617           0 :     if (! is_recursive_t(tx)) { pari_putc('\n'); x += lx; continue; }
    2618           0 :     a = (GEN*)x + lontyp[tx]; x += lx;
    2619           0 :     for (  ; a < (GEN*)x; a++)
    2620             :     {
    2621           0 :       if (*a == gen_0)
    2622           0 :         pari_puts("  gen_0");
    2623           0 :       else if (*a == gen_1)
    2624           0 :         pari_puts("  gen_1");
    2625           0 :       else if (*a == gen_m1)
    2626           0 :         pari_puts("  gen_m1");
    2627           0 :       else if (*a == gen_2)
    2628           0 :         pari_puts("  gen_2");
    2629           0 :       else if (*a == gen_m2)
    2630           0 :         pari_puts("  gen_m2");
    2631           0 :       else if (*a == ghalf)
    2632           0 :         pari_puts("  ghalf");
    2633           0 :       else if (isclone(*a))
    2634           0 :         pari_printf("  %Ps (clone)", *a);
    2635             :       else
    2636           0 :         pari_printf("  %Ps [%ld]", *a, *a - (GEN)avma);
    2637           0 :       if (a+1 < (GEN*)x) pari_putc(',');
    2638             :     }
    2639           0 :     pari_printf("\n");
    2640             :   }
    2641           0 : }
    2642             : void
    2643           0 : dbg_gerepileupto(GEN q)
    2644             : {
    2645           0 :   err_printf("%Ps:\n", q);
    2646           0 :   dbg_gerepile((pari_sp) (q+lg(q)));
    2647           0 : }
    2648             : 
    2649             : GEN
    2650  1212076050 : gerepile(pari_sp av, pari_sp tetpil, GEN q)
    2651             : {
    2652  1212076050 :   const size_t dec = av - tetpil;
    2653  1212076050 :   const pari_sp av0 = avma;
    2654             :   GEN x, a;
    2655             : 
    2656  1212076050 :   if (dec == 0) return q;
    2657  1016057825 :   if ((long)dec < 0) pari_err(e_MISC,"lbot>ltop in gerepile");
    2658             : 
    2659             :   /* dec_gerepile(&q, av0, av, tetpil, dec), saving 1 comparison */
    2660  1016487034 :   if (q >= (GEN)av0 && q < (GEN)tetpil)
    2661   688410993 :     q = (GEN) (((pari_sp)q) + dec);
    2662             : 
    2663 27763453738 :   for (x = (GEN)av, a = (GEN)tetpil; a > (GEN)av0; ) *--x = *--a;
    2664  1016487034 :   set_avma((pari_sp)x);
    2665  6205363272 :   while (x < (GEN)av)
    2666             :   {
    2667  5187977451 :     const long tx = typ(x), lx = lg(x);
    2668             : 
    2669  5187977451 :     if (! is_recursive_t(tx)) { x += lx; continue; }
    2670  1075310852 :     a = x + lontyp[tx]; x += lx;
    2671  5434898956 :     for (  ; a < x; a++) dec_gerepile((pari_sp*)a, av0, av, tetpil, dec);
    2672             :   }
    2673  1017385821 :   return q;
    2674             : }
    2675             : 
    2676             : void
    2677           0 : dbg_fill_stack(void)
    2678             : {
    2679             : #ifdef LONG_IS_64BIT
    2680           0 :   const long JUNK = 0xBADC0FFEE0DDF00D;
    2681             : #else
    2682           0 :   const long JUNK = 0xDEADBEEF;
    2683             : #endif
    2684           0 :   GEN x = ((GEN)pari_mainstack->bot);
    2685           0 :   while (x < (GEN)avma) *x++ = JUNK;
    2686           0 : }
    2687             : 
    2688             : void
    2689           0 : debug_stack(void)
    2690             : {
    2691           0 :   pari_sp top = pari_mainstack->top, bot = pari_mainstack->bot;
    2692             :   GEN z;
    2693           0 :   err_printf("bot=0x%lx\ttop=0x%lx\tavma=0x%lx\n", bot, top, avma);
    2694           0 :   for (z = ((GEN)top)-1; z >= (GEN)avma; z--)
    2695           0 :     err_printf("%p:\t0x%lx\t%lu\n",z,*z,*z);
    2696           0 : }
    2697             : 
    2698             : void
    2699      324248 : setdebugvar(long n) { DEBUGVAR=n; }
    2700             : 
    2701             : long
    2702      324508 : getdebugvar(void) { return DEBUGVAR; }
    2703             : 
    2704             : long
    2705           7 : getstack(void) { return pari_mainstack->top-avma; }
    2706             : 
    2707             : /*******************************************************************/
    2708             : /*                                                                 */
    2709             : /*                               timer_delay                             */
    2710             : /*                                                                 */
    2711             : /*******************************************************************/
    2712             : 
    2713             : #if defined(USE_CLOCK_GETTIME)
    2714             : #if defined(_POSIX_THREAD_CPUTIME)
    2715             : static THREAD clockid_t time_type = CLOCK_THREAD_CPUTIME_ID;
    2716             : #else
    2717             : static const THREAD clockid_t time_type = CLOCK_PROCESS_CPUTIME_ID;
    2718             : #endif
    2719             : static void
    2720             : pari_init_timer(void)
    2721             : {
    2722             : #if defined(_POSIX_THREAD_CPUTIME)
    2723             :   time_type = CLOCK_PROCESS_CPUTIME_ID;
    2724             : #endif
    2725             : }
    2726             : 
    2727             : void
    2728             : timer_start(pari_timer *T)
    2729             : {
    2730             :   struct timespec t;
    2731             :   clock_gettime(time_type,&t);
    2732             :   T->us = t.tv_nsec / 1000;
    2733             :   T->s  = t.tv_sec;
    2734             : }
    2735             : #elif defined(USE_GETRUSAGE)
    2736             : #ifdef RUSAGE_THREAD
    2737             : static THREAD int rusage_type = RUSAGE_THREAD;
    2738             : #else
    2739             : static const THREAD int rusage_type = RUSAGE_SELF;
    2740             : #endif /*RUSAGE_THREAD*/
    2741             : static void
    2742        1858 : pari_init_timer(void)
    2743             : {
    2744             : #ifdef RUSAGE_THREAD
    2745        1858 :   rusage_type = RUSAGE_SELF;
    2746             : #endif
    2747        1858 : }
    2748             : 
    2749             : void
    2750      348116 : timer_start(pari_timer *T)
    2751             : {
    2752             :   struct rusage r;
    2753      348116 :   getrusage(rusage_type,&r);
    2754      348120 :   T->us = r.ru_utime.tv_usec;
    2755      348120 :   T->s  = r.ru_utime.tv_sec;
    2756      348120 : }
    2757             : #elif defined(USE_FTIME)
    2758             : 
    2759             : static void
    2760             : pari_init_timer(void) { }
    2761             : 
    2762             : void
    2763             : timer_start(pari_timer *T)
    2764             : {
    2765             :   struct timeb t;
    2766             :   ftime(&t);
    2767             :   T->us = ((long)t.millitm) * 1000;
    2768             :   T->s  = t.time;
    2769             : }
    2770             : 
    2771             : #else
    2772             : 
    2773             : static void
    2774             : _get_time(pari_timer *T, long Ticks, long TickPerSecond)
    2775             : {
    2776             :   T->us = (long) ((Ticks % TickPerSecond) * (1000000. / TickPerSecond));
    2777             :   T->s  = Ticks / TickPerSecond;
    2778             : }
    2779             : 
    2780             : # ifdef USE_TIMES
    2781             : static void
    2782             : pari_init_timer(void) { }
    2783             : 
    2784             : void
    2785             : timer_start(pari_timer *T)
    2786             : {
    2787             : # ifdef _SC_CLK_TCK
    2788             :   long tck = sysconf(_SC_CLK_TCK);
    2789             : # else
    2790             :   long tck = CLK_TCK;
    2791             : # endif
    2792             :   struct tms t; times(&t);
    2793             :   _get_time(T, t.tms_utime, tck);
    2794             : }
    2795             : # elif defined(_WIN32)
    2796             : static void
    2797             : pari_init_timer(void) { }
    2798             : 
    2799             : void
    2800             : timer_start(pari_timer *T)
    2801             : { _get_time(T, win32_timer(), 1000); }
    2802             : # else
    2803             : #  include <time.h>
    2804             : #  ifndef CLOCKS_PER_SEC
    2805             : #   define CLOCKS_PER_SEC 1000000 /* may be false on YOUR system */
    2806             : #  endif
    2807             : static void
    2808             : pari_init_timer(void) { }
    2809             : 
    2810             : void
    2811             : timer_start(pari_timer *T)
    2812             : { _get_time(T, clock(), CLOCKS_PER_SEC); }
    2813             : # endif
    2814             : #endif
    2815             : 
    2816             : /* round microseconds to milliseconds */
    2817             : static long
    2818      212886 : rndus(long x) { return (x + 500) / 1000; }
    2819             : static long
    2820      212878 : timer_aux(pari_timer *T, pari_timer *U, void (*settime)(pari_timer *))
    2821             : {
    2822      212878 :   long s = T->s, us = T->us;
    2823      212878 :   settime(U); return 1000 * (U->s - s) + rndus(U->us - us);
    2824             : }
    2825             : 
    2826             : /* return delay, set timer checkpoint */
    2827             : long
    2828      107318 : timer_delay(pari_timer *T) { return timer_aux(T, T, &timer_start); }
    2829             : /* return delay, don't set checkpoint */
    2830             : long
    2831        1876 : timer_get(pari_timer *T) {pari_timer t; return timer_aux(T, &t, &timer_start);}
    2832             : 
    2833             : static void
    2834           0 : timer_vprintf(pari_timer *T, const char *format, va_list args)
    2835             : {
    2836           0 :   out_puts(pariErr, "Time ");
    2837           0 :   out_vprintf(pariErr, format,args);
    2838           0 :   out_printf(pariErr, ": %ld\n", timer_delay(T));
    2839           0 :   pariErr->flush();
    2840           0 : }
    2841             : void
    2842           0 : timer_printf(pari_timer *T, const char *format, ...)
    2843             : {
    2844           0 :   va_list args; va_start(args, format);
    2845           0 :   timer_vprintf(T, format, args);
    2846           0 :   va_end(args);
    2847           0 : }
    2848             : 
    2849             : long
    2850           0 : timer(void)  { static THREAD pari_timer T; return timer_delay(&T);}
    2851             : long
    2852        3634 : gettime(void)  { static THREAD pari_timer T; return timer_delay(&T);}
    2853             : 
    2854             : static THREAD pari_timer timer2_T, abstimer_T;
    2855             : long
    2856           0 : timer2(void) {  return timer_delay(&timer2_T);}
    2857             : void
    2858           0 : msgtimer(const char *format, ...)
    2859             : {
    2860           0 :   va_list args; va_start(args, format);
    2861           0 :   timer_vprintf(&timer2_T, format, args);
    2862           0 :   va_end(args);
    2863           0 : }
    2864             : long
    2865        1870 : getabstime(void)  { return timer_get(&abstimer_T);}
    2866             : 
    2867             : void
    2868      243493 : walltimer_start(pari_timer *ti)
    2869             : {
    2870             : #if defined(USE_CLOCK_GETTIME)
    2871             :   struct timespec t;
    2872             :   if (!clock_gettime(CLOCK_REALTIME,&t))
    2873             :   { ti->s = t.tv_sec; ti->us = rndus(t.tv_nsec); return; }
    2874             : #elif defined(USE_GETTIMEOFDAY)
    2875             :   struct timeval tv;
    2876      243493 :   if (!gettimeofday(&tv, NULL))
    2877      243493 :   {  ti->s = tv.tv_sec; ti->us = tv.tv_usec; return; }
    2878             : #elif defined(USE_FTIMEFORWALLTIME)
    2879             :   struct timeb tp;
    2880             :   if (!ftime(&tp))
    2881             :   { ti->s = tp.time; ti->us = tp.millitm*1000; return; }
    2882             : #endif
    2883           0 :   timer_start(ti);
    2884             : }
    2885             : /* return delay, set timer checkpoint */
    2886             : long
    2887      103684 : walltimer_delay(pari_timer *T) { return timer_aux(T, T, &walltimer_start); }
    2888             : /* return delay, don't set checkpoint */
    2889             : long
    2890           0 : walltimer_get(pari_timer *T)
    2891             : {
    2892             :   pari_timer t;
    2893           0 :   return timer_aux(T, &t, &walltimer_start);
    2894             : }
    2895             : 
    2896             : static GEN
    2897           8 : timetoi(ulong s, ulong m)
    2898             : {
    2899           8 :   pari_sp av = avma;
    2900           8 :   return gerepileuptoint(av, addiu(muluu(s, 1000), m));
    2901             : }
    2902             : GEN
    2903           8 : getwalltime(void)
    2904             : {
    2905             :   pari_timer ti;
    2906           8 :   walltimer_start(&ti);
    2907           8 :   return timetoi(ti.s, rndus(ti.us));
    2908             : }
    2909             : 
    2910             : /*******************************************************************/
    2911             : /*                                                                 */
    2912             : /*                   FUNCTIONS KNOWN TO THE ANALYZER               */
    2913             : /*                                                                 */
    2914             : /*******************************************************************/
    2915             : 
    2916             : GEN
    2917         127 : setdebug(const char *s, long n)
    2918             : {
    2919         127 :   long i, l = numberof(pari_DEBUGLEVEL_str);
    2920             :   GEN V, V1, V2;
    2921         127 :   if (s)
    2922             :   {
    2923         120 :     if (n > 20)
    2924           0 :       pari_err_DOMAIN("setdebug", "n", ">", utoipos(20), stoi(n));
    2925        2276 :     for (i = 0; i < l; i++)
    2926        2248 :       if (!strcmp(s, pari_DEBUGLEVEL_str[i])) break;
    2927         120 :     if (i == l)
    2928          28 :       pari_err_DOMAIN("setdebug", s, "not a valid",
    2929             :                       strtoGENstr("debug domain"), strtoGENstr(s));
    2930          92 :     if (n >= 0) { *pari_DEBUGLEVEL_ptr[i] = n; return gnil; }
    2931          42 :     return stoi(*pari_DEBUGLEVEL_ptr[i]);
    2932             :   }
    2933           7 :   V = cgetg(3,t_MAT);
    2934           7 :   V1 = gel(V,1) = cgetg(l+1, t_COL);
    2935           7 :   V2 = gel(V,2) = cgetg(l+1, t_COL);
    2936         427 :   for (i = 0; i < l; i++)
    2937             :   {
    2938         420 :     gel(V1, i+1) = strtoGENstr(pari_DEBUGLEVEL_str[i]);
    2939         420 :     gel(V2, i+1) = stoi(*pari_DEBUGLEVEL_ptr[i]);
    2940             :   }
    2941           7 :   return V;
    2942             : }
    2943             : 
    2944             : GEN
    2945           7 : pari_version(void)
    2946             : {
    2947           7 :   const ulong mask = (1UL<<PARI_VERSION_SHIFT) - 1;
    2948           7 :   ulong major, minor, patch, n = paricfg_version_code;
    2949           7 :   patch = n & mask; n >>= PARI_VERSION_SHIFT;
    2950           7 :   minor = n & mask; n >>= PARI_VERSION_SHIFT;
    2951           7 :   major = n;
    2952           7 :   if (*paricfg_vcsversion) {
    2953           7 :     const char *ver = paricfg_vcsversion;
    2954           7 :     const char *s = strchr(ver, '-');
    2955             :     char t[8];
    2956           7 :     const long len = s-ver;
    2957             :     GEN v;
    2958           7 :     if (!s || len > 6) pari_err_BUG("pari_version()"); /* paranoia */
    2959           7 :     memcpy(t, ver, len); t[len] = 0;
    2960           7 :     v = cgetg(6, t_VEC);
    2961           7 :     gel(v,1) = utoi(major);
    2962           7 :     gel(v,2) = utoi(minor);
    2963           7 :     gel(v,3) = utoi(patch);
    2964           7 :     gel(v,4) = stoi( atoi(t) );
    2965           7 :     gel(v,5) = strtoGENstr(s+1);
    2966           7 :     return v;
    2967             :   } else {
    2968           0 :     GEN v = cgetg(4, t_VEC);
    2969           0 :     gel(v,1) = utoi(major);
    2970           0 :     gel(v,2) = utoi(minor);
    2971           0 :     gel(v,3) = utoi(patch);
    2972           0 :     return v;
    2973             :   }
    2974             : }
    2975             : 
    2976             : /* List of GP functions: generated from the description system.
    2977             :  * Format (struct entree) :
    2978             :  *   char *name   : name (under GP).
    2979             :  *   ulong valence: (EpNEW, EpALIAS,EpVAR, EpINSTALL)|EpSTATIC
    2980             :  *   void *value  : For PREDEFINED FUNCTIONS: C function to call.
    2981             :  *                  For USER FUNCTIONS: pointer to defining data (block) =
    2982             :  *                   entree*: NULL, list of entree (arguments), NULL
    2983             :  *                   char*  : function text
    2984             :  *   long menu    : which help section do we belong to
    2985             :  *                   1: Standard monadic or dyadic OPERATORS
    2986             :  *                   2: CONVERSIONS and similar elementary functions
    2987             :  *                   3: functions related to COMBINATORICS
    2988             :  *                   4: TRANSCENDENTAL functions, etc.
    2989             :  *   char *code   : GP prototype, aka Parser Code (see libpari's manual)
    2990             :  *                  if NULL, use valence instead.
    2991             :  *   char *help   : short help text (init to NULL).
    2992             :  *   void *pvalue : push_val history.
    2993             :  *   long arity   : maximum number of arguments.
    2994             :  *   entree *next : next entree (init to NULL, used in hashing code). */
    2995             : #include "init.h"
    2996             : #include "default.h"

Generated by: LCOV version 1.16