Line data Source code
1 : /* Copyright (C) 2000-2003 The PARI group.
2 :
3 : This file is part of the PARI/GP package.
4 :
5 : PARI/GP is free software; you can redistribute it and/or modify it under the
6 : terms of the GNU General Public License as published by the Free Software
7 : Foundation; either version 2 of the License, or (at your option) any later
8 : version. It is distributed in the hope that it will be useful, but WITHOUT
9 : ANY WARRANTY WHATSOEVER.
10 :
11 : Check the License for details. You should have received a copy of it, along
12 : with the package; see the file 'COPYING'. If not, write to the Free Software
13 : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
14 :
15 : #include "pari.h"
16 : #include "paripriv.h"
17 :
18 : #define DEBUGLEVEL DEBUGLEVEL_galois
19 :
20 : /*************************************************************************/
21 : /** **/
22 : /** GALOIS CONJUGATES **/
23 : /** **/
24 : /*************************************************************************/
25 :
26 : static int
27 10927 : is2sparse(GEN x)
28 : {
29 10927 : long i, l = lg(x);
30 10927 : if (odd(l-3)) return 0;
31 32676 : for(i=3; i<l; i+=2)
32 21924 : if (signe(gel(x,i))) return 0;
33 10752 : return 1;
34 : }
35 :
36 : static GEN
37 35847 : galoisconj1(GEN nf)
38 : {
39 35847 : GEN x = get_nfpol(nf, &nf), f = nf? nf : x, y, z;
40 35847 : long i, lz, v = varn(x), nbmax;
41 35847 : pari_sp av = avma;
42 35847 : RgX_check_ZX(x, "nfgaloisconj");
43 35847 : nbmax = numberofconjugates(x, 2);
44 35846 : if (nbmax==1) retmkcol(pol_x(v));
45 11550 : if (nbmax==2 && is2sparse(x))
46 : {
47 10752 : GEN res = cgetg(3,t_COL);
48 10752 : gel(res,1) = deg1pol_shallow(gen_m1, gen_0, v);
49 10752 : gel(res,2) = pol_x(v);
50 10752 : return res;
51 : }
52 798 : x = leafcopy(x); setvarn(x, fetch_var_higher());
53 798 : z = nfroots(f, x); lz = lg(z);
54 798 : y = cgetg(lz, t_COL);
55 3857 : for (i = 1; i < lz; i++)
56 : {
57 3059 : GEN t = lift(gel(z,i));
58 3059 : if (typ(t) == t_POL) setvarn(t, v);
59 3059 : gel(y,i) = t;
60 : }
61 798 : (void)delete_var();
62 798 : return gerepileupto(av, y);
63 : }
64 :
65 : /*************************************************************************/
66 : /** **/
67 : /** GALOISCONJ4 **/
68 : /** **/
69 : /*************************************************************************/
70 : /*DEBUGLEVEL:
71 : 1: timing
72 : 2: outline
73 : 4: complete outline
74 : 6: detail
75 : 7: memory
76 : 9: complete detail
77 : */
78 : struct galois_borne {
79 : GEN l;
80 : long valsol;
81 : long valabs;
82 : GEN bornesol;
83 : GEN ladicsol;
84 : GEN ladicabs;
85 : GEN dis;
86 : };
87 : struct galois_lift {
88 : GEN T;
89 : GEN den;
90 : GEN p;
91 : GEN L;
92 : GEN Lden;
93 : long e;
94 : GEN Q;
95 : GEN TQ;
96 : struct galois_borne *gb;
97 : };
98 : struct galois_testlift {
99 : long n;
100 : long f;
101 : long g;
102 : GEN bezoutcoeff;
103 : GEN pauto;
104 : GEN C;
105 : GEN Cd;
106 : };
107 : struct galois_test { /* data for permutation test */
108 : GEN order; /* order of tests pour galois_test_perm */
109 : GEN borne, lborne; /* coefficient bounds */
110 : GEN ladic;
111 : GEN PV; /* NULL or vector of test matrices (Vmatrix) */
112 : GEN TM; /* transpose of M */
113 : GEN L; /* p-adic roots, known mod ladic */
114 : GEN M; /* vandermonde inverse */
115 : };
116 : /* result of the study of Frobenius degrees */
117 : enum ga_code {ga_all_normal=1,ga_ext_2=2,ga_non_wss=4,
118 : ga_all_nilpotent=8,ga_easy=16};
119 : struct galois_analysis {
120 : long p; /* prime to be lifted */
121 : long deg; /* degree of the lift */
122 : long ord;
123 : long l; /* l: prime number such that T is totally split mod l */
124 : long p4;
125 : long group;
126 : };
127 : struct galois_frobenius {
128 : long p;
129 : long fp;
130 : long deg;
131 : GEN Tmod;
132 : GEN psi;
133 : };
134 :
135 : /* #r = r1 + r2 */
136 : GEN
137 36 : embed_roots(GEN ro, long r1)
138 : {
139 36 : long r2 = lg(ro)-1-r1;
140 : GEN L;
141 36 : if (!r2) L = ro;
142 : else
143 : {
144 33 : long i,j, N = r1+2*r2;
145 33 : L = cgetg(N+1, t_VEC);
146 145 : for (i = 1; i <= r1; i++) gel(L,i) = gel(ro,i);
147 212 : for (j = i; j <= N; i++)
148 : {
149 179 : GEN z = gel(ro,i);
150 179 : gel(L,j++) = z;
151 179 : gel(L,j++) = mkcomplex(gel(z,1), gneg(gel(z,2)));
152 : }
153 : }
154 36 : return L;
155 : }
156 : GEN
157 135078 : embed_disc(GEN z, long r1, long prec)
158 : {
159 135078 : pari_sp av = avma;
160 135078 : GEN t = real_1(prec);
161 135083 : long i, j, n = lg(z)-1, r2 = n-r1;
162 230115 : for (i = 1; i < r1; i++)
163 : {
164 95032 : GEN zi = gel(z,i);
165 665042 : for (j = i+1; j <= r1; j++) t = gmul(t, gsub(zi, gel(z,j)));
166 : }
167 714104 : for (j = r1+1; j <= n; j++)
168 : {
169 579049 : GEN zj = gel(z,j), a = gel(zj,1), b = gel(zj,2), b2 = gsqr(b);
170 595992 : for (i = 1; i <= r1; i++)
171 : {
172 16982 : GEN zi = gel(z,i);
173 16982 : t = gmul(t, gadd(gsqr(gsub(zi, a)), b2));
174 : }
175 579010 : t = gmul(t, b);
176 : }
177 135055 : if (r2) t = gmul2n(t, r2);
178 135061 : if (r2 > 1)
179 : {
180 123305 : GEN T = real_1(prec);
181 577674 : for (i = r1+1; i < n; i++)
182 : {
183 454418 : GEN zi = gel(z,i), a = gel(zi,1), b = gel(zi,2);
184 1955320 : for (j = i+1; j <= n; j++)
185 : {
186 1500951 : GEN zj = gel(z,j), c = gel(zj,1), d = gel(zj,2);
187 1500951 : GEN f = gsqr(gsub(a,c)), g = gsqr(gsub(b,d)), h = gsqr(gadd(b,d));
188 1500904 : T = gmul(T, gmul(gadd(f,g), gadd(f,h)));
189 : }
190 : }
191 123256 : t = gmul(t, T);
192 : }
193 135052 : t = gsqr(t);
194 135065 : if (odd(r2)) t = gneg(t);
195 135074 : return gerepileupto(av, t);
196 : }
197 :
198 : /* Compute bound for the coefficients of automorphisms.
199 : * T a ZX, den a t_INT denominator or NULL */
200 : GEN
201 80522 : initgaloisborne(GEN T, GEN den, long prec, GEN *pL, GEN *pprep, GEN *pD)
202 : {
203 : GEN L, prep, nf, r;
204 : pari_timer ti;
205 :
206 80522 : if (DEBUGLEVEL>=4) timer_start(&ti);
207 80522 : T = get_nfpol(T, &nf);
208 80522 : r = nf ? nf_get_roots(nf) : NULL;
209 80522 : if (nf && precision(gel(r, 1)) >= prec)
210 36 : L = embed_roots(r, nf_get_r1(nf));
211 : else
212 80486 : L = QX_complex_roots(T, prec);
213 80521 : if (DEBUGLEVEL>=4) timer_printf(&ti,"roots");
214 80521 : prep = vandermondeinverseinit(L);
215 80521 : if (!den || pD)
216 : {
217 61904 : GEN res = RgV_prod(gabs(prep,prec));
218 61900 : GEN D = ZX_disc_all(T, 1 + gexpo(res)); /* +1 for safety */
219 61903 : if (pD) *pD = D;
220 61903 : if (!den) den = indexpartial(T,D);
221 : }
222 80524 : if (pprep) *pprep = prep;
223 80524 : *pL = L; return den;
224 : }
225 :
226 : /* ||| M ||| with respect to || x ||_oo, M t_MAT */
227 : GEN
228 114523 : matrixnorm(GEN M, long prec)
229 : {
230 114523 : long i,j,m, l = lg(M);
231 114523 : GEN B = real_0(prec);
232 :
233 114523 : if (l == 1) return B;
234 114523 : m = lgcols(M);
235 388292 : for (i = 1; i < m; i++)
236 : {
237 273772 : GEN z = gabs(gcoeff(M,i,1), prec);
238 2644706 : for (j = 2; j < l; j++) z = gadd(z, gabs(gcoeff(M,i,j), prec));
239 273758 : if (gcmp(z, B) > 0) B = z;
240 : }
241 114520 : return B;
242 : }
243 :
244 : static GEN
245 31248 : galoisborne(GEN T, GEN dn, struct galois_borne *gb, long d)
246 : {
247 : pari_sp ltop, av2;
248 : GEN borne, borneroots, bornetrace, borneabs;
249 : long prec;
250 : GEN L, M, prep, den;
251 : pari_timer ti;
252 31248 : const long step=3;
253 :
254 31248 : prec = nbits2prec(bit_accuracy(ZX_max_lg(T)));
255 31248 : den = initgaloisborne(T,dn,prec, &L,&prep,&gb->dis);
256 31247 : if (!dn) dn = den;
257 31247 : ltop = avma;
258 31247 : if (DEBUGLEVEL>=4) timer_start(&ti);
259 31247 : M = vandermondeinverse(L, RgX_gtofp(T, prec), den, prep);
260 31244 : if (DEBUGLEVEL>=4) timer_printf(&ti,"vandermondeinverse");
261 31244 : borne = matrixnorm(M, prec);
262 31243 : borneroots = gsupnorm(L, prec); /*t_REAL*/
263 31246 : bornetrace = mulur((2*step)*degpol(T)/d,
264 31246 : powru(borneroots, minss(degpol(T), step)));
265 31244 : borneroots = ceil_safe(gmul(borne, borneroots));
266 31242 : borneabs = ceil_safe(gmax_shallow(gmul(borne, bornetrace),
267 : powru(bornetrace, d)));
268 31244 : av2 = avma;
269 : /*We use d-1 test, so we must overlift to 2^BITS_IN_LONG*/
270 31244 : gb->valsol = logint(shifti(borneroots,2+BITS_IN_LONG), gb->l) + 1;
271 31247 : gb->valabs = logint(shifti(borneabs,2), gb->l) + 1;
272 31247 : gb->valabs = maxss(gb->valsol, gb->valabs);
273 31247 : if (DEBUGLEVEL >= 4)
274 0 : err_printf("GaloisConj: val1=%ld val2=%ld\n", gb->valsol, gb->valabs);
275 31247 : set_avma(av2);
276 31247 : gb->bornesol = gerepileuptoint(ltop, shifti(borneroots,1));
277 31245 : if (DEBUGLEVEL >= 9)
278 0 : err_printf("GaloisConj: Bound %Ps\n",borneroots);
279 31245 : gb->ladicsol = powiu(gb->l, gb->valsol);
280 31248 : gb->ladicabs = powiu(gb->l, gb->valabs);
281 31247 : return dn;
282 : }
283 :
284 : static GEN
285 29791 : makeLden(GEN L,GEN den, struct galois_borne *gb)
286 29791 : { return FpC_Fp_mul(L, den, gb->ladicsol); }
287 :
288 : /* Initialize the galois_lift structure */
289 : static void
290 29897 : initlift(GEN T, GEN den, ulong p, GEN L, GEN Lden, struct galois_borne *gb, struct galois_lift *gl)
291 : {
292 : pari_sp av;
293 : long e;
294 29897 : gl->gb = gb;
295 29897 : gl->T = T;
296 29897 : gl->den = is_pm1(den)? gen_1: den;
297 29897 : gl->p = utoipos(p);
298 29897 : gl->L = L;
299 29897 : gl->Lden = Lden;
300 29897 : av = avma;
301 29897 : e = logint(shifti(gb->bornesol, 2+BITS_IN_LONG), gl->p) + 1;
302 29896 : set_avma(av);
303 29896 : if (e < 2) e = 2;
304 29896 : gl->e = e;
305 29896 : gl->Q = powuu(p, e);
306 29897 : gl->TQ = FpX_red(T,gl->Q);
307 29895 : }
308 :
309 : /* Check whether f is (with high probability) a solution and compute its
310 : * permutation */
311 : static int
312 65844 : poltopermtest(GEN f, struct galois_lift *gl, GEN pf)
313 : {
314 : pari_sp av;
315 65844 : GEN fx, fp, B = gl->gb->bornesol;
316 : long i, j, ll;
317 283218 : for (i = 2; i < lg(f); i++)
318 225673 : if (abscmpii(gel(f,i),B) > 0)
319 : {
320 8299 : if (DEBUGLEVEL>=4) err_printf("GaloisConj: Solution too large.\n");
321 8299 : if (DEBUGLEVEL>=8) err_printf("f=%Ps\n borne=%Ps\n",f,B);
322 8299 : return 0;
323 : }
324 57545 : ll = lg(gl->L);
325 57545 : fp = const_vecsmall(ll-1, 1); /* left on stack */
326 57546 : av = avma;
327 283919 : for (i = 1; i < ll; i++, set_avma(av))
328 : {
329 226532 : fx = FpX_eval(f, gel(gl->L,i), gl->gb->ladicsol);
330 1245494 : for (j = 1; j < ll; j++)
331 1245330 : if (fp[j] && equalii(fx, gel(gl->Lden,j))) { pf[i]=j; fp[j]=0; break; }
332 226530 : if (j == ll) return 0;
333 : }
334 57389 : return 1;
335 : }
336 :
337 : static long
338 60494 : galoisfrobeniustest(GEN aut, struct galois_lift *gl, GEN frob)
339 : {
340 60494 : pari_sp av = avma;
341 60494 : GEN tlift = aut;
342 60494 : if (gl->den != gen_1) tlift = FpX_Fp_mul(tlift, gl->den, gl->Q);
343 60492 : tlift = FpX_center_i(tlift, gl->Q, shifti(gl->Q,-1));
344 60492 : return gc_long(av, poltopermtest(tlift, gl, frob));
345 : }
346 :
347 : static GEN
348 64299 : monoratlift(void *E, GEN S, GEN q)
349 : {
350 64299 : pari_sp ltop = avma;
351 64299 : struct galois_lift *gl = (struct galois_lift *) E;
352 64299 : GEN qm1 = sqrti(shifti(q,-2)), N = gl->Q;
353 64293 : GEN tlift = FpX_ratlift(S, q, qm1, qm1, gl->den);
354 64299 : if (tlift)
355 : {
356 27282 : pari_sp ltop = avma;
357 27282 : GEN frob = cgetg(lg(gl->L), t_VECSMALL);
358 27282 : if(DEBUGLEVEL>=4)
359 0 : err_printf("MonomorphismLift: trying early solution %Ps\n",tlift);
360 27282 : if (gl->den != gen_1)
361 23384 : tlift = FpX_Fp_mul(FpX_red(Q_muli_to_int(tlift, gl->den), N),
362 : Fp_inv(gl->den, N), N);
363 27280 : if (galoisfrobeniustest(tlift, gl, frob))
364 : {
365 27124 : if(DEBUGLEVEL>=4) err_printf("MonomorphismLift: true early solution.\n");
366 27124 : return gerepilecopy(ltop, tlift);
367 : }
368 157 : if(DEBUGLEVEL>=4) err_printf("MonomorphismLift: false early solution.\n");
369 : }
370 37174 : set_avma(ltop);
371 37175 : return NULL;
372 : }
373 :
374 : static GEN
375 31933 : monomorphismratlift(GEN P, GEN S, struct galois_lift *gl)
376 : {
377 : pari_timer ti;
378 31933 : if (DEBUGLEVEL >= 1) timer_start(&ti);
379 31933 : S = ZpX_ZpXQ_liftroot_ea(P,S,gl->T,gl->p, gl->e, (void*)gl, monoratlift);
380 31932 : if (DEBUGLEVEL >= 1) timer_printf(&ti, "monomorphismlift()");
381 31932 : return S;
382 : }
383 :
384 : /* Let T be a polynomial in Z[X] , p a prime number, S in Fp[X]/(T) so
385 : * that T(S)=0 [p,T]. Lift S in S_0 so that T(S_0)=0 [T,p^e]
386 : * Unclean stack */
387 : static GEN
388 31933 : automorphismlift(GEN S, struct galois_lift *gl)
389 : {
390 31933 : return monomorphismratlift(gl->T, S, gl);
391 : }
392 :
393 : static GEN
394 29895 : galoisdolift(struct galois_lift *gl)
395 : {
396 29895 : pari_sp av = avma;
397 29895 : GEN Tp = FpX_red(gl->T, gl->p);
398 29896 : GEN S = FpX_Frobenius(Tp, gl->p);
399 29896 : return gerepileupto(av, automorphismlift(S, gl));
400 : }
401 :
402 : static GEN
403 1449 : galoisdoliftn(struct galois_lift *gl, long e)
404 : {
405 1449 : pari_sp av = avma;
406 1449 : GEN Tp = FpX_red(gl->T, gl->p);
407 1449 : GEN S = FpXQ_autpow(FpX_Frobenius(Tp, gl->p), e, Tp, gl->p);
408 1449 : return gerepileupto(av, automorphismlift(S, gl));
409 : }
410 :
411 : static ulong
412 93 : findpsi(GEN D, ulong pstart, GEN P, GEN S, long o, GEN *Tmod, GEN *Tpsi)
413 : {
414 : forprime_t iter;
415 : ulong p;
416 93 : long n = degpol(P), i, j, g = n/o;
417 93 : GEN psi = cgetg(g+1, t_VECSMALL);
418 93 : u_forprime_init(&iter, pstart, ULONG_MAX);
419 2396 : while ((p = u_forprime_next(&iter)))
420 : {
421 : GEN F, Sp;
422 2396 : long gp = 0;
423 2396 : if (smodis(D, p) == 0)
424 235 : continue;
425 2161 : F = gel(Flx_factor(ZX_to_Flx(P, p), p), 1);
426 2161 : if (lg(F)-1 != g) continue;
427 738 : Sp = RgX_to_Flx(S, p);
428 1929 : for (j = 1; j <= g; j++)
429 : {
430 1765 : GEN Fj = gel(F, j);
431 1765 : GEN Sj = Flx_rem(Sp, Fj, p);
432 1765 : GEN A = Flxq_autpowers(Flx_Frobenius(Fj, p), o, Fj, p);
433 5945 : for (i = 1; i <= o; i++)
434 5371 : if (gequal(Sj, gel(A,i+1)))
435 : {
436 1191 : psi[j] = i; break;
437 : }
438 1765 : if (i > o) break;
439 1191 : if (gp==0 && i==1) gp=j;
440 : }
441 738 : if (gp && j > g)
442 : {
443 : /* Normalize result so that psi[l]=1 */
444 93 : if (gp!=1)
445 : {
446 15 : swap(gel(F,1),gel(F,gp));
447 15 : lswap(uel(psi,1),uel(psi,gp));
448 : }
449 93 : *Tpsi = Flv_Fl_div(psi,psi[g],o);
450 93 : *Tmod = FlxV_to_ZXV(F);
451 93 : return p;
452 : }
453 : }
454 0 : return 0;
455 : }
456 :
457 : static void
458 1757 : inittestlift(GEN plift, GEN Tmod, struct galois_lift *gl,
459 : struct galois_testlift *gt)
460 : {
461 : pari_timer ti;
462 1757 : gt->n = lg(gl->L) - 1;
463 1757 : gt->g = lg(Tmod) - 1;
464 1757 : gt->f = gt->n / gt->g;
465 1757 : gt->bezoutcoeff = bezout_lift_fact(gl->T, Tmod, gl->p, gl->e);
466 1757 : if (DEBUGLEVEL >= 2) timer_start(&ti);
467 1757 : gt->pauto = FpXQ_autpowers(plift, gt->f-1, gl->TQ, gl->Q);
468 1757 : if (DEBUGLEVEL >= 2) timer_printf(&ti, "Frobenius power");
469 1757 : }
470 :
471 : /* Explanation of the intheadlong technique:
472 : * Let C be a bound, B = BITS_IN_LONG, M > C*2^B a modulus and 0 <= a_i < M for
473 : * i=1,...,n where n < 2^B. We want to test if there exist k,l, |k| < C < M/2^B,
474 : * such that sum a_i = k + l*M
475 : * We write a_i*2^B/M = b_i+c_i with b_i integer and 0<=c_i<1, so that
476 : * sum b_i - l*2^B = k*2^B/M - sum c_i
477 : * Since -1 < k*2^B/M < 1 and 0<=c_i<1, it follows that
478 : * -n-1 < sum b_i - l*2^B < 1 i.e. -n <= sum b_i -l*2^B <= 0
479 : * So we compute z = - sum b_i [mod 2^B] and check if 0 <= z <= n. */
480 :
481 : /* Assume 0 <= x < mod. */
482 : static ulong
483 1360230 : intheadlong(GEN x, GEN mod)
484 : {
485 1360230 : pari_sp av = avma;
486 1360230 : long res = (long) itou(divii(shifti(x,BITS_IN_LONG),mod));
487 1360229 : return gc_long(av,res);
488 : }
489 : static GEN
490 58317 : vecheadlong(GEN W, GEN mod)
491 : {
492 58317 : long i, l = lg(W);
493 58317 : GEN V = cgetg(l, t_VECSMALL);
494 1379662 : for(i=1; i<l; i++) V[i] = intheadlong(gel(W,i), mod);
495 58316 : return V;
496 : }
497 : static GEN
498 5705 : matheadlong(GEN W, GEN mod)
499 : {
500 5705 : long i, l = lg(W);
501 5705 : GEN V = cgetg(l,t_MAT);
502 64022 : for(i=1; i<l; i++) gel(V,i) = vecheadlong(gel(W,i), mod);
503 5705 : return V;
504 : }
505 : static ulong
506 38884 : polheadlong(GEN P, long n, GEN mod)
507 : {
508 38884 : return (lg(P)>n+2)? intheadlong(gel(P,n+2),mod): 0;
509 : }
510 :
511 : #define headlongisint(Z,N) (-(ulong)(Z)<=(ulong)(N))
512 :
513 : static long
514 2002 : frobeniusliftall(GEN sg, long el, GEN *psi, struct galois_lift *gl,
515 : struct galois_testlift *gt, GEN frob)
516 : {
517 2002 : pari_sp av, ltop2, ltop = avma;
518 2002 : long i,j,k, c = lg(sg)-1, n = lg(gl->L)-1, m = gt->g, d = m / c;
519 : GEN pf, u, v, C, Cd, SG, cache;
520 2002 : long N1, N2, R1, Ni, ord = gt->f, c_idx = gt->g-1;
521 : ulong headcache;
522 2002 : long hop = 0;
523 : GEN NN, NQ;
524 : pari_timer ti;
525 :
526 2002 : *psi = pf = cgetg(m, t_VECSMALL);
527 2002 : ltop2 = avma;
528 2002 : NN = diviiexact(mpfact(m), mului(c, powiu(mpfact(d), c)));
529 2002 : if (DEBUGLEVEL >= 4)
530 0 : err_printf("GaloisConj: I will try %Ps permutations\n", NN);
531 2002 : N1=10000000;
532 2002 : NQ=divis_rem(NN,N1,&R1);
533 2002 : if (abscmpiu(NQ,1000000000)>0)
534 : {
535 0 : pari_warn(warner,"Combinatorics too hard : would need %Ps tests!\n"
536 : "I will skip it, but it may induce an infinite loop",NN);
537 0 : *psi = NULL; return gc_long(ltop,0);
538 : }
539 2002 : N2=itos(NQ); if(!N2) N1=R1;
540 2002 : if (DEBUGLEVEL>=4) timer_start(&ti);
541 2002 : set_avma(ltop2);
542 2002 : C = gt->C;
543 2002 : Cd= gt->Cd;
544 2002 : v = FpXQ_mul(gel(gt->pauto, 1+el%ord), gel(gt->bezoutcoeff, m),gl->TQ,gl->Q);
545 2002 : if (gl->den != gen_1) v = FpX_Fp_mul(v, gl->den, gl->Q);
546 2002 : SG = cgetg(lg(sg),t_VECSMALL);
547 6552 : for(i=1; i<lg(SG); i++) SG[i] = (el*sg[i])%ord + 1;
548 2002 : cache = cgetg(m+1,t_VECSMALL); cache[m] = polheadlong(v,1,gl->Q);
549 2002 : headcache = polheadlong(v,2,gl->Q);
550 5306 : for (i = 1; i < m; i++) pf[i] = 1 + i/d;
551 2002 : av = avma;
552 2002 : for (Ni = 0, i = 0; ;i++)
553 : {
554 301525 : for (j = c_idx ; j > 0; j--)
555 : {
556 236903 : long h = SG[pf[j]];
557 236903 : if (!mael(C,h,j))
558 : {
559 5025 : pari_sp av3 = avma;
560 5025 : GEN r = FpXQ_mul(gel(gt->pauto,h), gel(gt->bezoutcoeff,j),gl->TQ,gl->Q);
561 5025 : if (gl->den != gen_1) r = FpX_Fp_mul(r, gl->den, gl->Q);
562 5025 : gmael(C,h,j) = gclone(r);
563 5025 : mael(Cd,h,j) = polheadlong(r,1,gl->Q);
564 5025 : set_avma(av3);
565 : }
566 236903 : uel(cache,j) = uel(cache,j+1)+umael(Cd,h,j);
567 : }
568 64622 : if (headlongisint(uel(cache,1),n))
569 : {
570 3360 : ulong head = headcache;
571 33215 : for (j = 1; j < m; j++) head += polheadlong(gmael(C,SG[pf[j]],j),2,gl->Q);
572 3360 : if (headlongisint(head,n))
573 : {
574 1722 : u = v;
575 4193 : for (j = 1; j < m; j++) u = ZX_add(u, gmael(C,SG[pf[j]],j));
576 1722 : u = FpX_center_i(FpX_red(u, gl->Q), gl->Q, shifti(gl->Q,-1));
577 1722 : if (poltopermtest(u, gl, frob))
578 : {
579 1715 : if (DEBUGLEVEL >= 4)
580 : {
581 0 : timer_printf(&ti, "");
582 0 : err_printf("GaloisConj: %d hops on %Ps tests\n",hop,addis(mulss(Ni,N1),i));
583 : }
584 1715 : return gc_long(ltop2,1);
585 : }
586 7 : if (DEBUGLEVEL >= 4) err_printf("M");
587 : }
588 1638 : else hop++;
589 : }
590 62907 : if (DEBUGLEVEL >= 4 && i % maxss(N1/20, 1) == 0)
591 0 : timer_printf(&ti, "GaloisConj:Testing %Ps", addis(mulss(Ni,N1),i));
592 62907 : set_avma(av);
593 62907 : if (i == N1-1)
594 : {
595 287 : if (Ni==N2-1) N1 = R1;
596 287 : if (Ni==N2) break;
597 0 : Ni++; i = 0;
598 0 : if (DEBUGLEVEL>=4) timer_start(&ti);
599 : }
600 170979 : for (j = 2; j < m && pf[j-1] >= pf[j]; j++)
601 : /*empty*/; /* to kill clang Warning */
602 101749 : for (k = 1; k < j-k && pf[k] != pf[j-k]; k++) { lswap(pf[k], pf[j-k]); }
603 109730 : for (k = j - 1; pf[k] >= pf[j]; k--)
604 : /*empty*/;
605 62620 : lswap(pf[j], pf[k]); c_idx = j;
606 : }
607 287 : if (DEBUGLEVEL>=4) err_printf("GaloisConj: not found, %d hops \n",hop);
608 287 : *psi = NULL; return gc_long(ltop,0);
609 : }
610 :
611 : /* Compute the test matrix for the i-th line of V. Clone. */
612 : static GEN
613 5705 : Vmatrix(long i, struct galois_test *td)
614 : {
615 5705 : pari_sp av = avma;
616 5705 : GEN m = gclone( matheadlong(FpC_FpV_mul(td->L, row(td->M,i), td->ladic), td->ladic));
617 5705 : set_avma(av); return m;
618 : }
619 :
620 : /* Initialize galois_test */
621 : static void
622 5467 : inittest(GEN L, GEN M, GEN borne, GEN ladic, struct galois_test *td)
623 : {
624 5467 : long i, n = lg(L)-1;
625 5467 : GEN p = cgetg(n+1, t_VECSMALL);
626 5467 : if (DEBUGLEVEL >= 8) err_printf("GaloisConj: Init Test\n");
627 5467 : td->order = p;
628 44828 : for (i = 1; i <= n-2; i++) p[i] = i+2;
629 5467 : p[n-1] = 1; p[n] = 2;
630 5467 : td->borne = borne;
631 5467 : td->lborne = subii(ladic, borne);
632 5467 : td->ladic = ladic;
633 5467 : td->L = L;
634 5467 : td->M = M;
635 5467 : td->TM = shallowtrans(M);
636 5467 : td->PV = zero_zv(n);
637 5467 : gel(td->PV, 2) = Vmatrix(2, td);
638 5467 : }
639 :
640 : /* Free clones stored inside galois_test */
641 : static void
642 5467 : freetest(struct galois_test *td)
643 : {
644 : long i;
645 55762 : for (i = 1; i < lg(td->PV); i++)
646 50295 : if (td->PV[i]) { gunclone(gel(td->PV,i)); td->PV[i] = 0; }
647 5467 : }
648 :
649 : /* Check if the integer P seen as a p-adic number is close to an integer less
650 : * than td->borne in absolute value */
651 : static long
652 94570 : padicisint(GEN P, struct galois_test *td)
653 : {
654 94570 : pari_sp ltop = avma;
655 94570 : GEN U = modii(P, td->ladic);
656 94570 : long r = cmpii(U, td->borne) <= 0 || cmpii(U, td->lborne) >= 0;
657 94570 : return gc_long(ltop, r);
658 : }
659 :
660 : /* Check if the permutation pf is valid according to td.
661 : * If not, update td to make subsequent test faster (hopefully) */
662 : static long
663 117334 : galois_test_perm(struct galois_test *td, GEN pf)
664 : {
665 117334 : pari_sp av = avma;
666 117334 : long i, j, n = lg(td->L)-1;
667 117334 : GEN V, P = NULL;
668 213101 : for (i = 1; i < n; i++)
669 : {
670 206234 : long ord = td->order[i];
671 206234 : GEN PW = gel(td->PV, ord);
672 206234 : if (PW)
673 : {
674 111664 : ulong head = umael(PW,1,pf[1]);
675 7121072 : for (j = 2; j <= n; j++) head += umael(PW,j,pf[j]);
676 111664 : if (!headlongisint(head,n)) break;
677 : } else
678 : {
679 94570 : if (!P) P = vecpermute(td->L, pf);
680 94570 : V = FpV_dotproduct(gel(td->TM,ord), P, td->ladic);
681 94570 : if (!padicisint(V, td)) {
682 238 : gel(td->PV, ord) = Vmatrix(ord, td);
683 238 : if (DEBUGLEVEL >= 4) err_printf("M");
684 238 : break;
685 : }
686 : }
687 : }
688 117334 : if (i == n) return gc_long(av,1);
689 110467 : if (DEBUGLEVEL >= 4) err_printf("%d.", i);
690 110467 : if (i > 1)
691 : {
692 742 : long z = td->order[i];
693 1526 : for (j = i; j > 1; j--) td->order[j] = td->order[j-1];
694 742 : td->order[1] = z;
695 742 : if (DEBUGLEVEL >= 8) err_printf("%Ps", td->order);
696 : }
697 110467 : return gc_long(av,0);
698 : }
699 : /*Compute a*b/c when a*b will overflow*/
700 : static long
701 0 : muldiv(long a,long b,long c)
702 : {
703 0 : return (long)((double)a*(double)b/c);
704 : }
705 :
706 : /* F = cycle decomposition of sigma,
707 : * B = cycle decomposition of cl(tau).
708 : * Check all permutations pf who can possibly correspond to tau, such that
709 : * tau*sigma*tau^-1 = sigma^s and tau^d = sigma^t, where d = ord cl(tau)
710 : * x: vector of choices,
711 : * G: vector allowing linear access to elts of F.
712 : * Choices multiple of e are not changed. */
713 : static GEN
714 8925 : testpermutation(GEN F, GEN B, GEN x, long s, long e, long cut,
715 : struct galois_test *td)
716 : {
717 8925 : pari_sp av, avm = avma;
718 : long a, b, c, d, n, p1, p2, p3, p4, p5, p6, l1, l2, N1, N2, R1;
719 8925 : long i, j, cx, hop = 0, start = 0;
720 : GEN pf, ar, G, W, NN, NQ;
721 : pari_timer ti;
722 8925 : if (DEBUGLEVEL >= 1) timer_start(&ti);
723 8925 : a = lg(F)-1; b = lg(gel(F,1))-1;
724 8925 : c = lg(B)-1; d = lg(gel(B,1))-1;
725 8925 : n = a*b;
726 8925 : s = (b+s) % b;
727 8925 : pf = cgetg(n+1, t_VECSMALL);
728 8925 : av = avma;
729 8925 : ar = cgetg(a+2, t_VECSMALL); ar[a+1]=0;
730 8925 : G = cgetg(a+1, t_VECSMALL);
731 8925 : W = gel(td->PV, td->order[n]);
732 58744 : for (cx=1, i=1, j=1; cx <= a; cx++, i++)
733 : {
734 49819 : gel(G,cx) = gel(F, coeff(B,i,j));
735 49819 : if (i == d) { i = 0; j++; }
736 : }
737 8925 : NN = divis(powuu(b, c * (d - d/e)), cut);
738 8925 : if (DEBUGLEVEL>=4) err_printf("GaloisConj: I will try %Ps permutations\n", NN);
739 8925 : N1 = 1000000;
740 8925 : NQ = divis_rem(NN,N1,&R1);
741 8925 : if (abscmpiu(NQ,100000000)>0)
742 : {
743 0 : set_avma(avm);
744 0 : pari_warn(warner,"Combinatorics too hard: would need %Ps tests!\n"
745 : "I'll skip it but you will get a partial result...",NN);
746 0 : return identity_perm(n);
747 : }
748 8925 : N2 = itos(NQ);
749 11277 : for (l2 = 0; l2 <= N2; l2++)
750 : {
751 8925 : long nbiter = (l2<N2) ? N1: R1;
752 8925 : if (DEBUGLEVEL >= 2 && N2) err_printf("%d%% ", muldiv(l2,100,N2));
753 10148203 : for (l1 = 0; l1 < nbiter; l1++)
754 : {
755 10145851 : if (start)
756 : {
757 18145582 : for (i=1, j=e; i < a;)
758 : {
759 18145582 : if ((++(x[i])) != b) break;
760 8008656 : x[i++] = 0;
761 8008656 : if (i == j) { i++; j += e; }
762 : }
763 : }
764 8925 : else { start=1; i = a-1; }
765 : /* intheadlong test: overflow in + is OK, we compute mod 2^BIL */
766 46092454 : for (p1 = i+1, p5 = p1%d - 1 ; p1 >= 1; p1--, p5--) /* p5 = (p1%d) - 1 */
767 : {
768 : GEN G1, G6;
769 35946603 : ulong V = 0;
770 35946603 : if (p5 == - 1) { p5 = d - 1; p6 = p1 + 1 - d; } else p6 = p1 + 1;
771 35946603 : G1 = gel(G,p1); G6 = gel(G,p6);
772 35946603 : p4 = p5 ? x[p1-1] : 0;
773 109414000 : for (p2 = 1+p4, p3 = 1 + x[p1]; p2 <= b; p2++)
774 : {
775 73467397 : V += umael(W,uel(G6,p3),uel(G1,p2));
776 73467397 : p3 += s; if (p3 > b) p3 -= b;
777 : }
778 35946603 : p3 = 1 + x[p1] - s; if (p3 <= 0) p3 += b;
779 50260116 : for (p2 = p4; p2 >= 1; p2--)
780 : {
781 14313513 : V += umael(W,uel(G6,p3),uel(G1,p2));
782 14313513 : p3 -= s; if (p3 <= 0) p3 += b;
783 : }
784 35946603 : uel(ar,p1) = uel(ar,p1+1) + V;
785 : }
786 10145851 : if (!headlongisint(uel(ar,1),n)) continue;
787 :
788 : /* intheadlong succeeds. Full computation */
789 3510899 : for (p1=1, p5=d; p1 <= a; p1++, p5++)
790 : {
791 3393859 : if (p5 == d) { p5 = 0; p4 = 0; } else p4 = x[p1-1];
792 3393859 : if (p5 == d-1) p6 = p1+1-d; else p6 = p1+1;
793 9516612 : for (p2 = 1+p4, p3 = 1 + x[p1]; p2 <= b; p2++)
794 : {
795 6122753 : pf[mael(G,p1,p2)] = mael(G,p6,p3);
796 6122753 : p3 += s; if (p3 > b) p3 -= b;
797 : }
798 3393859 : p3 = 1 + x[p1] - s; if (p3 <= 0) p3 += b;
799 4417882 : for (p2 = p4; p2 >= 1; p2--)
800 : {
801 1024023 : pf[mael(G,p1,p2)] = mael(G,p6,p3);
802 1024023 : p3 -= s; if (p3 <= 0) p3 += b;
803 : }
804 : }
805 117040 : if (galois_test_perm(td, pf))
806 : {
807 6573 : if (DEBUGLEVEL >= 1)
808 : {
809 0 : GEN nb = addis(mulss(l2,N1),l1);
810 0 : timer_printf(&ti, "testpermutation(%Ps)", nb);
811 0 : if (DEBUGLEVEL >= 2 && hop)
812 0 : err_printf("GaloisConj: %d hop over %Ps iterations\n", hop, nb);
813 : }
814 6573 : set_avma(av); return pf;
815 : }
816 110467 : hop++;
817 : }
818 : }
819 2352 : if (DEBUGLEVEL >= 1)
820 : {
821 0 : timer_printf(&ti, "testpermutation(%Ps)", NN);
822 0 : if (DEBUGLEVEL >= 2 && hop)
823 0 : err_printf("GaloisConj: %d hop over %Ps iterations\n", hop, NN);
824 : }
825 2352 : return gc_NULL(avm);
826 : }
827 :
828 : /* List of subgroups of (Z/mZ)^* whose order divide o, and return the list
829 : * of their elements, sorted by increasing order */
830 : static GEN
831 1764 : listznstarelts(long m, long o)
832 : {
833 1764 : pari_sp av = avma;
834 : GEN L, zn, zns;
835 : long i, phi, ind, l;
836 1764 : if (m == 2) retmkvec(mkvecsmall(1));
837 1750 : zn = znstar(stoi(m));
838 1750 : phi = itos(gel(zn,1));
839 1750 : o = ugcd(o, phi); /* do we impose this on input ? */
840 1750 : zns = znstar_small(zn);
841 1750 : L = cgetg(o+1, t_VEC);
842 5782 : for (i=1,ind = phi; ind; ind -= phi/o, i++) /* by *decreasing* exact index */
843 4032 : gel(L,i) = subgrouplist(gel(zn,2), mkvec(utoipos(ind)));
844 1750 : L = shallowconcat1(L); l = lg(L);
845 5502 : for (i = 1; i < l; i++) gel(L,i) = znstar_hnf_elts(zns, gel(L,i));
846 1750 : return gerepilecopy(av, L);
847 : }
848 :
849 : /* A sympol is a symmetric polynomial
850 : *
851 : * Currently sympol are couple of t_VECSMALL [v,w]
852 : * v[1]...v[k], w[1]...w[k] represent the polynomial sum(i=1,k,v[i]*s_w[i])
853 : * where s_i(X_1,...,X_n) = sum(j=1,n,X_j^i) */
854 :
855 : static GEN
856 17732 : Flm_newtonsum(GEN M, ulong e, ulong p)
857 : {
858 17732 : long f = lg(M), g = lg(gel(M,1)), i, j;
859 17732 : GEN NS = cgetg(f, t_VECSMALL);
860 83847 : for(i=1; i<f; i++)
861 : {
862 66115 : ulong s = 0;
863 66115 : GEN Mi = gel(M,i);
864 294974 : for(j = 1; j < g; j++)
865 228859 : s = Fl_add(s, Fl_powu(uel(Mi,j), e, p), p);
866 66115 : uel(NS,i) = s;
867 : }
868 17732 : return NS;
869 : }
870 :
871 : static GEN
872 11200 : Flv_sympol_eval(GEN v, GEN NS, ulong p)
873 : {
874 11200 : pari_sp av = avma;
875 11200 : long i, l = lg(v);
876 11200 : GEN S = Flv_Fl_mul(gel(NS,1), uel(v,1), p);
877 11840 : for (i=2; i<l; i++)
878 640 : if (v[i]) S = Flv_add(S, Flv_Fl_mul(gel(NS,i), uel(v,i), p), p);
879 11200 : return gerepileuptoleaf(av, S);
880 : }
881 :
882 : static GEN
883 11200 : sympol_eval_newtonsum(long e, GEN O, GEN mod)
884 : {
885 11200 : long f = lg(O), g = lg(gel(O,1)), i, j;
886 11200 : GEN PL = cgetg(f, t_COL);
887 55039 : for(i=1; i<f; i++)
888 : {
889 43840 : pari_sp av = avma;
890 43840 : GEN s = gen_0;
891 179813 : for(j=1; j<g; j++) s = addii(s, Fp_powu(gmael(O,i,j), e, mod));
892 43840 : gel(PL,i) = gerepileuptoint(av, remii(s,mod));
893 : }
894 11199 : return PL;
895 : }
896 :
897 : static GEN
898 11137 : sympol_eval(GEN sym, GEN O, GEN mod)
899 : {
900 11137 : pari_sp av = avma;
901 : long i;
902 11137 : GEN v = gel(sym,1), w = gel(sym,2);
903 11137 : GEN S = gen_0;
904 22850 : for (i=1; i<lg(v); i++)
905 11714 : if (v[i]) S = gadd(S, gmulsg(v[i], sympol_eval_newtonsum(w[i], O, mod)));
906 11136 : return gerepileupto(av, S);
907 : }
908 :
909 : /* Let sigma be an automorphism of L (as a polynomial with rational coefs)
910 : * Let 'sym' be a symmetric polynomial defining alpha in L.
911 : * We have alpha = sym(x,sigma(x),,,sigma^(g-1)(x)). Compute alpha mod p */
912 : static GEN
913 5390 : sympol_aut_evalmod(GEN sym, long g, GEN sigma, GEN Tp, GEN p)
914 : {
915 5390 : pari_sp ltop=avma;
916 5390 : long i, j, npows = brent_kung_optpow(degpol(Tp)-1, g-1, 1);
917 5390 : GEN s, f, pows, v = zv_to_ZV(gel(sym,1)), w = zv_to_ZV(gel(sym,2));
918 5390 : sigma = RgX_to_FpX(sigma, p);
919 5390 : pows = FpXQ_powers(sigma,npows,Tp,p);
920 5390 : f = pol_x(varn(sigma));
921 5390 : s = pol_0(varn(sigma));
922 22162 : for(i=1; i<=g;i++)
923 : {
924 16772 : if (i > 1) f = FpX_FpXQV_eval(f,pows,Tp,p);
925 33922 : for(j=1; j<lg(v); j++)
926 17150 : s = FpX_add(s, FpX_Fp_mul(FpXQ_pow(f,gel(w,j),Tp,p),gel(v,j),p),p);
927 : }
928 5390 : return gerepileupto(ltop, s);
929 : }
930 :
931 : /* Let Sp be as computed with sympol_aut_evalmod
932 : * Let Tmod be the factorisation of T mod p.
933 : * Return the factorisation of the minimal polynomial of S mod p w.r.t. Tmod */
934 : static GEN
935 5390 : fixedfieldfactmod(GEN Sp, GEN p, GEN Tmod)
936 : {
937 5390 : long i, l = lg(Tmod);
938 5390 : GEN F = cgetg(l,t_VEC);
939 18172 : for(i=1; i<l; i++)
940 : {
941 12782 : GEN Ti = gel(Tmod,i);
942 12782 : gel(F,i) = FpXQ_minpoly(FpX_rem(Sp,Ti,p), Ti,p);
943 : }
944 5390 : return F;
945 : }
946 :
947 : static GEN
948 11137 : fixedfieldsurmer(ulong l, GEN NS, GEN W)
949 : {
950 11137 : const long step=3;
951 11137 : long i, j, n = lg(W)-1, m = 1L<<((n-1)<<1);
952 11137 : GEN sym = cgetg(n+1,t_VECSMALL);
953 11714 : for (j=1;j<n;j++) sym[j] = step;
954 11137 : sym[n] = 0;
955 11137 : if (DEBUGLEVEL>=4) err_printf("FixedField: Weight: %Ps\n",W);
956 11200 : for (i=0; i<m; i++)
957 : {
958 11200 : pari_sp av = avma;
959 : GEN L;
960 11777 : for (j=1; sym[j]==step; j++) sym[j]=0;
961 11200 : sym[j]++;
962 11200 : if (DEBUGLEVEL>=6) err_printf("FixedField: Sym: %Ps\n",sym);
963 11200 : L = Flv_sympol_eval(sym, NS, l);
964 11200 : if (!vecsmall_is1to1(L)) { set_avma(av); continue; }
965 11137 : return mkvec2(sym,W);
966 : }
967 0 : return NULL;
968 : }
969 :
970 : /*Check whether the line of NS are pair-wise distinct.*/
971 : static long
972 11714 : sympol_is1to1_lg(GEN NS, long n)
973 : {
974 11714 : long i, j, k, l = lgcols(NS);
975 54298 : for (i=1; i<l; i++)
976 172574 : for(j=i+1; j<l; j++)
977 : {
978 133797 : for(k=1; k<n; k++)
979 133220 : if (mael(NS,k,j)!=mael(NS,k,i)) break;
980 129990 : if (k>=n) return 0;
981 : }
982 11137 : return 1;
983 : }
984 :
985 : /* Let O a set of orbits of roots (see fixedfieldorbits) modulo mod,
986 : * l | mod and p two prime numbers. Return a vector [sym,s,P] where:
987 : * sym is a sympol, s is the set of images of sym on O and
988 : * P is the polynomial with roots s. */
989 : static GEN
990 11137 : fixedfieldsympol(GEN O, ulong l)
991 : {
992 11137 : pari_sp ltop=avma;
993 11137 : const long n=(BITS_IN_LONG>>1)-1;
994 11137 : GEN NS = cgetg(n+1,t_MAT), sym = NULL, W = cgetg(n+1,t_VECSMALL);
995 11137 : long i, e=1;
996 11137 : if (DEBUGLEVEL>=4)
997 0 : err_printf("FixedField: Size: %ldx%ld\n",lg(O)-1,lg(gel(O,1))-1);
998 11137 : O = ZM_to_Flm(O,l);
999 22851 : for (i=1; !sym && i<=n; i++)
1000 : {
1001 11714 : GEN L = Flm_newtonsum(O, e++, l);
1002 11714 : if (lg(O)>2)
1003 17634 : while (vecsmall_isconst(L)) L = Flm_newtonsum(O, e++, l);
1004 11714 : W[i] = e-1; gel(NS,i) = L;
1005 11714 : if (sympol_is1to1_lg(NS,i+1))
1006 11137 : sym = fixedfieldsurmer(l,NS,vecsmall_shorten(W,i));
1007 : }
1008 11137 : if (!sym) pari_err_BUG("fixedfieldsympol [p too small]");
1009 11137 : if (DEBUGLEVEL>=2) err_printf("FixedField: Found: %Ps\n",gel(sym,1));
1010 11137 : return gerepilecopy(ltop,sym);
1011 : }
1012 :
1013 : /* Let O a set of orbits as indices and L the corresponding roots.
1014 : * Return the set of orbits as roots. */
1015 : static GEN
1016 11137 : fixedfieldorbits(GEN O, GEN L)
1017 : {
1018 11137 : GEN S = cgetg(lg(O), t_MAT);
1019 : long i;
1020 53662 : for (i = 1; i < lg(O); i++) gel(S,i) = vecpermute(L, gel(O,i));
1021 11137 : return S;
1022 : }
1023 :
1024 : static GEN
1025 1043 : fixedfieldinclusion(GEN O, GEN PL)
1026 : {
1027 1043 : long i, j, f = lg(O)-1, g = lg(gel(O,1))-1;
1028 1043 : GEN S = cgetg(f*g + 1, t_COL);
1029 7042 : for (i = 1; i <= f; i++)
1030 : {
1031 5999 : GEN Oi = gel(O,i);
1032 24836 : for (j = 1; j <= g; j++) gel(S, Oi[j]) = gel(PL, i);
1033 : }
1034 1043 : return S;
1035 : }
1036 :
1037 : /* Polynomial attached to a vector of conjugates. Not stack clean */
1038 : static GEN
1039 51412 : vectopol(GEN v, GEN M, GEN den , GEN mod, GEN mod2, long x)
1040 : {
1041 51412 : long l = lg(v)+1, i;
1042 51412 : GEN z = cgetg(l,t_POL);
1043 51413 : z[1] = evalsigne(1)|evalvarn(x);
1044 394918 : for (i=2; i<l; i++)
1045 343505 : gel(z,i) = gdiv(centermodii(ZMrow_ZC_mul(M,v,i-1), mod, mod2), den);
1046 51413 : return normalizepol_lg(z, l);
1047 : }
1048 :
1049 : /* Polynomial associate to a permutation of the roots. Not stack clean */
1050 : static GEN
1051 49830 : permtopol(GEN p, GEN L, GEN M, GEN den, GEN mod, GEN mod2, long x)
1052 : {
1053 49830 : if (lg(p) != lg(L)) pari_err_TYPE("permtopol [permutation]", p);
1054 49830 : return vectopol(vecpermute(L,p), M, den, mod, mod2, x);
1055 : }
1056 :
1057 : static GEN
1058 8694 : galoisvecpermtopol(GEN gal, GEN vec, GEN mod, GEN mod2)
1059 : {
1060 8694 : long i, l = lg(vec);
1061 8694 : long v = varn(gal_get_pol(gal));
1062 8694 : GEN L = gal_get_roots(gal);
1063 8694 : GEN M = gal_get_invvdm(gal);
1064 8694 : GEN P = cgetg(l, t_MAT);
1065 45395 : for (i=1; i<l; i++)
1066 36701 : gel(P, i) = vecpermute(L,gel(vec,i));
1067 8694 : P = RgM_to_RgXV(FpM_center(FpM_mul(M, P, mod), mod, mod2), v);
1068 8694 : return gdiv(P, gal_get_den(gal));
1069 : }
1070 :
1071 : static void
1072 66715 : notgalois(long p, struct galois_analysis *ga)
1073 : {
1074 66715 : if (DEBUGLEVEL >= 2) err_printf("GaloisAnalysis:non Galois for p=%ld\n", p);
1075 66715 : ga->p = p;
1076 66715 : ga->deg = 0;
1077 66715 : }
1078 :
1079 : /*Gather information about the group*/
1080 : static long
1081 96604 : init_group(long n, long np, GEN Fp, GEN Fe, long *porder)
1082 : {
1083 96604 : const long prim_nonwss_orders[] = { 48,56,60,72,75,80,196,200,216 };
1084 96604 : long i, phi_order = 1, order = 1, group = 0;
1085 : ulong p;
1086 :
1087 : /* non-WSS groups of this order? */
1088 965855 : for (i=0; i < (long)numberof(prim_nonwss_orders); i++)
1089 869272 : if (n % prim_nonwss_orders[i] == 0) { group |= ga_non_wss; break; }
1090 96604 : if (np == 2 && Fp[2] == 3 && Fe[2] == 1 && Fe[1] > 2) group |= ga_ext_2;
1091 :
1092 143358 : for (i = np; i > 0; i--)
1093 : {
1094 102429 : long p = Fp[i];
1095 102429 : if (phi_order % p == 0) { group |= ga_all_normal; break; }
1096 96654 : order *= p; phi_order *= p-1;
1097 96654 : if (Fe[i] > 1) break;
1098 : }
1099 96604 : if (uisprimepower(n, &p) || n == 135) group |= ga_all_nilpotent;
1100 96604 : if (n <= 104) group |= ga_easy; /* no need to use polynomial algo */
1101 96604 : *porder = order; return group;
1102 : }
1103 :
1104 : /*is a "better" than b ? (if so, update karma) */
1105 : static int
1106 162671 : improves(long a, long b, long plift, long p, long n, long *karma)
1107 : {
1108 162671 : if (!plift || a > b) { *karma = ugcd(p-1,n); return 1; }
1109 159031 : if (a == b) {
1110 156322 : long k = ugcd(p-1,n);
1111 156325 : if (k > *karma) { *karma = k; return 1; }
1112 : }
1113 141067 : return 0; /* worse */
1114 : }
1115 :
1116 : /* return 0 if not galois or not wss */
1117 : static int
1118 96605 : galoisanalysis(GEN T, struct galois_analysis *ga, long calcul_l, GEN bad)
1119 : {
1120 96605 : pari_sp ltop = avma, av;
1121 96605 : long group, linf, n, p, i, karma = 0;
1122 : GEN F, Fp, Fe, Fpe, O;
1123 : long np, order, plift, nbmax, nbtest, deg;
1124 : pari_timer ti;
1125 : forprime_t S;
1126 96605 : if (DEBUGLEVEL >= 1) timer_start(&ti);
1127 96605 : n = degpol(T);
1128 96606 : O = zero_zv(n);
1129 96603 : F = factoru_pow(n);
1130 96604 : Fp = gel(F,1); np = lg(Fp)-1;
1131 96604 : Fe = gel(F,2);
1132 96604 : Fpe= gel(F,3);
1133 96604 : group = init_group(n, np, Fp, Fe, &order);
1134 :
1135 : /*Now we study the orders of the Frobenius elements*/
1136 96604 : deg = Fp[np]; /* largest prime | n */
1137 96604 : plift = 0;
1138 96604 : nbtest = 0;
1139 96604 : nbmax = 8+(n>>1);
1140 96604 : u_forprime_init(&S, n*maxss(expu(n)-3, 2), ULONG_MAX);
1141 96605 : av = avma;
1142 533647 : while (!plift || (nbtest < nbmax && (nbtest <=8 || order < (n>>1)))
1143 30163 : || ((n == 24 || n==36) && O[6] == 0 && O[4] == 0)
1144 30163 : || ((group&ga_non_wss) && order == Fp[np]))
1145 : {
1146 503483 : long d, o, norm_o = 1;
1147 : GEN D, Tp;
1148 :
1149 503483 : if ((group&ga_non_wss) && nbtest >= 3*nbmax) break; /* in all cases */
1150 503483 : nbtest++; set_avma(av);
1151 503482 : p = u_forprime_next(&S);
1152 503480 : if (!p) pari_err_OVERFLOW("galoisanalysis [ran out of primes]");
1153 543450 : if (bad && dvdiu(bad, p)) continue;
1154 503480 : Tp = ZX_to_Flx(T,p);
1155 503458 : if (!Flx_is_squarefree(Tp,p)) { if (!--nbtest) nbtest = 1; continue; }
1156 :
1157 463469 : D = Flx_nbfact_by_degree(Tp, &d, p);
1158 463495 : o = n / d; /* d factors, all should have degree o */
1159 463495 : if (D[o] != d) { notgalois(p, ga); return gc_bool(ltop,0); }
1160 :
1161 397089 : if (!O[o]) O[o] = p;
1162 397089 : if (o % deg) goto ga_end; /* NB: deg > 1 */
1163 245681 : if ((group&ga_all_normal) && o < order) goto ga_end;
1164 :
1165 : /*Frob_p has order o > 1, find a power which generates a normal subgroup*/
1166 245527 : if (o * Fp[1] >= n)
1167 229784 : norm_o = o; /*subgroups of smallest index are normal*/
1168 : else
1169 : {
1170 19202 : for (i = np; i > 0; i--)
1171 : {
1172 19203 : if (o % Fpe[i]) break;
1173 3459 : norm_o *= Fpe[i];
1174 : }
1175 : }
1176 : /* Frob_p^(o/norm_o) generates a normal subgroup of order norm_o */
1177 245527 : if (norm_o != 1)
1178 : {
1179 233244 : if (!(group&ga_all_normal) || o > order)
1180 82856 : karma = ugcd(p-1,n);
1181 150388 : else if (!improves(norm_o, deg, plift,p,n, &karma)) goto ga_end;
1182 : /* karma0=0, deg0<=norm_o -> the first improves() returns 1 */
1183 102140 : deg = norm_o; group |= ga_all_normal; /* STORE */
1184 : }
1185 12283 : else if (group&ga_all_normal) goto ga_end;
1186 12283 : else if (!improves(o, order, plift,p,n, &karma)) goto ga_end;
1187 :
1188 104464 : order = o; plift = p; /* STORE */
1189 397094 : ga_end:
1190 397094 : if (DEBUGLEVEL >= 5)
1191 0 : err_printf("GaloisAnalysis:Nbtest=%ld,p=%ld,o=%ld,n_o=%d,best p=%ld,ord=%ld,k=%ld\n", nbtest, p, o, norm_o, plift, order,karma);
1192 : }
1193 : /* To avoid looping on non-WSS group.
1194 : * TODO: check for large groups. Would it be better to disable this check if
1195 : * we are in a good case (ga_all_normal && !(ga_ext_2) (e.g. 60)) ?*/
1196 30164 : ga->p = plift;
1197 30164 : if (!plift || ((group&ga_non_wss) && order == Fp[np]))
1198 : {
1199 0 : if (DEBUGLEVEL)
1200 0 : pari_warn(warner,"Galois group probably not weakly super solvable");
1201 0 : return 0;
1202 : }
1203 30169 : linf = 2*n*usqrt(n);
1204 30169 : if (calcul_l && O[1] <= linf)
1205 : {
1206 : pari_sp av2;
1207 : forprime_t S2;
1208 : ulong p;
1209 6895 : u_forprime_init(&S2, linf+1,ULONG_MAX);
1210 6895 : av2 = avma;
1211 97622 : while ((p = u_forprime_next(&S2)))
1212 : { /*find a totally split prime l > linf*/
1213 97622 : GEN Tp = ZX_to_Flx(T, p);
1214 97622 : long nb = Flx_nbroots(Tp, p);
1215 97622 : if (nb == n) { O[1] = p; break; }
1216 91007 : if (nb && Flx_is_squarefree(Tp,p)) { notgalois(p,ga); return gc_bool(ltop,0); }
1217 90727 : set_avma(av2);
1218 : }
1219 6615 : if (!p) pari_err_OVERFLOW("galoisanalysis [ran out of primes]");
1220 : }
1221 29889 : ga->group = group;
1222 29889 : ga->deg = deg;
1223 29889 : ga->ord = order;
1224 29889 : ga->l = O[1];
1225 29889 : ga->p4 = n >= 4 ? O[4] : 0;
1226 29889 : if (DEBUGLEVEL >= 4)
1227 0 : err_printf("GaloisAnalysis:p=%ld l=%ld group=%ld deg=%ld ord=%ld\n",
1228 0 : plift, O[1], group, deg, order);
1229 29889 : if (DEBUGLEVEL >= 1) timer_printf(&ti, "galoisanalysis()");
1230 29889 : return gc_bool(ltop,1);
1231 : }
1232 :
1233 : static GEN
1234 98 : a4galoisgen(struct galois_test *td)
1235 : {
1236 98 : const long n = 12;
1237 98 : pari_sp ltop = avma, av, av2;
1238 98 : long i, j, k, N, hop = 0;
1239 : GEN MT, O,O1,O2,O3, ar, mt, t, u, res, orb, pft, pfu, pfv;
1240 :
1241 98 : res = cgetg(3, t_VEC);
1242 98 : pft = cgetg(n+1, t_VECSMALL);
1243 98 : pfu = cgetg(n+1, t_VECSMALL);
1244 98 : pfv = cgetg(n+1, t_VECSMALL);
1245 98 : gel(res,1) = mkvec3(pft,pfu,pfv);
1246 98 : gel(res,2) = mkvecsmall3(2,2,3);
1247 98 : av = avma;
1248 98 : ar = cgetg(5, t_VECSMALL);
1249 98 : mt = gel(td->PV, td->order[n]);
1250 98 : t = identity_perm(n) + 1; /* Sorry for this hack */
1251 98 : u = cgetg(n+1, t_VECSMALL) + 1; /* too lazy to correct */
1252 98 : MT = cgetg(n+1, t_MAT);
1253 1274 : for (j = 1; j <= n; j++) gel(MT,j) = cgetg(n+1, t_VECSMALL);
1254 1274 : for (j = 1; j <= n; j++)
1255 7644 : for (i = 1; i < j; i++)
1256 6468 : ucoeff(MT,i,j) = ucoeff(MT,j,i) = ucoeff(mt,i,j)+ucoeff(mt,j,i);
1257 : /* MT(i,i) unused */
1258 :
1259 98 : av2 = avma;
1260 : /* N = itos(gdiv(mpfact(n), mpfact(n >> 1))) >> (n >> 1); */
1261 : /* n = 2k = 12; N = (2k)! / (k! * 2^k) = 10395 */
1262 98 : N = 10395;
1263 98 : if (DEBUGLEVEL>=4) err_printf("A4GaloisConj: will test %ld permutations\n", N);
1264 98 : uel(ar,4) = umael(MT,11,12);
1265 98 : uel(ar,3) = uel(ar,4) + umael(MT,9,10);
1266 98 : uel(ar,2) = uel(ar,3) + umael(MT,7,8);
1267 98 : uel(ar,1) = uel(ar,2) + umael(MT,5,6);
1268 232152 : for (i = 0; i < N; i++)
1269 : {
1270 : long g;
1271 232152 : if (i)
1272 : {
1273 232054 : long a, x = i, y = 1;
1274 327152 : do { y += 2; a = x%y; x = x/y; } while (!a);
1275 232054 : switch (y)
1276 : {
1277 154744 : case 3:
1278 154744 : lswap(t[2], t[2-a]);
1279 154744 : break;
1280 61875 : case 5:
1281 61875 : x = t[0]; t[0] = t[2]; t[2] = t[1]; t[1] = x;
1282 61875 : lswap(t[4], t[4-a]);
1283 61875 : uel(ar,1) = uel(ar,2) + umael(MT,t[4],t[5]);
1284 61875 : break;
1285 13274 : case 7:
1286 13274 : x = t[0]; t[0] = t[4]; t[4] = t[3]; t[3] = t[1]; t[1] = t[2]; t[2] = x;
1287 13274 : lswap(t[6], t[6-a]);
1288 13274 : uel(ar,2) = uel(ar,3) + umael(MT,t[6],t[7]);
1289 13274 : uel(ar,1) = uel(ar,2) + umael(MT,t[4],t[5]);
1290 13274 : break;
1291 1969 : case 9:
1292 1969 : x = t[0]; t[0] = t[6]; t[6] = t[5]; t[5] = t[3]; t[3] = x;
1293 1969 : lswap(t[1], t[4]);
1294 1969 : lswap(t[8], t[8-a]);
1295 1969 : uel(ar,3) = uel(ar,4) + umael(MT,t[8],t[9]);
1296 1969 : uel(ar,2) = uel(ar,3) + umael(MT,t[6],t[7]);
1297 1969 : uel(ar,1) = uel(ar,2) + umael(MT,t[4],t[5]);
1298 1969 : break;
1299 192 : case 11:
1300 192 : x = t[0]; t[0] = t[8]; t[8] = t[7]; t[7] = t[5]; t[5] = t[1];
1301 192 : t[1] = t[6]; t[6] = t[3]; t[3] = t[2]; t[2] = t[4]; t[4] = x;
1302 192 : lswap(t[10], t[10-a]);
1303 192 : uel(ar,4) = umael(MT,t[10],t[11]);
1304 192 : uel(ar,3) = uel(ar,4) + umael(MT,t[8],t[9]);
1305 192 : uel(ar,2) = uel(ar,3) + umael(MT,t[6],t[7]);
1306 192 : uel(ar,1) = uel(ar,2) + umael(MT,t[4],t[5]);
1307 : }
1308 98 : }
1309 232152 : g = uel(ar,1)+umael(MT,t[0],t[1])+umael(MT,t[2],t[3]);
1310 232152 : if (headlongisint(g,n))
1311 : {
1312 686 : for (k = 0; k < n; k += 2)
1313 : {
1314 588 : pft[t[k]] = t[k+1];
1315 588 : pft[t[k+1]] = t[k];
1316 : }
1317 98 : if (galois_test_perm(td, pft)) break;
1318 0 : hop++;
1319 : }
1320 232054 : set_avma(av2);
1321 : }
1322 98 : if (DEBUGLEVEL >= 1 && hop)
1323 0 : err_printf("A4GaloisConj: %ld hop over %ld iterations\n", hop, N);
1324 98 : if (i == N) return gc_NULL(ltop);
1325 : /* N = itos(gdiv(mpfact(n >> 1), mpfact(n >> 2))) >> 1; */
1326 98 : N = 60;
1327 98 : if (DEBUGLEVEL >= 4) err_printf("A4GaloisConj: sigma=%Ps \n", pft);
1328 392 : for (k = 0; k < n; k += 4)
1329 : {
1330 294 : u[k+3] = t[k+3];
1331 294 : u[k+2] = t[k+1];
1332 294 : u[k+1] = t[k+2];
1333 294 : u[k] = t[k];
1334 : }
1335 4658 : for (i = 0; i < N; i++)
1336 : {
1337 4658 : ulong g = 0;
1338 4658 : if (i)
1339 : {
1340 4560 : long a, x = i, y = -2;
1341 7178 : do { y += 4; a = x%y; x = x/y; } while (!a);
1342 4560 : lswap(u[0],u[2]);
1343 4560 : switch (y)
1344 : {
1345 2280 : case 2:
1346 2280 : break;
1347 1942 : case 6:
1348 1942 : lswap(u[4],u[6]);
1349 1942 : if (!(a & 1))
1350 : {
1351 795 : a = 4 - (a>>1);
1352 795 : lswap(u[6], u[a]);
1353 795 : lswap(u[4], u[a-2]);
1354 : }
1355 1942 : break;
1356 338 : case 10:
1357 338 : x = u[6];
1358 338 : u[6] = u[3];
1359 338 : u[3] = u[2];
1360 338 : u[2] = u[4];
1361 338 : u[4] = u[1];
1362 338 : u[1] = u[0];
1363 338 : u[0] = x;
1364 338 : if (a >= 3) a += 2;
1365 338 : a = 8 - a;
1366 338 : lswap(u[10],u[a]);
1367 338 : lswap(u[8], u[a-2]);
1368 338 : break;
1369 : }
1370 98 : }
1371 32606 : for (k = 0; k < n; k += 2) g += mael(MT,u[k],u[k+1]);
1372 4658 : if (headlongisint(g,n))
1373 : {
1374 686 : for (k = 0; k < n; k += 2)
1375 : {
1376 588 : pfu[u[k]] = u[k+1];
1377 588 : pfu[u[k+1]] = u[k];
1378 : }
1379 98 : if (galois_test_perm(td, pfu)) break;
1380 0 : hop++;
1381 : }
1382 4560 : set_avma(av2);
1383 : }
1384 98 : if (i == N) return gc_NULL(ltop);
1385 98 : if (DEBUGLEVEL >= 1 && hop)
1386 0 : err_printf("A4GaloisConj: %ld hop over %ld iterations\n", hop, N);
1387 98 : if (DEBUGLEVEL >= 4) err_printf("A4GaloisConj: tau=%Ps \n", pfu);
1388 98 : set_avma(av2);
1389 98 : orb = mkvec2(pft,pfu);
1390 98 : O = vecperm_orbits(orb, 12);
1391 98 : if (DEBUGLEVEL >= 4) {
1392 0 : err_printf("A4GaloisConj: orb=%Ps\n", orb);
1393 0 : err_printf("A4GaloisConj: O=%Ps \n", O);
1394 : }
1395 98 : av2 = avma;
1396 98 : O1 = gel(O,1); O2 = gel(O,2); O3 = gel(O,3);
1397 140 : for (j = 0; j < 2; j++)
1398 : {
1399 140 : pfv[O1[1]] = O2[1];
1400 140 : pfv[O1[2]] = O2[3+j];
1401 140 : pfv[O1[3]] = O2[4 - (j << 1)];
1402 140 : pfv[O1[4]] = O2[2+j];
1403 469 : for (i = 0; i < 4; i++)
1404 : {
1405 427 : ulong g = 0;
1406 427 : switch (i)
1407 : {
1408 140 : case 0: break;
1409 119 : case 1: lswap(O3[1], O3[2]); lswap(O3[3], O3[4]); break;
1410 105 : case 2: lswap(O3[1], O3[4]); lswap(O3[2], O3[3]); break;
1411 63 : case 3: lswap(O3[1], O3[2]); lswap(O3[3], O3[4]); break;
1412 : }
1413 427 : pfv[O2[1]] = O3[1];
1414 427 : pfv[O2[3+j]] = O3[4-j];
1415 427 : pfv[O2[4 - (j<<1)]] = O3[2 + (j<<1)];
1416 427 : pfv[O2[2+j]] = O3[3-j];
1417 427 : pfv[O3[1]] = O1[1];
1418 427 : pfv[O3[4-j]] = O1[2];
1419 427 : pfv[O3[2 + (j<<1)]] = O1[3];
1420 427 : pfv[O3[3-j]] = O1[4];
1421 5551 : for (k = 1; k <= n; k++) g += mael(mt,k,pfv[k]);
1422 427 : if (headlongisint(g,n) && galois_test_perm(td, pfv))
1423 : {
1424 98 : set_avma(av);
1425 98 : if (DEBUGLEVEL >= 1)
1426 0 : err_printf("A4GaloisConj: %ld hop over %d iterations max\n",
1427 : hop, 10395 + 68);
1428 98 : return res;
1429 : }
1430 329 : hop++; set_avma(av2);
1431 : }
1432 : }
1433 0 : return gc_NULL(ltop);
1434 : }
1435 :
1436 : /* S4 */
1437 : static GEN
1438 1414 : s4makelift(GEN u, struct galois_lift *gl)
1439 1414 : { return FpXQ_powers(u, degpol(gl->T)-1, gl->TQ, gl->Q); }
1440 :
1441 : static long
1442 239230 : s4test(GEN u, GEN liftpow, struct galois_lift *gl, GEN phi)
1443 : {
1444 239230 : pari_sp av = avma;
1445 : GEN res, Q, Q2;
1446 239230 : long bl, i, d = lg(u)-2;
1447 : pari_timer ti;
1448 239230 : if (DEBUGLEVEL >= 6) timer_start(&ti);
1449 239230 : if (!d) return 0;
1450 239230 : Q = gl->Q; Q2 = shifti(Q,-1);
1451 239230 : res = gel(u,2);
1452 8526848 : for (i = 2; i <= d; i++)
1453 8287618 : if (lg(gel(liftpow,i))>2)
1454 8287618 : res = addii(res, mulii(gmael(liftpow,i,2), gel(u,i+1)));
1455 239230 : res = remii(res,Q);
1456 239230 : if (gl->den != gen_1) res = mulii(res, gl->den);
1457 239230 : res = centermodii(res, Q,Q2);
1458 239230 : if (abscmpii(res, gl->gb->bornesol) > 0) return gc_long(av,0);
1459 3630 : res = scalar_ZX_shallow(gel(u,2),varn(u));
1460 120958 : for (i = 2; i <= d ; i++)
1461 117328 : if (lg(gel(liftpow,i))>2)
1462 117328 : res = ZX_add(res, ZX_Z_mul(gel(liftpow,i), gel(u,i+1)));
1463 3630 : res = FpX_red(res, Q);
1464 3630 : if (gl->den != gen_1) res = FpX_Fp_mul(res, gl->den, Q);
1465 3630 : res = FpX_center_i(res, Q, shifti(Q,-1));
1466 3630 : bl = poltopermtest(res, gl, phi);
1467 3630 : if (DEBUGLEVEL >= 6) timer_printf(&ti, "s4test()");
1468 3630 : return gc_long(av,bl);
1469 : }
1470 :
1471 : static GEN
1472 469 : s4releveauto(GEN M, GEN B, GEN T, GEN p,long a1,long a2,long a3,long a4,long a5,long a6)
1473 : {
1474 469 : GEN F = ZX_mul(gel(M,a1),gel(B,a2));
1475 469 : F = ZX_add(F, ZX_mul(gel(M,a2),gel(B,a1)));
1476 469 : F = ZX_add(F, ZX_mul(gel(M,a3),gel(B,a4)));
1477 469 : F = ZX_add(F, ZX_mul(gel(M,a4),gel(B,a3)));
1478 469 : F = ZX_add(F, ZX_mul(gel(M,a5),gel(B,a6)));
1479 469 : F = ZX_add(F, ZX_mul(gel(M,a6),gel(B,a5)));
1480 469 : return FpXQ_red(F, T, p);
1481 : }
1482 :
1483 : static GEN
1484 319203 : lincomb(GEN B, long a, long b, long j)
1485 : {
1486 319203 : long k = (-j) & 3;
1487 319203 : return ZX_add(gmael(B,a,j+1), gmael(B,b,k+1));
1488 : }
1489 :
1490 : static GEN
1491 91 : FpXV_ffisom(GEN V, GEN p)
1492 : {
1493 91 : pari_sp av = avma;
1494 91 : long i, j, l = lg(V);
1495 91 : GEN S = cgetg(l, t_VEC), Si = cgetg(l, t_VEC), M;
1496 679 : for (i = 1; i < l; i++)
1497 : {
1498 588 : gel(S,i) = FpX_ffisom(gel(V,1), gel(V,i), p);
1499 588 : gel(Si,i) = FpXQ_ffisom_inv(gel(S,i), gel(V,i), p);
1500 : }
1501 91 : M = cgetg(l, t_MAT);
1502 679 : for (j = 1; j < l; j++)
1503 588 : gel(M,j) = FpXC_FpXQ_eval(Si, gel(S,j), gel(V,j), p);
1504 91 : return gerepileupto(av, M);
1505 : }
1506 :
1507 : static GEN
1508 91 : mkliftpow(GEN x, GEN T, GEN p, struct galois_lift *gl)
1509 679 : { pari_APPLY_same(automorphismlift(FpXV_chinese(gel(x,i), T, p, NULL), gl)) }
1510 :
1511 : #define rot3(x,y,z) {long _t=x; x=y; y=z; z=_t;}
1512 : #define rot4(x,y,z,t) {long _u=x; x=y; y=z; z=t; t=_u;}
1513 :
1514 : /* FIXME: could use the intheadlong technique */
1515 : static GEN
1516 77 : s4galoisgen(struct galois_lift *gl)
1517 : {
1518 77 : const long n = 24;
1519 : struct galois_testlift gt;
1520 77 : pari_sp av, ltop2, ltop = avma;
1521 : long i, j;
1522 77 : GEN sigma, tau, phi, res, r1,r2,r3,r4, pj, p = gl->p, Q = gl->Q, TQ = gl->TQ;
1523 : GEN sg, Tp, Tmod, misom, B, Bcoeff, liftpow, liftp, aut;
1524 :
1525 77 : res = cgetg(3, t_VEC);
1526 77 : r1 = cgetg(n+1, t_VECSMALL);
1527 77 : r2 = cgetg(n+1, t_VECSMALL);
1528 77 : r3 = cgetg(n+1, t_VECSMALL);
1529 77 : r4 = cgetg(n+1, t_VECSMALL);
1530 77 : gel(res,1)= mkvec4(r1,r2,r3,r4);
1531 77 : gel(res,2) = mkvecsmall4(2,2,3,2);
1532 77 : ltop2 = avma;
1533 77 : sg = identity_perm(6);
1534 77 : pj = zero_zv(6);
1535 77 : sigma = cgetg(n+1, t_VECSMALL);
1536 77 : tau = cgetg(n+1, t_VECSMALL);
1537 77 : phi = cgetg(n+1, t_VECSMALL);
1538 77 : Tp = FpX_red(gl->T,p);
1539 77 : Tmod = gel(FpX_factor(Tp,p), 1);
1540 77 : misom = FpXV_ffisom(Tmod, p);
1541 77 : aut = galoisdolift(gl);
1542 77 : inittestlift(aut, Tmod, gl, >);
1543 77 : B = FqC_FqV_mul(gt.pauto, gt.bezoutcoeff, gl->TQ, Q);
1544 77 : Bcoeff = gt.bezoutcoeff;
1545 77 : liftp = mkliftpow(shallowtrans(misom), Tmod, p, gl);
1546 77 : av = avma;
1547 119 : for (i = 0; i < 3; i++, set_avma(av))
1548 : {
1549 : pari_sp av1, av2, av3;
1550 : GEN u, u1, u2, u3;
1551 : long j1, j2, j3;
1552 119 : if (i)
1553 : {
1554 42 : if (i == 1) { lswap(sg[2],sg[3]); }
1555 0 : else { lswap(sg[1],sg[3]); }
1556 : }
1557 119 : u = s4releveauto(liftp,Bcoeff,TQ,Q,sg[1],sg[2],sg[3],sg[4],sg[5],sg[6]);
1558 119 : liftpow = s4makelift(u, gl);
1559 119 : av1 = avma;
1560 354 : for (j1 = 0; j1 < 4; j1++, set_avma(av1))
1561 : {
1562 312 : u1 = lincomb(B,sg[5],sg[6],j1);
1563 312 : av2 = avma;
1564 1377 : for (j2 = 0; j2 < 4; j2++, set_avma(av2))
1565 : {
1566 1142 : u2 = lincomb(B,sg[3],sg[4],j2);
1567 1142 : u2 = FpX_add(u1, u2, Q); av3 = avma;
1568 5508 : for (j3 = 0; j3 < 4; j3++, set_avma(av3))
1569 : {
1570 4443 : u3 = lincomb(B,sg[1],sg[2],j3);
1571 4443 : u3 = FpX_add(u2, u3, Q);
1572 4443 : if (DEBUGLEVEL >= 4)
1573 0 : err_printf("S4GaloisConj: Testing %d/3:%d/4:%d/4:%d/4:%Ps\n",
1574 : i,j1,j2,j3, sg);
1575 4443 : if (s4test(u3, liftpow, gl, sigma))
1576 : {
1577 77 : pj[1] = j3;
1578 77 : pj[2] = j2;
1579 77 : pj[3] = j1; goto suites4;
1580 : }
1581 : }
1582 : }
1583 : }
1584 : }
1585 0 : return gc_NULL(ltop);
1586 77 : suites4:
1587 77 : if (DEBUGLEVEL >= 4) err_printf("S4GaloisConj: sigma=%Ps\n", sigma);
1588 77 : if (DEBUGLEVEL >= 4) err_printf("S4GaloisConj: pj=%Ps\n", pj);
1589 77 : set_avma(av);
1590 154 : for (j = 1; j <= 3; j++)
1591 : {
1592 : pari_sp av2, av3;
1593 : GEN u;
1594 : long w, l;
1595 154 : rot3(sg[1], sg[3], sg[5])
1596 154 : rot3(sg[2], sg[4], sg[6])
1597 154 : rot3(pj[1], pj[2], pj[3])
1598 350 : for (l = 0; l < 2; l++, set_avma(av))
1599 : {
1600 273 : u = s4releveauto(liftp,Bcoeff,TQ,Q,sg[1],sg[3],sg[2],sg[4],sg[5],sg[6]);
1601 273 : liftpow = s4makelift(u, gl);
1602 273 : av2 = avma;
1603 742 : for (w = 0; w < 4; w += 2, set_avma(av2))
1604 : {
1605 : GEN uu;
1606 546 : pj[6] = (w + pj[3]) & 3;
1607 546 : uu = lincomb(B, sg[5], sg[6], pj[6]);
1608 546 : uu = FpX_red(uu, Q);
1609 546 : av3 = avma;
1610 2541 : for (i = 0; i < 4; i++, set_avma(av3))
1611 : {
1612 : GEN u;
1613 2072 : pj[4] = i;
1614 2072 : pj[5] = (i + pj[2] - pj[1]) & 3;
1615 2072 : if (DEBUGLEVEL >= 4)
1616 0 : err_printf("S4GaloisConj: Testing %d/3:%d/2:%d/2:%d/4:%Ps:%Ps\n",
1617 : j-1, w >> 1, l, i, sg, pj);
1618 2072 : u = ZX_add(lincomb(B,sg[1],sg[3],pj[4]),
1619 2072 : lincomb(B,sg[2],sg[4],pj[5]));
1620 2072 : u = FpX_add(uu,u,Q);
1621 2072 : if (s4test(u, liftpow, gl, tau)) goto suites4_2;
1622 : }
1623 : }
1624 196 : lswap(sg[3], sg[4]);
1625 196 : pj[2] = (-pj[2]) & 3;
1626 : }
1627 : }
1628 0 : return gc_NULL(ltop);
1629 77 : suites4_2:
1630 77 : set_avma(av);
1631 : {
1632 77 : long abc = (pj[1] + pj[2] + pj[3]) & 3;
1633 77 : long abcdef = ((abc + pj[4] + pj[5] - pj[6]) & 3) >> 1;
1634 : GEN u;
1635 : pari_sp av2;
1636 77 : u = s4releveauto(liftp,Bcoeff,TQ,Q,sg[1],sg[4],sg[2],sg[5],sg[3],sg[6]);
1637 77 : liftpow = s4makelift(u, gl);
1638 77 : av2 = avma;
1639 322 : for (j = 0; j < 8; j++)
1640 : {
1641 : long h, g, i;
1642 322 : h = j & 3;
1643 322 : g = (abcdef + ((j & 4) >> 1)) & 3;
1644 322 : i = (h + abc - g) & 3;
1645 322 : u = ZX_add(lincomb(B,sg[1],sg[4], g), lincomb(B,sg[2],sg[5], h));
1646 322 : u = FpX_add(u, lincomb(B,sg[3],sg[6], i),Q);
1647 322 : if (DEBUGLEVEL >= 4)
1648 0 : err_printf("S4GaloisConj: Testing %d/8 %d:%d:%d\n", j,g,h,i);
1649 322 : if (s4test(u, liftpow, gl, phi)) break;
1650 245 : set_avma(av2);
1651 : }
1652 : }
1653 77 : if (j == 8) return gc_NULL(ltop);
1654 1925 : for (i = 1; i <= n; i++)
1655 : {
1656 1848 : r1[i] = sigma[tau[i]];
1657 1848 : r2[i] = phi[sigma[tau[phi[i]]]];
1658 1848 : r3[i] = phi[sigma[i]];
1659 1848 : r4[i] = sigma[i];
1660 : }
1661 77 : set_avma(ltop2); return res;
1662 : }
1663 :
1664 : static GEN
1665 910 : f36releveauto2(GEN Bl, GEN T, GEN p,GEN a)
1666 : {
1667 910 : GEN F = gmael(Bl,a[1],a[1]);
1668 910 : F = ZX_add(F,gmael(Bl,a[2],a[3]));
1669 910 : F = ZX_add(F,gmael(Bl,a[3],a[2]));
1670 910 : F = ZX_add(F,gmael(Bl,a[4],a[5]));
1671 910 : F = ZX_add(F,gmael(Bl,a[5],a[4]));
1672 910 : F = ZX_add(F,gmael(Bl,a[6],a[7]));
1673 910 : F = ZX_add(F,gmael(Bl,a[7],a[6]));
1674 910 : F = ZX_add(F,gmael(Bl,a[8],a[9]));
1675 910 : F = ZX_add(F,gmael(Bl,a[9],a[8]));
1676 910 : return FpXQ_red(F, T, p);
1677 : }
1678 :
1679 : static GEN
1680 35 : f36releveauto4(GEN Bl, GEN T, GEN p,GEN a)
1681 : {
1682 35 : GEN F = gmael(Bl,a[1],a[1]);
1683 35 : F = ZX_add(F,gmael(Bl,a[2],a[3]));
1684 35 : F = ZX_add(F,gmael(Bl,a[3],a[4]));
1685 35 : F = ZX_add(F,gmael(Bl,a[4],a[5]));
1686 35 : F = ZX_add(F,gmael(Bl,a[5],a[2]));
1687 35 : F = ZX_add(F,gmael(Bl,a[6],a[7]));
1688 35 : F = ZX_add(F,gmael(Bl,a[7],a[8]));
1689 35 : F = ZX_add(F,gmael(Bl,a[8],a[9]));
1690 35 : F = ZX_add(F,gmael(Bl,a[9],a[6]));
1691 35 : return FpXQ_red(F, T, p);
1692 : }
1693 :
1694 : static GEN
1695 14 : f36galoisgen(struct galois_lift *gl)
1696 : {
1697 14 : const long n = 36;
1698 : struct galois_testlift gt;
1699 14 : pari_sp av, ltop2, ltop = avma;
1700 : long i;
1701 14 : GEN sigma, tau, rho, res, r1,r2,r3, pj, pk, p = gl->p, Q = gl->Q, TQ = gl->TQ;
1702 : GEN sg, s4, sp, Tp, Tmod, misom, Bcoeff, liftpow, aut, liftp, B, Bl, tam;
1703 14 : res = cgetg(3, t_VEC);
1704 14 : r1 = cgetg(n+1, t_VECSMALL);
1705 14 : r2 = cgetg(n+1, t_VECSMALL);
1706 14 : r3 = cgetg(n+1, t_VECSMALL);
1707 14 : gel(res,1)= mkvec3(r1,r2,r3);
1708 14 : gel(res,2) = mkvecsmall3(3,3,4);
1709 14 : ltop2 = avma;
1710 14 : sg = identity_perm(9);
1711 14 : s4 = identity_perm(9);
1712 14 : sp = identity_perm(9);
1713 14 : pj = zero_zv(4);
1714 14 : pk = zero_zv(2);
1715 14 : sigma = cgetg(n+1, t_VECSMALL);
1716 14 : tau = r3;
1717 14 : rho = cgetg(n+1, t_VECSMALL);
1718 14 : Tp = FpX_red(gl->T,p);
1719 14 : Tmod = gel(FpX_factor(Tp,p), 1);
1720 14 : misom = FpXV_ffisom(Tmod, p);
1721 14 : aut = galoisdolift(gl);
1722 14 : inittestlift(aut, Tmod, gl, >);
1723 14 : Bcoeff = gt.bezoutcoeff;
1724 14 : B = FqC_FqV_mul(gt.pauto, Bcoeff, gl->TQ, gl->Q);
1725 14 : liftp = mkliftpow(shallowtrans(misom), Tmod, p, gl);
1726 14 : Bl = FqC_FqV_mul(liftp,Bcoeff, gl->TQ, gl->Q);
1727 14 : av = avma;
1728 910 : for (i = 0; i < 105; i++, set_avma(av))
1729 : {
1730 : pari_sp av0, av1, av2, av3;
1731 : GEN u0, u1, u2, u3;
1732 : long j0, j1, j2, j3, s;
1733 910 : if (i)
1734 : {
1735 896 : rot3(sg[7],sg[8],sg[9])
1736 896 : if (i%3==0)
1737 : {
1738 294 : s=sg[5]; sg[5]=sg[6]; sg[6]=sg[7]; sg[7]=sg[8]; sg[8]=sg[9]; sg[9]=s;
1739 294 : if (i%15==0)
1740 : {
1741 49 : s=sg[3]; sg[3]=sg[4]; sg[4]=sg[5];
1742 49 : sg[5]=sg[6]; sg[6]=sg[7]; sg[7]=sg[8]; sg[8]=sg[9]; sg[9]=s;
1743 : }
1744 : }
1745 : }
1746 910 : liftpow = s4makelift(f36releveauto2(Bl, TQ, Q, sg), gl);
1747 910 : av0 = avma;
1748 4522 : for (j0 = 0; j0 < 4; j0++, set_avma(av0))
1749 : {
1750 3626 : u0 = lincomb(B,sg[8],sg[9],j0);
1751 3626 : u0 = FpX_add(u0, gmael(B,sg[1],3), Q); av1 = avma;
1752 18095 : for (j1 = 0; j1 < 4; j1++, set_avma(av1))
1753 : {
1754 14483 : u1 = lincomb(B,sg[6],sg[7],j1);
1755 14483 : u1 = FpX_add(u0, u1, Q); av2 = avma;
1756 72380 : for (j2 = 0; j2 < 4; j2++, set_avma(av2))
1757 : {
1758 57911 : u2 = lincomb(B,sg[4],sg[5],j2);
1759 57911 : u2 = FpX_add(u1, u2, Q); av3 = avma;
1760 289527 : for (j3 = 0; j3 < 4; j3++, set_avma(av3))
1761 : {
1762 231630 : u3 = lincomb(B,sg[2],sg[3],j3);
1763 231630 : u3 = FpX_add(u2, u3, Q);
1764 231630 : if (s4test(u3, liftpow, gl, sigma))
1765 : {
1766 14 : pj[1] = j3;
1767 14 : pj[2] = j2;
1768 14 : pj[3] = j1;
1769 14 : pj[4] = j0; goto suitef36;
1770 : }
1771 : }
1772 : }
1773 : }
1774 : }
1775 : }
1776 0 : return gc_NULL(ltop);
1777 14 : suitef36:
1778 14 : s4[1]=sg[1]; s4[2]=sg[2]; s4[4]=sg[3];
1779 14 : s4[3]=sg[4]; s4[5]=sg[5]; s4[6]=sg[6];
1780 14 : s4[8]=sg[7]; s4[7]=sg[8]; s4[9]=sg[9];
1781 14 : for (i = 0; i < 12; i++, set_avma(av))
1782 : {
1783 : pari_sp av0, av1;
1784 : GEN u0, u1;
1785 : long j0, j1;
1786 14 : if (i)
1787 : {
1788 0 : lswap(s4[3],s4[5]); pj[2] = (-pj[2])&3;
1789 0 : if (odd(i)) { lswap(s4[7],s4[9]); pj[4]=(-pj[4])&3; }
1790 0 : if (i%4==0)
1791 : {
1792 0 : rot3(s4[3],s4[6],s4[7]);
1793 0 : rot3(s4[5],s4[8],s4[9]);
1794 0 : rot3(pj[2],pj[3],pj[4]);
1795 : }
1796 : }
1797 14 : liftpow = s4makelift(f36releveauto4(Bl, TQ, Q, s4), gl);
1798 14 : av0 = avma;
1799 21 : for (j0 = 0; j0 < 4; j0++, set_avma(av0))
1800 : {
1801 21 : u0 = FpX_add(gmael(B,s4[1],2), gmael(B,s4[2],1+j0),Q);
1802 21 : u0 = FpX_add(u0, gmael(B,s4[3],1+smodss(pj[2]-j0,4)),Q);
1803 21 : u0 = FpX_add(u0, gmael(B,s4[4],1+smodss(j0-pj[1]-pj[2],4)),Q);
1804 21 : u0 = FpX_add(u0, gmael(B,s4[5],1+smodss(pj[1]-j0,4)),Q);
1805 21 : av1 = avma;
1806 84 : for (j1 = 0; j1 < 4; j1++, set_avma(av1))
1807 : {
1808 77 : u1 = FpX_add(u0, gmael(B,s4[6],1+j1),Q);
1809 77 : u1 = FpX_add(u1, gmael(B,s4[7],1+smodss(pj[4]-j1,4)),Q);
1810 77 : u1 = FpX_add(u1, gmael(B,s4[8],1+smodss(j1-pj[3]-pj[4],4)),Q);
1811 77 : u1 = FpX_add(u1, gmael(B,s4[9],1+smodss(pj[3]-j1,4)),Q);
1812 77 : if (s4test(u1, liftpow, gl, tau))
1813 : {
1814 14 : pk[1] = j0;
1815 14 : pk[2] = j1; goto suitef36_2;
1816 : }
1817 : }
1818 : }
1819 : }
1820 0 : return gc_NULL(ltop);
1821 14 : suitef36_2:
1822 14 : sp[1]=s4[9]; sp[2]=s4[1]; sp[3]=s4[2];
1823 14 : sp[4]=s4[7]; sp[5]=s4[3]; sp[6]=s4[8];
1824 14 : sp[8]=s4[4]; sp[7]=s4[5]; sp[9]=s4[6];
1825 21 : for (i = 0; i < 4; i++, set_avma(av))
1826 : {
1827 21 : const int w[4][6]={{0,0,1,3,0,2},{1,0,2,1,1,2},{3,3,2,0,3,1},{0,1,3,0,0,3}};
1828 : pari_sp av0, av1, av2;
1829 : GEN u0, u1, u2;
1830 : long j0, j1,j2,j3,j4,j5;
1831 21 : if (i)
1832 : {
1833 7 : rot4(sp[3],sp[5],sp[8],sp[7])
1834 7 : pk[1]=(-pk[1])&3;
1835 : }
1836 21 : liftpow = s4makelift(f36releveauto4(Bl,TQ,Q,sp), gl);
1837 21 : av0 = avma;
1838 56 : for (j0 = 0; j0 < 4; j0++, set_avma(av0))
1839 : {
1840 49 : u0 = FpX_add(gmael(B,sp[1],2), gmael(B,sp[2],1+j0),Q);
1841 49 : av1 = avma;
1842 210 : for (j1 = 0; j1 < 4; j1++, set_avma(av1))
1843 : {
1844 175 : u1 = FpX_add(u0, gmael(B,sp[3],1+j1),Q);
1845 175 : j3 = (-pk[1]-pj[3]+j0+j1-w[i][0]*pj[1]-w[i][3]*pj[2])&3;
1846 175 : u1 = FpX_add(u1, gmael(B,sp[6],1+j3),Q);
1847 175 : j5 = (-pk[1]+2*j0+2*j1-w[i][2]*pj[1]-w[i][5]*pj[2])&3;
1848 175 : u1 = FpX_add(u1, gmael(B,sp[8],1+j5),Q);
1849 175 : av2 = avma;
1850 847 : for (j2 = 0; j2 < 4; j2++, set_avma(av2))
1851 : {
1852 686 : u2 = FpX_add(u1, gmael(B,sp[4],1+j2),Q);
1853 686 : u2 = FpX_add(u2, gmael(B,sp[5],1+smodss(-j0-j1-j2,4)),Q);
1854 686 : j4 = (-pk[1]-pk[2]+pj[3]+pj[4]-j2-w[i][1]*pj[1]-w[i][4]*pj[2])&3;
1855 686 : u2 = FpX_add(u2, gmael(B,sp[7],1+j4),Q);
1856 686 : u2 = FpX_add(u2, gmael(B,sp[9],1+smodss(-j3-j4-j5,4)),Q);
1857 686 : if (s4test(u2, liftpow, gl, rho))
1858 14 : goto suitef36_3;
1859 : }
1860 : }
1861 : }
1862 : }
1863 0 : return gc_NULL(ltop);
1864 14 : suitef36_3:
1865 14 : tam = perm_inv(tau);
1866 518 : for (i = 1; i <= n; i++)
1867 : {
1868 504 : r1[tau[i]] = rho[i];
1869 504 : r2[i] = tam[rho[i]];
1870 : }
1871 14 : set_avma(ltop2); return res;
1872 : }
1873 :
1874 : static GEN
1875 98 : galoisfindgroups(GEN lo, GEN sg, long f)
1876 : {
1877 98 : pari_sp ltop = avma;
1878 : long i, j, k;
1879 98 : GEN V = cgetg(lg(lo), t_VEC);
1880 287 : for(j=1,i=1; i<lg(lo); i++)
1881 : {
1882 189 : pari_sp av = avma;
1883 189 : GEN loi = gel(lo,i), W = cgetg(lg(loi),t_VECSMALL);
1884 476 : for (k=1; k<lg(loi); k++) W[k] = loi[k] % f;
1885 189 : W = vecsmall_uniq(W);
1886 189 : if (zv_equal(W, sg)) gel(V,j++) = loi;
1887 189 : set_avma(av);
1888 : }
1889 98 : setlg(V,j); return gerepilecopy(ltop,V);
1890 : }
1891 :
1892 : static GEN
1893 1715 : galoismakepsi(long g, GEN sg, GEN pf)
1894 : {
1895 1715 : GEN psi=cgetg(g+1,t_VECSMALL);
1896 : long i;
1897 4172 : for (i = 1; i < g; i++) psi[i] = sg[pf[i]];
1898 1715 : psi[g] = sg[1]; return psi;
1899 : }
1900 :
1901 : static GEN
1902 27594 : galoisfrobeniuslift_nilp(GEN T, GEN den, GEN L, GEN Lden,
1903 : struct galois_frobenius *gf, struct galois_borne *gb)
1904 : {
1905 27594 : pari_sp ltop=avma, av2;
1906 : struct galois_lift gl;
1907 27594 : long i, k, deg = 1, g = lg(gf->Tmod)-1;
1908 27594 : GEN F,Fp,Fe, aut, frob, res = cgetg(lg(L), t_VECSMALL);
1909 27594 : gf->psi = const_vecsmall(g,1);
1910 27594 : av2 = avma;
1911 27594 : initlift(T, den, gf->p, L, Lden, gb, &gl);
1912 27592 : if (DEBUGLEVEL >= 4)
1913 0 : err_printf("GaloisConj: p=%ld e=%ld deg=%ld fp=%ld\n",
1914 : gf->p, gl.e, deg, gf->fp);
1915 27592 : aut = galoisdolift(&gl);
1916 27593 : if (galoisfrobeniustest(aut,&gl,res))
1917 : {
1918 26213 : set_avma(av2); gf->deg = gf->fp; return res;
1919 : }
1920 :
1921 1379 : F =factoru(gf->fp);
1922 1379 : Fp = gel(F,1);
1923 1379 : Fe = gel(F,2);
1924 1379 : frob = cgetg(lg(L), t_VECSMALL);
1925 2758 : for(k = lg(Fp)-1; k>=1; k--)
1926 : {
1927 1379 : pari_sp btop=avma;
1928 1379 : GEN fres=NULL;
1929 1379 : long el = gf->fp, dg = 1, dgf = 1, e, pr;
1930 2702 : for(e=1; e<=Fe[k]; e++)
1931 : {
1932 2702 : dg *= Fp[k]; el /= Fp[k];
1933 2702 : if (DEBUGLEVEL>=4) err_printf("Trying degre %d.\n",dg);
1934 2702 : if (el==1) break;
1935 1449 : aut = galoisdoliftn(&gl, el);
1936 1449 : if (!galoisfrobeniustest(aut,&gl,frob))
1937 126 : break;
1938 1323 : dgf = dg; fres = gcopy(frob);
1939 : }
1940 1379 : if (dgf == 1) { set_avma(btop); continue; }
1941 1260 : pr = deg*dgf;
1942 1260 : if (deg == 1)
1943 : {
1944 16072 : for(i=1;i<lg(res);i++) res[i]=fres[i];
1945 : }
1946 : else
1947 : {
1948 0 : GEN cp = perm_mul(res,fres);
1949 0 : for(i=1;i<lg(res);i++) res[i] = cp[i];
1950 : }
1951 1260 : deg = pr; set_avma(btop);
1952 : }
1953 1379 : if (DEBUGLEVEL>=4 && res) err_printf("Best lift: %d\n",deg);
1954 1379 : if (deg==1) return gc_NULL(ltop);
1955 : else
1956 : {
1957 1260 : set_avma(av2); gf->deg = deg; return res;
1958 : }
1959 : }
1960 :
1961 : static GEN
1962 2212 : galoisfrobeniuslift(GEN T, GEN den, GEN L, GEN Lden,
1963 : struct galois_frobenius *gf, struct galois_borne *gb)
1964 : {
1965 2212 : pari_sp ltop=avma, av2;
1966 : struct galois_testlift gt;
1967 : struct galois_lift gl;
1968 2212 : long i, j, k, n = lg(L)-1, deg = 1, g = lg(gf->Tmod)-1;
1969 2212 : GEN F,Fp,Fe, aut, frob, res = cgetg(lg(L), t_VECSMALL);
1970 2212 : gf->psi = const_vecsmall(g,1);
1971 2212 : av2 = avma;
1972 2212 : initlift(T, den, gf->p, L, Lden, gb, &gl);
1973 2212 : if (DEBUGLEVEL >= 4)
1974 0 : err_printf("GaloisConj: p=%ld e=%ld deg=%ld fp=%ld\n",
1975 : gf->p, gl.e, deg, gf->fp);
1976 2212 : aut = galoisdolift(&gl);
1977 2212 : if (galoisfrobeniustest(aut,&gl,res))
1978 : {
1979 546 : set_avma(av2); gf->deg = gf->fp; return res;
1980 : }
1981 1666 : inittestlift(aut,gf->Tmod, &gl, >);
1982 1666 : gt.C = cgetg(gf->fp+1,t_VEC);
1983 1666 : gt.Cd= cgetg(gf->fp+1,t_VEC);
1984 9303 : for (i = 1; i <= gf->fp; i++) {
1985 7637 : gel(gt.C,i) = zero_zv(gt.g);
1986 7637 : gel(gt.Cd,i) = zero_zv(gt.g);
1987 : }
1988 :
1989 1666 : F =factoru(gf->fp);
1990 1666 : Fp = gel(F,1);
1991 1666 : Fe = gel(F,2);
1992 1666 : frob = cgetg(lg(L), t_VECSMALL);
1993 3528 : for(k=lg(Fp)-1;k>=1;k--)
1994 : {
1995 1862 : pari_sp btop=avma;
1996 1862 : GEN psi=NULL, fres=NULL, sg = identity_perm(1);
1997 1862 : long el=gf->fp, dg=1, dgf=1, e, pr;
1998 3773 : for(e=1; e<=Fe[k]; e++)
1999 : {
2000 : GEN lo, pf;
2001 : long l;
2002 1960 : dg *= Fp[k]; el /= Fp[k];
2003 1960 : if (DEBUGLEVEL>=4) err_printf("Trying degre %d.\n",dg);
2004 1960 : if (galoisfrobeniustest(gel(gt.pauto,el+1),&gl,frob))
2005 : {
2006 196 : psi = const_vecsmall(g,1);
2007 196 : dgf = dg; fres = gcopy(frob); continue;
2008 : }
2009 1764 : lo = listznstarelts(dg, n / gf->fp);
2010 1764 : if (e!=1) lo = galoisfindgroups(lo, sg, dgf);
2011 1764 : if (DEBUGLEVEL>=4) err_printf("Galoisconj:Subgroups list:%Ps\n", lo);
2012 3745 : for (l = 1; l < lg(lo); l++)
2013 3696 : if (lg(gel(lo,l))>2 && frobeniusliftall(gel(lo,l), el, &pf,&gl,>, frob))
2014 : {
2015 1715 : sg = gcopy(gel(lo,l));
2016 1715 : psi = galoismakepsi(g,sg,pf);
2017 1715 : dgf = dg; fres = gcopy(frob); break;
2018 : }
2019 1764 : if (l == lg(lo)) break;
2020 : }
2021 1862 : if (dgf == 1) { set_avma(btop); continue; }
2022 1827 : pr = deg*dgf;
2023 1827 : if (deg == 1)
2024 : {
2025 20454 : for(i=1;i<lg(res);i++) res[i]=fres[i];
2026 5719 : for(i=1;i<lg(psi);i++) gf->psi[i]=psi[i];
2027 : }
2028 : else
2029 : {
2030 161 : GEN cp = perm_mul(res,fres);
2031 3059 : for(i=1;i<lg(res);i++) res[i] = cp[i];
2032 525 : for(i=1;i<lg(psi);i++) gf->psi[i] = (dgf*gf->psi[i] + deg*psi[i]) % pr;
2033 : }
2034 1827 : deg = pr; set_avma(btop);
2035 : }
2036 9303 : for (i = 1; i <= gf->fp; i++)
2037 26425 : for (j = 1; j <= gt.g; j++) guncloneNULL(gmael(gt.C,i,j));
2038 1666 : if (DEBUGLEVEL>=4 && res) err_printf("Best lift: %d\n",deg);
2039 1666 : if (deg==1) return gc_NULL(ltop);
2040 : else
2041 : {
2042 : /* Normalize result so that psi[g]=1 */
2043 1666 : ulong im = Fl_inv(gf->psi[g], deg);
2044 1666 : GEN cp = perm_powu(res, im);
2045 20454 : for(i=1;i<lg(res);i++) res[i] = cp[i];
2046 5719 : for(i=1;i<lg(gf->psi);i++) gf->psi[i] = Fl_mul(im, gf->psi[i], deg);
2047 1666 : set_avma(av2); gf->deg = deg; return res;
2048 : }
2049 : }
2050 :
2051 : /* return NULL if not Galois */
2052 : static GEN
2053 29701 : galoisfindfrobenius(GEN T, GEN L, GEN den, GEN bad, struct galois_frobenius *gf,
2054 : struct galois_borne *gb, const struct galois_analysis *ga)
2055 : {
2056 29701 : pari_sp ltop = avma, av;
2057 29701 : long Try = 0, n = degpol(T), deg, gmask = (ga->group&ga_ext_2)? 3: 1;
2058 29700 : GEN frob, Lden = makeLden(L,den,gb);
2059 29700 : long is_nilpotent = ga->group&ga_all_nilpotent;
2060 : forprime_t S;
2061 :
2062 29700 : u_forprime_init(&S, ga->p, ULONG_MAX);
2063 29700 : av = avma;
2064 29700 : deg = gf->deg = ga->deg;
2065 29833 : while ((gf->p = u_forprime_next(&S)))
2066 : {
2067 : pari_sp lbot;
2068 : GEN Ti, Tp;
2069 : long nb, d;
2070 29833 : set_avma(av);
2071 29833 : Tp = ZX_to_Flx(T, gf->p);
2072 29832 : if (!Flx_is_squarefree(Tp, gf->p)) continue;
2073 29832 : if (bad && dvdiu(bad, gf->p)) continue;
2074 29832 : Ti = gel(Flx_factor(Tp, gf->p), 1);
2075 29834 : nb = lg(Ti)-1; d = degpol(gel(Ti,1));
2076 29834 : if (nb > 1 && degpol(gel(Ti,nb)) != d) return gc_NULL(ltop);
2077 29820 : if (((gmask&1)==0 || d % deg) && ((gmask&2)==0 || odd(d))) continue;
2078 29806 : if (DEBUGLEVEL >= 1) err_printf("GaloisConj: Trying p=%ld\n", gf->p);
2079 29806 : FlxV_to_ZXV_inplace(Ti);
2080 29806 : gf->fp = d;
2081 29806 : gf->Tmod = Ti; lbot = avma;
2082 29806 : if (is_nilpotent)
2083 27594 : frob = galoisfrobeniuslift_nilp(T, den, L, Lden, gf, gb);
2084 : else
2085 2212 : frob = galoisfrobeniuslift(T, den, L, Lden, gf, gb);
2086 29804 : if (frob)
2087 : {
2088 : GEN *gptr[3];
2089 29685 : gf->Tmod = gcopy(Ti);
2090 29686 : gptr[0]=&gf->Tmod; gptr[1]=&gf->psi; gptr[2]=&frob;
2091 29686 : gerepilemanysp(ltop,lbot,gptr,3); return frob;
2092 : }
2093 119 : if (is_nilpotent) continue;
2094 0 : if ((ga->group&ga_all_normal) && d % deg == 0) gmask &= ~1;
2095 : /* The first prime degree is always divisible by deg, so we don't
2096 : * have to worry about ext_2 being used before regular supersolvable*/
2097 0 : if (!gmask) return gc_NULL(ltop);
2098 0 : if ((ga->group&ga_non_wss) && ++Try > ((3*n)>>1))
2099 : {
2100 0 : if (DEBUGLEVEL)
2101 0 : pari_warn(warner,"Galois group probably not weakly super solvable");
2102 0 : return NULL;
2103 : }
2104 : }
2105 0 : if (!gf->p) pari_err_OVERFLOW("galoisfindfrobenius [ran out of primes]");
2106 0 : return NULL;
2107 : }
2108 :
2109 : /* compute g such that tau(Pmod[#])= tau(Pmod[g]) */
2110 :
2111 : static long
2112 6594 : get_image(GEN tau, GEN P, GEN Pmod, GEN p)
2113 : {
2114 6594 : pari_sp av = avma;
2115 6594 : long g, gp = lg(Pmod)-1;
2116 6594 : tau = RgX_to_FpX(tau, p);
2117 6594 : tau = FpX_FpXQ_eval(gel(Pmod, gp), tau, P, p);
2118 6594 : tau = FpX_normalize(FpX_gcd(P, tau, p), p);
2119 10695 : for (g = 1; g <= gp; g++)
2120 10695 : if (ZX_equal(tau, gel(Pmod,g))) return gc_long(av,g);
2121 0 : return gc_long(av,0);
2122 : }
2123 :
2124 : static GEN
2125 33928 : gg_get_std(GEN G)
2126 : {
2127 33928 : return !G ? NULL: lg(G)==3 ? G: mkvec2(gel(G,1),gmael(G,5,1));
2128 : }
2129 :
2130 : static GEN galoisgen(GEN T, GEN L, GEN M, GEN den, GEN bad, struct galois_borne *gb,
2131 : const struct galois_analysis *ga);
2132 :
2133 : static GEN
2134 5390 : galoisgenfixedfield(GEN Tp, GEN Pmod, GEN PL, GEN P, GEN ip, GEN bad, struct galois_borne *gb)
2135 : {
2136 : GEN Pden, PM;
2137 : GEN tau, PG, Pg;
2138 : long g, lP;
2139 5390 : long x = varn(Tp);
2140 5390 : GEN Pp = FpX_red(P, ip);
2141 5390 : if (DEBUGLEVEL>=6)
2142 0 : err_printf("GaloisConj: Fixed field %Ps\n",P);
2143 5390 : if (degpol(P)==2 && !bad)
2144 : {
2145 4046 : PG=cgetg(3,t_VEC);
2146 4046 : gel(PG,1) = mkvec( mkvecsmall2(2,1) );
2147 4046 : gel(PG,2) = mkvecsmall(2);
2148 4046 : tau = deg1pol_shallow(gen_m1, negi(gel(P,3)), x);
2149 4046 : g = get_image(tau, Pp, Pmod, ip);
2150 4046 : if (!g) return NULL;
2151 4046 : Pg = mkvecsmall(g);
2152 : }
2153 : else
2154 : {
2155 : struct galois_analysis Pga;
2156 : struct galois_borne Pgb;
2157 : GEN mod, mod2;
2158 : long j;
2159 1351 : if (!galoisanalysis(P, &Pga, 0, NULL)) return NULL;
2160 1330 : if (bad) Pga.group &= ~ga_easy;
2161 1330 : Pgb.l = gb->l;
2162 1330 : Pden = galoisborne(P, NULL, &Pgb, degpol(P));
2163 :
2164 1330 : if (Pgb.valabs > gb->valabs)
2165 : {
2166 125 : if (DEBUGLEVEL>=4)
2167 0 : err_printf("GaloisConj: increase prec of p-adic roots of %ld.\n"
2168 0 : ,Pgb.valabs-gb->valabs);
2169 125 : PL = ZpX_liftroots(P,PL,gb->l,Pgb.valabs);
2170 : }
2171 1205 : else if (Pgb.valabs < gb->valabs)
2172 1123 : PL = FpC_red(PL, Pgb.ladicabs);
2173 1330 : PM = FpV_invVandermonde(PL, Pden, Pgb.ladicabs);
2174 1330 : PG = galoisgen(P, PL, PM, Pden, bad ? lcmii(Pgb.dis, bad): NULL, &Pgb, &Pga);
2175 1330 : if (!PG) return NULL;
2176 1323 : lP = lg(gel(PG,1));
2177 1323 : mod = Pgb.ladicabs; mod2 = shifti(mod, -1);
2178 1323 : Pg = cgetg(lP, t_VECSMALL);
2179 3871 : for (j = 1; j < lP; j++)
2180 : {
2181 2548 : pari_sp btop=avma;
2182 2548 : tau = permtopol(gmael(PG,1,j), PL, PM, Pden, mod, mod2, x);
2183 2548 : g = get_image(tau, Pp, Pmod, ip);
2184 2548 : if (!g) return NULL;
2185 2548 : Pg[j] = g;
2186 2548 : set_avma(btop);
2187 : }
2188 : }
2189 5369 : return mkvec2(PG,Pg);
2190 : }
2191 :
2192 : static GEN
2193 5390 : galoisgenfixedfield0(GEN O, GEN L, GEN sigma, GEN T, GEN bad, GEN *pt_V,
2194 : struct galois_frobenius *gf, struct galois_borne *gb)
2195 : {
2196 5390 : pari_sp btop = avma;
2197 5390 : long vT = varn(T);
2198 5390 : GEN mod = gb->ladicabs, mod2 = shifti(gb->ladicabs,-1);
2199 : GEN OL, sym, P, PL, p, Tp, Sp, Pmod, PG;
2200 5390 : OL = fixedfieldorbits(O,L);
2201 5390 : sym = fixedfieldsympol(OL, itou(gb->l));
2202 5390 : PL = sympol_eval(sym, OL, mod);
2203 5390 : P = FpX_center_i(FpV_roots_to_pol(PL, mod, vT), mod, mod2);
2204 5390 : if (!FpX_is_squarefree(P,utoipos(gf->p)))
2205 : {
2206 93 : GEN badp = lcmii(bad? bad: gb->dis, ZX_disc(P));
2207 93 : gf->p = findpsi(badp, gf->p, T, sigma, gf->deg, &gf->Tmod, &gf->psi);
2208 : }
2209 5390 : p = utoipos(gf->p);
2210 5390 : Tp = FpX_red(T,p);
2211 5390 : Sp = sympol_aut_evalmod(sym, gf->deg, sigma, Tp, p);
2212 5390 : Pmod = fixedfieldfactmod(Sp, p, gf->Tmod);
2213 5390 : PG = galoisgenfixedfield(Tp, Pmod, PL, P, p, bad, gb);
2214 5390 : if (PG == NULL) return NULL;
2215 5369 : if (DEBUGLEVEL >= 4)
2216 0 : err_printf("GaloisConj: Back to Earth:%Ps\n", gg_get_std(gel(PG,1)));
2217 5369 : if (pt_V) *pt_V = mkvec3(sym, PL, P);
2218 5369 : return gc_all(btop, pt_V ? 4: 3, &PG, &gf->Tmod, &gf->psi, pt_V);
2219 : }
2220 :
2221 : /* Let sigma^m=1, tau*sigma*tau^-1=sigma^s. Return n = sum_{0<=k<e,0} s^k mod m
2222 : * so that (sigma*tau)^e = sigma^n*tau^e. N.B. n*(1-s) = 1-s^e mod m,
2223 : * unfortunately (1-s) may not invertible mod m */
2224 : static long
2225 14672 : stpow(long s, long e, long m)
2226 : {
2227 14672 : long i, n = 1;
2228 22946 : for (i = 1; i < e; i++) n = (1 + n * s) % m;
2229 14672 : return n;
2230 : }
2231 :
2232 : static GEN
2233 6594 : wpow(long s, long m, long e, long n)
2234 : {
2235 6594 : GEN w = cgetg(n+1,t_VECSMALL);
2236 6594 : long si = s;
2237 : long i;
2238 6594 : w[1] = 1;
2239 7336 : for(i=2; i<=n; i++) w[i] = w[i-1]*e;
2240 13930 : for(i=n; i>=1; i--)
2241 : {
2242 7336 : si = Fl_powu(si,e,m);
2243 7336 : w[i] = Fl_mul(s-1, stpow(si, w[i], m), m);
2244 : }
2245 6594 : return w;
2246 : }
2247 :
2248 : static GEN
2249 6594 : galoisgenliftauto(GEN O, GEN gj, long s, long n, struct galois_test *td)
2250 : {
2251 6594 : pari_sp av = avma;
2252 : long sr, k;
2253 6594 : long deg = lg(gel(O,1))-1;
2254 6594 : GEN X = cgetg(lg(O), t_VECSMALL);
2255 6594 : GEN oX = cgetg(lg(O), t_VECSMALL);
2256 6594 : GEN B = perm_cycles(gj);
2257 6594 : long oj = lg(gel(B,1)) - 1;
2258 6594 : GEN F = factoru(oj);
2259 6594 : GEN Fp = gel(F,1);
2260 6594 : GEN Fe = gel(F,2);
2261 6594 : GEN pf = identity_perm(n);
2262 6594 : if (DEBUGLEVEL >= 6)
2263 0 : err_printf("GaloisConj: %Ps of relative order %d\n", gj, oj);
2264 12425 : for (k=lg(Fp)-1; k>=1; k--)
2265 : {
2266 6594 : long f, dg = 1, el = oj, osel = 1, a = 0;
2267 6594 : long p = Fp[k], e = Fe[k], op = oj / upowuu(p,e);
2268 : long i;
2269 6594 : GEN pf1 = NULL, w, wg, Be = cgetg(e+1,t_VEC);
2270 6594 : gel(Be,e) = cyc_pow(B, op);
2271 7336 : for(i=e-1; i>=1; i--) gel(Be,i) = cyc_pow(gel(Be,i+1), p);
2272 6594 : w = wpow(Fl_powu(s,op,deg),deg,p,e);
2273 6594 : wg = cgetg(e+2,t_VECSMALL);
2274 6594 : wg[e+1] = deg;
2275 13930 : for (i=e; i>=1; i--) wg[i] = ugcd(wg[i+1], w[i]);
2276 36645 : for (i=1; i<lg(O); i++) oX[i] = 0;
2277 13167 : for (f=1; f<=e; f++)
2278 : {
2279 : long sel, t;
2280 7336 : GEN Bel = gel(Be,f);
2281 7336 : dg *= p; el /= p;
2282 7336 : sel = Fl_powu(s,el,deg);
2283 7336 : if (DEBUGLEVEL >= 6) err_printf("GaloisConj: B=%Ps\n", Bel);
2284 7336 : sr = ugcd(stpow(sel,p,deg),deg);
2285 7336 : if (DEBUGLEVEL >= 6)
2286 0 : err_printf("GaloisConj: exp %d: s=%ld [%ld] a=%ld w=%ld wg=%ld sr=%ld\n",
2287 0 : dg, sel, deg, a, w[f], wg[f+1], sr);
2288 9758 : for (t = 0; t < sr; t++)
2289 8995 : if ((a+t*w[f])%wg[f+1]==0)
2290 : {
2291 : long i, j, k, st;
2292 58744 : for (i = 1; i < lg(X); i++) X[i] = 0;
2293 30870 : for (i = 0; i < lg(X)-1; i+=dg)
2294 46732 : for (j = 1, k = p, st = t; k <= dg; j++, k += p)
2295 : {
2296 24787 : X[k+i] = (oX[j+i] + st)%deg;
2297 24787 : st = (t + st*osel)%deg;
2298 : }
2299 8925 : pf1 = testpermutation(O, Bel, X, sel, p, sr, td);
2300 8925 : if (pf1) break;
2301 : }
2302 7336 : if (!pf1) return NULL;
2303 42882 : for (i=1; i<lg(O); i++) oX[i] = X[i];
2304 6573 : osel = sel; a = (a+t*w[f])%deg;
2305 : }
2306 5831 : pf = perm_mul(pf, perm_powu(pf1, el));
2307 : }
2308 5831 : return gerepileuptoleaf(av, pf);
2309 : }
2310 :
2311 : static GEN
2312 0 : FlxV_Flx_gcd(GEN x, GEN T, ulong p)
2313 0 : { pari_APPLY_same(Flx_normalize(Flx_gcd(gel(x,i),T,p),p)) }
2314 :
2315 : static GEN
2316 0 : Flx_FlxV_minpolymod(GEN y, GEN x, ulong p)
2317 0 : { pari_APPLY_same(Flxq_minpoly(Flx_rem(y, gel(x,i), p), gel(x,i), p)) }
2318 :
2319 : static GEN
2320 0 : FlxV_minpolymod(GEN x, GEN y, ulong p)
2321 0 : { pari_APPLY_same(Flx_FlxV_minpolymod(gel(x,i), y, p)) }
2322 :
2323 : static GEN
2324 0 : factperm(GEN x)
2325 : {
2326 0 : pari_APPLY_same(gen_indexsort(gel(x,i), (void*)cmp_Flx, cmp_nodata))
2327 : }
2328 :
2329 : /* compute (prod p_i^e_i)(1) */
2330 :
2331 : static long
2332 0 : permprodeval(GEN p, GEN e, long s)
2333 : {
2334 0 : long i, j, l = lg(p);
2335 0 : for (i=l-1; i>=1; i--)
2336 : {
2337 0 : GEN pi = gel(p,i);
2338 0 : long ei = uel(e,i);
2339 0 : for(j = 1; j <= ei; j++)
2340 0 : s = uel(pi, s);
2341 : }
2342 0 : return s;
2343 : }
2344 :
2345 : static GEN
2346 0 : pc_to_perm(GEN pc, GEN gen, long n)
2347 : {
2348 0 : long i, l = lg(pc);
2349 0 : GEN s = identity_perm(n);
2350 0 : for (i=1; i<l; i++)
2351 0 : s = perm_mul(gel(gen,pc[i]),s);
2352 0 : return s;
2353 : }
2354 :
2355 : static GEN
2356 0 : genorbit(GEN ordH, GEN permfact_Hp, long fr, long n, long k, long j)
2357 : {
2358 0 : pari_sp av = avma;
2359 0 : long l = lg(gel(permfact_Hp,1))-1, no = 1, b, i;
2360 0 : GEN W = zero_zv(l);
2361 0 : GEN orb = cgetg(l+1, t_VECSMALL);
2362 0 : GEN gen = cgetg(l+1, t_VEC);
2363 0 : GEN E = cgetg(k+1, t_VECSMALL);
2364 0 : for(b = 0; b < n; b++)
2365 : {
2366 0 : long bb = b, s;
2367 0 : for(i = 1; i <= k; i++)
2368 : {
2369 0 : uel(E,i) = bb % uel(ordH,i);
2370 0 : bb /= uel(ordH,i);
2371 : }
2372 0 : if (E[j]) continue;
2373 0 : s = permprodeval(permfact_Hp, E, fr);
2374 0 : if (s>lg(W)-1) pari_err_BUG("W1");
2375 0 : if (W[s]) continue;
2376 0 : W[s] = 1;
2377 0 : if (no > l) pari_err_BUG("genorbit");
2378 0 : uel(orb,no) = s;
2379 0 : gel(gen,no) = zv_copy(E);
2380 0 : no++;
2381 : }
2382 0 : if(no<l) pari_err_BUG("genorbit");
2383 0 : return gerepilecopy(av, mkvec2(orb,gen));
2384 : }
2385 :
2386 0 : INLINE GEN br_get(GEN br, long i, long j) { return gmael(br,j,i-j); }
2387 0 : static GEN pcgrp_get_ord(GEN G) { return gel(G,1); }
2388 0 : static GEN pcgrp_get_pow(GEN G) { return gel(G,2); }
2389 0 : static GEN pcgrp_get_br(GEN G) { return gel(G,3); }
2390 :
2391 : static GEN
2392 24296 : cyclic_pc(long n)
2393 : {
2394 24296 : return mkvec3(mkvecsmall(n),mkvec(cgetg(1,t_VECSMALL)), mkvec(cgetg(1,t_VEC)));
2395 : }
2396 :
2397 : static GEN
2398 0 : pc_normalize(GEN g, GEN G)
2399 : {
2400 0 : long i, l = lg(g)-1, o = 1;
2401 0 : GEN ord = pcgrp_get_ord(G), pw = pcgrp_get_pow(G), br = pcgrp_get_br(G);
2402 0 : for (i = 1; i < l; i++)
2403 : {
2404 0 : if (g[i] == g[i+1])
2405 : {
2406 0 : if (++o == ord[g[i]])
2407 : {
2408 0 : GEN v = vecsmall_concat(vecslice(g,1,i-o+1),gel(pw,g[i]));
2409 0 : GEN w = vecsmall_concat(v,vecslice(g,i+2,l));
2410 0 : return pc_normalize(w, G);
2411 : }
2412 : }
2413 0 : else if (g[i] > g[i+1])
2414 : {
2415 0 : GEN v = vecsmall_concat(vecslice(g,1,i-1), br_get(br,g[i],g[i+1]));
2416 0 : GEN w = vecsmall_concat(mkvecsmall2(g[i+1],g[i]),vecslice(g,i+2,l));
2417 0 : v = vecsmall_concat(v, w);
2418 0 : return pc_normalize(v, G);
2419 : }
2420 0 : else o = 1;
2421 : }
2422 0 : return g;
2423 : }
2424 :
2425 : static GEN
2426 0 : pc_inv(GEN g, GEN G)
2427 : {
2428 0 : long i, l = lg(g);
2429 0 : GEN ord = pcgrp_get_ord(G), pw = pcgrp_get_pow(G);
2430 0 : GEN v = cgetg(l, t_VEC);
2431 0 : if (l==1) return v;
2432 0 : for(i = 1; i < l; i++)
2433 : {
2434 0 : ulong gi = uel(g,i);
2435 0 : gel(v,l-i) = vecsmall_concat(pc_inv(gel(pw, gi), G),
2436 0 : const_vecsmall(uel(ord,gi)-1,gi));
2437 : }
2438 0 : return pc_normalize(shallowconcat1(v), G);
2439 : }
2440 :
2441 : static GEN
2442 0 : pc_mul(GEN g, GEN h, GEN G)
2443 : {
2444 0 : return pc_normalize(vecsmall_concat(g,h), G);
2445 : }
2446 :
2447 : static GEN
2448 0 : pc_bracket(GEN g, GEN h, GEN G)
2449 : {
2450 0 : GEN gh = pc_mul(g, h, G);
2451 0 : GEN hg = pc_mul(h, g, G);
2452 0 : long i, l1 = lg(gh), l2 = lg(hg), lm = minss(l1,l2);
2453 0 : for (i = 1; i < lm; i++)
2454 0 : if (gh[l1-i] != hg[l2-i]) break;
2455 0 : return pc_mul(vecsmall_shorten(gh,l1-i), pc_inv(vecsmall_shorten(hg,l2-i), G), G);
2456 : }
2457 :
2458 : static GEN
2459 0 : pc_exp(GEN v)
2460 : {
2461 0 : long i, l = lg(v);
2462 0 : GEN w = cgetg(l, t_VEC);
2463 0 : if (l==1) return w;
2464 0 : for (i = 1; i < l; i++)
2465 0 : gel(w,i) = const_vecsmall(v[i], i+1);
2466 0 : return shallowconcat1(w);
2467 : }
2468 : static GEN
2469 0 : vecsmall_increase(GEN x)
2470 0 : { pari_APPLY_ulong(x[i]+1) }
2471 :
2472 : static GEN
2473 0 : vecvecsmall_increase(GEN x)
2474 0 : { pari_APPLY_same(vecsmall_increase(gel(x,i))) }
2475 :
2476 : static GEN
2477 0 : pcgrp_lift(GEN G, long deg)
2478 : {
2479 0 : GEN ord = pcgrp_get_ord(G), pw = pcgrp_get_pow(G), br = pcgrp_get_br(G);
2480 0 : long i, l = lg(br);
2481 0 : GEN Ord = vecsmall_prepend(ord, deg);
2482 0 : GEN Pw = vec_prepend(vecvecsmall_increase(pw), cgetg(1,t_VECSMALL));
2483 0 : GEN Br = cgetg(l+1, t_VEC);
2484 0 : gel(Br,1) = const_vec(l-1, cgetg(1, t_VECSMALL));
2485 0 : for (i = 1; i < l; i++)
2486 0 : gel(Br,i+1) = vecvecsmall_increase(gel(br, i));
2487 0 : return mkvec3(Ord, Pw, Br);
2488 : }
2489 :
2490 : static GEN
2491 0 : brl_add(GEN x, GEN a)
2492 : {
2493 0 : pari_APPLY_same(vecsmall_concat(const_vecsmall(uel(a,i),1),gel(x,i)))
2494 : }
2495 :
2496 : static void
2497 0 : pcgrp_insert(GEN G, long j, GEN a)
2498 : {
2499 0 : GEN pw = pcgrp_get_pow(G), br = pcgrp_get_br(G);
2500 0 : gel(pw,j) = vecsmall_concat(gel(a,1),gel(pw, j));
2501 0 : gel(br,j) = brl_add(gel(br, j), gel(a,2));
2502 0 : }
2503 :
2504 : static long
2505 0 : getfr(GEN f, GEN h)
2506 : {
2507 0 : long i, l = lg(f);
2508 0 : for (i = 1; i < l; i++)
2509 0 : if (zv_equal(gel(f,i), h)) return i;
2510 0 : pari_err_BUG("galoisinit");
2511 0 : return 0;
2512 : }
2513 :
2514 : static long
2515 0 : get_pow(GEN pf, ulong o, GEN pw, GEN gen)
2516 : {
2517 0 : long i, n = lg(pf)-1;
2518 0 : GEN p1 = perm_powu(pf, o);
2519 0 : GEN p2 = pc_to_perm(pw, gen, n);
2520 0 : for(i = 0; ; i++)
2521 : {
2522 0 : if (zv_equal(p1, p2)) break;
2523 0 : p2 = perm_mul(gel(gen,1), p2);
2524 : }
2525 0 : return i;
2526 : }
2527 :
2528 : struct galois_perm
2529 : {
2530 : GEN L;
2531 : GEN M;
2532 : GEN den;
2533 : GEN mod, mod2;
2534 : long x;
2535 : GEN cache;
2536 : };
2537 :
2538 : static void
2539 0 : galoisperm_init(struct galois_perm *gp, GEN L, GEN M, GEN den, GEN mod, GEN mod2, long x)
2540 : {
2541 0 : gp->L = L;
2542 0 : gp->M = M;
2543 0 : gp->den = den;
2544 0 : gp->mod = mod;
2545 0 : gp->mod2 = mod2;
2546 0 : gp->x = x;
2547 0 : gp->cache = zerovec(lg(L)-1);
2548 0 : }
2549 :
2550 : static void
2551 0 : galoisperm_free(struct galois_perm *gp)
2552 : {
2553 0 : long i, l = lg(gp->cache);
2554 0 : for (i=1; i<l; i++)
2555 0 : if (!isintzero(gel(gp->cache,i)))
2556 0 : gunclone(gel(gp->cache,i));
2557 0 : }
2558 :
2559 : static GEN
2560 0 : permtoaut(GEN p, struct galois_perm *gp)
2561 : {
2562 0 : pari_sp av = avma;
2563 0 : if (isintzero(gel(gp->cache,p[1])))
2564 : {
2565 0 : GEN pol = permtopol(p, gp->L, gp->M, gp->den, gp->mod, gp->mod2, gp->x);
2566 0 : gel(gp->cache,p[1]) = gclone(pol);
2567 : }
2568 0 : set_avma(av);
2569 0 : return gel(gp->cache,p[1]);
2570 : }
2571 :
2572 : static GEN
2573 0 : pc_evalcache(GEN W, GEN u, GEN sp, GEN T, GEN p, struct galois_perm *gp)
2574 : {
2575 : GEN v;
2576 0 : long ns = sp[1];
2577 0 : if (!isintzero(gel(W,ns))) return gel(W,ns);
2578 0 : v = RgX_to_FpX(permtoaut(sp, gp), p);
2579 0 : gel(W,ns) = FpX_FpXQV_eval(v, u, T, p);
2580 0 : return gel(W,ns);
2581 : }
2582 :
2583 : static ulong
2584 0 : findp(GEN D, GEN P, GEN S, long o, GEN *Tmod)
2585 : {
2586 : forprime_t iter;
2587 : ulong p;
2588 0 : long n = degpol(P);
2589 0 : u_forprime_init(&iter, n*maxss(expu(n)-3, 2), ULONG_MAX);
2590 0 : while ((p = u_forprime_next(&iter)))
2591 : {
2592 : GEN F, F1, Sp;
2593 0 : if (smodis(D, p) == 0)
2594 0 : continue;
2595 0 : F = gel(Flx_factor(ZX_to_Flx(P, p), p), 1);
2596 0 : F1 = gel(F,1);
2597 0 : if (degpol(F1) != o)
2598 0 : continue;
2599 0 : Sp = RgX_to_Flx(S, p);
2600 0 : if (gequal(Flx_rem(Sp, F1, p), Flx_Frobenius(F1, p)))
2601 : {
2602 0 : *Tmod = FlxV_to_ZXV(F);
2603 0 : return p;
2604 : }
2605 : }
2606 0 : return 0;
2607 : }
2608 :
2609 : static GEN
2610 0 : nilp_froblift(GEN genG, GEN autH, long j, GEN pcgrp,
2611 : GEN idp, GEN incl, GEN H, struct galois_lift *gl, struct galois_perm *gp)
2612 : {
2613 0 : pari_sp av = avma;
2614 0 : GEN T = gl->T, p = gl->p, pe = gl->Q;
2615 0 : ulong pp = itou(p);
2616 0 : long e = gl->e;
2617 0 : GEN pf = cgetg(lg(gl->L), t_VECSMALL);
2618 0 : GEN Tp = ZX_to_Flx(T, pp);
2619 0 : GEN Hp = ZX_to_Flx(H, pp);
2620 0 : GEN ord = pcgrp_get_ord(pcgrp);
2621 0 : GEN pcp = gel(pcgrp_get_pow(pcgrp),j+1);
2622 0 : long o = uel(ord,1);
2623 0 : GEN ordH = vecslice(ord,2,lg(ord)-1);
2624 0 : long n = zv_prod(ordH), k = lg(ordH)-1, l = k-j, m = upowuu(o, l), v = varn(T);
2625 0 : GEN factTp = gel(Flx_factor(Tp, pp), 1);
2626 0 : long fp = degpol(gel(factTp, 1));
2627 0 : GEN frobp = Flxq_autpow(Flx_Frobenius(Tp, pp), fp-1, Tp, pp);
2628 0 : GEN frob = ZpX_ZpXQ_liftroot(T, Flx_to_ZX(frobp), T, p, e);
2629 0 : if (galoisfrobeniustest(frob, gl, pf))
2630 : {
2631 0 : GEN pfi = perm_inv(pf);
2632 0 : long d = get_pow(pfi, uel(ord,j+1), pcp, genG);
2633 0 : return mkvec3(pfi, mkvec2(const_vecsmall(d,1),zero_zv(l+1)), gel(factTp, 1));
2634 : }
2635 : else
2636 : {
2637 0 : GEN frobG = FpXQ_powers(frob, usqrt(degpol(T)), T, pe);
2638 0 : GEN autHp = RgXV_to_FlxV(autH,pp);
2639 0 : GEN inclp = RgX_to_Flx(incl,pp);
2640 0 : GEN factHp = gel(Flx_factor(Hp, pp),1);
2641 0 : long fr = getfr(factHp, idp);
2642 0 : GEN minHp = FlxV_minpolymod(autHp, factHp, pp);
2643 0 : GEN permfact_Hp = factperm(minHp);
2644 0 : GEN permfact_Gp = FlxV_Flx_gcd(FlxC_Flxq_eval(factHp, inclp, Tp, pp), Tp, pp);
2645 0 : GEN bezout_Gpe = bezout_lift_fact(T, FlxV_to_ZXV(permfact_Gp), p, e);
2646 0 : GEN id = gmael(Flx_factor(gel(permfact_Gp, fr),pp),1,1);
2647 0 : GEN orbgen = genorbit(ordH, permfact_Hp, fr, n, k, j);
2648 0 : GEN orb = gel(orbgen,1), gen = gel(orbgen,2);
2649 0 : long nborb = lg(orb)-1;
2650 0 : GEN A = cgetg(l+1, t_VECSMALL);
2651 0 : GEN W = zerovec(lg(gl->L)-1);
2652 0 : GEN U = zeromatcopy(nborb,degpol(T));
2653 0 : GEN br = pcgrp_get_br(pcgrp), brj = gcopy(gel(br, j+1));
2654 0 : GEN Ui = cgetg(nborb+1, t_VEC);
2655 : long a, b, i;
2656 0 : for(a = 0; a < m; a++)
2657 : {
2658 : pari_timer ti;
2659 : pari_sp av2;
2660 0 : GEN B = pol_0(v);
2661 0 : long aa = a;
2662 0 : if (DEBUGLEVEL>=4) timer_start(&ti);
2663 0 : for(i = 1; i <= l; i++)
2664 : {
2665 0 : uel(A,i) = aa % o;
2666 0 : aa /= o;
2667 : }
2668 0 : gel(br,j+1) = brl_add(brj, A);
2669 0 : for(b = 1; b <= nborb; b++)
2670 : {
2671 0 : GEN br = pc_bracket(pc_exp(gel(gen,b)), mkvecsmall(j+1), pcgrp);
2672 0 : GEN sp = pc_to_perm(br, genG, degpol(T));
2673 0 : long u = sp[1];
2674 0 : long s = permprodeval(permfact_Hp, gel(gen,b), fr);
2675 0 : if (isintzero(gmael(U,u,s)))
2676 : {
2677 0 : GEN Ub = pc_evalcache(W, frobG, sp, T, pe, gp);
2678 0 : gmael(U,u,s) = FpXQ_mul(Ub, gel(bezout_Gpe,s), T, pe);
2679 : }
2680 0 : gel(Ui, b) = gmael(U,u,s);
2681 : }
2682 0 : av2 = avma;
2683 0 : for(b = 1; b <= nborb; b++)
2684 0 : B = FpX_add(B, gel(Ui,b), pe);
2685 0 : if (DEBUGLEVEL >= 4) timer_printf(&ti,"Testing candidate %ld",a);
2686 0 : if (galoisfrobeniustest(B, gl, pf))
2687 : {
2688 0 : GEN pfi = perm_inv(pf);
2689 0 : long d = get_pow(pfi, uel(ord,j+1), pcp, genG);
2690 0 : gel(br,j+1) = brj;
2691 0 : return gerepilecopy(av,mkvec3(pfi,mkvec2(const_vecsmall(d,1),A),id));
2692 : }
2693 0 : set_avma(av2);
2694 : }
2695 0 : return gc_NULL(av);
2696 : }
2697 : }
2698 :
2699 : static GEN
2700 0 : galoisgenlift_nilp(GEN PG, GEN O, GEN V, GEN T, GEN frob, GEN sigma,
2701 : struct galois_borne *gb, struct galois_frobenius *gf, struct galois_perm *gp)
2702 : {
2703 0 : long j, n = degpol(T), deg = gf->deg;
2704 0 : ulong p = gf->p;
2705 0 : GEN L = gp->L, M = gp->M, den = gp->den;
2706 0 : GEN S = fixedfieldinclusion(O, gel(V,2));
2707 0 : GEN incl = vectopol(S, M, den, gb->ladicabs, shifti(gb->ladicabs,-1), varn(T));
2708 0 : GEN H = gel(V,3);
2709 0 : GEN PG1 = gmael(PG, 1, 1);
2710 0 : GEN PG2 = gmael(PG, 1, 2);
2711 0 : GEN PG3 = gmael(PG, 1, 3);
2712 0 : GEN PG4 = gmael(PG, 1, 4);
2713 0 : long lP = lg(PG1);
2714 0 : GEN PG5 = pcgrp_lift(gmael(PG, 1, 5), deg);
2715 0 : GEN res = cgetg(6, t_VEC), res1, res2, res3;
2716 0 : gel(res,1) = res1 = cgetg(lP + 1, t_VEC);
2717 0 : gel(res,2) = res2 = cgetg(lP + 1, t_VEC);
2718 0 : gel(res,3) = res3 = cgetg(lP + 1, t_VEC);
2719 0 : gel(res,4) = vecsmall_prepend(PG4, p);
2720 0 : gel(res,5) = PG5;
2721 0 : gel(res1, 1) = frob;
2722 0 : gel(res2, 1) = ZX_to_Flx(gel(gf->Tmod,1), p);
2723 0 : gel(res3, 1) = sigma;
2724 0 : for (j = 1; j < lP; j++)
2725 : {
2726 : struct galois_lift gl;
2727 0 : GEN Lden = makeLden(L,den,gb);
2728 : GEN pf;
2729 0 : initlift(T, den, uel(PG4,j), L, Lden, gb, &gl);
2730 0 : pf = nilp_froblift(vecslice(res1,1,j), PG3, j, PG5, gel(PG2,j), incl, H, &gl, gp);
2731 0 : if (!pf) return NULL;
2732 0 : if (DEBUGLEVEL>=2)
2733 0 : err_printf("found: %ld/%ld: %Ps: %Ps\n", n, j+1, gel(pf,2),gel(pf,1));
2734 0 : pcgrp_insert(PG5, j+1, gel(pf,2));
2735 0 : gel(res1, j+1) = gel(pf,1);
2736 0 : gel(res2, j+1) = gel(pf,3);
2737 0 : gel(res3, j+1) = gcopy(permtoaut(gel(pf,1), gp));
2738 : }
2739 0 : if (DEBUGLEVEL >= 4) err_printf("GaloisConj: Fini!\n");
2740 0 : return res;
2741 : }
2742 :
2743 : static GEN
2744 5369 : galoisgenlift(GEN PG, GEN Pg, GEN O, GEN L, GEN M, GEN frob,
2745 : struct galois_borne *gb, struct galois_frobenius *gf)
2746 : {
2747 : struct galois_test td;
2748 : GEN res, res1;
2749 5369 : GEN PG1 = gel(PG, 1), PG2 = gel(PG, 2);
2750 5369 : long lP = lg(PG1), j, n = lg(L)-1;
2751 5369 : inittest(L, M, gb->bornesol, gb->ladicsol, &td);
2752 5369 : res = cgetg(3, t_VEC);
2753 5369 : gel(res,1) = res1 = cgetg(lP + 1, t_VEC);
2754 5369 : gel(res,2) = vecsmall_prepend(PG2, gf->deg);
2755 5369 : gel(res1, 1) = vecsmall_copy(frob);
2756 11200 : for (j = 1; j < lP; j++)
2757 : {
2758 6594 : GEN pf = galoisgenliftauto(O, gel(PG1, j), gf->psi[Pg[j]], n, &td);
2759 6594 : if (!pf) { freetest(&td); return NULL; }
2760 5831 : gel(res1, j+1) = pf;
2761 : }
2762 4606 : if (DEBUGLEVEL >= 4) err_printf("GaloisConj: Fini!\n");
2763 4606 : freetest(&td);
2764 4606 : return res;
2765 : }
2766 :
2767 : static ulong
2768 29685 : psi_order(GEN psi, ulong d)
2769 : {
2770 29685 : long i, l = lg(psi);
2771 29685 : ulong s = 1;
2772 66636 : for (i=1; i<l; i++)
2773 36951 : s = clcm(s, d/cgcd(uel(psi, i)-1, d));
2774 29685 : return s;
2775 : }
2776 :
2777 : static GEN
2778 29890 : galoisgen(GEN T, GEN L, GEN M, GEN den, GEN bad, struct galois_borne *gb,
2779 : const struct galois_analysis *ga)
2780 : {
2781 : struct galois_test td;
2782 : struct galois_frobenius gf, ogf;
2783 29890 : pari_sp ltop = avma;
2784 29890 : long x, n = degpol(T), is_central;
2785 : ulong po;
2786 29890 : GEN sigma, res, frob, O, PG, V, ofrob = NULL;
2787 :
2788 29890 : if (!ga->deg) return NULL;
2789 29890 : x = varn(T);
2790 29890 : if (DEBUGLEVEL >= 9) err_printf("GaloisConj: denominator:%Ps\n", den);
2791 29890 : if (n == 12 && ga->ord==3 && !ga->p4)
2792 : { /* A4 is very probable: test it first */
2793 98 : pari_sp av = avma;
2794 98 : if (DEBUGLEVEL >= 4) err_printf("GaloisConj: Testing A4 first\n");
2795 98 : inittest(L, M, gb->bornesol, gb->ladicsol, &td);
2796 98 : PG = a4galoisgen(&td);
2797 98 : freetest(&td);
2798 98 : if (PG) return gerepileupto(ltop, PG);
2799 0 : set_avma(av);
2800 : }
2801 29792 : if (n == 24 && ga->ord==3 && ga->p4)
2802 : { /* S4 is very probable: test it first */
2803 77 : pari_sp av = avma;
2804 : struct galois_lift gl;
2805 77 : if (DEBUGLEVEL >= 4) err_printf("GaloisConj: Testing S4 first\n");
2806 77 : initlift(T, den, ga->p4, L, makeLden(L,den,gb), gb, &gl);
2807 77 : PG = s4galoisgen(&gl);
2808 77 : if (PG) return gerepileupto(ltop, PG);
2809 0 : set_avma(av);
2810 : }
2811 29715 : if (n == 36 && ga->ord==3 && ga->p4)
2812 : { /* F36 is very probable: test it first */
2813 14 : pari_sp av = avma;
2814 : struct galois_lift gl;
2815 14 : if (DEBUGLEVEL >= 4) err_printf("GaloisConj: Testing 3x3:4 first (p=%ld)\n",ga->p4);
2816 14 : initlift(T, den, ga->p4, L, makeLden(L,den,gb), gb, &gl);
2817 14 : PG = f36galoisgen(&gl);
2818 14 : if (PG) return gerepileupto(ltop, PG);
2819 0 : set_avma(av);
2820 : }
2821 29701 : frob = galoisfindfrobenius(T, L, den, bad, &gf, gb, ga);
2822 29699 : if (!frob) return gc_NULL(ltop);
2823 29685 : po = psi_order(gf.psi, gf.deg);
2824 29685 : if (!(ga->group&ga_easy) && po < (ulong) gf.deg && gf.deg/radicalu(gf.deg)%po == 0)
2825 : {
2826 0 : is_central = 1;
2827 0 : if (!bad) bad = gb->dis;
2828 0 : if (po > 1)
2829 : {
2830 0 : ofrob = frob; ogf = gf;
2831 0 : frob = perm_powu(frob, po);
2832 0 : gf.deg /= po;
2833 : }
2834 29685 : } else is_central = 0;
2835 29685 : sigma = permtopol(frob, L, M, den, gb->ladicabs, shifti(gb->ladicabs,-1), x);
2836 29686 : if (is_central && gf.fp != gf.deg)
2837 0 : { gf.p = findp(bad, T, sigma, gf.deg, &gf.Tmod); gf.fp = gf.deg;
2838 0 : gf.psi = const_vecsmall(lg(gf.Tmod)-1, 1);
2839 : }
2840 29686 : if (gf.deg == n) /* cyclic */
2841 : {
2842 24296 : GEN Tp = ZX_to_Flx(gel(gf.Tmod,1), gf.p);
2843 24296 : res = mkvec5(mkvec(frob), mkvec(Tp), mkvec(sigma), mkvecsmall(gf.p), cyclic_pc(n));
2844 24296 : return gerepilecopy(ltop, res);
2845 : }
2846 5390 : O = perm_cycles(frob);
2847 5390 : if (DEBUGLEVEL >= 9) err_printf("GaloisConj: Frobenius:%Ps\n", sigma);
2848 5390 : PG = galoisgenfixedfield0(O, L, sigma, T, is_central ? bad: NULL,
2849 : is_central ? &V: NULL, &gf, gb);
2850 5390 : if (PG == NULL) return gc_NULL(ltop);
2851 5369 : if (is_central && lg(gel(PG,1))!=3)
2852 0 : {
2853 : struct galois_perm gp;
2854 0 : galoisperm_init(&gp, L, M, den, gb->ladicabs, shifti(gb->ladicabs,-1), varn(T));
2855 0 : res = galoisgenlift_nilp(PG, O, V, T, frob, sigma, gb, &gf, &gp);
2856 0 : galoisperm_free(&gp);
2857 : }
2858 : else
2859 : {
2860 5369 : if (is_central && po > 1)
2861 : { /* backtrack powering of frob */
2862 0 : frob = ofrob; gf = ogf;
2863 0 : O = perm_cycles(ofrob);
2864 0 : sigma = permtopol(ofrob, L, M, den, gb->ladicabs, shifti(gb->ladicabs,-1), x);
2865 0 : PG = galoisgenfixedfield0(O, L, sigma, T, NULL, NULL, &gf, gb);
2866 0 : if (PG == NULL) return gc_NULL(ltop);
2867 : }
2868 5369 : res = galoisgenlift(gg_get_std(gel(PG,1)), gel(PG,2), O, L, M, frob, gb, &gf);
2869 : }
2870 5369 : if (!res) return gc_NULL(ltop);
2871 4606 : return gerepilecopy(ltop, res);
2872 : }
2873 :
2874 : /* T = polcyclo(N) */
2875 : static GEN
2876 966 : conjcyclo(GEN T, long N)
2877 : {
2878 966 : pari_sp av = avma;
2879 966 : long i, k = 1, d = eulerphiu(N), v = varn(T);
2880 966 : GEN L = cgetg(d+1,t_COL);
2881 14546 : for (i=1; i<=N; i++)
2882 13580 : if (ugcd(i, N)==1)
2883 : {
2884 6356 : GEN s = pol_xn(i, v);
2885 6356 : if (i >= d) s = ZX_rem(s, T);
2886 6356 : gel(L,k++) = s;
2887 : }
2888 966 : return gerepileupto(av, gen_sort(L, (void*)&gcmp, &gen_cmp_RgX));
2889 : }
2890 :
2891 : static GEN
2892 1246 : aut_to_groupelts(GEN aut, GEN L, ulong p)
2893 : {
2894 1246 : pari_sp av = avma;
2895 1246 : long i, d = lg(aut)-1;
2896 1246 : GEN P = ZV_to_Flv(L, p);
2897 1246 : GEN N = FlxV_Flv_multieval(aut, P, p);
2898 1246 : GEN q = perm_inv(vecsmall_indexsort(P));
2899 1246 : GEN G = cgetg(d+1, t_VEC);
2900 35945 : for (i=1; i<=d; i++)
2901 34699 : gel(G,i) = perm_mul(vecsmall_indexsort(gel(N,i)), q);
2902 1246 : return gerepilecopy(av, vecvecsmall_sort_shallow(G));
2903 : }
2904 :
2905 : static ulong
2906 7 : galois_find_totally_split(GEN P, GEN Q)
2907 : {
2908 7 : pari_sp av = avma;
2909 : forprime_t iter;
2910 : ulong p;
2911 7 : long n = degpol(P);
2912 7 : u_forprime_init(&iter, n*maxss(expu(n)-3, 2), ULONG_MAX);
2913 714 : while ((p = u_forprime_next(&iter)))
2914 : {
2915 714 : if (Flx_is_totally_split(ZX_to_Flx(P, p), p)
2916 7 : && (!Q || Flx_is_squarefree(ZX_to_Flx(Q, p), p)))
2917 7 : return gc_ulong(av, p);
2918 707 : set_avma(av);
2919 : }
2920 0 : return 0;
2921 : }
2922 :
2923 : GEN
2924 1253 : galoisinitfromaut(GEN T, GEN aut, ulong l)
2925 : {
2926 1253 : pari_sp ltop = avma;
2927 1253 : GEN nf, A, G, L, M, grp, den=NULL;
2928 : struct galois_analysis ga;
2929 : struct galois_borne gb;
2930 : long n;
2931 : pari_timer ti;
2932 :
2933 1253 : T = get_nfpol(T, &nf);
2934 1253 : n = degpol(T);
2935 1253 : if (nf)
2936 0 : { if (!den) den = nf_get_zkden(nf); }
2937 : else
2938 : {
2939 1253 : if (n <= 0) pari_err_IRREDPOL("galoisinit",T);
2940 1253 : RgX_check_ZX(T, "galoisinit");
2941 1253 : if (!ZX_is_squarefree(T))
2942 0 : pari_err_DOMAIN("galoisinit","issquarefree(pol)","=",gen_0,T);
2943 1253 : if (!gequal1(gel(T,n+2))) pari_err_IMPL("galoisinit(nonmonic)");
2944 : }
2945 1253 : if (lg(aut)-1 != n)
2946 7 : return gen_0;
2947 1246 : ga.l = l? l: galois_find_totally_split(T, NULL);
2948 1246 : if (!l) aut = RgXV_to_FlxV(aut, ga.l);
2949 1246 : gb.l = utoipos(ga.l);
2950 1246 : if (DEBUGLEVEL >= 1) timer_start(&ti);
2951 1246 : den = galoisborne(T, den, &gb, degpol(T));
2952 1246 : if (DEBUGLEVEL >= 1) timer_printf(&ti, "galoisborne()");
2953 1246 : L = ZpX_roots(T, gb.l, gb.valabs);
2954 1246 : if (DEBUGLEVEL >= 1) timer_printf(&ti, "ZpX_roots");
2955 1246 : M = FpV_invVandermonde(L, den, gb.ladicabs);
2956 1246 : if (DEBUGLEVEL >= 1) timer_printf(&ti, "FpV_invVandermonde()");
2957 1246 : A = aut_to_groupelts(aut, L, ga.l);
2958 1246 : G = groupelts_to_group(A);
2959 1246 : if (!G) G = trivialgroup();
2960 1239 : else A = group_elts(G,n);
2961 1246 : grp = cgetg(9, t_VEC);
2962 1246 : gel(grp,1) = T;
2963 1246 : gel(grp,2) = mkvec3(utoipos(ga.l), utoipos(gb.valabs), gb.ladicabs);
2964 1246 : gel(grp,3) = L;
2965 1246 : gel(grp,4) = M;
2966 1246 : gel(grp,5) = den;
2967 1246 : gel(grp,6) = A;
2968 1246 : gel(grp,7) = gel(G,1);
2969 1246 : gel(grp,8) = gel(G,2);
2970 1246 : return gerepilecopy(ltop, grp);
2971 : }
2972 :
2973 : GEN
2974 1239 : galoissplittinginit(GEN T, GEN D)
2975 : {
2976 1239 : pari_sp av = avma;
2977 1239 : GEN R = nfsplitting0(T, D, 3), P = gel(R,1), aut = gel(R,2);
2978 1232 : ulong p = itou(gel(R,3));
2979 1232 : return gerepileupto(av, galoisinitfromaut(P, aut, p));
2980 : }
2981 :
2982 : /* T: polynomial or nf, den multiple of common denominator of solutions or
2983 : * NULL (unknown). If T is nf, and den unknown, use den = denom(nf.zk) */
2984 : static GEN
2985 96714 : galoisconj4_main(GEN T, GEN den, long flag)
2986 : {
2987 96714 : pari_sp ltop = avma;
2988 : GEN nf, G, L, M, aut, grp;
2989 : struct galois_analysis ga;
2990 : struct galois_borne gb;
2991 : long n;
2992 : pari_timer ti;
2993 :
2994 96714 : T = get_nfpol(T, &nf);
2995 96715 : n = poliscyclo(T);
2996 96713 : if (n) return flag? galoiscyclo(n, varn(T)): conjcyclo(T, n);
2997 95292 : n = degpol(T);
2998 95291 : if (nf)
2999 53984 : { if (!den) den = nf_get_zkden(nf); }
3000 : else
3001 : {
3002 41307 : if (n <= 0) pari_err_IRREDPOL("galoisinit",T);
3003 41307 : RgX_check_ZX(T, "galoisinit");
3004 41308 : if (!ZX_is_squarefree(T))
3005 7 : pari_err_DOMAIN("galoisinit","issquarefree(pol)","=",gen_0,T);
3006 41306 : if (!gequal1(gel(T,n+2))) pari_err_IMPL("galoisinit(nonmonic)");
3007 : }
3008 95282 : if (n == 1)
3009 : {
3010 21 : if (!flag) { G = cgetg(2, t_COL); gel(G,1) = pol_x(varn(T)); return G;}
3011 21 : ga.l = 3;
3012 21 : ga.deg = 1;
3013 21 : den = gen_1;
3014 : }
3015 95261 : else if (!galoisanalysis(T, &ga, 1, NULL)) return gc_NULL(ltop);
3016 :
3017 28581 : if (den)
3018 : {
3019 17983 : if (typ(den) != t_INT) pari_err_TYPE("galoisconj4 [2nd arg integer]", den);
3020 17983 : den = absi_shallow(den);
3021 : }
3022 28581 : gb.l = utoipos(ga.l);
3023 28581 : if (DEBUGLEVEL >= 1) timer_start(&ti);
3024 28581 : den = galoisborne(T, den, &gb, degpol(T));
3025 28580 : if (DEBUGLEVEL >= 1) timer_printf(&ti, "galoisborne()");
3026 28580 : L = ZpX_roots(T, gb.l, gb.valabs);
3027 28579 : if (DEBUGLEVEL >= 1) timer_printf(&ti, "ZpX_roots");
3028 28579 : M = FpV_invVandermonde(L, den, gb.ladicabs);
3029 28581 : if (DEBUGLEVEL >= 1) timer_printf(&ti, "FpV_invVandermonde()");
3030 28581 : if (n == 1)
3031 : {
3032 21 : G = cgetg(3, t_VEC);
3033 21 : gel(G,1) = cgetg(1, t_VEC);
3034 21 : gel(G,2) = cgetg(1, t_VECSMALL);
3035 : }
3036 : else
3037 28560 : G = gg_get_std(galoisgen(T, L, M, den, NULL, &gb, &ga));
3038 28580 : if (DEBUGLEVEL >= 6) err_printf("GaloisConj: %Ps\n", G);
3039 28580 : if (!G) return gc_NULL(ltop);
3040 27789 : if (DEBUGLEVEL >= 1) timer_start(&ti);
3041 27789 : grp = cgetg(9, t_VEC);
3042 27789 : gel(grp,1) = T;
3043 27789 : gel(grp,2) = mkvec3(utoipos(ga.l), utoipos(gb.valabs), gb.ladicabs);
3044 27789 : gel(grp,3) = L;
3045 27789 : gel(grp,4) = M;
3046 27789 : gel(grp,5) = den;
3047 27789 : gel(grp,6) = group_elts(G,n);
3048 27790 : gel(grp,7) = gel(G,1);
3049 27790 : gel(grp,8) = gel(G,2);
3050 27790 : if (flag) return gerepilecopy(ltop, grp);
3051 8456 : aut = galoisvecpermtopol(grp, gal_get_group(grp), gb.ladicabs, shifti(gb.ladicabs,-1));
3052 8456 : settyp(aut, t_COL);
3053 8456 : if (DEBUGLEVEL >= 1) timer_printf(&ti, "Computation of polynomials");
3054 8456 : return gerepileupto(ltop, gen_sort(aut, (void*)&gcmp, &gen_cmp_RgX));
3055 : }
3056 :
3057 : /* Heuristic computation of #Aut(T), pinit = first prime to be tested */
3058 : long
3059 35847 : numberofconjugates(GEN T, long pinit)
3060 : {
3061 35847 : pari_sp av = avma;
3062 35847 : long c, nbtest, nbmax, n = degpol(T);
3063 : ulong p;
3064 : forprime_t S;
3065 :
3066 35847 : if (n == 1) return 1;
3067 35847 : nbmax = (n < 10)? 20: (n<<1) + 1;
3068 35847 : nbtest = 0;
3069 : #if 0
3070 : c = ZX_sturm(T); c = ugcd(c, n-c); /* too costly: finite primes are cheaper */
3071 : #else
3072 35847 : c = n;
3073 : #endif
3074 35847 : u_forprime_init(&S, pinit, ULONG_MAX);
3075 336987 : while((p = u_forprime_next(&S)))
3076 : {
3077 336976 : GEN L, Tp = ZX_to_Flx(T,p);
3078 : long i, nb;
3079 336945 : if (!Flx_is_squarefree(Tp, p)) continue;
3080 : /* unramified */
3081 279218 : nbtest++;
3082 279218 : L = Flx_nbfact_by_degree(Tp, &nb, p); /* L[i] = #factors of degree i */
3083 279239 : if (L[n/nb] == nb) {
3084 233316 : if (c == n && nbtest > 10) break; /* probably Galois */
3085 : }
3086 : else
3087 : {
3088 81713 : c = ugcd(c, L[1]);
3089 280587 : for (i = 2; i <= n; i++)
3090 223158 : if (L[i]) { c = ugcd(c, L[i]*i); if (c == 1) break; }
3091 81725 : if (c == 1) break;
3092 : }
3093 254899 : if (nbtest == nbmax) break;
3094 243405 : if (DEBUGLEVEL >= 6)
3095 0 : err_printf("NumberOfConjugates [%ld]:c=%ld,p=%ld\n", nbtest,c,p);
3096 243405 : set_avma(av);
3097 : }
3098 35846 : if (DEBUGLEVEL >= 2) err_printf("NumberOfConjugates:c=%ld,p=%ld\n", c, p);
3099 35846 : return gc_long(av,c);
3100 : }
3101 : static GEN
3102 0 : galoisconj4(GEN nf, GEN d)
3103 : {
3104 0 : pari_sp av = avma;
3105 : GEN G, T;
3106 0 : G = galoisconj4_main(nf, d, 0);
3107 0 : if (G) return G; /* Success */
3108 0 : set_avma(av); T = get_nfpol(nf, &nf);
3109 0 : G = cgetg(2, t_COL); gel(G,1) = pol_x(varn(T)); return G; /* Fail */
3110 :
3111 : }
3112 :
3113 : /* d multiplicative bound for the automorphism's denominators */
3114 : GEN
3115 69738 : galoisconj(GEN nf, GEN d)
3116 : {
3117 69738 : pari_sp av = avma;
3118 69738 : GEN G, NF, T = get_nfpol(nf,&NF);
3119 69738 : if (degpol(T) == 2)
3120 : { /* fast shortcut */
3121 24477 : GEN a = gel(T,4), b = gel(T,3);
3122 24477 : long v = varn(T);
3123 24477 : RgX_check_ZX(T, "galoisconj");
3124 24477 : if (!gequal1(a)) pari_err_IMPL("galoisconj(nonmonic)");
3125 24477 : G = cgetg(3, t_COL);
3126 24477 : gel(G,1) = deg1pol_shallow(gen_m1, negi(b), v);
3127 24476 : gel(G,2) = pol_x(v);
3128 24477 : return G;
3129 : }
3130 45261 : G = galoisconj4_main(nf, d, 0);
3131 45262 : if (G) return G; /* Success */
3132 35840 : set_avma(av); return galoisconj1(nf);
3133 : }
3134 :
3135 : /* FIXME: obsolete, use galoisconj(nf, d) directly */
3136 : GEN
3137 63 : galoisconj0(GEN nf, long flag, GEN d, long prec)
3138 : {
3139 : (void)prec;
3140 63 : switch(flag) {
3141 56 : case 2:
3142 56 : case 0: return galoisconj(nf, d);
3143 7 : case 1: return galoisconj1(nf);
3144 0 : case 4: return galoisconj4(nf, d);
3145 : }
3146 0 : pari_err_FLAG("nfgaloisconj");
3147 : return NULL; /*LCOV_EXCL_LINE*/
3148 : }
3149 :
3150 : /******************************************************************************/
3151 : /* Galois theory related algorithms */
3152 : /******************************************************************************/
3153 : GEN
3154 30330 : checkgal(GEN gal)
3155 : {
3156 30330 : if (typ(gal) == t_POL) pari_err_TYPE("checkgal [apply galoisinit first]",gal);
3157 30330 : if (typ(gal) != t_VEC || lg(gal) != 9) pari_err_TYPE("checkgal",gal);
3158 30323 : return gal;
3159 : }
3160 :
3161 : GEN
3162 51467 : galoisinit(GEN nf, GEN den)
3163 : {
3164 : GEN G;
3165 51467 : if (is_vec_t(typ(nf)) && lg(nf)==3 && is_vec_t(typ(gel(nf,2))))
3166 14 : return galoisinitfromaut(gel(nf,1), gel(nf,2), 0);
3167 51453 : G = galoisconj4_main(nf, den, 1);
3168 51442 : return G? G: gen_0;
3169 : }
3170 :
3171 : static GEN
3172 17835 : galoispermtopol_i(GEN gal, GEN perm, GEN mod, GEN mod2)
3173 : {
3174 17835 : switch (typ(perm))
3175 : {
3176 17597 : case t_VECSMALL:
3177 17597 : return permtopol(perm, gal_get_roots(gal), gal_get_invvdm(gal),
3178 : gal_get_den(gal), mod, mod2,
3179 17597 : varn(gal_get_pol(gal)));
3180 238 : case t_VEC: case t_COL: case t_MAT:
3181 238 : return galoisvecpermtopol(gal, perm, mod, mod2);
3182 : }
3183 0 : pari_err_TYPE("galoispermtopol", perm);
3184 : return NULL; /* LCOV_EXCL_LINE */
3185 : }
3186 :
3187 : GEN
3188 17835 : galoispermtopol(GEN gal, GEN perm)
3189 : {
3190 17835 : pari_sp av = avma;
3191 : GEN mod, mod2;
3192 17835 : gal = checkgal(gal);
3193 17835 : mod = gal_get_mod(gal);
3194 17835 : mod2 = shifti(mod,-1);
3195 17835 : return gerepilecopy(av, galoispermtopol_i(gal, perm, mod, mod2));
3196 : }
3197 :
3198 : GEN
3199 91 : galoiscosets(GEN O, GEN perm)
3200 : {
3201 91 : long i, j, k, u, f, l = lg(O);
3202 91 : GEN RC, C = cgetg(l,t_VECSMALL), o = gel(O,1);
3203 91 : pari_sp av = avma;
3204 91 : f = lg(o); u = o[1]; RC = zero_zv(lg(perm)-1);
3205 371 : for(i=1,j=1; j<l; i++)
3206 : {
3207 280 : GEN p = gel(perm,i);
3208 280 : if (RC[ p[u] ]) continue;
3209 763 : for(k=1; k<f; k++) RC[ p[ o[k] ] ] = 1;
3210 224 : C[j++] = i;
3211 : }
3212 91 : set_avma(av); return C;
3213 : }
3214 :
3215 : static GEN
3216 91 : fixedfieldfactor(GEN L, GEN O, GEN perm, GEN M, GEN den, GEN mod, GEN mod2,
3217 : long x,long y)
3218 : {
3219 91 : pari_sp ltop = avma;
3220 91 : long i, j, k, l = lg(O), lo = lg(gel(O,1));
3221 91 : GEN V, res, cosets = galoiscosets(O,perm), F = cgetg(lo+1,t_COL);
3222 :
3223 91 : gel(F, lo) = gen_1;
3224 91 : if (DEBUGLEVEL>=4) err_printf("GaloisFixedField:cosets=%Ps \n",cosets);
3225 91 : if (DEBUGLEVEL>=6) err_printf("GaloisFixedField:den=%Ps mod=%Ps \n",den,mod);
3226 91 : V = cgetg(l,t_COL); res = cgetg(l,t_VEC);
3227 315 : for (i = 1; i < l; i++)
3228 : {
3229 224 : pari_sp av = avma;
3230 224 : GEN G = cgetg(l,t_VEC), Lp = vecpermute(L, gel(perm, cosets[i]));
3231 938 : for (k = 1; k < l; k++)
3232 714 : gel(G,k) = FpV_roots_to_pol(vecpermute(Lp, gel(O,k)), mod, x);
3233 763 : for (j = 1; j < lo; j++)
3234 : {
3235 1834 : for(k = 1; k < l; k++) gel(V,k) = gmael(G,k,j+1);
3236 539 : gel(F,j) = vectopol(V, M, den, mod, mod2, y);
3237 : }
3238 224 : gel(res,i) = gerepileupto(av,gtopolyrev(F,x));
3239 : }
3240 91 : return gerepileupto(ltop,res);
3241 : }
3242 :
3243 : static void
3244 7406 : chk_perm(GEN perm, long n)
3245 : {
3246 7406 : if (typ(perm) != t_VECSMALL || lg(perm)!=n+1)
3247 0 : pari_err_TYPE("galoisfixedfield", perm);
3248 7406 : }
3249 :
3250 : static int
3251 12285 : is_group(GEN g)
3252 : {
3253 12285 : if (typ(g) == t_VEC && lg(g) == 3)
3254 : {
3255 1974 : GEN a = gel(g,1), o = gel(g,2);
3256 1974 : return typ(a)==t_VEC && typ(o)==t_VECSMALL && lg(a) == lg(o);
3257 : }
3258 10311 : return 0;
3259 : }
3260 :
3261 : GEN
3262 5747 : galoisfixedfield(GEN gal, GEN perm, long flag, long y)
3263 : {
3264 5747 : pari_sp ltop = avma;
3265 : GEN T, L, P, S, PL, O, res, mod, mod2, OL, sym;
3266 : long vT, n, i;
3267 5747 : if (flag<0 || flag>2) pari_err_FLAG("galoisfixedfield");
3268 5747 : gal = checkgal(gal); T = gal_get_pol(gal);
3269 5747 : vT = varn(T);
3270 5747 : L = gal_get_roots(gal); n = lg(L)-1;
3271 5747 : mod = gal_get_mod(gal);
3272 5747 : if (typ(perm) == t_VEC)
3273 : {
3274 4634 : if (is_group(perm)) perm = gel(perm, 1);
3275 10927 : for (i = 1; i < lg(perm); i++) chk_perm(gel(perm,i), n);
3276 4634 : O = vecperm_orbits(perm, n);
3277 : }
3278 : else
3279 : {
3280 1113 : chk_perm(perm, n);
3281 1113 : O = perm_cycles(perm);
3282 : }
3283 5747 : mod2 = shifti(mod,-1);
3284 5747 : OL = fixedfieldorbits(O, L);
3285 5747 : sym = fixedfieldsympol(OL, itou(gal_get_p(gal)));
3286 5747 : PL = sympol_eval(sym, OL, mod);
3287 5747 : P = FpX_center_i(FpV_roots_to_pol(PL, mod, vT), mod, mod2);
3288 5746 : if (flag==1) return gerepilecopy(ltop,P);
3289 1043 : S = fixedfieldinclusion(O, PL);
3290 1043 : S = vectopol(S, gal_get_invvdm(gal), gal_get_den(gal), mod, mod2, vT);
3291 1043 : if (flag==0)
3292 952 : res = cgetg(3, t_VEC);
3293 : else
3294 : {
3295 : GEN PM, Pden;
3296 : struct galois_borne Pgb;
3297 91 : long val = itos(gal_get_e(gal));
3298 91 : Pgb.l = gal_get_p(gal);
3299 91 : Pden = galoisborne(P, NULL, &Pgb, degpol(T)/degpol(P));
3300 91 : if (Pgb.valabs > val)
3301 : {
3302 7 : if (DEBUGLEVEL>=4)
3303 0 : err_printf("GaloisConj: increase p-adic prec by %ld.\n", Pgb.valabs-val);
3304 7 : PL = ZpX_liftroots(P, PL, Pgb.l, Pgb.valabs);
3305 7 : L = ZpX_liftroots(T, L, Pgb.l, Pgb.valabs);
3306 7 : mod = Pgb.ladicabs; mod2 = shifti(mod,-1);
3307 : }
3308 91 : PM = FpV_invVandermonde(PL, Pden, mod);
3309 91 : if (y < 0) y = 1;
3310 91 : if (varncmp(y, vT) <= 0)
3311 0 : pari_err_PRIORITY("galoisfixedfield", T, "<=", y);
3312 91 : setvarn(P, y);
3313 91 : res = cgetg(4, t_VEC);
3314 91 : gel(res,3) = fixedfieldfactor(L,O,gal_get_group(gal), PM,Pden,mod,mod2,vT,y);
3315 : }
3316 1043 : gel(res,1) = gcopy(P);
3317 1043 : gel(res,2) = gmodulo(S, T);
3318 1043 : return gerepileupto(ltop, res);
3319 : }
3320 :
3321 : /* gal a galois group output the underlying wss group */
3322 : GEN
3323 3556 : galois_group(GEN gal) { return mkvec2(gal_get_gen(gal), gal_get_orders(gal)); }
3324 :
3325 : GEN
3326 3038 : checkgroup(GEN g, GEN *S)
3327 : {
3328 3038 : if (is_group(g)) { *S = NULL; return g; }
3329 2681 : g = checkgal(g);
3330 2674 : *S = gal_get_group(g); return galois_group(g);
3331 : }
3332 :
3333 : GEN
3334 4627 : checkgroupelts(GEN G)
3335 : {
3336 : long i, n;
3337 4627 : if (typ(G)!=t_VEC) pari_err_TYPE("checkgroupelts", G);
3338 4613 : if (is_group(G))
3339 : { /* subgroup of S_n */
3340 371 : if (lg(gel(G,1))==1) return mkvec(mkvecsmall(1));
3341 371 : return group_elts(G, group_domain(G));
3342 : }
3343 4242 : if (lg(G)==9 && typ(gel(G,1))==t_POL)
3344 3913 : return gal_get_group(G); /* galoisinit */
3345 : /* vector of permutations ? */
3346 329 : n = lg(G)-1;
3347 329 : if (n==0) pari_err_DIM("checkgroupelts");
3348 5418 : for (i = 1; i <= n; i++)
3349 : {
3350 5131 : if (typ(gel(G,i)) != t_VECSMALL)
3351 28 : pari_err_TYPE("checkgroupelts (element)", gel(G,i));
3352 5103 : if (lg(gel(G,i)) != lg(gel(G,1)))
3353 14 : pari_err_DIM("checkgroupelts [length of permutations]");
3354 : }
3355 287 : return G;
3356 : }
3357 :
3358 : GEN
3359 224 : galoisisabelian(GEN gal, long flag)
3360 : {
3361 224 : pari_sp av = avma;
3362 224 : GEN S, G = checkgroup(gal,&S);
3363 224 : if (!group_isabelian(G)) { set_avma(av); return gen_0; }
3364 203 : switch(flag)
3365 : {
3366 49 : case 0: return gerepileupto(av, group_abelianHNF(G,S));
3367 49 : case 1: set_avma(av); return gen_1;
3368 105 : case 2: return gerepileupto(av, group_abelianSNF(G,S));
3369 0 : default: pari_err_FLAG("galoisisabelian");
3370 : }
3371 : return NULL; /* LCOV_EXCL_LINE */
3372 : }
3373 :
3374 : long
3375 56 : galoisisnormal(GEN gal, GEN sub)
3376 : {
3377 56 : pari_sp av = avma;
3378 56 : GEN S, G = checkgroup(gal, &S), H = checkgroup(sub, &S);
3379 56 : long res = group_subgroup_isnormal(G, H);
3380 56 : set_avma(av);
3381 56 : return res;
3382 : }
3383 :
3384 : static GEN
3385 308 : conjclasses_count(GEN conj, long nb)
3386 : {
3387 308 : long i, l = lg(conj);
3388 308 : GEN c = zero_zv(nb);
3389 4039 : for (i = 1; i < l; i++) c[conj[i]]++;
3390 308 : return c;
3391 : }
3392 : GEN
3393 308 : galoisconjclasses(GEN G)
3394 : {
3395 308 : pari_sp av = avma;
3396 308 : GEN c, e, cc = group_to_cc(G);
3397 308 : GEN elts = gel(cc,1), conj = gel(cc,2), repr = gel(cc,3);
3398 308 : long i, l = lg(conj), lc = lg(repr);
3399 308 : c = conjclasses_count(conj, lc-1);
3400 308 : e = cgetg(lc, t_VEC);
3401 3143 : for (i = 1; i < lc; i++) gel(e,i) = cgetg(c[i]+1, t_VEC);
3402 4039 : for (i = 1; i < l; i++)
3403 : {
3404 3731 : long ci = conj[i];
3405 3731 : gmael(e, ci, c[ci]) = gel(elts, i);
3406 3731 : c[ci]--;
3407 : }
3408 308 : return gerepilecopy(av, e);
3409 : }
3410 :
3411 : static GEN
3412 406 : groupelts_to_group_or_elts(GEN elts)
3413 : {
3414 406 : GEN G = groupelts_to_group(elts);
3415 406 : return G ? G: gcopy(elts);
3416 : }
3417 :
3418 : static GEN
3419 7 : vec_groupelts_to_group_or_elts(GEN x)
3420 413 : { pari_APPLY_same(groupelts_to_group_or_elts(gel(x,i))) }
3421 :
3422 : GEN
3423 1981 : galoissubgroups(GEN gal)
3424 : {
3425 1981 : pari_sp av = avma;
3426 1981 : GEN S, G = checkgroup(gal,&S);
3427 1981 : if (lg(gel(G,1))==1 && lg(S)>2)
3428 7 : return gerepileupto(av,
3429 : vec_groupelts_to_group_or_elts(groupelts_solvablesubgroups(S)));
3430 1974 : return gerepileupto(av, group_subgroups(G));
3431 : }
3432 :
3433 : GEN
3434 84 : galoissubfields(GEN G, long flag, long v)
3435 : {
3436 84 : pari_sp av = avma;
3437 84 : GEN L = galoissubgroups(G);
3438 84 : long i, l = lg(L);
3439 84 : GEN S = cgetg(l, t_VEC);
3440 1309 : for (i = 1; i < l; ++i) gel(S,i) = galoisfixedfield(G, gmael(L,i,1), flag, v);
3441 84 : return gerepileupto(av, S);
3442 : }
3443 :
3444 : GEN
3445 28 : galoisexport(GEN gal, long format)
3446 : {
3447 28 : pari_sp av = avma;
3448 28 : GEN S, G = checkgroup(gal,&S);
3449 28 : return gerepileupto(av, group_export(G,format));
3450 : }
3451 :
3452 : GEN
3453 497 : galoisidentify(GEN gal)
3454 : {
3455 497 : pari_sp av = avma;
3456 497 : GEN S, G = checkgroup(gal,&S);
3457 490 : long idx = group_ident(G,S), card = S ? lg(S)-1: group_order(G);
3458 490 : set_avma(av); return mkvec2s(card, idx);
3459 : }
3460 :
3461 : /* index of conjugacy class containing g */
3462 : static long
3463 36939 : cc_id(GEN cc, GEN g)
3464 : {
3465 36939 : GEN conj = gel(cc,2);
3466 36939 : long k = signe(gel(cc,4))? g[1]: vecvecsmall_search(gel(cc,1), g);
3467 36939 : return conj[k];
3468 : }
3469 :
3470 : static GEN
3471 4186 : Qevproj_RgX(GEN c, long d, GEN pro)
3472 4186 : { return RgV_to_RgX(Qevproj_down(RgX_to_RgC(c,d), pro), varn(c)); }
3473 : /* c in Z[X] / (X^o-1), To = polcyclo(o), T = polcyclo(expo), e = expo/o
3474 : * return c(X^e) mod T as an element of Z[X] / (To) */
3475 : static GEN
3476 3920 : chival(GEN c, GEN T, GEN To, long e, GEN pro, long phie)
3477 : {
3478 3920 : c = ZX_rem(c, To);
3479 3920 : if (e != 1) c = ZX_rem(RgX_inflate(c,e), T);
3480 3920 : if (pro) c = Qevproj_RgX(c, phie, pro);
3481 3920 : return c;
3482 : }
3483 : /* chi(g^l) = sum_{k=0}^{o-1} a_k zeta_o^{l*k} for all l;
3484 : * => a_k = 1/o sum_{l=0}^{o-1} chi(g^l) zeta_o^{-k*l}. Assume o > 1 */
3485 : static GEN
3486 861 : chiFT(GEN cp, GEN jg, GEN vze, long e, long o, ulong p, ulong pov2)
3487 : {
3488 861 : const long var = 1;
3489 861 : ulong oinv = Fl_inv(o,p);
3490 : long k, l;
3491 861 : GEN c = cgetg(o+2, t_POL);
3492 5642 : for (k = 0; k < o; k++)
3493 : {
3494 4781 : ulong a = 0;
3495 51478 : for (l=0; l<o; l++)
3496 : {
3497 46697 : ulong z = vze[Fl_mul(k,l,o)*e + 1];/* zeta_o^{-k*l} */
3498 46697 : a = Fl_add(a, Fl_mul(uel(cp,jg[l+1]), z, p), p);
3499 : }
3500 4781 : gel(c,k+2) = stoi(Fl_center(Fl_mul(a,oinv,p), p, pov2)); /* a_k */
3501 : }
3502 861 : c[1] = evalvarn(var) | evalsigne(1); return ZX_renormalize(c,o+2);
3503 : }
3504 : static GEN
3505 546 : cc_chartable(GEN cc)
3506 : {
3507 : GEN al, elts, rep, ctp, ct, dec, id, vjg, H, vord, operm;
3508 : long i, j, k, f, l, expo, lcl, n;
3509 : ulong p, pov2;
3510 :
3511 546 : elts = gel(cc,1); n = lg(elts)-1;
3512 546 : if (n == 1) return mkvec2(mkmat(mkcol(gen_1)), gen_1);
3513 532 : rep = gel(cc,3);
3514 532 : lcl = lg(rep);
3515 532 : vjg = cgetg(lcl, t_VEC);
3516 532 : vord = cgetg(lcl,t_VECSMALL);
3517 532 : id = identity_perm(lg(gel(elts,1))-1);
3518 532 : expo = 1;
3519 4879 : for(j=1;j<lcl;j++)
3520 : {
3521 4347 : GEN jg, h = id, g = gel(elts,rep[j]);
3522 : long o;
3523 4347 : vord[j] = o = perm_orderu(g);
3524 4347 : expo = ulcm(expo, o);
3525 4347 : gel(vjg,j) = jg = cgetg(o+1,t_VECSMALL);
3526 27671 : for (l=1; l<=o; l++)
3527 : {
3528 23324 : jg[l] = cc_id(cc, h); /* index of conjugacy class of g^(l-1) */
3529 23324 : if (l < o) h = perm_mul(h, g);
3530 : }
3531 : }
3532 : /* would sort conjugacy classes by inc. order */
3533 532 : operm = vecsmall_indexsort(vord);
3534 :
3535 : /* expo > 1, exponent of G */
3536 532 : p = unextprime(2*n+1);
3537 1043 : while (p%expo != 1) p = unextprime(p+1);
3538 : /* compute character table modulo p: idempotents of Z(KG) */
3539 532 : al = conjclasses_algcenter(cc, utoipos(p));
3540 532 : dec = algsimpledec_ss(al,1);
3541 532 : ctp = cgetg(lcl,t_VEC);
3542 4879 : for(i=1; i<lcl; i++)
3543 : {
3544 4347 : GEN e = ZV_to_Flv(gmael3(dec,i,3,1), p); /*(1/n)[(dim chi)chi(g): g in G]*/
3545 4347 : ulong d = usqrt(Fl_mul(e[1], n, p)); /* = chi(1) <= sqrt(n) < sqrt(p) */
3546 4347 : gel(ctp,i) = Flv_Fl_mul(e,Fl_div(n,d,p), p); /*[chi(g): g in G]*/
3547 : }
3548 : /* Find minimal f such that table is defined over Q(zeta(f)): the conductor
3549 : * of the class field Q(\zeta_e)^H defined by subgroup
3550 : * H = { k in (Z/e)^*: g^k ~ g, for all g } */
3551 532 : H = coprimes_zv(expo);
3552 3458 : for (k = 2; k < expo; k++)
3553 : {
3554 2926 : if (!H[k]) continue;
3555 2548 : for (j = 2; j < lcl; j++) /* skip g ~ 1 */
3556 2366 : if (umael(vjg,j,(k % vord[j])+1) != umael(vjg,j,2)) { H[k] = 0; break; }
3557 : }
3558 532 : f = znstar_conductor_bits(Flv_to_F2v(H));
3559 : /* lift character table to Z[zeta_f] */
3560 532 : pov2 = p>>1;
3561 532 : ct = cgetg(lcl, t_MAT);
3562 532 : if (f == 1)
3563 : { /* rational representation */
3564 938 : for (j=1; j<lcl; j++) gel(ct,j) = cgetg(lcl,t_COL);
3565 938 : for(j=1; j<lcl; j++)
3566 : {
3567 791 : GEN jg = gel(vjg,j); /* jg[l+1] = class of g^l */
3568 791 : long t = lg(jg) > 2? jg[2]: jg[1];
3569 6706 : for(i=1; i<lcl; i++)
3570 : {
3571 5915 : GEN cp = gel(ctp,i); /* cp[i] = chi(g_i) mod \P */
3572 5915 : gcoeff(ct,j,i) = stoi(Fl_center(cp[t], p, pov2));
3573 : }
3574 : }
3575 : }
3576 : else
3577 : {
3578 385 : const long var = 1;
3579 385 : ulong ze = Fl_powu(pgener_Fl(p),(p-1)/expo, p); /* seen as zeta_e^(-1) */
3580 385 : GEN vze = Fl_powers(ze, expo-1, p); /* vze[i] = ze^(i-1) */
3581 385 : GEN vzeZX = const_vec(p, gen_0);
3582 385 : GEN T = polcyclo(expo, var), vT = const_vec(expo,NULL), pro = NULL;
3583 385 : long phie = degpol(T), id1 = gel(vjg,1)[1]; /* index of 1_G, always 1 ? */
3584 385 : gel(vT, expo) = T;
3585 385 : if (f != expo)
3586 : {
3587 147 : long phif = eulerphiu(f);
3588 147 : GEN zf = ZX_rem(pol_xn(expo/f,var), T), zfj = zf;
3589 147 : GEN M = cgetg(phif+1, t_MAT);
3590 147 : gel(M,1) = col_ei(phie,1);
3591 518 : for (j = 2; j <= phif; j++)
3592 : {
3593 371 : gel(M,j) = RgX_to_RgC(zfj, phie);
3594 371 : if (j < phif) zfj = ZX_rem(ZX_mul(zfj, zf), T);
3595 : }
3596 147 : pro = Qevproj_init(M);
3597 : }
3598 385 : gel(vzeZX,1) = pol_1(var);
3599 3416 : for (i = 2; i <= expo; i++)
3600 : {
3601 3031 : GEN t = ZX_rem(pol_xn(expo-(i-1), var), T);
3602 3031 : if (pro) t = Qevproj_RgX(t, phie, pro);
3603 3031 : gel(vzeZX, vze[i]) = t;
3604 : }
3605 3941 : for(i=1; i<lcl; i++)
3606 : { /* loop over characters */
3607 3556 : GEN cp = gel(ctp,i), C, cj; /* cp[j] = chi(g_j) mod \P */
3608 3556 : long dim = cp[id1];
3609 3556 : gel(ct, i) = C = const_col(lcl-1, NULL);
3610 3556 : gel(C,operm[1]) = utoi(dim); /* chi(1_G) */
3611 40978 : for (j=lcl-1; j > 1; j--)
3612 : { /* loop over conjugacy classes, decreasing order: skip 1_G */
3613 37422 : long e, jperm = operm[j], o = vord[jperm];
3614 37422 : GEN To, jg = gel(vjg,jperm); /* jg[l+1] = class of g^l */
3615 :
3616 37422 : if (gel(C, jperm)) continue; /* done already */
3617 35903 : if (dim == 1) { gel(C, jperm) = gel(vzeZX, cp[jg[2]]); continue; }
3618 861 : e = expo / o;
3619 861 : cj = chiFT(cp, jg, vze, e, o, p, pov2);
3620 861 : To = gel(vT, o); if (!To) To = gel(vT,o) = polcyclo(o, var);
3621 861 : gel(C, jperm) = chival(cj, T, To, e, pro, phie);
3622 3920 : for (k = 2; k < o; k++)
3623 : {
3624 3059 : GEN ck = RgX_inflate(cj, k); /* chi(g^k) */
3625 3059 : gel(C, jg[k+1]) = chival(ck, T, To, e, pro, phie);
3626 : }
3627 : }
3628 : }
3629 : }
3630 532 : ct = gen_sort_shallow(ct,(void*)cmp_universal,cmp_nodata);
3631 1736 : i = 1; while (!vec_isconst(gel(ct,i))) i++;
3632 532 : if (i > 1) swap(gel(ct,1), gel(ct,i));
3633 532 : return mkvec2(ct, utoipos(f));
3634 : }
3635 : GEN
3636 553 : galoischartable(GEN gal)
3637 : {
3638 553 : pari_sp av = avma;
3639 553 : GEN cc = group_to_cc(gal);
3640 546 : return gerepilecopy(av, cc_chartable(cc));
3641 : }
3642 :
3643 : static void
3644 1491 : checkgaloischar(GEN ch, GEN repr)
3645 : {
3646 1491 : if (gvar(ch) == 0) pari_err_PRIORITY("galoischarpoly",ch,"=",0);
3647 1491 : if (!is_vec_t(typ(ch))) pari_err_TYPE("galoischarpoly", ch);
3648 1491 : if (lg(repr) != lg(ch)) pari_err_DIM("galoischarpoly");
3649 1491 : }
3650 :
3651 : static long
3652 1547 : galoischar_dim(GEN ch)
3653 : {
3654 1547 : pari_sp av = avma;
3655 1547 : long d = gtos(simplify_shallow(lift_shallow(gel(ch,1))));
3656 1547 : return gc_long(av,d);
3657 : }
3658 :
3659 : static GEN
3660 12355 : galoischar_aut_charpoly(GEN cc, GEN ch, GEN p, long d)
3661 : {
3662 12355 : GEN q = p, V = cgetg(d+2, t_POL);
3663 : long i;
3664 12355 : V[1] = evalsigne(1)|evalvarn(0);
3665 25970 : for (i = 1; i <= d; i++)
3666 : {
3667 13615 : gel(V,i+1) = gel(ch, cc_id(cc,q));
3668 13615 : if (i < d) q = perm_mul(q, p);
3669 : }
3670 12355 : return liftpol_shallow(RgXn_expint(RgX_neg(V),d+1));
3671 : }
3672 :
3673 : static GEN
3674 1491 : galoischar_charpoly(GEN cc, GEN ch, long o)
3675 : {
3676 1491 : GEN chm, V, elts = gel(cc,1), repr = gel(cc,3);
3677 1491 : long i, d, l = lg(ch), v = gvar(ch);
3678 1491 : checkgaloischar(ch, repr);
3679 1491 : chm = v < 0 ? ch: gmodulo(ch, polcyclo(o, v));
3680 1491 : V = cgetg(l, t_COL); d = galoischar_dim(ch);
3681 13846 : for (i = 1; i < l; i++)
3682 12355 : gel(V,i) = galoischar_aut_charpoly(cc, chm, gel(elts,repr[i]), d);
3683 1491 : return V;
3684 : }
3685 :
3686 : GEN
3687 1435 : galoischarpoly(GEN gal, GEN ch, long o)
3688 : {
3689 1435 : pari_sp av = avma;
3690 1435 : GEN cc = group_to_cc(gal);
3691 1435 : return gerepilecopy(av, galoischar_charpoly(cc, ch, o));
3692 : }
3693 :
3694 : static GEN
3695 56 : cc_char_det(GEN cc, GEN ch, long o)
3696 : {
3697 56 : long i, l = lg(ch), d = galoischar_dim(ch);
3698 56 : GEN V = galoischar_charpoly(cc, ch, o);
3699 280 : for (i = 1; i < l; i++) gel(V,i) = leading_coeff(gel(V,i));
3700 56 : return odd(d)? gneg(V): V;
3701 : }
3702 :
3703 : GEN
3704 56 : galoischardet(GEN gal, GEN ch, long o)
3705 : {
3706 56 : pari_sp av = avma;
3707 56 : GEN cc = group_to_cc(gal);
3708 56 : return gerepilecopy(av, cc_char_det(cc, ch, o));
3709 : }
|