Line data Source code
1 : /* Copyright (C) 2000-2004 The PARI group.
2 :
3 : This file is part of the PARI/GP package.
4 :
5 : PARI/GP is free software; you can redistribute it and/or modify it under the
6 : terms of the GNU General Public License as published by the Free Software
7 : Foundation; either version 2 of the License, or (at your option) any later
8 : version. It is distributed in the hope that it will be useful, but WITHOUT
9 : ANY WARRANTY WHATSOEVER.
10 :
11 : Check the License for details. You should have received a copy of it, along
12 : with the package; see the file 'COPYING'. If not, write to the Free Software
13 : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
14 :
15 : /***********************************************************************/
16 : /** **/
17 : /** GENERIC ALGORITHMS ON BLACKBOX GROUP **/
18 : /** **/
19 : /***********************************************************************/
20 : #include "pari.h"
21 : #include "paripriv.h"
22 : #undef pow /* AIX: pow(a,b) is a macro, wrongly expanded on grp->pow(a,b,c) */
23 :
24 : #define DEBUGLEVEL DEBUGLEVEL_bb_group
25 :
26 : /***********************************************************************/
27 : /** **/
28 : /** POWERING **/
29 : /** **/
30 : /***********************************************************************/
31 :
32 : /* return (n>>(i+1-l)) & ((1<<l)-1) */
33 : static ulong
34 8041038 : int_block(GEN n, long i, long l)
35 : {
36 8041038 : long q = divsBIL(i), r = remsBIL(i)+1, lr;
37 8040448 : GEN nw = int_W(n, q);
38 8040448 : ulong w = (ulong) *nw, w2;
39 8040448 : if (r>=l) return (w>>(r-l))&((1UL<<l)-1);
40 542286 : w &= (1UL<<r)-1; lr = l-r;
41 542286 : w2 = (ulong) *int_precW(nw); w2 >>= (BITS_IN_LONG-lr);
42 542286 : return (w<<lr)|w2;
43 : }
44 :
45 : /* assume n != 0, t_INT. Compute x^|n| using sliding window powering */
46 : static GEN
47 9139819 : sliding_window_powu(GEN x, ulong n, long e, void *E, GEN (*sqr)(void*,GEN),
48 : GEN (*mul)(void*,GEN,GEN))
49 : {
50 9139819 : long i, l = expu(n), u = (1UL<<(e-1));
51 9139660 : GEN tab = cgetg(1+u, t_VEC), x2 = sqr(E, x), z = NULL;
52 :
53 9142772 : gel(tab, 1) = x;
54 24333063 : for (i = 2; i <= u; i++) gel(tab,i) = mul(E, gel(tab,i-1), x2);
55 67580004 : while (l >= 0)
56 : {
57 : long w, v;
58 : GEN tw;
59 58579212 : if (e > l+1) e = l+1;
60 58579212 : w = (n>>(l+1-e)) & ((1UL<<e)-1); v = vals(w); l-=e;
61 58583402 : tw = gel(tab, 1 + (w>>(v+1)));
62 58583402 : if (!z) z = tw;
63 : else
64 : {
65 134880588 : for (i = 1; i <= e-v; i++) z = sqr(E, z);
66 49321833 : z = mul(E, z, tw);
67 : }
68 95604511 : for (i = 1; i <= v; i++) z = sqr(E, z);
69 107545106 : while (l >= 0)
70 : {
71 98641865 : if (n&(1UL<<l)) break;
72 49095982 : z = sqr(E, z); l--;
73 : }
74 : }
75 9000792 : return z;
76 : }
77 :
78 : /* assume n != 0, t_INT. Compute x^|n| using sliding window powering */
79 : static GEN
80 220811 : sliding_window_pow(GEN x, GEN n, long e, void *E, GEN (*sqr)(void*,GEN),
81 : GEN (*mul)(void*,GEN,GEN))
82 : {
83 : pari_sp av;
84 220811 : long i, l = expi(n), u = (1UL<<(e-1));
85 220808 : GEN tab = cgetg(1+u, t_VEC);
86 220809 : GEN x2 = sqr(E, x), z = NULL;
87 :
88 220753 : gel(tab, 1) = x;
89 2856852 : for (i=2; i<=u; i++) gel(tab,i) = mul(E, gel(tab,i-1), x2);
90 220823 : av = avma;
91 7681901 : while (l >= 0)
92 : {
93 : long w, v;
94 : GEN tw;
95 7435302 : if (e > l+1) e = l+1;
96 7435302 : w = int_block(n,l,e); v = vals(w); l-=e;
97 7436102 : tw = gel(tab, 1+(w>>(v+1)));
98 7436102 : if (!z) z = tw;
99 : else
100 : {
101 39133437 : for (i = 1; i <= e-v; i++) z = sqr(E, z);
102 7208350 : z = mul(E, z, tw);
103 : }
104 14162460 : for (i = 1; i <= v; i++) z = sqr(E, z);
105 20055633 : while (l >= 0)
106 : {
107 19809276 : if (gc_needed(av,1))
108 : {
109 3138 : if (DEBUGMEM>1) pari_warn(warnmem,"sliding_window_pow (%ld)", l);
110 3138 : z = gc_GEN(av, z);
111 : }
112 19809276 : if (int_bit(n,l)) break;
113 12580925 : z = sqr(E, z); l--;
114 : }
115 : }
116 246599 : return z;
117 : }
118 :
119 : /* assume n != 0, t_INT. Compute x^|n| using leftright binary powering */
120 : static GEN
121 96661065 : leftright_binary_powu(GEN x, ulong n, void *E, GEN (*sqr)(void*,GEN),
122 : GEN (*mul)(void*,GEN,GEN))
123 : {
124 96661065 : pari_sp av = avma;
125 : GEN y;
126 : int j;
127 :
128 96661065 : if (n == 1) return x;
129 96661065 : y = x; j = 1+bfffo(n);
130 : /* normalize, i.e set highest bit to 1 (we know n != 0) */
131 96661065 : n<<=j; j = BITS_IN_LONG-j;
132 : /* first bit is now implicit */
133 317870149 : for (; j; n<<=1,j--)
134 : {
135 221234273 : y = sqr(E,y);
136 221210469 : if (n & HIGHBIT) y = mul(E,y,x); /* first bit set: multiply by base */
137 221204022 : if (gc_needed(av,1))
138 : {
139 0 : if (DEBUGMEM>1) pari_warn(warnmem,"leftright_powu (%d)", j);
140 0 : y = gc_GEN(av, y);
141 : }
142 : }
143 96635876 : return y;
144 : }
145 :
146 : GEN
147 108836613 : gen_powu_i(GEN x, ulong n, void *E, GEN (*sqr)(void*,GEN),
148 : GEN (*mul)(void*,GEN,GEN))
149 : {
150 108836613 : if (n == 1) return x;
151 105799593 : if (n < 512)
152 96661276 : return leftright_binary_powu(x, n, E, sqr, mul);
153 : else
154 9138317 : return sliding_window_powu(x, n, n < (1UL<<25)? 2: 3, E, sqr, mul);
155 : }
156 :
157 : GEN
158 215637 : gen_powu(GEN x, ulong n, void *E, GEN (*sqr)(void*,GEN),
159 : GEN (*mul)(void*,GEN,GEN))
160 : {
161 215637 : pari_sp av = avma;
162 215637 : if (n == 1) return gcopy(x);
163 181315 : return gc_GEN(av, gen_powu_i(x,n,E,sqr,mul));
164 : }
165 :
166 : GEN
167 37948626 : gen_pow_i(GEN x, GEN n, void *E, GEN (*sqr)(void*,GEN),
168 : GEN (*mul)(void*,GEN,GEN))
169 : {
170 : long l, e;
171 37948626 : if (lgefint(n)==3) return gen_powu_i(x, uel(n,2), E, sqr, mul);
172 220812 : l = expi(n);
173 220809 : if (l<=64) e = 3;
174 161883 : else if (l<=160) e = 4;
175 81971 : else if (l<=384) e = 5;
176 20617 : else if (l<=896) e = 6;
177 10840 : else e = 7;
178 220809 : return sliding_window_pow(x, n, e, E, sqr, mul);
179 : }
180 :
181 : GEN
182 11454284 : gen_pow(GEN x, GEN n, void *E, GEN (*sqr)(void*,GEN),
183 : GEN (*mul)(void*,GEN,GEN))
184 : {
185 11454284 : pari_sp av = avma;
186 11454284 : return gc_GEN(av, gen_pow_i(x,n,E,sqr,mul));
187 : }
188 :
189 : /* assume n > 0. Compute x^n using left-right binary powering */
190 : GEN
191 1244038 : gen_powu_fold_i(GEN x, ulong n, void *E, GEN (*sqr)(void*,GEN),
192 : GEN (*msqr)(void*,GEN))
193 : {
194 1244038 : pari_sp av = avma;
195 : GEN y;
196 : int j;
197 :
198 1244038 : if (n == 1) return x;
199 1244038 : y = x; j = 1+bfffo(n);
200 : /* normalize, i.e set highest bit to 1 (we know n != 0) */
201 1244038 : n<<=j; j = BITS_IN_LONG-j;
202 : /* first bit is now implicit */
203 5734650 : for (; j; n<<=1,j--)
204 : {
205 4490647 : if (n & HIGHBIT) y = msqr(E,y); /* first bit set: multiply by base */
206 3224893 : else y = sqr(E,y);
207 4490626 : if (gc_needed(av,1))
208 : {
209 0 : if (DEBUGMEM>1) pari_warn(warnmem,"gen_powu_fold (%d)", j);
210 0 : y = gc_GEN(av, y);
211 : }
212 : }
213 1244003 : return y;
214 : }
215 : GEN
216 216184 : gen_powu_fold(GEN x, ulong n, void *E, GEN (*sqr)(void*,GEN),
217 : GEN (*msqr)(void*,GEN))
218 : {
219 216184 : pari_sp av = avma;
220 216184 : if (n == 1) return gcopy(x);
221 201442 : return gc_GEN(av, gen_powu_fold_i(x,n,E,sqr,msqr));
222 : }
223 :
224 : /* assume N != 0, t_INT. Compute x^|N| using left-right binary powering */
225 : GEN
226 687526 : gen_pow_fold_i(GEN x, GEN N, void *E, GEN (*sqr)(void*,GEN),
227 : GEN (*msqr)(void*,GEN))
228 : {
229 687526 : long ln = lgefint(N);
230 687526 : if (ln == 3) return gen_powu_fold_i(x, N[2], E, sqr, msqr);
231 : else
232 : {
233 152582 : GEN nd = int_MSW(N), y = x;
234 152582 : ulong n = *nd;
235 : long i;
236 : int j;
237 152582 : pari_sp av = avma;
238 :
239 152582 : if (n == 1)
240 7168 : j = 0;
241 : else
242 : {
243 145414 : j = 1+bfffo(n); /* < BIL */
244 : /* normalize, i.e set highest bit to 1 (we know n != 0) */
245 145414 : n <<= j; j = BITS_IN_LONG - j;
246 : }
247 : /* first bit is now implicit */
248 152582 : for (i=ln-2;;)
249 : {
250 50928687 : for (; j; n<<=1,j--)
251 : {
252 50143516 : if (n & HIGHBIT) y = msqr(E,y); /* first bit set: multiply by base */
253 35637811 : else y = sqr(E,y);
254 50155277 : if (gc_needed(av,1))
255 : {
256 0 : if (DEBUGMEM>1) pari_warn(warnmem,"gen_pow_fold (%ld,%d)", i,j);
257 0 : y = gc_GEN(av, y);
258 : }
259 : }
260 803587 : if (--i == 0) return y;
261 827239 : nd = int_precW(nd);
262 827239 : n = *nd; j = BITS_IN_LONG;
263 : }
264 : }
265 : }
266 : GEN
267 534944 : gen_pow_fold(GEN x, GEN n, void *E, GEN (*sqr)(void*,GEN),
268 : GEN (*msqr)(void*,GEN))
269 : {
270 534944 : pari_sp av = avma;
271 534944 : return gc_GEN(av, gen_pow_fold_i(x,n,E,sqr,msqr));
272 : }
273 :
274 : GEN
275 147 : gen_pow_init(GEN x, GEN n, long k, void *E, GEN (*sqr)(void*,GEN), GEN (*mul)(void*,GEN,GEN))
276 : {
277 147 : long i, j, l = expi(n);
278 147 : long m = 1UL<<(k-1);
279 147 : GEN R = cgetg(m+1, t_VEC);
280 147 : GEN x2 = sqr(E, x), y = gcopy(x);
281 644 : for(i = 1; i <= m; i++)
282 : {
283 497 : GEN C = cgetg(l+1, t_VEC);
284 497 : gel(C,1) = y;
285 27118 : for(j = 2; j <= l; j++)
286 26621 : gel(C,j) = sqr(E, gel(C,j-1));
287 497 : gel(R,i) = C;
288 497 : y = mul(E, y, x2);
289 : }
290 147 : return R;
291 : }
292 :
293 : GEN
294 53329 : gen_pow_table(GEN R, GEN n, void *E, GEN (*one)(void*), GEN (*mul)(void*,GEN,GEN))
295 : {
296 53329 : long e = expu(lg(R)-1) + 1;
297 53329 : long l = expi(n);
298 : long i, w;
299 53329 : GEN z = one(E), tw;
300 1241239 : for(i=0; i<=l; )
301 : {
302 1187910 : if (int_bit(n, i)==0) { i++; continue; }
303 605920 : if (i+e-1>l) e = l+1-i;
304 605920 : w = int_block(n,i+e-1,e);
305 605920 : tw = gmael(R, 1+(w>>1), i+1);
306 605920 : z = mul(E, z, tw);
307 605920 : i += e;
308 : }
309 53329 : return z;
310 : }
311 :
312 : GEN
313 28219502 : gen_powers(GEN x, long l, int use_sqr, void *E, GEN (*sqr)(void*,GEN),
314 : GEN (*mul)(void*,GEN,GEN), GEN (*one)(void*))
315 : {
316 : long i;
317 28219502 : GEN V = cgetg(l+2,t_VEC);
318 28219425 : gel(V,1) = one(E); if (l==0) return V;
319 28193597 : gel(V,2) = gcopy(x); if (l==1) return V;
320 12306214 : gel(V,3) = sqr(E,x);
321 12308741 : if (use_sqr)
322 41207207 : for(i = 4; i < l+2; i++)
323 31370721 : gel(V,i) = (i&1)? sqr(E,gel(V, (i+1)>>1))
324 31371334 : : mul(E,gel(V, i-1),x);
325 : else
326 6841818 : for(i = 4; i < l+2; i++)
327 4369579 : gel(V,i) = mul(E,gel(V,i-1),x);
328 12308112 : return V;
329 : }
330 :
331 : GEN
332 59319990 : producttree_scheme(long n)
333 : {
334 : GEN v, w;
335 : long i, j, k, u, l;
336 59319990 : if (n<=2) return mkvecsmall(n);
337 49542008 : u = expu(n-1);
338 49542014 : v = cgetg(n+1,t_VECSMALL);
339 49542021 : w = cgetg(n+1,t_VECSMALL);
340 49542376 : v[1] = n; l = 1;
341 160132660 : for (i=1; i<=u; i++)
342 : {
343 422632621 : for(j=1, k=1; j<=l; j++, k+=2)
344 : {
345 312042337 : long vj = v[j], v2 = vj>>1;
346 312042337 : w[k] = vj-v2;
347 312042337 : w[k+1] = v2;
348 : }
349 110590284 : swap(v,w); l<<=1;
350 : }
351 49542376 : fixlg(v, l+1); set_avma((pari_sp)v); return v;
352 : }
353 :
354 : GEN
355 61690816 : gen_product(GEN x, void *E, GEN (*mul)(void *,GEN,GEN))
356 : {
357 : pari_sp av;
358 61690816 : long i, k, l = lg(x);
359 : pari_timer ti;
360 : GEN y, v;
361 :
362 61690816 : if (l <= 2) return l == 1? gen_1: gcopy(gel(x,1));
363 57541669 : y = cgetg(l, t_VEC); av = avma;
364 57541794 : v = producttree_scheme(l-1);
365 57542304 : l = lg(v);
366 57542304 : if (DEBUGLEVEL>7) timer_start(&ti);
367 418914331 : for (k = i = 1; k < l; i += v[k++])
368 361374681 : gel(y,k) = v[k]==1? gel(x,i): mul(E, gel(x,i), gel(x,i+1));
369 165073432 : while (k > 2)
370 : {
371 107531396 : long n = k - 1;
372 107531396 : if (DEBUGLEVEL>7) timer_printf(&ti,"gen_product: remaining objects %ld",n);
373 411351818 : for (k = i = 1; i < n; i += 2) gel(y,k++) = mul(E, gel(y,i), gel(y,i+1));
374 107532589 : if (gc_needed(av,1)) gc_slice(av, y+1, k-1);
375 : }
376 57542036 : return gel(y,1);
377 : }
378 :
379 : /***********************************************************************/
380 : /** **/
381 : /** DISCRETE LOGARITHM **/
382 : /** **/
383 : /***********************************************************************/
384 : static GEN
385 63432204 : iter_rho(GEN x, GEN g, GEN q, GEN A, ulong h, void *E, const struct bb_group *grp)
386 : {
387 63432204 : GEN a = gel(A,1), b = gel(A,2), c = gel(A,3);
388 63432204 : switch((h | grp->hash(a)) % 3UL)
389 : {
390 21165806 : case 0: return mkvec3(grp->pow(E,a,gen_2), Fp_mulu(b,2,q), Fp_mulu(c,2,q));
391 21141335 : case 1: return mkvec3(grp->mul(E,a,x), addiu(b,1), c);
392 21125063 : case 2: return mkvec3(grp->mul(E,a,g), b, addiu(c,1));
393 : }
394 : return NULL; /* LCOV_EXCL_LINE */
395 : }
396 :
397 : /*Generic Pollard rho discrete log algorithm*/
398 : static GEN
399 49 : gen_Pollard_log(GEN x, GEN g, GEN q, void *E, const struct bb_group *grp)
400 : {
401 49 : pari_sp av=avma;
402 49 : GEN A, B, l, sqrt4q = sqrti(shifti(q,4));
403 49 : ulong i, h = 0, imax = itou_or_0(sqrt4q);
404 49 : if (!imax) imax = ULONG_MAX;
405 : do {
406 49 : rho_restart:
407 49 : A = B = mkvec3(x,gen_1,gen_0);
408 49 : i=0;
409 : do {
410 21144068 : if (i>imax)
411 : {
412 0 : h++;
413 0 : if (DEBUGLEVEL)
414 0 : pari_warn(warner,"changing Pollard rho hash seed to %ld",h);
415 0 : goto rho_restart;
416 : }
417 21144068 : A = iter_rho(x, g, q, A, h, E, grp);
418 21144068 : B = iter_rho(x, g, q, B, h, E, grp);
419 21144068 : B = iter_rho(x, g, q, B, h, E, grp);
420 21144068 : if (gc_needed(av,2))
421 : {
422 1474 : if(DEBUGMEM>1) pari_warn(warnmem,"gen_Pollard_log");
423 1474 : (void)gc_all(av, 2, &A, &B);
424 : }
425 21144068 : i++;
426 21144068 : } while (!grp->equal(gel(A,1), gel(B,1)));
427 49 : gel(A,2) = modii(gel(A,2), q);
428 49 : gel(B,2) = modii(gel(B,2), q);
429 49 : h++;
430 49 : } while (equalii(gel(A,2), gel(B,2)));
431 49 : l = Fp_div(Fp_sub(gel(B,3), gel(A,3),q),Fp_sub(gel(A,2), gel(B,2), q), q);
432 49 : return gc_INT(av, l);
433 : }
434 :
435 : /* compute a hash of g^(i-1), 1<=i<=n. Return [sorted hash, perm, g^-n] */
436 : GEN
437 750693 : gen_Shanks_init(GEN g, long n, void *E, const struct bb_group *grp)
438 : {
439 750693 : GEN p1 = g, G, perm, table = cgetg(n+1,t_VECSMALL);
440 750693 : pari_sp av=avma;
441 : long i;
442 750693 : table[1] = grp->hash(grp->pow(E,g,gen_0));
443 4300999 : for (i=2; i<=n; i++)
444 : {
445 3550306 : table[i] = grp->hash(p1);
446 3550306 : p1 = grp->mul(E,p1,g);
447 3550306 : if (gc_needed(av,2))
448 : {
449 0 : if(DEBUGMEM>1) pari_warn(warnmem,"gen_Shanks_log, baby = %ld", i);
450 0 : p1 = gc_upto(av, p1);
451 : }
452 : }
453 750693 : G = gc_upto(av, grp->pow(E,p1,gen_m1)); /* g^-n */
454 750693 : perm = vecsmall_indexsort(table);
455 750693 : table = vecsmallpermute(table,perm);
456 750693 : return mkvec4(table,perm,g,G);
457 : }
458 : /* T from gen_Shanks_init(g,n). Return v < n*N such that x = g^v or NULL */
459 : GEN
460 756069 : gen_Shanks(GEN T, GEN x, ulong N, void *E, const struct bb_group *grp)
461 : {
462 756069 : pari_sp av=avma;
463 756069 : GEN table = gel(T,1), perm = gel(T,2), g = gel(T,3), G = gel(T,4);
464 756069 : GEN p1 = x;
465 756069 : long n = lg(table)-1;
466 : ulong k;
467 4541246 : for (k=0; k<N; k++)
468 : { /* p1 = x G^k, G = g^-n */
469 4196032 : long h = grp->hash(p1), i = zv_search(table, h);
470 4196032 : if (i)
471 : {
472 411763 : do i--; while (i && table[i] == h);
473 410855 : for (i++; i <= n && table[i] == h; i++)
474 : {
475 410855 : GEN v = addiu(muluu(n,k), perm[i]-1);
476 410855 : if (grp->equal(grp->pow(E,g,v),x)) return gc_INT(av,v);
477 0 : if (DEBUGLEVEL)
478 0 : err_printf("gen_Shanks_log: false positive %lu, %lu\n", k,h);
479 : }
480 : }
481 3785177 : p1 = grp->mul(E,p1,G);
482 3785177 : if (gc_needed(av,2))
483 : {
484 0 : if(DEBUGMEM>1) pari_warn(warnmem,"gen_Shanks_log, k = %lu", k);
485 0 : p1 = gc_upto(av, p1);
486 : }
487 : }
488 345214 : return NULL;
489 : }
490 : /* Generic Shanks baby-step/giant-step algorithm. Return log_g(x), ord g = q.
491 : * One-shot: use gen_Shanks_init/log if many logs are desired; early abort
492 : * if log < sqrt(q) */
493 : static GEN
494 1548991 : gen_Shanks_log(GEN x, GEN g, GEN q, void *E, const struct bb_group *grp)
495 : {
496 1548991 : pari_sp av=avma, av1;
497 : long lbaby, i, k;
498 : GEN p1, table, giant, perm, ginv;
499 1548991 : p1 = sqrti(q);
500 1548991 : if (abscmpiu(p1,LGBITS) >= 0)
501 0 : pari_err_OVERFLOW("gen_Shanks_log [order too large]");
502 1548991 : lbaby = itos(p1)+1; table = cgetg(lbaby+1,t_VECSMALL);
503 1548991 : ginv = grp->pow(E,g,gen_m1);
504 1548991 : av1 = avma;
505 5116017 : for (p1=x, i=1;;i++)
506 : {
507 5116017 : if (grp->equal1(p1)) return gc_stoi(av, i-1);
508 4960217 : table[i] = grp->hash(p1); if (i==lbaby) break;
509 3567026 : p1 = grp->mul(E,p1,ginv);
510 3567026 : if (gc_needed(av1,2))
511 : {
512 6 : if(DEBUGMEM>1) pari_warn(warnmem,"gen_Shanks_log, baby = %ld", i);
513 6 : p1 = gc_upto(av1, p1);
514 : }
515 : }
516 1393191 : p1 = giant = gc_upto(av1, grp->mul(E,x,grp->pow(E, p1, gen_m1)));
517 1393191 : perm = vecsmall_indexsort(table);
518 1393191 : table = vecsmallpermute(table,perm);
519 1393191 : av1 = avma;
520 2180575 : for (k=1; k<= lbaby; k++)
521 : {
522 2180575 : long h = grp->hash(p1), i = zv_search(table, h);
523 2180575 : if (i)
524 : {
525 2786382 : while (table[i] == h && i) i--;
526 1393191 : for (i++; i <= lbaby && table[i] == h; i++)
527 : {
528 1393191 : GEN v = addiu(mulss(lbaby-1,k),perm[i]-1);
529 1393190 : if (grp->equal(grp->pow(E,g,v),x)) return gc_INT(av,v);
530 0 : if (DEBUGLEVEL)
531 0 : err_printf("gen_Shanks_log: false positive %ld, %lu\n", k,h);
532 : }
533 : }
534 787384 : p1 = grp->mul(E,p1,giant);
535 787384 : if (gc_needed(av1,2))
536 : {
537 0 : if(DEBUGMEM>1) pari_warn(warnmem,"gen_Shanks_log, k = %ld", k);
538 0 : p1 = gc_upto(av1, p1);
539 : }
540 : }
541 0 : retgc_const(av, cgetg(1, t_VEC)); /* no solution */
542 : }
543 :
544 : /*Generic discrete logarithme in a group of prime order p*/
545 : GEN
546 6088530 : gen_plog(GEN x, GEN g, GEN p, void *E, const struct bb_group *grp)
547 : {
548 6088530 : if (grp->easylog)
549 : {
550 6010568 : GEN e = grp->easylog(E, x, g, p);
551 6010568 : if (e) return e;
552 : }
553 2799002 : if (grp->equal1(x)) return gen_0;
554 2792688 : if (grp->equal(x,g)) return gen_1;
555 1549040 : if (expi(p)<32) return gen_Shanks_log(x,g,p,E,grp);
556 49 : return gen_Pollard_log(x, g, p, E, grp);
557 : }
558 :
559 : GEN
560 11614220 : get_arith_ZZM(GEN o)
561 : {
562 11614220 : if (!o) return NULL;
563 11614220 : switch(typ(o))
564 : {
565 3831212 : case t_INT:
566 3831212 : if (signe(o) > 0) return mkvec2(o, Z_factor(o));
567 8 : break;
568 1420543 : case t_MAT:
569 1420543 : if (is_Z_factorpos(o)) return mkvec2(factorback(o), o);
570 14 : break;
571 6362463 : case t_VEC:
572 6362463 : if (lg(o) == 3 && signe(gel(o,1)) > 0 && is_Z_factorpos(gel(o,2))) return o;
573 0 : break;
574 : }
575 24 : pari_err_TYPE("generic discrete logarithm (order factorization)",o);
576 : return NULL; /* LCOV_EXCL_LINE */
577 : }
578 : GEN
579 1425831 : get_arith_Z(GEN o)
580 : {
581 1425831 : if (!o) return NULL;
582 1425831 : switch(typ(o))
583 : {
584 1160243 : case t_INT:
585 1160243 : if (signe(o) > 0) return o;
586 7 : break;
587 14 : case t_MAT:
588 14 : o = factorback(o);
589 0 : if (typ(o) == t_INT && signe(o) > 0) return o;
590 0 : break;
591 265567 : case t_VEC:
592 265567 : if (lg(o) != 3) break;
593 265567 : o = gel(o,1);
594 265567 : if (typ(o) == t_INT && signe(o) > 0) return o;
595 0 : break;
596 : }
597 14 : pari_err_TYPE("generic discrete logarithm (order factorization)",o);
598 : return NULL; /* LCOV_EXCL_LINE */
599 : }
600 :
601 : /* Generic Pohlig-Hellman discrete logarithm: smallest integer n >= 0 such that
602 : * g^n=a. Assume ord(g) | ord; grp->easylog() is an optional trapdoor
603 : * function that catches easy logarithms */
604 : GEN
605 5022061 : gen_PH_log(GEN a, GEN g, GEN ord, void *E, const struct bb_group *grp)
606 : {
607 5022061 : pari_sp av = avma;
608 : GEN v, ginv, fa, ex;
609 : long i, j, l;
610 :
611 5022061 : if (grp->equal(g, a)) /* frequent special case */
612 1083849 : return grp->equal1(g)? gen_0: gen_1;
613 3938210 : if (grp->easylog)
614 : {
615 3888382 : GEN e = grp->easylog(E, a, g, ord);
616 3888353 : if (e) return e;
617 : }
618 2199398 : v = get_arith_ZZM(ord);
619 2199424 : ord= gel(v,1);
620 2199424 : fa = gel(v,2);
621 2199424 : ex = gel(fa,2);
622 2199424 : fa = gel(fa,1); l = lg(fa);
623 2199424 : ginv = grp->pow(E,g,gen_m1);
624 2199424 : v = cgetg(l, t_VEC);
625 6400392 : for (i = 1; i < l; i++)
626 : {
627 4200976 : GEN q = gel(fa,i), qj, gq, nq, ginv0, a0, t0;
628 4200976 : long e = itos(gel(ex,i));
629 4200976 : if (DEBUGLEVEL>5)
630 0 : err_printf("Pohlig-Hellman: DL mod %Ps^%ld\n",q,e);
631 4200976 : qj = new_chunk(e+1);
632 4200976 : gel(qj,0) = gen_1;
633 4200976 : gel(qj,1) = q;
634 6057751 : for (j = 2; j <= e; j++) gel(qj,j) = mulii(gel(qj,j-1), q);
635 4200976 : t0 = diviiexact(ord, gel(qj,e));
636 4200975 : a0 = grp->pow(E, a, t0);
637 4200976 : ginv0 = grp->pow(E, ginv, t0); /* ord(ginv0) | q^e */
638 4200976 : if (grp->equal1(ginv0)) { gel(v,i) = mkintmod(gen_0, gen_1); continue; }
639 4200969 : do gq = grp->pow(E,g, mulii(t0, gel(qj,--e))); while (grp->equal1(gq));
640 : /* ord(gq) = q */
641 4200962 : nq = gen_0;
642 4200962 : for (j = 0;; j++)
643 1856740 : { /* nq = sum_{i<j} b_i q^i */
644 6057702 : GEN b = grp->pow(E,a0, gel(qj,e-j));
645 : /* cheap early abort: wrong local order */
646 6057702 : if (j == 0 && !grp->equal1(grp->pow(E,b,q))) {
647 7 : retgc_const(av, cgetg(1, t_VEC));
648 : }
649 6057695 : b = gen_plog(b, gq, q, E, grp);
650 6057694 : if (typ(b) != t_INT) retgc_const(av, cgetg(1, t_VEC));
651 6057694 : nq = addii(nq, mulii(b, gel(qj,j)));
652 6057694 : if (j == e) break;
653 :
654 1856740 : a0 = grp->mul(E,a0, grp->pow(E,ginv0, b));
655 1856739 : ginv0 = grp->pow(E,ginv0, q);
656 : }
657 4200954 : gel(v,i) = mkintmod(nq, gel(qj,e+1));
658 : }
659 2199416 : return gc_INT(av, lift(chinese1_coprime_Z(v)));
660 : }
661 :
662 : /***********************************************************************/
663 : /** **/
664 : /** ORDER OF AN ELEMENT **/
665 : /** **/
666 : /***********************************************************************/
667 :
668 : static GEN
669 12023272 : rec_order(GEN a, GEN o, GEN m,
670 : void *E, const struct bb_group *grp, long x, long y)
671 : {
672 12023272 : pari_sp av = avma;
673 12023272 : if (grp->equal1(a)) return gen_1;
674 9849849 : if (x == y)
675 : {
676 5928088 : GEN p = gcoeff(m, x, 1);
677 5928088 : long i, e = itos(gcoeff(m, x, 2));
678 5928087 : if (e == 1) return icopy(p); /* a != 1 */
679 4039933 : for (i = 1; i < e; i++)
680 : {
681 2931314 : a = grp->pow(E, a, p);
682 2931301 : if (grp->equal1(a)) break;
683 : }
684 1911806 : return gc_INT(av, powiu(p, i));
685 : }
686 : else
687 : {
688 3921761 : GEN b, o1, o2, cof = gen_1;
689 3921761 : long i, z = (x+y)/2;
690 8902872 : for (i = x; i <= z; i++)
691 4980017 : cof = mulii(cof, powii(gcoeff(m, i, 1), gcoeff(m, i, 2)));
692 3922855 : b = grp->pow(E, a, cof);
693 3922453 : o1 = rec_order(b, diviiexact(o, cof), m, E, grp, z+1, y);
694 3922776 : b = grp->pow(E, a, o1);
695 3922401 : o2 = rec_order(b, diviiexact(o, o1), m, E, grp, x, z);
696 3922871 : return gc_INT(av, mulii(o1, o2));
697 : }
698 : }
699 :
700 : /*Find the exact order of a assuming a^o==1*/
701 : GEN
702 4187136 : gen_order(GEN a, GEN o, void *E, const struct bb_group *grp)
703 : {
704 4187136 : pari_sp av = avma;
705 : long l;
706 : GEN m;
707 :
708 4187136 : m = get_arith_ZZM(o);
709 4187131 : if (!m) pari_err_TYPE("gen_order [missing order]",a);
710 4187131 : o = gel(m,1);
711 4187131 : if (is_pm1(o)) return gc_const(av, gen_1);
712 4179615 : m = gel(m,2); l = lgcols(m);
713 4179615 : return gc_INT(av, rec_order(a, o, m, E, grp, 1, l-1));
714 : }
715 :
716 : /*Find the exact order of a assuming a^o==1, return [order,factor(order)] */
717 : GEN
718 5488 : gen_factored_order(GEN a, GEN o, void *E, const struct bb_group *grp)
719 : {
720 5488 : pari_sp av = avma;
721 : long i, l, ind;
722 : GEN m, F, P;
723 :
724 5488 : m = get_arith_ZZM(o);
725 5488 : if (!m) pari_err_TYPE("gen_factored_order [missing order]",a);
726 5488 : o = gel(m,1);
727 5488 : m = gel(m,2); l = lgcols(m);
728 5488 : P = cgetg(l, t_COL); ind = 1;
729 5488 : F = cgetg(l, t_COL);
730 18935 : for (i = l-1; i; i--)
731 : {
732 13447 : GEN t, y, p = gcoeff(m,i,1);
733 13447 : long j, e = itos(gcoeff(m,i,2));
734 13447 : if (l == 2) {
735 679 : t = gen_1;
736 679 : y = a;
737 : } else {
738 12768 : t = diviiexact(o, powiu(p,e));
739 12768 : y = grp->pow(E, a, t);
740 : }
741 13447 : if (grp->equal1(y)) o = t;
742 : else {
743 16233 : for (j = 1; j < e; j++)
744 : {
745 4655 : y = grp->pow(E, y, p);
746 4655 : if (grp->equal1(y)) break;
747 : }
748 12285 : gel(P,ind) = p;
749 12285 : gel(F,ind) = utoipos(j);
750 12285 : if (j < e) {
751 707 : if (j > 1) p = powiu(p, j);
752 707 : o = mulii(t, p);
753 : }
754 12285 : ind++;
755 : }
756 : }
757 5488 : setlg(P, ind); P = vecreverse(P);
758 5488 : setlg(F, ind); F = vecreverse(F);
759 5488 : return gc_GEN(av, mkvec2(o, mkmat2(P,F)));
760 : }
761 :
762 : /* E has order o[1], ..., or o[#o], draw random points until all solutions
763 : * but one are eliminated */
764 : GEN
765 980 : gen_select_order(GEN o, void *E, const struct bb_group *grp)
766 : {
767 980 : pari_sp ltop = avma, btop;
768 : GEN lastgood, so, vo;
769 980 : long lo = lg(o), nbo=lo-1;
770 980 : if (nbo == 1) return icopy(gel(o,1));
771 441 : so = ZV_indexsort(o); /* minimize max( o[i+1] - o[i] ) */
772 441 : vo = zero_zv(lo);
773 441 : lastgood = gel(o, so[nbo]);
774 441 : btop = avma;
775 : for(;;)
776 0 : {
777 441 : GEN lasto = gen_0;
778 441 : GEN P = grp->rand(E), t = mkvec(gen_0);
779 : long i;
780 567 : for (i = 1; i < lo; i++)
781 : {
782 567 : GEN newo = gel(o, so[i]);
783 567 : if (vo[i]) continue;
784 567 : t = grp->mul(E,t, grp->pow(E, P, subii(newo,lasto)));/*P^o[i]*/
785 567 : lasto = newo;
786 567 : if (!grp->equal1(t))
787 : {
788 483 : if (--nbo == 1) { set_avma(ltop); return icopy(lastgood); }
789 42 : vo[i] = 1;
790 : }
791 : else
792 84 : lastgood = lasto;
793 : }
794 0 : set_avma(btop);
795 : }
796 : }
797 :
798 : /*******************************************************************/
799 : /* */
800 : /* n-th ROOT */
801 : /* */
802 : /*******************************************************************/
803 : /* Assume l is prime. Return a generator of the l-th Sylow and set *zeta to an element
804 : * of order l.
805 : *
806 : * q = l^e*r, e>=1, (r,l)=1
807 : * UNCLEAN */
808 : static GEN
809 123205 : gen_lgener(GEN l, long e, GEN r,GEN *zeta, void *E, const struct bb_group *grp)
810 : {
811 123205 : const pari_sp av1 = avma;
812 : GEN m, m1;
813 : long i;
814 58002 : for (;; set_avma(av1))
815 : {
816 181207 : m1 = m = grp->pow(E, grp->rand(E), r);
817 181205 : if (grp->equal1(m)) continue;
818 217043 : for (i=1; i<e; i++)
819 : {
820 93839 : m = grp->pow(E,m,l);
821 93839 : if (grp->equal1(m)) break;
822 : }
823 140730 : if (i==e) break;
824 : }
825 123204 : *zeta = m; return m1;
826 : }
827 :
828 : /* Let G be a cyclic group of order o>1. Returns a (random) generator */
829 :
830 : GEN
831 15883 : gen_gener(GEN o, void *E, const struct bb_group *grp)
832 : {
833 15883 : pari_sp ltop = avma, av;
834 : long i, lpr;
835 15883 : GEN F, N, pr, z=NULL;
836 15883 : F = get_arith_ZZM(o);
837 15883 : N = gel(F,1); pr = gel(F,2); lpr = lgcols(pr);
838 15883 : av = avma;
839 :
840 51562 : for (i = 1; i < lpr; i++)
841 : {
842 35679 : GEN l = gcoeff(pr,i,1);
843 35679 : long e = itos(gcoeff(pr,i,2));
844 35679 : GEN r = diviiexact(N,powis(l,e));
845 35679 : GEN zetan, zl = gen_lgener(l,e,r,&zetan,E,grp);
846 35679 : z = i==1 ? zl: grp->mul(E,z,zl);
847 35679 : if (gc_needed(av,2))
848 : { /* n can have lots of prime factors*/
849 0 : if(DEBUGMEM>1) pari_warn(warnmem,"gen_gener");
850 0 : z = gc_upto(av, z);
851 : }
852 : }
853 15883 : return gc_upto(ltop, z);
854 : }
855 :
856 : /* solve x^l = a , l prime in G of order q.
857 : *
858 : * q = (l^e)*r, e >= 1, (r,l) = 1
859 : * y is not an l-th power, hence generates the l-Sylow of G
860 : * m = y^(q/l) != 1 */
861 : static GEN
862 87581 : gen_Shanks_sqrtl(GEN a, GEN l, long e, GEN r, GEN y, GEN m,void *E,
863 : const struct bb_group *grp)
864 : {
865 87581 : pari_sp av = avma;
866 : long k;
867 : GEN p1, u1, u2, v, w, z, dl;
868 :
869 87581 : (void)bezout(r,l,&u1,&u2);
870 87582 : v = grp->pow(E,a,u2);
871 87581 : w = grp->pow(E,v,l);
872 87583 : w = grp->mul(E,w,grp->pow(E,a,gen_m1));
873 118417 : while (!grp->equal1(w))
874 : {
875 30841 : k = 0;
876 30841 : p1 = w;
877 : do
878 : {
879 45770 : z = p1; p1 = grp->pow(E,p1,l);
880 45771 : k++;
881 45771 : } while(!grp->equal1(p1));
882 30842 : if (k==e) return gc_NULL(av);
883 30835 : dl = gen_plog(z,m,l,E,grp);
884 30835 : if (typ(dl) != t_INT) return gc_NULL(av);
885 30835 : dl = negi(dl);
886 30835 : p1 = grp->pow(E, grp->pow(E,y, dl), powiu(l,e-k-1));
887 30835 : m = grp->pow(E,m,dl);
888 30835 : e = k;
889 30835 : v = grp->mul(E,p1,v);
890 30834 : y = grp->pow(E,p1,l);
891 30834 : w = grp->mul(E,y,w);
892 30835 : if (gc_needed(av,1))
893 : {
894 0 : if(DEBUGMEM>1) pari_warn(warnmem,"gen_Shanks_sqrtl");
895 0 : (void)gc_all(av,4, &y,&v,&w,&m);
896 : }
897 : }
898 87572 : return gc_GEN(av, v);
899 : }
900 : /* Return one solution of x^n = a in a cyclic group of order q
901 : *
902 : * 1) If there is no solution, return NULL.
903 : *
904 : * 2) If there is a solution, there are exactly m of them [m = gcd(q-1,n)].
905 : * If zetan!=NULL, *zetan is set to a primitive m-th root of unity so that
906 : * the set of solutions is { x*zetan^k; k=0..m-1 }
907 : */
908 : GEN
909 90202 : gen_Shanks_sqrtn(GEN a, GEN n, GEN q, GEN *zetan, void *E, const struct bb_group *grp)
910 : {
911 90202 : pari_sp ltop = avma;
912 : GEN m, u1, u2, z;
913 : int is_1;
914 :
915 90202 : if (is_pm1(n))
916 : {
917 28 : if (zetan) *zetan = grp->pow(E,a,gen_0);
918 28 : return signe(n) < 0? grp->pow(E,a,gen_m1): gcopy(a);
919 : }
920 90174 : is_1 = grp->equal1(a);
921 90174 : if (is_1 && !zetan) return gcopy(a);
922 :
923 88286 : m = bezout(n,q,&u1,&u2);
924 88284 : z = grp->pow(E,a,gen_0);
925 88282 : if (!is_pm1(m))
926 : {
927 87512 : GEN F = Z_factor(m);
928 : long i, j, e;
929 : GEN r, zeta, y, l;
930 87517 : pari_sp av1 = avma;
931 175037 : for (i = nbrows(F); i; i--)
932 : {
933 87524 : l = gcoeff(F,i,1);
934 87524 : j = itos(gcoeff(F,i,2));
935 87524 : e = Z_pvalrem(q,l,&r);
936 87526 : y = gen_lgener(l,e,r,&zeta,E,grp);
937 87525 : if (zetan) z = grp->mul(E,z, grp->pow(E,y,powiu(l,e-j)));
938 87525 : if (!is_1) {
939 : do
940 : {
941 87581 : a = gen_Shanks_sqrtl(a,l,e,r,y,zeta,E,grp);
942 87583 : if (!a) return gc_NULL(ltop);
943 87576 : } while (--j);
944 : }
945 87520 : if (gc_needed(ltop,1))
946 : { /* n can have lots of prime factors*/
947 0 : if(DEBUGMEM>1) pari_warn(warnmem,"gen_Shanks_sqrtn");
948 0 : (void)gc_all(av1, zetan? 2: 1, &a, &z);
949 : }
950 : }
951 : }
952 88283 : if (!equalii(m, n))
953 784 : a = grp->pow(E,a,modii(u1,q));
954 88283 : if (zetan)
955 : {
956 1600 : *zetan = z;
957 1600 : (void)gc_all(ltop,2,&a,zetan);
958 : }
959 : else /* is_1 is 0: a was modified above -> gc_upto valid */
960 86683 : a = gc_upto(ltop, a);
961 88283 : return a;
962 : }
963 :
964 : /*******************************************************************/
965 : /* */
966 : /* structure of groups with pairing */
967 : /* */
968 : /*******************************************************************/
969 : /* return c = \prod_{p^2 | (N,d^2)} p^{v_p(N)} and factor(c); multiple of d2 */
970 : static GEN
971 146671 : d2_multiple(GEN N, GEN d)
972 : {
973 146671 : GEN P, E, Q = gel(Z_factor(gcdii(N,d)), 1);
974 146671 : long i, j, l = lg(Q);
975 146671 : P = cgetg(l, t_COL);
976 146671 : E = cgetg(l, t_COL);
977 271747 : for (i = 1, j = 1; i < l; i++)
978 : {
979 125076 : long v = Z_pval(N, gel(Q,i));
980 125076 : if (v <= 1) continue;
981 68558 : gel(P, j) = gel(Q,i);
982 68558 : gel(E, j) = utoipos(v); j++;
983 : }
984 146671 : if (j == 1) return NULL;
985 65100 : setlg(P,j); setlg(E,j);
986 65100 : return mkvec2(factorback2(P, E), mkmat2(P, E));
987 : }
988 :
989 : /* Return elementary divisors [d1, d2], d2 | d1, for group of order N.
990 : * We have d2 | d */
991 : GEN
992 146832 : gen_ellgroup(GEN N, GEN d, GEN *pm, void *E, const struct bb_group *grp,
993 : GEN pairorder(void *E, GEN P, GEN Q, GEN m, GEN F))
994 : {
995 146832 : pari_sp av = avma;
996 146832 : GEN N0, N1, F, fa0, L0, E0, g1 = gen_1, g2 = gen_1, mm = gen_1;
997 146832 : long n = 0, n0, j;
998 :
999 146832 : if (pm) *pm = gen_1;
1000 146832 : if (is_pm1(N)) return cgetg(1,t_VEC);
1001 146671 : F = d2_multiple(N, d); if (!F) { set_avma(av); return mkveccopy(N); }
1002 65100 : N0 = gel(F,1); fa0 = gel(F,2); /* N0 a multiple of d2, fa0 = factor(N0) */
1003 65100 : N1 = diviiexact(N, N0); /* N1 | d1 */
1004 65100 : L0 = gel(fa0, 1); n0 = lg(L0)-1; /* primes dividing N0 */
1005 65100 : E0 = ZV_to_nv(gel(fa0, 2)); /* ... and their exponents */
1006 : while (1)
1007 51705 : { /* g1 | (d1/N1), g2 | d2 */
1008 116805 : pari_sp av2 = avma;
1009 : GEN P, Q, s, t, m, mo;
1010 :
1011 116805 : P = grp->pow(E,grp->rand(E), N1);
1012 116805 : s = gen_order(P, F, E, grp); /* s | N0 */
1013 116805 : if (equalii(s, N0)) { set_avma(av); return mkveccopy(N); }
1014 :
1015 93404 : Q = grp->pow(E,grp->rand(E), N1);
1016 93404 : t = gen_order(Q, F, E, grp); /* t | N0 */
1017 93404 : if (equalii(t, N0)) { set_avma(av); return mkveccopy(N); }
1018 :
1019 82190 : m = lcmii(s, t); /* m | N0 */
1020 82190 : mo = mulii(m, pairorder(E, P, Q, m, F));
1021 :
1022 : /* For each prime l dividing N0, check whether P and Q
1023 : * generate all rational points of order a power of l */
1024 139077 : for (j = 1; j <= n0; j++)
1025 : {
1026 : GEN l;
1027 87372 : ulong e = uel(E0, j);
1028 87372 : if (e == 0) continue;
1029 85403 : l = gel(L0, j);
1030 85403 : if ((ulong)Z_pval(mo, l) == e)
1031 : {
1032 33341 : long vm = Z_pval(m, l);
1033 33341 : GEN l1 = powiu(l, vm);
1034 33341 : GEN l2 = powiu(l, e - vm);
1035 33341 : g1 = mulii(g1, l1);
1036 33341 : g2 = mulii(g2, l2);
1037 33341 : if (pm)
1038 33341 : mm = mulii(mm, cmpii(l1, l2) > 0 ? l1: l2);
1039 33341 : if (++n == n0)
1040 : { /* done with all primes l */
1041 : GEN g;
1042 30485 : if (equali1(g2)) { set_avma(av); return mkveccopy(N); }
1043 30296 : if (pm) *pm = g1;
1044 30296 : g = mkvec2(mulii(g1, N1), g2);
1045 30296 : if (!pm) return gc_GEN(av, g);
1046 30296 : *pm = mm; return gc_all(av, 2, &g, pm);
1047 : }
1048 2856 : uel(E0, j) = 0; /* done with this prime l */
1049 : }
1050 : }
1051 51705 : (void)gc_all(av2, pm ? 3: 2, &g1, &g2, &mm);
1052 : }
1053 : }
1054 :
1055 : GEN
1056 2716 : gen_ellgens(GEN D1, GEN d2, GEN m, void *E, const struct bb_group *grp,
1057 : GEN pairorder(void *E, GEN P, GEN Q, GEN m, GEN F))
1058 : {
1059 2716 : pari_sp ltop = avma, av;
1060 : GEN F, d1, dm;
1061 : GEN P, Q, d, s;
1062 2716 : F = get_arith_ZZM(D1);
1063 2716 : d1 = gel(F, 1), dm = diviiexact(d1,m);
1064 2716 : av = avma;
1065 : do
1066 : {
1067 6621 : set_avma(av);
1068 6621 : P = grp->rand(E);
1069 6621 : s = gen_order(P, F, E, grp);
1070 6621 : } while (!equalii(s, d1));
1071 2716 : av = avma;
1072 : do
1073 : {
1074 5480 : set_avma(av);
1075 5480 : Q = grp->rand(E);
1076 5480 : d = pairorder(E, grp->pow(E, P, dm), grp->pow(E, Q, dm), m, F);
1077 5480 : } while (!equalii(d, d2));
1078 2716 : return gc_GEN(ltop, mkvec2(P,Q));
1079 : }
|