1 /* Copyright (C) 2004  The PARI group.
2 
3 This file is part of the PARI/GP package.
4 
5 PARI/GP is free software; you can redistribute it and/or modify it under the
6 terms of the GNU General Public License as published by the Free Software
7 Foundation; either version 2 of the License, or (at your option) any later
8 version. It is distributed in the hope that it will be useful, but WITHOUT
9 ANY WARRANTY WHATSOEVER.
10 
11 Check the License for details. You should have received a copy of it, along
12 with the package; see the file 'COPYING'. If not, write to the Free Software
13 Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
14 
15 BEGINEXTERN
16 
17 /* for qsort */
18 typedef int (*QSCOMP)(const void *, const void *);
19 
20 #define uel(a,i)            (((ulong*)(a))[i])
21 #define ucoeff(a,i,j)       (((ulong**)(a))[j][i])
22 #define umael(a,i,j)        (((ulong**)(a))[i][j])
23 #define umael2(a,i,j)       (((ulong**)(a))[i][j])
24 #define umael3(a,i,j,k)     (((ulong***)(a))[i][j][k])
25 #define umael4(a,i,j,k,l)   (((ulong****)(a))[i][j][k][l])
26 #define umael5(a,i,j,k,l,m) (((ulong*****)(a))[i][j][k][l][m])
27 
28 #define numberof(x) (sizeof(x) / sizeof((x)[0]))
29 
30 /* to manipulate 'blocs' */
31 #define BL_HEAD 8
32 #define bl_base(x) (void*)((x) - BL_HEAD)
33 #define bl_height(x) (((GEN)x)[-8])
34 #define bl_left(x)   (((GEN*)x)[-7])
35 #define bl_right(x)  (((GEN*)x)[-6])
36 #define bl_size(x)   (((GEN)x)[-5])
37 #define bl_refc(x)   (((GEN)x)[-4])
38 #define bl_next(x)   (((GEN*)x)[-3])
39 #define bl_prev(x)   (((GEN*)x)[-2])
40 #define bl_num(x)    (((GEN)x)[-1])
41 
42 void clone_lock(GEN C);
43 void clone_unlock(GEN C);
44 void clone_unlock_deep(GEN C);
45 
46 /* swap */
47 #define lswap(x,y) {long _z=x; x=y; y=_z;}
48 #define pswap(x,y) {GEN *_z=x; x=y; y=_z;}
49 #define swap(x,y)  {GEN  _z=x; x=y; y=_z;}
50 #define dswap(x,y) { double _t=x; x=y; y=_t; }
51 #define pdswap(x,y) { double* _t=x; x=y; y=_t; }
52 #define swapspec(x,y, nx,ny) {swap(x,y); lswap(nx,ny);}
53 
54 /* loops */
55 GEN incloop(GEN a);
56 GEN resetloop(GEN a, GEN b);
57 GEN setloop(GEN a);
58 
59 /* parser */
60 
61 /* GP control structures */
62 #define EXPR_WRAP(code, call) \
63 { GEN z; GEN __E = code; \
64   push_lex(gen_0, __E); z = call; pop_lex(1); return z; }
65 #define EXPRVOID_WRAP(code, call) \
66 { GEN __E = code; \
67   push_lex(gen_0, __E); call; pop_lex(1); }
68 #define EXPR_ARG __E, &gp_eval
69 #define EXPR_ARGPREC __E, &gp_evalprec
70 #define EXPR_ARGUPTO __E, &gp_evalupto
71 #define EXPR_ARGBOOL __E, &gp_evalbool
72 #define EXPR_ARGVOID __E, &gp_evalvoid
73 
74 GEN  iferrpari(GEN a, GEN b, GEN c);
75 void forfactored(GEN a, GEN b, GEN code);
76 void forpari(GEN a, GEN b, GEN node);
77 void foreachpari(GEN a, GEN node);
78 void forsquarefree(GEN a, GEN b, GEN code);
79 void untilpari(GEN a, GEN b);
80 void whilepari(GEN a, GEN b);
81 GEN  ifpari(GEN g, GEN a, GEN b);
82 GEN  andpari(GEN a, GEN b);
83 GEN  orpari(GEN a, GEN b);
84 void ifpari_void(GEN g, GEN a, GEN b);
85 GEN  ifpari_multi(GEN g, GEN a);
86 GEN  geval_gp(GEN x, GEN t);
87 
88 GEN  gadde(GEN *x, GEN y);
89 GEN  gadd1e(GEN *x);
90 GEN  gdive(GEN *x, GEN y);
91 GEN  gdivente(GEN *x, GEN y);
92 GEN  gdivrounde(GEN *x, GEN y);
93 GEN  gmode(GEN *x, GEN y);
94 GEN  gmule(GEN *x, GEN y);
95 GEN  gshiftle(GEN *x, long n);
96 GEN  gshiftre(GEN *x, long n);
97 GEN  gsube(GEN *x, GEN y);
98 GEN  gsub1e(GEN *x);
99 GEN  gshift_right(GEN x, long n);
100 
101 GEN  asympnum0(GEN u, GEN alpha, long prec);
102 GEN  asympnumraw0(GEN u, long LIM, GEN alpha, long prec);
103 GEN  derivnum0(GEN a, GEN code, GEN ind, long prec);
104 GEN  derivfun0(GEN args, GEN def, GEN code, long k, long prec);
105 GEN  direuler0(GEN a, GEN b, GEN code, GEN c);
106 GEN  direuler_bad(void *E, GEN (*eval)(void *, GEN, long), GEN a, GEN b, GEN c, GEN Sbad);
107 void forcomposite(GEN a, GEN b, GEN code);
108 void fordiv(GEN a, GEN code);
109 void fordivfactored(GEN a, GEN code);
110 void forell0(long a, long b, GEN code, long flag);
111 void forperm0(GEN k, GEN code);
112 void forprime(GEN a, GEN b, GEN code);
113 void forprimestep(GEN a, GEN b, GEN q, GEN code);
114 void forstep(GEN a, GEN b, GEN s, GEN code);
115 void forsubgroup0(GEN cyc, GEN bound, GEN code);
116 void forsubset0(GEN nk, GEN code);
117 void forvec(GEN x, GEN code, long flag);
118 void forpart0(GEN k, GEN code , GEN nbound, GEN abound);
119 GEN  intcirc0(GEN a, GEN R, GEN code, GEN tab, long prec);
120 GEN  intfuncinit0(GEN a, GEN b, GEN code, long m, long prec);
121 GEN  intnum0(GEN a, GEN b, GEN code, GEN tab, long prec);
122 GEN  intnumgauss0(GEN a, GEN b, GEN code, GEN tab, long prec);
123 GEN  intnumromb0_bitprec(GEN a, GEN b, GEN code, long flag, long bit);
124 GEN  laurentseries0(GEN f, long M, long v, long prec);
125 GEN  limitnum0(GEN u, GEN alpha, long prec);
126 GEN  matrice(GEN nlig, GEN ncol, GEN code);
127 void pariplot0(GEN a, GEN b, GEN code, GEN ysmlu, GEN ybigu, long prec);
128 GEN  prodeuler0(GEN a, GEN b, GEN code, long prec);
129 GEN  prodinf0(GEN a, GEN code, long flag, long prec);
130 GEN  produit(GEN a, GEN b, GEN code, GEN x);
131 GEN  somme(GEN a, GEN b, GEN code, GEN x);
132 GEN  sumalt0(GEN a, GEN code,long flag, long prec);
133 GEN  sumdivexpr(GEN num, GEN code);
134 GEN  sumdivmultexpr0(GEN num, GEN code);
135 GEN  suminf0_bitprec(GEN a, GEN code, long bit);
136 GEN  sumnum0(GEN a, GEN code, GEN tab, long prec);
137 GEN  sumnumap0(GEN a, GEN code, GEN tab, long prec);
138 GEN  sumnumlagrange0(GEN a, GEN code, GEN tab, long prec);
139 GEN  sumnummonien0(GEN a, GEN code, GEN tab, long prec);
140 GEN  sumpos0(GEN a, GEN code, long flag,long prec);
141 GEN  vecexpr0(GEN nmax, GEN code, GEN pred);
142 GEN  vecexpr1(GEN nmax, GEN code, GEN pred);
143 GEN  vecteursmall(GEN nmax, GEN code);
144 GEN  vecteur(GEN nmax, GEN n);
145 GEN  vvecteur(GEN nmax, GEN n);
146 GEN  zbrent0(GEN a, GEN b, GEN code, long prec);
147 GEN  solvestep0(GEN a, GEN b, GEN step, GEN code, long flag, long prec);
148 
149 GEN  ploth0(GEN a, GEN b, GEN code, long flag, long n, long prec);
150 GEN  plothexport0(GEN fmt, GEN a, GEN b, GEN code, long flags, long n, long prec);
151 GEN  psploth0(GEN a,GEN b,GEN code,long flag,long n,long prec);
152 GEN  plotrecth0(long ne,GEN a,GEN b,GEN code,ulong flags,long n,long prec);
153 
154 GEN  listcreate_gp(long n);
155 
156 /* mt */
157 void mt_sigint(void);
158 void mt_err_recover(long er);
159 void mt_export_add(const char *str, GEN val);
160 void mt_export_del(const char *str);
161 void mt_init_stack(size_t s);
162 int  mt_is_thread(void);
163 
164 GEN  eisker_worker(GEN Ei, GEN M, GEN D, GEN co, GEN CD);
165 GEN  pareval_worker(GEN code);
166 GEN  parselect_worker(GEN d, GEN code);
167 void parfor0(GEN a, GEN b, GEN code, GEN code2);
168 GEN  parfor_worker(GEN i, GEN C);
169 void parforeach0(GEN x, GEN code, GEN code2);
170 void parforprime0(GEN a, GEN b, GEN code, GEN code2);
171 void parforprimestep0(GEN a, GEN b, GEN q, GEN code, GEN code2);
172 void parforvec0(GEN a, GEN code, GEN code2, long flag);
173 GEN  parvector_worker(GEN i, GEN C);
174 GEN  polmodular_worker(GEN pt, ulong L, GEN hilb, GEN factu,
175        GEN vne, GEN vinfo, long compute_derivs, GEN j_powers, GEN fdb);
176 GEN  nmV_polint_center_tree_worker(GEN Va, GEN T, GEN R, GEN xa, GEN m2);
177 GEN  nmV_chinese_center_tree_seq(GEN A, GEN P, GEN T, GEN R);
178 GEN  nxMV_polint_center_tree_worker(GEN Va, GEN T, GEN R, GEN xa, GEN m2);
179 GEN  nxMV_chinese_center_tree_seq(GEN A, GEN P, GEN T, GEN R);
180 GEN  F2xq_log_Coppersmith_worker(GEN u, long i, GEN V, GEN R);
181 GEN  Flxq_log_Coppersmith_worker(GEN u, long i, GEN V, GEN R);
182 GEN  Fp_log_sieve_worker(long a, long prmax, GEN C, GEN c, GEN Ci, GEN ci, GEN pr, GEN sz);
183 GEN  QM_charpoly_ZX_worker(GEN P, GEN M, GEN dM);
184 GEN  QXQ_div_worker(GEN P, GEN A, GEN B, GEN C);
185 GEN  QXQ_inv_worker(GEN P, GEN A, GEN B);
186 GEN  ZX_resultant_worker(GEN P, GEN A, GEN B, GEN dB);
187 GEN  ZXQX_resultant_worker(GEN P, GEN A, GEN B, GEN T, GEN dB);
188 GEN  ZX_ZXY_resultant_worker(GEN P, GEN A, GEN B, GEN dB, GEN v);
189 GEN  ZX_direct_compositum_worker(GEN P, GEN A, GEN B);
190 GEN  ZXQX_direct_compositum_worker(GEN P, GEN A, GEN B, GEN C);
191 GEN  ZX_gcd_worker(GEN P, GEN A, GEN B, GEN g);
192 GEN  ZXQ_minpoly_worker(GEN P, GEN A, GEN B, long d);
193 GEN  ZM_det_worker(GEN P, GEN A);
194 GEN  ZM_inv_worker(GEN P, GEN A);
195 GEN  ZM_ker_worker(GEN P, GEN A);
196 GEN  ZM_mul_worker(GEN P, GEN A, GEN B);
197 GEN  ZabM_inv_worker(GEN P, GEN A, GEN Q);
198 GEN  aprcl_step4_worker(ulong q, GEN pC, GEN N, GEN v);
199 GEN  aprcl_step6_worker(GEN r, long t, GEN N, GEN N1, GEN et);
200 GEN  ecpp_sqrt_worker(GEN g, GEN N, GEN p);
201 GEN  ecpp_ispsp_worker(GEN N);
202 GEN  ecpp_step2_worker(GEN S, GEN HD, GEN primelist);
203 GEN  primecertisvalid_ecpp_worker(GEN certi);
204 GEN  lfuninit_worker(long r, GEN K, GEN L, GEN peh2d, GEN vroots, GEN dr, GEN di, GEN an, GEN bn);
205 GEN  lfuninit_theta2_worker(long r, GEN L, GEN qk, GEN a, GEN di, GEN an, GEN bn);
206 GEN  gen_parapply(GEN worker, GEN D);
207 GEN  parapply_slice_worker(GEN worker, GEN D);
208 GEN  gen_parapply_slice(GEN worker, GEN D, long mmin);
209 GEN  gen_crt(const char *str, GEN worker, forprime_t *S, GEN dB, ulong bound, long mmin, GEN *pt_mod,
210              GEN crt(GEN, GEN, GEN*), GEN center(GEN, GEN, GEN));
211 void gen_inccrt(const char *str, GEN worker, GEN dB, long n, long mmin,
212            forprime_t *S, GEN *pt_H, GEN *pt_mod, GEN crt(GEN, GEN, GEN*),
213            GEN center(GEN, GEN, GEN));
214 void gen_inccrt_i(const char *str, GEN worker, GEN dB, long n, long mmin,
215            forprime_t *S, GEN *pH, GEN *pmod, GEN crt(GEN, GEN, GEN*),
216            GEN center(GEN, GEN, GEN));
217 GEN  direllnf_worker(GEN P, ulong X, GEN E);
218 GEN  dirartin_worker(GEN P, ulong X, GEN nf, GEN G, GEN V, GEN aut);
219 GEN  direllsympow_worker(GEN P, ulong X, GEN E, ulong m);
220 GEN  dirgenus2_worker(GEN P, ulong X, GEN Q);
221 GEN  pardireuler(GEN worker, GEN a, GEN b, GEN c, GEN Sbad);
222 GEN  FpM_ratlift_worker(GEN A, GEN mod, GEN B);
223 GEN  ellQ_factorback_worker(GEN A, GEN P, GEN L, GEN c4);
224 GEN  chinese_unit_worker(GEN P, GEN A, GEN U, GEN B, GEN D, GEN C);
225 
226 /* Relative number fields */
227 enum { rnf_NFABS = 1, rnf_MAPS };
228 
229 /* Finite fields */
230 enum { t_FF_FpXQ = 0, t_FF_Flxq = 1, t_FF_F2xq = 2 };
231 GEN FF_ellinit(GEN E, GEN fg);
232 GEN FF_elldata(GEN E, GEN fg);
233 
234 /* L functions */
235 enum { t_LFUN_GENERIC, t_LFUN_ZETA, t_LFUN_NF, t_LFUN_ELL, t_LFUN_KRONECKER,
236        t_LFUN_CHIZ, t_LFUN_CHIGEN, t_LFUN_ETA,
237        t_LFUN_DIV, t_LFUN_MUL, t_LFUN_CONJ,
238        t_LFUN_SYMPOW_ELL, t_LFUN_QF, t_LFUN_ARTIN, t_LFUN_MFCLOS,
239        t_LFUN_GENUS2, t_LFUN_TWIST, t_LFUN_CLOSURE0, t_LFUN_SHIFT};
240 enum { t_LDESC_INIT, t_LDESC_THETA, t_LDESC_PRODUCT };
241 
242 /* Elliptic curves */
243 /* common to Q and Rg */
244 enum { R_PERIODS = 1, R_ETA, R_ROOTS, R_AB };
245 
246 enum { Qp_ROOT = 1, Qp_TATE };
247 enum { Q_GROUPGEN = 5, Q_GLOBALRED, Q_ROOTNO, Q_MINIMALMODEL };
248 enum { NF_MINIMALMODEL = 1, NF_GLOBALRED, NF_MINIMALPRIMES, NF_ROOTNO, NF_NF };
249 
250 /* common to Fp and Fq */
251 enum { FF_CARD = 1, FF_GROUP, FF_GROUPGEN, FF_O };
252 
253 /* for Buchall_param */
254 enum { fupb_NONE = 0, fupb_RELAT, fupb_LARGE, fupb_PRECI };
255 
256 /* Polycyclic presentation for the classgroup of discriminant D */
257 typedef struct {
258   long D; /* Negative discriminant */
259   long h; /* Size of classgroup */
260   long enum_cnt; /* Either h or h/2 (if L0 is set) */
261   /* If nonzero, L0=L[0] and n[0]=2 and classpoly is a perfect square
262    * (and we enumerate each double root just once), default is 0 */
263   long L0;
264   /* Product of primes L that are prohibited as norms of generators or
265    * auxilliary prime forms (by default, primes that make enumeration hard) */
266   long Lfilter;
267   /* Norms of implicit generators (primeforms a=(L*x^2+b*x*y+c*y^2)
268    * with norm L and b >=0) */
269   long *L;
270   long *m; /* products of relative orders: m[i] is the order of <g_1,...,g_i> */
271   long *n; /* Relative orders */
272   long *o; /* Absolute orders */
273   /* Power relations (a[i]^n[i] = a[0]^e[0]*...*a[i-1]^e[i-1], where e
274    * is an exponent vector of length i stored at offset binom(i,2) of r) */
275   long *r;
276   long *orient_p; /* Optional list of norms of orienting primes p ... */
277   long *orient_q; /* or product of primes p*q (q=1 when only p is needed) */
278   long *orient_reps; /* Representation of orienting norm p*q in terms of Ls */
279   long inv; /* Attached invariant */
280   long k; /* Number of generators */
281   GEN _data; /* Storage space for the above arrays */
282 } classgp_pcp_struct;
283 typedef classgp_pcp_struct classgp_pcp_t[1];
284 
285 /* Represents the data in the equation(s)
286  *   4p = t^2 - v^2 D = t^2 - v^2 u^2 D_K = w^2 D_K.
287  * t is the absolute trace, so always > 0.
288  * T is a twisting parameter, which satisfies (T|p) == -1. */
289 typedef struct {
290   long D, t, u, v;
291   ulong p, pi, s2, T;
292 } norm_eqn_struct;
293 typedef norm_eqn_struct norm_eqn_t[1];
294 
295 #define zv_to_longptr(v) (&((v)[1]))
296 #define zv_to_ulongptr(v) ((ulong *)&((v)[1]))
297 
298 /* Modular invariants */
299 #define INV_J       0
300 #define INV_F       1
301 #define INV_F2      2
302 #define INV_F3      3
303 #define INV_F4      4
304 #define INV_G2      5
305 #define INV_W2W3    6
306 #define INV_F8      8
307 #define INV_W3W3    9
308 #define INV_W2W5    10
309 #define INV_W2W7    14
310 #define INV_W3W5    15
311 #define INV_W3W7    21
312 #define INV_W2W3E2  23
313 #define INV_W2W5E2  24
314 #define INV_W2W13   26
315 #define INV_W2W7E2  27
316 #define INV_W3W3E2  28
317 #define INV_W5W7    35
318 #define INV_W3W13   39
319 
320 /* Get coefficient of x^d in f, assuming f is nonzero. */
Flx_coeff(GEN f,long d)321 INLINE ulong Flx_coeff(GEN f, long d) { return f[d + 2]; }
322 /* Return the root of f, assuming deg(f) = 1. */
Flx_deg1_root(GEN f,ulong p)323 INLINE ulong Flx_deg1_root(GEN f, ulong p) {
324   if (degpol(f) != 1) pari_err_BUG("Flx_deg1_root");
325   return Fl_div(Fl_neg(Flx_coeff(f, 0), p), Flx_coeff(f, 1), p);
326 }
327 
328 /* Allocation / gerepile */
329 long   getdebugvar(void);
330 void   setdebugvar(long n);
331 void   debug_stack(void);
332 void   fill_stack(void);
333 void   minim_alloc(long n, double ***q, GEN *x, double **y,  double **z, double **v);
334 int    pop_entree_block(entree *ep, long loc);
335 int    pop_val_if_newer(entree *ep, long loc);
336 
337 /* general printing */
338 void print_errcontext(PariOUT *out, const char *msg, const char *s, const char *entry);
339 void print_prefixed_text(PariOUT *out, const char *s, const char *prefix, const char *str);
340 INLINE void
print_text(const char * s)341 print_text(const char *s) { print_prefixed_text(pariOut, s,NULL,NULL); }
342 INLINE void
out_print_text(PariOUT * out,const char * s)343 out_print_text(PariOUT *out, const char *s) { print_prefixed_text(out, s,NULL,NULL); }
344 INLINE long
is_keyword_char(char c)345 is_keyword_char(char c) { return (isalnum((int)c) || c=='_'); }
346 
347 /* Interfaces (GP, etc.) */
348 hashtable *hash_from_link(GEN e, GEN names, int use_stack);
349 void gen_relink(GEN x, hashtable *table);
350 entree* do_alias(entree *ep);
351 char* get_sep(const char *t);
352 long get_int(const char *s, long dflt);
353 ulong get_uint(const char *s);
354 void gp_initrc(pari_stack *p_A);
355 
356 void pari_sigint(const char *s);
357 void* get_stack(double fraction, long min);
358 void  free_graph(void);
359 void  initout(int initerr);
360 void  resetout(int initerr);
361 void  init_linewrap(long w);
362 void  print_functions_hash(const char *s);
363 GEN   readbin(const char *name, FILE *f, int *vector);
364 int   term_height(void);
365 int   term_width(void);
366 /* gp_colors */
367 void decode_color(long n, long *c);
368 
369 /* defaults */
370 extern long precreal;
371 
372 void lim_lines_output(char *s, long n, long max);
373 int tex2mail_output(GEN z, long n);
374 void gen_output(GEN x);
375 void fputGEN_pariout(GEN x, pariout_t *T, FILE *out);
376 
377 void parsestate_reset(void);
378 void parsestate_save(struct pari_parsestate *state);
379 void parsestate_restore(struct pari_parsestate *state);
380 
381 void compilestate_reset(void);
382 void compilestate_save(struct pari_compilestate *comp);
383 void compilestate_restore(struct pari_compilestate *comp);
384 
385 void filestate_save(struct pari_filestate *file);
386 void filestate_restore(struct pari_filestate *file);
387 void tmp_restore(pariFILE *F);
388 
389 long evalstate_get_trace(void);
390 void evalstate_set_trace(long lvl);
391 void evalstate_clone(void);
392 void evalstate_reset(void);
393 void evalstate_restore(struct pari_evalstate *state);
394 GEN  evalstate_restore_err(struct pari_evalstate *state);
395 void evalstate_save(struct pari_evalstate *state);
396 void varstate_save(struct pari_varstate *s);
397 void varstate_restore(struct pari_varstate *s);
398 
399 void mtstate_save(struct pari_mtstate *s);
400 void mtstate_reset(void);
401 void mtstate_restore(struct pari_mtstate *s);
402 
403 void debug_context(void);
404 
405 typedef struct {
406   const char *s;
407   size_t ls;
408   char **dir;
409 } forpath_t;
410 void forpath_init(forpath_t *T, gp_path *path, const char *s);
411 char *forpath_next(forpath_t *T);
412 
413 /* GP output && output format */
414 void gpwritebin(const char *s, GEN x);
415 extern char *current_logfile;
416 
417 /* colors */
418 extern long    gp_colors[];
419 extern int     disable_color;
420 
421 /* entrees */
422 #define EpVALENCE(ep) ((ep)->valence & 0xFF)
423 #define EpSTATIC(ep) ((ep)->valence & 0x100)
424 #define EpSETSTATIC(ep) ((ep)->valence |= 0x100)
425 enum { EpNEW = 100, EpALIAS, EpVAR, EpINSTALL };
426 #define initial_value(ep) ((ep)+1)
427 
428 /* functions lists */
429 extern const long functions_tblsz;  /* hashcodes table size */
430 extern entree **functions_hash;   /* functions hashtable */
431 extern entree **defaults_hash;    /* defaults hashtable */
432 
433 /* buffers */
434 typedef struct Buffer {
435   char *buf;
436   ulong len;
437   jmp_buf env;
438 } Buffer;
439 Buffer *new_buffer(void);
440 void delete_buffer(Buffer *b);
441 void fix_buffer(Buffer *b, long newlbuf);
442 
443 typedef struct {
444   const char *s; /* source */
445   char *t, *end; /* target, last char read */
446   int in_string, in_comment, more_input, wait_for_brace;
447   Buffer *buf;
448 } filtre_t;
449 void init_filtre(filtre_t *F, Buffer *buf);
450 Buffer *filtered_buffer(filtre_t *F);
451 void kill_buffers_upto_including(Buffer *B);
452 void pop_buffer(void);
453 void kill_buffers_upto(Buffer *B);
454 int gp_read_line(filtre_t *F, const char *PROMPT);
455 void parse_key_val(char *src, char **ps, char **pt);
456 extern int (*cb_pari_get_line_interactive)(const char*, const char*, filtre_t *F);
457 extern char *(*cb_pari_fgets_interactive)(char *s, int n, FILE *f);
458 int get_line_from_file(const char *prompt, filtre_t *F, FILE *file);
459 void pari_skip_space(char **s);
460 void pari_skip_alpha(char **s);
461 char *pari_translate_string(const char *src, char *s, char *entry);
462 
463 gp_data *default_gp_data(void);
464 
465 typedef char *(*fgets_t)(char *, int, void*);
466 
467 typedef struct input_method {
468 /* optional */
469   fgets_t myfgets;  /* like libc fgets() but last argument is (void*) */
470 /* mandatory */
471   char * (*getline)(char**, int f, struct input_method*, filtre_t *F);
472   int free; /* boolean: must we free the output of getline() ? */
473 /* optional */
474   const char *prompt, *prompt_cont;
475   void *file;  /* can be used as last argument for fgets() */
476 } input_method;
477 
478 int input_loop(filtre_t *F, input_method *IM);
479 char *file_input(char **s0, int junk, input_method *IM, filtre_t *F);
480 char *file_getline(Buffer *b, char **s0, input_method *IM);
481 
482 /* readline */
483 typedef struct {
484   /* pointers to readline variables/functions */
485   char **line_buffer;
486   int *point;
487   int *end;
488   char **(*completion_matches)(const char *, char *(*)(const char*, int));
489   char *(*filename_completion_function)(const char *, int);
490   char *(*username_completion_function)(const char *, int);
491   int (*insert)(int, int);
492   int *completion_append_character;
493 
494   /* PARI-specific */
495   int back;  /* rewind the cursor by this number of chars */
496 } pari_rl_interface;
497 
498 /* Code which wants to use readline needs to do the following:
499 
500 #include <readline.h>
501 #include <paripriv.h>
502 pari_rl_interface pari_rl;
503 pari_use_readline(pari_rl);
504 
505 This will initialize the pari_rl structure. A pointer to this structure
506 must be given as first argument to all PARI readline functions. */
507 
508 /* IMPLEMENTATION NOTE: this really must be a macro (not a function),
509  * since we refer to readline symbols. */
510 #define pari_use_readline(pari_rl) do {\
511     (pari_rl).line_buffer = &rl_line_buffer; \
512     (pari_rl).point = &rl_point; \
513     (pari_rl).end = &rl_end; \
514     (pari_rl).completion_matches = &rl_completion_matches; \
515     (pari_rl).filename_completion_function = &rl_filename_completion_function; \
516     (pari_rl).username_completion_function = &rl_username_completion_function; \
517     (pari_rl).insert = &rl_insert; \
518     (pari_rl).completion_append_character = &rl_completion_append_character; \
519     (pari_rl).back = 0; } while(0)
520 
521 /* FIXME: EXPORT AND DOCUMENT THE FOLLOWING */
522 
523 /* PROBABLY NOT IN THE RIGHT FILE, SORT BY THEME */
524 
525 /* multiprecision */
526 GEN   adduispec_offset(ulong s, GEN x, long offset, long nx);
527 int   lgcdii(ulong* d, ulong* d1, ulong* u, ulong* u1, ulong* v, ulong* v1, ulong vmax);
528 ulong rgcduu(ulong d, ulong d1, ulong vmax, ulong* u, ulong* u1, ulong* v, ulong* v1, long *s);
529 ulong xgcduu(ulong d, ulong d1, int f, ulong* v, ulong* v1, long *s);
530 ulong xxgcduu(ulong d, ulong d1, int f, ulong* u, ulong* u1, ulong* v, ulong* v1, long *s);
531 GEN   divgunu(GEN x, ulong i);
532 GEN   divrunu(GEN x, ulong i);
533 GEN   muliispec(GEN x, GEN y, long nx, long ny);
534 GEN   red_montgomery(GEN T, GEN N, ulong inv);
535 GEN   sqrispec(GEN x, long nx);
536 ulong *convi(GEN x, long *l);
537 
538 int approx_0(GEN x, GEN y);
539 
540 /* powers */
541 GEN    rpowuu(ulong a, ulong n, long prec);
542 
543 /* floats */
544 double dabs(double s, double t);
545 double darg(double s, double t);
546 void   dcxlog(double s, double t, double *a, double *b);
547 double dnorm(double s, double t);
548 double dbllog2(GEN z);
549 
550 /* hnf */
551 GEN hnfadd(GEN m,GEN p,GEN* ptdep,GEN* ptA,GEN* ptC,GEN extramat,GEN extraC);
552 GEN hnfadd_i(GEN m,GEN p,GEN* ptdep,GEN* ptA,GEN* ptC,GEN extramat,GEN extraC);
553 GEN hnfspec_i(GEN m,GEN p,GEN* ptdep,GEN* ptA,GEN* ptC,long k0);
554 GEN hnfspec(GEN m,GEN p,GEN* ptdep,GEN* ptA,GEN* ptC,long k0);
555 GEN mathnfspec(GEN x, GEN *ptperm, GEN *ptdep, GEN *ptB, GEN *ptC);
556 GEN ZM_hnfmodall_i(GEN x, GEN dm, long flag);
557 
558 GEN LLL_check_progress(GEN Bnorm, long n0, GEN m, int final, long *ti_LLL);
559 
560 /* integer factorization / discrete log */
561 ulong is_kth_power(GEN x, ulong p, GEN *pt);
562 GEN   mpqs(GEN N);
563 
564 /* Polynomials */
565 /* a) Arithmetic/conversions */
566 GEN  lift_if_rational(GEN x);
567 GEN  monomial(GEN a, long degpol, long v);
568 GEN  monomialcopy(GEN a, long degpol, long v);
569 GEN  ser2pol_i(GEN x, long lx);
570 GEN  ser2rfrac_i(GEN x);
571 GEN  swap_vars(GEN b0, long v);
572 GEN  RgX_recipspec_shallow(GEN x, long l, long n);
573 
574 /* b) Modular */
575 GEN  bezout_lift_fact(GEN T, GEN Tmod, GEN p, long e);
576 GEN  polsym_gen(GEN P, GEN y0, long n, GEN T, GEN N);
577 GEN  ZXQ_charpoly_sqf(GEN A, GEN B, long *lambda, long v);
578 GEN  ZX_disc_all(GEN,ulong);
579 GEN  ZX_resultant_all(GEN A, GEN B, GEN dB, ulong bound);
580 GEN  ZX_ZXY_resultant_all(GEN A, GEN B, long *lambda, GEN *LPRS);
581 
582 GEN FlxqM_mul_Kronecker(GEN A, GEN B, GEN T, ulong p);
583 GEN FqM_mul_Kronecker(GEN x, GEN y, GEN T, GEN p);
584 
585 /* c) factorization */
586 GEN chk_factors_get(GEN lt, GEN famod, GEN c, GEN T, GEN N);
587 long cmbf_maxK(long nb);
588 GEN ZX_DDF(GEN x);
589 GEN initgaloisborne(GEN T, GEN dn, long prec, GEN *pL, GEN *pprep, GEN *pdis);
590 GEN quicktrace(GEN x, GEN sym);
591 
592 /* pari_init / pari_close */
593 void pari_close_compiler(void);
594 void pari_close_evaluator(void);
595 void pari_close_files(void);
596 void pari_close_floats(void);
597 void pari_close_homedir(void);
598 void pari_close_parser(void);
599 void pari_close_paths(void);
600 void pari_close_primes(void);
601 void pari_init_buffers(void);
602 void pari_init_compiler(void);
603 void pari_init_defaults(void);
604 void pari_init_evaluator(void);
605 void pari_init_files(void);
606 void pari_init_floats(void);
607 void pari_init_homedir(void);
608 void pari_init_graphics(void);
609 void pari_init_parser(void);
610 void pari_init_rand(void);
611 void pari_init_paths(void);
612 void pari_init_primetab(void);
613 void pari_init_seadata(void);
614 GEN pari_get_seadata(void);
615 void pari_set_primetab(GEN global_primetab);
616 void pari_set_seadata(GEN seadata);
617 void pari_set_varstate(long *vp, struct pari_varstate *vs);
618 void pari_thread_close_files(void);
619 
620 void export_add(const char *str, GEN val);
621 void export_del(const char *str);
622 GEN  export_get(const char *str);
623 void exportall(void);
624 void unexportall(void);
625 
626 /* BY FILES */
627 
628 /* parinf.h */
629 
630 GEN fincke_pohst(GEN a,GEN BOUND,long stockmax,long PREC, FP_chk_fun *CHECK);
631 void init_zlog(zlog_S *S, GEN bid);
632 GEN  log_gen_arch(zlog_S *S, long index);
633 GEN  log_gen_pr(zlog_S *S, long index, GEN nf, long e);
634 GEN  sprk_log_gen_pr(GEN nf, GEN sprk, long e);
635 GEN  sprk_log_prk1(GEN nf, GEN a, GEN sprk);
636 GEN    poltobasis(GEN nf,GEN x);
637 GEN    coltoalg(GEN nf,GEN x);
638 
639 GEN    rnfdisc_get_T(GEN nf, GEN P, GEN *lim);
640 GEN    make_integral(GEN nf, GEN L0, GEN f, GEN listpr);
641 GEN    rnfallbase(GEN nf, GEN pol, GEN lim, GEN rnfeq, GEN *pD, GEN *pfi,
642                   GEN *pdKP);
643 GEN    subgroupcondlist(GEN cyc, GEN bound, GEN listKer);
644 
645 /* Qfb.c */
646 
647 GEN     redimagsl2(GEN q, GEN *U);
648 GEN     redrealsl2(GEN V, GEN d, GEN rd);
649 GEN     redrealsl2step(GEN A, GEN d, GEN rd);
650 
651 /* alglin1.c */
652 
653 typedef long (*pivot_fun)(GEN,GEN,long,GEN);
654 GEN ZM_pivots(GEN x0, long *rr);
655 GEN RgM_pivots(GEN x0, GEN data, long *rr, pivot_fun pivot);
656 void RgMs_structelim_col(GEN M, long nbcol, long nbrow, GEN A, GEN *p_col, GEN *p_lin);
657 
658 /* arith1.c */
659 
660 int     is_gener_Fp(GEN x, GEN p, GEN p_1, GEN L);
661 int     is_gener_Fl(ulong x, ulong p, ulong p_1, GEN L);
662 
663 /* arith2.c */
664 
665 int     divisors_init(GEN n, GEN *pP, GEN *pE);
666 long    set_optimize(long what, GEN g);
667 
668 /* base1.c */
669 
670 GEN zk_galoisapplymod(GEN nf, GEN z, GEN S, GEN p);
671 
672 /* base2.c */
673 
674 GEN     dim1proj(GEN prh);
675 GEN     gen_if_principal(GEN bnf, GEN x);
676 
677 /* base3.c */
678 
679 void    check_nfelt(GEN x, GEN *den);
680 GEN     zk_ei_mul(GEN nf, GEN x, long i);
681 GEN     log_prk(GEN nf, GEN a, GEN sprk, GEN mod);
682 GEN     log_prk_units(GEN nf, GEN D, GEN sprk);
683 GEN     log_prk_units_init(GEN bnf);
684 GEN     veclog_prk(GEN nf, GEN v, GEN sprk);
685 GEN     log_prk_init(GEN nf, GEN pr, long k, GEN mod);
686 
687 /* base4.c */
688 
689 GEN     factorbackprime(GEN nf, GEN L, GEN e);
690 
691 /* bb_group.c */
692 
693 GEN     producttree_scheme(long n);
694 
695 /* bern.c */
696 long bernbitprec(long N);
697 
698 /* bibli2.c */
699 
700 GEN sort_factor_pol(GEN y, int (*cmp)(GEN,GEN));
701 
702 /* buch1.c */
703 
704 long   bnf_increase_LIMC(long LIMC, long LIMCMAX);
705 
706 /* buch2.c */
707 
708 typedef struct GRHprime_t { ulong p; double logp; GEN dec; } GRHprime_t;
709 typedef struct GRHcheck_t { double cD, cN; GRHprime_t *primes; long clone, nprimes, maxprimes; ulong limp; forprime_t P; } GRHcheck_t;
710 void    free_GRHcheck(GRHcheck_t *S);
711 void    init_GRHcheck(GRHcheck_t *S, long N, long R1, double LOGD);
712 void    GRH_ensure(GRHcheck_t *S, long nb);
713 ulong   GRH_last_prime(GRHcheck_t *S);
714 int     GRHok(GRHcheck_t *S, double L, double SA, double SB);
715 GEN     extract_full_lattice(GEN x);
716 GEN     init_red_mod_units(GEN bnf, long prec);
717 GEN     isprincipalarch(GEN bnf, GEN col, GEN kNx, GEN e, GEN dx, long *pe);
718 GEN     red_mod_units(GEN col, GEN z);
719 
720 /* buch3.c */
721 
722 GEN     minkowski_bound(GEN D, long N, long r2, long prec);
723 int     subgroup_conductor_ok(GEN H, GEN L);
724 GEN     subgrouplist_cond_sub(GEN bnr, GEN C, GEN bound);
725 
726 /* buch4.c */
727 
728 GEN     nf_quadchar_modpr(GEN nf, GEN z, GEN modpr, GEN pstar);
729 
730 /* crvwtors.c */
731 
732 void random_curves_with_m_torsion(ulong *a4, ulong *a6, ulong *tx, ulong *ty, long ncurves, long m, ulong p);
733 
734 /* dirichlet.c */
735 GEN direuler_factor(GEN s, long n);
736 
737 /* elliptic.c */
738 
739 void ellprint(GEN e);
740 GEN  elltors_psylow(GEN e, ulong p);
741 GEN  ellintegralbmodel(GEN e, GEN *pv);
742 GEN  ellQ_genreduce(GEN E, GEN G, long prec);
743 
744 /* es.c */
745 
746 void    killallfiles(void);
747 pariFILE* newfile(FILE *f, const char *name, int type);
748 int     popinfile(void);
749 pariFILE* try_pipe(const char *cmd, int flag);
750 
751 /* F2m.c */
752 
753 GEN     F2m_gauss_pivot(GEN x, long *rr);
754 GEN     F2m_gauss_sp(GEN a, GEN b);
755 GEN     F2m_invimage_i(GEN A, GEN B);
756 
757 /* Fle.c */
758 
759 void    FleV_add_pre_inplace(GEN P, GEN Q, GEN a4, ulong p, ulong pi);
760 void    FleV_dbl_pre_inplace(GEN P, GEN a4, ulong p, ulong pi);
761 void    FleV_mulu_pre_inplace(GEN P, ulong n, GEN a4, ulong p, ulong pi);
762 void    FleV_sub_pre_inplace(GEN P, GEN Q, GEN a4, ulong p, ulong pi);
763 
764 /* Flv.c */
765 
766 GEN     Flm_gauss_sp(GEN a, GEN b, ulong *detp, ulong p);
767 GEN     Flm_invimage_i(GEN A, GEN B, ulong p);
768 GEN     Flm_inv_sp(GEN a, ulong *detp, ulong p);
769 GEN     Flm_pivots(GEN x, ulong p, long *rr, long inplace);
770 
771 /* Flxq_log.c */
772 
773 GEN Flxq_log_index(GEN a0, GEN b0, GEN m, GEN T0, ulong p);
774 int Flxq_log_use_index(GEN m, GEN T0, ulong p);
775 
776 /* FlxqE.c */
777 
778 GEN     ZpXQ_norm_pcyc(GEN x, GEN T, GEN q, GEN p);
779 long    zx_is_pcyc(GEN T);
780 
781 /* FpV.c */
782 
783 GEN FpMs_leftkernel_elt_col(GEN M, long nbcol, long nbrow, GEN p);
784 GEN FpX_to_mod_raw(GEN z, GEN p);
785 
786 /* FpX.c */
787 
788 GEN     ZlXQXn_expint(GEN h, long e, GEN T, GEN p, ulong pp);
789 
790 /* FpX_factor.c */
791 
792 GEN     ddf_to_ddf2(GEN V);
793 long    ddf_to_nbfact(GEN D);
794 GEN     vddf_to_simplefact(GEN V, long d);
795 
796 /* FpXQX_factor.c */
797 
798 GEN     FpXQX_factor_Berlekamp(GEN x, GEN T, GEN p);
799 
800 /* forprime.c*/
801 
802 void    init_modular_big(forprime_t *S);
803 void    init_modular_small(forprime_t *S);
804 
805 /* galconj.c */
806 
807 GEN     galoiscosets(GEN O, GEN perm);
808 GEN     matrixnorm(GEN M, long prec);
809 
810 /* gen1.c */
811 
812 GEN     gred_rfrac_simple(GEN n, GEN d);
813 GEN     sqr_ser_part(GEN x, long l1, long l2);
814 
815 /* hash.c */
816 
817 hashtable *hashstr_import_static(hashentry *e, ulong size);
818 
819 /* hyperell.c */
820 
821 GEN     ZlXQX_hyperellpadicfrobenius(GEN H, GEN T, ulong p, long n);
822 GEN     hyperell_redsl2(GEN Q);
823 
824 /* ifactor1.c */
825 
826 ulong snextpr(ulong p, byteptr *d, long *rcn, long *q, int (*ispsp)(ulong));
827 
828 /* intnum.c */
829 
830 GEN     contfraceval_inv(GEN CF, GEN tinv, long nlim);
831 
832 /* mftrace.c */
833 
834 void pari_close_mf(void);
835 long polishomogeneous(GEN P);
836 GEN sertocol(GEN S);
837 
838 /* prime.c */
839 
840 long    BPSW_psp_nosmalldiv(GEN N);
841 int     MR_Jaeschke(GEN n);
842 long    isanypower_nosmalldiv(GEN N, GEN *px);
843 void    prime_table_next_p(ulong a, byteptr *pd, ulong *pp, ulong *pn);
844 
845 /* perm.c */
846 
847 long    cosets_perm_search(GEN C, GEN p);
848 GEN     perm_generate(GEN S, GEN H, long o);
849 long    perm_relorder(GEN p, GEN S);
850 GEN     vecperm_extendschreier(GEN C, GEN v, long n);
851 
852 /* polclass.c */
853 
854 GEN polclass0(long D, long inv, long xvar, GEN *db);
855 
856 /* polmodular.c */
857 
858 GEN polmodular0_ZM(long L, long inv, GEN J, GEN Q, int compute_derivs, GEN *db);
859 GEN Flm_Fl_polmodular_evalx(GEN phi, long L, ulong j, ulong p, ulong pi);
860 GEN polmodular_db_init(long inv);
861 void polmodular_db_clear(GEN db);
862 void polmodular_db_add_level(GEN *db, long L, long inv);
863 void polmodular_db_add_levels(GEN *db, long *levels, long k, long inv);
864 GEN polmodular_db_for_inv(GEN db, long inv);
865 GEN polmodular_db_getp(GEN fdb, long L, ulong p);
866 
867 long modinv_level(long inv);
868 long modinv_degree(long *p1, long *p2, long inv);
869 long modinv_ramified(long D, long inv);
870 long modinv_j_from_2double_eta(GEN F, long inv, ulong x0, ulong x1, ulong p, ulong pi);
871 GEN double_eta_raw(long inv);
872 ulong modfn_root(ulong j, norm_eqn_t ne, long inv);
873 long modfn_unambiguous_root(ulong *r, long inv, ulong j0, norm_eqn_t ne, GEN jdb);
874 GEN qfb_nform(long D, long n);
875 
876 /* Fle.c */
877 
878 ulong   Flj_order_ufact(GEN P, ulong n, GEN F, ulong a4, ulong p, ulong pi);
879 
880 /* polarit3.c */
881 
882 GEN     Flm_Frobenius_pow(GEN M, long d, GEN T, ulong p);
883 GEN     FpM_Frobenius_pow(GEN M, long d, GEN T, GEN p);
884 GEN     FpX_compositum(GEN A, GEN B, GEN p);
885 GEN     Flx_direct_compositum(GEN A, GEN B, ulong p);
886 GEN     FlxV_direct_compositum(GEN V, ulong p);
887 GEN     FlxqX_direct_compositum(GEN P, GEN Q, GEN T, ulong p);
888 GEN     FpX_direct_compositum(GEN A, GEN B, GEN p);
889 GEN     FpXV_direct_compositum(GEN V, GEN p);
890 GEN     nf_direct_compositum(GEN nf, GEN A, GEN B);
891 ulong   ZX_ZXY_ResBound(GEN A, GEN B, GEN dB);
892 GEN     ffinit_Artin_Schreier(ulong p, long l);
893 GEN     ffinit_rand(GEN p, long n);
894 
895 /* readline.c */
896 
897 char**  pari_completion(pari_rl_interface *pari_rl, char *text, int START, int END);
898 char**  pari_completion_matches(pari_rl_interface *pari_rl, const char *s, long pos, long *wordpos);
899 
900 /* RgX.c */
901 
902 GEN     RgX_homogenous_evalpow(GEN P, GEN A, GEN B);
903 GEN     QXQX_homogenous_evalpow(GEN P, GEN A, GEN B, GEN T);
904 
905 /* subcyclo.c */
906 
907 GEN     galoiscyclo(long n, long v);
908 GEN     znstar_bits(long n, GEN H);
909 long    znstar_conductor(GEN H);
910 long    znstar_conductor_bits(GEN bits);
911 GEN     znstar_cosets(long n, long phi_n, GEN H);
912 GEN     znstar_elts(long n, GEN H);
913 GEN     znstar_generate(long n, GEN V);
914 GEN     znstar_hnf(GEN Z, GEN M);
915 GEN     znstar_hnf_elts(GEN Z, GEN H);
916 GEN     znstar_hnf_generators(GEN Z, GEN M);
917 GEN     znstar_reduce_modulus(GEN H, long n);
918 GEN     znstar_small(GEN zn);
919 
920 /* trans1.c */
921 
922 struct abpq { GEN *a, *b, *p, *q; };
923 struct abpq_res { GEN P, Q, B, T; };
924 void    abpq_init(struct abpq *A, long n);
925 void    abpq_sum(struct abpq_res *r, long n1, long n2, struct abpq *A);
926 GEN     logagmcx(GEN q, long prec);
927 GEN     zellagmcx(GEN a0, GEN b0, GEN r, GEN t, long prec);
928 
929 /* trans2.c */
930 
931 GEN     trans_fix_arg(long *prec, GEN *s0, GEN *sig, GEN *tau, pari_sp *av, GEN *res);
932 
933 /* trans3.c */
934 
935 GEN     double_eta_quotient(GEN a, GEN w, GEN D, long p, long q, GEN pq, GEN sqrtD);
936 GEN     inv_szeta_euler(long n, long prec);
937 
938 /* volcano.c */
939 
940 long j_level_in_volcano(GEN phi, ulong j, ulong p, ulong pi, long L, long depth);
941 ulong ascend_volcano(GEN phi, ulong j, ulong p, ulong pi, long level, long L, long depth, long steps);
942 ulong descend_volcano(GEN phi, ulong j, ulong p, ulong pi, long level, long L, long depth, long steps);
943 long next_surface_nbr(ulong *nJ, GEN phi, long L, long h, ulong J, const ulong *pJ, ulong p, ulong pi);
944 GEN enum_roots(ulong j, norm_eqn_t ne, GEN fdb, classgp_pcp_t G);
945 
946 ENDEXTERN
947