Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - language - init.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 20079-8a65571) Lines: 832 1168 71.2 %
Date: 2017-01-18 05:50:33 Functions: 100 132 75.8 %
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. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : 
      14             : /*******************************************************************/
      15             : /*                                                                 */
      16             : /*        INITIALIZING THE SYSTEM, ERRORS, STACK MANAGEMENT        */
      17             : /*                                                                 */
      18             : /*******************************************************************/
      19             : /* _GNU_SOURCE is needed before first include to get RUSAGE_THREAD */
      20             : #undef _GNU_SOURCE /* avoid warning */
      21             : #define _GNU_SOURCE
      22             : #include <string.h>
      23             : #if defined(_WIN32) || defined(__CYGWIN32__)
      24             : #  include "../systems/mingw/mingw.h"
      25             : #  include <process.h>
      26             : #endif
      27             : #include "paricfg.h"
      28             : #if defined(STACK_CHECK) && !defined(__EMX__)
      29             : #  include <sys/types.h>
      30             : #  include <sys/time.h>
      31             : #  include <sys/resource.h>
      32             : #endif
      33             : #if defined(HAS_WAITPID) && defined(HAS_SETSID)
      34             : #  include <sys/wait.h>
      35             : #endif
      36             : #ifdef HAS_MMAP
      37             : #  include <sys/mman.h>
      38             : #endif
      39             : #if defined(USE_GETTIMEOFDAY) || defined(USE_GETRUSAGE) || defined(USE_TIMES)
      40             : #  include <sys/time.h>
      41             : #endif
      42             : #if defined(USE_GETRUSAGE)
      43             : #  include <sys/resource.h>
      44             : #endif
      45             : #if defined(USE_FTIME) || defined(USE_FTIMEFORWALLTIME)
      46             : #  include <sys/timeb.h>
      47             : #endif
      48             : #if defined(USE_CLOCK_GETTIME) || defined(USE_TIMES)
      49             : #  include <time.h>
      50             : #endif
      51             : #if defined(USE_TIMES)
      52             : #  include <sys/times.h>
      53             : #endif
      54             : #include "pari.h"
      55             : #include "paripriv.h"
      56             : #include "anal.h"
      57             : 
      58             : const double LOG2    = 0.6931471805599453; /* log(2) */
      59             : const double LOG10_2 = 0.3010299956639812; /* log_10(2) */
      60             : const double LOG2_10 = 3.321928094887362;  /* log_2(10) */
      61             : 
      62             : GEN gnil, gen_0, gen_1, gen_m1, gen_2, gen_m2, ghalf, err_e_STACK;
      63             : 
      64             : static const ulong readonly_constants[] = {
      65             :   evaltyp(t_INT) | _evallg(2),  /* gen_0 */
      66             :   evallgefint(2),
      67             :   evaltyp(t_INT) | _evallg(2),  /* gnil */
      68             :   evallgefint(2),
      69             :   evaltyp(t_INT) | _evallg(3),  /* gen_1 */
      70             :   evalsigne(1) | evallgefint(3),
      71             :   1,
      72             :   evaltyp(t_INT) | _evallg(3),  /* gen_2 */
      73             :   evalsigne(1) | evallgefint(3),
      74             :   2,
      75             :   evaltyp(t_INT) | _evallg(3),  /* gen_m1 */
      76             :   evalsigne(-1) | evallgefint(3),
      77             :   1,
      78             :   evaltyp(t_INT) | _evallg(3),  /* gen_m2 */
      79             :   evalsigne(-1) | evallgefint(3),
      80             :   2,
      81             :   evaltyp(t_FRAC) | _evallg(3), /* ghalf */
      82             :   (ulong)(readonly_constants+4),
      83             :   (ulong)(readonly_constants+7)
      84             : };
      85             : static const ulong readonly_err_STACK[] = {
      86             :   evaltyp(t_ERROR) | _evallg(2),
      87             :   e_STACK
      88             : };
      89             : THREAD GEN    bernzone;
      90             : GEN     primetab; /* private primetable */
      91             : byteptr diffptr;
      92             : FILE    *pari_outfile, *pari_errfile, *pari_logfile, *pari_infile;
      93             : char    *current_logfile, *current_psfile, *pari_datadir;
      94             : long    gp_colors[c_LAST];
      95             : int     disable_color;
      96             : ulong   DEBUGFILES, DEBUGLEVEL, DEBUGMEM;
      97             : long    DEBUGVAR;
      98             : ulong   pari_mt_nbthreads;
      99             : long    precreal;
     100             : ulong   precdl, logstyle;
     101             : gp_data *GP_DATA;
     102             : 
     103             : GEN colormap, pari_graphcolors;
     104             : 
     105             : entree  **varentries;
     106             : THREAD long *varpriority;
     107             : 
     108             : THREAD pari_sp avma;
     109             : THREAD struct pari_mainstack *pari_mainstack;
     110             : 
     111             : static void ** MODULES;
     112             : static pari_stack s_MODULES;
     113             : const long functions_tblsz = 135; /* size of functions_hash */
     114             : entree **functions_hash, **defaults_hash;
     115             : 
     116             : char *(*cb_pari_fgets_interactive)(char *s, int n, FILE *f);
     117             : int (*cb_pari_get_line_interactive)(const char*, const char*, filtre_t *F);
     118             : void (*cb_pari_quit)(long);
     119             : void (*cb_pari_init_histfile)(void);
     120             : void (*cb_pari_ask_confirm)(const char *);
     121             : int  (*cb_pari_handle_exception)(long);
     122             : int  (*cb_pari_err_handle)(GEN);
     123             : int  (*cb_pari_whatnow)(PariOUT *out, const char *, int);
     124             : void (*cb_pari_sigint)(void);
     125             : void (*cb_pari_pre_recover)(long);
     126             : void (*cb_pari_err_recover)(long);
     127             : int (*cb_pari_break_loop)(int);
     128             : int (*cb_pari_is_interactive)(void);
     129             : void (*cb_pari_start_output)();
     130             : 
     131             : const char * pari_library_path = NULL;
     132             : 
     133             : static THREAD GEN global_err_data;
     134             : THREAD jmp_buf *iferr_env;
     135             : const long CATCH_ALL = -1;
     136             : 
     137             : static void pari_init_timer(void);
     138             : 
     139             : /*********************************************************************/
     140             : /*                                                                   */
     141             : /*                       BLOCKS & CLONES                             */
     142             : /*                                                                   */
     143             : /*********************************************************************/
     144             : /*#define DEBUG*/
     145             : static THREAD long next_block;
     146             : static THREAD GEN cur_block; /* current block in block list */
     147             : #ifdef DEBUG
     148             : static THREAD long NUM;
     149             : #endif
     150             : 
     151             : static void
     152       63050 : pari_init_blocks(void)
     153             : {
     154       63050 :   next_block = 0; cur_block = NULL;
     155             : #ifdef DEBUG
     156             :   NUM = 0;
     157             : #endif
     158       63050 : }
     159             : 
     160             : static void
     161       64345 : pari_close_blocks(void)
     162             : {
     163       64345 :   while (cur_block) killblock(cur_block);
     164       64237 : }
     165             : 
     166             : /* Return x, where:
     167             :  * x[-4]: reference count
     168             :  * x[-3]: adress of next block
     169             :  * x[-2]: adress of preceding block.
     170             :  * x[-1]: number of allocated blocs.
     171             :  * x[0..n-1]: malloc-ed memory. */
     172             : GEN
     173   137670266 : newblock(size_t n)
     174             : {
     175   137670266 :   long *x = (long *) pari_malloc((n + BL_HEAD)*sizeof(long)) + BL_HEAD;
     176             : 
     177   137671304 :   bl_refc(x) = 1;
     178   137671304 :   bl_next(x) = NULL;
     179   137671304 :   bl_prev(x) = cur_block;
     180   137671304 :   bl_num(x)  = next_block++;
     181   137671304 :   if (cur_block) bl_next(cur_block) = x;
     182             : #ifdef DEBUG
     183             :   err_printf("+ %ld\n", ++NUM);
     184             : #endif
     185   137671304 :   if (DEBUGMEM)
     186             :   {
     187           0 :     if (!n) pari_warn(warner,"mallocing NULL object in newblock");
     188           0 :     if (DEBUGMEM > 2)
     189           0 :       err_printf("new block, size %6lu (no %ld): %08lx\n", n, next_block-1, x);
     190             :   }
     191   137671178 :   return cur_block = x;
     192             : }
     193             : 
     194             : GEN
     195         303 : gcloneref(GEN x)
     196             : {
     197         303 :   if (isclone(x)) { ++bl_refc(x); return x; }
     198         303 :   else return gclone(x);
     199             : }
     200             : 
     201             : void
     202           0 : gclone_refc(GEN x) { ++bl_refc(x); }
     203             : 
     204             : void
     205   179967366 : gunclone(GEN x)
     206             : {
     207   359934723 :   if (--bl_refc(x) > 0) return;
     208   137678638 :   BLOCK_SIGINT_START;
     209   137678870 :   if (bl_next(x)) bl_prev(bl_next(x)) = bl_prev(x);
     210             :   else
     211             :   {
     212    24593385 :     cur_block = bl_prev(x);
     213    24593385 :     next_block = bl_num(x);
     214             :   }
     215   137678870 :   if (bl_prev(x)) bl_next(bl_prev(x)) = bl_next(x);
     216   137678870 :   if (DEBUGMEM > 2)
     217           0 :     err_printf("killing block (no %ld): %08lx\n", bl_num(x), x);
     218   137678747 :   free((void*)bl_base(x)); /* pari_free not needed: we already block */
     219   137678747 :   BLOCK_SIGINT_END;
     220             : #ifdef DEBUG
     221             :   err_printf("- %ld\n", NUM--);
     222             : #endif
     223             : }
     224             : 
     225             : /* Recursively look for clones in the container and kill them. Then kill
     226             :  * container if clone. SIGINT could be blocked until it returns */
     227             : void
     228  2540429907 : gunclone_deep(GEN x)
     229             : {
     230             :   long i, lx;
     231             :   GEN v;
     232  5080859814 :   if (isclone(x) && bl_refc(x) > 1) { --bl_refc(x); return; }
     233  2518851476 :   BLOCK_SIGINT_START;
     234  2518851476 :   switch(typ(x))
     235             :   {
     236             :     case t_VEC: case t_COL: case t_MAT:
     237    33385214 :       lx = lg(x);
     238    33385214 :       for (i=1;i<lx;i++) gunclone_deep(gel(x,i));
     239    33385214 :       break;
     240             :     case t_LIST:
     241         270 :       v = list_data(x); lx = v? lg(v): 1;
     242         270 :       for (i=1;i<lx;i++) gunclone_deep(gel(v,i));
     243         270 :       if (v) killblock(v);
     244         270 :       break;
     245             :   }
     246  2518851476 :   if (isclone(x)) gunclone(x);
     247  2518851476 :   BLOCK_SIGINT_END;
     248             : }
     249             : 
     250             : int
     251      138572 : pop_entree_block(entree *ep, long loc)
     252             : {
     253      138572 :   GEN x = (GEN)ep->value;
     254      138572 :   if (bl_num(x) < loc) return 0; /* older */
     255         112 :   if (DEBUGMEM>2)
     256           0 :     err_printf("popping %s (block no %ld)\n", ep->name, bl_num(x));
     257         112 :   gunclone_deep(x); return 1;
     258             : }
     259             : 
     260             : /*********************************************************************/
     261             : /*                                                                   */
     262             : /*                       C STACK SIZE CONTROL                        */
     263             : /*                                                                   */
     264             : /*********************************************************************/
     265             : /* Avoid core dump on deep recursion. Adapted Perl code by Dominic Dunlop */
     266             : THREAD void *PARI_stack_limit = NULL;
     267             : 
     268             : #ifdef STACK_CHECK
     269             : 
     270             : #  ifdef __EMX__                                /* Emulate */
     271             : void
     272             : pari_stackcheck_init(void *pari_stack_base)
     273             : {
     274             :   (void) pari_stack_base;
     275             :   if (!pari_stack_base) { PARI_stack_limit = NULL; return; }
     276             :   PARI_stack_limit = get_stack(1./16, 32*1024);
     277             : }
     278             : #  else /* !__EMX__ */
     279             : /* Set PARI_stack_limit to (a little above) the lowest safe address that can be
     280             :  * used on the stack. Leave PARI_stack_limit at its initial value (NULL) to
     281             :  * show no check should be made [init failed]. Assume stack grows downward. */
     282             : void
     283        1399 : pari_stackcheck_init(void *pari_stack_base)
     284             : {
     285             :   struct rlimit rip;
     286             :   ulong size;
     287        1399 :   if (!pari_stack_base) { PARI_stack_limit = NULL; return; }
     288        1399 :   if (getrlimit(RLIMIT_STACK, &rip)) return;
     289        1399 :   size = rip.rlim_cur;
     290        1399 :   if (size == (ulong)RLIM_INFINITY || size > (ulong)pari_stack_base)
     291           0 :     PARI_stack_limit = (void*)(((ulong)pari_stack_base) / 16);
     292             :   else
     293        1399 :     PARI_stack_limit = (void*)((ulong)pari_stack_base - (size/16)*15);
     294             : }
     295             : #  endif /* !__EMX__ */
     296             : 
     297             : #else
     298             : void
     299             : pari_stackcheck_init(void *pari_stack_base)
     300             : {
     301             :   (void) pari_stack_base; PARI_stack_limit = NULL;
     302             : }
     303             : #endif /* STACK_CHECK */
     304             : 
     305             : /*******************************************************************/
     306             : /*                         HEAP TRAVERSAL                          */
     307             : /*******************************************************************/
     308             : struct getheap_t { long n, l; };
     309             : static void
     310        6615 : f_getheap(GEN x, void *D)
     311             : {
     312        6615 :   struct getheap_t *T = (struct getheap_t*)D;
     313        6615 :   T->n++;
     314        6615 :   T->l += gsizeword(x);
     315        6615 : }
     316             : GEN
     317          84 : getheap(void)
     318             : {
     319          84 :   struct getheap_t T = { 0, 0 };
     320          84 :   traverseheap(&f_getheap, &T);
     321          84 :   return mkvec2s(T.n, T.l + BL_HEAD * T.n);
     322             : }
     323             : 
     324             : void
     325          84 : traverseheap( void(*f)(GEN, void *), void *data )
     326             : {
     327             :   GEN x;
     328          84 :   for (x = cur_block; x; x = bl_prev(x)) f(x, data);
     329          84 : }
     330             : 
     331             : /*********************************************************************/
     332             : /*                          DAEMON / FORK                            */
     333             : /*********************************************************************/
     334             : #if defined(HAS_WAITPID) && defined(HAS_SETSID)
     335             : /* Properly fork a process, detaching from main process group without creating
     336             :  * zombies on exit. Parent returns 1, son returns 0 */
     337             : int
     338         508 : pari_daemon(void)
     339             : {
     340         508 :   pid_t pid = fork();
     341        1016 :   switch(pid) {
     342           0 :       case -1: return 1; /* father, fork failed */
     343             :       case 0:
     344         508 :         (void)setsid(); /* son becomes process group leader */
     345         508 :         if (fork()) exit(0); /* now son exits, also when fork fails */
     346         508 :         break; /* grandson: its father is the son, which exited,
     347             :                 * hence father becomes 'init', that'll take care of it */
     348             :       default: /* father, fork succeeded */
     349         508 :         (void)waitpid(pid,NULL,0); /* wait for son to exit, immediate */
     350         508 :         return 1;
     351             :   }
     352             :   /* grandson */
     353         508 :   return 0;
     354             : }
     355             : #else
     356             : int
     357             : pari_daemon(void)
     358             : {
     359             :   pari_err_IMPL("pari_daemon without waitpid & setsid");
     360             :   return 0;
     361             : }
     362             : #endif
     363             : 
     364             : /*********************************************************************/
     365             : /*                                                                   */
     366             : /*                       SYSTEM INITIALIZATION                       */
     367             : /*                                                                   */
     368             : /*********************************************************************/
     369             : static int try_to_recover = 0;
     370             : THREAD VOLATILE int PARI_SIGINT_block = 0, PARI_SIGINT_pending = 0;
     371             : 
     372             : /*********************************************************************/
     373             : /*                         SIGNAL HANDLERS                           */
     374             : /*********************************************************************/
     375             : static void
     376           0 : dflt_sigint_fun(void) { pari_err(e_MISC, "user interrupt"); }
     377             : 
     378             : #if defined(_WIN32) || defined(__CYGWIN32__)
     379             : int win32ctrlc = 0, win32alrm = 0;
     380             : void
     381             : dowin32ctrlc(void)
     382             : {
     383             :   win32ctrlc = 0;
     384             :   cb_pari_sigint();
     385             : }
     386             : #endif
     387             : 
     388             : static void
     389           0 : pari_handle_SIGINT(void)
     390             : {
     391             : #ifdef _WIN32
     392             :   if (++win32ctrlc >= 5) _exit(3);
     393             : #else
     394           0 :   cb_pari_sigint();
     395             : #endif
     396           0 : }
     397             : 
     398             : void
     399           0 : pari_sighandler(int sig)
     400             : {
     401             :   const char *msg;
     402             : #ifndef HAS_SIGACTION
     403             :   /*SYSV reset the signal handler in the handler*/
     404             :   (void)os_signal(sig,pari_sighandler);
     405             : #endif
     406           0 :   switch(sig)
     407             :   {
     408             : #ifdef SIGBREAK
     409             :     case SIGBREAK:
     410             :       if (PARI_SIGINT_block==1)
     411             :       {
     412             :         PARI_SIGINT_pending=SIGBREAK;
     413             :         mt_sigint();
     414             :       }
     415             :       else pari_handle_SIGINT();
     416             :       return;
     417             : #endif
     418             : 
     419             : #ifdef SIGINT
     420             :     case SIGINT:
     421           0 :       if (PARI_SIGINT_block==1)
     422             :       {
     423           0 :         PARI_SIGINT_pending=SIGINT;
     424           0 :         mt_sigint();
     425             :       }
     426           0 :       else pari_handle_SIGINT();
     427           0 :       return;
     428             : #endif
     429             : 
     430             : #ifdef SIGSEGV
     431             :     case SIGSEGV:
     432           0 :       msg="PARI/GP (Segmentation Fault)"; break;
     433             : #endif
     434             : #ifdef SIGBUS
     435             :     case SIGBUS:
     436           0 :       msg="PARI/GP (Bus Error)"; break;
     437             : #endif
     438             : #ifdef SIGFPE
     439             :     case SIGFPE:
     440           0 :       msg="PARI/GP (Floating Point Exception)"; break;
     441             : #endif
     442             : 
     443             : #ifdef SIGPIPE
     444             :     case SIGPIPE:
     445             :     {
     446           0 :       pariFILE *f = GP_DATA->pp->file;
     447           0 :       if (f && pari_outfile == f->file)
     448             :       {
     449           0 :         GP_DATA->pp->file = NULL; /* to avoid oo recursion on error */
     450           0 :         pari_outfile = stdout; pari_fclose(f);
     451           0 :         pari_err(e_MISC, "Broken Pipe, resetting file stack...");
     452             :       }
     453             :       return; /* LCOV_EXCL_LINE */
     454             :     }
     455             : #endif
     456             : 
     457           0 :     default: msg="signal handling"; break;
     458             :   }
     459           0 :   pari_err_BUG(msg);
     460             : }
     461             : 
     462             : void
     463        3305 : pari_sig_init(void (*f)(int))
     464             : {
     465             : #ifdef SIGBUS
     466        3305 :   (void)os_signal(SIGBUS,f);
     467             : #endif
     468             : #ifdef SIGFPE
     469        3305 :   (void)os_signal(SIGFPE,f);
     470             : #endif
     471             : #ifdef SIGINT
     472        3305 :   (void)os_signal(SIGINT,f);
     473             : #endif
     474             : #ifdef SIGBREAK
     475             :   (void)os_signal(SIGBREAK,f);
     476             : #endif
     477             : #ifdef SIGPIPE
     478        3305 :   (void)os_signal(SIGPIPE,f);
     479             : #endif
     480             : #ifdef SIGSEGV
     481        3305 :   (void)os_signal(SIGSEGV,f);
     482             : #endif
     483        3305 : }
     484             : 
     485             : /*********************************************************************/
     486             : /*                      STACK AND UNIVERSAL CONSTANTS                */
     487             : /*********************************************************************/
     488             : static void
     489        1399 : init_universal_constants(void)
     490             : {
     491        1399 :   gen_0  = (GEN)readonly_constants;
     492        1399 :   gnil   = (GEN)readonly_constants+2;
     493        1399 :   gen_1  = (GEN)readonly_constants+4;
     494        1399 :   gen_2  = (GEN)readonly_constants+7;
     495        1399 :   gen_m1 = (GEN)readonly_constants+10;
     496        1399 :   gen_m2 = (GEN)readonly_constants+13;
     497        1399 :   ghalf  = (GEN)readonly_constants+16;
     498        1399 :   err_e_STACK = (GEN)readonly_err_STACK;
     499        1399 : }
     500             : 
     501             : static void
     502       63219 : pari_init_errcatch(void)
     503             : {
     504       63219 :   iferr_env = NULL;
     505       63219 :   global_err_data = NULL;
     506       63219 : }
     507             : 
     508             : /*********************************************************************/
     509             : /*                           INIT DEFAULTS                           */
     510             : /*********************************************************************/
     511             : void
     512        1399 : pari_init_defaults(void)
     513             : {
     514             :   long i;
     515        1399 :   initout(1);
     516             : 
     517             : #ifdef LONG_IS_64BIT
     518        1201 :   precreal = 128;
     519             : #else
     520         198 :   precreal = 96;
     521             : #endif
     522             : 
     523        1399 :   precdl = 16;
     524        1399 :   DEBUGFILES = DEBUGLEVEL = DEBUGMEM = 0;
     525        1399 :   disable_color = 1;
     526        1399 :   logstyle = logstyle_none;
     527             : 
     528        1399 :   current_psfile = pari_strdup("pari.ps");
     529        1399 :   current_logfile= pari_strdup("pari.log");
     530        1399 :   pari_logfile = NULL;
     531             : 
     532        1399 :   pari_datadir = os_getenv("GP_DATA_DIR");
     533        1399 :   if (!pari_datadir)
     534             :   {
     535             : #if defined(_WIN32) || defined(__CYGWIN32__)
     536             :     if (paricfg_datadir[0]=='@' && paricfg_datadir[1]==0)
     537             :       pari_datadir = win32_datadir();
     538             :     else
     539             : #endif
     540        1399 :       pari_datadir = pari_strdup(paricfg_datadir);
     541             :   }
     542           0 :   else pari_datadir= pari_strdup(pari_datadir);
     543        1399 :   for (i=0; i<c_LAST; i++) gp_colors[i] = c_NONE;
     544        1399 :   colormap = NULL; pari_graphcolors = NULL;
     545        1399 : }
     546             : 
     547             : /*********************************************************************/
     548             : /*                   FUNCTION HASHTABLES, MODULES                    */
     549             : /*********************************************************************/
     550             : extern entree functions_basic[], functions_default[];
     551             : static void
     552        1399 : pari_init_functions(void)
     553             : {
     554        1399 :   pari_stack_init(&s_MODULES, sizeof(*MODULES),(void**)&MODULES);
     555        1399 :   pari_stack_pushp(&s_MODULES,functions_basic);
     556        1399 :   functions_hash = (entree**) pari_calloc(sizeof(entree*)*functions_tblsz);
     557        1399 :   pari_fill_hashtable(functions_hash, functions_basic);
     558        1399 :   defaults_hash = (entree**) pari_calloc(sizeof(entree*)*functions_tblsz);
     559        1399 :   pari_add_defaults_module(functions_default);
     560        1399 : }
     561             : 
     562             : void
     563        2796 : pari_add_module(entree *ep)
     564             : {
     565        2796 :   pari_fill_hashtable(functions_hash, ep);
     566        2796 :   pari_stack_pushp(&s_MODULES, ep);
     567        2796 : }
     568             : 
     569             : void
     570        1399 : pari_add_defaults_module(entree *ep)
     571        1399 : { pari_fill_hashtable(defaults_hash, ep); }
     572             : 
     573             : /*********************************************************************/
     574             : /*                       PARI MAIN STACK                             */
     575             : /*********************************************************************/
     576             : 
     577             : #ifdef HAS_MMAP
     578             : #define PARI_STACK_ALIGN (sysconf(_SC_PAGE_SIZE))
     579             : #ifndef MAP_ANONYMOUS
     580             : #define MAP_ANONYMOUS MAP_ANON
     581             : #endif
     582             : static void *
     583       64439 : pari_mainstack_malloc(size_t size)
     584             : {
     585       64439 :   void *b = mmap(NULL, size, PROT_READ|PROT_WRITE,
     586             :                              MAP_PRIVATE|MAP_ANONYMOUS, -1, 0);
     587       64439 :   return (b == MAP_FAILED) ? NULL: b;
     588             : }
     589             : 
     590             : static void
     591       64946 : pari_mainstack_mfree(void *s, size_t size)
     592             : {
     593       64946 :   munmap(s, size);
     594       64946 : }
     595             : 
     596             : /* Completely discard the memory mapped between the addresses "from"
     597             :  * and "to" (which must be page-aligned).
     598             :  *
     599             :  * We use mmap() with PROT_NONE, which means that the underlying memory
     600             :  * is freed and that the kernel should not commit memory for it. We
     601             :  * still keep the mapping such that we can change the flags to
     602             :  * PROT_READ|PROT_WRITE later.
     603             :  *
     604             :  * NOTE: remapping with MAP_FIXED and PROT_NONE is not the same as
     605             :  * calling mprotect(..., PROT_NONE) because the latter will keep the
     606             :  * memory committed (this is in particular relevant on Linux with
     607             :  * vm.overcommit = 2). This remains true even when calling
     608             :  * madvise(..., MADV_DONTNEED). */
     609             : static void
     610      134843 : pari_mainstack_mreset(pari_sp from, pari_sp to)
     611             : {
     612      134843 :   size_t s = to - from;
     613      134843 :   mmap((void*)from, s, PROT_NONE, MAP_FIXED|MAP_PRIVATE|MAP_ANONYMOUS, -1, 0);
     614      134843 : }
     615             : 
     616             : /* Commit (make available) the virtual memory mapped between the
     617             :  * addresses "from" and "to" (which must be page-aligned).
     618             :  * Return 0 if successful, -1 if failed. */
     619             : static int
     620      134843 : pari_mainstack_mextend(pari_sp from, pari_sp to)
     621             : {
     622      134843 :   size_t s = to - from;
     623      134843 :   return mprotect((void*)from, s, PROT_READ|PROT_WRITE);
     624             : }
     625             : 
     626             : /* Set actual stack size to the given size. This sets st->size and
     627             :  * st->bot. If not enough system memory is available, this can fail.
     628             :  * Return 1 if successful, 0 if failed (in that case, st->size is not
     629             :  * changed) */
     630             : static int
     631      134843 : pari_mainstack_setsize(struct pari_mainstack *st, size_t size)
     632             : {
     633      134843 :   pari_sp newbot = st->top - size;
     634             :   /* Align newbot to pagesize */
     635      134843 :   pari_sp alignbot = newbot & ~(pari_sp)(PARI_STACK_ALIGN - 1);
     636      134843 :   if (pari_mainstack_mextend(alignbot, st->top))
     637             :   {
     638             :     /* Making the memory available did not work: limit vsize to the
     639             :      * current actual stack size. */
     640           0 :     st->vsize = st->size;
     641           0 :     pari_warn(warnstack, st->vsize);
     642           0 :     return 0;
     643             :   }
     644      134843 :   pari_mainstack_mreset(st->vbot, alignbot);
     645      134843 :   st->bot = newbot;
     646      134843 :   st->size = size;
     647      134843 :   return 1;
     648             : }
     649             : 
     650             : #else
     651             : #define PARI_STACK_ALIGN (0x40UL)
     652             : static void *
     653             : pari_mainstack_malloc(size_t s)
     654             : {
     655             :   return malloc(s); /* NOT pari_malloc, e_MEM would be deadly */
     656             : }
     657             : 
     658             : static void
     659             : pari_mainstack_mfree(void *s, size_t size) { (void) size; free(s); }
     660             : 
     661             : static int
     662             : pari_mainstack_setsize(struct pari_mainstack *st, size_t size)
     663             : {
     664             :   st->bot = st->top - size;
     665             :   st->size = size;
     666             :   return 1;
     667             : }
     668             : 
     669             : #endif
     670             : 
     671             : static const size_t MIN_STACK = 500032UL;
     672             : static size_t
     673      129385 : fix_size(size_t a)
     674             : {
     675      129385 :   size_t ps = PARI_STACK_ALIGN;
     676      129385 :   size_t b = a & ~(ps - 1); /* Align */
     677      129385 :   if (b < a && b < ~(ps - 1)) b += ps;
     678      129385 :   if (b < MIN_STACK) b = MIN_STACK;
     679      129385 :   return b;
     680             : }
     681             : 
     682             : static void
     683       64439 : pari_mainstack_alloc(struct pari_mainstack *st, size_t rsize, size_t vsize)
     684             : {
     685       64439 :   size_t sizemax = vsize ? vsize: rsize, s = fix_size(sizemax);
     686             :   for (;;)
     687             :   {
     688       64439 :     st->vbot = (pari_sp)pari_mainstack_malloc(s);
     689       64439 :     if (st->vbot) break;
     690           0 :     if (s == MIN_STACK) pari_err(e_MEM); /* no way out. Die */
     691           0 :     s = fix_size(s >> 1);
     692           0 :     pari_warn(warnstack, s);
     693           0 :   }
     694       64439 :   st->vsize = vsize ? s: 0;
     695       64439 :   st->rsize = minuu(rsize, s);
     696       64439 :   st->top = st->vbot+s;
     697       64439 :   if (!pari_mainstack_setsize(st, st->rsize))
     698             :   {
     699             :     /* This should never happen since we only decrease the allocated space */
     700           0 :     pari_err(e_MEM);
     701             :   }
     702       64439 :   st->memused = 0;
     703       64439 : }
     704             : 
     705             : static void
     706       64946 : pari_mainstack_free(struct pari_mainstack *st)
     707             : {
     708       64946 :   pari_mainstack_mfree((void*)st->vbot, st->vsize ? st->vsize : fix_size(st->rsize));
     709       64946 :   st->top = st->bot = st->vbot = 0;
     710       64946 :   st->size = st->vsize = 0;
     711       64946 : }
     712             : 
     713             : static void
     714         176 : pari_mainstack_resize(struct pari_mainstack *st, size_t rsize, size_t vsize)
     715             : {
     716         176 :   BLOCK_SIGINT_START;
     717         176 :   pari_mainstack_free(st);
     718         176 :   pari_mainstack_alloc(st, rsize, vsize);
     719         176 :   BLOCK_SIGINT_END;
     720         176 : }
     721             : 
     722             : static void
     723       63087 : pari_mainstack_use(struct pari_mainstack *st)
     724             : {
     725       63087 :   pari_mainstack = st;
     726       63087 :   avma = st->top;
     727       63087 : }
     728             : 
     729             : static void
     730        1399 : paristack_alloc(size_t rsize, size_t vsize)
     731             : {
     732        1399 :   pari_mainstack_alloc(pari_mainstack, rsize, vsize);
     733        1399 :   pari_mainstack_use(pari_mainstack);
     734        1399 : }
     735             : 
     736             : void
     737           0 : paristack_setsize(size_t rsize, size_t vsize)
     738             : {
     739           0 :   pari_mainstack_resize(pari_mainstack, rsize, vsize);
     740           0 :   pari_mainstack_use(pari_mainstack);
     741           0 : }
     742             : 
     743             : void
     744           0 : parivstack_resize(ulong newsize)
     745             : {
     746             :   size_t s;
     747           0 :   if (newsize && newsize < pari_mainstack->rsize)
     748           0 :     pari_err_DIM("stack sizes [parisizemax < parisize]");
     749           0 :   if (newsize == pari_mainstack->vsize) return;
     750           0 :   evalstate_reset();
     751           0 :   paristack_setsize(pari_mainstack->rsize, newsize);
     752           0 :   s = pari_mainstack->vsize ? pari_mainstack->vsize : pari_mainstack->rsize;
     753           0 :   pari_warn(warner,"new maximum stack size = %lu (%.3f Mbytes)", s, s/1048576.);
     754           0 :   pari_init_errcatch();
     755           0 :   cb_pari_err_recover(-1);
     756             : }
     757             : 
     758             : void
     759         182 : paristack_newrsize(ulong newsize)
     760             : {
     761         182 :   size_t s, vsize = pari_mainstack->vsize;
     762         182 :   if (!newsize) newsize = pari_mainstack->rsize << 1;
     763         182 :   if (newsize != pari_mainstack->rsize)
     764         176 :     pari_mainstack_resize(pari_mainstack, newsize, vsize);
     765         182 :   evalstate_reset();
     766         182 :   s = pari_mainstack->rsize;
     767         182 :   pari_warn(warner,"new stack size = %lu (%.3f Mbytes)", s, s/1048576.);
     768         182 :   pari_init_errcatch();
     769         182 :   cb_pari_err_recover(-1);
     770           0 : }
     771             : 
     772             : void
     773           7 : paristack_resize(ulong newsize)
     774             : {
     775           7 :   if (!newsize)
     776           0 :     newsize = 2 * pari_mainstack->size;
     777           7 :   newsize = minuu(newsize, pari_mainstack->vsize);
     778          14 :   if (newsize <= pari_mainstack->size) return;
     779           0 :   if (pari_mainstack_setsize(pari_mainstack, newsize))
     780             :   {
     781           0 :     pari_warn(warner, "increasing stack size to %lu", pari_mainstack->size);
     782             :   }
     783             : }
     784             : 
     785             : void
     786       70397 : parivstack_reset(void)
     787             : {
     788       70397 :   pari_mainstack_setsize(pari_mainstack, pari_mainstack->rsize);
     789       70397 : }
     790             : 
     791             : /* Enlarge the stack if needed such that the unused portion of the stack
     792             :  * (between bot and avma) is large enough to contain x longs. */
     793             : void
     794           7 : new_chunk_resize(size_t x)
     795             : {
     796             :   ulong size, newsize, avail;
     797           7 :   avail = (avma - pari_mainstack->bot) / sizeof(long);
     798           7 :   if (avail >= x) return;
     799             : 
     800             :   /* We need to enlarge the stack. We try to at least double the
     801             :    * stack, to avoid increasing the stack a lot of times by a small
     802             :    * amount. */
     803           7 :   size = pari_mainstack->size;
     804           7 :   newsize = size + maxuu((x - avail) * sizeof(long), size);
     805           7 :   paristack_resize(newsize);
     806             : 
     807             :   /* Verify that we have enough space. Using a division here instead
     808             :    * of a multiplication is safe against overflow. */
     809           7 :   avail = (avma - pari_mainstack->bot) / sizeof(long);
     810           7 :   if (avail < x)
     811             :   {
     812             :     /* Restore old size and error out */
     813           7 :     pari_mainstack_setsize(pari_mainstack, size);
     814           7 :     pari_err(e_STACK);
     815             :   }
     816             : }
     817             : 
     818             : /*********************************************************************/
     819             : /*                       PARI THREAD                                 */
     820             : /*********************************************************************/
     821             : 
     822             : /* Initial PARI thread structure t with a stack of size s and virtual size v
     823             :  * and argument arg */
     824             : 
     825             : void
     826           0 : pari_thread_valloc(struct pari_thread *t, size_t s, size_t v, GEN arg)
     827             : {
     828           0 :   pari_mainstack_alloc(&t->st,s,v);
     829           0 :   t->data = arg;
     830           0 : }
     831             : 
     832             : /* Initial PARI thread structure t with a stack of size s and
     833             :  * argument arg */
     834             : 
     835             : void
     836       62864 : pari_thread_alloc(struct pari_thread *t, size_t s, GEN arg)
     837             : {
     838       62864 :   pari_mainstack_alloc(&t->st,s,0);
     839       62864 :   t->data = arg;
     840       62864 : }
     841             : 
     842             : void
     843       62864 : pari_thread_free(struct pari_thread *t)
     844             : {
     845       62864 :   pari_mainstack_free(&t->st);
     846       62864 : }
     847             : 
     848             : void
     849       62997 : pari_thread_init(void)
     850             : {
     851       62997 :   pari_init_blocks();
     852       63191 :   pari_init_errcatch();
     853       63046 :   pari_init_rand();
     854       64090 :   pari_init_floats();
     855       64110 :   pari_init_parser();
     856       64146 :   pari_init_compiler();
     857       64107 :   pari_init_evaluator();
     858       63986 :   pari_init_files();
     859       63839 :   pari_thread_init_seadata();
     860       63789 : }
     861             : 
     862             : void
     863        3929 : pari_thread_sync(void)
     864             : {
     865        3929 :   pari_pthread_init_varstate();
     866        3929 :   pari_pthread_init_seadata();
     867        3929 : }
     868             : 
     869             : void
     870       64622 : pari_thread_close(void)
     871             : {
     872       64622 :   pari_thread_close_files();
     873       64618 :   pari_close_evaluator();
     874       64630 :   pari_close_compiler();
     875       64538 :   pari_close_parser();
     876       64453 :   pari_close_floats();
     877       64450 :   pari_close_blocks();
     878       64050 : }
     879             : 
     880             : GEN
     881       61684 : pari_thread_start(struct pari_thread *t)
     882             : {
     883       61684 :   pari_mainstack_use(&t->st);
     884       61601 :   pari_thread_init();
     885       62365 :   pari_thread_init_varstate();
     886       62829 :   return t->data;
     887             : }
     888             : 
     889             : /*********************************************************************/
     890             : /*                       LIBPARI INIT / CLOSE                        */
     891             : /*********************************************************************/
     892             : 
     893             : static void
     894           0 : pari_exit(void)
     895             : {
     896           0 :   err_printf("  ***   Error in the PARI system. End of program.\n");
     897           0 :   exit(1);
     898             : }
     899             : 
     900             : static void
     901           0 : dflt_err_recover(long errnum) { (void) errnum; pari_exit(); }
     902             : 
     903             : static void
     904           0 : dflt_pari_quit(long err) { (void)err; /*do nothing*/; }
     905             : 
     906             : static int pari_err_display(GEN err);
     907             : 
     908             : /* initialize PARI data. Initialize [new|old]fun to NULL for default set. */
     909             : void
     910        1399 : pari_init_opts(size_t parisize, ulong maxprime, ulong init_opts)
     911             : {
     912             :   ulong u;
     913             : 
     914        1399 :   pari_mt_nbthreads = 0;
     915        1399 :   cb_pari_quit = dflt_pari_quit;
     916        1399 :   cb_pari_init_histfile = NULL;
     917        1399 :   cb_pari_get_line_interactive = NULL;
     918        1399 :   cb_pari_fgets_interactive = NULL;
     919        1399 :   cb_pari_whatnow = NULL;
     920        1399 :   cb_pari_handle_exception = NULL;
     921        1399 :   cb_pari_err_handle = pari_err_display;
     922        1399 :   cb_pari_pre_recover = NULL;
     923        1399 :   cb_pari_break_loop = NULL;
     924        1399 :   cb_pari_is_interactive = NULL;
     925        1399 :   cb_pari_start_output = NULL;
     926        1399 :   cb_pari_sigint = dflt_sigint_fun;
     927        1399 :   if (init_opts&INIT_JMPm) cb_pari_err_recover = dflt_err_recover;
     928             : 
     929        1399 :   pari_stackcheck_init(&u);
     930        1399 :   pari_init_homedir();
     931        1399 :   if (init_opts&INIT_DFTm) {
     932           0 :     pari_init_defaults();
     933           0 :     GP_DATA = default_gp_data();
     934           0 :     gp_expand_path(GP_DATA->path);
     935             :   }
     936             : 
     937        1399 :   pari_mainstack = (struct pari_mainstack *) malloc(sizeof(*pari_mainstack));
     938        1399 :   paristack_alloc(parisize, 0);
     939        1399 :   init_universal_constants();
     940        1399 :   diffptr = NULL;
     941        1399 :   if (!(init_opts&INIT_noPRIMEm))  pari_init_primes(maxprime);
     942        1399 :   if (!(init_opts&INIT_noINTGMPm)) pari_kernel_init();
     943             : 
     944        1399 :   primetab = cgetalloc(t_VEC, 1);
     945        1399 :   pari_init_seadata();
     946        1399 :   pari_thread_init();
     947        1399 :   pari_init_functions();
     948        1399 :   pari_var_init();
     949        1399 :   pari_init_timer();
     950        1399 :   pari_init_buffers();
     951        1399 :   (void)getabstime();
     952        1399 :   try_to_recover = 1;
     953        1399 :   if (!(init_opts&INIT_noIMTm)) pari_mt_init();
     954        1399 :   if ((init_opts&INIT_SIGm)) pari_sig_init(pari_sighandler);
     955        1399 : }
     956             : 
     957             : void
     958           0 : pari_init(size_t parisize, ulong maxprime)
     959           0 : { pari_init_opts(parisize, maxprime, INIT_JMPm | INIT_SIGm | INIT_DFTm); }
     960             : 
     961             : void
     962        1906 : pari_close_opts(ulong init_opts)
     963             : {
     964             :   long i;
     965             : 
     966        1906 :   BLOCK_SIGINT_START;
     967        1906 :   if ((init_opts&INIT_SIGm)) pari_sig_init(SIG_DFL);
     968        1906 :   if (!(init_opts&INIT_noIMTm)) pari_mt_close();
     969             : 
     970      259216 :   for (i = 0; i < functions_tblsz; i++)
     971             :   {
     972      257310 :     entree *ep = functions_hash[i];
     973     2446746 :     while (ep) {
     974     1932126 :       entree *EP = ep->next;
     975     1932126 :       if (!EpSTATIC(ep)) { freeep(ep); free(ep); }
     976     1932126 :       ep = EP;
     977             :     }
     978             :   }
     979        1906 :   pari_var_close();
     980        1906 :   free((void*)primetab);
     981        1906 :   pari_thread_close();
     982        1906 :   pari_close_files();
     983        1906 :   pari_close_homedir();
     984        1906 :   pari_kernel_close();
     985             : 
     986        1906 :   free((void*)functions_hash);
     987        1906 :   free((void*)defaults_hash);
     988        1906 :   if (diffptr) pari_close_primes();
     989        1906 :   free(current_logfile);
     990        1906 :   free(current_psfile);
     991        1906 :   pari_mainstack_free(pari_mainstack);
     992        1906 :   free((void*)pari_mainstack);
     993        1906 :   pari_stack_delete(&s_MODULES);
     994        1906 :   if (pari_datadir) free(pari_datadir);
     995        1906 :   if (init_opts&INIT_DFTm)
     996             :   { /* delete GP_DATA */
     997        1906 :     if (GP_DATA->hist->v) free((void*)GP_DATA->hist->v);
     998        1906 :     if (GP_DATA->pp->cmd) free((void*)GP_DATA->pp->cmd);
     999        1906 :     delete_dirs(GP_DATA->path);
    1000        1906 :     free((void*)GP_DATA->path->PATH);
    1001        1906 :     delete_dirs(GP_DATA->sopath);
    1002        1906 :     free((void*)GP_DATA->sopath->PATH);
    1003        1906 :     if (GP_DATA->help) free((void*)GP_DATA->help);
    1004        1906 :     free((void*)GP_DATA->prompt);
    1005        1906 :     free((void*)GP_DATA->prompt_cont);
    1006        1906 :     free((void*)GP_DATA->histfile);
    1007             :   }
    1008        1906 :   BLOCK_SIGINT_END;
    1009        1906 : }
    1010             : 
    1011             : void
    1012        1906 : pari_close(void)
    1013        1906 : { pari_close_opts(INIT_JMPm | INIT_SIGm | INIT_DFTm); }
    1014             : 
    1015             : /*******************************************************************/
    1016             : /*                                                                 */
    1017             : /*                         ERROR RECOVERY                          */
    1018             : /*                                                                 */
    1019             : /*******************************************************************/
    1020             : void
    1021       84971 : gp_context_save(struct gp_context* rec)
    1022             : {
    1023       84971 :   rec->file = pari_last_tmp_file();
    1024       84971 :   if (DEBUGFILES>1)
    1025           0 :     err_printf("gp_context_save: %s\n", rec->file ? rec->file->name: "NULL");
    1026       84971 :   rec->prettyp = GP_DATA->fmt->prettyp;
    1027       84971 :   rec->listloc = next_block;
    1028       84971 :   rec->iferr_env = iferr_env;
    1029       84971 :   rec->err_data  = global_err_data;
    1030       84971 :   varstate_save(&rec->var);
    1031       84971 :   evalstate_save(&rec->eval);
    1032       84971 :   parsestate_save(&rec->parse);
    1033       84971 : }
    1034             : 
    1035             : void
    1036       16365 : gp_context_restore(struct gp_context* rec)
    1037             : {
    1038             :   long i;
    1039             : 
    1040       32730 :   if (!try_to_recover) return;
    1041             :   /* disable gp_context_restore() and SIGINT */
    1042       16365 :   try_to_recover = 0;
    1043       16365 :   BLOCK_SIGINT_START
    1044       16365 :   if (DEBUGMEM>2) err_printf("entering recover(), loc = %ld\n", rec->listloc);
    1045       16365 :   evalstate_restore(&rec->eval);
    1046       16365 :   parsestate_restore(&rec->parse);
    1047       16365 :   filestate_restore(rec->file);
    1048       16365 :   global_err_data = rec->err_data;
    1049       16365 :   iferr_env = rec->iferr_env;
    1050       16365 :   GP_DATA->fmt->prettyp = rec->prettyp;
    1051             : 
    1052     2225640 :   for (i = 0; i < functions_tblsz; i++)
    1053             :   {
    1054     2209275 :     entree *ep = functions_hash[i];
    1055    21618176 :     while (ep)
    1056             :     {
    1057    17199626 :       entree *EP = ep->next;
    1058    17199626 :       switch(EpVALENCE(ep))
    1059             :       {
    1060             :         case EpVAR:
    1061      176129 :           while (pop_val_if_newer(ep,rec->listloc)) /* empty */;
    1062      176129 :           break;
    1063      642132 :         case EpNEW: break;
    1064             :       }
    1065    17199626 :       ep = EP;
    1066             :     }
    1067             :   }
    1068       16365 :   varstate_restore(&rec->var);
    1069       16365 :   if (DEBUGMEM>2) err_printf("leaving recover()\n");
    1070       16365 :   BLOCK_SIGINT_END
    1071       16365 :   try_to_recover = 1;
    1072             : }
    1073             : 
    1074             : static void
    1075       16295 : err_recover(long numerr)
    1076             : {
    1077       16295 :   if (cb_pari_pre_recover)
    1078       16295 :     cb_pari_pre_recover(numerr);
    1079           0 :   evalstate_reset();
    1080           0 :   killallfiles();
    1081           0 :   pari_init_errcatch();
    1082           0 :   cb_pari_err_recover(numerr);
    1083           0 : }
    1084             : 
    1085             : static void
    1086       16743 : err_init(void)
    1087             : {
    1088             :   /* make sure pari_err msg starts at the beginning of line */
    1089       16743 :   if (!pari_last_was_newline()) pari_putc('\n');
    1090       16743 :   pariOut->flush();
    1091       16743 :   pariErr->flush();
    1092       16743 :   out_term_color(pariErr, c_ERR);
    1093       16743 : }
    1094             : 
    1095             : static void
    1096       16708 : err_init_msg(int user)
    1097             : {
    1098             :   const char *gp_function_name;
    1099       16708 :   out_puts(pariErr, "  *** ");
    1100       16708 :   if (!user && (gp_function_name = closure_func_err()))
    1101       11857 :     out_printf(pariErr, "%s: ", gp_function_name);
    1102             :   else
    1103        4851 :     out_puts(pariErr, "  ");
    1104       16708 : }
    1105             : 
    1106             : void
    1107         427 : pari_warn(int numerr, ...)
    1108             : {
    1109             :   char *ch1;
    1110             :   va_list ap;
    1111             : 
    1112         427 :   va_start(ap,numerr);
    1113             : 
    1114         427 :   err_init();
    1115         427 :   err_init_msg(numerr==warnuser || numerr==warnstack);
    1116         427 :   switch (numerr)
    1117             :   {
    1118             :     case warnuser:
    1119           7 :       out_puts(pariErr, "user warning: ");
    1120           7 :       out_print0(pariErr, NULL, va_arg(ap, GEN), f_RAW);
    1121           7 :       break;
    1122             : 
    1123             :     case warnmem:
    1124           0 :       out_puts(pariErr, "collecting garbage in "); ch1=va_arg(ap, char*);
    1125           0 :       out_vprintf(pariErr, ch1,ap); out_putc(pariErr, '.');
    1126           0 :       break;
    1127             : 
    1128             :     case warner:
    1129         420 :       out_puts(pariErr, "Warning: "); ch1=va_arg(ap, char*);
    1130         420 :       out_vprintf(pariErr, ch1,ap); out_putc(pariErr, '.');
    1131         420 :       break;
    1132             : 
    1133             :     case warnprec:
    1134           0 :       out_vprintf(pariErr, "Warning: increasing prec in %s; new prec = %ld",
    1135             :                       ap);
    1136           0 :       break;
    1137             : 
    1138             :     case warnfile:
    1139           0 :       out_puts(pariErr, "Warning: failed to "),
    1140           0 :       ch1 = va_arg(ap, char*);
    1141           0 :       out_printf(pariErr, "%s: %s", ch1, va_arg(ap, char*));
    1142           0 :       break;
    1143             : 
    1144             :     case warnstack:
    1145             :     {
    1146           0 :       ulong  s = va_arg(ap, ulong);
    1147             :       char buf[128];
    1148           0 :       sprintf(buf,"Warning: not enough memory, new stack %lu", (ulong)s);
    1149           0 :       out_puts(pariErr,buf);
    1150           0 :       break;
    1151             :     }
    1152             : 
    1153             :   }
    1154         427 :   va_end(ap);
    1155         427 :   out_term_color(pariErr, c_NONE);
    1156         427 :   out_putc(pariErr, '\n');
    1157         427 :   pariErr->flush();
    1158         427 : }
    1159             : void
    1160           0 : pari_sigint(const char *time_s)
    1161             : {
    1162           0 :   int recover=0;
    1163           0 :   BLOCK_SIGALRM_START
    1164           0 :   err_init();
    1165           0 :   closure_err(0);
    1166           0 :   err_init_msg(0);
    1167           0 :   out_puts(pariErr, "user interrupt after ");
    1168           0 :   out_puts(pariErr, time_s);
    1169           0 :   out_term_color(pariErr, c_NONE);
    1170           0 :   pariErr->flush();
    1171           0 :   if (cb_pari_handle_exception)
    1172           0 :     recover = cb_pari_handle_exception(-1);
    1173           0 :   if (!recover && !block)
    1174           0 :     PARI_SIGINT_pending = 0;
    1175           0 :   BLOCK_SIGINT_END
    1176           0 :   if (!recover) err_recover(e_MISC);
    1177           0 : }
    1178             : 
    1179             : #define retmkerr2(x,y)\
    1180             :   do { GEN _v = cgetg(3, t_ERROR);\
    1181             :        _v[1] = (x);\
    1182             :        gel(_v,2) = (y); return _v; } while(0)
    1183             : #define retmkerr3(x,y,z)\
    1184             :   do { GEN _v = cgetg(4, t_ERROR);\
    1185             :        _v[1] = (x);\
    1186             :        gel(_v,2) = (y);\
    1187             :        gel(_v,3) = (z); return _v; } while(0)
    1188             : #define retmkerr4(x,y,z,t)\
    1189             :   do { GEN _v = cgetg(5, t_ERROR);\
    1190             :        _v[1] = (x);\
    1191             :        gel(_v,2) = (y);\
    1192             :        gel(_v,3) = (z);\
    1193             :        gel(_v,4) = (t); return _v; } while(0)
    1194             : #define retmkerr5(x,y,z,t,u)\
    1195             :   do { GEN _v = cgetg(6, t_ERROR);\
    1196             :        _v[1] = (x);\
    1197             :        gel(_v,2) = (y);\
    1198             :        gel(_v,3) = (z);\
    1199             :        gel(_v,4) = (t);\
    1200             :        gel(_v,5) = (u); return _v; } while(0)
    1201             : #define retmkerr6(x,y,z,t,u,v)\
    1202             :   do { GEN _v = cgetg(7, t_ERROR);\
    1203             :        _v[1] = (x);\
    1204             :        gel(_v,2) = (y);\
    1205             :        gel(_v,3) = (z);\
    1206             :        gel(_v,4) = (t);\
    1207             :        gel(_v,5) = (u);\
    1208             :        gel(_v,6) = (v); return _v; } while(0)
    1209             : 
    1210             : static GEN
    1211       34957 : pari_err2GEN(long numerr, va_list ap)
    1212             : {
    1213       34957 :   switch ((enum err_list) numerr)
    1214             :   {
    1215             :   case e_SYNTAX:
    1216             :     {
    1217          35 :       const char *msg = va_arg(ap, char*);
    1218          35 :       const char *s = va_arg(ap,char *);
    1219          35 :       const char *entry = va_arg(ap,char *);
    1220          35 :       retmkerr3(numerr,strtoGENstr(msg), mkvecsmall2((long)s,(long)entry));
    1221             :     }
    1222             :   case e_MISC: case e_ALARM:
    1223             :     {
    1224        2692 :       const char *ch1 = va_arg(ap, char*);
    1225        2692 :       retmkerr2(numerr, gvsprintf(ch1,ap));
    1226             :     }
    1227             :   case e_NOTFUNC:
    1228             :   case e_USER:
    1229        2737 :     retmkerr2(numerr,va_arg(ap, GEN));
    1230             :   case e_FILE:
    1231             :     {
    1232           0 :       const char *f = va_arg(ap, const char*);
    1233           0 :       retmkerr3(numerr, strtoGENstr(f), strtoGENstr(va_arg(ap, char*)));
    1234             :     }
    1235             :   case e_OVERFLOW:
    1236             :   case e_IMPL:
    1237             :   case e_DIM:
    1238             :   case e_CONSTPOL:
    1239             :   case e_ROOTS0:
    1240             :   case e_FLAG:
    1241             :   case e_PREC:
    1242             :   case e_BUG:
    1243             :   case e_ARCH:
    1244             :   case e_PACKAGE:
    1245        1633 :     retmkerr2(numerr, strtoGENstr(va_arg(ap, char*)));
    1246             :   case e_MODULUS:
    1247             :   case e_VAR:
    1248             :     {
    1249         938 :       const char *f = va_arg(ap, const char*);
    1250         938 :       GEN x = va_arg(ap, GEN);
    1251         938 :       GEN y = va_arg(ap, GEN);
    1252         938 :       retmkerr4(numerr, strtoGENstr(f), x,y);
    1253             :     }
    1254             :   case e_INV:
    1255             :   case e_IRREDPOL:
    1256             :   case e_PRIME:
    1257             :   case e_SQRTN:
    1258             :   case e_TYPE:
    1259             :     {
    1260       17962 :       const char *f = va_arg(ap, const char*);
    1261       17962 :       GEN x = va_arg(ap, GEN);
    1262       17962 :       retmkerr3(numerr, strtoGENstr(f), x);
    1263             :     }
    1264             :   case e_COPRIME: case e_OP: case e_TYPE2:
    1265             :     {
    1266        3423 :       const char *f = va_arg(ap, const char*);
    1267        3423 :       GEN x = va_arg(ap, GEN);
    1268        3423 :       GEN y = va_arg(ap, GEN);
    1269        3423 :       retmkerr4(numerr,strtoGENstr(f),x,y);
    1270             :     }
    1271             :   case e_COMPONENT:
    1272             :     {
    1273         224 :       const char *f= va_arg(ap, const char *);
    1274         224 :       const char *op = va_arg(ap, const char *);
    1275         224 :       GEN l = va_arg(ap, GEN);
    1276         224 :       GEN x = va_arg(ap, GEN);
    1277         224 :       retmkerr5(numerr,strtoGENstr(f),strtoGENstr(op),l,x);
    1278             :     }
    1279             :   case e_DOMAIN:
    1280             :     {
    1281        5187 :       const char *f = va_arg(ap, const char*);
    1282        5187 :       const char *v = va_arg(ap, const char *);
    1283        5187 :       const char *op = va_arg(ap, const char *);
    1284        5187 :       GEN l = va_arg(ap, GEN);
    1285        5187 :       GEN x = va_arg(ap, GEN);
    1286        5187 :       retmkerr6(numerr,strtoGENstr(f),strtoGENstr(v),strtoGENstr(op),l,x);
    1287             :     }
    1288             :   case e_PRIORITY:
    1289             :     {
    1290         119 :       const char *f = va_arg(ap, const char*);
    1291         119 :       GEN x = va_arg(ap, GEN);
    1292         119 :       const char *op = va_arg(ap, const char *);
    1293         119 :       long v = va_arg(ap, long);
    1294         119 :       retmkerr5(numerr,strtoGENstr(f),x,strtoGENstr(op),stoi(v));
    1295             :     }
    1296             :   case e_MAXPRIME:
    1297           0 :     retmkerr2(numerr, utoi(va_arg(ap, ulong)));
    1298             :   case e_STACK:
    1299           7 :     return err_e_STACK;
    1300             :   case e_STACKTHREAD:
    1301           0 :     retmkerr3(numerr, utoi(va_arg(ap, ulong)), utoi(va_arg(ap, ulong)));
    1302             :   default:
    1303           0 :     return mkerr(numerr);
    1304             :   }
    1305             : }
    1306             : 
    1307             : static char *
    1308        6258 : type_dim(GEN x)
    1309             : {
    1310        6258 :   char *v = stack_malloc(64);
    1311        6258 :   switch(typ(x))
    1312             :   {
    1313             :     case t_MAT:
    1314             :     {
    1315          63 :       long l = lg(x), r = (l == 1)? 1: lgcols(x);
    1316          63 :       sprintf(v, "t_MAT (%ldx%ld)", r-1,l-1);
    1317          63 :       break;
    1318             :     }
    1319             :     case t_COL:
    1320          91 :       sprintf(v, "t_COL (%ld elts)", lg(x)-1);
    1321          91 :       break;
    1322             :     case t_VEC:
    1323         133 :       sprintf(v, "t_VEC (%ld elts)", lg(x)-1);
    1324         133 :       break;
    1325             :     default:
    1326        5971 :       v = (char*)type_name(typ(x));
    1327             :   }
    1328        6258 :   return v;
    1329             : }
    1330             : 
    1331             : static char *
    1332        2240 : gdisplay(GEN x)
    1333             : {
    1334        2240 :   char *s = GENtostr_raw(x);
    1335        2240 :   if (strlen(s) < 1600) return s;
    1336          21 :   if (! GP_DATA->breakloop) return (char*)"(...)";
    1337           0 :   return stack_sprintf("\n  ***  (...) Huge %s omitted; you can access it via dbg_err()", type_name(typ(x)));
    1338             : }
    1339             : 
    1340             : char *
    1341       23344 : pari_err2str(GEN e)
    1342             : {
    1343       23344 :   long numerr = err_get_num(e);
    1344       23344 :   switch ((enum err_list) numerr)
    1345             :   {
    1346             :   case e_ALARM:
    1347           0 :     return pari_sprintf("alarm interrupt after %Ps.",gel(e,2));
    1348             :   case e_MISC:
    1349        2690 :     return pari_sprintf("%Ps.",gel(e,2));
    1350             : 
    1351             :   case e_ARCH:
    1352           0 :     return pari_sprintf("sorry, '%Ps' not available on this system.",gel(e,2));
    1353             :   case e_BUG:
    1354          14 :     return pari_sprintf("bug in %Ps, please report.",gel(e,2));
    1355             :   case e_CONSTPOL:
    1356           7 :     return pari_sprintf("constant polynomial in %Ps.", gel(e,2));
    1357             :   case e_COPRIME:
    1358         189 :     return pari_sprintf("elements not coprime in %Ps:\n    %s\n    %s",
    1359         189 :                         gel(e,2), gdisplay(gel(e,3)), gdisplay(gel(e,4)));
    1360             :   case e_DIM:
    1361         925 :     return pari_sprintf("inconsistent dimensions in %Ps.", gel(e,2));
    1362             :   case e_FILE:
    1363           0 :     return pari_sprintf("error opening %Ps: `%Ps'.", gel(e,2), gel(e,3));
    1364             :   case e_FLAG:
    1365          21 :     return pari_sprintf("invalid flag in %Ps.", gel(e,2));
    1366             :   case e_IMPL:
    1367         294 :     return pari_sprintf("sorry, %Ps is not yet implemented.", gel(e,2));
    1368             :   case e_PACKAGE:
    1369           0 :     return pari_sprintf("package %Ps is required, please install it.", gel(e,2));
    1370             :   case e_INV:
    1371         490 :     return pari_sprintf("impossible inverse in %Ps: %s.", gel(e,2),
    1372         490 :                         gdisplay(gel(e,3)));
    1373             :   case e_IRREDPOL:
    1374          42 :     return pari_sprintf("not an irreducible polynomial in %Ps: %s.",
    1375          42 :                         gel(e,2), gdisplay(gel(e,3)));
    1376             :   case e_MAXPRIME:
    1377             :     {
    1378           0 :       const char * msg = "not enough precomputed primes";
    1379           0 :       ulong c = itou(gel(e,2));
    1380           0 :       if (c) return pari_sprintf("%s, need primelimit ~ %lu.",msg, c);
    1381           0 :       else   return pari_strdup(msg);
    1382             :     }
    1383             :   case e_MEM:
    1384           0 :     return pari_strdup("not enough memory");
    1385             :   case e_MODULUS:
    1386             :     {
    1387         735 :       GEN x = gel(e,3), y = gel(e,4);
    1388        1470 :       return pari_sprintf("inconsistent moduli in %Ps: %s != %s",
    1389         735 :                           gel(e,2), gdisplay(x), gdisplay(y));
    1390             :     }
    1391           0 :   case e_NONE: return NULL;
    1392             :   case e_NOTFUNC:
    1393        2723 :     return pari_strdup("not a function in function call");
    1394             :   case e_OP: case e_TYPE2:
    1395             :     {
    1396        3129 :       pari_sp av = avma;
    1397             :       char *v;
    1398        3129 :       const char *f, *op = GSTR(gel(e,2));
    1399        3129 :       const char *what = numerr == e_OP? "inconsistent": "forbidden";
    1400        3129 :       GEN x = gel(e,3);
    1401        3129 :       GEN y = gel(e,4);
    1402        3129 :       switch(*op)
    1403             :       {
    1404          14 :       case '+': f = "addition"; break;
    1405         105 :       case '*': f = "multiplication"; break;
    1406        2387 :       case '/': case '%': case '\\': f = "division"; break;
    1407           0 :       case '=': op = "-->"; f = "assignment"; break;
    1408         623 :       default:  f = op; op = ","; break;
    1409             :       }
    1410        3129 :       v = pari_sprintf("%s %s %s %s %s.", what,f,type_dim(x),op,type_dim(y));
    1411        3129 :       avma = av; return v;
    1412             :     }
    1413             :   case e_COMPONENT:
    1414             :     {
    1415         224 :       const char *f= GSTR(gel(e,2));
    1416         224 :       const char *op= GSTR(gel(e,3));
    1417         224 :       GEN l = gel(e,4);
    1418         224 :       if (!*f)
    1419         126 :         return pari_sprintf("non-existent component: index %s %Ps",op,l);
    1420          98 :       return pari_sprintf("non-existent component in %s: index %s %Ps",f,op,l);
    1421             :     }
    1422             :   case e_DOMAIN:
    1423             :     {
    1424        5089 :       const char *f = GSTR(gel(e,2));
    1425        5089 :       const char *v = GSTR(gel(e,3));
    1426        5089 :       const char *op= GSTR(gel(e,4));
    1427        5089 :       GEN l = gel(e,5);
    1428        5089 :       if (!*op)
    1429          28 :         return pari_sprintf("domain error in %s: %s out of range",f,v);
    1430        5061 :       return pari_sprintf("domain error in %s: %s %s %Ps",f,v,op,l);
    1431             :     }
    1432             :   case e_PRIORITY:
    1433             :     {
    1434          70 :       const char *f = GSTR(gel(e,2));
    1435          70 :       long vx = gvar(gel(e,3));
    1436          70 :       const char *op= GSTR(gel(e,4));
    1437          70 :       long v = itos(gel(e,5));
    1438          70 :       return pari_sprintf("incorrect priority in %s: variable %Ps %s %Ps",f,
    1439             :              pol_x(vx), op, pol_x(v));
    1440             :     }
    1441             :   case e_OVERFLOW:
    1442          77 :     return pari_sprintf("overflow in %Ps.", gel(e,2));
    1443             :   case e_PREC:
    1444         210 :     return pari_sprintf("precision too low in %Ps.", gel(e,2));
    1445             :   case e_PRIME:
    1446         112 :     return pari_sprintf("not a prime number in %Ps: %s.",
    1447         112 :                         gel(e,2), gdisplay(gel(e,3)));
    1448             :   case e_ROOTS0:
    1449          42 :     return pari_sprintf("zero polynomial in %Ps.", gel(e,2));
    1450             :   case e_SQRTN:
    1451         154 :     return pari_sprintf("not an n-th power residue in %Ps: %s.",
    1452         154 :                         gel(e,2), gdisplay(gel(e,3)));
    1453             :   case e_STACK:
    1454             :   case e_STACKTHREAD:
    1455             :     {
    1456           7 :       const char *stack = numerr == e_STACK? "PARI": "thread";
    1457           7 :       const char *var = numerr == e_STACK? "parisizemax": "threadsizemax";
    1458           8 :       size_t rsize = numerr == e_STACKTHREAD && GP_DATA->threadsize ?
    1459           7 :                                 GP_DATA->threadsize: pari_mainstack->rsize;
    1460           7 :       size_t vsize = numerr == e_STACK? pari_mainstack->vsize:
    1461           0 :                                         GP_DATA->threadsizemax;
    1462           7 :       char *buf = (char *) pari_malloc(512*sizeof(char));
    1463           7 :       if (vsize)
    1464             :       {
    1465           0 :         sprintf(buf, "the %s stack overflows !\n"
    1466             :             "  current stack size: %lu (%.3f Mbytes)\n"
    1467             :             "  [hint] you can increase '%s' using default()\n",
    1468           0 :             stack, (ulong)vsize, (double)vsize/1048576., var);
    1469             :       }
    1470             :       else
    1471             :       {
    1472           7 :         sprintf(buf, "the %s stack overflows !\n"
    1473             :             "  current stack size: %lu (%.3f Mbytes)\n"
    1474             :             "  [hint] set '%s' to a non-zero value in your GPRC\n",
    1475           7 :             stack, (ulong)rsize, (double)rsize/1048576., var);
    1476             :       }
    1477           7 :       return buf;
    1478             :     }
    1479             :   case e_SYNTAX:
    1480           0 :     return pari_strdup(GSTR(gel(e,2)));
    1481             :   case e_TYPE:
    1482       12326 :     return pari_sprintf("incorrect type in %Ps (%s).",
    1483       12326 :                         gel(e,2), type_name(typ(gel(e,3))));
    1484             :   case e_USER:
    1485          14 :     return pari_sprint0("user error: ", gel(e,2), f_RAW);
    1486             :   case e_VAR:
    1487             :     {
    1488         203 :       GEN x = gel(e,3), y = gel(e,4);
    1489         609 :       return pari_sprintf("inconsistent variables in %Ps, %Ps != %Ps.",
    1490         609 :                           gel(e,2), pol_x(varn(x)), pol_x(varn(y)));
    1491             :     }
    1492             :   }
    1493             :   return NULL; /*LCOV_EXCL_LINE*/
    1494             : }
    1495             : 
    1496             : static int
    1497       16316 : pari_err_display(GEN err)
    1498             : {
    1499       16316 :   long numerr=err_get_num(err);
    1500       16316 :   err_init();
    1501       16316 :   if (numerr==e_SYNTAX)
    1502             :   {
    1503          35 :     const char *msg = GSTR(gel(err,2));
    1504          35 :     const char *s     = (const char *) gmael(err,3,1);
    1505          35 :     const char *entry = (const char *) gmael(err,3,2);
    1506          35 :     print_errcontext(pariErr, msg, s, entry);
    1507             :   }
    1508             :   else
    1509             :   {
    1510       16281 :     char *s = pari_err2str(err);
    1511       16281 :     closure_err(0);
    1512       16281 :     err_init_msg(numerr==e_USER);
    1513       16281 :     pariErr->puts(s);
    1514       16281 :     if (numerr==e_NOTFUNC)
    1515             :     {
    1516        2723 :       GEN fun = gel(err,2);
    1517        2723 :       if (gequalX(fun))
    1518             :       {
    1519        2723 :         entree *ep = varentries[varn(fun)];
    1520        2723 :         const char *s = ep->name;
    1521        2723 :         if (cb_pari_whatnow) cb_pari_whatnow(pariErr,s,1);
    1522             :       }
    1523             :     }
    1524       16267 :     pari_free(s);
    1525             :   }
    1526       16302 :   out_term_color(pariErr, c_NONE);
    1527       16302 :   pariErr->flush(); return 0;
    1528             : }
    1529             : 
    1530             : void
    1531       34973 : pari_err(int numerr, ...)
    1532             : {
    1533             :   va_list ap;
    1534             :   GEN E;
    1535             : 
    1536       34973 :   va_start(ap,numerr);
    1537             : 
    1538       34973 :   if (numerr)
    1539       34957 :     E = pari_err2GEN(numerr,ap);
    1540             :   else
    1541             :   {
    1542          16 :     E = va_arg(ap,GEN);
    1543          16 :     numerr = err_get_num(E);
    1544             :   }
    1545       34973 :   global_err_data = E;
    1546       34973 :   if (*iferr_env) longjmp(*iferr_env, numerr);
    1547       16318 :   mt_err_recover(numerr);
    1548       16316 :   va_end(ap);
    1549       32618 :   if (cb_pari_err_handle &&
    1550       16316 :       cb_pari_err_handle(E)) return;
    1551       32597 :   if (cb_pari_handle_exception &&
    1552       16302 :       cb_pari_handle_exception(numerr)) return;
    1553       16295 :   err_recover(numerr);
    1554             : }
    1555             : 
    1556             : GEN
    1557       37299 : pari_err_last(void) { return global_err_data; }
    1558             : 
    1559             : const char *
    1560        7340 : numerr_name(long numerr)
    1561             : {
    1562        7340 :   switch ((enum err_list) numerr)
    1563             :   {
    1564           0 :   case e_ALARM:    return "e_ALARM";
    1565           0 :   case e_ARCH:     return "e_ARCH";
    1566           0 :   case e_BUG:      return "e_BUG";
    1567           0 :   case e_COMPONENT: return "e_COMPONENT";
    1568           0 :   case e_CONSTPOL: return "e_CONSTPOL";
    1569           0 :   case e_COPRIME:  return "e_COPRIME";
    1570           0 :   case e_DIM:      return "e_DIM";
    1571          56 :   case e_DOMAIN:   return "e_DOMAIN";
    1572           0 :   case e_FILE:     return "e_FILE";
    1573           7 :   case e_FLAG:     return "e_FLAG";
    1574          28 :   case e_IMPL:     return "e_IMPL";
    1575          95 :   case e_INV:      return "e_INV";
    1576           0 :   case e_IRREDPOL: return "e_IRREDPOL";
    1577           0 :   case e_MAXPRIME: return "e_MAXPRIME";
    1578           0 :   case e_MEM:      return "e_MEM";
    1579           0 :   case e_MISC:     return "e_MISC";
    1580           0 :   case e_MODULUS:  return "e_MODULUS";
    1581           0 :   case e_NONE:     return "e_NONE";
    1582           0 :   case e_NOTFUNC:  return "e_NOTFUNC";
    1583           0 :   case e_OP:       return "e_OP";
    1584           0 :   case e_OVERFLOW: return "e_OVERFLOW";
    1585           0 :   case e_PACKAGE:  return "e_PACKAGE";
    1586           0 :   case e_PREC:     return "e_PREC";
    1587           0 :   case e_PRIME:    return "e_PRIME";
    1588          49 :   case e_PRIORITY: return "e_PRIORITY";
    1589           0 :   case e_ROOTS0:   return "e_ROOTS0";
    1590           0 :   case e_SQRTN:    return "e_SQRTN";
    1591           0 :   case e_STACK:    return "e_STACK";
    1592           0 :   case e_SYNTAX:   return "e_SYNTAX";
    1593           0 :   case e_STACKTHREAD:   return "e_STACKTHREAD";
    1594           0 :   case e_TYPE2:    return "e_TYPE2";
    1595        7105 :   case e_TYPE:     return "e_TYPE";
    1596           0 :   case e_USER:     return "e_USER";
    1597           0 :   case e_VAR:      return "e_VAR";
    1598             :   }
    1599           0 :   return "invalid error number";
    1600             : }
    1601             : 
    1602             : long
    1603          21 : name_numerr(const char *s)
    1604             : {
    1605          21 :   if (!strcmp(s,"e_ALARM"))    return e_ALARM;
    1606          21 :   if (!strcmp(s,"e_ARCH"))     return e_ARCH;
    1607          21 :   if (!strcmp(s,"e_BUG"))      return e_BUG;
    1608          21 :   if (!strcmp(s,"e_COMPONENT")) return e_COMPONENT;
    1609          21 :   if (!strcmp(s,"e_CONSTPOL")) return e_CONSTPOL;
    1610          21 :   if (!strcmp(s,"e_COPRIME"))  return e_COPRIME;
    1611          21 :   if (!strcmp(s,"e_DIM"))      return e_DIM;
    1612          21 :   if (!strcmp(s,"e_DOMAIN"))   return e_DOMAIN;
    1613          21 :   if (!strcmp(s,"e_FILE"))     return e_FILE;
    1614          21 :   if (!strcmp(s,"e_FLAG"))     return e_FLAG;
    1615          21 :   if (!strcmp(s,"e_IMPL"))     return e_IMPL;
    1616          21 :   if (!strcmp(s,"e_INV"))      return e_INV;
    1617           0 :   if (!strcmp(s,"e_IRREDPOL")) return e_IRREDPOL;
    1618           0 :   if (!strcmp(s,"e_MAXPRIME")) return e_MAXPRIME;
    1619           0 :   if (!strcmp(s,"e_MEM"))      return e_MEM;
    1620           0 :   if (!strcmp(s,"e_MISC"))     return e_MISC;
    1621           0 :   if (!strcmp(s,"e_MODULUS"))  return e_MODULUS;
    1622           0 :   if (!strcmp(s,"e_NONE"))     return e_NONE;
    1623           0 :   if (!strcmp(s,"e_NOTFUNC"))  return e_NOTFUNC;
    1624           0 :   if (!strcmp(s,"e_OP"))       return e_OP;
    1625           0 :   if (!strcmp(s,"e_OVERFLOW")) return e_OVERFLOW;
    1626           0 :   if (!strcmp(s,"e_PACKAGE"))  return e_PACKAGE;
    1627           0 :   if (!strcmp(s,"e_PREC"))     return e_PREC;
    1628           0 :   if (!strcmp(s,"e_PRIME"))    return e_PRIME;
    1629           0 :   if (!strcmp(s,"e_PRIORITY")) return e_PRIORITY;
    1630           0 :   if (!strcmp(s,"e_ROOTS0"))   return e_ROOTS0;
    1631           0 :   if (!strcmp(s,"e_SQRTN"))    return e_SQRTN;
    1632           0 :   if (!strcmp(s,"e_STACK"))    return e_STACK;
    1633           0 :   if (!strcmp(s,"e_SYNTAX"))   return e_SYNTAX;
    1634           0 :   if (!strcmp(s,"e_TYPE"))     return e_TYPE;
    1635           0 :   if (!strcmp(s,"e_TYPE2"))    return e_TYPE2;
    1636           0 :   if (!strcmp(s,"e_USER"))     return e_USER;
    1637           0 :   if (!strcmp(s,"e_VAR"))      return e_VAR;
    1638           0 :   pari_err(e_MISC,"unknown error name");
    1639             :   return -1; /* LCOV_EXCL_LINE */
    1640             : }
    1641             : 
    1642             : GEN
    1643        7340 : errname(GEN err)
    1644             : {
    1645        7340 :   if (typ(err)!=t_ERROR) pari_err_TYPE("errname",err);
    1646        7340 :   return strtoGENstr(numerr_name(err_get_num(err)));
    1647             : }
    1648             : 
    1649             : /* Try f (trapping error e), recover using r (break_loop, if NULL) */
    1650             : GEN
    1651          21 : trap0(const char *e, GEN r, GEN f)
    1652             : {
    1653          21 :   long numerr = CATCH_ALL;
    1654             :   GEN x;
    1655          21 :   if (!e || !*e) numerr = CATCH_ALL;
    1656          21 :   else numerr = name_numerr(e);
    1657          21 :   if (!f) {
    1658           0 :     pari_warn(warner,"default handlers are no longer supported --> ignored");
    1659           0 :     return gnil;
    1660             :   }
    1661          21 :   x = closure_trapgen(f, numerr);
    1662          14 :   if (x == (GEN)1L) x = r? closure_evalgen(r): gnil;
    1663          14 :   return x;
    1664             : }
    1665             : 
    1666             : /*******************************************************************/
    1667             : /*                                                                */
    1668             : /*                       CLONING & COPY                            */
    1669             : /*                  Replicate an existing GEN                      */
    1670             : /*                                                                 */
    1671             : /*******************************************************************/
    1672             : /* lontyp[tx] = 0 (non recursive type) or number of codewords for type tx */
    1673             : 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 };
    1674             : 
    1675             : static GEN
    1676         690 : list_internal_copy(GEN z, long nmax)
    1677             : {
    1678             :   long i, l;
    1679             :   GEN a;
    1680         690 :   if (!z) return NULL;
    1681         560 :   l = lg(z);
    1682         560 :   a = newblock(nmax+1);
    1683         560 :   for (i = 1; i < l; i++) gel(a,i) = gel(z,i)? gclone(gel(z,i)): gen_0;
    1684         560 :   a[0] = z[0]; return a;
    1685             : }
    1686             : 
    1687             : static void
    1688         690 : listassign(GEN x, GEN y)
    1689             : {
    1690         690 :   long nmax = list_nmax(x);
    1691         690 :   GEN L = list_data(x);
    1692         690 :   if (!nmax && L) nmax = lg(L) + 32; /* not malloc'ed yet */
    1693         690 :   y[1] = evaltyp(list_typ(x))|evallg(nmax);
    1694         690 :   list_data(y) = list_internal_copy(L, nmax);
    1695         690 : }
    1696             : 
    1697             : /* copy list on the PARI stack */
    1698             : GEN
    1699         165 : listcopy(GEN x)
    1700             : {
    1701         165 :   GEN y = mklist(), L = list_data(x);
    1702         165 :   if (L) list_data(y) = gcopy(L);
    1703         165 :   y[1] = evaltyp(list_typ(x));
    1704         165 :   return y;
    1705             : }
    1706             : 
    1707             : GEN
    1708  3099346933 : gcopy(GEN x)
    1709             : {
    1710  3099346933 :   long tx = typ(x), lx, i;
    1711             :   GEN y;
    1712  3099346933 :   switch(tx)
    1713             :   { /* non recursive types */
    1714  2877473814 :     case t_INT: return signe(x)? icopy(x): gen_0;
    1715             :     case t_REAL:
    1716             :     case t_STR:
    1717   123142072 :     case t_VECSMALL: return leafcopy(x);
    1718             :     /* one more special case */
    1719         165 :     case t_LIST: return listcopy(x);
    1720             :   }
    1721    98730882 :   y = cgetg_copy(x, &lx);
    1722    98730881 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    1723    98730881 :   for (; i<lx; i++) gel(y,i) = gcopy(gel(x,i));
    1724    98730882 :   return y;
    1725             : }
    1726             : 
    1727             : /* as gcopy, but truncate to the first lx components if recursive type
    1728             :  * [ leaves use their own lg ]. No checks. */
    1729             : GEN
    1730         588 : gcopy_lg(GEN x, long lx)
    1731             : {
    1732         588 :   long tx = typ(x), i;
    1733             :   GEN y;
    1734         588 :   switch(tx)
    1735             :   { /* non recursive types */
    1736           0 :     case t_INT: return signe(x)? icopy(x): gen_0;
    1737             :     case t_REAL:
    1738             :     case t_STR:
    1739           0 :     case t_VECSMALL: return leafcopy(x);
    1740             :     /* one more special case */
    1741           0 :     case t_LIST: return listcopy(x);
    1742             :   }
    1743         588 :   y = cgetg(lx, tx);
    1744         588 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    1745         588 :   for (; i<lx; i++) gel(y,i) = gcopy(gel(x,i));
    1746         588 :   return y;
    1747             : }
    1748             : 
    1749             : /* cf cgetg_copy: "allocate" (by updating first codeword only) for subsequent
    1750             :  * copy of x, as if avma = *AVMA */
    1751             : INLINE GEN
    1752   189076296 : cgetg_copy_avma(GEN x, long *plx, pari_sp *AVMA) {
    1753             :   GEN z;
    1754   189076296 :   *plx = lg(x);
    1755   189076296 :   z = ((GEN)*AVMA) - *plx;
    1756   189076296 :   z[0] = x[0] & (TYPBITS|LGBITS);
    1757   189076296 :   *AVMA = (pari_sp)z; return z;
    1758             : }
    1759             : INLINE GEN
    1760         133 : cgetlist_avma(pari_sp *AVMA)
    1761             : {
    1762         133 :   GEN y = ((GEN)*AVMA) - 3;
    1763         133 :   y[0] = _evallg(3) | evaltyp(t_LIST);
    1764         133 :   *AVMA = (pari_sp)y; return y;
    1765             : }
    1766             : 
    1767             : /* copy x as if avma = *AVMA, update *AVMA */
    1768             : GEN
    1769  2537758712 : gcopy_avma(GEN x, pari_sp *AVMA)
    1770             : {
    1771  2537758712 :   long i, lx, tx = typ(x);
    1772             :   GEN y;
    1773             : 
    1774  2537758712 :   switch(typ(x))
    1775             :   { /* non recursive types */
    1776             :     case t_INT:
    1777  2480200550 :       *AVMA = (pari_sp)icopy_avma(x, *AVMA);
    1778  2480200504 :       return (GEN)*AVMA;
    1779             :     case t_REAL: case t_STR: case t_VECSMALL:
    1780    13988406 :       *AVMA = (pari_sp)leafcopy_avma(x, *AVMA);
    1781    13988409 :       return (GEN)*AVMA;
    1782             : 
    1783             :     /* one more special case */
    1784             :     case t_LIST:
    1785         133 :       y = cgetlist_avma(AVMA);
    1786         133 :       listassign(x, y); return y;
    1787             : 
    1788             :   }
    1789    43569623 :   y = cgetg_copy_avma(x, &lx, AVMA);
    1790    43569621 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    1791    43569621 :   for (; i<lx; i++) gel(y,i) = gcopy_avma(gel(x,i), AVMA);
    1792    43569622 :   return y;
    1793             : }
    1794             : 
    1795             : /* [copy_bin/bin_copy:] same as gcopy_avma but use NULL to code an exact 0, and
    1796             :  * make shallow copies of finalized t_LISTs */
    1797             : static GEN
    1798   777247986 : gcopy_av0(GEN x, pari_sp *AVMA)
    1799             : {
    1800   777247986 :   long i, lx, tx = typ(x);
    1801             :   GEN y;
    1802             : 
    1803   777247986 :   switch(tx)
    1804             :   { /* non recursive types */
    1805             :     case t_INT:
    1806   546975587 :       if (!signe(x)) return NULL; /* special marker */
    1807   284480659 :       *AVMA = (pari_sp)icopy_avma(x, *AVMA);
    1808   284480682 :       return (GEN)*AVMA;
    1809             :     case t_LIST:
    1810          49 :       if (list_data(x) && !list_nmax(x)) break; /* not finalized, need copy */
    1811             :       /* else finalized: shallow copy */
    1812             :     case t_REAL: case t_STR: case t_VECSMALL:
    1813    84765823 :       *AVMA = (pari_sp)leafcopy_avma(x, *AVMA);
    1814    84765860 :       return (GEN)*AVMA;
    1815             :   }
    1816   145506576 :   y = cgetg_copy_avma(x, &lx, AVMA);
    1817   145506578 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    1818   145506578 :   for (; i<lx; i++) gel(y,i) = gcopy_av0(gel(x,i), AVMA);
    1819   145506690 :   return y;
    1820             : }
    1821             : 
    1822             : INLINE GEN
    1823           0 : icopy_avma_canon(GEN x, pari_sp AVMA)
    1824             : {
    1825           0 :   long i, lx = lgefint(x);
    1826           0 :   GEN y = ((GEN)AVMA) - lx;
    1827           0 :   y[0] = evaltyp(t_INT)|evallg(lx); /* kills isclone */
    1828           0 :   y[1] = x[1]; x = int_MSW(x);
    1829           0 :   for (i=2; i<lx; i++, x = int_precW(x)) y[i] = *x;
    1830           0 :   return y;
    1831             : }
    1832             : 
    1833             : /* [copy_bin_canon:] same as gcopy_av0, but copy integers in
    1834             :  * canonical (native kernel) form and make a full copy of t_LISTs */
    1835             : static GEN
    1836           0 : gcopy_av0_canon(GEN x, pari_sp *AVMA)
    1837             : {
    1838           0 :   long i, lx, tx = typ(x);
    1839             :   GEN y;
    1840             : 
    1841           0 :   switch(tx)
    1842             :   { /* non recursive types */
    1843             :     case t_INT:
    1844           0 :       if (!signe(x)) return NULL; /* special marker */
    1845           0 :       *AVMA = (pari_sp)icopy_avma_canon(x, *AVMA);
    1846           0 :       return (GEN)*AVMA;
    1847             :     case t_REAL: case t_STR: case t_VECSMALL:
    1848           0 :       *AVMA = (pari_sp)leafcopy_avma(x, *AVMA);
    1849           0 :       return (GEN)*AVMA;
    1850             : 
    1851             :     /* one more special case */
    1852             :     case t_LIST:
    1853             :     {
    1854           0 :       long t = list_typ(x);
    1855           0 :       GEN y = cgetlist_avma(AVMA), z = list_data(x);
    1856           0 :       if (z) {
    1857           0 :         list_data(y) = gcopy_av0_canon(z, AVMA);
    1858           0 :         y[1] = evaltyp(t)|evallg(lg(z)-1);
    1859             :       } else {
    1860           0 :         list_data(y) = NULL;
    1861           0 :         y[1] = evaltyp(t);
    1862             :       }
    1863           0 :       return y;
    1864             :     }
    1865             :   }
    1866           0 :   y = cgetg_copy_avma(x, &lx, AVMA);
    1867           0 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    1868           0 :   for (; i<lx; i++) gel(y,i) = gcopy_av0_canon(gel(x,i), AVMA);
    1869           0 :   return y;
    1870             : }
    1871             : 
    1872             : /* [copy_bin/bin_copy:] size (number of words) required for
    1873             :  * gcopy_av0_canon(x) */
    1874             : static long
    1875           0 : taille0_canon(GEN x)
    1876             : {
    1877           0 :   long i,n,lx, tx = typ(x);
    1878           0 :   switch(tx)
    1879             :   { /* non recursive types */
    1880           0 :     case t_INT: return signe(x)? lgefint(x): 0;
    1881             :     case t_REAL:
    1882             :     case t_STR:
    1883           0 :     case t_VECSMALL: return lg(x);
    1884             : 
    1885             :     /* one more special case */
    1886             :     case t_LIST:
    1887             :     {
    1888           0 :       GEN L = list_data(x);
    1889           0 :       return L? 3 + taille0_canon(L): 3;
    1890             :     }
    1891             :   }
    1892           0 :   n = lx = lg(x);
    1893           0 :   for (i=lontyp[tx]; i<lx; i++) n += taille0_canon(gel(x,i));
    1894           0 :   return n;
    1895             : }
    1896             : 
    1897             : /* [copy_bin/bin_copy:] size (number of words) required for gcopy_av0(x) */
    1898             : static long
    1899   777248492 : taille0(GEN x)
    1900             : {
    1901   777248492 :   long i,n,lx, tx = typ(x);
    1902   777248492 :   switch(tx)
    1903             :   { /* non recursive types */
    1904             :     case t_INT:
    1905   546975930 :       lx = lgefint(x);
    1906   546975930 :       return lx == 2? 0: lx;
    1907             :     case t_LIST:
    1908             :     {
    1909          49 :       GEN L = list_data(x);
    1910          49 :       if (L && !list_nmax(x)) break; /* not finalized, deep copy */
    1911             :     }
    1912             :     /* else finalized: shallow */
    1913             :     case t_REAL:
    1914             :     case t_STR:
    1915             :     case t_VECSMALL:
    1916    84765851 :       return lg(x);
    1917             :   }
    1918   145506711 :   n = lx = lg(x);
    1919   145506711 :   for (i=lontyp[tx]; i<lx; i++) n += taille0(gel(x,i));
    1920   145506718 :   return n;
    1921             : }
    1922             : 
    1923             : /* How many words do we need to allocate to copy x ? t_LIST is a special case
    1924             :  * since list_data() is malloc'ed later, in list_internal_copy() */
    1925             : static long
    1926  2675412063 : words_to_allocate(GEN x)
    1927             : {
    1928  2675412063 :   long i,n,lx, tx = typ(x);
    1929  2675412063 :   switch(tx)
    1930             :   { /* non recursive types */
    1931  2575556515 :     case t_INT: return lgefint(x);
    1932             :     case t_REAL:
    1933             :     case t_STR:
    1934    19770734 :     case t_VECSMALL: return lg(x);
    1935             : 
    1936         697 :     case t_LIST: return 3;
    1937             :     default:
    1938    80084117 :       n = lx = lg(x);
    1939    80084117 :       for (i=lontyp[tx]; i<lx; i++) n += words_to_allocate(gel(x,i));
    1940    80084068 :       return n;
    1941             :   }
    1942             : }
    1943             : 
    1944             : long
    1945        6636 : gsizeword(GEN x)
    1946             : {
    1947             :   GEN L;
    1948        6636 :   if (typ(x) != t_LIST) return words_to_allocate(x);
    1949             :   /* For t_LIST, return the actual list size, words_to_allocate() is always 3 */
    1950         147 :   L = list_data(x);
    1951         147 :   return L? 3 + words_to_allocate(L): 3;
    1952             : }
    1953             : long
    1954          21 : gsizebyte(GEN x) { return gsizeword(x) * sizeof(long); }
    1955             : 
    1956             : /* return a clone of x structured as a gcopy */
    1957             : GENbin*
    1958    40457210 : copy_bin(GEN x)
    1959             : {
    1960    40457210 :   long t = taille0(x);
    1961    40457213 :   GENbin *p = (GENbin*)pari_malloc(sizeof(GENbin) + t*sizeof(long));
    1962    40457244 :   pari_sp AVMA = (pari_sp)(GENbinbase(p) + t);
    1963    40457230 :   p->rebase = &shiftaddress;
    1964    40457230 :   p->len = t;
    1965    40457230 :   p->x   = gcopy_av0(x, &AVMA);
    1966    40457213 :   p->base= (GEN)AVMA; return p;
    1967             : }
    1968             : 
    1969             : /* same, writing t_INT in canonical native form */
    1970             : GENbin*
    1971           0 : copy_bin_canon(GEN x)
    1972             : {
    1973           0 :   long t = taille0_canon(x);
    1974           0 :   GENbin *p = (GENbin*)pari_malloc(sizeof(GENbin) + t*sizeof(long));
    1975           0 :   pari_sp AVMA = (pari_sp)(GENbinbase(p) + t);
    1976           0 :   p->rebase = &shiftaddress_canon;
    1977           0 :   p->len = t;
    1978           0 :   p->x   = gcopy_av0_canon(x, &AVMA);
    1979           0 :   p->base= (GEN)AVMA; return p;
    1980             : }
    1981             : 
    1982             : GEN
    1983   137566223 : gclone(GEN x)
    1984             : {
    1985   137566223 :   long i,lx,tx = typ(x), t = words_to_allocate(x);
    1986   137566223 :   GEN y = newblock(t);
    1987   137566225 :   switch(tx)
    1988             :   { /* non recursive types */
    1989             :     case t_INT:
    1990    95219280 :       lx = lgefint(x);
    1991    95219280 :       y[0] = evaltyp(t_INT)|evallg(lx);
    1992    95219280 :       for (i=1; i<lx; i++) y[i] = x[i];
    1993    95219280 :       break;
    1994             :     case t_REAL:
    1995             :     case t_STR:
    1996             :     case t_VECSMALL:
    1997     5869812 :       lx = lg(x);
    1998     5869812 :       for (i=0; i<lx; i++) y[i] = x[i];
    1999     5869812 :       break;
    2000             : 
    2001             :     /* one more special case */
    2002             :     case t_LIST:
    2003         557 :       y[0] = evaltyp(t_LIST)|_evallg(3);
    2004         557 :       listassign(x, y);
    2005         557 :       break;
    2006             :     default: {
    2007    36476576 :       pari_sp AVMA = (pari_sp)(y + t);
    2008    36476576 :       lx = lg(x);
    2009    36476576 :       y[0] = x[0];
    2010    36476576 :       if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    2011    36476576 :       for (; i<lx; i++) gel(y,i) = gcopy_avma(gel(x,i), &AVMA);
    2012             :     }
    2013             :   }
    2014   137566223 :   setisclone(y); return y;
    2015             : }
    2016             : 
    2017             : void
    2018   514753121 : shiftaddress(GEN x, long dec)
    2019             : {
    2020   514753121 :   long i, lx, tx = typ(x);
    2021   514753121 :   if (is_recursive_t(tx))
    2022             :   {
    2023   145506760 :     if (tx == t_LIST)
    2024             :     {
    2025   514753235 :       if (!list_data(x) || list_nmax(x)) return; /* empty or finalized */
    2026             :       /* not finalized, update pointers  */
    2027             :     }
    2028   145506718 :     lx = lg(x);
    2029   882297828 :     for (i=lontyp[tx]; i<lx; i++) {
    2030   736791051 :       if (!x[i]) gel(x,i) = gen_0;
    2031             :       else
    2032             :       {
    2033   474298535 :         x[i] += dec;
    2034   474298535 :         shiftaddress(gel(x,i), dec);
    2035             :       }
    2036             :     }
    2037             :   }
    2038             : }
    2039             : 
    2040             : void
    2041           0 : shiftaddress_canon(GEN x, long dec)
    2042             : {
    2043           0 :   long i, lx, tx = typ(x);
    2044           0 :   switch(tx)
    2045             :   { /* non recursive types */
    2046             :     case t_INT: {
    2047             :       GEN y;
    2048           0 :       lx = lgefint(x); if (lx <= 3) return;
    2049           0 :       y = x + 2;
    2050           0 :       x = int_MSW(x);  if (x == y) return;
    2051           0 :       while (x > y) { lswap(*x, *y); x = int_precW(x); y++; }
    2052           0 :       break;
    2053             :     }
    2054             :     case t_REAL:
    2055             :     case t_STR:
    2056             :     case t_VECSMALL:
    2057           0 :       break;
    2058             : 
    2059             :     /* one more special case */
    2060             :     case t_LIST: {
    2061           0 :       GEN Lx = list_data(x);
    2062           0 :       if (Lx) {
    2063           0 :         pari_sp av = avma;
    2064           0 :         GEN L = (GEN)((long)Lx+dec);
    2065           0 :         shiftaddress_canon(L, dec);
    2066           0 :         list_data(x) = list_internal_copy(L, lg(L)); avma = av;
    2067             :       }
    2068             :     }
    2069             :     default:
    2070           0 :       lx = lg(x);
    2071           0 :       for (i=lontyp[tx]; i<lx; i++) {
    2072           0 :         if (!x[i]) gel(x,i) = gen_0;
    2073             :         else
    2074             :         {
    2075           0 :           x[i] += dec;
    2076           0 :           shiftaddress_canon(gel(x,i), dec);
    2077             :         }
    2078             :       }
    2079             :   }
    2080             : }
    2081             : 
    2082             : /********************************************************************/
    2083             : /**                                                                **/
    2084             : /**                INSERT DYNAMIC OBJECT IN STRUCTURE              **/
    2085             : /**                                                                **/
    2086             : /********************************************************************/
    2087             : GEN
    2088          21 : obj_reinit(GEN S)
    2089             : {
    2090          21 :   GEN s, T = leafcopy(S);
    2091          21 :   long a = lg(T)-1;
    2092          21 :   s = gel(T,a);
    2093          21 :   gel(T,a) = zerovec(lg(s)-1);
    2094          21 :   return T;
    2095             : }
    2096             : 
    2097             : GEN
    2098     1834228 : obj_init(long d, long n)
    2099             : {
    2100     1834228 :   GEN S = cgetg(d+2, t_VEC);
    2101     1834228 :   gel(S, d+1) = zerovec(n);
    2102     1834229 :   return S;
    2103             : }
    2104             : /* insert O in S [last position] at position K, return it */
    2105             : GEN
    2106     1086312 : obj_insert(GEN S, long K, GEN O)
    2107     1086312 : { return obj_insert_shallow(S, K, gclone(O)); }
    2108             : /* as obj_insert. WITHOUT cloning (for libpari, when creating a *new* obj
    2109             :  * from an existing one) */
    2110             : GEN
    2111     1090141 : obj_insert_shallow(GEN S, long K, GEN O)
    2112             : {
    2113     1090141 :   GEN o, v = gel(S, lg(S)-1);
    2114     1090141 :   if (typ(v) != t_VEC) pari_err_TYPE("obj_insert", S);
    2115     1090141 :   o = gel(v,K);
    2116     1090141 :   gel(v,K) = O; /*SIGINT: before unclone(o)*/
    2117     1090141 :   if (isclone(o)) gunclone(o);
    2118     1090141 :   return gel(v,K);
    2119             : }
    2120             : 
    2121             : /* Does S [last position] contain data at position K ? Return it, or NULL */
    2122             : GEN
    2123     2394844 : obj_check(GEN S, long K)
    2124             : {
    2125     2394844 :   GEN O, v = gel(S, lg(S)-1);
    2126     2394844 :   if (typ(v) != t_VEC || K >= lg(v)) pari_err_TYPE("obj_check", S);
    2127     2394845 :   O = gel(v,K); return isintzero(O)? NULL: O;
    2128             : }
    2129             : 
    2130             : GEN
    2131      723088 : obj_checkbuild(GEN S, long tag, GEN (*build)(GEN))
    2132             : {
    2133      723088 :   GEN O = obj_check(S, tag);
    2134      723089 :   if (!O)
    2135      628635 :   { pari_sp av = avma; O = obj_insert(S, tag, build(S)); avma = av; }
    2136      723089 :   return O;
    2137             : }
    2138             : 
    2139             : GEN
    2140       27112 : obj_checkbuild_prec(GEN S, long tag, GEN (*build)(GEN,long),
    2141             :   long (*pr)(GEN), long prec)
    2142             : {
    2143       27112 :   pari_sp av = avma;
    2144       27112 :   GEN w = obj_check(S, tag);
    2145       27112 :   if (!w || pr(w) < prec) w = obj_insert(S, tag, build(S, prec));
    2146       27112 :   avma = av; return gcopy(w);
    2147             : }
    2148             : GEN
    2149        4341 : obj_checkbuild_realprec(GEN S, long tag, GEN (*build)(GEN,long), long prec)
    2150        4341 : { return obj_checkbuild_prec(S,tag,build,gprecision,prec); }
    2151             : GEN
    2152         455 : obj_checkbuild_padicprec(GEN S, long tag, GEN (*build)(GEN,long), long prec)
    2153         455 : { return obj_checkbuild_prec(S,tag,build,padicprec_relative,prec); }
    2154             : 
    2155             : /* Reset S [last position], freeing all clones */
    2156             : void
    2157        5096 : obj_free(GEN S)
    2158             : {
    2159        5096 :   GEN v = gel(S, lg(S)-1);
    2160             :   long i;
    2161        5096 :   if (typ(v) != t_VEC) pari_err_TYPE("obj_free", S);
    2162       28980 :   for (i = 1; i < lg(v); i++)
    2163             :   {
    2164       23884 :     GEN o = gel(v,i);
    2165       23884 :     if (isclone(o)) gunclone(o);
    2166       23884 :     gel(v,i) = gen_0;
    2167             :   }
    2168        5096 : }
    2169             : 
    2170             : /*******************************************************************/
    2171             : /*                                                                 */
    2172             : /*                         STACK MANAGEMENT                        */
    2173             : /*                                                                 */
    2174             : /*******************************************************************/
    2175             : INLINE void
    2176  2019792709 : dec_gerepile(pari_sp *x, pari_sp av0, pari_sp av, pari_sp tetpil, size_t dec)
    2177             : {
    2178  2019792709 :   if (*x < av && *x >= av0)
    2179             :   { /* update address if in stack */
    2180  1652313936 :     if (*x < tetpil) *x += dec;
    2181           0 :     else pari_err_BUG("gerepile, significant pointers lost");
    2182             :   }
    2183  2019792709 : }
    2184             : 
    2185             : void
    2186      187894 : gerepileallsp(pari_sp av, pari_sp tetpil, int n, ...)
    2187             : {
    2188      187894 :   const pari_sp av0 = avma;
    2189      187894 :   const size_t dec = av-tetpil;
    2190             :   int i;
    2191      187894 :   va_list a; va_start(a, n);
    2192      187894 :   (void)gerepile(av,tetpil,NULL);
    2193      187894 :   for (i=0; i<n; i++) dec_gerepile((pari_sp*)va_arg(a,GEN*), av0,av,tetpil,dec);
    2194      187894 :   va_end(a);
    2195      187894 : }
    2196             : 
    2197             : /* Takes an array of pointers to GENs, of length n.
    2198             :  * Cleans up the stack between av and tetpil, updating those GENs. */
    2199             : void
    2200     5263924 : gerepilemanysp(pari_sp av, pari_sp tetpil, GEN* gptr[], int n)
    2201             : {
    2202     5263924 :   const pari_sp av0 = avma;
    2203     5263924 :   const size_t dec = av-tetpil;
    2204             :   int i;
    2205     5263924 :   (void)gerepile(av,tetpil,NULL);
    2206     5263924 :   for (i=0; i<n; i++) dec_gerepile((pari_sp*)gptr[i], av0, av, tetpil, dec);
    2207     5263924 : }
    2208             : 
    2209             : /* Takes an array of GENs (cast to longs), of length n.
    2210             :  * Cleans up the stack between av and tetpil, updating those GENs. */
    2211             : void
    2212   124010319 : gerepilecoeffssp(pari_sp av, pari_sp tetpil, long *g, int n)
    2213             : {
    2214   124010319 :   const pari_sp av0 = avma;
    2215   124010319 :   const size_t dec = av-tetpil;
    2216             :   int i;
    2217   124010319 :   (void)gerepile(av,tetpil,NULL);
    2218   124010319 :   for (i=0; i<n; i++,g++) dec_gerepile((pari_sp*)g, av0, av, tetpil, dec);
    2219   124010319 : }
    2220             : 
    2221             : static int
    2222           0 : dochk_gerepileupto(GEN av, GEN x)
    2223             : {
    2224             :   long i,lx,tx;
    2225           0 :   if (!isonstack(x)) return 1;
    2226           0 :   if (x > av)
    2227             :   {
    2228           0 :     pari_warn(warner,"bad object %Ps",x);
    2229           0 :     return 0;
    2230             :   }
    2231           0 :   tx = typ(x);
    2232           0 :   if (! is_recursive_t(tx)) return 1;
    2233             : 
    2234           0 :   lx = lg(x);
    2235           0 :   for (i=lontyp[tx]; i<lx; i++)
    2236           0 :     if (!dochk_gerepileupto(av, gel(x,i)))
    2237             :     {
    2238           0 :       pari_warn(warner,"bad component %ld in object %Ps",i,x);
    2239           0 :       return 0;
    2240             :     }
    2241           0 :   return 1;
    2242             : }
    2243             : /* check that x and all its components are out of stack, or have been
    2244             :  * created after av */
    2245             : int
    2246           0 : chk_gerepileupto(GEN x) { return dochk_gerepileupto(x, x); }
    2247             : 
    2248             : /* print stack between avma & av */
    2249             : void
    2250           0 : dbg_gerepile(pari_sp av)
    2251             : {
    2252           0 :   GEN x = (GEN)avma;
    2253           0 :   while (x < (GEN)av)
    2254             :   {
    2255           0 :     const long tx = typ(x), lx = lg(x);
    2256             :     GEN *a;
    2257             : 
    2258           0 :     pari_printf(" [%ld] %Ps:", x - (GEN)avma, x);
    2259           0 :     if (! is_recursive_t(tx)) { pari_putc('\n'); x += lx; continue; }
    2260           0 :     a = (GEN*)x + lontyp[tx]; x += lx;
    2261           0 :     for (  ; a < (GEN*)x; a++)
    2262             :     {
    2263           0 :       if (*a == gen_0)
    2264           0 :         pari_puts("  gen_0");
    2265           0 :       else if (*a == gen_1)
    2266           0 :         pari_puts("  gen_1");
    2267           0 :       else if (*a == gen_m1)
    2268           0 :         pari_puts("  gen_m1");
    2269           0 :       else if (*a == gen_2)
    2270           0 :         pari_puts("  gen_2");
    2271           0 :       else if (*a == gen_m2)
    2272           0 :         pari_puts("  gen_m2");
    2273           0 :       else if (*a == ghalf)
    2274           0 :         pari_puts("  ghalf");
    2275           0 :       else if (isclone(*a))
    2276           0 :         pari_printf("  %Ps (clone)", *a);
    2277             :       else
    2278           0 :         pari_printf("  %Ps [%ld]", *a, *a - (GEN)avma);
    2279           0 :       if (a+1 < (GEN*)x) pari_putc(',');
    2280             :     }
    2281           0 :     pari_printf("\n");
    2282             :   }
    2283           0 : }
    2284             : void
    2285           0 : dbg_gerepileupto(GEN q)
    2286             : {
    2287           0 :   err_printf("%Ps:\n", q);
    2288           0 :   dbg_gerepile((pari_sp) (q+lg(q)));
    2289           0 : }
    2290             : 
    2291             : GEN
    2292   427880207 : gerepile(pari_sp av, pari_sp tetpil, GEN q)
    2293             : {
    2294   427880207 :   const size_t dec = av - tetpil;
    2295   427880207 :   const pari_sp av0 = avma;
    2296             :   GEN x, a;
    2297             : 
    2298   427880207 :   if (dec == 0) return q;
    2299   364888773 :   if ((long)dec < 0) pari_err(e_MISC,"lbot>ltop in gerepile");
    2300             : 
    2301             :   /* dec_gerepile(&q, av0, av, tetpil, dec), saving 1 comparison */
    2302   364862690 :   if (q >= (GEN)av0 && q < (GEN)tetpil)
    2303   239009267 :     q = (GEN) (((pari_sp)q) + dec);
    2304             : 
    2305   364862690 :   for (x = (GEN)av, a = (GEN)tetpil; a > (GEN)av0; ) *--x = *--a;
    2306   364862690 :   avma = (pari_sp)x;
    2307  2705074796 :   while (x < (GEN)av)
    2308             :   {
    2309  1975304458 :     const long tx = typ(x), lx = lg(x);
    2310             : 
    2311  1975304458 :     if (! is_recursive_t(tx)) { x += lx; continue; }
    2312   414673798 :     a = x + lontyp[tx]; x += lx;
    2313   414673798 :     for (  ; a < x; a++) dec_gerepile((pari_sp*)a, av0, av, tetpil, dec);
    2314             :   }
    2315   364907648 :   return q;
    2316             : }
    2317             : 
    2318             : void
    2319           0 : fill_stack(void)
    2320             : {
    2321           0 :   GEN x = ((GEN)pari_mainstack->bot);
    2322           0 :   while (x < (GEN)avma) *x++ = 0xfefefefeUL;
    2323           0 : }
    2324             : 
    2325             : void
    2326           0 : debug_stack(void)
    2327             : {
    2328           0 :   pari_sp top = pari_mainstack->top, bot = pari_mainstack->bot;
    2329             :   GEN z;
    2330           0 :   err_printf("bot=0x%lx\ttop=0x%lx\tavma=0x%lx\n", bot, top, avma);
    2331           0 :   for (z = ((GEN)top)-1; z >= (GEN)avma; z--)
    2332           0 :     err_printf("%p:\t0x%lx\t%lu\n",z,*z,*z);
    2333           0 : }
    2334             : 
    2335             : void
    2336           0 : setdebugvar(long n) { DEBUGVAR=n; }
    2337             : 
    2338             : long
    2339           0 : getdebugvar(void) { return DEBUGVAR; }
    2340             : 
    2341             : long
    2342           7 : getstack(void) { return pari_mainstack->top-avma; }
    2343             : 
    2344             : /*******************************************************************/
    2345             : /*                                                                 */
    2346             : /*                               timer_delay                             */
    2347             : /*                                                                 */
    2348             : /*******************************************************************/
    2349             : 
    2350             : #if defined(USE_CLOCK_GETTIME)
    2351             : #if defined(_POSIX_THREAD_CPUTIME)
    2352             : static THREAD clockid_t time_type = CLOCK_THREAD_CPUTIME_ID;
    2353             : #else
    2354             : static const THREAD clockid_t time_type = CLOCK_PROCESS_CPUTIME_ID;
    2355             : #endif
    2356             : static void
    2357             : pari_init_timer(void)
    2358             : {
    2359             : #if defined(_POSIX_THREAD_CPUTIME)
    2360             :   time_type = CLOCK_PROCESS_CPUTIME_ID;
    2361             : #endif
    2362             : }
    2363             : 
    2364             : void
    2365             : timer_start(pari_timer *T)
    2366             : {
    2367             :   struct timespec t;
    2368             :   clock_gettime(time_type,&t);
    2369             :   T->us = t.tv_nsec / 1000;
    2370             :   T->s  = t.tv_sec;
    2371             : }
    2372             : #elif defined(USE_GETRUSAGE)
    2373             : #ifdef RUSAGE_THREAD
    2374             : static THREAD int rusage_type = RUSAGE_THREAD;
    2375             : #else
    2376             : static const THREAD int rusage_type = RUSAGE_SELF;
    2377             : #endif /*RUSAGE_THREAD*/
    2378             : static void
    2379        1399 : pari_init_timer(void)
    2380             : {
    2381             : #ifdef RUSAGE_THREAD
    2382        1399 :   rusage_type = RUSAGE_SELF;
    2383             : #endif
    2384        1399 : }
    2385             : 
    2386             : void
    2387      462346 : timer_start(pari_timer *T)
    2388             : {
    2389             :   struct rusage r;
    2390      462346 :   getrusage(rusage_type,&r);
    2391      462347 :   T->us = r.ru_utime.tv_usec;
    2392      462347 :   T->s  = r.ru_utime.tv_sec;
    2393      462347 : }
    2394             : #elif defined(USE_FTIME)
    2395             : 
    2396             : static void
    2397             : pari_init_timer(void) { }
    2398             : 
    2399             : void
    2400             : timer_start(pari_timer *T)
    2401             : {
    2402             :   struct timeb t;
    2403             :   ftime(&t);
    2404             :   T->us = ((long)t.millitm) * 1000;
    2405             :   T->s  = t.time;
    2406             : }
    2407             : 
    2408             : #else
    2409             : 
    2410             : static void
    2411             : _get_time(pari_timer *T, long Ticks, long TickPerSecond)
    2412             : {
    2413             :   T->us = (long) ((Ticks % TickPerSecond) * (1000000. / TickPerSecond));
    2414             :   T->s  = Ticks / TickPerSecond;
    2415             : }
    2416             : 
    2417             : # ifdef USE_TIMES
    2418             : static void
    2419             : pari_init_timer(void) { }
    2420             : 
    2421             : void
    2422             : timer_start(pari_timer *T)
    2423             : {
    2424             : # ifdef _SC_CLK_TCK
    2425             :   long tck = sysconf(_SC_CLK_TCK);
    2426             : # else
    2427             :   long tck = CLK_TCK;
    2428             : # endif
    2429             :   struct tms t; times(&t);
    2430             :   _get_time(T, t.tms_utime, tck);
    2431             : }
    2432             : # elif defined(_WIN32)
    2433             : static void
    2434             : pari_init_timer(void) { }
    2435             : 
    2436             : void
    2437             : timer_start(pari_timer *T)
    2438             : { _get_time(T, win32_timer(), 1000); }
    2439             : # else
    2440             : #  include <time.h>
    2441             : #  ifndef CLOCKS_PER_SEC
    2442             : #   define CLOCKS_PER_SEC 1000000 /* may be false on YOUR system */
    2443             : #  endif
    2444             : static void
    2445             : pari_init_timer(void) { }
    2446             : 
    2447             : void
    2448             : timer_start(pari_timer *T)
    2449             : { _get_time(T, clock(), CLOCKS_PER_SEC); }
    2450             : # endif
    2451             : #endif
    2452             : 
    2453             : static long
    2454       58606 : timer_aux(pari_timer *T, pari_timer *U)
    2455             : {
    2456       58606 :   long s = T->s, us = T->us; timer_start(U);
    2457       58606 :   return 1000 * (U->s - s) + (U->us - us + 500) / 1000;
    2458             : }
    2459             : /* return delay, reset timer */
    2460             : long
    2461       57205 : timer_delay(pari_timer *T) { return timer_aux(T, T); }
    2462             : /* return delay, don't reset timer */
    2463             : long
    2464        1401 : timer_get(pari_timer *T) { pari_timer t; return timer_aux(T, &t); }
    2465             : 
    2466             : static void
    2467           0 : timer_vprintf(pari_timer *T, const char *format, va_list args)
    2468             : {
    2469           0 :   out_puts(pariErr, "Time ");
    2470           0 :   out_vprintf(pariErr, format,args);
    2471           0 :   out_printf(pariErr, ": %ld\n", timer_delay(T));
    2472           0 :   pariErr->flush();
    2473           0 : }
    2474             : void
    2475           0 : timer_printf(pari_timer *T, const char *format, ...)
    2476             : {
    2477           0 :   va_list args; va_start(args, format);
    2478           0 :   timer_vprintf(T, format, args);
    2479           0 :   va_end(args);
    2480           0 : }
    2481             : 
    2482             : long
    2483           0 : timer(void)  { static THREAD pari_timer T; return timer_delay(&T);}
    2484             : long
    2485        3278 : gettime(void)  { static THREAD pari_timer T; return timer_delay(&T);}
    2486             : 
    2487             : static THREAD pari_timer timer2_T, abstimer_T;
    2488             : long
    2489           0 : timer2(void) {  return timer_delay(&timer2_T);}
    2490             : void
    2491           0 : msgtimer(const char *format, ...)
    2492             : {
    2493           0 :   va_list args; va_start(args, format);
    2494           0 :   timer_vprintf(&timer2_T, format, args);
    2495           0 :   va_end(args);
    2496           0 : }
    2497             : long
    2498        1399 : getabstime(void)  { return timer_get(&abstimer_T);}
    2499             : #if defined(USE_CLOCK_GETTIME) || defined(USE_GETTIMEOFDAY) \
    2500             :  || defined(USE_FTIMEFORWALLTIME)
    2501             : static GEN
    2502           0 : timetoi(ulong s, ulong m)
    2503             : {
    2504           0 :   pari_sp av = avma;
    2505           0 :   GEN r = addiu(muliu(utoi(s), 1000), m);
    2506           0 :   return gerepileuptoint(av, r);
    2507             : }
    2508             : #endif
    2509             : GEN
    2510           0 : getwalltime(void)
    2511             : {
    2512             : #if defined(USE_CLOCK_GETTIME)
    2513             :   struct timespec t;
    2514             :   if (!clock_gettime(CLOCK_REALTIME,&t))
    2515             :     return timetoi(t.tv_sec, (t.tv_nsec + 500000)/1000000);
    2516             : #elif defined(USE_GETTIMEOFDAY)
    2517             :   struct timeval tv;
    2518           0 :   if (!gettimeofday(&tv, NULL))
    2519           0 :     return timetoi(tv.tv_sec, (tv.tv_usec + 500)/1000);
    2520             : #elif defined(USE_FTIMEFORWALLTIME)
    2521             :   struct timeb tp;
    2522             :   ftime(&tp); return timetoi(tp.time, tp.millitm);
    2523             : #endif
    2524           0 :   return utoi(getabstime());
    2525             : }
    2526             : 
    2527             : /*******************************************************************/
    2528             : /*                                                                 */
    2529             : /*                   FUNCTIONS KNOWN TO THE ANALYZER               */
    2530             : /*                                                                 */
    2531             : /*******************************************************************/
    2532             : GEN
    2533           7 : pari_version(void)
    2534             : {
    2535           7 :   const ulong mask = (1UL<<PARI_VERSION_SHIFT) - 1;
    2536           7 :   ulong major, minor, patch, n = paricfg_version_code;
    2537           7 :   patch = n & mask; n >>= PARI_VERSION_SHIFT;
    2538           7 :   minor = n & mask; n >>= PARI_VERSION_SHIFT;
    2539           7 :   major = n;
    2540           7 :   if (*paricfg_vcsversion) {
    2541           7 :     const char *ver = paricfg_vcsversion;
    2542           7 :     const char *s = strchr(ver, '-');
    2543             :     char t[8];
    2544           7 :     const long len = s-ver;
    2545             :     GEN v;
    2546           7 :     if (!s || len > 6) pari_err_BUG("pari_version()"); /* paranoia */
    2547           7 :     memcpy(t, ver, len); t[len] = 0;
    2548           7 :     v = cgetg(6, t_VEC);
    2549           7 :     gel(v,1) = utoi(major);
    2550           7 :     gel(v,2) = utoi(minor);
    2551           7 :     gel(v,3) = utoi(patch);
    2552           7 :     gel(v,4) = stoi( atoi(t) );
    2553           7 :     gel(v,5) = strtoGENstr(s+1);
    2554           7 :     return v;
    2555             :   } else {
    2556           0 :     GEN v = cgetg(4, t_VEC);
    2557           0 :     gel(v,1) = utoi(major);
    2558           0 :     gel(v,2) = utoi(minor);
    2559           0 :     gel(v,3) = utoi(patch);
    2560           0 :     return v;
    2561             :   }
    2562             : }
    2563             : 
    2564             : /* List of GP functions: generated from the description system.
    2565             :  * Format (struct entree) :
    2566             :  *   char *name   : name (under GP).
    2567             :  *   ulong valence: (EpNEW, EpALIAS,EpVAR, EpINSTALL)|EpSTATIC
    2568             :  *   void *value  : For PREDEFINED FUNCTIONS: C function to call.
    2569             :  *                  For USER FUNCTIONS: pointer to defining data (block) =
    2570             :  *                   entree*: NULL, list of entree (arguments), NULL
    2571             :  *                   char*  : function text
    2572             :  *   long menu    : which help section do we belong to
    2573             :  *                   1: Standard monadic or dyadic OPERATORS
    2574             :  *                   2: CONVERSIONS and similar elementary functions
    2575             :  *                   3: TRANSCENDENTAL functions, etc.
    2576             :  *   char *code   : GP prototype, aka Parser Code (see libpari's manual)
    2577             :  *                  if NULL, use valence instead.
    2578             :  *   char *help   : short help text (init to NULL).
    2579             :  *   void *pvalue : push_val history.
    2580             :  *   long arity   : maximum number of arguments.
    2581             :  *   entree *next : next entree (init to NULL, used in hashing code). */
    2582             : #include "init.h"
    2583             : #include "default.h"

Generated by: LCOV version 1.11