Line data Source code
1 : #line 2 "../src/kernel/none/mp.c"
2 : /* Copyright (C) 2000-2003 The PARI group.
3 :
4 : This file is part of the PARI/GP package.
5 :
6 : PARI/GP is free software; you can redistribute it and/or modify it under the
7 : terms of the GNU General Public License as published by the Free Software
8 : Foundation; either version 2 of the License, or (at your option) any later
9 : version. It is distributed in the hope that it will be useful, but WITHOUT
10 : ANY WARRANTY WHATSOEVER.
11 :
12 : Check the License for details. You should have received a copy of it, along
13 : with the package; see the file 'COPYING'. If not, write to the Free Software
14 : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
15 :
16 : /***********************************************************************/
17 : /** **/
18 : /** MULTIPRECISION KERNEL **/
19 : /** **/
20 : /***********************************************************************/
21 : #include "pari.h"
22 : #include "paripriv.h"
23 : #include "../src/kernel/none/tune-gen.h"
24 :
25 : void
26 788 : pari_kernel_init(void) { }
27 : void
28 786 : pari_kernel_close(void) { }
29 : const char *
30 2 : pari_kernel_version(void) { return ""; }
31 :
32 : /* NOTE: arguments of "spec" routines (muliispec, addiispec, etc.) aren't
33 : * GENs but pairs (long *a, long na) representing a list of digits (in basis
34 : * BITS_IN_LONG) : a[0], ..., a[na-1]. In order to facilitate splitting: no
35 : * need to reintroduce codewords. */
36 :
37 : #define LIMBS(x) ((x)+2)
38 : #define NLIMBS(x) (lgefint(x)-2)
39 :
40 : /* Normalize a nonnegative integer */
41 : GEN
42 851292963 : int_normalize(GEN x, long known_zero_words)
43 : {
44 851292963 : long i, lx = lgefint(x);
45 : GEN x0;
46 851292963 : if (lx == 2) { x[1] = evalsigne(0) | evallgefint(2); return x; }
47 851292963 : if (!known_zero_words && x[2]) return x;
48 3533401425 : for (i = 2+known_zero_words; i < lx; i++)
49 3465822531 : if (x[i]) break;
50 332086341 : x0 = x; i -= 2; x += i;
51 332086341 : if (x0 == (GEN)avma) set_avma((pari_sp)x);
52 199207680 : else stackdummy((pari_sp)(x0+i), (pari_sp)x0);
53 332086341 : lx -= i;
54 332086341 : x[0] = evaltyp(t_INT) | evallg(lx);
55 332086341 : if (lx == 2) x[1] = evalsigne(0) | evallgefint(lx);
56 264507447 : else x[1] = evalsigne(1) | evallgefint(lx);
57 332086341 : return x;
58 : }
59 :
60 : /***********************************************************************/
61 : /** **/
62 : /** ADDITION / SUBTRACTION **/
63 : /** **/
64 : /***********************************************************************/
65 :
66 : GEN
67 2247126 : setloop(GEN a)
68 : {
69 2247126 : pari_sp av = avma;
70 2247126 : (void)cgetg(lgefint(a) + 3, t_VECSMALL);
71 2247126 : return icopy_avma(a, av); /* two cells of extra space before a */
72 : }
73 :
74 : /* we had a = setloop(?), then some incloops. Reset a to b */
75 : GEN
76 130656 : resetloop(GEN a, GEN b) {
77 130656 : long lb = lgefint(b);
78 130656 : a += lgefint(a) - lb;
79 130656 : a[0] = evaltyp(t_INT) | evallg(lb);
80 130656 : affii(b, a); return a;
81 : }
82 :
83 : /* assume a > 0, initialized by setloop. Do a++ */
84 : static GEN
85 34954581 : incpos(GEN a)
86 : {
87 34954581 : long i, l = lgefint(a);
88 34954584 : for (i=l-1; i>1; i--)
89 34954581 : if (++a[i]) return a;
90 3 : l++; a--; /* use extra cell */
91 3 : a[0]=evaltyp(t_INT) | _evallg(l);
92 3 : a[1]=evalsigne(1) | evallgefint(l);
93 3 : a[2]=1; return a;
94 : }
95 :
96 : /* assume a < 0, initialized by setloop. Do a++ */
97 : static GEN
98 50013 : incneg(GEN a)
99 : {
100 50013 : long i, l = lgefint(a)-1;
101 50013 : if (uel(a,l)--)
102 : {
103 50010 : if (l == 2 && !a[2])
104 : {
105 1485 : a++; /* save one cell */
106 1485 : a[0] = evaltyp(t_INT) | _evallg(2);
107 1485 : a[1] = evalsigne(0) | evallgefint(2);
108 : }
109 50010 : return a;
110 : }
111 3 : for (i = l-1;; i--) /* finishes since a[2] != 0 */
112 3 : if (uel(a,i)--) break;
113 3 : if (!a[2])
114 : {
115 3 : a++; /* save one cell */
116 3 : a[0] = evaltyp(t_INT) | _evallg(l);
117 3 : a[1] = evalsigne(-1) | evallgefint(l);
118 : }
119 3 : return a;
120 : }
121 :
122 : /* assume a initialized by setloop. Do a++ */
123 : GEN
124 35257065 : incloop(GEN a)
125 : {
126 35257065 : switch(signe(a))
127 : {
128 252471 : case 0: a--; /* use extra cell */
129 252471 : a[0]=evaltyp(t_INT) | _evallg(3);
130 252471 : a[1]=evalsigne(1) | evallgefint(3);
131 252471 : a[2]=1; return a;
132 50013 : case -1: return incneg(a);
133 34954581 : default: return incpos(a);
134 : }
135 : }
136 :
137 : INLINE GEN
138 2389385580 : adduispec(ulong s, GEN x, long nx)
139 : {
140 2389385580 : GEN xd, zd = (GEN)avma;
141 : long lz;
142 :
143 2389385580 : if (nx == 1) return adduu(s, uel(x,0));
144 877466025 : lz = nx+3; (void)new_chunk(lz);
145 877466025 : xd = x + nx;
146 877466025 : *--zd = (ulong)*--xd + s;
147 877466025 : if ((ulong)*zd < s)
148 : for(;;)
149 : {
150 263402418 : if (xd == x) { *--zd = 1; break; } /* enlarge z */
151 259692009 : *--zd = ((ulong)*--xd) + 1;
152 259692009 : if (*zd) { lz--; break; }
153 : }
154 620752422 : else lz--;
155 2175454917 : while (xd > x) *--zd = *--xd;
156 877466025 : *--zd = evalsigne(1) | evallgefint(lz);
157 877466025 : *--zd = evaltyp(t_INT) | evallg(lz);
158 877466025 : return gc_const((pari_sp)zd, zd);
159 : }
160 :
161 : GEN
162 491859555 : adduispec_offset(ulong s, GEN x, long offset, long nx)
163 : {
164 491859555 : GEN xd = x+lgefint(x)-nx-offset;
165 606147528 : while (nx && *xd==0) {xd++; nx--;}
166 491859555 : if (!nx) return utoi(s);
167 459733728 : return adduispec(s,xd,nx);
168 : }
169 :
170 : static GEN
171 4572364206 : addiispec(GEN x, GEN y, long nx, long ny)
172 : {
173 : GEN xd, yd, zd;
174 4572364206 : long lz, i = -2;
175 : LOCAL_OVERFLOW;
176 :
177 4572364206 : if (nx < ny) swapspec(x,y, nx,ny);
178 4572364206 : if (ny == 1) return adduispec(*y,x,nx);
179 2707839393 : zd = (GEN)avma;
180 2707839393 : lz = nx+3; (void)new_chunk(lz);
181 2707839393 : xd = x + nx;
182 2707839393 : yd = y + ny;
183 2707839393 : zd[-1] = addll(xd[-1], yd[-1]);
184 : #ifdef addllx8
185 2411298922 : for ( ; i-8 > -ny; i-=8)
186 1508685791 : addllx8(xd+i, yd+i, zd+i, overflow);
187 : #endif
188 38136247604 : for ( ; i >= -ny; i--) zd[i] = addllx(xd[i], yd[i]);
189 2707839393 : if (overflow)
190 : for(;;)
191 : {
192 550422789 : if (i < -nx) { zd[i] = 1; i--; break; } /* enlarge z */
193 367193094 : zd[i] = uel(xd,i) + 1;
194 367193094 : if (zd[i]) { i--; lz--; break; }
195 63481803 : i--;
196 : }
197 2220898407 : else lz--;
198 20542356666 : for (; i >= -nx; i--) zd[i] = xd[i];
199 2707839393 : zd += i+1;
200 2707839393 : *--zd = evalsigne(1) | evallgefint(lz);
201 2707839393 : *--zd = evaltyp(t_INT) | evallg(lz);
202 2707839393 : return gc_const((pari_sp)zd, zd);
203 : }
204 :
205 : /* assume x >= s */
206 : INLINE GEN
207 1587892914 : subiuspec(GEN x, ulong s, long nx)
208 : {
209 1587892914 : GEN xd, zd = (GEN)avma;
210 : long lz;
211 : LOCAL_OVERFLOW;
212 :
213 1587892914 : if (nx == 1) return utoi(x[0] - s);
214 :
215 364081998 : lz = nx+2; (void)new_chunk(lz);
216 364081998 : xd = x + nx;
217 364081998 : *--zd = subll(*--xd, s);
218 364081998 : if (overflow)
219 : for(;;)
220 : {
221 161370519 : *--zd = ((ulong)*--xd) - 1;
222 161370519 : if (*xd) break;
223 : }
224 364081998 : if (xd == x)
225 266834175 : while (*zd == 0) { zd++; lz--; } /* shorten z */
226 : else
227 4743060762 : do *--zd = *--xd; while (xd > x);
228 364081998 : *--zd = evalsigne(1) | evallgefint(lz);
229 364081998 : *--zd = evaltyp(t_INT) | evallg(lz);
230 364081998 : return gc_const((pari_sp)zd, zd);
231 : }
232 :
233 : /* assume x > y */
234 : static GEN
235 3382659858 : subiispec(GEN x, GEN y, long nx, long ny)
236 : {
237 : GEN xd,yd,zd;
238 3382659858 : long lz, i = -2;
239 : LOCAL_OVERFLOW;
240 :
241 3382659858 : if (ny==1) return subiuspec(x,*y,nx);
242 1954395003 : zd = (GEN)avma;
243 1954395003 : lz = nx+2; (void)new_chunk(lz);
244 1954395003 : xd = x + nx;
245 1954395003 : yd = y + ny;
246 1954395003 : zd[-1] = subll(xd[-1], yd[-1]);
247 : #ifdef subllx8
248 2184298345 : for ( ; i-8 > -ny; i-=8)
249 1532833344 : subllx8(xd+i, yd+i, zd+i, overflow);
250 : #endif
251 34142554398 : for ( ; i >= -ny; i--) zd[i] = subllx(xd[i], yd[i]);
252 1954395003 : if (overflow)
253 : for(;;)
254 : {
255 991641735 : zd[i] = uel(xd,i) - 1;
256 991641735 : if (xd[i--]) break;
257 : }
258 1954395003 : if (i>=-nx)
259 4545381054 : for (; i >= -nx; i--) zd[i] = xd[i];
260 : else
261 2344893591 : while (zd[i+1] == 0) { i++; lz--; } /* shorten z */
262 1954395003 : zd += i+1;
263 1954395003 : *--zd = evalsigne(1) | evallgefint(lz);
264 1954395003 : *--zd = evaltyp(t_INT) | evallg(lz);
265 1954395003 : return gc_const((pari_sp)zd, zd);
266 : }
267 :
268 : static void
269 444314239 : roundr_up_ip(GEN x, long l)
270 : {
271 444314239 : long i = l;
272 : for(;;)
273 : {
274 445274980 : if (++uel(x,--i)) break;
275 1287942 : if (i == 2) { x[2] = (long)HIGHBIT; shiftr_inplace(x, 1); break; }
276 : }
277 444314239 : }
278 :
279 : void
280 320525793 : affir(GEN x, GEN y)
281 : {
282 320525793 : const long s = signe(x), ly = lg(y);
283 : long lx, sh, i;
284 :
285 320525793 : if (!s)
286 : {
287 31185198 : y[1] = evalexpo(-bit_accuracy(ly));
288 31185198 : return;
289 : }
290 :
291 289340595 : lx = lgefint(x); sh = bfffo(x[2]);
292 289340595 : y[1] = evalsigne(s) | evalexpo(bit_accuracy(lx)-sh-1);
293 289340595 : if (sh) {
294 284292597 : if (lx <= ly)
295 : {
296 748939722 : for (i=lx; i<ly; i++) y[i]=0;
297 203830758 : shift_left(y,x,2,lx-1, 0,sh);
298 203830758 : return;
299 : }
300 80461839 : shift_left(y,x,2,ly-1, x[ly],sh);
301 : /* lx > ly: round properly */
302 80461839 : if ((uel(x,ly)<<sh) & HIGHBIT) roundr_up_ip(y, ly);
303 : }
304 : else {
305 5047998 : if (lx <= ly)
306 : {
307 4681614 : for (i=2; i<lx; i++) y[i]=x[i];
308 3732465 : for ( ; i<ly; i++) y[i]=0;
309 1270461 : return;
310 : }
311 9519327 : for (i=2; i<ly; i++) y[i]=x[i];
312 : /* lx > ly: round properly */
313 3777537 : if (uel(x,ly) & HIGHBIT) roundr_up_ip(y, ly);
314 : }
315 : }
316 :
317 : INLINE GEN
318 1275317085 : shiftispec(GEN x, long nx, long n)
319 : {
320 : long ny, i, m;
321 : GEN y, yd;
322 1275317085 : if (!n) return icopyspec(x, nx);
323 :
324 1188610017 : if (n > 0)
325 : {
326 731149005 : GEN z = (GEN)avma;
327 731149005 : long d = dvmdsBIL(n, &m);
328 :
329 731149005 : ny = nx+d; y = new_chunk(ny + 2); yd = y + 2;
330 6731159784 : for ( ; d; d--) *--z = 0;
331 1908144075 : if (!m) for (i=0; i<nx; i++) yd[i]=x[i];
332 : else
333 : {
334 710137992 : const ulong sh = BITS_IN_LONG - m;
335 710137992 : shift_left(yd,x, 0,nx-1, 0,m);
336 710137992 : i = uel(x,0) >> sh;
337 : /* Extend y on the left? */
338 710137992 : if (i) { ny++; y = new_chunk(1); y[2] = i; }
339 : }
340 : }
341 : else
342 : {
343 457461012 : ny = nx - dvmdsBIL(-n, &m);
344 457461012 : if (ny<1) return gen_0;
345 456198552 : y = new_chunk(ny + 2); yd = y + 2;
346 456198552 : if (m) {
347 277898688 : shift_right(yd,x, 0,ny, 0,m);
348 277898688 : if (yd[0] == 0)
349 : {
350 33237783 : if (ny==1) return gc_const((pari_sp)(y+3), gen_0);
351 25583100 : ny--; set_avma((pari_sp)(++y));
352 : }
353 : } else {
354 7431945717 : for (i=0; i<ny; i++) yd[i]=x[i];
355 : }
356 : }
357 1179692874 : y[1] = evalsigne(1)|evallgefint(ny + 2);
358 1179692874 : y[0] = evaltyp(t_INT)|evallg(ny + 2); return y;
359 : }
360 :
361 : GEN
362 48718818 : mantissa2nr(GEN x, long n)
363 : { /*This is a kludge since x is not an integer*/
364 48718818 : long s = signe(x);
365 : GEN y;
366 :
367 48718818 : if(s == 0) return gen_0;
368 48717900 : y = shiftispec(x + 2, lg(x) - 2, n);
369 48717900 : if (signe(y)) setsigne(y, s);
370 48717900 : return y;
371 : }
372 :
373 : GEN
374 2621976 : truncr(GEN x)
375 : {
376 : long d,e,i,s,m;
377 : GEN y;
378 :
379 2621976 : if ((s=signe(x)) == 0 || (e=expo(x)) < 0) return gen_0;
380 1099497 : d = nbits2lg(e+1); m = remsBIL(e);
381 1099497 : if (d > lg(x)) pari_err_PREC( "truncr (precision loss in truncation)");
382 :
383 1099494 : y=cgeti(d); y[1] = evalsigne(s) | evallgefint(d);
384 1099494 : if (++m == BITS_IN_LONG)
385 834 : for (i=2; i<d; i++) y[i]=x[i];
386 : else
387 1099140 : shift_right(y,x, 2,d,0, BITS_IN_LONG - m);
388 1099494 : return y;
389 : }
390 :
391 : /* integral part */
392 : GEN
393 5246184 : floorr(GEN x)
394 : {
395 : long d,e,i,lx,m;
396 : GEN y;
397 :
398 5246184 : if (signe(x) >= 0) return truncr(x);
399 3155790 : if ((e=expo(x)) < 0) return gen_m1;
400 2634549 : d = nbits2lg(e+1); m = remsBIL(e);
401 2634549 : lx=lg(x); if (d>lx) pari_err_PREC( "floorr (precision loss in truncation)");
402 2634549 : y = new_chunk(d);
403 2634549 : if (++m == BITS_IN_LONG)
404 : {
405 501 : for (i=2; i<d; i++) y[i]=x[i];
406 210 : i=d; while (i<lx && !x[i]) i++;
407 174 : if (i==lx) goto END;
408 : }
409 : else
410 : {
411 2634375 : shift_right(y,x, 2,d,0, BITS_IN_LONG - m);
412 2634375 : if (uel(x,d-1)<<m == 0)
413 : {
414 316512 : i=d; while (i<lx && !x[i]) i++;
415 85230 : if (i==lx) goto END;
416 : }
417 : }
418 : /* set y:=y+1 */
419 2577030 : for (i=d-1; i>=2; i--) { uel(y,i)++; if (y[i]) goto END; }
420 0 : y=new_chunk(1); y[2]=1; d++;
421 2634549 : END:
422 2634549 : y[1] = evalsigne(-1) | evallgefint(d);
423 2634549 : y[0] = evaltyp(t_INT) | evallg(d); return y;
424 : }
425 :
426 : INLINE int
427 4098786696 : cmpiispec(GEN x, GEN y, long lx, long ly)
428 : {
429 : long i;
430 4098786696 : if (lx < ly) return -1;
431 3807874386 : if (lx > ly) return 1;
432 3793356909 : i = 0; while (i<lx && x[i]==y[i]) i++;
433 3294138363 : if (i==lx) return 0;
434 3067286304 : return (uel(x,i) > uel(y,i))? 1: -1;
435 : }
436 :
437 : INLINE int
438 198584547 : equaliispec(GEN x, GEN y, long lx, long ly)
439 : {
440 : long i;
441 198584547 : if (lx != ly) return 0;
442 363317703 : i = ly-1; while (i>=0 && x[i]==y[i]) i--;
443 198502095 : return i < 0;
444 : }
445 :
446 : /***********************************************************************/
447 : /** **/
448 : /** MULTIPLICATION **/
449 : /** **/
450 : /***********************************************************************/
451 : /* assume ny > 0 */
452 : INLINE GEN
453 4809036840 : muluispec(ulong x, GEN y, long ny)
454 : {
455 4809036840 : GEN yd, z = (GEN)avma;
456 4809036840 : long lz = ny+3;
457 : LOCAL_HIREMAINDER;
458 :
459 4809036840 : (void)new_chunk(lz);
460 4809036840 : yd = y + ny; *--z = mulll(x, *--yd);
461 15389186991 : while (yd > y) *--z = addmul(x,*--yd);
462 4809036840 : if (hiremainder) *--z = hiremainder; else lz--;
463 4809036840 : *--z = evalsigne(1) | evallgefint(lz);
464 4809036840 : *--z = evaltyp(t_INT) | evallg(lz);
465 4809036840 : return gc_const((pari_sp)z, z);
466 : }
467 :
468 : /* a + b*|Y| */
469 : GEN
470 0 : addumului(ulong a, ulong b, GEN Y)
471 : {
472 : GEN yd,y,z;
473 : long ny,lz;
474 : LOCAL_HIREMAINDER;
475 : LOCAL_OVERFLOW;
476 :
477 0 : if (!b || !signe(Y)) return utoi(a);
478 :
479 0 : y = LIMBS(Y); z = (GEN)avma;
480 0 : ny = NLIMBS(Y);
481 0 : lz = ny+3;
482 :
483 0 : (void)new_chunk(lz);
484 0 : yd = y + ny; *--z = addll(a, mulll(b, *--yd));
485 0 : if (overflow) hiremainder++; /* can't overflow */
486 0 : while (yd > y) *--z = addmul(b,*--yd);
487 0 : if (hiremainder) *--z = hiremainder; else lz--;
488 0 : *--z = evalsigne(1) | evallgefint(lz);
489 0 : *--z = evaltyp(t_INT) | evallg(lz);
490 0 : return gc_const((pari_sp)z, z);
491 : }
492 :
493 : /***********************************************************************/
494 : /** **/
495 : /** DIVISION **/
496 : /** **/
497 : /***********************************************************************/
498 :
499 : ulong
500 1407658350 : umodiu(GEN y, ulong x)
501 : {
502 1407658350 : long sy=signe(y),ly,i;
503 : ulong xi;
504 : LOCAL_HIREMAINDER;
505 :
506 1407658350 : if (!x) pari_err_INV("umodiu",gen_0);
507 1407658350 : if (!sy) return 0;
508 1103753172 : ly = lgefint(y);
509 1103753172 : if (x <= uel(y,2))
510 : {
511 332838717 : hiremainder=0;
512 332838717 : if (ly==3)
513 : {
514 302591793 : hiremainder=uel(y,2)%x;
515 302591793 : if (!hiremainder) return 0;
516 254063850 : return (sy > 0)? hiremainder: x - hiremainder;
517 : }
518 : }
519 : else
520 : {
521 770914455 : if (ly==3) return (sy > 0)? uel(y,2): x - uel(y,2);
522 97581375 : hiremainder=uel(y,2); ly--; y++;
523 : }
524 127828299 : xi = get_Fl_red(x);
525 905732916 : for (i=2; i<ly; i++) (void)divll_pre(y[i],x,xi);
526 127828299 : if (!hiremainder) return 0;
527 121864986 : return (sy > 0)? hiremainder: x - hiremainder;
528 : }
529 :
530 : /* return |y| \/ x */
531 : GEN
532 278806557 : absdiviu_rem(GEN y, ulong x, ulong *rem)
533 : {
534 : long ly,i;
535 : GEN z;
536 : ulong xi;
537 : LOCAL_HIREMAINDER;
538 :
539 278806557 : if (!x) pari_err_INV("absdiviu_rem",gen_0);
540 278806557 : if (!signe(y)) { *rem = 0; return gen_0; }
541 :
542 258800718 : ly = lgefint(y);
543 258800718 : if (x <= uel(y,2))
544 : {
545 231511236 : hiremainder=0;
546 231511236 : if (ly==3)
547 : {
548 208036485 : z = cgetipos(3);
549 208036485 : z[2] = divll(uel(y,2),x);
550 208036485 : *rem = hiremainder; return z;
551 : }
552 : }
553 : else
554 : {
555 27289482 : if (ly==3) { *rem = uel(y,2); return gen_0; }
556 6813600 : hiremainder = uel(y,2); ly--; y++;
557 : }
558 30288351 : xi = get_Fl_red(x);
559 30288351 : z = cgetipos(ly);
560 164540673 : for (i=2; i<ly; i++) z[i]=divll_pre(y[i],x,xi);
561 30288351 : *rem = hiremainder; return z;
562 : }
563 :
564 : GEN
565 65060814 : divis_rem(GEN y, long x, long *rem)
566 : {
567 65060814 : long sy=signe(y),ly,s,i;
568 : GEN z;
569 : ulong xi;
570 : LOCAL_HIREMAINDER;
571 :
572 65060814 : if (!x) pari_err_INV("divis_rem",gen_0);
573 65060814 : if (!sy) { *rem=0; return gen_0; }
574 45955152 : if (x<0) { s = -sy; x = -x; } else s = sy;
575 :
576 45955152 : ly = lgefint(y);
577 45955152 : if ((ulong)x <= uel(y,2))
578 : {
579 31723251 : hiremainder=0;
580 31723251 : if (ly==3)
581 : {
582 31421574 : z = cgeti(3); z[1] = evallgefint(3) | evalsigne(s);
583 31421574 : z[2] = divll(uel(y,2),x);
584 31421574 : if (sy<0) hiremainder = - ((long)hiremainder);
585 31421574 : *rem = (long)hiremainder; return z;
586 : }
587 : }
588 : else
589 : {
590 14231901 : if (ly==3) { *rem = itos(y); return gen_0; }
591 256056 : hiremainder = uel(y,2); ly--; y++;
592 : }
593 557733 : xi = get_Fl_red(x);
594 557733 : z = cgeti(ly); z[1] = evallgefint(ly) | evalsigne(s);
595 2947848 : for (i=2; i<ly; i++) z[i]=divll_pre(y[i],x,xi);
596 557733 : if (sy<0) hiremainder = - ((long)hiremainder);
597 557733 : *rem = (long)hiremainder; return z;
598 : }
599 :
600 : GEN
601 721647 : divis(GEN y, long x)
602 : {
603 721647 : long sy=signe(y),ly,s,i;
604 : ulong xi;
605 : GEN z;
606 : LOCAL_HIREMAINDER;
607 :
608 721647 : if (!x) pari_err_INV("divis",gen_0);
609 721647 : if (!sy) return gen_0;
610 721611 : if (x<0) { s = -sy; x = -x; } else s = sy;
611 :
612 721611 : ly = lgefint(y);
613 721611 : if ((ulong)x <= uel(y,2))
614 : {
615 712572 : hiremainder=0;
616 712572 : if (ly==3)
617 : {
618 643170 : z = cgeti(3); z[1] = evallgefint(3) | evalsigne(s);
619 643170 : z[2] = divll(y[2],x);
620 643170 : return z;
621 : }
622 : }
623 : else
624 : {
625 9039 : if (ly==3) return gen_0;
626 8805 : hiremainder=y[2]; ly--; y++;
627 : }
628 78207 : xi = get_Fl_red(x);
629 78207 : z = cgeti(ly); z[1] = evallgefint(ly) | evalsigne(s);
630 593091 : for (i=2; i<ly; i++) z[i]=divll_pre(y[i],x, xi);
631 78207 : return z;
632 : }
633 :
634 : GEN
635 129485079 : divrr(GEN x, GEN y)
636 : {
637 129485079 : long sx=signe(x), sy=signe(y), lx,ly,lr,e,i,j;
638 : ulong y0,y1;
639 : GEN r, r1;
640 :
641 129485079 : if (!sy) pari_err_INV("divrr",y);
642 129485079 : e = expo(x) - expo(y);
643 129485079 : if (!sx) return real_0_bit(e);
644 129127245 : if (sy<0) sx = -sx;
645 :
646 129127245 : lx=lg(x); ly=lg(y);
647 129127245 : if (ly==3)
648 : {
649 23516181 : ulong k = x[2], l = (lx>3)? x[3]: 0;
650 : LOCAL_HIREMAINDER;
651 23516181 : if (k < uel(y,2)) e--;
652 : else
653 : {
654 6883794 : l >>= 1; if (k&1) l |= HIGHBIT;
655 6883794 : k >>= 1;
656 : }
657 23516181 : hiremainder = k; k = divll(l,y[2]);
658 23516181 : if (hiremainder > (uel(y,2) >> 1) && !++k) { k = HIGHBIT; e++; }
659 23516181 : r = cgetg(3, t_REAL);
660 23516181 : r[1] = evalsigne(sx) | evalexpo(e);
661 23516181 : r[2] = k; return r;
662 : }
663 :
664 105611064 : lr = minss(lx,ly); r = new_chunk(lr);
665 105611064 : r1 = r-1;
666 748975704 : r1[1] = 0; for (i=2; i<lr; i++) r1[i]=x[i];
667 105611064 : r1[lr] = (lx>ly)? x[lr]: 0;
668 105611064 : y0 = y[2]; y1 = y[3];
669 854586768 : for (i=0; i<lr-1; i++)
670 : { /* r1 = r + (i-1), OK up to r1[2] (accesses at most r[lr]) */
671 : ulong k, qp;
672 : LOCAL_HIREMAINDER;
673 : LOCAL_OVERFLOW;
674 :
675 748975704 : if (uel(r1,1) == y0) { qp = ULONG_MAX; k = addll(y0,r1[2]); }
676 : else
677 : {
678 747399921 : if (uel(r1,1) > y0) /* can't happen if i=0 */
679 : {
680 0 : GEN y1 = y+1;
681 0 : j = lr-i; r1[j] = subll(r1[j],y1[j]);
682 0 : for (j--; j>0; j--) r1[j] = subllx(r1[j],y1[j]);
683 0 : j=i; do uel(r,--j)++; while (j && !uel(r,j));
684 : }
685 747399921 : hiremainder = r1[1]; overflow = 0;
686 747399921 : qp = divll(r1[2],y0); k = hiremainder;
687 : }
688 748975704 : j = lr-i+1;
689 748975704 : if (!overflow)
690 : {
691 : long k3, k4;
692 747722946 : k3 = mulll(qp,y1);
693 747722946 : if (j == 3) /* i = lr - 2 maximal, r1[3] undefined -> 0 */
694 105543561 : k4 = subll(hiremainder,k);
695 : else
696 : {
697 642179385 : k3 = subll(k3, r1[3]);
698 642179385 : k4 = subllx(hiremainder,k);
699 : }
700 990730706 : while (!overflow && k4) { qp--; k3 = subll(k3,y1); k4 = subllx(k4,y0); }
701 : }
702 748975704 : if (j<ly) (void)mulll(qp,y[j]); else { hiremainder = 0 ; j = ly; }
703 4962159117 : for (j--; j>1; j--)
704 : {
705 4213183413 : r1[j] = subll(r1[j], addmul(qp,y[j]));
706 4213183413 : hiremainder += overflow;
707 : }
708 748975704 : if (uel(r1,1) != hiremainder)
709 : {
710 596505 : if (uel(r1,1) < hiremainder)
711 : {
712 596505 : qp--;
713 596505 : j = lr-i-(lr-i>=ly); r1[j] = addll(r1[j], y[j]);
714 3336804 : for (j--; j>1; j--) r1[j] = addllx(r1[j], y[j]);
715 : }
716 : else
717 : {
718 0 : r1[1] -= hiremainder;
719 0 : while (r1[1])
720 : {
721 0 : qp++; if (!qp) { j=i; do uel(r,--j)++; while (j && !r[j]); }
722 0 : j = lr-i-(lr-i>=ly); r1[j] = subll(r1[j],y[j]);
723 0 : for (j--; j>1; j--) r1[j] = subllx(r1[j],y[j]);
724 0 : r1[1] -= overflow;
725 : }
726 : }
727 : }
728 748975704 : *++r1 = qp;
729 : }
730 : /* i = lr-1 */
731 : /* round correctly */
732 105611064 : if (uel(r1,1) > (y0>>1))
733 : {
734 51864660 : j=i; do uel(r,--j)++; while (j && !r[j]);
735 : }
736 748975704 : r1 = r-1; for (j=i; j>=2; j--) r[j]=r1[j];
737 105611064 : if (r[0] == 0) e--;
738 45948471 : else if (r[0] == 1) { shift_right(r,r, 2,lr, 1,1); }
739 : else { /* possible only when rounding up to 0x2 0x0 ... */
740 6 : r[2] = (long)HIGHBIT; e++;
741 : }
742 105611064 : r[0] = evaltyp(t_REAL)|evallg(lr);
743 105611064 : r[1] = evalsigne(sx) | evalexpo(e);
744 105611064 : return r;
745 : }
746 :
747 : GEN
748 113500974 : divri(GEN x, GEN y)
749 : {
750 113500974 : long lx, s = signe(y);
751 : pari_sp av;
752 : GEN z;
753 :
754 113500974 : if (!s) pari_err_INV("divri",y);
755 113500974 : if (!signe(x)) return real_0_bit(expo(x) - expi(y));
756 113331228 : if (!is_bigint(y)) {
757 89113542 : GEN z = divru(x, y[2]);
758 89113542 : if (s < 0) togglesign(z);
759 89113542 : return z;
760 : }
761 24217686 : lx = lg(x); z = cgetg(lx, t_REAL); av = avma;
762 24217686 : affrr(divrr(x, itor(y, lg2prec(lx+1))), z);
763 24217686 : return gc_const(av, z);
764 : }
765 :
766 : /* Integer division x / y: such that sign(r) = sign(x)
767 : * if z = ONLY_REM return remainder, otherwise return quotient
768 : * if z != NULL set *z to remainder
769 : * *z is the last object on stack (and can be disposed of with cgiv
770 : * If *z is zero, we put gen_0 here and no copy.
771 : * space needed: lx + ly */
772 : GEN
773 1671314826 : dvmdii(GEN x, GEN y, GEN *z)
774 : {
775 1671314826 : long sx = signe(x), sy = signe(y);
776 1671314826 : long lx, ly = lgefint(y), lz, i, j, sh, lq, lr;
777 : pari_sp av;
778 : ulong y0,y0i,y1, *xd,*rd,*qd;
779 : GEN q, r, r1;
780 :
781 1671314826 : if (!sx)
782 : {
783 52961784 : if (ly < 3) pari_err_INV("dvmdii",gen_0);
784 52961781 : if (!z || z == ONLY_REM) return gen_0;
785 32205201 : *z=gen_0; return gen_0;
786 : }
787 1618353042 : if (ly <= 3)
788 : {
789 : ulong rem;
790 651531069 : if (ly < 3) pari_err_INV("dvmdii",gen_0);
791 651531069 : if (z == ONLY_REM)
792 : {
793 445374519 : rem = umodiu(x,uel(y,2));
794 445374519 : if (!rem) return gen_0;
795 402964671 : return (sx < 0)? utoineg(uel(y,2) - rem): utoipos(rem);
796 : }
797 206156550 : q = absdiviu_rem(x, uel(y,2), &rem);
798 206156550 : if (sx != sy) togglesign(q);
799 206156550 : if (!z) return q;
800 202931553 : if (!rem) *z = gen_0;
801 57281508 : else *z = sx < 0? utoineg(rem): utoipos(rem);
802 202931553 : return q;
803 : }
804 966821973 : lx=lgefint(x);
805 966821973 : lz=lx-ly;
806 966821973 : if (lz <= 0)
807 : {
808 440366862 : if (lz == 0)
809 : {
810 334413666 : for (i=2; i<lx; i++)
811 333787131 : if (x[i] != y[i])
812 : {
813 317471697 : if (uel(x,i) > uel(y,i)) goto DIVIDE;
814 44701197 : goto TRIVIAL;
815 : }
816 626535 : if (z == ONLY_REM) return gen_0;
817 65814 : if (z) *z = gen_0;
818 65814 : if (sx < 0) sy = -sy;
819 65814 : return stoi(sy);
820 : }
821 122268630 : TRIVIAL:
822 166969827 : if (z == ONLY_REM) return icopy(x);
823 2163063 : if (z) *z = icopy(x);
824 2163063 : return gen_0;
825 : }
826 526455111 : DIVIDE: /* quotient is nonzero */
827 799225611 : av=avma; if (sx<0) sy = -sy;
828 799225611 : r1 = new_chunk(lx); sh = bfffo(y[2]);
829 799225611 : if (sh)
830 : { /* normalize so that highbit(y) = 1 (shift left x and y by sh bits)*/
831 790306080 : const ulong m = BITS_IN_LONG - sh;
832 790306080 : r = new_chunk(ly);
833 790306080 : shift_left(r, y,2,ly-1, 0,sh); y = r;
834 790306080 : shift_left(r1,x,2,lx-1, 0,sh);
835 790306080 : r1[1] = uel(x,2) >> m;
836 : }
837 : else
838 : {
839 95127831 : r1[1] = 0; for (j=2; j<lx; j++) r1[j] = x[j];
840 : }
841 799225611 : x = r1;
842 799225611 : y0 = y[2]; y0i = get_Fl_red(y0);
843 799225611 : y1 = y[3];
844 3264224913 : for (i=0; i<=lz; i++)
845 : { /* r1 = x + i */
846 : ulong k, qp;
847 : LOCAL_HIREMAINDER;
848 : LOCAL_OVERFLOW;
849 :
850 2464999302 : if (uel(r1,1) == y0)
851 : {
852 49749 : qp = ULONG_MAX; k = addll(y0,r1[2]);
853 : }
854 : else
855 : {
856 2464949553 : hiremainder = r1[1]; overflow = 0;
857 2464949553 : qp = divll_pre(r1[2],y0,y0i); k = hiremainder;
858 : }
859 2464999302 : if (!overflow)
860 : {
861 2464949145 : long k3 = subll(mulll(qp,y1), r1[3]);
862 2464949145 : long k4 = subllx(hiremainder,k);
863 3006125908 : while (!overflow && k4) { qp--; k3 = subll(k3,y1); k4 = subllx(k4,y0); }
864 : }
865 2464999302 : hiremainder = 0; j = ly;
866 63799441413 : for (j--; j>1; j--)
867 : {
868 61334442111 : r1[j] = subll(r1[j], addmul(qp,y[j]));
869 61334442111 : hiremainder += overflow;
870 : }
871 2464999302 : if (uel(r1,1) < hiremainder)
872 : {
873 5931406 : qp--;
874 5931406 : j = ly-1; r1[j] = addll(r1[j],y[j]);
875 31280652 : for (j--; j>1; j--) r1[j] = addllx(r1[j],y[j]);
876 : }
877 2464999302 : *++r1 = qp;
878 : }
879 :
880 799225611 : lq = lz+2;
881 799225611 : if (!z)
882 : {
883 2808330 : qd = (ulong*)av;
884 2808330 : xd = (ulong*)(x + lq);
885 2808330 : if (x[1]) { lz++; lq++; }
886 34968915 : while (lz--) *--qd = *--xd;
887 2808330 : *--qd = evalsigne(sy) | evallgefint(lq);
888 2808330 : *--qd = evaltyp(t_INT) | evallg(lq);
889 2808330 : return gc_const((pari_sp)qd, (GEN)qd);
890 : }
891 :
892 888438585 : j=lq; while (j<lx && !x[j]) j++;
893 796417281 : lz = lx-j;
894 796417281 : if (z == ONLY_REM)
895 : {
896 516687504 : if (lz==0) return gc_const(av, gen_0);
897 507369852 : rd = (ulong*)av; lr = lz+2;
898 507369852 : xd = (ulong*)(x + lx);
899 541081473 : if (!sh) while (lz--) *--rd = *--xd;
900 : else
901 : { /* shift remainder right by sh bits */
902 499221753 : const ulong shl = BITS_IN_LONG - sh;
903 : ulong l;
904 499221753 : xd--;
905 1497803343 : while (--lz) /* fill r[3..] */
906 : {
907 998581590 : l = *xd >> sh;
908 998581590 : *--rd = l | (*--xd << shl);
909 : }
910 499221753 : l = *xd >> sh;
911 499221753 : if (l) *--rd = l; else lr--;
912 : }
913 507369852 : *--rd = evalsigne(sx) | evallgefint(lr);
914 507369852 : *--rd = evaltyp(t_INT) | evallg(lr);
915 507369852 : return gc_const((pari_sp)rd, (GEN)rd);
916 : }
917 :
918 279729777 : lr = lz+2;
919 279729777 : rd = NULL; /* gcc -Wall */
920 279729777 : if (lz)
921 : { /* non zero remainder: initialize rd */
922 275068578 : xd = (ulong*)(x + lx);
923 275068578 : if (!sh)
924 : {
925 572136 : rd = (ulong*)avma; (void)new_chunk(lr);
926 5761080 : while (lz--) *--rd = *--xd;
927 : }
928 : else
929 : { /* shift remainder right by sh bits */
930 274496442 : const ulong shl = BITS_IN_LONG - sh;
931 : ulong l;
932 274496442 : rd = (ulong*)x; /* overwrite shifted y */
933 274496442 : xd--;
934 1223230170 : while (--lz)
935 : {
936 948733728 : l = *xd >> sh;
937 948733728 : *--rd = l | (*--xd << shl);
938 : }
939 274496442 : l = *xd >> sh;
940 274496442 : if (l) *--rd = l; else lr--;
941 : }
942 275068578 : *--rd = evalsigne(sx) | evallgefint(lr);
943 275068578 : *--rd = evaltyp(t_INT) | evallg(lr);
944 275068578 : rd += lr;
945 : }
946 279729777 : qd = (ulong*)av;
947 279729777 : xd = (ulong*)(x + lq);
948 279729777 : if (x[1]) lq++;
949 874962522 : j = lq-2; while (j--) *--qd = *--xd;
950 279729777 : *--qd = evalsigne(sy) | evallgefint(lq);
951 279729777 : *--qd = evaltyp(t_INT) | evallg(lq);
952 279729777 : q = (GEN)qd;
953 279729777 : if (lr==2) *z = gen_0;
954 : else
955 : { /* rd has been properly initialized: we had lz > 0 */
956 1873444296 : while (lr--) *--qd = *--rd;
957 275068578 : *z = (GEN)qd;
958 : }
959 279729777 : return gc_const((pari_sp)qd, q);
960 : }
961 :
962 : /* Montgomery reduction.
963 : * N has k words, assume T >= 0 has less than 2k.
964 : * Return res := T / B^k mod N, where B = 2^BIL
965 : * such that 0 <= res < T/B^k + N and res has less than k words */
966 : GEN
967 36737337 : red_montgomery(GEN T, GEN N, ulong inv)
968 : {
969 : pari_sp av;
970 : GEN Te, Td, Ne, Nd, scratch;
971 36737337 : ulong i, j, m, t, d, k = NLIMBS(N);
972 : int carry;
973 : LOCAL_HIREMAINDER;
974 : LOCAL_OVERFLOW;
975 :
976 36737337 : if (k == 0) return gen_0;
977 36737337 : d = NLIMBS(T); /* <= 2*k */
978 36737337 : if (d == 0) return gen_0;
979 : #ifdef DEBUG
980 : if (d > 2*k) pari_err_BUG("red_montgomery");
981 : #endif
982 36737328 : if (k == 1)
983 : { /* as below, special cased for efficiency */
984 163341 : ulong n = uel(N,2);
985 163341 : if (d == 1) {
986 163194 : hiremainder = uel(T,2);
987 163194 : m = hiremainder * inv;
988 163194 : (void)addmul(m, n); /* t + m*n = 0 */
989 163194 : return utoi(hiremainder);
990 : } else { /* d = 2 */
991 147 : hiremainder = uel(T,3);
992 147 : m = hiremainder * inv;
993 147 : (void)addmul(m, n); /* t + m*n = 0 */
994 147 : t = addll(hiremainder, uel(T,2));
995 147 : if (overflow) t -= n; /* t > n doesn't fit in 1 word */
996 147 : return utoi(t);
997 : }
998 : }
999 : /* assume k >= 2 */
1000 36573987 : av = avma; scratch = new_chunk(k<<1); /* >= k + 2: result fits */
1001 :
1002 : /* copy T to scratch space (pad with zeroes to 2k words) */
1003 36573987 : Td = (GEN)av;
1004 36573987 : Te = T + (d+2);
1005 810094839 : for (i=0; i < d ; i++) *--Td = *--Te;
1006 62086605 : for ( ; i < (k<<1); i++) *--Td = 0;
1007 :
1008 36573987 : Te = (GEN)av; /* 1 beyond end of current T mantissa (in scratch) */
1009 36573987 : Ne = N + k+2; /* 1 beyond end of N mantissa */
1010 :
1011 36573987 : carry = 0;
1012 436090722 : for (i=0; i<k; i++) /* set T := T/B nod N, k times */
1013 : {
1014 399516735 : Td = Te; /* one beyond end of (new) T mantissa */
1015 399516735 : Nd = Ne;
1016 399516735 : hiremainder = *--Td;
1017 399516735 : m = hiremainder * inv; /* solve T + m N = O(B) */
1018 :
1019 : /* set T := (T + mN) / B */
1020 399516735 : Te = Td;
1021 399516735 : (void)addmul(m, *--Nd); /* = 0 */
1022 6662926071 : for (j=1; j<k; j++)
1023 : {
1024 6263409336 : t = addll(addmul(m, *--Nd), *--Td);
1025 6263409336 : *Td = t;
1026 6263409336 : hiremainder += overflow;
1027 : }
1028 399516735 : t = addll(hiremainder, *--Td); *Td = t + carry;
1029 399516735 : carry = (overflow || (carry && *Td == 0));
1030 : }
1031 36573987 : if (carry)
1032 : { /* Td > N overflows (k+1 words), set Td := Td - N */
1033 384435 : Td = Te;
1034 384435 : Nd = Ne;
1035 384435 : t = subll(*--Td, *--Nd); *Td = t;
1036 7176597 : while (Td > scratch) { t = subllx(*--Td, *--Nd); *Td = t; }
1037 : }
1038 :
1039 : /* copy result */
1040 36573987 : Td = (GEN)av;
1041 40372947 : while (*scratch == 0 && Te > scratch) scratch++; /* strip leading 0s */
1042 432291762 : while (Te > scratch) *--Td = *--Te;
1043 36573987 : k = (GEN)av - Td; if (!k) return gc_const(av, gen_0);
1044 36573987 : k += 2;
1045 36573987 : *--Td = evalsigne(1) | evallgefint(k);
1046 36573987 : *--Td = evaltyp(t_INT) | evallg(k);
1047 : #ifdef DEBUG
1048 : {
1049 : long l = NLIMBS(N), s = BITS_IN_LONG*l;
1050 : GEN R = int2n(s);
1051 : GEN res = remii(mulii(T, Fp_inv(R, N)), N);
1052 : if (k > lgefint(N)
1053 : || !equalii(remii(Td,N),res)
1054 : || cmpii(Td, addii(shifti(T, -s), N)) >= 0) pari_err_BUG("red_montgomery");
1055 : }
1056 : #endif
1057 36573987 : return gc_const((pari_sp)Td, Td);
1058 : }
1059 :
1060 : /* EXACT INTEGER DIVISION */
1061 :
1062 : /* assume xy>0, the division is exact and y is odd. Destroy x */
1063 : static GEN
1064 29512044 : diviuexact_i(GEN x, ulong y)
1065 : {
1066 : long i, lz, lx;
1067 : ulong q, yinv;
1068 : GEN z, z0, x0, x0min;
1069 :
1070 29512044 : if (y == 1) return icopy(x);
1071 23753835 : lx = lgefint(x);
1072 23753835 : if (lx == 3)
1073 : {
1074 849078 : q = uel(x,2) / y;
1075 849078 : if (!q) pari_err_OP("exact division", x, utoi(y));
1076 849078 : return utoipos(q);
1077 : }
1078 22904757 : yinv = invmod2BIL(y);
1079 22904757 : lz = (y <= uel(x,2)) ? lx : lx-1;
1080 22904757 : z = new_chunk(lz);
1081 22904757 : z0 = z + lz;
1082 22904757 : x0 = x + lx; x0min = x + lx-lz+2;
1083 :
1084 83397501 : while (x0 > x0min)
1085 : {
1086 60492744 : *--z0 = q = yinv*uel(--x0,0); /* i-th quotient */
1087 60492744 : if (!q) continue;
1088 : /* x := x - q * y */
1089 : { /* update neither lowest word (could set it to 0) nor highest ones */
1090 59963271 : GEN x1 = x0 - 1;
1091 : LOCAL_HIREMAINDER;
1092 59963271 : (void)mulll(q,y);
1093 59963271 : if (hiremainder)
1094 : {
1095 48044196 : if (uel(x1,0) < hiremainder)
1096 : {
1097 140112 : uel(x1,0) -= hiremainder;
1098 148449 : do uel(--x1,0)--; while (uel(x1,0) == ULONG_MAX);
1099 : }
1100 : else
1101 47904084 : uel(x1,0) -= hiremainder;
1102 : }
1103 : }
1104 : }
1105 22904757 : i=2; while(!z[i]) i++;
1106 22904757 : z += i-2; lz -= i-2;
1107 22904757 : z[0] = evaltyp(t_INT)|evallg(lz);
1108 22904757 : z[1] = evalsigne(1)|evallg(lz);
1109 22904757 : if (lz == 2) pari_err_OP("exact division", x, utoi(y));
1110 22904757 : return gc_const((pari_sp)z, z);
1111 : }
1112 :
1113 : /* assume y != 0 and the division is exact */
1114 : GEN
1115 28476120 : diviuexact(GEN x, ulong y)
1116 : {
1117 : pari_sp av;
1118 28476120 : long lx, vy, s = signe(x);
1119 : GEN z;
1120 :
1121 28476120 : if (!s) return gen_0;
1122 27629661 : if (y == 1) return icopy(x);
1123 24635454 : lx = lgefint(x);
1124 24635454 : if (lx == 3) {
1125 20344452 : ulong q = uel(x,2) / y;
1126 20344452 : if (!q) pari_err_OP("exact division", x, utoi(y));
1127 20344452 : return (s > 0)? utoipos(q): utoineg(q);
1128 : }
1129 4291002 : av = avma; (void)new_chunk(lx); vy = vals(y);
1130 4291002 : if (vy) {
1131 1713636 : y >>= vy;
1132 1713636 : if (y == 1) { set_avma(av); return shifti(x, -vy); }
1133 805233 : x = shifti(x, -vy);
1134 805233 : if (lx == 3) {
1135 0 : ulong q = uel(x,2) / y;
1136 0 : set_avma(av);
1137 0 : if (!q) pari_err_OP("exact division", x, utoi(y));
1138 0 : return (s > 0)? utoipos(q): utoineg(q);
1139 : }
1140 2577366 : } else x = icopy(x);
1141 3382599 : set_avma(av);
1142 3382599 : z = diviuexact_i(x, y);
1143 3382599 : setsigne(z, s); return z;
1144 : }
1145 :
1146 : /* Find z such that x=y*z, knowing that y | x (unchecked)
1147 : * Method: y0 z0 = x0 mod B = 2^BITS_IN_LONG ==> z0 = 1/y0 mod B.
1148 : * Set x := (x - z0 y) / B, updating only relevant words, and repeat */
1149 : GEN
1150 390661989 : diviiexact(GEN x, GEN y)
1151 : {
1152 390661989 : long lx, ly, lz, vy, i, ii, sx = signe(x), sy = signe(y);
1153 : pari_sp av;
1154 : ulong y0inv,q;
1155 : GEN z;
1156 :
1157 390661989 : if (!sy) pari_err_INV("diviiexact",gen_0);
1158 390661989 : if (!sx) return gen_0;
1159 322543005 : lx = lgefint(x);
1160 322543005 : if (lx == 3) {
1161 259448136 : q = uel(x,2) / uel(y,2);
1162 259448136 : if (!q) pari_err_OP("exact division", x, y);
1163 259448136 : return (sx+sy) ? utoipos(q): utoineg(q);
1164 : }
1165 63094869 : vy = vali(y); av = avma;
1166 63094869 : (void)new_chunk(lx); /* enough room for z */
1167 63094869 : if (vy)
1168 : { /* make y odd */
1169 32246580 : y = shifti(y,-vy);
1170 32246580 : x = shifti(x,-vy); lx = lgefint(x);
1171 : }
1172 30848289 : else x = icopy(x); /* necessary because we destroy x */
1173 63094869 : set_avma(av); /* will erase our x,y when exiting */
1174 : /* now y is odd */
1175 63094869 : ly = lgefint(y);
1176 63094869 : if (ly == 3)
1177 : {
1178 26129445 : z = diviuexact_i(x,uel(y,2)); /* x != 0 */
1179 26129445 : setsigne(z, (sx+sy)? 1: -1); return z;
1180 : }
1181 36965424 : y0inv = invmod2BIL(y[ly-1]);
1182 58326495 : i=2; while (i<ly && y[i]==x[i]) i++;
1183 36965424 : lz = (i==ly || uel(y,i) < uel(x,i)) ? lx-ly+3 : lx-ly+2;
1184 36965424 : z = new_chunk(lz);
1185 :
1186 36965424 : y += ly - 1; /* now y[-i] = i-th word of y */
1187 173149464 : for (ii=lx-1,i=lz-1; i>=2; i--,ii--)
1188 : {
1189 : long limj;
1190 : LOCAL_HIREMAINDER;
1191 : LOCAL_OVERFLOW;
1192 :
1193 136184040 : z[i] = q = y0inv*uel(x,ii); /* i-th quotient */
1194 136184040 : if (!q) continue;
1195 :
1196 : /* x := x - q * y */
1197 136052262 : (void)mulll(q,y[0]); limj = maxss(lx - lz, ii+3-ly);
1198 : { /* update neither lowest word (could set it to 0) nor highest ones */
1199 136052262 : GEN x0 = x + (ii - 1), y0 = y - 1, xlim = x + limj;
1200 2310051297 : for (; x0 >= xlim; x0--, y0--)
1201 : {
1202 2173999035 : *x0 = subll(*x0, addmul(q,*y0));
1203 2173999035 : hiremainder += overflow;
1204 : }
1205 136052262 : if (hiremainder && limj != lx - lz)
1206 : {
1207 72092880 : if ((ulong)*x0 < hiremainder)
1208 : {
1209 831369 : *x0 -= hiremainder;
1210 831387 : do (*--x0)--; while ((ulong)*x0 == ULONG_MAX);
1211 : }
1212 : else
1213 71261511 : *x0 -= hiremainder;
1214 : }
1215 : }
1216 : }
1217 36965424 : i=2; while(!z[i]) i++;
1218 36965424 : z += i-2; lz -= (i-2);
1219 36965424 : z[0] = evaltyp(t_INT)|evallg(lz);
1220 36965424 : z[1] = evalsigne((sx+sy)? 1: -1) | evallg(lz);
1221 36965424 : if (lz == 2) pari_err_OP("exact division", x, y);
1222 36965424 : return gc_const((pari_sp)z, z);
1223 : }
1224 :
1225 : /* assume yz != and yz | x */
1226 : GEN
1227 149898 : diviuuexact(GEN x, ulong y, ulong z)
1228 : {
1229 : long tmp[4];
1230 : ulong t;
1231 : LOCAL_HIREMAINDER;
1232 149898 : t = mulll(y, z);
1233 149898 : if (!hiremainder) return diviuexact(x, t);
1234 0 : tmp[0] = evaltyp(t_INT)|_evallg(4);
1235 0 : tmp[1] = evalsigne(1)|evallgefint(4);
1236 0 : tmp[2] = hiremainder;
1237 0 : tmp[3] = t;
1238 0 : return diviiexact(x, tmp);
1239 : }
1240 :
1241 : /********************************************************************/
1242 : /** **/
1243 : /** INTEGER MULTIPLICATION (BASECASE) **/
1244 : /** **/
1245 : /********************************************************************/
1246 : /* nx >= ny = num. of digits of x, y (not GEN, see mulii) */
1247 : INLINE GEN
1248 5199216126 : muliispec_basecase(GEN x, GEN y, long nx, long ny)
1249 : {
1250 : GEN z2e,z2d,yd,xd,ye,zd;
1251 : long p1,lz;
1252 : LOCAL_HIREMAINDER;
1253 :
1254 5199216126 : if (ny == 1) return muluispec((ulong)*y, x, nx);
1255 1129314264 : if (ny == 0) return gen_0;
1256 1128096465 : zd = (GEN)avma; lz = nx+ny+2;
1257 1128096465 : (void)new_chunk(lz);
1258 1128096465 : xd = x + nx;
1259 1128096465 : yd = y + ny;
1260 1128096465 : ye = yd; p1 = *--xd;
1261 :
1262 1128096465 : *--zd = mulll(p1, *--yd); z2e = zd;
1263 9455354517 : while (yd > y) *--zd = addmul(p1, *--yd);
1264 1128096465 : *--zd = hiremainder;
1265 :
1266 10810737375 : while (xd > x)
1267 : {
1268 : LOCAL_OVERFLOW;
1269 9682640910 : yd = ye; p1 = *--xd;
1270 :
1271 9682640910 : z2d = --z2e;
1272 9682640910 : *z2d = addll(mulll(p1, *--yd), *z2d); z2d--;
1273 >12306*10^7 : while (yd > y)
1274 : {
1275 >11338*10^7 : hiremainder += overflow;
1276 >11338*10^7 : *z2d = addll(addmul(p1, *--yd), *z2d); z2d--;
1277 : }
1278 9682640910 : *--zd = hiremainder + overflow;
1279 : }
1280 1128096465 : if (*zd == 0) { zd++; lz--; } /* normalize */
1281 1128096465 : *--zd = evalsigne(1) | evallgefint(lz);
1282 1128096465 : *--zd = evaltyp(t_INT) | evallg(lz);
1283 1128096465 : return gc_const((pari_sp)zd, zd);
1284 : }
1285 :
1286 : INLINE GEN
1287 928113927 : sqrispec_basecase(GEN x, long nx)
1288 : {
1289 : GEN z2e,z2d,yd,xd,zd,x0,z0;
1290 : long p1,lz;
1291 : LOCAL_HIREMAINDER;
1292 : LOCAL_OVERFLOW;
1293 :
1294 928113927 : if (nx == 1) return sqru((ulong)*x);
1295 627678324 : if (nx == 0) return gen_0;
1296 223614711 : zd = (GEN)avma; lz = (nx+1) << 1;
1297 223614711 : z0 = new_chunk(lz);
1298 223614711 : if (nx == 1)
1299 : {
1300 0 : *--zd = mulll(*x, *x);
1301 0 : *--zd = hiremainder; goto END;
1302 : }
1303 223614711 : xd = x + nx;
1304 :
1305 : /* compute double products --> zd */
1306 223614711 : p1 = *--xd; yd = xd; --zd;
1307 223614711 : *--zd = mulll(p1, *--yd); z2e = zd;
1308 1355624118 : while (yd > x) *--zd = addmul(p1, *--yd);
1309 223614711 : *--zd = hiremainder;
1310 :
1311 223614711 : x0 = x+1;
1312 1355624118 : while (xd > x0)
1313 : {
1314 : LOCAL_OVERFLOW;
1315 1132009407 : p1 = *--xd; yd = xd;
1316 :
1317 1132009407 : z2e -= 2; z2d = z2e;
1318 1132009407 : *z2d = addll(mulll(p1, *--yd), *z2d); z2d--;
1319 10555383444 : while (yd > x)
1320 : {
1321 9423374037 : hiremainder += overflow;
1322 9423374037 : *z2d = addll(addmul(p1, *--yd), *z2d); z2d--;
1323 : }
1324 1132009407 : *--zd = hiremainder + overflow;
1325 : }
1326 : /* multiply zd by 2 (put result in zd - 1) */
1327 223614711 : zd[-1] = ((*zd & HIGHBIT) != 0);
1328 223614711 : shift_left(zd, zd, 0, (nx<<1)-3, 0, 1);
1329 :
1330 : /* add the squares */
1331 223614711 : xd = x + nx; zd = z0 + lz;
1332 223614711 : p1 = *--xd;
1333 223614711 : zd--; *zd = mulll(p1,p1);
1334 223614711 : zd--; *zd = addll(hiremainder, *zd);
1335 1579238829 : while (xd > x)
1336 : {
1337 1355624118 : p1 = *--xd;
1338 1355624118 : zd--; *zd = addll(mulll(p1,p1)+ overflow, *zd);
1339 1355624118 : zd--; *zd = addll(hiremainder + overflow, *zd);
1340 : }
1341 :
1342 223614711 : END:
1343 223614711 : if (*zd == 0) { zd++; lz--; } /* normalize */
1344 223614711 : *--zd = evalsigne(1) | evallgefint(lz);
1345 223614711 : *--zd = evaltyp(t_INT) | evallg(lz);
1346 223614711 : return gc_const((pari_sp)zd, zd);
1347 : }
1348 :
1349 : /********************************************************************/
1350 : /** **/
1351 : /** INTEGER MULTIPLICATION (FFT) **/
1352 : /** **/
1353 : /********************************************************************/
1354 :
1355 : /*
1356 : Compute parameters for FFT:
1357 : len: result length
1358 : k: FFT depth.
1359 : n: number of blocks (2^k)
1360 : bs: block size
1361 : mod: Modulus is M=2^(BIL*mod)+1
1362 : ord: order of 2 in Z/MZ.
1363 : We must have:
1364 : bs*n >= l
1365 : 2^(BIL*mod) > nb*2^(2*BIL*bs)
1366 : 2^k | 2*BIL*mod
1367 : */
1368 : static void
1369 85290 : mulliifft_params(long len, long *k, long *mod, long *bs, long *n, ulong *ord)
1370 : {
1371 : long r;
1372 85290 : *k = expu((3*len)>>2)-3;
1373 : do {
1374 85293 : (*k)--;
1375 85293 : r = *k-(TWOPOTBITS_IN_LONG+2);
1376 85293 : *n = 1L<<*k;
1377 85293 : *bs = (len+*n-1)>>*k;
1378 85293 : *mod= 2**bs+1;
1379 85293 : if (r>0)
1380 5181 : *mod=((*mod+(1L<<r)-1)>>r)<<r;
1381 85293 : } while(*mod>=3**bs);
1382 85290 : *ord= 4**mod*BITS_IN_LONG;
1383 85290 : }
1384 :
1385 : /* Zf_: arithmetic in ring Z/MZ where M= 2^(BITS_IN_LONG*mod)+1
1386 : * for some mod.
1387 : * Do not garbage collect.
1388 : */
1389 :
1390 : static GEN
1391 187658496 : Zf_add(GEN a, GEN b, GEN M)
1392 : {
1393 187658496 : GEN y, z = addii(a,b);
1394 187658496 : long mod = lgefint(M)-3;
1395 187658496 : long l = NLIMBS(z);
1396 187658496 : if (l<=mod) return z;
1397 72714114 : y = subiu(z, 1);
1398 72714114 : if (NLIMBS(y)<=mod) return z;
1399 72714114 : return int_normalize(y,1);
1400 : }
1401 :
1402 : static GEN
1403 191070978 : Zf_sub(GEN a, GEN b, GEN M)
1404 : {
1405 191070978 : GEN z = subii(a,b);
1406 191070978 : return signe(z)>=0? z: addii(M,z);
1407 : }
1408 :
1409 : /* destroy z */
1410 : static GEN
1411 399271983 : Zf_red_destroy(GEN z, GEN M)
1412 : {
1413 399271983 : long mod = lgefint(M)-3;
1414 399271983 : long l = NLIMBS(z);
1415 : GEN y;
1416 399271983 : if (l<=mod) return z;
1417 177220563 : y = shifti(z, -mod*BITS_IN_LONG);
1418 177220563 : z = int_normalize(z, NLIMBS(y));
1419 177220563 : y = Zf_red_destroy(y, M);
1420 177220563 : z = subii(z, y);
1421 177220563 : if (signe(z)<0) z = addii(z, M);
1422 177220563 : return z;
1423 : }
1424 :
1425 : INLINE GEN
1426 206264028 : Zf_shift(GEN a, ulong s, GEN M) { return Zf_red_destroy(shifti(a, s), M); }
1427 :
1428 : /*
1429 : Multiply by sqrt(2)^s
1430 : We use the formula sqrt(2)=z_8*(1-z_4)) && z_8=2^(ord/16) [2^(ord/4)+1]
1431 : */
1432 :
1433 : static GEN
1434 187658496 : Zf_mulsqrt2(GEN a, ulong s, ulong ord, GEN M)
1435 : {
1436 187658496 : ulong hord = ord>>1;
1437 187658496 : if (!signe(a)) return gen_0;
1438 183651672 : if (odd(s)) /* Multiply by 2^(s/2) */
1439 : {
1440 3412482 : GEN az8 = Zf_shift(a, ord>>4, M);
1441 3412482 : GEN az83 = Zf_shift(az8, ord>>3, M);
1442 3412482 : a = Zf_sub(az8, az83, M);
1443 3412482 : s--;
1444 : }
1445 183651672 : if (s < hord)
1446 136430355 : return Zf_shift(a, s>>1, M);
1447 : else
1448 47221317 : return subii(M,Zf_shift(a, (s-hord)>>1, M));
1449 : }
1450 :
1451 : INLINE GEN
1452 448896 : Zf_sqr(GEN a, GEN M) { return Zf_red_destroy(sqri(a), M); }
1453 :
1454 : INLINE GEN
1455 15338496 : Zf_mul(GEN a, GEN b, GEN M) { return Zf_red_destroy(mulii(a,b), M); }
1456 :
1457 : /* In place, bit reversing FFT */
1458 : static void
1459 30956967 : muliifft_dit(ulong o, ulong ord, GEN M, GEN FFT, long d, long step)
1460 : {
1461 30956967 : pari_sp av = avma;
1462 : long i;
1463 30956967 : ulong j, no = (o<<1)%ord;
1464 30956967 : long hstep=step>>1;
1465 155404263 : for (i = d+1, j = 0; i <= d+hstep; ++i, j =(j+o)%ord)
1466 : {
1467 124447296 : GEN a = Zf_add(gel(FFT,i), gel(FFT,i+hstep), M);
1468 124447296 : GEN b = Zf_mulsqrt2(Zf_sub(gel(FFT,i), gel(FFT,i+hstep), M), j, ord, M);
1469 124447296 : affii(a,gel(FFT,i));
1470 124447296 : affii(b,gel(FFT,i+hstep));
1471 124447296 : set_avma(av);
1472 : }
1473 30956967 : if (hstep>1)
1474 : {
1475 15394023 : muliifft_dit(no, ord, M, FFT, d, hstep);
1476 15394023 : muliifft_dit(no, ord, M, FFT, d+hstep, hstep);
1477 : }
1478 30956967 : }
1479 :
1480 : /* In place, bit reversed FFT, inverse of muliifft_dit */
1481 : static void
1482 15702102 : muliifft_dis(ulong o, ulong ord, GEN M, GEN FFT, long d, long step)
1483 : {
1484 15702102 : pari_sp av = avma;
1485 : long i;
1486 15702102 : ulong j, no = (o<<1)%ord;
1487 15702102 : long hstep=step>>1;
1488 15702102 : if (hstep>1)
1489 : {
1490 7808406 : muliifft_dis(no, ord, M, FFT, d, hstep);
1491 7808406 : muliifft_dis(no, ord, M, FFT, d+hstep, hstep);
1492 : }
1493 78913302 : for (i = d+1, j = 0; i <= d+hstep; ++i, j =(j+o)%ord)
1494 : {
1495 63211200 : GEN z = Zf_mulsqrt2(gel(FFT,i+hstep), j, ord, M);
1496 63211200 : GEN a = Zf_add(gel(FFT,i), z, M);
1497 63211200 : GEN b = Zf_sub(gel(FFT,i), z, M);
1498 63211200 : affii(a,gel(FFT,i));
1499 63211200 : affii(b,gel(FFT,i+hstep));
1500 63211200 : set_avma(av);
1501 : }
1502 15702102 : }
1503 :
1504 : static GEN
1505 168921 : muliifft_spliti(GEN a, long na, long bs, long n, long mod)
1506 : {
1507 168921 : GEN ap = a+na-1;
1508 168921 : GEN c = cgetg(n+1, t_VEC);
1509 : long i,j;
1510 31294809 : for(i=1;i<=n;i++)
1511 : {
1512 31125888 : GEN z = cgeti(mod+3);
1513 31125888 : if (na)
1514 : {
1515 15326127 : long m = minss(bs, na), v=0;
1516 15326127 : GEN zp, aa=ap-m+1;
1517 83251410 : while (!*aa && v<m) {aa++; v++;}
1518 15326127 : zp = z+m-v+1;
1519 381514206 : for (j=v; j < m; j++)
1520 366188079 : *zp-- = *ap--;
1521 15326127 : ap -= v; na -= m;
1522 15326127 : z[1] = evalsigne(m!=v) | evallgefint(m-v+2);
1523 : } else
1524 15799761 : z[1] = evalsigne(0) | evallgefint(2);
1525 31125888 : gel(c, i) = z;
1526 : }
1527 168921 : return c;
1528 : }
1529 :
1530 : static GEN
1531 85290 : muliifft_unspliti(GEN V, long bs, long len)
1532 : {
1533 85290 : long s, i, j, l = lg(V);
1534 85290 : GEN a = cgeti(len);
1535 85290 : a[1] = evalsigne(1)|evallgefint(len);
1536 440997075 : for(i=2;i<len;i++)
1537 440911785 : a[i] = 0;
1538 15872682 : for(i=1, s=0; i<l; i++, s+=bs)
1539 : {
1540 15787392 : GEN u = gel(V,i);
1541 15787392 : if (signe(u))
1542 : {
1543 15211731 : GEN ap = int_W(a,s);
1544 15211731 : GEN up = int_LSW(u);
1545 15211731 : long lu = NLIMBS(u);
1546 : LOCAL_OVERFLOW;
1547 15211731 : *ap = addll(*ap, *up--); ap--;
1548 862324044 : for (j=1; j<lu; j++)
1549 847112313 : { *ap = addllx(*ap, *up--); ap--; }
1550 15214287 : while (overflow)
1551 2556 : { *ap = addll(*ap, 1); ap--; }
1552 : }
1553 : }
1554 85290 : return int_normalize(a,0);
1555 : }
1556 :
1557 : static GEN
1558 1659 : sqrispec_fft(GEN a, long na)
1559 : {
1560 1659 : pari_sp av, ltop = avma;
1561 1659 : long len = 2*na;
1562 : long k, mod, bs, n;
1563 : GEN FFT, M;
1564 : long i;
1565 : ulong o, ord;
1566 :
1567 1659 : mulliifft_params(len,&k,&mod,&bs,&n,&ord);
1568 1659 : o = ord>>k;
1569 1659 : M = int2n(mod*BITS_IN_LONG);
1570 1659 : M[2+mod] = 1;
1571 1659 : FFT = muliifft_spliti(a, na, bs, n, mod);
1572 1659 : muliifft_dit(o, ord, M, FFT, 0, n);
1573 1659 : av = avma;
1574 450555 : for(i=1; i<=n; i++)
1575 : {
1576 448896 : affii(Zf_sqr(gel(FFT,i), M), gel(FFT,i));
1577 448896 : set_avma(av);
1578 : }
1579 1659 : muliifft_dis(ord-o, ord, M, FFT, 0, n);
1580 450555 : for(i=1; i<=n; i++)
1581 : {
1582 448896 : affii(Zf_shift(gel(FFT,i), (ord>>1)-k, M), gel(FFT,i));
1583 448896 : set_avma(av);
1584 : }
1585 1659 : return gc_INT(ltop, muliifft_unspliti(FFT,bs,2+len));
1586 : }
1587 :
1588 : static GEN
1589 83631 : muliispec_fft(GEN a, GEN b, long na, long nb)
1590 : {
1591 83631 : pari_sp av, av2, ltop = avma;
1592 83631 : long len = na+nb;
1593 : long k, mod, bs, n;
1594 : GEN FFT, FFTb, M;
1595 : long i;
1596 : ulong o, ord;
1597 :
1598 83631 : mulliifft_params(len,&k,&mod,&bs,&n,&ord);
1599 83631 : o = ord>>k;
1600 83631 : M = int2n(mod*BITS_IN_LONG);
1601 83631 : M[2+mod] = 1;
1602 83631 : FFT = muliifft_spliti(a, na, bs, n, mod);
1603 83631 : av=avma;
1604 83631 : muliifft_dit(o, ord, M, FFT, 0, n);
1605 83631 : FFTb = muliifft_spliti(b, nb, bs, n, mod);
1606 83631 : av2 = avma;
1607 83631 : muliifft_dit(o, ord, M, FFTb, 0, n);
1608 15422127 : for(i=1; i<=n; i++)
1609 : {
1610 15338496 : affii(Zf_mul(gel(FFT,i), gel(FFTb,i), M), gel(FFT,i));
1611 15338496 : set_avma(av2);
1612 : }
1613 83631 : set_avma(av);
1614 83631 : muliifft_dis(ord-o, ord, M, FFT, 0, n);
1615 15422127 : for(i=1; i<=n; i++)
1616 : {
1617 15338496 : affii(Zf_shift(gel(FFT,i),(ord>>1)-k,M), gel(FFT,i));
1618 15338496 : set_avma(av);
1619 : }
1620 83631 : return gc_INT(ltop, muliifft_unspliti(FFT,bs,2+len));
1621 : }
1622 :
1623 : /********************************************************************/
1624 : /** **/
1625 : /** INTEGER MULTIPLICATION (KARATSUBA) **/
1626 : /** **/
1627 : /********************************************************************/
1628 :
1629 : /* return (x shifted left d words) + y. Assume d > 0, x > 0 and y >= 0 */
1630 : static GEN
1631 704270652 : addshiftw(GEN x, GEN y, long d)
1632 : {
1633 704270652 : GEN z,z0,y0,yd, zd = (GEN)avma;
1634 704270652 : long a,lz,ly = lgefint(y);
1635 :
1636 704270652 : z0 = new_chunk(d);
1637 704270652 : a = ly-2; yd = y+ly;
1638 704270652 : if (a >= d)
1639 : {
1640 12707316654 : y0 = yd-d; while (yd > y0) *--zd = *--yd; /* copy last d words of y */
1641 699720489 : a -= d;
1642 699720489 : if (a)
1643 467340426 : z = addiispec(LIMBS(x), LIMBS(y), NLIMBS(x), a);
1644 : else
1645 232380063 : z = icopy(x);
1646 : }
1647 : else
1648 : {
1649 16648863 : y0 = yd-a; while (yd > y0) *--zd = *--yd; /* copy last a words of y */
1650 69657453 : while (zd > z0) *--zd = 0; /* complete with 0s */
1651 4550163 : z = icopy(x);
1652 : }
1653 704270652 : lz = lgefint(z)+d;
1654 704270652 : z[1] = evalsigne(1) | evallgefint(lz);
1655 704270652 : z[0] = evaltyp(t_INT) | evallg(lz); return z;
1656 : }
1657 :
1658 : /* Fast product (Karatsuba) of integers. a and b are "special" GENs
1659 : * c,c0,c1,c2 are genuine GENs.
1660 : */
1661 : GEN
1662 5426141082 : muliispec(GEN a, GEN b, long na, long nb)
1663 : {
1664 : GEN a0,c,c0;
1665 : long n0, n0a, i;
1666 : pari_sp av;
1667 :
1668 5426141082 : if (na < nb) swapspec(a,b, na,nb);
1669 5426141082 : if (nb < MULII_KARATSUBA_LIMIT) return muliispec_basecase(a,b,na,nb);
1670 226924956 : if (nb >= MULII_FFT_LIMIT) return muliispec_fft(a,b,na,nb);
1671 226841325 : i=(na>>1); n0=na-i; na=i;
1672 226841325 : av=avma; a0=a+na; n0a=n0;
1673 340243788 : while (n0a && !*a0) { a0++; n0a--; }
1674 :
1675 226841325 : if (n0a && nb > n0)
1676 223351743 : { /* nb <= na <= n0 */
1677 : GEN b0,c1,c2;
1678 : long n0b;
1679 :
1680 223351743 : nb -= n0;
1681 223351743 : c = muliispec(a,b,na,nb);
1682 223351743 : b0 = b+nb; n0b = n0;
1683 321639213 : while (n0b && !*b0) { b0++; n0b--; }
1684 223351743 : if (n0b)
1685 : {
1686 222586539 : c0 = muliispec(a0,b0, n0a,n0b);
1687 :
1688 222586539 : c2 = addiispec(a0,a, n0a,na);
1689 222586539 : c1 = addiispec(b0,b, n0b,nb);
1690 222586539 : c1 = muliispec(LIMBS(c1),LIMBS(c2), NLIMBS(c1),NLIMBS(c2));
1691 222586539 : c2 = addiispec(LIMBS(c0),LIMBS(c), NLIMBS(c0),NLIMBS(c));
1692 :
1693 222586539 : c1 = subiispec(LIMBS(c1),LIMBS(c2), NLIMBS(c1),NLIMBS(c2));
1694 : }
1695 : else
1696 : {
1697 765204 : c0 = gen_0;
1698 765204 : c1 = muliispec(a0,b, n0a,nb);
1699 : }
1700 223351743 : c = addshiftw(c,c1, n0);
1701 : }
1702 : else
1703 : {
1704 3489582 : c = muliispec(a,b,na,nb);
1705 3489582 : c0 = muliispec(a0,b,n0a,nb);
1706 : }
1707 226841325 : return gc_INT(av, addshiftw(c,c0, n0));
1708 : }
1709 : GEN
1710 166566 : muluui(ulong x, ulong y, GEN z)
1711 : {
1712 166566 : long t, s = signe(z);
1713 : GEN r;
1714 : LOCAL_HIREMAINDER;
1715 :
1716 166566 : if (!x || !y || !signe(z)) return gen_0;
1717 166269 : t = mulll(x,y);
1718 166269 : if (!hiremainder)
1719 166269 : r = muluispec(t, z+2, lgefint(z)-2);
1720 : else
1721 : {
1722 : long tmp[2];
1723 0 : tmp[0] = hiremainder;
1724 0 : tmp[1] = t;
1725 0 : r = muliispec(z+2,tmp,lgefint(z)-2,2);
1726 : }
1727 166269 : setsigne(r,s); return r;
1728 : }
1729 :
1730 : #define sqrispec_mirror sqrispec
1731 : #define muliispec_mirror muliispec
1732 :
1733 : /* x % (2^n), assuming n >= 0 */
1734 : GEN
1735 19702815 : remi2n(GEN x, long n)
1736 : {
1737 19702815 : long hi,l,k,lx,ly, sx = signe(x);
1738 : GEN z, xd, zd;
1739 :
1740 19702815 : if (!sx || !n) return gen_0;
1741 :
1742 19679310 : k = dvmdsBIL(n, &l);
1743 19679310 : lx = lgefint(x);
1744 19679310 : if (lx < k+3) return icopy(x);
1745 :
1746 19292469 : xd = x + (lx-k-1);
1747 : /* x = |_|...|#|1|...|k| : copy the last l bits of # and the last k words
1748 : * ^--- initial xd */
1749 19292469 : hi = ((ulong)*xd) & ((1UL<<l)-1); /* last l bits of # = top bits of result */
1750 19292469 : if (!hi)
1751 : { /* strip leading zeroes from result */
1752 1391364 : xd++; while (k && !*xd) { k--; xd++; }
1753 1361142 : if (!k) return gen_0;
1754 600975 : ly = k+2; xd--;
1755 : }
1756 : else
1757 17931327 : ly = k+3;
1758 :
1759 18532302 : zd = z = cgeti(ly);
1760 18532302 : *++zd = evalsigne(sx) | evallgefint(ly);
1761 18532302 : if (hi) *++zd = hi;
1762 102722385 : for ( ;k; k--) *++zd = *++xd;
1763 18532302 : return z;
1764 : }
1765 :
1766 : GEN
1767 938053002 : sqrispec(GEN a, long na)
1768 : {
1769 : GEN a0,c;
1770 : long n0, n0a, i;
1771 : pari_sp av;
1772 :
1773 938053002 : if (na < SQRI_KARATSUBA_LIMIT) return sqrispec_basecase(a,na);
1774 9939075 : if (na >= SQRI_FFT_LIMIT) return sqrispec_fft(a,na);
1775 9937416 : i=(na>>1); n0=na-i; na=i;
1776 9937416 : av=avma; a0=a+na; n0a=n0;
1777 14932608 : while (n0a && !*a0) { a0++; n0a--; }
1778 9937416 : c = sqrispec(a,na);
1779 9937416 : if (n0a)
1780 : {
1781 9928518 : GEN t, c1, c0 = sqrispec(a0,n0a);
1782 : #if 0
1783 : c1 = shifti(muliispec(a0,a, n0a,na),1);
1784 : #else /* faster */
1785 9928518 : t = addiispec(a0,a,n0a,na);
1786 9928518 : t = sqrispec(LIMBS(t),NLIMBS(t));
1787 9928518 : c1= addiispec(LIMBS(c0),LIMBS(c), NLIMBS(c0), NLIMBS(c));
1788 9928518 : c1= subiispec(LIMBS(t),LIMBS(c1), NLIMBS(t), NLIMBS(c1));
1789 : #endif
1790 9928518 : c = addshiftw(c,c1, n0);
1791 9928518 : c = addshiftw(c,c0, n0);
1792 : }
1793 : else
1794 8898 : c = addshiftw(c,gen_0,n0<<1);
1795 9937416 : return gc_INT(av, c);
1796 : }
1797 :
1798 : /********************************************************************/
1799 : /** **/
1800 : /** KARATSUBA SQUARE ROOT **/
1801 : /** adapted from Paul Zimmermann's implementation of **/
1802 : /** his algorithm in GMP (mpn_sqrtrem) **/
1803 : /** **/
1804 : /********************************************************************/
1805 :
1806 : /* Square roots table */
1807 : static const unsigned char approx_tab[192] = {
1808 : 128,128,129,130,131,132,133,134,135,136,137,138,139,140,141,142,
1809 : 143,144,144,145,146,147,148,149,150,150,151,152,153,154,155,155,
1810 : 156,157,158,159,160,160,161,162,163,163,164,165,166,167,167,168,
1811 : 169,170,170,171,172,173,173,174,175,176,176,177,178,178,179,180,
1812 : 181,181,182,183,183,184,185,185,186,187,187,188,189,189,190,191,
1813 : 192,192,193,193,194,195,195,196,197,197,198,199,199,200,201,201,
1814 : 202,203,203,204,204,205,206,206,207,208,208,209,209,210,211,211,
1815 : 212,212,213,214,214,215,215,216,217,217,218,218,219,219,220,221,
1816 : 221,222,222,223,224,224,225,225,226,226,227,227,228,229,229,230,
1817 : 230,231,231,232,232,233,234,234,235,235,236,236,237,237,238,238,
1818 : 239,240,240,241,241,242,242,243,243,244,244,245,245,246,246,247,
1819 : 247,248,248,249,249,250,250,251,251,252,252,253,253,254,254,255
1820 : };
1821 :
1822 : /* N[0], assume N[0] >= 2^(BIL-2).
1823 : * Return r,s such that s^2 + r = N, 0 <= r <= 2s */
1824 : static void
1825 94687737 : p_sqrtu1(ulong *N, ulong *ps, ulong *pr)
1826 : {
1827 94687737 : ulong prec, r, s, q, u, n0 = N[0];
1828 :
1829 94687737 : q = n0 >> (BITS_IN_LONG - 8);
1830 : /* 2^6 = 64 <= q < 256 = 2^8 */
1831 94687737 : s = approx_tab[q - 64]; /* 128 <= s < 255 */
1832 94687737 : r = (n0 >> (BITS_IN_LONG - 16)) - s * s; /* r <= 2*s */
1833 94687737 : if (r > (s << 1)) { r -= (s << 1) | 1; s++; }
1834 :
1835 : /* 8-bit approximation from the high 8-bits of N[0] */
1836 94687737 : prec = 8;
1837 94687737 : n0 <<= 2 * prec;
1838 284063211 : while (2 * prec < BITS_IN_LONG)
1839 : { /* invariant: s has prec bits, and r <= 2*s */
1840 189375474 : r = (r << prec) + (n0 >> (BITS_IN_LONG - prec));
1841 189375474 : n0 <<= prec;
1842 189375474 : u = 2 * s;
1843 189375474 : q = r / u; u = r - q * u;
1844 189375474 : s = (s << prec) + q;
1845 189375474 : u = (u << prec) + (n0 >> (BITS_IN_LONG - prec));
1846 189375474 : q = q * q;
1847 189375474 : r = u - q;
1848 189375474 : if (u < q) { s--; r += (s << 1) | 1; }
1849 189375474 : n0 <<= prec;
1850 189375474 : prec = 2 * prec;
1851 : }
1852 94687737 : *ps = s;
1853 94687737 : *pr = r;
1854 94687737 : }
1855 :
1856 : /* N[0..1], assume N[0] >= 2^(BIL-2).
1857 : * Return 1 if remainder overflows, 0 otherwise */
1858 : static int
1859 91985496 : p_sqrtu2(ulong *N, ulong *ps, ulong *pr)
1860 : {
1861 91985496 : ulong cc, qhl, r, s, q, u, n1 = N[1];
1862 : LOCAL_OVERFLOW;
1863 :
1864 91985496 : p_sqrtu1(N, &s, &r); /* r <= 2s */
1865 137605458 : qhl = 0; while (r >= s) { qhl++; r -= s; }
1866 : /* now r < s < 2^(BIL/2) */
1867 91985496 : r = (r << BITS_IN_HALFULONG) | (n1 >> BITS_IN_HALFULONG);
1868 91985496 : u = s << 1;
1869 91985496 : q = r / u; u = r - q * u;
1870 91985496 : q += (qhl & 1) << (BITS_IN_HALFULONG - 1);
1871 91985496 : qhl >>= 1;
1872 : /* (initial r)<<(BIL/2) + n1>>(BIL/2) = (qhl<<(BIL/2) + q) * 2s + u */
1873 91985496 : s = ((s + qhl) << BITS_IN_HALFULONG) + q;
1874 91985496 : cc = u >> BITS_IN_HALFULONG;
1875 91985496 : r = (u << BITS_IN_HALFULONG) | (n1 & LOWMASK);
1876 91985496 : r = subll(r, q * q);
1877 91985496 : cc -= overflow + qhl;
1878 : /* now subtract 2*q*2^(BIL/2) + 2^BIL if qhl is set */
1879 91985496 : if ((long)cc < 0)
1880 : {
1881 23472189 : if (s) {
1882 23419161 : r = addll(r, s);
1883 23419161 : cc += overflow;
1884 23419161 : s--;
1885 : } else {
1886 53028 : cc++;
1887 53028 : s = ~0UL;
1888 : }
1889 23472189 : r = addll(r, s);
1890 23472189 : cc += overflow;
1891 : }
1892 91985496 : *ps = s;
1893 91985496 : *pr = r; return cc;
1894 : }
1895 :
1896 : static void
1897 90924714 : xmpn_zero(GEN x, long n)
1898 : {
1899 695325138 : while (--n >= 0) x[n]=0;
1900 90924714 : }
1901 : static void
1902 1068383172 : xmpn_copy(GEN z, GEN x, long n)
1903 : {
1904 1068383172 : long k = n;
1905 4198937577 : while (--k >= 0) z[k] = x[k];
1906 1068383172 : }
1907 : /* a[0..la-1] * 2^(lb BIL) | b[0..lb-1] */
1908 : static GEN
1909 468423300 : catii(GEN a, long la, GEN b, long lb)
1910 : {
1911 468423300 : long l = la + lb + 2;
1912 468423300 : GEN z = cgetipos(l);
1913 468423300 : xmpn_copy(LIMBS(z), a, la);
1914 468423300 : xmpn_copy(LIMBS(z) + la, b, lb);
1915 468423300 : return int_normalize(z, 0);
1916 : }
1917 :
1918 : /* sqrt n[0..1], assume n normalized */
1919 : static GEN
1920 91710672 : sqrtispec2(GEN n, GEN *pr)
1921 : {
1922 : ulong s, r;
1923 91710672 : int hi = p_sqrtu2((ulong*)n, &s, &r);
1924 91710672 : GEN S = utoi(s);
1925 91710672 : *pr = hi? uutoi(1,r): utoi(r);
1926 91710672 : return S;
1927 : }
1928 :
1929 : /* sqrt n[0], _dont_ assume n normalized */
1930 : static GEN
1931 2702241 : sqrtispec1_sh(GEN n, GEN *pr)
1932 : {
1933 : GEN S;
1934 2702241 : ulong r, s, u0 = uel(n,0);
1935 2702241 : int sh = bfffo(u0) & ~1UL;
1936 2702241 : if (sh) u0 <<= sh;
1937 2702241 : p_sqrtu1(&u0, &s, &r);
1938 : /* s^2 + r = u0, s < 2^(BIL/2). Rescale back:
1939 : * 2^(2k) n = S^2 + R
1940 : * so 2^(2k) n = (S - s0)^2 + (2*S*s0 - s0^2 + R), s0 = S mod 2^k. */
1941 2702241 : if (sh) {
1942 1630353 : int k = sh >> 1;
1943 1630353 : ulong s0 = s & ((1L<<k) - 1);
1944 1630353 : r += s * (s0<<1);
1945 1630353 : s >>= k;
1946 1630353 : r >>= sh;
1947 : }
1948 2702241 : S = utoi(s);
1949 2702241 : if (pr) *pr = utoi(r);
1950 2702241 : return S;
1951 : }
1952 :
1953 : /* sqrt n[0..1], _dont_ assume n normalized */
1954 : static GEN
1955 274824 : sqrtispec2_sh(GEN n, GEN *pr)
1956 : {
1957 : GEN S;
1958 274824 : ulong U[2], r, s, u0 = uel(n,0), u1 = uel(n,1);
1959 274824 : int hi, sh = bfffo(u0) & ~1UL;
1960 274824 : if (sh) {
1961 246855 : u0 = (u0 << sh) | (u1 >> (BITS_IN_LONG-sh));
1962 246855 : u1 <<= sh;
1963 : }
1964 274824 : U[0] = u0;
1965 274824 : U[1] = u1; hi = p_sqrtu2(U, &s, &r);
1966 : /* s^2 + R = u0|u1. Rescale back:
1967 : * 2^(2k) n = S^2 + R
1968 : * so 2^(2k) n = (S - s0)^2 + (2*S*s0 - s0^2 + R), s0 = S mod 2^k. */
1969 274824 : if (sh) {
1970 246855 : int k = sh >> 1;
1971 246855 : ulong s0 = s & ((1L<<k) - 1);
1972 : LOCAL_HIREMAINDER;
1973 : LOCAL_OVERFLOW;
1974 246855 : r = addll(r, mulll(s, (s0<<1)));
1975 246855 : if (overflow) hiremainder++;
1976 246855 : hiremainder += hi; /* + 0 or 1 */
1977 246855 : s >>= k;
1978 246855 : r = (r>>sh) | (hiremainder << (BITS_IN_LONG-sh));
1979 246855 : hi = (hiremainder & (1L<<sh));
1980 : }
1981 274824 : S = utoi(s);
1982 274824 : if (pr) *pr = hi? uutoi(1,r): utoi(r);
1983 274824 : return S;
1984 : }
1985 :
1986 : /* Let N = N[0..2n-1]. Return S (and set R) s.t S^2 + R = N, 0 <= R <= 2S
1987 : * Assume N normalized */
1988 : static GEN
1989 325922322 : sqrtispec(GEN N, long n, GEN *r)
1990 : {
1991 : GEN S, R, q, z, u;
1992 : long l, h;
1993 :
1994 325922322 : if (n == 1) return sqrtispec2(N, r);
1995 234211650 : l = n >> 1;
1996 234211650 : h = n - l; /* N = a3(h) | a2(h) | a1(l) | a0(l words) */
1997 234211650 : S = sqrtispec(N, h, &R); /* S^2 + R = a3|a2 */
1998 :
1999 234211650 : z = catii(LIMBS(R), NLIMBS(R), N + 2*h, l); /* = R | a1(l) */
2000 234211650 : q = dvmdii(z, shifti(S,1), &u);
2001 234211650 : z = catii(LIMBS(u), NLIMBS(u), N + n + h, l); /* = u | a0(l) */
2002 :
2003 234211650 : S = addshiftw(S, q, l);
2004 234211650 : R = subii(z, sqri(q));
2005 234211650 : if (signe(R) < 0)
2006 : {
2007 40027875 : GEN S2 = shifti(S,1);
2008 40027875 : R = addis(subiispec(LIMBS(S2),LIMBS(R), NLIMBS(S2),NLIMBS(R)), -1);
2009 40027875 : S = addis(S, -1);
2010 : }
2011 234211650 : *r = R; return S;
2012 : }
2013 :
2014 : /* Return S (and set R) s.t S^2 + R = N, 0 <= R <= 2S.
2015 : * As for dvmdii, R is last on stack and guaranteed to be gen_0 in case the
2016 : * remainder is 0. R = NULL is allowed. */
2017 : GEN
2018 3763578 : sqrtremi(GEN N, GEN *r)
2019 : {
2020 : pari_sp av;
2021 3763578 : GEN S, R, n = N+2;
2022 3763578 : long k, l2, ln = NLIMBS(N);
2023 : int sh;
2024 :
2025 3763578 : if (ln <= 2)
2026 : {
2027 2977620 : if (ln == 2) return sqrtispec2_sh(n, r);
2028 2702796 : if (ln == 1) return sqrtispec1_sh(n, r);
2029 555 : if (r) *r = gen_0;
2030 555 : return gen_0;
2031 : }
2032 785958 : av = avma;
2033 785958 : sh = bfffo(n[0]) >> 1;
2034 785958 : l2 = (ln + 1) >> 1;
2035 785958 : if (sh || (ln & 1)) { /* normalize n, so that n[0] >= 2^BIL / 4 */
2036 785280 : GEN s0, t = new_chunk(ln + 1);
2037 785280 : t[ln] = 0;
2038 785280 : if (sh)
2039 783450 : shift_left(t, n, 0,ln-1, 0, sh << 1);
2040 : else
2041 1830 : xmpn_copy(t, n, ln);
2042 785280 : S = sqrtispec(t, l2, &R); /* t normalized, 2 * l2 words */
2043 : /* Rescale back:
2044 : * 2^(2k) n = S^2 + R, k = sh + (ln & 1)*BIL/2
2045 : * so 2^(2k) n = (S - s0)^2 + (2*S*s0 - s0^2 + R), s0 = S mod 2^k. */
2046 785280 : k = sh + (ln & 1) * (BITS_IN_LONG/2);
2047 785280 : s0 = remi2n(S, k);
2048 785280 : R = addii(shifti(R,-1), mulii(s0, S));
2049 785280 : R = shifti(R, 1 - (k<<1));
2050 785280 : S = shifti(S, -k);
2051 : }
2052 : else
2053 678 : S = sqrtispec(n, l2, &R);
2054 :
2055 785958 : if (!r) { set_avma((pari_sp)S); return gc_INT(av, S); }
2056 722325 : *r = R; return gc_all(av, 2, &S, r);
2057 : }
2058 :
2059 : /* compute sqrt(|a|), assuming a != 0 */
2060 :
2061 : #if 1
2062 : GEN
2063 90924714 : sqrtr_abs(GEN x)
2064 : {
2065 90924714 : long l = lg(x) - 2, e = expo(x), er = e>>1;
2066 90924714 : GEN b, c, res = cgetg(2 + l, t_REAL);
2067 90924714 : res[1] = evalsigne(1) | evalexpo(er);
2068 90924714 : if (e&1) {
2069 40610028 : b = new_chunk(l << 1);
2070 40610028 : xmpn_copy(b, x+2, l);
2071 40610028 : xmpn_zero(b + l,l);
2072 40610028 : b = sqrtispec(b, l, &c);
2073 40610028 : xmpn_copy(res+2, b+2, l);
2074 40610028 : if (cmpii(c, b) > 0) roundr_up_ip(res, l+2);
2075 : } else {
2076 : ulong u;
2077 50314686 : b = new_chunk(2 + (l << 1));
2078 50314686 : shift_left(b+1, x+2, 0,l-1, 0, BITS_IN_LONG-1);
2079 50314686 : b[0] = uel(x,2)>>1;
2080 50314686 : xmpn_zero(b + l+1,l+1);
2081 50314686 : b = sqrtispec(b, l+1, &c);
2082 50314686 : xmpn_copy(res+2, b+2, l);
2083 50314686 : u = uel(b,l+2);
2084 50314686 : if ( u&HIGHBIT || (u == ~HIGHBIT && cmpii(c,b) > 0))
2085 24823002 : roundr_up_ip(res, l+2);
2086 : }
2087 90924714 : return gc_const((pari_sp)res, res);
2088 : }
2089 :
2090 : #else /* use t_REAL: currently much slower (quadratic division) */
2091 :
2092 : #ifdef LONG_IS_64BIT
2093 : /* 64 bits of b = sqrt(a[0] * 2^64 + a[1]) [ up to 1ulp ] */
2094 : static ulong
2095 : sqrtu2(ulong *a)
2096 : {
2097 : ulong c, b = dblmantissa( sqrt((double)a[0]) );
2098 : LOCAL_HIREMAINDER;
2099 : LOCAL_OVERFLOW;
2100 :
2101 : /* > 32 correct bits, 1 Newton iteration to reach 64 */
2102 : if (b <= a[0]) return HIGHBIT | (a[0] >> 1);
2103 : hiremainder = a[0]; c = divll(a[1], b);
2104 : return (addll(c, b) >> 1) | HIGHBIT;
2105 : }
2106 : /* 64 bits of sqrt(a[0] * 2^63) */
2107 : static ulong
2108 : sqrtu2_1(ulong *a)
2109 : {
2110 : ulong t[2];
2111 : t[0] = (a[0] >> 1);
2112 : t[1] = (a[0] << (BITS_IN_LONG-1)) | (a[1] >> 1);
2113 : return sqrtu2(t);
2114 : }
2115 : #else
2116 : /* 32 bits of sqrt(a[0] * 2^32) */
2117 : static ulong
2118 : sqrtu2(ulong *a) { return dblmantissa( sqrt((double)a[0]) ); }
2119 : /* 32 bits of sqrt(a[0] * 2^31) */
2120 : static ulong
2121 : sqrtu2_1(ulong *a) { return dblmantissa( sqrt(2. * a[0]) ); }
2122 : #endif
2123 :
2124 : GEN
2125 : sqrtr_abs(GEN x)
2126 : {
2127 : long l1, i, l = lg(x), ex = expo(x);
2128 : GEN a, t, y = cgetg(l, t_REAL);
2129 : pari_sp av, av0 = avma;
2130 :
2131 : a = rtor(x, lg2prec(l+1));
2132 : t = cgetg(l+1, t_REAL);
2133 : if (ex & 1) { /* odd exponent */
2134 : a[1] = evalsigne(1) | _evalexpo(1);
2135 : t[2] = (long)sqrtu2((ulong*)a + 2);
2136 : } else { /* even exponent */
2137 : a[1] = evalsigne(1) | _evalexpo(0);
2138 : t[2] = (long)sqrtu2_1((ulong*)a + 2);
2139 : }
2140 : t[1] = evalsigne(1) | _evalexpo(0);
2141 : for (i = 3; i <= l; i++) t[i] = 0;
2142 :
2143 : /* |x| = 2^(ex/2) a, t ~ sqrt(a) */
2144 : l--; l1 = 1; av = avma;
2145 : while (l1 < l) { /* let t := (t + a/t)/2 */
2146 : l1 <<= 1; if (l1 > l) l1 = l;
2147 : setlg(a, l1 + 2);
2148 : setlg(t, l1 + 2);
2149 : affrr(addrr(t, divrr(a,t)), t); shiftr_inplace(t, -1);
2150 : set_avma(av);
2151 : }
2152 : affrr(t,y); shiftr_inplace(y, (ex>>1));
2153 : return gc_const(av0, y);
2154 : }
2155 :
2156 : #endif
2157 :
2158 : /*******************************************************************
2159 : * *
2160 : * Base Conversion *
2161 : * *
2162 : *******************************************************************/
2163 :
2164 : static void
2165 736128 : convi_dac(GEN x, ulong l, ulong *res)
2166 : {
2167 736128 : pari_sp ltop=avma;
2168 : ulong m;
2169 : GEN x1,x2;
2170 736128 : if (l==1) { *res=itou(x); return; }
2171 349119 : m=l>>1;
2172 349119 : x1=dvmdii(x,powuu(1000000000UL,m),&x2);
2173 349119 : convi_dac(x1,l-m,res+m);
2174 349119 : convi_dac(x2,m,res);
2175 349119 : set_avma(ltop);
2176 : }
2177 :
2178 : /* convert integer --> base 10^9 [not memory clean] */
2179 : ulong *
2180 317011 : convi(GEN x, long *l)
2181 : {
2182 317011 : long lz, lx = lgefint(x);
2183 : ulong *z;
2184 317011 : if (lx == 3 && uel(x,2) < 1000000000UL) {
2185 279121 : z = (ulong*)new_chunk(1);
2186 279121 : *z = x[2];
2187 279121 : *l = 1; return z+1;
2188 : }
2189 37890 : lz = 1 + (long)bit_accuracy_mul(lx, LOG10_2/9);
2190 37890 : z = (ulong*)new_chunk(lz);
2191 37890 : convi_dac(x,(ulong)lz,z);
2192 67953 : while (z[lz-1]==0) lz--;
2193 37890 : *l=lz; return z+lz;
2194 : }
2195 :
|