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

Generated by: LCOV version 1.14