Line data Source code
1 : /* Copyright (C) 2000 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 : #include "pari.h"
16 : #include "paripriv.h"
17 :
18 : GEN
19 1278237 : iferrpari(GEN a, GEN b, GEN c)
20 : {
21 : GEN res;
22 : struct pari_evalstate state;
23 1278237 : evalstate_save(&state);
24 1278237 : pari_CATCH(CATCH_ALL)
25 : {
26 : GEN E;
27 72460 : if (!b&&!c) return gnil;
28 36237 : E = evalstate_restore_err(&state);
29 36237 : if (c)
30 : {
31 292 : push_lex(E,c);
32 292 : res = closure_evalnobrk(c);
33 285 : pop_lex(1);
34 285 : if (gequal0(res))
35 7 : pari_err(0, E);
36 : }
37 36223 : if (!b) return gnil;
38 36223 : push_lex(E,b);
39 36223 : res = closure_evalgen(b);
40 36223 : pop_lex(1);
41 36223 : return res;
42 : } pari_TRY {
43 1278237 : res = closure_evalgen(a);
44 1242000 : } pari_ENDCATCH;
45 1242000 : return res;
46 : }
47 :
48 : /********************************************************************/
49 : /** **/
50 : /** ITERATIONS **/
51 : /** **/
52 : /********************************************************************/
53 :
54 : static void
55 5112657 : forparii(GEN a, GEN b, GEN code)
56 : {
57 5112657 : pari_sp av, av0 = avma;
58 : GEN aa;
59 5112657 : if (gcmp(b,a) < 0) return;
60 5030848 : if (typ(b) != t_INFINITY) b = gfloor(b);
61 5030848 : aa = a = setloop(a);
62 5030848 : av=avma;
63 5030848 : push_lex(a,code);
64 103437667 : while (gcmp(a,b) <= 0)
65 : {
66 98534109 : closure_evalvoid(code); if (loop_break()) break;
67 98426645 : a = get_lex(-1);
68 98410721 : if (a == aa)
69 : {
70 98411163 : a = incloop(a);
71 98406791 : if (a != aa) { set_lex(-1,a); aa = a; }
72 : }
73 : else
74 : { /* 'code' modified a ! Be careful (and slow) from now on */
75 24 : a = gaddgs(a,1);
76 28 : if (gc_needed(av,1))
77 : {
78 0 : if (DEBUGMEM>1) pari_warn(warnmem,"forparii");
79 0 : a = gerepileupto(av,a);
80 : }
81 28 : set_lex(-1,a);
82 : }
83 : }
84 5027239 : pop_lex(1); set_avma(av0);
85 : }
86 :
87 : void
88 5112664 : forpari(GEN a, GEN b, GEN code)
89 : {
90 5112664 : pari_sp ltop=avma, av;
91 5112664 : if (typ(a) == t_INT) { forparii(a,b,code); return; }
92 7 : b = gcopy(b); /* Kludge to work-around the a+(a=2) bug */
93 7 : av=avma;
94 7 : push_lex(a,code);
95 28 : while (gcmp(a,b) <= 0)
96 : {
97 21 : closure_evalvoid(code); if (loop_break()) break;
98 21 : a = get_lex(-1); a = gaddgs(a,1);
99 21 : if (gc_needed(av,1))
100 : {
101 0 : if (DEBUGMEM>1) pari_warn(warnmem,"forpari");
102 0 : a = gerepileupto(av,a);
103 : }
104 21 : set_lex(-1, a);
105 : }
106 7 : pop_lex(1); set_avma(ltop);
107 : }
108 :
109 : void
110 980 : foreachpari(GEN x, GEN code)
111 : {
112 : long i, l;
113 980 : switch(typ(x))
114 : {
115 14 : case t_LIST:
116 14 : x = list_data(x); /* FALL THROUGH */
117 14 : if (!x) return;
118 : case t_MAT: case t_VEC: case t_COL:
119 966 : break;
120 7 : default:
121 7 : pari_err_TYPE("foreach",x);
122 : return; /*LCOV_EXCL_LINE*/
123 : }
124 966 : clone_lock(x); l = lg(x);
125 966 : push_lex(gen_0,code);
126 4704 : for (i = 1; i < l; i++)
127 : {
128 3738 : set_lex(-1, gel(x,i));
129 3738 : closure_evalvoid(code); if (loop_break()) break;
130 : }
131 966 : pop_lex(1); clone_unlock_deep(x);
132 : }
133 :
134 : /* is it better to sieve [a,b] or to factor individually ? */
135 : static int
136 189 : no_sieve(ulong a, ulong b)
137 189 : { return b - a < usqrt(b) / tridiv_boundu(b); }
138 :
139 : /* 0 < a <= b. Using small consecutive chunks to 1) limit memory use, 2) allow
140 : * cheap early abort */
141 : static int
142 63 : forfactoredpos(ulong a, ulong b, GEN code)
143 : {
144 63 : ulong x1, step = maxuu(2 * usqrt(b), 1024);
145 63 : pari_sp av = avma;
146 63 : if (no_sieve(a, b))
147 : {
148 : ulong n;
149 0 : for (n = a; n <= b; n++, set_avma(av))
150 : {
151 0 : GEN m = factoru(n);
152 0 : set_lex(-1, mkvec2(utoipos(n), Flm_to_ZM(m)));
153 0 : closure_evalvoid(code); if (loop_break()) return 1;
154 : }
155 0 : return 0;
156 : }
157 3549 : for(x1 = a;; x1 += step, set_avma(av))
158 3486 : { /* beware overflow, fuse last two bins (avoid a tiny remainder) */
159 3549 : ulong j, lv, x2 = (b >= 2*step && b - 2*step >= x1)? x1-1 + step: b;
160 3549 : GEN v = vecfactoru_i(x1, x2);
161 3549 : lv = lg(v);
162 7005082 : for (j = 1; j < lv; j++)
163 : {
164 7001547 : ulong n = x1-1 + j;
165 7001547 : set_lex(-1, mkvec2(utoipos(n), Flm_to_ZM(gel(v,j))));
166 7001547 : closure_evalvoid(code);
167 7001547 : if (loop_break()) return 1;
168 : }
169 3535 : if (x2 == b) break;
170 3486 : set_lex(-1, gen_0);
171 : }
172 49 : return 0;
173 : }
174 :
175 : /* vector of primes to squarefree factorization */
176 : static GEN
177 4255559 : zv_to_ZM(GEN v)
178 4255559 : { return mkmat2(zc_to_ZC(v), const_col(lg(v)-1,gen_1)); }
179 : /* vector of primes to negative squarefree factorization */
180 : static GEN
181 4255559 : zv_to_mZM(GEN v)
182 : {
183 4255559 : long i, l = lg(v);
184 4255559 : GEN w = cgetg(l+1, t_COL);
185 15388443 : gel(w,1) = gen_m1; for (i = 1; i < l; i++) gel(w,i+1) = utoipos(v[i]);
186 4255559 : return mkmat2(w, const_col(l,gen_1));
187 : }
188 : /* 0 <= a <= b. Using small consecutive chunks to 1) limit memory use, 2) allow
189 : * cheap early abort */
190 : static void
191 21 : forsquarefreepos(ulong a, ulong b, GEN code)
192 : {
193 21 : const ulong step = maxuu(1024, 2 * usqrt(b));
194 21 : pari_sp av = avma;
195 : ulong x1;
196 21 : if (no_sieve(a, b))
197 : {
198 : ulong n;
199 0 : for (n = a; n <= b; n++, set_avma(av))
200 : {
201 0 : GEN m = factoru(n);
202 0 : if (!uissquarefree_fact(m)) continue;
203 0 : set_lex(-1, mkvec2(utoipos(n), Flm_to_ZM(m)));
204 0 : closure_evalvoid(code); if (loop_break()) return;
205 : }
206 0 : return;
207 : }
208 3507 : for(x1 = a;; x1 += step, set_avma(av))
209 3486 : { /* beware overflow, fuse last two bins (avoid a tiny remainder) */
210 3507 : ulong j, lv, x2 = (b >= 2*step && b - 2*step >= x1)? x1-1 + step: b;
211 3507 : GEN v = vecfactorsquarefreeu(x1, x2);
212 3507 : lv = lg(v);
213 7003619 : for (j = 1; j < lv; j++) if (gel(v,j))
214 : {
215 4255559 : ulong n = x1-1 + j;
216 4255559 : set_lex(-1, mkvec2(utoipos(n), zv_to_ZM(gel(v,j))));
217 4255559 : closure_evalvoid(code); if (loop_break()) return;
218 : }
219 3507 : if (x2 == b) break;
220 3486 : set_lex(-1, gen_0);
221 : }
222 : }
223 : /* 0 <= a <= b. Loop from -b, ... -a through squarefree integers */
224 : static void
225 21 : forsquarefreeneg(ulong a, ulong b, GEN code)
226 : {
227 21 : const ulong step = maxuu(1024, 2 * usqrt(b));
228 21 : pari_sp av = avma;
229 : ulong x2;
230 21 : if (no_sieve(a, b))
231 : {
232 : ulong n;
233 0 : for (n = b; n >= a; n--, set_avma(av))
234 : {
235 0 : GEN m = factoru(n);
236 0 : if (!uissquarefree_fact(m)) continue;
237 0 : set_lex(-1, mkvec2(utoineg(n), zv_to_mZM(gel(m,1))));
238 0 : closure_evalvoid(code); if (loop_break()) return;
239 : }
240 0 : return;
241 : }
242 3507 : for(x2 = b;; x2 -= step, set_avma(av))
243 3486 : { /* beware overflow, fuse last two bins (avoid a tiny remainder) */
244 3507 : ulong j, x1 = (x2 >= 2*step && x2-2*step >= a)? x2+1 - step: a;
245 3507 : GEN v = vecfactorsquarefreeu(x1, x2);
246 7003619 : for (j = lg(v)-1; j > 0; j--) if (gel(v,j))
247 : {
248 4255559 : ulong n = x1-1 + j;
249 4255559 : set_lex(-1, mkvec2(utoineg(n), zv_to_mZM(gel(v,j))));
250 4255559 : closure_evalvoid(code); if (loop_break()) return;
251 : }
252 3507 : if (x1 == a) break;
253 3486 : set_lex(-1, gen_0);
254 : }
255 : }
256 : void
257 35 : forsquarefree(GEN a, GEN b, GEN code)
258 : {
259 35 : pari_sp av = avma;
260 : long s;
261 35 : if (typ(a) != t_INT) pari_err_TYPE("forsquarefree", a);
262 35 : if (typ(b) != t_INT) pari_err_TYPE("forsquarefree", b);
263 35 : if (cmpii(a,b) > 0) return;
264 35 : s = signe(a); push_lex(NULL,code);
265 35 : if (s < 0)
266 : {
267 21 : if (signe(b) <= 0)
268 14 : forsquarefreeneg(itou(b), itou(a), code);
269 : else
270 : {
271 7 : forsquarefreeneg(1, itou(a), code);
272 7 : forsquarefreepos(1, itou(b), code);
273 : }
274 : }
275 : else
276 14 : forsquarefreepos(itou(a), itou(b), code);
277 35 : pop_lex(1); set_avma(av);
278 : }
279 :
280 : /* convert factoru(n) to factor(-n); M pre-allocated factorization matrix
281 : * with (-1)^1 already set */
282 : static void
283 7001582 : Flm2negfact(GEN v, GEN M)
284 : {
285 7001582 : GEN p = gel(v,1), e = gel(v,2), P = gel(M,1), E = gel(M,2);
286 7001582 : long i, l = lg(p);
287 26980058 : for (i = 1; i < l; i++)
288 : {
289 19978476 : gel(P,i+1) = utoipos(p[i]);
290 19978476 : gel(E,i+1) = utoipos(e[i]);
291 : }
292 7001582 : setlg(P,l+1);
293 7001582 : setlg(E,l+1);
294 7001582 : }
295 : /* 0 < a <= b, from -b to -a */
296 : static int
297 84 : forfactoredneg(ulong a, ulong b, GEN code)
298 : {
299 84 : ulong x2, step = maxuu(2 * usqrt(b), 1024);
300 : GEN P, E, M;
301 : pari_sp av;
302 :
303 84 : P = cgetg(18, t_COL); gel(P,1) = gen_m1;
304 84 : E = cgetg(18, t_COL); gel(E,1) = gen_1;
305 84 : M = mkmat2(P,E);
306 84 : av = avma;
307 84 : if (no_sieve(a, b))
308 : {
309 : ulong n;
310 0 : for (n = b; n >= a; n--, set_avma(av))
311 : {
312 0 : GEN m = factoru(n);
313 0 : Flm2negfact(m, M);
314 0 : set_lex(-1, mkvec2(utoineg(n), M));
315 0 : closure_evalvoid(code); if (loop_break()) return 1;
316 : }
317 0 : return 0;
318 : }
319 3570 : for (x2 = b;; x2 -= step, set_avma(av))
320 3486 : { /* beware overflow, fuse last two bins (avoid a tiny remainder) */
321 3570 : ulong j, x1 = (x2 >= 2*step && x2-2*step >= a)? x2+1 - step: a;
322 3570 : GEN v = vecfactoru_i(x1, x2);
323 7005131 : for (j = lg(v)-1; j; j--)
324 : { /* run backward: from factor(x1..x2) to factor(-x2..-x1) */
325 7001582 : ulong n = x1-1 + j;
326 7001582 : Flm2negfact(gel(v,j), M);
327 7001582 : set_lex(-1, mkvec2(utoineg(n), M));
328 7001582 : closure_evalvoid(code); if (loop_break()) return 1;
329 : }
330 3549 : if (x1 == a) break;
331 3486 : set_lex(-1, gen_0);
332 : }
333 63 : return 0;
334 : }
335 : static int
336 70 : eval0(GEN code)
337 : {
338 70 : pari_sp av = avma;
339 70 : set_lex(-1, mkvec2(gen_0, mkmat2(mkcol(gen_0),mkcol(gen_1))));
340 70 : closure_evalvoid(code); set_avma(av);
341 70 : return loop_break();
342 : }
343 : void
344 140 : forfactored(GEN a, GEN b, GEN code)
345 : {
346 140 : pari_sp av = avma;
347 140 : long sa, sb, stop = 0;
348 140 : if (typ(a) != t_INT) pari_err_TYPE("forfactored", a);
349 140 : if (typ(b) != t_INT) pari_err_TYPE("forfactored", b);
350 140 : if (cmpii(a,b) > 0) return;
351 133 : push_lex(NULL,code);
352 133 : sa = signe(a);
353 133 : sb = signe(b);
354 133 : if (sa < 0)
355 : {
356 84 : stop = forfactoredneg((sb < 0)? uel(b,2): 1UL, itou(a), code);
357 84 : if (!stop && sb >= 0) stop = eval0(code);
358 84 : if (!stop && sb > 0) forfactoredpos(1UL, b[2], code);
359 : }
360 : else
361 : {
362 49 : if (!sa) stop = eval0(code);
363 49 : if (!stop && sb) forfactoredpos(sa? uel(a,2): 1UL, itou(b), code);
364 : }
365 133 : pop_lex(1); set_avma(av);
366 : }
367 : void
368 1794728 : whilepari(GEN a, GEN b)
369 : {
370 1794728 : pari_sp av = avma;
371 : for(;;)
372 17595068 : {
373 19389796 : GEN res = closure_evalnobrk(a);
374 19389796 : if (gequal0(res)) break;
375 17595117 : set_avma(av);
376 17595117 : closure_evalvoid(b); if (loop_break()) break;
377 : }
378 1794728 : set_avma(av);
379 1794728 : }
380 :
381 : void
382 222074 : untilpari(GEN a, GEN b)
383 : {
384 222074 : pari_sp av = avma;
385 : for(;;)
386 1456780 : {
387 : GEN res;
388 1678854 : closure_evalvoid(b); if (loop_break()) break;
389 1678854 : res = closure_evalnobrk(a);
390 1678854 : if (!gequal0(res)) break;
391 1456780 : set_avma(av);
392 : }
393 222074 : set_avma(av);
394 222074 : }
395 :
396 28 : static int negcmp(GEN x, GEN y) { return gcmp(y,x); }
397 :
398 : void
399 1638 : forstep(GEN a, GEN b, GEN s, GEN code)
400 : {
401 : long ss, i;
402 1638 : pari_sp av, av0 = avma;
403 1638 : GEN v = NULL;
404 : int (*cmp)(GEN,GEN);
405 :
406 1638 : b = gcopy(b);
407 1638 : s = gcopy(s); av = avma;
408 1638 : switch(typ(s))
409 : {
410 14 : case t_VEC: case t_COL: ss = gsigne(vecsum(s)); v = s; break;
411 21 : case t_INTMOD:
412 21 : if (typ(a) != t_INT) a = gceil(a);
413 21 : a = addii(a, modii(subii(gel(s,2),a), gel(s,1)));
414 21 : s = gel(s,1); /* FALL THROUGH */
415 1624 : default: ss = gsigne(s);
416 : }
417 1638 : if (!ss) pari_err_DOMAIN("forstep","step","=",gen_0,s);
418 1631 : cmp = (ss > 0)? &gcmp: &negcmp;
419 1631 : i = 0;
420 1631 : push_lex(a,code);
421 49756 : while (cmp(a,b) <= 0)
422 : {
423 48125 : closure_evalvoid(code); if (loop_break()) break;
424 48125 : if (v)
425 : {
426 98 : if (++i >= lg(v)) i = 1;
427 98 : s = gel(v,i);
428 : }
429 48125 : a = get_lex(-1); a = gadd(a,s);
430 :
431 48125 : if (gc_needed(av,1))
432 : {
433 0 : if (DEBUGMEM>1) pari_warn(warnmem,"forstep");
434 0 : a = gerepileupto(av,a);
435 : }
436 48125 : set_lex(-1,a);
437 : }
438 1631 : pop_lex(1); set_avma(av0);
439 1631 : }
440 :
441 : static void
442 28 : _fordiv(GEN a, GEN code, GEN (*D)(GEN))
443 : {
444 28 : pari_sp av = avma;
445 : long i, l;
446 28 : GEN t = D(a);
447 28 : push_lex(gen_0,code); l = lg(t);
448 231 : for (i=1; i<l; i++)
449 : {
450 203 : set_lex(-1,gel(t,i));
451 203 : closure_evalvoid(code); if (loop_break()) break;
452 : }
453 28 : pop_lex(1); set_avma(av);
454 28 : }
455 : void
456 14 : fordiv(GEN a, GEN code) { return _fordiv(a, code, &divisors); }
457 : void
458 14 : fordivfactored(GEN a, GEN code) { return _fordiv(a, code, &divisors_factored); }
459 :
460 : /* Embedded for loops:
461 : * fl = 0: execute ch (a), where a = (ai) runs through all n-uplets in
462 : * [m1,M1] x ... x [mn,Mn]
463 : * fl = 1: impose a1 <= ... <= an
464 : * fl = 2: a1 < ... < an
465 : */
466 : /* increment and return d->a [over integers]*/
467 : static GEN
468 183656 : _next_i(forvec_t *d)
469 : {
470 183656 : long i = d->n;
471 183656 : if (d->first) { d->first = 0; return (GEN)d->a; }
472 : for (;;) {
473 235762 : if (cmpii(d->a[i], d->M[i]) < 0) {
474 183280 : d->a[i] = incloop(d->a[i]);
475 183280 : return (GEN)d->a;
476 : }
477 52482 : d->a[i] = resetloop(d->a[i], d->m[i]);
478 52482 : if (--i <= 0) return NULL;
479 : }
480 : }
481 : /* increment and return d->a [generic]*/
482 : static GEN
483 63 : _next(forvec_t *d)
484 : {
485 63 : long i = d->n;
486 63 : if (d->first) { d->first = 0; return (GEN)d->a; }
487 : for (;;) {
488 98 : d->a[i] = gaddgs(d->a[i], 1);
489 98 : if (gcmp(d->a[i], d->M[i]) <= 0) return (GEN)d->a;
490 49 : d->a[i] = d->m[i];
491 49 : if (--i <= 0) return NULL;
492 : }
493 : }
494 :
495 : /* nondecreasing order [over integers] */
496 : static GEN
497 206 : _next_le_i(forvec_t *d)
498 : {
499 206 : long i = d->n;
500 206 : if (d->first) { d->first = 0; return (GEN)d->a; }
501 : for (;;) {
502 294 : if (cmpii(d->a[i], d->M[i]) < 0)
503 : {
504 152 : d->a[i] = incloop(d->a[i]);
505 : /* m[i] < a[i] <= M[i] <= M[i+1] */
506 233 : while (i < d->n)
507 : {
508 : GEN t;
509 81 : i++;
510 81 : if (cmpii(d->a[i-1], d->a[i]) <= 0) continue;
511 : /* a[i] < a[i-1] <= M[i-1] <= M[i] */
512 81 : t = d->a[i-1]; if (cmpii(t, d->m[i]) < 0) t = d->m[i];
513 81 : d->a[i] = resetloop(d->a[i], t);/*a[i]:=max(a[i-1],m[i])*/
514 : }
515 152 : return (GEN)d->a;
516 : }
517 142 : d->a[i] = resetloop(d->a[i], d->m[i]);
518 142 : if (--i <= 0) return NULL;
519 : }
520 : }
521 : /* nondecreasing order [generic] */
522 : static GEN
523 154 : _next_le(forvec_t *d)
524 : {
525 154 : long i = d->n;
526 154 : if (d->first) { d->first = 0; return (GEN)d->a; }
527 : for (;;) {
528 266 : d->a[i] = gaddgs(d->a[i], 1);
529 266 : if (gcmp(d->a[i], d->M[i]) <= 0)
530 : {
531 224 : while (i < d->n)
532 : {
533 : GEN c;
534 98 : i++;
535 98 : if (gcmp(d->a[i-1], d->a[i]) <= 0) continue;
536 : /* M[i] >= M[i-1] >= a[i-1] > a[i] */
537 98 : c = gceil(gsub(d->a[i-1], d->a[i]));
538 98 : d->a[i] = gadd(d->a[i], c);
539 : /* a[i-1] <= a[i] < M[i-1] + 1 => a[i] < M[i]+1 => a[i] <= M[i] */
540 : }
541 126 : return (GEN)d->a;
542 : }
543 140 : d->a[i] = d->m[i];
544 140 : if (--i <= 0) return NULL;
545 : }
546 : }
547 : /* strictly increasing order [over integers] */
548 : static GEN
549 1173574 : _next_lt_i(forvec_t *d)
550 : {
551 1173574 : long i = d->n;
552 1173574 : if (d->first) { d->first = 0; return (GEN)d->a; }
553 : for (;;) {
554 1290100 : if (cmpii(d->a[i], d->M[i]) < 0)
555 : {
556 1159954 : d->a[i] = incloop(d->a[i]);
557 : /* m[i] < a[i] <= M[i] < M[i+1] */
558 1276466 : while (i < d->n)
559 : {
560 : pari_sp av;
561 : GEN t;
562 116512 : i++;
563 116512 : if (cmpii(d->a[i-1], d->a[i]) < 0) continue;
564 116512 : av = avma;
565 : /* M[i] > M[i-1] >= a[i-1] */
566 116512 : t = addiu(d->a[i-1],1); if (cmpii(t, d->m[i]) < 0) t = d->m[i];
567 116512 : d->a[i] = resetloop(d->a[i], t);/*a[i]:=max(a[i-1]+1,m[i]) <= M[i]*/
568 116512 : set_avma(av);
569 : }
570 1159954 : return (GEN)d->a;
571 : }
572 130146 : d->a[i] = resetloop(d->a[i], d->m[i]);
573 130146 : if (--i <= 0) return NULL;
574 : }
575 : }
576 : /* strictly increasing order [generic] */
577 : static GEN
578 84 : _next_lt(forvec_t *d)
579 : {
580 84 : long i = d->n;
581 84 : if (d->first) { d->first = 0; return (GEN)d->a; }
582 : for (;;) {
583 133 : d->a[i] = gaddgs(d->a[i], 1);
584 133 : if (gcmp(d->a[i], d->M[i]) <= 0)
585 : {
586 91 : while (i < d->n)
587 : {
588 : GEN c;
589 35 : i++;
590 35 : if (gcmp(d->a[i-1], d->a[i]) < 0) continue;
591 : /* M[i] > M[i-1] >= a[i-1] >= a[i] */
592 35 : c = addiu(gfloor(gsub(d->a[i-1], d->a[i])), 1); /* > a[i-1] - a[i] */
593 35 : d->a[i] = gadd(d->a[i], c);
594 : /* a[i-1] < a[i] <= M[i-1] + 1 => a[i] < M[i]+1 => a[i] <= M[i] */
595 : }
596 56 : return (GEN)d->a;
597 : }
598 77 : d->a[i] = d->m[i];
599 77 : if (--i <= 0) return NULL;
600 : }
601 : }
602 :
603 : /* on Z^n /(cyc Z^n) [over integers]
604 : * torsion (cyc>0) and free (cyc=0) components may be interleaved */
605 : static GEN
606 8463 : _next_mod_cyc(forvec_t *d)
607 : { /* keep free components indices t1 < t2 last nonzero < t3 */
608 8463 : long t, t1 = 0, t2 = 0, t3 = 0;
609 8463 : if (d->first) { d->first = 0; return (GEN)d->a; }
610 27293 : for (t = d->n; t > 0; t--)
611 : {
612 24332 : if (signe(d->M[t]) > 0)
613 : { /* torsion component */
614 10738 : d->a[t] = incloop(d->a[t]);
615 10738 : if (cmpii(d->a[t], d->M[t]) < 0) return (GEN)d->a;
616 5278 : d->a[t] = resetloop(d->a[t], gen_0);
617 : }
618 : else
619 : { /* set or update t1,t2,t3 */
620 13594 : if (t2 && !t1) t1 = t;
621 13594 : if (!t2 && signe(d->a[t])) t2 = t;
622 13594 : if (!t2) t3 = t;
623 : }
624 : }
625 2961 : if (!t3 && !t2) return NULL; /* no free component, stop */
626 2947 : if (!t2) d->a[t3] = resetloop(d->a[t3], gen_m1);
627 2919 : else if (!t3 && signe(d->a[t2]) < 0) togglesign(d->a[t2]);
628 1757 : else if (signe(d->a[t2]) < 0)
629 : {
630 315 : d->a[t2] = incloop(d->a[t2]);
631 315 : d->a[t3] = resetloop(d->a[t3], gen_m1);
632 : }
633 1442 : else if (!t1) { d->a[t2] = incloop(d->a[t2]); togglesign(d->a[t2]); }
634 : else
635 : {
636 1197 : if (signe(d->a[t1]) < 0)
637 490 : { d->a[t2] = incloop(d->a[t2]); togglesign(d->a[t2]); }
638 : else
639 707 : { togglesign(d->a[t2]); d->a[t2] = incloop(d->a[t2]); }
640 1197 : d->a[t1] = incloop(d->a[t1]);
641 : }
642 2947 : return (GEN)d->a;
643 : }
644 : /* for forvec(v=[],) */
645 : static GEN
646 14 : _next_void(forvec_t *d)
647 : {
648 14 : if (d->first) { d->first = 0; return (GEN)d->a; }
649 7 : return NULL;
650 : }
651 : static int
652 7144 : RgV_is_ZV_nonneg(GEN x)
653 : {
654 : long i;
655 7298 : for (i = lg(x)-1; i > 0; i--)
656 7256 : if (typ(gel(x,i)) != t_INT || signe(gel(x, i)) < 0) return 0;
657 42 : return 1;
658 : }
659 : /* x assumed to be cyc vector, l>1 */
660 : static int
661 42 : forvec_mod_cyc_init(forvec_t *d, GEN x)
662 : {
663 42 : long i, tx = typ(x), l = lg(x);
664 42 : d->a = (GEN*)cgetg(l,tx); /* current */
665 42 : d->M = (GEN*)cgetg(l,tx); /* cyc */
666 175 : for (i = 1; i < l; i++)
667 : {
668 133 : d->a[i] = setloop(gen_0);
669 133 : d->M[i] = setloop(gel(x, i));
670 : }
671 42 : d->first = 1;
672 42 : d->n = l-1;
673 42 : d->m = NULL;
674 42 : d->next = &_next_mod_cyc;
675 42 : return 1;
676 : }
677 :
678 : /* Initialize minima (m) and maxima (M); guarantee M[i] - m[i] integer and
679 : * if flag = 1: m[i-1] <= m[i] <= M[i] <= M[i+1]
680 : * if flag = 2: m[i-1] < m[i] <= M[i] < M[i+1],
681 : * for all i */
682 : int
683 7151 : forvec_init(forvec_t *d, GEN x, long flag)
684 : {
685 7151 : long i, tx = typ(x), l = lg(x), t = t_INT;
686 7151 : if (!is_vec_t(tx)) pari_err_TYPE("forvec [not a vector]", x);
687 7151 : if (l > 1 && RgV_is_ZV_nonneg(x))
688 42 : return forvec_mod_cyc_init(d, x);
689 7109 : d->first = 1;
690 7109 : d->n = l - 1;
691 7109 : d->a = (GEN*)cgetg(l,tx);
692 7109 : d->m = (GEN*)cgetg(l,tx);
693 7109 : d->M = (GEN*)cgetg(l,tx);
694 7109 : if (l == 1) { d->next = &_next_void; return 1; }
695 21537 : for (i = 1; i < l; i++)
696 : {
697 14470 : GEN a, e = gel(x,i), m = gel(e,1), M = gel(e,2);
698 14470 : tx = typ(e);
699 14470 : if (! is_vec_t(tx) || lg(e)!=3)
700 21 : pari_err_TYPE("forvec [expected vector not of type [min,MAX]]",e);
701 14449 : if (typ(m) != t_INT) t = t_REAL;
702 14449 : if (i > 1) switch(flag)
703 : {
704 62 : case 1: /* a >= m[i-1] - m */
705 62 : a = gceil(gsub(d->m[i-1], m));
706 62 : if (typ(a) != t_INT) pari_err_TYPE("forvec",a);
707 62 : if (signe(a) > 0) m = gadd(m, a); else m = gcopy(m);
708 62 : break;
709 6859 : case 2: /* a > m[i-1] - m */
710 6859 : a = gfloor(gsub(d->m[i-1], m));
711 6859 : if (typ(a) != t_INT) pari_err_TYPE("forvec",a);
712 6859 : a = addiu(a, 1);
713 6859 : if (signe(a) > 0) m = gadd(m, a); else m = gcopy(m);
714 6859 : break;
715 440 : default: m = gcopy(m);
716 440 : break;
717 : }
718 14449 : M = gadd(m, gfloor(gsub(M,m))); /* ensure M-m is an integer */
719 14442 : if (gcmp(m,M) > 0) { d->a = NULL; d->next = &_next; return 0; }
720 14435 : d->m[i] = m;
721 14435 : d->M[i] = M;
722 : }
723 7129 : if (flag == 1) for (i = l-2; i >= 1; i--)
724 : {
725 62 : GEN M = d->M[i], a = gfloor(gsub(d->M[i+1], M));
726 62 : if (typ(a) != t_INT) pari_err_TYPE("forvec",a);
727 : /* M[i]+a <= M[i+1] */
728 62 : if (signe(a) < 0) d->M[i] = gadd(M, a);
729 : }
730 13878 : else if (flag == 2) for (i = l-2; i >= 1; i--)
731 : {
732 6852 : GEN M = d->M[i], a = gceil(gsub(d->M[i+1], M));
733 6852 : if (typ(a) != t_INT) pari_err_TYPE("forvec",a);
734 6852 : a = subiu(a, 1);
735 : /* M[i]+a < M[i+1] */
736 6852 : if (signe(a) < 0) d->M[i] = gadd(M, a);
737 : }
738 7067 : if (t == t_INT) {
739 21348 : for (i = 1; i < l; i++) {
740 14316 : d->a[i] = setloop(d->m[i]);
741 14316 : if (typ(d->M[i]) != t_INT) d->M[i] = gfloor(d->M[i]);
742 : }
743 : } else {
744 140 : for (i = 1; i < l; i++) d->a[i] = d->m[i];
745 : }
746 7067 : switch(flag)
747 : {
748 195 : case 0: d->next = t==t_INT? &_next_i: &_next; break;
749 41 : case 1: d->next = t==t_INT? &_next_le_i: &_next_le; break;
750 6824 : case 2: d->next = t==t_INT? &_next_lt_i: &_next_lt; break;
751 7 : default: pari_err_FLAG("forvec");
752 : }
753 7060 : return 1;
754 : }
755 : GEN
756 1366214 : forvec_next(forvec_t *d) { return d->next(d); }
757 :
758 : void
759 7070 : forvec(GEN x, GEN code, long flag)
760 : {
761 7070 : pari_sp av = avma;
762 : forvec_t T;
763 : GEN v;
764 7070 : if (!forvec_init(&T, x, flag)) { set_avma(av); return; }
765 7028 : push_lex((GEN)T.a, code);
766 1365602 : while ((v = forvec_next(&T)))
767 : {
768 1358602 : closure_evalvoid(code);
769 1358602 : if (loop_break()) break;
770 : }
771 7028 : pop_lex(1); set_avma(av);
772 : }
773 :
774 : /********************************************************************/
775 : /** **/
776 : /** SUMS **/
777 : /** **/
778 : /********************************************************************/
779 :
780 : GEN
781 70231 : somme(GEN a, GEN b, GEN code, GEN x)
782 : {
783 70231 : pari_sp av, av0 = avma;
784 : GEN p1;
785 :
786 70231 : if (typ(a) != t_INT) pari_err_TYPE("sum",a);
787 70231 : if (!x) x = gen_0;
788 70231 : if (gcmp(b,a) < 0) return gcopy(x);
789 :
790 70231 : b = gfloor(b);
791 70231 : a = setloop(a);
792 70231 : av=avma;
793 70231 : push_lex(a,code);
794 : for(;;)
795 : {
796 1870456 : p1 = closure_evalnobrk(code);
797 1870456 : x=gadd(x,p1); if (cmpii(a,b) >= 0) break;
798 1800225 : a = incloop(a);
799 1800225 : if (gc_needed(av,1))
800 : {
801 0 : if (DEBUGMEM>1) pari_warn(warnmem,"sum");
802 0 : x = gerepileupto(av,x);
803 : }
804 1800225 : set_lex(-1,a);
805 : }
806 70231 : pop_lex(1); return gerepileupto(av0,x);
807 : }
808 :
809 : static GEN
810 28 : sum_init(GEN x0, GEN t)
811 : {
812 28 : long tp = typ(t);
813 : GEN x;
814 28 : if (is_vec_t(tp))
815 : {
816 7 : x = const_vec(lg(t)-1, x0);
817 7 : settyp(x, tp);
818 : }
819 : else
820 21 : x = x0;
821 28 : return x;
822 : }
823 :
824 : GEN
825 28 : suminf(void *E, GEN (*eval)(void *, GEN), GEN a, long bit)
826 : {
827 28 : long fl = 0, G = bit + 1;
828 28 : pari_sp av0 = avma, av;
829 28 : GEN x = NULL, _1;
830 :
831 28 : if (typ(a) != t_INT) pari_err_TYPE("suminf",a);
832 28 : a = setloop(a); av = avma;
833 : for(;;)
834 15617 : {
835 15645 : GEN t = eval(E, a);
836 15645 : if (!x) _1 = x = sum_init(real_1_bit(bit), t);
837 :
838 15645 : x = gadd(x,t);
839 15645 : if (!gequal0(t) && gexpo(t) > gexpo(x)-G)
840 15449 : fl = 0;
841 196 : else if (++fl == 3)
842 28 : break;
843 15617 : a = incloop(a);
844 15617 : if (gc_needed(av,1))
845 : {
846 0 : if (DEBUGMEM>1) pari_warn(warnmem,"suminf");
847 0 : gerepileall(av,2, &x, &_1);
848 : }
849 : }
850 28 : return gerepileupto(av0, gsub(x, _1));
851 : }
852 : GEN
853 28 : suminf0(GEN a, GEN code, long bit)
854 28 : { EXPR_WRAP(code, suminf(EXPR_ARG, a, bit)); }
855 :
856 : GEN
857 56 : sumdivexpr(GEN num, GEN code)
858 : {
859 56 : pari_sp av = avma;
860 56 : GEN y = gen_0, t = divisors(num);
861 56 : long i, l = lg(t);
862 :
863 56 : push_lex(gen_0, code);
864 9352 : for (i=1; i<l; i++)
865 : {
866 9296 : set_lex(-1,gel(t,i));
867 9296 : y = gadd(y, closure_evalnobrk(code));
868 : }
869 56 : pop_lex(1); return gerepileupto(av,y);
870 : }
871 :
872 : GEN
873 49 : sumdivmultexpr(void *D, GEN (*fun)(void*, GEN), GEN num)
874 : {
875 49 : pari_sp av = avma;
876 49 : GEN y = gen_1, P,E;
877 49 : int isint = divisors_init(num, &P,&E);
878 49 : long i, l = lg(P);
879 : GEN (*mul)(GEN,GEN);
880 :
881 49 : if (l == 1) return gc_const(av, gen_1);
882 49 : mul = isint? mulii: gmul;
883 224 : for (i=1; i<l; i++)
884 : {
885 175 : GEN p = gel(P,i), q = p, z = gen_1;
886 175 : long j, e = E[i];
887 581 : for (j = 1; j <= e; j++, q = mul(q, p))
888 : {
889 581 : z = gadd(z, fun(D, q));
890 581 : if (j == e) break;
891 : }
892 175 : y = gmul(y, z);
893 : }
894 49 : return gerepileupto(av,y);
895 : }
896 :
897 : GEN
898 49 : sumdivmultexpr0(GEN num, GEN code)
899 49 : { EXPR_WRAP(code, sumdivmultexpr(EXPR_ARG, num)) }
900 :
901 : /********************************************************************/
902 : /** **/
903 : /** PRODUCTS **/
904 : /** **/
905 : /********************************************************************/
906 :
907 : GEN
908 120694 : produit(GEN a, GEN b, GEN code, GEN x)
909 : {
910 120694 : pari_sp av, av0 = avma;
911 : GEN p1;
912 :
913 120694 : if (typ(a) != t_INT) pari_err_TYPE("prod",a);
914 120694 : if (!x) x = gen_1;
915 120694 : if (gcmp(b,a) < 0) return gcopy(x);
916 :
917 115416 : b = gfloor(b);
918 115416 : a = setloop(a);
919 115416 : av=avma;
920 115416 : push_lex(a,code);
921 : for(;;)
922 : {
923 349804 : p1 = closure_evalnobrk(code);
924 349804 : x = gmul(x,p1); if (cmpii(a,b) >= 0) break;
925 234388 : a = incloop(a);
926 234388 : if (gc_needed(av,1))
927 : {
928 0 : if (DEBUGMEM>1) pari_warn(warnmem,"prod");
929 0 : x = gerepileupto(av,x);
930 : }
931 234388 : set_lex(-1,a);
932 : }
933 115416 : pop_lex(1); return gerepileupto(av0,x);
934 : }
935 :
936 : GEN
937 14 : prodinf(void *E, GEN (*eval)(void *, GEN), GEN a, long prec)
938 : {
939 14 : pari_sp av0 = avma, av;
940 : long fl,G;
941 14 : GEN p1,x = real_1(prec);
942 :
943 14 : if (typ(a) != t_INT) pari_err_TYPE("prodinf",a);
944 14 : a = setloop(a);
945 14 : av = avma;
946 14 : fl=0; G = -prec-5;
947 : for(;;)
948 : {
949 1897 : p1 = eval(E, a); if (gequal0(p1)) { x = p1; break; }
950 1897 : x = gmul(x,p1); a = incloop(a);
951 1897 : p1 = gsubgs(p1, 1);
952 1897 : if (gequal0(p1) || gexpo(p1) <= G) { if (++fl==3) break; } else fl=0;
953 1883 : if (gc_needed(av,1))
954 : {
955 0 : if (DEBUGMEM>1) pari_warn(warnmem,"prodinf");
956 0 : x = gerepileupto(av,x);
957 : }
958 : }
959 14 : return gerepilecopy(av0,x);
960 : }
961 : GEN
962 7 : prodinf1(void *E, GEN (*eval)(void *, GEN), GEN a, long prec)
963 : {
964 7 : pari_sp av0 = avma, av;
965 : long fl,G;
966 7 : GEN p1,p2,x = real_1(prec);
967 :
968 7 : if (typ(a) != t_INT) pari_err_TYPE("prodinf1",a);
969 7 : a = setloop(a);
970 7 : av = avma;
971 7 : fl=0; G = -prec-5;
972 : for(;;)
973 : {
974 952 : p2 = eval(E, a); p1 = gaddgs(p2,1);
975 952 : if (gequal0(p1)) { x = p1; break; }
976 952 : x = gmul(x,p1); a = incloop(a);
977 952 : if (gequal0(p2) || gexpo(p2) <= G) { if (++fl==3) break; } else fl=0;
978 945 : if (gc_needed(av,1))
979 : {
980 0 : if (DEBUGMEM>1) pari_warn(warnmem,"prodinf1");
981 0 : x = gerepileupto(av,x);
982 : }
983 : }
984 7 : return gerepilecopy(av0,x);
985 : }
986 : GEN
987 28 : prodinf0(GEN a, GEN code, long flag, long prec)
988 : {
989 28 : switch(flag)
990 : {
991 14 : case 0: EXPR_WRAP(code, prodinf (EXPR_ARG, a, prec));
992 7 : case 1: EXPR_WRAP(code, prodinf1(EXPR_ARG, a, prec));
993 : }
994 7 : pari_err_FLAG("prodinf");
995 : return NULL; /* LCOV_EXCL_LINE */
996 : }
997 :
998 : GEN
999 14 : prodeuler(void *E, GEN (*eval)(void *, GEN), GEN a, GEN b, long prec)
1000 : {
1001 14 : pari_sp av, av0 = avma;
1002 14 : GEN x = real_1(prec), prime;
1003 : forprime_t T;
1004 :
1005 14 : av = avma;
1006 14 : if (!forprime_init(&T, a,b)) return gc_const(av, x);
1007 :
1008 14 : av = avma;
1009 8645 : while ( (prime = forprime_next(&T)) )
1010 : {
1011 8631 : x = gmul(x, eval(E, prime));
1012 8631 : if (gc_needed(av,1))
1013 : {
1014 0 : if (DEBUGMEM>1) pari_warn(warnmem,"prodeuler");
1015 0 : x = gerepilecopy(av, x);
1016 : }
1017 : }
1018 14 : return gerepilecopy(av0,x);
1019 : }
1020 : GEN
1021 14 : prodeuler0(GEN a, GEN b, GEN code, long prec)
1022 14 : { EXPR_WRAP(code, prodeuler(EXPR_ARG, a, b, prec)); }
1023 : GEN
1024 133 : direuler0(GEN a, GEN b, GEN code, GEN c)
1025 133 : { EXPR_WRAP(code, direuler(EXPR_ARG, a, b, c)); }
1026 :
1027 : /********************************************************************/
1028 : /** **/
1029 : /** VECTORS & MATRICES **/
1030 : /** **/
1031 : /********************************************************************/
1032 :
1033 : INLINE GEN
1034 2841666 : copyupto(GEN z, GEN t)
1035 : {
1036 2841666 : if (is_universal_constant(z) || (z>(GEN)pari_mainstack->bot && z<=t))
1037 2841659 : return z;
1038 : else
1039 7 : return gcopy(z);
1040 : }
1041 :
1042 : GEN
1043 132724 : vecexpr0(GEN vec, GEN code, GEN pred)
1044 : {
1045 132724 : switch(typ(vec))
1046 : {
1047 21 : case t_LIST:
1048 : {
1049 21 : if (list_typ(vec)==t_LIST_MAP)
1050 7 : vec = mapdomain_shallow(vec);
1051 : else
1052 14 : vec = list_data(vec);
1053 21 : if (!vec) return cgetg(1, t_VEC);
1054 14 : break;
1055 : }
1056 7 : case t_VECSMALL:
1057 7 : vec = vecsmall_to_vec(vec);
1058 7 : break;
1059 132696 : case t_VEC: case t_COL: case t_MAT: break;
1060 0 : default: pari_err_TYPE("[_|_<-_,_]",vec);
1061 : }
1062 132717 : if (pred && code)
1063 448 : EXPR_WRAP(code,vecselapply((void*)pred,&gp_evalbool,EXPR_ARGUPTO,vec))
1064 132269 : else if (code)
1065 132269 : EXPR_WRAP(code,vecapply(EXPR_ARGUPTO,vec))
1066 : else
1067 0 : EXPR_WRAP(pred,vecselect(EXPR_ARGBOOL,vec))
1068 : }
1069 :
1070 : GEN
1071 2114 : vecexpr1(GEN vec, GEN code, GEN pred)
1072 : {
1073 2114 : GEN v = vecexpr0(vec, code, pred);
1074 2114 : return lg(v) == 1? v: shallowconcat1(v);
1075 : }
1076 :
1077 : GEN
1078 2369639 : vecteur(GEN nmax, GEN code)
1079 : {
1080 : GEN y, c;
1081 2369639 : long i, m = gtos(nmax);
1082 :
1083 2369639 : if (m < 0) pari_err_DOMAIN("vector", "dimension", "<", gen_0, stoi(m));
1084 2369625 : if (!code) return zerovec(m);
1085 23449 : c = cgetipos(3); /* left on stack */
1086 23449 : y = cgetg(m+1,t_VEC); push_lex(c, code);
1087 899284 : for (i=1; i<=m; i++)
1088 : {
1089 875849 : c[2] = i;
1090 875849 : gel(y,i) = copyupto(closure_evalnobrk(code), y);
1091 875835 : set_lex(-1,c);
1092 : }
1093 23435 : pop_lex(1); return y;
1094 : }
1095 :
1096 : GEN
1097 791 : vecteursmall(GEN nmax, GEN code)
1098 : {
1099 : pari_sp av;
1100 : GEN y, c;
1101 791 : long i, m = gtos(nmax);
1102 :
1103 791 : if (m < 0) pari_err_DOMAIN("vectorsmall", "dimension", "<", gen_0, stoi(m));
1104 784 : if (!code) return zero_zv(m);
1105 763 : c = cgetipos(3); /* left on stack */
1106 763 : y = cgetg(m+1,t_VECSMALL); push_lex(c,code);
1107 763 : av = avma;
1108 10186883 : for (i = 1; i <= m; i++)
1109 : {
1110 10186127 : c[2] = i;
1111 10186127 : y[i] = gtos(closure_evalnobrk(code));
1112 10186120 : set_avma(av);
1113 10186120 : set_lex(-1,c);
1114 : }
1115 756 : pop_lex(1); return y;
1116 : }
1117 :
1118 : GEN
1119 1981 : vvecteur(GEN nmax, GEN n)
1120 : {
1121 1981 : GEN y = vecteur(nmax,n);
1122 1974 : settyp(y,t_COL); return y;
1123 : }
1124 :
1125 : GEN
1126 161168 : matrice(GEN nlig, GEN ncol, GEN code)
1127 : {
1128 : GEN c1, c2, y;
1129 : long i, m, n;
1130 :
1131 161168 : n = gtos(nlig);
1132 161168 : m = ncol? gtos(ncol): n;
1133 161168 : if (m < 0) pari_err_DOMAIN("matrix", "nbcols", "<", gen_0, stoi(m));
1134 161161 : if (n < 0) pari_err_DOMAIN("matrix", "nbrows", "<", gen_0, stoi(n));
1135 161154 : if (!m) return cgetg(1,t_MAT);
1136 161084 : if (!code || !n) return zeromatcopy(n, m);
1137 158606 : c1 = cgetipos(3); push_lex(c1,code);
1138 158606 : c2 = cgetipos(3); push_lex(c2,NULL); /* c1,c2 left on stack */
1139 158606 : y = cgetg(m+1,t_MAT);
1140 596862 : for (i = 1; i <= m; i++)
1141 : {
1142 438256 : GEN z = cgetg(n+1,t_COL);
1143 : long j;
1144 438256 : c2[2] = i; gel(y,i) = z;
1145 2404080 : for (j = 1; j <= n; j++)
1146 : {
1147 1965824 : c1[2] = j;
1148 1965824 : gel(z,j) = copyupto(closure_evalnobrk(code), y);
1149 1965824 : set_lex(-2,c1);
1150 1965824 : set_lex(-1,c2);
1151 : }
1152 : }
1153 158606 : pop_lex(2); return y;
1154 : }
1155 :
1156 : /********************************************************************/
1157 : /** **/
1158 : /** SUMMING SERIES **/
1159 : /** **/
1160 : /********************************************************************/
1161 : /* h = (2+2x)g'- g; g has t_INT coeffs */
1162 : static GEN
1163 1295 : delt(GEN g, long n)
1164 : {
1165 1295 : GEN h = cgetg(n+3,t_POL);
1166 : long k;
1167 1295 : h[1] = g[1];
1168 1295 : gel(h,2) = gel(g,2);
1169 359954 : for (k=1; k<n; k++)
1170 358659 : gel(h,k+2) = addii(mului(k+k+1,gel(g,k+2)), mului(k<<1,gel(g,k+1)));
1171 1295 : gel(h,n+2) = mului(n<<1, gel(g,n+1)); return h;
1172 : }
1173 :
1174 : #ifdef _MSC_VER /* Bill Daly: work around a MSVC bug */
1175 : #pragma optimize("g",off)
1176 : #endif
1177 : /* P = polzagier(n,m)(-X), unnormalized (P(0) != 1) */
1178 : static GEN
1179 84 : polzag1(long n, long m)
1180 : {
1181 84 : long d = n - m, i, k, d2, r, D;
1182 84 : pari_sp av = avma;
1183 : GEN g, T;
1184 :
1185 84 : if (d <= 0 || m < 0) return pol_0(0);
1186 77 : d2 = d << 1; r = (m+1) >> 1, D = (d+1) >> 1;
1187 77 : g = cgetg(d+2, t_POL);
1188 77 : g[1] = evalsigne(1)|evalvarn(0);
1189 77 : T = cgetg(d+1,t_VEC);
1190 : /* T[k+1] = binomial(2d,2k+1), 0 <= k < d */
1191 77 : gel(T,1) = utoipos(d2);
1192 1344 : for (k = 1; k < D; k++)
1193 : {
1194 1267 : long k2 = k<<1;
1195 1267 : gel(T,k+1) = diviiexact(mulii(gel(T,k), muluu(d2-k2+1, d2-k2)),
1196 1267 : muluu(k2,k2+1));
1197 : }
1198 1365 : for (; k < d; k++) gel(T,k+1) = gel(T,d-k);
1199 77 : gel(g,2) = gel(T,d); /* binomial(2d, 2(d-1)+1) */
1200 2632 : for (i = 1; i < d; i++)
1201 : {
1202 2555 : pari_sp av2 = avma;
1203 2555 : GEN s, t = gel(T,d-i); /* binomial(2d, 2(d-1-i)+1) */
1204 2555 : s = t;
1205 180635 : for (k = d-i; k < d; k++)
1206 : {
1207 178080 : long k2 = k<<1;
1208 178080 : t = diviiexact(mulii(t, muluu(d2-k2+1, d-k)), muluu(k2+1,k-(d-i)+1));
1209 178080 : s = addii(s, t);
1210 : }
1211 : /* g_i = sum_{d-1-i <= k < d}, binomial(2*d, 2*k+1)*binomial(k,d-1-i) */
1212 2555 : gel(g,i+2) = gerepileuptoint(av2, s);
1213 : }
1214 : /* sum_{0 <= i < d} g_i x^i * (x+x^2)^r */
1215 77 : g = RgX_mulXn(gmul(g, gpowgs(deg1pol(gen_1,gen_1,0),r)), r);
1216 77 : if (!odd(m)) g = delt(g, n);
1217 1337 : for (i = 1; i <= r; i++)
1218 : {
1219 1260 : g = delt(ZX_deriv(g), n);
1220 1260 : if (gc_needed(av,4))
1221 : {
1222 0 : if (DEBUGMEM>1) pari_warn(warnmem,"polzag, i = %ld/%ld", i,r);
1223 0 : g = gerepilecopy(av, g);
1224 : }
1225 : }
1226 77 : return g;
1227 : }
1228 : GEN
1229 35 : polzag(long n, long m)
1230 : {
1231 35 : pari_sp av = avma;
1232 35 : GEN g = polzag1(n,m);
1233 35 : if (lg(g) == 2) return g;
1234 28 : g = ZX_z_unscale(polzag1(n,m), -1);
1235 28 : return gerepileupto(av, RgX_Rg_div(g,gel(g,2)));
1236 : }
1237 :
1238 : /*0.39322 > 1/log_2(3+sqrt(8))*/
1239 : static ulong
1240 154 : sumalt_N(long prec)
1241 154 : { return (ulong)(0.39322*(prec + 7)); }
1242 :
1243 : GEN
1244 84 : sumalt(void *E, GEN (*eval)(void *, GEN), GEN a, long prec)
1245 : {
1246 : ulong k, N;
1247 84 : pari_sp av = avma, av2;
1248 : GEN s, az, c, d;
1249 :
1250 84 : if (typ(a) != t_INT) pari_err_TYPE("sumalt",a);
1251 84 : N = sumalt_N(prec);
1252 84 : d = powru(addsr(3, sqrtr(utor(8,prec))), N);
1253 84 : d = shiftr(addrr(d, invr(d)),-1);
1254 84 : a = setloop(a);
1255 84 : az = gen_m1; c = d;
1256 84 : s = gen_0;
1257 84 : av2 = avma;
1258 10752 : for (k=0; ; k++) /* k < N */
1259 : {
1260 10752 : c = addir(az,c); s = gadd(s, gmul(c, eval(E, a)));
1261 10752 : if (k==N-1) break;
1262 10668 : az = diviuuexact(muluui((N-k)<<1,N+k,az), k+1, (k<<1)+1);
1263 10668 : a = incloop(a); /* in place! */
1264 10668 : if (gc_needed(av,4))
1265 : {
1266 0 : if (DEBUGMEM>1) pari_warn(warnmem,"sumalt, k = %ld/%ld", k,N-1);
1267 0 : gerepileall(av2, 3, &az,&c,&s);
1268 : }
1269 : }
1270 84 : return gerepileupto(av, gdiv(s,d));
1271 : }
1272 :
1273 : GEN
1274 7 : sumalt2(void *E, GEN (*eval)(void *, GEN), GEN a, long prec)
1275 : {
1276 : long k, N;
1277 7 : pari_sp av = avma, av2;
1278 : GEN s, dn, pol;
1279 :
1280 7 : if (typ(a) != t_INT) pari_err_TYPE("sumalt",a);
1281 7 : N = (long)(0.307073*(prec + 5)); /*0.307073 > 1/log_2(\beta_B)*/
1282 7 : pol = ZX_div_by_X_1(polzag1(N,N>>1), &dn);
1283 7 : a = setloop(a);
1284 7 : N = degpol(pol);
1285 7 : s = gen_0;
1286 7 : av2 = avma;
1287 280 : for (k=0; k<=N; k++)
1288 : {
1289 280 : GEN t = itor(gel(pol,k+2), prec+EXTRAPREC64);
1290 280 : s = gadd(s, gmul(t, eval(E, a)));
1291 280 : if (k == N) break;
1292 273 : a = incloop(a); /* in place! */
1293 273 : if (gc_needed(av,4))
1294 : {
1295 0 : if (DEBUGMEM>1) pari_warn(warnmem,"sumalt2, k = %ld/%ld", k,N-1);
1296 0 : s = gerepileupto(av2, s);
1297 : }
1298 : }
1299 7 : return gerepileupto(av, gdiv(s,dn));
1300 : }
1301 :
1302 : GEN
1303 28 : sumalt0(GEN a, GEN code, long flag, long prec)
1304 : {
1305 28 : switch(flag)
1306 : {
1307 14 : case 0: EXPR_WRAP(code, sumalt (EXPR_ARG,a,prec));
1308 7 : case 1: EXPR_WRAP(code, sumalt2(EXPR_ARG,a,prec));
1309 7 : default: pari_err_FLAG("sumalt");
1310 : }
1311 : return NULL; /* LCOV_EXCL_LINE */
1312 : }
1313 :
1314 : /* For k > 0, set S[k*2^i] <- g(k*2^i), k*2^i <= N = #S.
1315 : * Only needed with k odd (but also works for g even). */
1316 : static void
1317 8953 : binsum(GEN S, ulong k, void *E, GEN (*f)(void *, GEN), GEN a,
1318 : long G, long prec)
1319 : {
1320 8953 : long e, i, N = lg(S)-1, l = expu(N / k); /* k 2^l <= N < k 2^(l+1) */
1321 8953 : pari_sp av = avma;
1322 8953 : GEN t = real_0(prec); /* unused unless f(a + k <<l) = 0 */
1323 :
1324 8953 : G -= l;
1325 8953 : if (!signe(a)) a = NULL;
1326 8953 : for (e = 0;; e++)
1327 5389657 : { /* compute g(k 2^l) with absolute error ~ 2^(G-l) */
1328 5398610 : GEN u, r = shifti(utoipos(k), l+e);
1329 5398610 : if (a) r = addii(r, a);
1330 5398610 : u = gtofp(f(E, r), prec);
1331 5398610 : if (typ(u) != t_REAL) pari_err_TYPE("sumpos",u);
1332 5398610 : if (!signe(u)) break;
1333 5398421 : if (!e)
1334 8764 : t = u;
1335 : else {
1336 5389657 : shiftr_inplace(u, e);
1337 5389657 : t = addrr(t,u); if (expo(u) < G) break;
1338 5380893 : if ((e & 0x1ff) == 0) t = gerepileuptoleaf(av, t);
1339 : }
1340 : }
1341 8953 : gel(S, k << l) = t = gerepileuptoleaf(av, t);
1342 : /* g(j) = 2g(2j) + f(a+j) for all j > 0 */
1343 17906 : for(i = l-1; i >= 0; i--)
1344 : { /* t ~ g(2 * k*2^i) with error ~ 2^(G-i-1) */
1345 : GEN u;
1346 8953 : av = avma; u = gtofp(f(E, a? addiu(a, k << i): utoipos(k << i)), prec);
1347 8953 : if (typ(u) != t_REAL) pari_err_TYPE("sumpos",u);
1348 8953 : t = addrr(gtofp(u,prec), mpshift(t,1)); /* ~ g(k*2^i) */
1349 8953 : gel(S, k << i) = t = gerepileuptoleaf(av, t);
1350 : }
1351 8953 : }
1352 : /* For k > 0, let g(k) := \sum_{e >= 0} 2^e f(a + k*2^e).
1353 : * Return [g(k), 1 <= k <= N] */
1354 : static GEN
1355 84 : sumpos_init(void *E, GEN (*f)(void *, GEN), GEN a, long N, long prec)
1356 : {
1357 84 : GEN S = cgetg(N+1,t_VEC);
1358 84 : long k, G = -prec - 5;
1359 9037 : for (k=1; k<=N; k+=2) binsum(S,k, E,f, a,G,prec);
1360 84 : return S;
1361 : }
1362 :
1363 : GEN
1364 70 : sumpos(void *E, GEN (*eval)(void *, GEN), GEN a, long prec)
1365 : {
1366 : ulong k, N;
1367 70 : pari_sp av = avma;
1368 : GEN s, az, c, d, S;
1369 :
1370 70 : if (typ(a) != t_INT) pari_err_TYPE("sumpos",a);
1371 70 : a = subiu(a,1);
1372 70 : N = sumalt_N(prec);
1373 70 : if (odd(N)) N++; /* extra precision for free */
1374 70 : d = powru(addsr(3, sqrtr(utor(8,prec))), N);
1375 70 : d = shiftr(addrr(d, invr(d)),-1);
1376 70 : az = gen_m1; c = d;
1377 :
1378 70 : S = sumpos_init(E, eval, a, N, prec);
1379 70 : s = gen_0;
1380 13454 : for (k=0; k<N; k++)
1381 : {
1382 : GEN t;
1383 13454 : c = addir(az,c);
1384 13454 : t = mulrr(gel(S,k+1), c);
1385 13454 : s = odd(k)? mpsub(s, t): mpadd(s, t);
1386 13454 : if (k == N-1) break;
1387 13384 : az = diviuuexact(muluui((N-k)<<1,N+k,az), k+1, (k<<1)+1);
1388 : }
1389 70 : return gerepileupto(av, gdiv(s,d));
1390 : }
1391 :
1392 : GEN
1393 14 : sumpos2(void *E, GEN (*eval)(void *, GEN), GEN a, long prec)
1394 : {
1395 : ulong k, N;
1396 14 : pari_sp av = avma;
1397 : GEN s, pol, dn, S;
1398 :
1399 14 : if (typ(a) != t_INT) pari_err_TYPE("sumpos2",a);
1400 14 : a = subiu(a,1);
1401 14 : N = (ulong)(0.31*(prec + 5));
1402 :
1403 14 : if (odd(N)) N++; /* extra precision for free */
1404 14 : S = sumpos_init(E, eval, a, N, prec);
1405 14 : pol = ZX_div_by_X_1(polzag1(N,N>>1), &dn);
1406 14 : s = gen_0;
1407 4466 : for (k=0; k<N; k++)
1408 : {
1409 4452 : GEN t = mulri(gel(S,k+1), gel(pol,k+2));
1410 4452 : s = odd(k)? mpsub(s,t): mpadd(s,t);
1411 : }
1412 14 : return gerepileupto(av, gdiv(s,dn));
1413 : }
1414 :
1415 : GEN
1416 91 : sumpos0(GEN a, GEN code, long flag, long prec)
1417 : {
1418 91 : switch(flag)
1419 : {
1420 70 : case 0: EXPR_WRAP(code, sumpos (EXPR_ARG,a,prec));
1421 14 : case 1: EXPR_WRAP(code, sumpos2(EXPR_ARG,a,prec));
1422 7 : default: pari_err_FLAG("sumpos");
1423 : }
1424 : return NULL; /* LCOV_EXCL_LINE */
1425 : }
1426 :
1427 : /********************************************************************/
1428 : /** **/
1429 : /** SEARCH FOR REAL ZEROS of an expression **/
1430 : /** **/
1431 : /********************************************************************/
1432 : /* Brent's method, [a,b] bracketing interval */
1433 : GEN
1434 23422 : zbrent(void *E, GEN (*eval)(void *, GEN), GEN a, GEN b, long prec)
1435 : {
1436 : long sig, iter, itmax, bit, bit0;
1437 23422 : pari_sp av = avma;
1438 : GEN c, d, e, fa, fb, fc;
1439 :
1440 23422 : if (typ(a) == t_INFINITY && typ(b) != t_INFINITY) swap(a,b);
1441 23422 : if (typ(a) == t_INFINITY && typ(b) == t_INFINITY)
1442 : {
1443 7 : long s = gsigne(eval(E, real_0(prec))), r = 0;
1444 7 : if (gidentical(gel(a,1), gel(b,1)))
1445 0 : pari_err_DOMAIN("solve", "a and b", "=", a, mkvec2(a, b));
1446 7 : a = real_m1(prec); /* domain = R */
1447 7 : b = real_1(prec);
1448 : for(;;)
1449 : {
1450 7 : fa = eval(E, a);
1451 7 : fb = eval(E, b);
1452 7 : if (gsigne(fa) != s)
1453 : {
1454 0 : if (r) b[1] = evalsigne(-1) | _evalexpo(r-1); else b = real_0(prec);
1455 0 : break;
1456 : }
1457 7 : if (gsigne(fb) != s)
1458 : {
1459 7 : if (r) a[1] = evalsigne(1) | _evalexpo(r-1); else a = real_0(prec);
1460 7 : break;
1461 : }
1462 0 : r++; setexpo(a, r); setexpo(b, r);
1463 : }
1464 7 : c = b;
1465 7 : goto SOLVE;
1466 : }
1467 23415 : if (typ(b) == t_INFINITY)
1468 : { /* a real, b == [+-]oo */
1469 28 : long s, r, minf = inf_get_sign(b) < 0;
1470 : GEN inc;
1471 28 : if (typ(a) != t_REAL || realprec(a) < prec) a = gtofp(a, prec);
1472 28 : fa = eval(E, a);
1473 28 : s = gsigne(fa);
1474 28 : inc = minf ? real_m1(prec) : real_1(prec);
1475 28 : r = gsigne(a) ? expo(a) : 0;
1476 : for(;;)
1477 : {
1478 570 : setexpo(inc, r);
1479 570 : b = addrr(a, inc); fb = eval(E, b);
1480 556 : if (gsigne(fb) != s) break;
1481 542 : a = b; fa = fb; r++;
1482 : }
1483 14 : if (minf) { c = a; swap(a, b); swap(fa, fb);} else c = b;
1484 14 : goto SOLVE;
1485 : }
1486 23387 : if (typ(a) != t_REAL || realprec(a) < prec) a = gtofp(a, prec);
1487 23387 : if (typ(b) != t_REAL || realprec(b) < prec) b = gtofp(b, prec);
1488 23387 : sig = cmprr(b, a);
1489 23387 : if (!sig) return gerepileupto(av, a);
1490 23387 : if (sig < 0) { c = a; swap(a, b); } else c = b;
1491 23387 : fa = eval(E, a);
1492 23387 : fb = eval(E, b);
1493 23387 : if (gsigne(fa)*gsigne(fb) > 0)
1494 7 : pari_err_DOMAIN("solve", "f(a)f(b)", ">", gen_0, mkvec2(fa, fb));
1495 23380 : SOLVE:
1496 23401 : bit0 = -prec; bit = 3+bit0; itmax = 1 - 2*bit0;
1497 23401 : fc = fb; e = d = NULL;
1498 225812 : for (iter = 1; iter <= itmax; ++iter)
1499 : { /* b = current best guess, a and c auxiliary points */
1500 : long bit2, exb;
1501 : GEN m;
1502 225812 : if (gsigne(fb)*gsigne(fc) > 0) { c = a; fc = fa; e = d = subrr(b, a); }
1503 225812 : if (gcmp(gabs(fc, 0), gabs(fb, 0)) < 0)
1504 50992 : { a = b; b = c; c = a; fa = fb; fb = fc; fc = fa; }
1505 225812 : m = subrr(c, b); shiftr_inplace(m, -1);
1506 225812 : exb = expo(b);
1507 225812 : if (bit < exb)
1508 : {
1509 225777 : bit2 = bit + exb - 1;
1510 225777 : if (expo(m) <= exb + bit0) break; /*SUCCESS*/
1511 : }
1512 : else
1513 : { /* b ~ 0 */
1514 35 : bit2 = 2*bit - 1;
1515 35 : if (expo(m) <= bit2) break; /*SUCCESS*/
1516 : }
1517 202476 : if (gequal0(fb)) break; /*SUCCESS*/
1518 :
1519 202411 : if (expo(e) > bit2 && gexpo(fa) > gexpo(fb))
1520 159658 : { /* interpolation, m != 0, |f(c)| >= |f(b)|, f(b)f(c) < 0 */
1521 159658 : GEN min1, min2, p, q, s = gdiv(fb, fa);
1522 159658 : if (a == c || equalrr(a,c))
1523 : {
1524 134067 : p = gmul2n(gmul(m, s), 1);
1525 134067 : q = gsubsg(1, s);
1526 : }
1527 : else
1528 : {
1529 25591 : GEN r = gdiv(fb, fc);
1530 25591 : q = gdiv(fa, fc);
1531 25591 : p = gmul2n(gmul(gsub(q, r), gmul(m, q)), 1);
1532 25591 : p = gmul(s, gsub(p, gmul(subrr(b, a), gsubgs(r, 1))));
1533 25591 : q = gmul(gmul(gsubgs(q, 1), gsubgs(r, 1)), gsubgs(s, 1));
1534 : }
1535 159658 : if (gsigne(p) > 0) q = gneg_i(q); else p = gneg_i(p);
1536 159658 : min1 = gsub(gmulsg(3, gmul(m,q)), gmul2n(gabs(q,0), bit2));
1537 159658 : min2 = gabs(gmul(e, q), 0);
1538 159658 : if (gcmp(gmul2n(p, 1), gmin_shallow(min1, min2)) < 0)
1539 158605 : { e = d; d = gdiv(p, q); } /* interpolation OK */
1540 : else
1541 1053 : e = d = m; /* failed, use bisection */
1542 : }
1543 42753 : else e = d = m; /* bound decreasing too slowly, use bisection */
1544 202411 : a = b; fa = fb;
1545 202411 : if (d == m) { b = addrr(c, b); shiftr_inplace(b,-1); }
1546 158605 : else if (gexpo(d) > bit2) b = gadd(b, d);
1547 24021 : else if (gsigne(m) > 0) b = addrr(b, real2n(bit2, LOWDEFAULTPREC));
1548 11859 : else b = subrr(b, real2n(bit2, LOWDEFAULTPREC));
1549 202411 : if (equalrr(a, b)) break;
1550 202411 : if (realprec(b) < prec) b = rtor(b, prec);
1551 202411 : fb = eval(E, b);
1552 : }
1553 23401 : if (iter > itmax) pari_err_IMPL("solve recovery [too many iterations]");
1554 23401 : return gerepileuptoleaf(av, rcopy(b));
1555 : }
1556 :
1557 : GEN
1558 77 : zbrent0(GEN a, GEN b, GEN code, long prec)
1559 77 : { EXPR_WRAP(code, zbrent(EXPR_ARG, a, b, prec)); }
1560 :
1561 : /* Find zeros of a function in the real interval [a,b] by interval splitting */
1562 : GEN
1563 119 : solvestep(void *E, GEN (*f)(void *,GEN), GEN a, GEN b, GEN step, long flag, long prec)
1564 : {
1565 119 : const long ITMAX = 10;
1566 119 : pari_sp av = avma;
1567 : GEN fa, a0, b0;
1568 119 : long sa0, it, bit = prec / 2, ct = 0, s = gcmp(a,b);
1569 :
1570 119 : if (!s) return gequal0(f(E, a)) ? gcopy(mkvec(a)): cgetg(1,t_VEC);
1571 119 : if (s > 0) swap(a, b);
1572 119 : if (flag&4)
1573 : {
1574 84 : if (gcmpgs(step,1)<=0) pari_err_DOMAIN("solvestep","step","<=",gen_1,step);
1575 84 : if (gsigne(a) <= 0) pari_err_DOMAIN("solvestep","a","<=",gen_0,a);
1576 : }
1577 35 : else if (gsigne(step) <= 0)
1578 7 : pari_err_DOMAIN("solvestep","step","<=",gen_0,step);
1579 112 : a0 = a = gtofp(a, prec); fa = f(E, a);
1580 112 : b0 = b = gtofp(b, prec); step = gtofp(step, prec);
1581 112 : sa0 = gsigne(fa);
1582 112 : if (gexpo(fa) < -bit) sa0 = 0;
1583 119 : for (it = 0; it < ITMAX; it++)
1584 : {
1585 119 : pari_sp av2 = avma;
1586 119 : GEN v = cgetg(1, t_VEC);
1587 119 : long sa = sa0;
1588 119 : a = a0; b = b0;
1589 37520 : while (gcmp(a,b) < 0)
1590 : {
1591 37401 : GEN fc, c = (flag&4)? gmul(a, step): gadd(a, step);
1592 : long sc;
1593 37401 : if (gcmp(c,b) > 0) c = b;
1594 37401 : fc = f(E, c); sc = gsigne(fc);
1595 37401 : if (gexpo(fc) < -bit) sc = 0;
1596 37401 : if (!sc || sa*sc < 0)
1597 : {
1598 22813 : GEN z = sc? zbrent(E, f, a, c, prec): c;
1599 : long e;
1600 22813 : (void)grndtoi(z, &e);
1601 22813 : if (e <= -bit) ct = 1;
1602 22813 : if ((flag&1) && ((!(flag&8)) || ct)) return gerepileupto(av, z);
1603 22813 : v = shallowconcat(v, z);
1604 : }
1605 37401 : a = c; fa = fc; sa = sc;
1606 37401 : if (gc_needed(av2,1))
1607 : {
1608 65 : if (DEBUGMEM>1) pari_warn(warnmem,"solvestep");
1609 65 : gerepileall(av2, 4, &a ,&fa, &v, &step);
1610 : }
1611 : }
1612 119 : if ((!(flag&2) || lg(v) > 1) && (!(flag&8) || ct))
1613 112 : return gerepilecopy(av, v);
1614 7 : step = (flag&4)? sqrtnr(step,4): gmul2n(step, -2);
1615 7 : gerepileall(av2, 2, &fa, &step);
1616 : }
1617 0 : pari_err_IMPL("solvestep recovery [too many iterations]");
1618 : return NULL;/*LCOV_EXCL_LINE*/
1619 : }
1620 :
1621 : GEN
1622 35 : solvestep0(GEN a, GEN b, GEN step, GEN code, long flag, long prec)
1623 35 : { EXPR_WRAP(code, solvestep(EXPR_ARG, a,b, step, flag, prec)); }
1624 :
1625 : /********************************************************************/
1626 : /** Numerical derivation **/
1627 : /********************************************************************/
1628 :
1629 : struct deriv_data
1630 : {
1631 : GEN code;
1632 : GEN args;
1633 : GEN def;
1634 : };
1635 :
1636 322 : static GEN deriv_eval(void *E, GEN x, long prec)
1637 : {
1638 322 : struct deriv_data *data=(struct deriv_data *)E;
1639 322 : gel(data->args,1)=x;
1640 322 : uel(data->def,1)=1;
1641 322 : return closure_callgenvecdefprec(data->code, data->args, data->def, prec);
1642 : }
1643 :
1644 : /* Rationale: (f(2^-e) - f(-2^-e) + O(2^-b)) / (2 * 2^-e) = f'(0) + O(2^-2e)
1645 : * since 2nd derivatives cancel.
1646 : * prec(LHS) = b - e
1647 : * prec(RHS) = 2e, equal when b = 3e = 3/2 b0 (b0 = required final bitprec)
1648 : *
1649 : * For f'(x), x far from 0: prec(LHS) = b - e - expo(x)
1650 : * --> pr = 3/2 b0 + expo(x) */
1651 : GEN
1652 966 : derivnum(void *E, GEN (*eval)(void *, GEN, long), GEN x, long prec)
1653 : {
1654 966 : long newprec, e, ex = gexpo(x), p = precision(x);
1655 966 : long b0 = prec2nbits(p? p: prec), b = (long)ceil(b0 * 1.5 + maxss(0,ex));
1656 : GEN eps, u, v, y;
1657 966 : pari_sp av = avma;
1658 966 : newprec = nbits2prec(b + BITS_IN_LONG);
1659 966 : switch(typ(x))
1660 : {
1661 385 : case t_REAL:
1662 : case t_COMPLEX:
1663 385 : x = gprec_w(x, newprec);
1664 : }
1665 966 : e = b0/2; /* 1/2 required prec (in sig. bits) */
1666 966 : b -= e; /* >= b0 */
1667 966 : eps = real2n(-e, ex < -e? newprec: nbits2prec(b));
1668 966 : u = eval(E, gsub(x, eps), newprec);
1669 966 : v = eval(E, gadd(x, eps), newprec);
1670 966 : y = gmul2n(gsub(v,u), e-1);
1671 966 : return gerepilecopy(av, gprec_wtrunc(y, nbits2prec(b0)));
1672 : }
1673 :
1674 : /* Fornberg interpolation algorithm for finite differences coefficients
1675 : * using 2N+1 equidistant grid points around 0 [ assume 2N even >= M ].
1676 : * Compute \delta[m]_{N,i} for all derivation orders m = 0..M such that
1677 : * h^m * f^{(m)}(0) = \sum_{i = 0}^n delta[m]_{N,i} f(a_i) + O(h^{N-m+1}),
1678 : * for step size h.
1679 : * Return a = [0,-1,1...,-N,N] and vector of vectors d: d[m+1][i+1]
1680 : * = w'(a_i) delta[m]_{2N,i}, i = 0..2N */
1681 : static void
1682 140 : FD(long M, long N2, GEN *pd, GEN *pa)
1683 : {
1684 : GEN d, a, b, W, F;
1685 140 : long N = N2>>1, m, i;
1686 :
1687 140 : F = cgetg(N2+2, t_VEC);
1688 140 : a = cgetg(N2+2, t_VEC);
1689 140 : b = cgetg(N+1, t_VEC);
1690 140 : gel(a,1) = gen_0;
1691 735 : for (i = 1; i <= N; i++)
1692 : {
1693 595 : gel(a,2*i) = utoineg(i);
1694 595 : gel(a,2*i+1) = utoipos(i);
1695 595 : gel(b,i) = sqru(i);
1696 : }
1697 : /* w = \prod (X - a[i]) = x W(x^2) */
1698 140 : W = roots_to_pol(b, 0);
1699 140 : gel(F,1) = RgX_inflate(W,2);
1700 735 : for (i = 1; i <= N; i++)
1701 : {
1702 595 : pari_sp av = avma;
1703 : GEN r, U, S;
1704 595 : U = RgX_inflate(RgX_div_by_X_x(W, gel(b,i), &r), 2);
1705 595 : U = RgXn_red_shallow(U, M); /* higher terms not needed */
1706 595 : U = RgX_shift_shallow(U,1); /* w(X) / (X^2-a[i]^2) mod X^(M+1) */
1707 595 : S = ZX_sub(RgX_shift_shallow(U,1),
1708 595 : ZX_Z_mul(U, gel(a,2*i+1)));
1709 595 : S = gerepileupto(av, S);
1710 595 : gel(F,2*i) = S;
1711 595 : gel(F,2*i+1) = ZX_z_unscale(S, -1);
1712 : }
1713 : /* F[i] = w(X) / (X-a[i]) + O(X^(M+1)) in Z[X] */
1714 140 : d = cgetg(M+2, t_VEC);
1715 693 : for (m = 0; m <= M; m++)
1716 : {
1717 553 : GEN v = cgetg(N2+2, t_VEC); /* coeff(F[i],X^m) */
1718 12222 : for (i = 0; i <= N2; i++) gel(v, i+1) = gmael(F, i+1, m+2);
1719 553 : gel(d,m+1) = v;
1720 : }
1721 140 : *pd = d;
1722 140 : *pa = a;
1723 140 : }
1724 :
1725 : static void
1726 385 : chk_ord(long m)
1727 : {
1728 385 : if (m < 0)
1729 14 : pari_err_DOMAIN("derivnumk", "derivation order", "<", gen_0, stoi(m));
1730 371 : }
1731 : /* m! / N! for m in ind; vecmax(ind) <= N. Result not a GEN if ind contains 0. */
1732 : static GEN
1733 140 : vfact(GEN ind, long N, long prec)
1734 : {
1735 : GEN v, iN;
1736 : long i, l;
1737 140 : ind = vecsmall_uniq(ind); chk_ord(ind[1]); l = lg(ind);
1738 133 : iN = invr(itor(mulu_interval(ind[1] + 1, N), prec));
1739 133 : v = const_vec(ind[l-1], NULL); gel(v, ind[1]) = iN;
1740 224 : for (i = 2; i < l; i++)
1741 91 : gel(v, ind[i]) = iN = mulri(iN, mulu_interval(ind[i-1] + 1, ind[i]));
1742 133 : return v;
1743 : }
1744 :
1745 : static GEN
1746 203 : chk_ind(GEN ind, long *M)
1747 : {
1748 203 : *M = 0;
1749 203 : switch(typ(ind))
1750 : {
1751 91 : case t_INT: ind = mkvecsmall(itos(ind)); break;
1752 0 : case t_VECSMALL:
1753 0 : if (lg(ind) == 1) return NULL;
1754 0 : break;
1755 105 : case t_VEC: case t_COL:
1756 105 : if (lg(ind) == 1) return NULL;
1757 98 : if (RgV_is_ZV(ind)) { ind = ZV_to_zv(ind); break; }
1758 : /* fall through */
1759 : default:
1760 7 : pari_err_TYPE("derivnum", ind);
1761 : return NULL; /*LCOV_EXCL_LINE*/
1762 : }
1763 189 : *M = vecsmall_max(ind); chk_ord(*M); return ind;
1764 : }
1765 : GEN
1766 168 : derivnumk(void *E, GEN (*eval)(void *, GEN, long), GEN x, GEN ind0, long prec)
1767 : {
1768 : GEN A, C, D, DM, T, X, F, v, ind, t;
1769 : long M, N, N2, fpr, p, i, pr, l, lA, e, ex, emin, emax, newprec;
1770 168 : pari_sp av = avma;
1771 168 : int allodd = 1;
1772 :
1773 168 : ind = chk_ind(ind0, &M); if (!ind) return cgetg(1, t_VEC);
1774 154 : l = lg(ind); F = cgetg(l, t_VEC);
1775 154 : if (!M) /* silly degenerate case */
1776 : {
1777 14 : X = eval(E, x, prec);
1778 28 : for (i = 1; i < l; i++) { chk_ord(ind[i]); gel(F,i) = X; }
1779 7 : if (typ(ind0) == t_INT) F = gel(F,1);
1780 7 : return gerepilecopy(av, F);
1781 : }
1782 140 : N2 = 3*M - 1; if (odd(N2)) N2++;
1783 140 : N = N2 >> 1;
1784 140 : FD(M, N2, &D,&A); /* optimal if 'eval' uses quadratic time */
1785 140 : C = vecbinomial(N2); DM = gel(D,M);
1786 140 : T = cgetg(N2+2, t_VEC);
1787 : /* (2N)! / w'(i) = (2N)! / w'(-i) = (-1)^(N-i) binom(2*N, N-i) */
1788 140 : t = gel(C, N+1);
1789 140 : gel(T,1) = odd(N)? negi(t): t;
1790 735 : for (i = 1; i <= N; i++)
1791 : {
1792 595 : t = gel(C, N-i+1);
1793 595 : gel(T,2*i) = gel(T,2*i+1) = odd(N-i)? negi(t): t;
1794 : }
1795 140 : N = N2 >> 1; emin = LONG_MAX; emax = 0;
1796 735 : for (i = 1; i <= N; i++)
1797 : {
1798 595 : e = expi(gel(DM,i)) + expi(gel(T,i));
1799 595 : if (e < 0) continue; /* 0 */
1800 504 : if (e < emin) emin = e;
1801 280 : else if (e > emax) emax = e;
1802 : }
1803 :
1804 140 : p = precision(x);
1805 140 : fpr = p ? p: prec;
1806 140 : e = (fpr + 3*M*log2((double)M)) / (2*M);
1807 140 : ex = gexpo(real_i(x));
1808 140 : if (ex < 0) ex = 0; /* near 0 */
1809 140 : pr = (long)ceil(fpr + e * M); /* ~ 3fpr/2 */
1810 140 : newprec = nbits2prec(pr + (emax - emin) + ex + BITS_IN_LONG);
1811 140 : switch(typ(x))
1812 : {
1813 28 : case t_REAL:
1814 : case t_COMPLEX:
1815 28 : x = gprec_w(x, newprec);
1816 : }
1817 140 : lA = lg(A); X = cgetg(lA, t_VEC);
1818 189 : for (i = 1; i < l; i++)
1819 147 : if (!odd(ind[i])) { allodd = 0; break; }
1820 : /* if only odd derivation orders, the value at 0 (A[1]) is not needed */
1821 140 : gel(X, 1) = gen_0;
1822 1428 : for (i = allodd? 2: 1; i < lA; i++)
1823 : {
1824 1288 : GEN t = eval(E, gadd(x, gmul2n(gel(A,i), -e)), newprec);
1825 1288 : t = gmul(t, gel(T,i));
1826 1288 : if (!gprecision(t))
1827 224 : t = is_scalar_t(typ(t))? gtofp(t, newprec): gmul(t, real_1(newprec));
1828 1288 : gel(X,i) = t;
1829 : }
1830 :
1831 140 : v = vfact(ind, N2, nbits2prec(fpr + 32));
1832 357 : for (i = 1; i < l; i++)
1833 : {
1834 224 : long m = ind[i];
1835 224 : GEN t = RgV_dotproduct(gel(D,m+1), X);
1836 224 : gel(F,i) = gmul(t, gmul2n(gel(v, m), e*m));
1837 : }
1838 133 : if (typ(ind0) == t_INT) F = gel(F,1);
1839 133 : return gerepilecopy(av, F);
1840 : }
1841 : /* v(t') */
1842 : static long
1843 14 : rfrac_val_deriv(GEN t)
1844 : {
1845 14 : long v = varn(gel(t,2));
1846 14 : return gvaluation(deriv(t, v), pol_x(v));
1847 : }
1848 :
1849 : GEN
1850 1190 : derivfunk(void *E, GEN (*eval)(void *, GEN, long), GEN x, GEN ind0, long prec)
1851 : {
1852 : pari_sp av;
1853 : GEN ind, xp, ixp, F, G;
1854 : long i, l, vx, M;
1855 1190 : if (!ind0) return derivfun(E, eval, x, prec);
1856 203 : switch(typ(x))
1857 : {
1858 140 : case t_REAL: case t_INT: case t_FRAC: case t_COMPLEX:
1859 140 : return derivnumk(E,eval, x, ind0, prec);
1860 21 : case t_POL:
1861 21 : ind = chk_ind(ind0,&M); if (!ind) return cgetg(1,t_VEC);
1862 21 : xp = RgX_deriv(x);
1863 21 : x = RgX_to_ser(x, precdl+2 + M * (1+RgX_val(xp)));
1864 21 : break;
1865 7 : case t_RFRAC:
1866 7 : ind = chk_ind(ind0,&M); if (!ind) return cgetg(1,t_VEC);
1867 7 : x = rfrac_to_ser_i(x, precdl+2 + M * (1+rfrac_val_deriv(x)));
1868 7 : xp = derivser(x);
1869 7 : break;
1870 7 : case t_SER:
1871 7 : ind = chk_ind(ind0,&M); if (!ind) return cgetg(1,t_VEC);
1872 7 : xp = derivser(x);
1873 7 : break;
1874 28 : default: pari_err_TYPE("numerical derivation",x);
1875 : return NULL; /*LCOV_EXCL_LINE*/
1876 : }
1877 35 : av = avma; vx = varn(x);
1878 35 : ixp = M? ginv(xp): NULL;
1879 35 : F = cgetg(M+2, t_VEC);
1880 35 : gel(F,1) = eval(E, x, prec);
1881 126 : for (i = 1; i <= M; i++) gel(F,i+1) = gmul(deriv(gel(F,i),vx), ixp);
1882 35 : l = lg(ind); G = cgetg(l, t_VEC);
1883 70 : for (i = 1; i < l; i++)
1884 : {
1885 35 : long m = ind[i]; chk_ord(m);
1886 35 : gel(G,i) = gel(F,m+1);
1887 : }
1888 35 : if (typ(ind0) == t_INT) G = gel(G,1);
1889 35 : return gerepilecopy(av, G);
1890 : }
1891 :
1892 : GEN
1893 987 : derivfun(void *E, GEN (*eval)(void *, GEN, long), GEN x, long prec)
1894 : {
1895 987 : pari_sp av = avma;
1896 : GEN xp;
1897 : long vx;
1898 987 : switch(typ(x))
1899 : {
1900 966 : case t_REAL: case t_INT: case t_FRAC: case t_COMPLEX:
1901 966 : return derivnum(E,eval, x, prec);
1902 7 : case t_POL:
1903 7 : xp = RgX_deriv(x);
1904 7 : x = RgX_to_ser(x, precdl+2+ (1 + RgX_val(xp)));
1905 7 : break;
1906 7 : case t_RFRAC:
1907 7 : x = rfrac_to_ser_i(x, precdl+2+ (1 + rfrac_val_deriv(x)));
1908 : /* fall through */
1909 14 : case t_SER:
1910 14 : xp = derivser(x);
1911 14 : break;
1912 0 : default: pari_err_TYPE("formal derivation",x);
1913 : return NULL; /*LCOV_EXCL_LINE*/
1914 : }
1915 21 : vx = varn(x);
1916 21 : return gerepileupto(av, gdiv(deriv(eval(E, x, prec),vx), xp));
1917 : }
1918 :
1919 : GEN
1920 21 : laurentseries(void *E, GEN (*f)(void*,GEN x, long), long M, long v, long prec)
1921 : {
1922 21 : pari_sp av = avma;
1923 : long d;
1924 :
1925 21 : if (v < 0) v = 0;
1926 21 : d = maxss(M+1,1);
1927 : for (;;)
1928 14 : {
1929 : long i, dr, vr;
1930 : GEN s;
1931 35 : s = cgetg(d+2, t_SER); s[1] = evalsigne(1) | evalvalser(1) | evalvarn(v);
1932 245 : gel(s, 2) = gen_1; for (i = 3; i <= d+1; i++) gel(s, i) = gen_0;
1933 35 : s = f(E, s, prec);
1934 35 : if (typ(s) != t_SER || varn(s) != v) pari_err_TYPE("laurentseries", s);
1935 35 : vr = valser(s);
1936 35 : if (M < vr) { set_avma(av); return zeroser(v, M); }
1937 35 : dr = lg(s) + vr - 3 - M;
1938 35 : if (dr >= 0) return gerepileupto(av, s);
1939 14 : set_avma(av); d -= dr;
1940 : }
1941 : }
1942 : static GEN
1943 35 : _evalclosprec(void *E, GEN x, long prec)
1944 : {
1945 : GEN s;
1946 35 : push_localprec(prec); s = closure_callgen1((GEN)E, x);
1947 35 : pop_localprec(); return s;
1948 : }
1949 : #define CLOS_ARGPREC __E, &_evalclosprec
1950 : GEN
1951 35 : laurentseries0(GEN f, long M, long v, long prec)
1952 : {
1953 35 : if (typ(f) != t_CLOSURE || closure_arity(f) != 1 || closure_is_variadic(f))
1954 14 : pari_err_TYPE("laurentseries",f);
1955 21 : EXPR_WRAP(f, laurentseries(CLOS_ARGPREC,M,v,prec));
1956 : }
1957 :
1958 : GEN
1959 1085 : derivnum0(GEN a, GEN code, GEN ind, long prec)
1960 1085 : { EXPR_WRAP(code, derivfunk(EXPR_ARGPREC,a,ind,prec)); }
1961 :
1962 : GEN
1963 105 : derivfun0(GEN args, GEN def, GEN code, long k, long prec)
1964 : {
1965 105 : pari_sp av = avma;
1966 : struct deriv_data E;
1967 : GEN z;
1968 105 : E.code=code; E.args=args; E.def=def;
1969 105 : z = gel(derivfunk((void*)&E, deriv_eval, gel(args,1), mkvecs(k), prec),1);
1970 77 : return gerepilecopy(av, z);
1971 : }
1972 :
1973 : /********************************************************************/
1974 : /** Numerical extrapolation **/
1975 : /********************************************************************/
1976 : /* [u(n), u <= N] */
1977 : static GEN
1978 140 : get_u(void *E, GEN (*f)(void *, GEN, long), long N, long prec)
1979 : {
1980 : long n;
1981 : GEN u;
1982 140 : if (f)
1983 : {
1984 126 : GEN v = f(E, utoipos(N), prec);
1985 126 : u = cgetg(N+1, t_VEC);
1986 126 : if (typ(v) != t_VEC || lg(v) != N+1) { gel(u,N) = v; v = NULL; }
1987 : else
1988 : {
1989 14 : GEN w = f(E, gen_1, LOWDEFAULTPREC);
1990 14 : if (typ(w) != t_VEC || lg(w) != 2) { gel(u,N) = v; v = NULL; }
1991 : }
1992 126 : if (v) u = v;
1993 : else
1994 9702 : for (n = 1; n < N; n++) gel(u,n) = f(E, utoipos(n), prec);
1995 : }
1996 : else
1997 : {
1998 14 : GEN v = (GEN)E;
1999 14 : long t = lg(v)-1;
2000 14 : if (t < N) pari_err_COMPONENT("limitnum","<",stoi(N), stoi(t));
2001 14 : u = vecslice(v, 1, N);
2002 : }
2003 12236 : for (n = 1; n <= N; n++)
2004 : {
2005 12096 : GEN un = gel(u,n);
2006 12096 : if (is_rational_t(typ(un))) gel(u,n) = gtofp(un, prec);
2007 : }
2008 140 : return u;
2009 : }
2010 :
2011 : struct limit
2012 : {
2013 : long prec; /* working accuracy */
2014 : long N; /* number of terms */
2015 : GEN na; /* [n^alpha, n <= N] */
2016 : GEN coef; /* or NULL (alpha != 1) */
2017 : };
2018 :
2019 : static GEN
2020 20822 : _gi(void *E, GEN x)
2021 : {
2022 20822 : GEN A = (GEN)E, y = gsubgs(x, 1);
2023 20822 : if (gequal0(y)) return A;
2024 20808 : return gdiv(gsubgs(gpow(x, A, LOWDEFAULTPREC), 1), y);
2025 : }
2026 : static GEN
2027 166 : _g(void *E, GEN x)
2028 : {
2029 166 : GEN D = (GEN)E, A = gel(D,1), T = gel(D,2);
2030 166 : const long prec = LOWDEFAULTPREC;
2031 166 : return gadd(glog(x,prec), intnum((void*)A, _gi, gen_0, gaddgs(x,1), T, prec));
2032 : }
2033 :
2034 : /* solve log(b) + int_0^{b+1} (x^(1/a)-1) / (x-1) dx = 0, b in [0,1]
2035 : * return -log_2(b), rounded up */
2036 : static double
2037 140 : get_accu(GEN a)
2038 : {
2039 140 : pari_sp av = avma;
2040 140 : const long prec = LOWDEFAULTPREC;
2041 140 : const double We2 = 1.844434455794; /* (W(1/e) + 1) / log(2) */
2042 : GEN b, T;
2043 140 : if (!a) return We2;
2044 49 : if (typ(a) == t_INT) switch(itos_or_0(a))
2045 : {
2046 0 : case 1: return We2;
2047 21 : case 2: return 1.186955309668;
2048 0 : case 3: return 0.883182331990;
2049 : }
2050 28 : else if (typ(a) == t_FRAC && equali1(gel(a,1))) switch(itos_or_0(gel(a,2)))
2051 : {
2052 14 : case 2: return 2.644090500290;
2053 0 : case 3: return 3.157759214459;
2054 0 : case 4: return 3.536383237500;
2055 : }
2056 14 : T = intnuminit(gen_0, gen_1, 0, prec);
2057 14 : b = zbrent((void*)mkvec2(ginv(a), T), &_g, dbltor(1E-5), gen_1, prec);
2058 14 : return gc_double(av, -dbllog2r(b));
2059 : }
2060 :
2061 : static double
2062 147 : get_c(GEN a)
2063 : {
2064 147 : double A = a? gtodouble(a): 1.0;
2065 147 : if (A <= 0) pari_err_DOMAIN("limitnum","alpha","<=",gen_0, a);
2066 140 : if (A >= 2) return 0.2270;
2067 105 : if (A >= 1) return 0.3318;
2068 14 : if (A >= 0.5) return 0.6212;
2069 0 : if (A >= 0.3333) return 1.2;
2070 0 : return 3; /* only tested for A >= 0.25 */
2071 : }
2072 : static void
2073 133 : limit_Nprec(struct limit *L, GEN alpha, long prec)
2074 : {
2075 133 : L->N = ceil(get_c(alpha) * prec);
2076 126 : L->prec = nbits2prec(prec + (long)ceil(get_accu(alpha) * L->N));
2077 126 : }
2078 : /* solve x - a log(x) = b; a, b >= 3 */
2079 : static double
2080 14 : solvedivlog(double a, double b) { return dbllemma526(a,1,1,b); }
2081 :
2082 : /* #u > 1, prod_{j != i} u[i] - u[j] */
2083 : static GEN
2084 3003 : proddiff(GEN u, long i)
2085 : {
2086 3003 : pari_sp av = avma;
2087 3003 : long l = lg(u), j;
2088 3003 : GEN p = NULL;
2089 3003 : if (i == 1)
2090 : {
2091 28 : p = gsub(gel(u,1), gel(u,2));
2092 2975 : for (j = 3; j < l; j++)
2093 2947 : p = gmul(p, gsub(gel(u,i), gel(u,j)));
2094 : }
2095 : else
2096 : {
2097 2975 : p = gsub(gel(u,i), gel(u,1));
2098 367346 : for (j = 2; j < l; j++)
2099 364371 : if (j != i) p = gmul(p, gsub(gel(u,i), gel(u,j)));
2100 : }
2101 3003 : return gerepileupto(av, p);
2102 : }
2103 : static GEN
2104 1883 : vecpows(GEN x, long N) { pari_APPLY_same(gpowgs(gel(x,i), N)); }
2105 :
2106 : static void
2107 140 : limit_init(struct limit *L, GEN alpha, int asymp)
2108 : {
2109 140 : long n, N = L->N, a = 0;
2110 140 : GEN c, v, T = NULL;
2111 :
2112 140 : if (!alpha) a = 1;
2113 49 : else if (typ(alpha) == t_INT)
2114 : {
2115 21 : a = itos_or_0(alpha);
2116 21 : if (a > 2) a = 0;
2117 : }
2118 28 : else if (typ(alpha) == t_FRAC)
2119 : {
2120 14 : long na = itos_or_0(gel(alpha,1)), da = itos_or_0(gel(alpha,2));
2121 14 : if (da && na && da <= 4 && na <= 4)
2122 : { /* don't bother with other cases */
2123 14 : long e = (N-1) % da, k = (N-1) / da;
2124 14 : if (e) { N += da - e; k++; } /* N = 1 (mod d) => simpler ^ (n/d)(N-1) */
2125 14 : L->N = N;
2126 14 : T = vecpowuu(N, na * k);
2127 : }
2128 : }
2129 140 : L->coef = v = cgetg(N+1, t_VEC);
2130 140 : if (!a)
2131 : {
2132 28 : long prec2 = gprecision(alpha);
2133 : GEN u;
2134 28 : if (prec2 && prec2 < L->prec) alpha = gprec_w(alpha, L->prec);
2135 28 : L->na = u = vecpowug(N, alpha, L->prec);
2136 28 : if (!T) T = vecpows(u, N-1);
2137 3031 : for (n = 1; n <= N; n++) gel(v,n) = gdiv(gel(T,n), proddiff(u,n));
2138 28 : return;
2139 : }
2140 112 : L->na = asymp? vecpowuu(N, a): NULL;
2141 112 : c = mpfactr(N-1, L->prec);
2142 112 : if (a == 1)
2143 : {
2144 91 : c = invr(c);
2145 91 : gel(v,1) = c; if (!odd(N)) togglesign(c);
2146 7651 : for (n = 2; n <= N; n++) gel(v,n) = divru(mulrs(gel(v,n-1), n-1-N), n);
2147 : }
2148 : else
2149 : { /* a = 2 */
2150 21 : c = invr(mulru(sqrr(c), (N*(N+1)) >> 1));
2151 21 : gel(v,1) = c; if (!odd(N)) togglesign(c);
2152 1442 : for (n = 2; n <= N; n++) gel(v,n) = divru(mulrs(gel(v,n-1), n-1-N), N+n);
2153 : }
2154 112 : T = vecpowuu(N, a*N);
2155 9093 : for (n = 2; n <= N; n++) gel(v,n) = mulri(gel(v,n), gel(T,n));
2156 : }
2157 :
2158 : /* Zagier/Lagrange extrapolation */
2159 : static GEN
2160 983 : limitnum_i(struct limit *L, GEN u, long prec)
2161 983 : { return gprec_w(RgV_dotproduct(u,L->coef), prec); }
2162 : GEN
2163 84 : limitnum(void *E, GEN (*f)(void *, GEN, long), GEN alpha, long prec)
2164 : {
2165 84 : pari_sp av = avma;
2166 : struct limit L;
2167 : GEN u;
2168 84 : limit_Nprec(&L, alpha, prec);
2169 77 : limit_init(&L, alpha, 0);
2170 77 : u = get_u(E, f, L.N, L.prec);
2171 77 : return gerepilecopy(av, limitnum_i(&L, u, prec));
2172 : }
2173 : typedef GEN (*LIMIT_FUN)(void*,GEN,long);
2174 161 : static LIMIT_FUN get_fun(GEN u, const char *s)
2175 : {
2176 161 : switch(typ(u))
2177 : {
2178 14 : case t_COL: case t_VEC: break;
2179 133 : case t_CLOSURE: return gp_callprec;
2180 14 : default: pari_err_TYPE(s, u);
2181 : }
2182 14 : return NULL;
2183 : }
2184 : GEN
2185 91 : limitnum0(GEN u, GEN alpha, long prec)
2186 91 : { return limitnum((void*)u,get_fun(u, "limitnum"), alpha, prec); }
2187 :
2188 : GEN
2189 49 : asympnum(void *E, GEN (*f)(void *, GEN, long), GEN alpha, long prec)
2190 : {
2191 49 : const long MAX = 100;
2192 49 : pari_sp av = avma;
2193 49 : GEN u, A = cgetg(MAX+1, t_VEC);
2194 49 : long i, B = prec2nbits(prec);
2195 49 : double LB = 0.9*expu(B); /* 0.9 and 0.95 below are heuristic */
2196 : struct limit L;
2197 49 : limit_Nprec(&L, alpha, prec);
2198 49 : if (alpha) LB *= gtodouble(alpha);
2199 49 : limit_init(&L, alpha, 1);
2200 49 : u = get_u(E, f, L.N, L.prec);
2201 906 : for(i = 1; i <= MAX; i++)
2202 : {
2203 906 : GEN a, v, q, s = limitnum_i(&L, u, prec);
2204 : long n;
2205 : /* NOT bestappr: lindep properly ignores the lower bits */
2206 906 : v = lindep_bit(mkvec2(gen_1, s), maxss((long)(0.95*floor(B - i*LB)), 32));
2207 906 : if (lg(v) == 1) break;
2208 899 : q = gel(v,2); if (!signe(q)) break;
2209 899 : a = gdiv(negi(gel(v,1)), q);
2210 899 : s = gsub(s, a);
2211 : /* |s|q^2 > eps */
2212 899 : if (!gequal0(s) && gexpo(s) + 2*expi(q) > -17) break;
2213 857 : gel(A,i) = a;
2214 82293 : for (n = 1; n <= L.N; n++) gel(u,n) = gmul(gsub(gel(u,n), a), gel(L.na,n));
2215 : }
2216 49 : setlg(A,i); return gerepilecopy(av, A);
2217 : }
2218 : GEN
2219 14 : asympnumraw(void *E, GEN (*f)(void *,GEN,long), long LIM, GEN alpha, long prec)
2220 : {
2221 14 : pari_sp av = avma;
2222 : double c, d, al;
2223 : long i, B;
2224 : GEN u, A;
2225 : struct limit L;
2226 :
2227 14 : if (LIM < 0) return cgetg(1, t_VEC);
2228 14 : c = get_c(alpha);
2229 14 : d = get_accu(alpha);
2230 14 : al = alpha? gtodouble(alpha): 1.0;
2231 14 : B = prec2nbits(prec);
2232 14 : L.N = ceil(solvedivlog(c * al * LIM / M_LN2, c * B));
2233 14 : L.prec = nbits2prec(ceil(B + L.N / c + d * L.N));
2234 14 : limit_init(&L, alpha, 1);
2235 14 : u = get_u(E, f, L.N, L.prec);
2236 14 : A = cgetg(LIM+2, t_VEC);
2237 217 : for(i = 0; i <= LIM; i++)
2238 : {
2239 203 : GEN a = RgV_dotproduct(u,L.coef);
2240 : long n;
2241 34461 : for (n = 1; n <= L.N; n++)
2242 34258 : gel(u,n) = gprec_wensure(gmul(gsub(gel(u,n), a), gel(L.na,n)), L.prec);
2243 203 : gel(A,i+1) = gprec_wtrunc(a, prec);
2244 : }
2245 14 : return gerepilecopy(av, A);
2246 : }
2247 : GEN
2248 56 : asympnum0(GEN u, GEN alpha, long prec)
2249 56 : { return asympnum((void*)u,get_fun(u, "asympnum"), alpha, prec); }
2250 : GEN
2251 14 : asympnumraw0(GEN u, long LIM, GEN alpha, long prec)
2252 14 : { return asympnumraw((void*)u,get_fun(u, "asympnumraw"), LIM, alpha, prec); }
|