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 : /*************************************************************************/
19 : /** **/
20 : /** Routines for handling VEC/COL **/
21 : /** **/
22 : /*************************************************************************/
23 : int
24 1841 : vec_isconst(GEN v)
25 : {
26 1841 : long i, l = lg(v);
27 : GEN w;
28 1841 : if (l==1) return 1;
29 1841 : w = gel(v,1);
30 6342 : for(i=2; i<l; i++)
31 5768 : if (!gequal(gel(v,i), w)) return 0;
32 574 : return 1;
33 : }
34 :
35 : int
36 17465 : vecsmall_isconst(GEN v)
37 : {
38 17465 : long i, l = lg(v);
39 : ulong w;
40 17465 : if (l==1) return 1;
41 17465 : w = uel(v,1);
42 30134 : for(i=2; i<l; i++)
43 24308 : if (uel(v,i) != w) return 0;
44 5826 : return 1;
45 : }
46 :
47 : /* Check if all the elements of v are different.
48 : * Use a quadratic algorithm. Could be done in n*log(n) by sorting. */
49 : int
50 0 : vec_is1to1(GEN v)
51 : {
52 0 : long i, j, l = lg(v);
53 0 : for (i=1; i<l; i++)
54 : {
55 0 : GEN w = gel(v,i);
56 0 : for(j=i+1; j<l; j++)
57 0 : if (gequal(gel(v,j), w)) return 0;
58 : }
59 0 : return 1;
60 : }
61 :
62 : GEN
63 98084 : vec_insert(GEN v, long n, GEN x)
64 : {
65 98084 : long i, l=lg(v);
66 98084 : GEN V = cgetg(l+1,t_VEC);
67 711340 : for(i=1; i<n; i++) gel(V,i) = gel(v,i);
68 98084 : gel(V,n) = x;
69 471646 : for(i=n+1; i<=l; i++) gel(V,i) = gel(v,i-1);
70 98084 : return V;
71 : }
72 : /*************************************************************************/
73 : /** **/
74 : /** Routines for handling VECSMALL **/
75 : /** **/
76 : /*************************************************************************/
77 : /* Sort v[0]...v[n-1] and put result in w[0]...w[n-1].
78 : * We accept v==w. w must be allocated. */
79 : static void
80 135247146 : vecsmall_sortspec(GEN v, long n, GEN w)
81 : {
82 135247146 : pari_sp ltop=avma;
83 135247146 : long nx=n>>1, ny=n-nx;
84 : long m, ix, iy;
85 : GEN x, y;
86 135247146 : if (n<=2)
87 : {
88 77865244 : if (n==1)
89 17067365 : w[0]=v[0];
90 60797879 : else if (n==2)
91 : {
92 62979872 : long v0=v[0], v1=v[1];
93 62979872 : if (v0<=v1) { w[0]=v0; w[1]=v1; }
94 3067736 : else { w[0]=v1; w[1]=v0; }
95 : }
96 77865244 : return;
97 : }
98 57381902 : x=new_chunk(nx); y=new_chunk(ny);
99 60778657 : vecsmall_sortspec(v,nx,x);
100 60622707 : vecsmall_sortspec(v+nx,ny,y);
101 270472463 : for (m=0, ix=0, iy=0; ix<nx && iy<ny; )
102 209043096 : if (x[ix]<=y[iy])
103 174939569 : w[m++]=x[ix++];
104 : else
105 34103527 : w[m++]=y[iy++];
106 63639715 : for(;ix<nx;) w[m++]=x[ix++];
107 228911185 : for(;iy<ny;) w[m++]=y[iy++];
108 61429367 : set_avma(ltop);
109 : }
110 :
111 : static long
112 25638083 : vecsmall_sort_max(GEN v)
113 : {
114 25638083 : long i, l = lg(v), max = -1;
115 91018830 : for (i = 1; i < l; i++)
116 89205467 : if (v[i] > max) { max = v[i]; if (max >= l) return -1; }
117 19262087 : else if (v[i] < 0) return -1;
118 1813363 : return max;
119 : }
120 : /* assume 0 <= v[i] <= M. In place. */
121 : void
122 1912906 : vecsmall_counting_sort(GEN v, long M)
123 : {
124 : pari_sp av;
125 : long i, j, k, l;
126 : GEN T;
127 1912906 : if (M == 0) return;
128 1912906 : av = avma; T = new_chunk(M + 1); l = lg(v);
129 9321451 : for (i = 0; i <= M; i++) T[i] = 0;
130 7466442 : for (i = 1; i < l; i++) T[v[i]]++; /* T[j] is # keys = j */
131 9321446 : for (j = 0, k = 1; j <= M; j++)
132 12962074 : for (i = 1; i <= T[j]; i++) v[k++] = j;
133 1912906 : set_avma(av);
134 : }
135 : /* not GC-clean, suitable for gerepileupto */
136 : GEN
137 16115 : vecsmall_counting_uniq(GEN v, long M)
138 : {
139 16115 : long i, k, l = lg(v);
140 : GEN T, U;
141 16115 : if (l == 1) return cgetg(1, t_VECSMALL);
142 16115 : if (M == 0) return mkvecsmall(0);
143 16115 : if (l == 2) return leafcopy(v);
144 15982 : U = new_chunk(M + 2);
145 15982 : T = U+1; /* allows to rewrite result over T also if T[0] = 1 */
146 124248 : for (i = 0; i <= M; i++) T[i] = 0;
147 201760 : for (i = 1; i < l; i++) T[v[i]] = 1;
148 124248 : for (i = 0, k = 1; i <= M; i++)
149 108266 : if (T[i]) U[k++] = i;
150 15982 : U[0] = evaltyp(t_VECSMALL) | _evallg(k); return U;
151 : }
152 : GEN
153 12432 : vecsmall_counting_indexsort(GEN v, long M)
154 : {
155 : pari_sp av;
156 12432 : long i, l = lg(v);
157 : GEN T, p;
158 12432 : if (M == 0 || l <= 2) return identity_zv(l - 1);
159 12418 : p = cgetg(l, t_VECSMALL); av = avma; T = new_chunk(M + 1);
160 55489 : for (i = 0; i <= M; i++) T[i] = 0;
161 14432025 : for (i = 1; i < l; i++) T[v[i]]++; /* T[j] is # keys = j */
162 43071 : for (i = 1; i <= M; i++) T[i] += T[i-1]; /* T[j] is # keys <= j */
163 14432025 : for (i = l-1; i > 0; i--) { p[T[v[i]]] = i; T[v[i]]--; }
164 12418 : return gc_const(av, p);
165 : }
166 :
167 : /* in place sort */
168 : void
169 30024283 : vecsmall_sort(GEN v)
170 : {
171 30024283 : long n = lg(v) - 1, max;
172 30024283 : if (n <= 1) return;
173 23318458 : if ((max = vecsmall_sort_max(v)) >= 0)
174 1912906 : vecsmall_counting_sort(v, max);
175 : else
176 21351967 : vecsmall_sortspec(v+1, n, v+1);
177 : }
178 :
179 : /* cf gen_sortspec */
180 : static GEN
181 7954213 : vecsmall_indexsortspec(GEN v, long n)
182 : {
183 : long nx, ny, m, ix, iy;
184 : GEN x, y, w;
185 7954213 : switch(n)
186 : {
187 152058 : case 1: return mkvecsmall(1);
188 3617621 : case 2: return (v[1] <= v[2])? mkvecsmall2(1,2): mkvecsmall2(2,1);
189 1406044 : case 3:
190 1406044 : if (v[1] <= v[2]) {
191 670654 : if (v[2] <= v[3]) return mkvecsmall3(1,2,3);
192 203056 : return (v[1] <= v[3])? mkvecsmall3(1,3,2)
193 613444 : : mkvecsmall3(3,1,2);
194 : } else {
195 735390 : if (v[1] <= v[3]) return mkvecsmall3(2,1,3);
196 276559 : return (v[2] <= v[3])? mkvecsmall3(2,3,1)
197 836323 : : mkvecsmall3(3,2,1);
198 : }
199 : }
200 2778490 : nx = n>>1; ny = n-nx;
201 2778490 : w = cgetg(n+1,t_VECSMALL);
202 2778495 : x = vecsmall_indexsortspec(v,nx);
203 2778495 : y = vecsmall_indexsortspec(v+nx,ny);
204 33633501 : for (m=1, ix=1, iy=1; ix<=nx && iy<=ny; )
205 30855005 : if (v[x[ix]] <= v[y[iy]+nx])
206 14947139 : w[m++] = x[ix++];
207 : else
208 15907866 : w[m++] = y[iy++]+nx;
209 5302152 : for(;ix<=nx;) w[m++] = x[ix++];
210 5308556 : for(;iy<=ny;) w[m++] = y[iy++]+nx;
211 2778496 : set_avma((pari_sp)w); return w;
212 : }
213 :
214 : /*indirect sort.*/
215 : GEN
216 2409723 : vecsmall_indexsort(GEN v)
217 : {
218 2409723 : long n = lg(v) - 1, max;
219 2409723 : if (n==0) return cgetg(1, t_VECSMALL);
220 2409660 : if ((max = vecsmall_sort_max(v)) >= 0)
221 12432 : return vecsmall_counting_indexsort(v, max);
222 : else
223 2397228 : return vecsmall_indexsortspec(v,n);
224 : }
225 :
226 : /* assume V sorted */
227 : GEN
228 31390 : vecsmall_uniq_sorted(GEN v)
229 : {
230 : long i, j, l;
231 31390 : GEN w = cgetg_copy(v, &l);
232 31390 : if (l == 1) return w;
233 31336 : w[1] = v[1];
234 34309 : for(i = j = 2; i < l; i++)
235 2973 : if (v[i] != w[j-1]) w[j++] = v[i];
236 31336 : stackdummy((pari_sp)(w + l), (pari_sp)(w + j));
237 31336 : setlg(w, j); return w;
238 : }
239 :
240 : GEN
241 16474 : vecsmall_uniq(GEN v)
242 : {
243 16474 : pari_sp av = avma;
244 : long max;
245 16474 : if ((max = vecsmall_sort_max(v)) >= 0)
246 16115 : v = vecsmall_counting_uniq(v, max);
247 : else
248 359 : { v = zv_copy(v); vecsmall_sort(v); v = vecsmall_uniq_sorted(v); }
249 16474 : return gerepileuptoleaf(av, v);
250 : }
251 :
252 : /* assume x sorted */
253 : long
254 0 : vecsmall_duplicate_sorted(GEN x)
255 : {
256 0 : long i,k,l=lg(x);
257 0 : if (l==1) return 0;
258 0 : for (k=x[1],i=2; i<l; k=x[i++])
259 0 : if (x[i] == k) return i;
260 0 : return 0;
261 : }
262 :
263 : long
264 20179 : vecsmall_duplicate(GEN x)
265 : {
266 20179 : pari_sp av=avma;
267 20179 : GEN p=vecsmall_indexsort(x);
268 20179 : long k,i,r=0,l=lg(x);
269 20179 : if (l==1) return 0;
270 28120 : for (k=x[p[1]],i=2; i<l; k=x[p[i++]])
271 7941 : if (x[p[i]] == k) { r=p[i]; break; }
272 20179 : set_avma(av);
273 20179 : return r;
274 : }
275 :
276 : static int
277 54275 : vecsmall_is1to1spec(GEN v, long n, GEN w)
278 : {
279 54275 : pari_sp ltop=avma;
280 54275 : long nx=n>>1, ny=n-nx;
281 : long m, ix, iy;
282 : GEN x, y;
283 54275 : if (n<=2)
284 : {
285 32660 : if (n==1)
286 11185 : w[0]=v[0];
287 21475 : else if (n==2)
288 : {
289 21475 : long v0=v[0], v1=v[1];
290 21475 : if (v0==v1) return 0;
291 21447 : else if (v0<v1) { w[0]=v0; w[1]=v1; }
292 4643 : else { w[0]=v1; w[1]=v0; }
293 : }
294 32632 : return 1;
295 : }
296 21615 : x = new_chunk(nx);
297 21615 : if (!vecsmall_is1to1spec(v,nx,x)) return 0;
298 21517 : y = new_chunk(ny);
299 21517 : if (!vecsmall_is1to1spec(v+nx,ny,y)) return 0;
300 84194 : for (m=0, ix=0, iy=0; ix<nx && iy<ny; )
301 62775 : if (x[ix]==y[iy]) return 0;
302 62726 : else if (x[ix]<y[iy])
303 38053 : w[m++]=x[ix++];
304 : else
305 24673 : w[m++]=y[iy++];
306 23524 : for(;ix<nx;) w[m++]=x[ix++];
307 53024 : for(;iy<ny;) w[m++]=y[iy++];
308 21419 : set_avma(ltop);
309 21419 : return 1;
310 : }
311 :
312 : int
313 11241 : vecsmall_is1to1(GEN V)
314 : {
315 11241 : pari_sp av = avma;
316 : long l;
317 11241 : GEN W = cgetg_copy(V, &l);
318 11241 : if (l <= 2) return 1;
319 11143 : return gc_bool(av, vecsmall_is1to1spec(V+1,l,W+1));
320 : }
321 :
322 : /*************************************************************************/
323 : /** **/
324 : /** Routines for handling vectors of VECSMALL **/
325 : /** **/
326 : /*************************************************************************/
327 :
328 : GEN
329 14 : vecvecsmall_sort(GEN x)
330 14 : { return gen_sort(x, (void*)&vecsmall_lexcmp, cmp_nodata); }
331 : GEN
332 365848 : vecvecsmall_sort_shallow(GEN x)
333 365848 : { return gen_sort_shallow(x, (void*)&vecsmall_lexcmp, cmp_nodata); }
334 :
335 : void
336 126 : vecvecsmall_sort_inplace(GEN x, GEN *perm)
337 126 : { gen_sort_inplace(x, (void*)&vecsmall_lexcmp, cmp_nodata, perm); }
338 :
339 : GEN
340 462 : vecvecsmall_sort_uniq(GEN x)
341 462 : { return gen_sort_uniq(x, (void*)&vecsmall_lexcmp, cmp_nodata); }
342 :
343 : GEN
344 861 : vecvecsmall_indexsort(GEN x)
345 861 : { return gen_indexsort(x, (void*)&vecsmall_lexcmp, cmp_nodata); }
346 :
347 : long
348 22934359 : vecvecsmall_search(GEN x, GEN y)
349 22934359 : { return gen_search(x,y,(void*)vecsmall_prefixcmp, cmp_nodata); }
350 :
351 : /* assume x non empty */
352 : long
353 0 : vecvecsmall_max(GEN x)
354 : {
355 0 : long i, l = lg(x), m = vecsmall_max(gel(x,1));
356 0 : for (i = 2; i < l; i++)
357 : {
358 0 : long t = vecsmall_max(gel(x,i));
359 0 : if (t > m) m = t;
360 : }
361 0 : return m;
362 : }
363 :
364 : /*************************************************************************/
365 : /** **/
366 : /** Routines for handling permutations **/
367 : /** **/
368 : /*************************************************************************/
369 :
370 : /* Permutations may be given by
371 : * perm (VECSMALL): a bijection from 1...n to 1...n i-->perm[i]
372 : * cyc (VEC of VECSMALL): a product of disjoint cycles. */
373 :
374 : /* Multiply (compose) two permutations, putting the result in the second one. */
375 : static void
376 21 : perm_mul_inplace2(GEN s, GEN t)
377 : {
378 21 : long i, l = lg(s);
379 525 : for (i = 1; i < l; i++) t[i] = s[t[i]];
380 21 : }
381 :
382 : GEN
383 0 : vecperm_extendschreier(GEN C, GEN v, long n)
384 : {
385 0 : pari_sp av = avma;
386 0 : long mj, lv = lg(v), m = 1, mtested = 1;
387 0 : GEN bit = const_vecsmall(n, 0);
388 0 : GEN cy = cgetg(n+1, t_VECSMALL);
389 0 : GEN sh = const_vec(n, gen_0);
390 0 : for(mj=1; mj<=n; mj++)
391 : {
392 0 : if (isintzero(gel(C,mj))) continue;
393 0 : gel(sh,mj) = gcopy(gel(C,mj));
394 0 : if (bit[mj]) continue;
395 0 : cy[m++] = mj;
396 0 : bit[mj] = 1;
397 : for(;;)
398 0 : {
399 0 : long o, mold = m;
400 0 : for (o = 1; o < lv; o++)
401 : {
402 0 : GEN vo = gel(v,o);
403 : long p;
404 0 : for (p = mtested; p < mold; p++) /* m increases! */
405 : {
406 0 : long j = vo[ cy[p] ];
407 0 : if (!bit[j])
408 : {
409 0 : gel(sh,j) = perm_mul(vo, gel(sh, cy[p]));
410 0 : cy[m++] = j;
411 : }
412 0 : bit[j] = 1;
413 : }
414 : }
415 0 : mtested = mold;
416 0 : if (m == mold) break;
417 : }
418 : }
419 0 : return gerepileupto(av, sh);
420 : }
421 :
422 : /* Orbits of the subgroup generated by v on {1,..,n} */
423 : static GEN
424 1410961 : vecperm_orbits_i(GEN v, long n)
425 : {
426 1410961 : long mj = 1, lv = lg(v), k, l;
427 1410961 : GEN cycle = cgetg(n+1, t_VEC), bit = const_vecsmall(n, 0);
428 8780287 : for (k = 1, l = 1; k <= n;)
429 : {
430 7369298 : pari_sp ltop = avma;
431 7369298 : long m = 1;
432 7369298 : GEN cy = cgetg(n+1, t_VECSMALL);
433 8762626 : for ( ; bit[mj]; mj++) /*empty*/;
434 7369240 : k++; cy[m++] = mj;
435 7369240 : bit[mj++] = 1;
436 : for(;;)
437 2405057 : {
438 9774297 : long o, mold = m;
439 19562901 : for (o = 1; o < lv; o++)
440 : {
441 9788604 : GEN vo = gel(v,o);
442 : long p;
443 30731909 : for (p = 1; p < m; p++) /* m increases! */
444 : {
445 20943305 : long j = vo[ cy[p] ];
446 20943305 : if (!bit[j]) cy[m++] = j;
447 20943305 : bit[j] = 1;
448 : }
449 : }
450 9774297 : if (m == mold) break;
451 2405057 : k += m - mold;
452 : }
453 7369240 : setlg(cy, m);
454 7369229 : gel(cycle,l++) = gerepileuptoleaf(ltop, cy);
455 : }
456 1410989 : setlg(cycle, l); return cycle;
457 : }
458 : /* memory clean version */
459 : GEN
460 4781 : vecperm_orbits(GEN v, long n)
461 : {
462 4781 : pari_sp av = avma;
463 4781 : return gerepilecopy(av, vecperm_orbits_i(v, n));
464 : }
465 :
466 : static int
467 2667 : isperm(GEN v)
468 : {
469 2667 : pari_sp av = avma;
470 2667 : long i, n = lg(v)-1;
471 : GEN w;
472 2667 : if (typ(v) != t_VECSMALL) return 0;
473 2667 : w = zero_zv(n);
474 26411 : for (i=1; i<=n; i++)
475 : {
476 23779 : long d = v[i];
477 23779 : if (d < 1 || d > n || w[d]) return gc_bool(av,0);
478 23744 : w[d] = 1;
479 : }
480 2632 : return gc_bool(av,1);
481 : }
482 :
483 : /* Compute the cyclic decomposition of a permutation */
484 : GEN
485 13312 : perm_cycles(GEN v)
486 : {
487 13312 : pari_sp av = avma;
488 13312 : return gerepilecopy(av, vecperm_orbits_i(mkvec(v), lg(v)-1));
489 : }
490 :
491 : GEN
492 259 : permcycles(GEN v)
493 : {
494 259 : if (!isperm(v)) pari_err_TYPE("permcycles",v);
495 252 : return perm_cycles(v);
496 : }
497 :
498 : /* Output the order of p */
499 : ulong
500 436352 : perm_orderu(GEN v)
501 : {
502 436352 : pari_sp av = avma;
503 436352 : GEN c = vecperm_orbits_i(mkvec(v), lg(v)-1);
504 : long i, d;
505 3169161 : for(i=1, d=1; i<lg(c); i++) d = ulcm(d, lg(gel(c,i))-1);
506 436362 : return gc_ulong(av,d);
507 : }
508 :
509 : static GEN
510 2002 : _domul(void *data, GEN x, GEN y)
511 : {
512 2002 : GEN (*mul)(GEN,GEN)=(GEN (*)(GEN,GEN)) data;
513 2002 : return mul(x,y);
514 : }
515 :
516 : /* Output the order of p */
517 : GEN
518 427 : perm_order(GEN v)
519 : {
520 427 : pari_sp av = avma;
521 427 : GEN c = vecperm_orbits_i(mkvec(v), lg(v)-1);
522 427 : long i, l = lg(c);
523 427 : GEN V = cgetg(l, t_VEC);
524 2856 : for (i = 1; i < l; i++)
525 2429 : gel(V,i) = utoi(lg(gel(c,i))-1);
526 427 : return gerepileuptoint(av, gen_product(V, (void *)lcmii, _domul));
527 : }
528 :
529 : GEN
530 434 : permorder(GEN v)
531 : {
532 434 : if (!isperm(v)) pari_err_TYPE("permorder",v);
533 427 : return perm_order(v);
534 : }
535 :
536 : /* sign of a permutation */
537 : long
538 956104 : perm_sign(GEN v)
539 : {
540 956104 : pari_sp av = avma;
541 956104 : GEN c = vecperm_orbits_i(mkvec(v), lg(v)-1);
542 956103 : long i, l = lg(c), s = 1;
543 5533854 : for (i = 1; i < l; i++)
544 4577751 : if (odd(lg(gel(c, i)))) s = -s;
545 956103 : return gc_long(av,s);
546 : }
547 :
548 : long
549 273 : permsign(GEN v)
550 : {
551 273 : if (!isperm(v)) pari_err_TYPE("permsign",v);
552 259 : return perm_sign(v);
553 : }
554 :
555 : GEN
556 5915 : Z_to_perm(long n, GEN x)
557 : {
558 : pari_sp av;
559 : ulong i, r;
560 5915 : GEN v = cgetg(n+1, t_VECSMALL);
561 5915 : if (n==0) return v;
562 5908 : uel(v,n) = 1; av = avma;
563 5908 : if (signe(x) <= 0) x = modii(x, mpfact(n));
564 27146 : for (r=n-1; r>=1; r--)
565 : {
566 : ulong a;
567 21238 : x = absdiviu_rem(x, n+1-r, &a);
568 71687 : for (i=r+1; i<=(ulong)n; i++)
569 50449 : if (uel(v,i) > a) uel(v,i)++;
570 21238 : uel(v,r) = a+1;
571 : }
572 5908 : set_avma(av); return v;
573 : }
574 : GEN
575 5915 : numtoperm(long n, GEN x)
576 : {
577 5915 : if (n < 0) pari_err_DOMAIN("numtoperm", "n", "<", gen_0, stoi(n));
578 5915 : if (typ(x) != t_INT) pari_err_TYPE("numtoperm",x);
579 5915 : return Z_to_perm(n, x);
580 : }
581 :
582 : /* destroys v */
583 : static GEN
584 1701 : perm_to_Z_inplace(GEN v)
585 : {
586 1701 : long l = lg(v), i, r;
587 1701 : GEN x = gen_0;
588 1701 : if (!isperm(v)) return NULL;
589 10143 : for (i = 1; i < l; i++)
590 : {
591 8449 : long vi = v[i];
592 8449 : if (vi <= 0) return NULL;
593 8449 : x = i==1 ? utoi(vi-1): addiu(muliu(x,l-i), vi-1);
594 25396 : for (r = i+1; r < l; r++)
595 16947 : if (v[r] > vi) v[r]--;
596 : }
597 1694 : return x;
598 : }
599 : GEN
600 1680 : perm_to_Z(GEN v)
601 : {
602 1680 : pari_sp av = avma;
603 1680 : GEN x = perm_to_Z_inplace(leafcopy(v));
604 1680 : if (!x) pari_err_TYPE("permtonum",v);
605 1680 : return gerepileuptoint(av, x);
606 : }
607 : GEN
608 1708 : permtonum(GEN p)
609 : {
610 1708 : pari_sp av = avma;
611 : GEN v, x;
612 1708 : switch(typ(p))
613 : {
614 1680 : case t_VECSMALL: return perm_to_Z(p);
615 21 : case t_VEC: case t_COL:
616 21 : if (RgV_is_ZV(p)) { v = ZV_to_zv(p); break; }
617 7 : default: pari_err_TYPE("permtonum",p);
618 : return NULL;/*LCOV_EXCL_LINE*/
619 : }
620 21 : x = perm_to_Z_inplace(v);
621 21 : if (!x) pari_err_TYPE("permtonum",p);
622 14 : return gerepileuptoint(av, x);
623 : }
624 :
625 : GEN
626 7300 : cyc_pow(GEN cyc, long exp)
627 : {
628 : long i, j, k, l, r;
629 : GEN c;
630 22167 : for (r = j = 1; j < lg(cyc); j++)
631 : {
632 14867 : long n = lg(gel(cyc,j)) - 1;
633 14867 : r += cgcd(n, exp);
634 : }
635 7300 : c = cgetg(r, t_VEC);
636 22167 : for (r = j = 1; j < lg(cyc); j++)
637 : {
638 14867 : GEN v = gel(cyc,j);
639 14867 : long n = lg(v) - 1, e = umodsu(exp,n), g = (long)ugcd(n, e), m = n / g;
640 31652 : for (i = 0; i < g; i++)
641 : {
642 16785 : GEN p = cgetg(m+1, t_VECSMALL);
643 16785 : gel(c,r++) = p;
644 54716 : for (k = 1, l = i; k <= m; k++)
645 : {
646 37931 : p[k] = v[l+1];
647 37931 : l += e; if (l >= n) l -= n;
648 : }
649 : }
650 : }
651 7300 : return c;
652 : }
653 :
654 : /* Compute the power of a permutation given by product of cycles
655 : * Ouput a perm, not a cyc */
656 : GEN
657 0 : cyc_pow_perm(GEN cyc, long exp)
658 : {
659 : long e, j, k, l, n;
660 : GEN p;
661 0 : for (n = 0, j = 1; j < lg(cyc); j++) n += lg(gel(cyc,j))-1;
662 0 : p = cgetg(n + 1, t_VECSMALL);
663 0 : for (j = 1; j < lg(cyc); j++)
664 : {
665 0 : GEN v = gel(cyc,j);
666 0 : n = lg(v) - 1; e = umodsu(exp, n);
667 0 : for (k = 1, l = e; k <= n; k++)
668 : {
669 0 : p[v[k]] = v[l+1];
670 0 : if (++l == n) l = 0;
671 : }
672 : }
673 0 : return p;
674 : }
675 :
676 : GEN
677 56 : perm_pow(GEN perm, GEN exp)
678 : {
679 56 : long i, r = lg(perm)-1;
680 56 : GEN p = zero_zv(r);
681 56 : pari_sp av = avma;
682 56 : GEN v = cgetg(r+1, t_VECSMALL);
683 196 : for (i=1; i<=r; i++)
684 : {
685 : long e, n, k, l;
686 140 : if (p[i]) continue;
687 56 : v[1] = i;
688 140 : for (n=1, k=perm[i]; k!=i; k=perm[k], n++) v[n+1] = k;
689 56 : e = umodiu(exp, n);
690 196 : for (k = 1, l = e; k <= n; k++)
691 : {
692 140 : p[v[k]] = v[l+1];
693 140 : if (++l == n) l = 0;
694 : }
695 : }
696 56 : set_avma(av); return p;
697 : }
698 :
699 : GEN
700 18690 : perm_powu(GEN perm, ulong exp)
701 : {
702 18690 : ulong i, r = lg(perm)-1;
703 18690 : GEN p = zero_zv(r);
704 18690 : pari_sp av = avma;
705 18690 : GEN v = cgetg(r+1, t_VECSMALL);
706 246540 : for (i=1; i<=r; i++)
707 : {
708 : ulong e, n, k, l;
709 227850 : if (p[i]) continue;
710 84707 : v[1] = i;
711 227850 : for (n=1, k=perm[i]; k!=i; k=perm[k], n++) v[n+1] = k;
712 84707 : e = exp % n;
713 312557 : for (k = 1, l = e; k <= n; k++)
714 : {
715 227850 : p[v[k]] = v[l+1];
716 227850 : if (++l == n) l = 0;
717 : }
718 : }
719 18690 : set_avma(av); return p;
720 : }
721 :
722 : GEN
723 21 : perm_to_GAP(GEN p)
724 : {
725 21 : pari_sp ltop=avma;
726 : GEN gap;
727 : GEN x;
728 : long i;
729 21 : long nb, c=0;
730 : char *s;
731 : long sz;
732 21 : long lp=lg(p)-1;
733 21 : if (typ(p) != t_VECSMALL) pari_err_TYPE("perm_to_GAP",p);
734 21 : x = perm_cycles(p);
735 21 : sz = (long) ((bfffo(lp)+1) * LOG10_2 + 1);
736 : /*Dry run*/
737 133 : for (i = 1, nb = 1; i < lg(x); ++i)
738 : {
739 112 : GEN z = gel(x,i);
740 112 : long lz = lg(z)-1;
741 112 : nb += 1+lz*(sz+2);
742 : }
743 21 : nb++;
744 : /*Real run*/
745 21 : gap = cgetg(nchar2nlong(nb) + 1, t_STR);
746 21 : s = GSTR(gap);
747 133 : for (i = 1; i < lg(x); ++i)
748 : {
749 : long j;
750 112 : GEN z = gel(x,i);
751 112 : if (lg(z) > 2)
752 : {
753 112 : s[c++] = '(';
754 364 : for (j = 1; j < lg(z); ++j)
755 : {
756 252 : if (j > 1)
757 : {
758 140 : s[c++] = ','; s[c++] = ' ';
759 : }
760 252 : sprintf(s+c,"%ld",z[j]);
761 567 : while(s[c++]) /* empty */;
762 252 : c--;
763 : }
764 112 : s[c++] = ')';
765 : }
766 : }
767 21 : if (!c) { s[c++]='('; s[c++]=')'; }
768 21 : s[c] = '\0';
769 21 : return gerepileupto(ltop,gap);
770 : }
771 :
772 : int
773 572495 : perm_commute(GEN s, GEN t)
774 : {
775 572495 : long i, l = lg(t);
776 40373487 : for (i = 1; i < l; i++)
777 39820382 : if (t[ s[i] ] != s[ t[i] ]) return 0;
778 553105 : return 1;
779 : }
780 :
781 : /*************************************************************************/
782 : /** **/
783 : /** Routines for handling groups **/
784 : /** **/
785 : /*************************************************************************/
786 : /* A Group is a t_VEC [gen,orders]
787 : * gen (vecvecsmall): list of generators given by permutations
788 : * orders (vecsmall): relatives orders of generators. */
789 942585 : INLINE GEN grp_get_gen(GEN G) { return gel(G,1); }
790 1597470 : INLINE GEN grp_get_ord(GEN G) { return gel(G,2); }
791 :
792 : /* A Quotient Group is a t_VEC [gen,coset]
793 : * gen (vecvecsmall): coset generators
794 : * coset (vecsmall): gen[coset[p[1]]] generate the p-coset.
795 : */
796 141820 : INLINE GEN quo_get_gen(GEN C) { return gel(C,1); }
797 30058 : INLINE GEN quo_get_coset(GEN C) { return gel(C,2); }
798 :
799 : static GEN
800 52458 : trivialsubgroups(void)
801 52458 : { GEN L = cgetg(2, t_VEC); gel(L,1) = trivialgroup(); return L; }
802 :
803 : /* Compute the order of p modulo the group given by a set */
804 : long
805 220220 : perm_relorder(GEN p, GEN set)
806 : {
807 220220 : pari_sp ltop = avma;
808 220220 : long n = 1, q = p[1];
809 654171 : while (!F2v_coeff(set,q)) { q = p[q]; n++; }
810 220220 : return gc_long(ltop,n);
811 : }
812 :
813 : GEN
814 13076 : perm_generate(GEN S, GEN H, long o)
815 : {
816 13076 : long i, n = lg(H)-1;
817 13076 : GEN L = cgetg(n*o + 1, t_VEC);
818 45885 : for(i=1; i<=n; i++) gel(L,i) = vecsmall_copy(gel(H,i));
819 50673 : for( ; i <= n*o; i++) gel(L,i) = perm_mul(gel(L,i-n), S);
820 13076 : return L;
821 : }
822 :
823 : /*Return the order (cardinality) of a group */
824 : long
825 711039 : group_order(GEN G)
826 : {
827 711039 : return zv_prod(grp_get_ord(G));
828 : }
829 :
830 : /* G being a subgroup of S_n, output n */
831 : long
832 26684 : group_domain(GEN G)
833 : {
834 26684 : GEN gen = grp_get_gen(G);
835 26684 : if (lg(gen) < 2) pari_err_DOMAIN("group_domain", "#G", "=", gen_1,G);
836 26684 : return lg(gel(gen,1)) - 1;
837 : }
838 :
839 : /*Left coset of g mod G: gG*/
840 : GEN
841 304430 : group_leftcoset(GEN G, GEN g)
842 : {
843 304430 : GEN gen = grp_get_gen(G), ord = grp_get_ord(G);
844 304430 : GEN res = cgetg(group_order(G)+1, t_VEC);
845 : long i, j, k;
846 304430 : gel(res,1) = vecsmall_copy(g);
847 304430 : k = 1;
848 560259 : for (i = 1; i < lg(gen); i++)
849 : {
850 255829 : long c = k * (ord[i] - 1);
851 704683 : for (j = 1; j <= c; j++) gel(res,++k) = perm_mul(gel(res,j), gel(gen,i));
852 : }
853 304430 : return res;
854 : }
855 : /*Right coset of g mod G: Gg*/
856 : GEN
857 182245 : group_rightcoset(GEN G, GEN g)
858 : {
859 182245 : GEN gen = grp_get_gen(G), ord = grp_get_ord(G);
860 182245 : GEN res = cgetg(group_order(G)+1, t_VEC);
861 : long i, j, k;
862 182245 : gel(res,1) = vecsmall_copy(g);
863 182245 : k = 1;
864 315553 : for (i = 1; i < lg(gen); i++)
865 : {
866 133308 : long c = k * (ord[i] - 1);
867 419517 : for (j = 1; j <= c; j++) gel(res,++k) = perm_mul(gel(gen,i), gel(res,j));
868 : }
869 182245 : return res;
870 : }
871 : /*Elements of a group from the generators, cf group_leftcoset*/
872 : GEN
873 140812 : group_elts(GEN G, long n)
874 : {
875 140812 : if (lg(G)==3 && typ(gel(G,1))==t_VEC)
876 : {
877 140812 : GEN gen = grp_get_gen(G), ord = grp_get_ord(G);
878 140812 : GEN res = cgetg(group_order(G)+1, t_VEC);
879 : long i, j, k;
880 140812 : gel(res,1) = identity_perm(n);
881 140812 : k = 1;
882 285397 : for (i = 1; i < lg(gen); i++)
883 : {
884 144585 : long c = k * (ord[i] - 1);
885 : /* j = 1, use res[1] = identity */
886 144585 : gel(res,++k) = vecsmall_copy(gel(gen,i));
887 384776 : for (j = 2; j <= c; j++) gel(res,++k) = perm_mul(gel(res,j), gel(gen,i));
888 : }
889 140812 : return res;
890 0 : } else return gcopy(G);
891 : }
892 :
893 : GEN
894 14448 : groupelts_conj_set(GEN elts, GEN p)
895 : {
896 14448 : long i, j, l = lg(elts), n = lg(p)-1;
897 14448 : GEN res = zero_F2v(n);
898 241465 : for(j = 1; j < n; j++)
899 241465 : if (p[j]==1) break;
900 101136 : for(i = 1; i < l; i++)
901 86688 : F2v_set(res, p[mael(elts,i,j)]);
902 14448 : return res;
903 : }
904 :
905 : GEN
906 28098 : groupelts_set(GEN elts, long n)
907 : {
908 28098 : GEN res = zero_F2v(n);
909 28098 : long i, l = lg(elts);
910 137606 : for(i=1; i<l; i++)
911 109508 : F2v_set(res,mael(elts,i,1));
912 28098 : return res;
913 : }
914 :
915 : /*Elements of a group from the generators, returned as a set (bitmap)*/
916 : GEN
917 90699 : group_set(GEN G, long n)
918 : {
919 90699 : GEN res = zero_F2v(n);
920 90699 : pari_sp av = avma;
921 90699 : GEN elts = group_elts(G, n);
922 90699 : long i, l = lg(elts);
923 284270 : for(i=1; i<l; i++)
924 193571 : F2v_set(res,mael(elts,i,1));
925 90699 : set_avma(av);
926 90699 : return res;
927 : }
928 :
929 : static int
930 17353 : sgcmp(GEN a, GEN b) { return vecsmall_lexcmp(gel(a,1),gel(b,1)); }
931 :
932 : GEN
933 497 : subgroups_tableset(GEN S, long n)
934 : {
935 497 : long i, l = lg(S);
936 497 : GEN v = cgetg(l, t_VEC);
937 5411 : for(i=1; i<l; i++)
938 4914 : gel(v,i) = mkvec2(group_set(gel(S,i), n), mkvecsmall(i));
939 497 : gen_sort_inplace(v,(void*)sgcmp,cmp_nodata, NULL);
940 497 : return v;
941 : }
942 :
943 : long
944 2002 : tableset_find_index(GEN tbl, GEN set)
945 : {
946 2002 : long i = tablesearch(tbl,mkvec2(set,mkvecsmall(0)),sgcmp);
947 2002 : if (!i) return 0;
948 2002 : return mael3(tbl,i,2,1);
949 : }
950 :
951 : GEN
952 52486 : trivialgroup(void) { retmkvec2(cgetg(1,t_VEC), cgetg(1,t_VECSMALL)); }
953 :
954 : /*Cyclic group generated by g of order s*/
955 : GEN
956 27916 : cyclicgroup(GEN g, long s)
957 27916 : { retmkvec2(mkvec( vecsmall_copy(g) ), mkvecsmall(s)); }
958 :
959 : /*Return the group generated by g1,g2 of relative orders s1,s2*/
960 : GEN
961 1085 : dicyclicgroup(GEN g1, GEN g2, long s1, long s2)
962 1085 : { retmkvec2( mkvec2(vecsmall_copy(g1), vecsmall_copy(g2)),
963 : mkvecsmall2(s1, s2) ); }
964 :
965 : /* return the quotient map G --> G/H */
966 : /*The ouput is [gen,hash]*/
967 : /* gen (vecvecsmall): coset generators
968 : * coset (vecsmall): vecsmall of coset number) */
969 : GEN
970 11725 : groupelts_quotient(GEN elt, GEN H)
971 : {
972 11725 : pari_sp ltop = avma;
973 : GEN p2, p3;
974 11725 : long i, j, a = 1;
975 11725 : long n = lg(gel(elt,1))-1, o = group_order(H);
976 : GEN el;
977 11725 : long le = lg(elt)-1;
978 11725 : GEN used = zero_F2v(le+1);
979 11725 : long l = le/o;
980 11725 : p2 = cgetg(l+1, t_VEC);
981 11725 : p3 = zero_zv(n);
982 11725 : el = zero_zv(n);
983 150073 : for (i = 1; i<=le; i++)
984 138348 : el[mael(elt,i,1)]=i;
985 68656 : for (i = 1; i <= l; ++i)
986 : {
987 : GEN V;
988 150619 : while(F2v_coeff(used,a)) a++;
989 56938 : V = group_leftcoset(H,gel(elt,a));
990 56938 : gel(p2,i) = gel(V,1);
991 195181 : for(j=1;j<lg(V);j++)
992 : {
993 138250 : long b = el[mael(V,j,1)];
994 138250 : if (b==0) pari_err_IMPL("group_quotient for a non-WSS group");
995 138243 : F2v_set(used,b);
996 : }
997 195167 : for (j = 1; j <= o; j++)
998 138236 : p3[mael(V, j, 1)] = i;
999 : }
1000 11718 : return gerepilecopy(ltop,mkvec2(p2,p3));
1001 : }
1002 :
1003 : GEN
1004 10213 : group_quotient(GEN G, GEN H)
1005 : {
1006 10213 : return groupelts_quotient(group_elts(G, group_domain(G)), H);
1007 : }
1008 :
1009 : /*Compute the image of a permutation by a quotient map.*/
1010 : GEN
1011 30058 : quotient_perm(GEN C, GEN p)
1012 : {
1013 30058 : GEN gen = quo_get_gen(C);
1014 30058 : GEN coset = quo_get_coset(C);
1015 30058 : long j, l = lg(gen);
1016 30058 : GEN p3 = cgetg(l, t_VECSMALL);
1017 283185 : for (j = 1; j < l; ++j)
1018 : {
1019 253127 : p3[j] = coset[p[mael(gen,j,1)]];
1020 253127 : if (p3[j]==0) pari_err_IMPL("quotient_perm for a non-WSS group");
1021 : }
1022 30058 : return p3;
1023 : }
1024 :
1025 : /* H is a subgroup of G, C is the quotient map G --> G/H
1026 : *
1027 : * Lift a subgroup S of G/H to a subgroup of G containing H */
1028 : GEN
1029 50778 : quotient_subgroup_lift(GEN C, GEN H, GEN S)
1030 : {
1031 50778 : GEN genH = grp_get_gen(H);
1032 50778 : GEN genS = grp_get_gen(S);
1033 50778 : GEN genC = quo_get_gen(C);
1034 50778 : long l1 = lg(genH)-1;
1035 50778 : long l2 = lg(genS)-1, j;
1036 50778 : GEN p1 = cgetg(3, t_VEC), L = cgetg(l1+l2+1, t_VEC);
1037 101724 : for (j = 1; j <= l1; ++j) gel(L,j) = gel(genH,j);
1038 118125 : for (j = 1; j <= l2; ++j) gel(L,l1+j) = gel(genC, mael(genS,j,1));
1039 50778 : gel(p1,1) = L;
1040 50778 : gel(p1,2) = vecsmall_concat(grp_get_ord(H), grp_get_ord(S));
1041 50778 : return p1;
1042 : }
1043 :
1044 : /* Let G a group and C a quotient map G --> G/H
1045 : * Assume H is normal, return the group G/H */
1046 : GEN
1047 10206 : quotient_group(GEN C, GEN G)
1048 : {
1049 10206 : pari_sp ltop = avma;
1050 : GEN Qgen, Qord, Qelt, Qset, Q;
1051 10206 : GEN Cgen = quo_get_gen(C);
1052 10206 : GEN Ggen = grp_get_gen(G);
1053 10206 : long i,j, n = lg(Cgen)-1, l = lg(Ggen);
1054 10206 : Qord = cgetg(l, t_VECSMALL);
1055 10206 : Qgen = cgetg(l, t_VEC);
1056 10206 : Qelt = mkvec(identity_perm(n));
1057 10206 : Qset = groupelts_set(Qelt, n);
1058 31164 : for (i = 1, j = 1; i < l; ++i)
1059 : {
1060 20958 : GEN g = quotient_perm(C, gel(Ggen,i));
1061 20958 : long o = perm_relorder(g, Qset);
1062 20958 : gel(Qgen,j) = g;
1063 20958 : Qord[j] = o;
1064 20958 : if (o != 1)
1065 : {
1066 13076 : Qelt = perm_generate(g, Qelt, o);
1067 13076 : Qset = groupelts_set(Qelt, n);
1068 13076 : j++;
1069 : }
1070 : }
1071 10206 : setlg(Qgen,j);
1072 10206 : setlg(Qord,j); Q = mkvec2(Qgen, Qord);
1073 10206 : return gerepilecopy(ltop,Q);
1074 : }
1075 :
1076 : GEN
1077 1512 : quotient_groupelts(GEN C)
1078 : {
1079 1512 : GEN G = quo_get_gen(C);
1080 1512 : long i, l = lg(G);
1081 1512 : GEN Q = cgetg(l, t_VEC);
1082 10612 : for (i = 1; i < l; ++i)
1083 9100 : gel(Q,i) = quotient_perm(C, gel(G,i));
1084 1512 : return Q;
1085 : }
1086 :
1087 : /* Return 1 if g normalizes N, 0 otherwise */
1088 : long
1089 182245 : group_perm_normalize(GEN N, GEN g)
1090 : {
1091 182245 : pari_sp ltop = avma;
1092 182245 : long r = gequal(vecvecsmall_sort_shallow(group_leftcoset(N, g)),
1093 : vecvecsmall_sort_shallow(group_rightcoset(N, g)));
1094 182245 : return gc_long(ltop, r);
1095 : }
1096 :
1097 : /* L is a list of subgroups, C is a coset and r a relative order.*/
1098 : static GEN
1099 65247 : liftlistsubgroups(GEN L, GEN C, long r)
1100 : {
1101 65247 : pari_sp ltop = avma;
1102 65247 : long c = lg(C)-1, l = lg(L)-1, n = lg(gel(C,1))-1, i, k;
1103 : GEN R;
1104 65247 : if (!l) return cgetg(1,t_VEC);
1105 58555 : R = cgetg(l*c+1, t_VEC);
1106 143024 : for (i = 1, k = 1; i <= l; ++i)
1107 : {
1108 84469 : GEN S = gel(L,i), Selt = group_set(S,n);
1109 84469 : GEN gen = grp_get_gen(S);
1110 84469 : GEN ord = grp_get_ord(S);
1111 : long j;
1112 279286 : for (j = 1; j <= c; ++j)
1113 : {
1114 194817 : GEN p = gel(C,j);
1115 194817 : if (perm_relorder(p, Selt) == r && group_perm_normalize(S, p))
1116 108731 : gel(R,k++) = mkvec2(vec_append(gen, p),
1117 : vecsmall_append(ord, r));
1118 : }
1119 : }
1120 58555 : setlg(R, k);
1121 58555 : return gerepilecopy(ltop, R);
1122 : }
1123 :
1124 : /* H is a normal subgroup, C is the quotient map G -->G/H,
1125 : * S is a subgroup of G/H, and G is embedded in Sym(l)
1126 : * Return all the subgroups K of G such that
1127 : * S= K mod H and K inter H={1} */
1128 : static GEN
1129 49266 : liftsubgroup(GEN C, GEN H, GEN S)
1130 : {
1131 49266 : pari_sp ltop = avma;
1132 49266 : GEN V = trivialsubgroups();
1133 49266 : GEN Sgen = grp_get_gen(S);
1134 49266 : GEN Sord = grp_get_ord(S);
1135 49266 : GEN Cgen = quo_get_gen(C);
1136 49266 : long n = lg(Sgen), i;
1137 114513 : for (i = 1; i < n; ++i)
1138 : { /*loop over generators of S*/
1139 65247 : GEN W = group_leftcoset(H, gel(Cgen, mael(Sgen, i, 1)));
1140 65247 : V = liftlistsubgroups(V, W, Sord[i]);
1141 : }
1142 49266 : return gerepilecopy(ltop,V);
1143 : }
1144 :
1145 : /* 1:A4, 2:S4, 3:F36, 0: other */
1146 : long
1147 10038 : group_isA4S4(GEN G)
1148 : {
1149 10038 : GEN elt = grp_get_gen(G);
1150 10038 : GEN ord = grp_get_ord(G);
1151 10038 : long n = lg(ord);
1152 10038 : if (n != 4 && n != 5) return 0;
1153 2219 : if (n==4 && ord[1]==3 && ord[2]==3 && ord[3]==4)
1154 : {
1155 : long i;
1156 7 : GEN p = gel(elt,1), q = gel(elt,2), r = gel(elt,3);
1157 259 : for(i=1; i<=36; i++)
1158 252 : if (p[r[i]]!=r[q[i]]) return 0;
1159 7 : return 3;
1160 : }
1161 2212 : if (ord[1]!=2 || ord[2]!=2 || ord[3]!=3) return 0;
1162 42 : if (perm_commute(gel(elt,1),gel(elt,3))) return 0;
1163 42 : if (n==4) return 1;
1164 21 : if (ord[4]!=2) return 0;
1165 21 : if (perm_commute(gel(elt,3),gel(elt,4))) return 0;
1166 21 : return 2;
1167 : }
1168 : /* compute all the subgroups of a group G */
1169 : GEN
1170 13230 : group_subgroups(GEN G)
1171 : {
1172 13230 : pari_sp ltop = avma;
1173 : GEN p1, H, C, Q, M, sg1, sg2, sg3;
1174 13230 : GEN gen = grp_get_gen(G);
1175 13230 : GEN ord = grp_get_ord(G);
1176 13230 : long lM, i, j, n = lg(gen);
1177 : long t;
1178 13230 : if (n == 1) return trivialsubgroups();
1179 10038 : t = group_isA4S4(G);
1180 10038 : if (t == 3)
1181 : {
1182 7 : GEN H = mkvec2(mkvec3(gel(gen,1), gel(gen,2), perm_sqr(gel(gen,3))),
1183 : mkvecsmall3(3, 3, 2));
1184 7 : GEN S = group_subgroups(H);
1185 7 : GEN V = cgetg(11,t_VEC);
1186 7 : gel(V,1) = cyclicgroup(gel(gen,3),4);
1187 63 : for (i=2; i<10; i++)
1188 56 : gel(V,i) = cyclicgroup(perm_mul(gmael3(V,i-1,1,1),gel(gen,i%3==1 ? 2:1)),4);
1189 7 : gel(V,10) = G;
1190 7 : return gerepilecopy(ltop,shallowconcat(S,V));
1191 : }
1192 10031 : else if (t)
1193 : {
1194 42 : GEN s = gel(gen,1); /*s = (1,2)(3,4) */
1195 42 : GEN t = gel(gen,2); /*t = (1,3)(2,4) */
1196 42 : GEN st = perm_mul(s, t); /*st = (1,4)(2,3) */
1197 42 : H = dicyclicgroup(s, t, 2, 2);
1198 : /* sg3 is the list of subgroups intersecting only partially with H*/
1199 42 : sg3 = cgetg((n==4)?4: 10, t_VEC);
1200 42 : gel(sg3,1) = cyclicgroup(s, 2);
1201 42 : gel(sg3,2) = cyclicgroup(t, 2);
1202 42 : gel(sg3,3) = cyclicgroup(st, 2);
1203 42 : if (n==5)
1204 : {
1205 21 : GEN u = gel(gen,3);
1206 21 : GEN v = gel(gen,4), w, u2;
1207 21 : if (zv_equal(perm_conj(u,s), t)) /*u=(2,3,4)*/
1208 21 : u2 = perm_sqr(u);
1209 : else
1210 : {
1211 0 : u2 = u;
1212 0 : u = perm_sqr(u);
1213 : }
1214 21 : if (perm_orderu(v)==2)
1215 : {
1216 21 : if (!perm_commute(s,v)) /*v=(1,2)*/
1217 : {
1218 0 : v = perm_conj(u,v);
1219 0 : if (!perm_commute(s,v)) v = perm_conj(u,v);
1220 : }
1221 21 : w = perm_mul(v,t); /*w=(1,4,2,3)*/
1222 : }
1223 : else
1224 : {
1225 0 : w = v;
1226 0 : if (!zv_equal(perm_sqr(w), s)) /*w=(1,4,2,3)*/
1227 : {
1228 0 : w = perm_conj(u,w);
1229 0 : if (!zv_equal(perm_sqr(w), s)) w = perm_conj(u,w);
1230 : }
1231 0 : v = perm_mul(w,t); /*v=(1,2)*/
1232 : }
1233 21 : gel(sg3,4) = dicyclicgroup(s,v,2,2);
1234 21 : gel(sg3,5) = dicyclicgroup(t,perm_conj(u,v),2,2);
1235 21 : gel(sg3,6) = dicyclicgroup(st,perm_conj(u2,v),2,2);
1236 21 : gel(sg3,7) = dicyclicgroup(s,w,2,2);
1237 21 : gel(sg3,8) = dicyclicgroup(t,perm_conj(u,w),2,2);
1238 21 : gel(sg3,9) = dicyclicgroup(st,perm_conj(u2,w),2,2);
1239 : }
1240 : }
1241 : else
1242 : {
1243 9989 : ulong osig = mael(factoru(ord[1]), 1, 1);
1244 9989 : GEN sig = perm_powu(gel(gen,1), ord[1]/osig);
1245 9989 : H = cyclicgroup(sig,osig);
1246 9989 : sg3 = NULL;
1247 : }
1248 10031 : C = group_quotient(G,H);
1249 10024 : Q = quotient_group(C,G);
1250 10024 : M = group_subgroups(Q); lM = lg(M);
1251 : /* sg1 is the list of subgroups containing H*/
1252 10017 : sg1 = cgetg(lM, t_VEC);
1253 59283 : for (i = 1; i < lM; ++i) gel(sg1,i) = quotient_subgroup_lift(C,H,gel(M,i));
1254 : /*sg2 is a list of lists of subgroups not intersecting with H*/
1255 10017 : sg2 = cgetg(lM, t_VEC);
1256 : /* Loop over all subgroups of G/H */
1257 59283 : for (j = 1; j < lM; ++j) gel(sg2,j) = liftsubgroup(C, H, gel(M,j));
1258 10017 : p1 = gconcat(sg1, shallowconcat1(sg2));
1259 10017 : if (sg3)
1260 : {
1261 42 : p1 = gconcat(p1, sg3);
1262 42 : if (n==5) /*ensure that the D4 subgroups of S4 are in supersolvable format*/
1263 84 : for(j = 3; j <= 5; j++)
1264 : {
1265 63 : GEN c = gmael(p1,j,1);
1266 63 : if (!perm_commute(gel(c,1),gel(c,3)))
1267 : {
1268 42 : if (perm_commute(gel(c,2),gel(c,3))) { swap(gel(c,1), gel(c,2)); }
1269 : else
1270 21 : perm_mul_inplace2(gel(c,2), gel(c,1));
1271 : }
1272 : }
1273 : }
1274 10017 : return gerepileupto(ltop,p1);
1275 : }
1276 :
1277 : /*return 1 if G is abelian, else 0*/
1278 : long
1279 8932 : group_isabelian(GEN G)
1280 : {
1281 8932 : GEN g = grp_get_gen(G);
1282 8932 : long i, j, n = lg(g);
1283 12852 : for(i=2; i<n; i++)
1284 13034 : for(j=1; j<i; j++)
1285 9114 : if (!perm_commute(gel(g,i), gel(g,j))) return 0;
1286 3990 : return 1;
1287 : }
1288 :
1289 : /*If G is abelian, return its HNF matrix*/
1290 : GEN
1291 385 : group_abelianHNF(GEN G, GEN S)
1292 : {
1293 385 : GEN M, g = grp_get_gen(G), o = grp_get_ord(G);
1294 385 : long i, j, k, n = lg(g);
1295 385 : if (!group_isabelian(G)) return NULL;
1296 315 : if (n==1) return cgetg(1,t_MAT);
1297 301 : if (!S) S = group_elts(G, group_domain(G));
1298 301 : M = cgetg(n,t_MAT);
1299 980 : for(i=1; i<n; i++)
1300 : {
1301 679 : GEN P, C = cgetg(n,t_COL);
1302 679 : pari_sp av = avma;
1303 679 : gel(M,i) = C;
1304 679 : P = perm_inv(perm_powu(gel(g,i), o[i]));
1305 959 : for(j=1; j<lg(S); j++)
1306 959 : if (zv_equal(P, gel(S,j))) break;
1307 679 : set_avma(av);
1308 679 : if (j==lg(S)) pari_err_BUG("galoisisabelian [inconsistent group]");
1309 679 : j--;
1310 1218 : for(k=1; k<i; k++)
1311 : {
1312 539 : long q = j / o[k];
1313 539 : gel(C,k) = stoi(j - q*o[k]);
1314 539 : j = q;
1315 : }
1316 679 : gel(C,k) = stoi(o[i]);
1317 1218 : for (k++; k<n; k++) gel(C,k) = gen_0;
1318 : }
1319 301 : return M;
1320 : }
1321 :
1322 : /*If G is abelian, return its abstract SNF matrix*/
1323 : GEN
1324 336 : group_abelianSNF(GEN G, GEN L)
1325 : {
1326 336 : pari_sp ltop = avma;
1327 336 : GEN H = group_abelianHNF(G,L);
1328 336 : if (!H) return NULL;
1329 266 : return gerepileupto(ltop, smithclean( ZM_snf(H) ));
1330 : }
1331 :
1332 : GEN
1333 434 : abelian_group(GEN v)
1334 : {
1335 434 : long card = zv_prod(v), i, d = 1, l = lg(v);
1336 434 : GEN G = cgetg(3,t_VEC), gen = cgetg(l,t_VEC);
1337 434 : gel(G,1) = gen;
1338 434 : gel(G,2) = vecsmall_copy(v);
1339 882 : for(i=1; i<l; i++)
1340 : {
1341 448 : GEN p = cgetg(card+1, t_VECSMALL);
1342 448 : long o = v[i], u = d*(o-1), j, k, l;
1343 448 : gel(gen, i) = p;
1344 : /* The following loop is over-optimized. Remember that I wrote it for
1345 : * testpermutation. Something has survived... BA */
1346 1036 : for(j=1;j<=card;)
1347 : {
1348 2296 : for(k=1;k<o;k++)
1349 4543 : for(l=1;l<=d; l++,j++) p[j] = j+d;
1350 1995 : for (l=1; l<=d; l++,j++) p[j] = j-u;
1351 : }
1352 448 : d += u;
1353 : }
1354 434 : return G;
1355 : }
1356 :
1357 : static long
1358 14609 : groupelts_subgroup_isnormal(GEN G, GEN H)
1359 : {
1360 14609 : long i, n = lg(G);
1361 64470 : for(i = 1; i < n; i++)
1362 62895 : if (!group_perm_normalize(H, gel(G,i))) return 0;
1363 1575 : return 1;
1364 : }
1365 :
1366 : /*return 1 if H is a normal subgroup of G*/
1367 : long
1368 336 : group_subgroup_isnormal(GEN G, GEN H)
1369 : {
1370 336 : if (lg(grp_get_gen(H)) > 1 && group_domain(G) != group_domain(H))
1371 0 : pari_err_DOMAIN("group_subgroup_isnormal","domain(H)","!=",
1372 : strtoGENstr("domain(G)"), H);
1373 336 : return groupelts_subgroup_isnormal(grp_get_gen(G), H);
1374 : }
1375 :
1376 : static GEN
1377 4816 : group_subgroup_kernel_set(GEN G, GEN H)
1378 : {
1379 : pari_sp av;
1380 4816 : GEN g = grp_get_gen(G);
1381 4816 : long i, n = lg(g);
1382 : GEN S, elts;
1383 4816 : long d = group_domain(G);
1384 4816 : if (lg(grp_get_gen(H)) > 1 && group_domain(G) != group_domain(H))
1385 0 : pari_err_DOMAIN("group_subgroup_isnormal","domain(H)","!=",
1386 : strtoGENstr("domain(G)"), H);
1387 4816 : elts = group_elts(H,d);
1388 4816 : S = groupelts_set(elts, d);
1389 4816 : av = avma;
1390 19264 : for(i=1; i<n; i++)
1391 : {
1392 14448 : F2v_and_inplace(S, groupelts_conj_set(elts,gel(g,i)));
1393 14448 : set_avma(av);
1394 : }
1395 4816 : return S;
1396 : }
1397 :
1398 : int
1399 4816 : group_subgroup_is_faithful(GEN G, GEN H)
1400 : {
1401 4816 : pari_sp av = avma;
1402 4816 : GEN K = group_subgroup_kernel_set(G,H);
1403 4816 : F2v_clear(K,1);
1404 4816 : return gc_long(av, F2v_equal0(K));
1405 : }
1406 :
1407 : long
1408 0 : groupelts_exponent(GEN elts)
1409 : {
1410 0 : long i, n = lg(elts)-1, expo = 1;
1411 0 : for(i=1; i<=n; i++) expo = ulcm(expo, perm_orderu(gel(elts,i)));
1412 0 : return expo;
1413 : }
1414 :
1415 : GEN
1416 700 : groupelts_center(GEN S)
1417 : {
1418 700 : pari_sp ltop = avma;
1419 700 : long i, j, n = lg(S)-1, l = n;
1420 700 : GEN V, elts = zero_F2v(n+1);
1421 25732 : for(i=1; i<=n; i++)
1422 : {
1423 25032 : if (F2v_coeff(elts,i)) { l--; continue; }
1424 573384 : for(j=1; j<=n; j++)
1425 563192 : if (!perm_commute(gel(S,i),gel(S,j)))
1426 : {
1427 14322 : F2v_set(elts,i);
1428 14322 : F2v_set(elts,j); l--; break;
1429 : }
1430 : }
1431 700 : V = cgetg(l+1,t_VEC);
1432 25732 : for (i=1, j=1; i<=n ;i++)
1433 25032 : if (!F2v_coeff(elts,i)) gel(V,j++) = vecsmall_copy(gel(S,i));
1434 700 : return gerepileupto(ltop,V);
1435 : }
1436 :
1437 : GEN
1438 4270 : groupelts_conjclasses(GEN elts, long *pnbcl)
1439 : {
1440 4270 : long i, j, cl = 0, n = lg(elts)-1;
1441 4270 : GEN c = const_vecsmall(n,0);
1442 4270 : pari_sp av = avma;
1443 52850 : for (i=1; i<=n; i++)
1444 : {
1445 48580 : GEN g = gel(elts,i);
1446 48580 : if (c[i]) continue;
1447 34965 : c[i] = ++cl;
1448 486871 : for(j=1; j<=n; j++)
1449 451906 : if (j != i)
1450 : {
1451 416941 : GEN h = perm_conj(gel(elts,j), g);
1452 416941 : long i2 = gen_search(elts,h,(void*)&vecsmall_lexcmp,&cmp_nodata);
1453 416941 : c[i2] = cl; set_avma(av);
1454 : }
1455 : }
1456 4270 : if (pnbcl) *pnbcl = cl;
1457 4270 : return c;
1458 : }
1459 :
1460 : GEN
1461 4270 : conjclasses_repr(GEN conj, long nb)
1462 : {
1463 4270 : long i, l = lg(conj);
1464 4270 : GEN e = const_vecsmall(nb, 0);
1465 52850 : for(i=1; i<l; i++)
1466 : {
1467 48580 : long ci = conj[i];
1468 48580 : if (!e[ci]) e[ci] = i;
1469 : }
1470 4270 : return e;
1471 : }
1472 :
1473 : /* elts of G sorted wrt vecsmall_lexcmp order: g in G is determined by g[1]
1474 : * so sort by increasing g[1] */
1475 : static GEN
1476 3885 : galois_elts_sorted(GEN gal)
1477 : {
1478 : long i, l;
1479 3885 : GEN elts = gal_get_group(gal), v = cgetg_copy(elts, &l);
1480 43141 : for (i = 1; i < l; i++) { GEN g = gel(elts,i); gel(v, g[1]) = g; }
1481 3885 : return v;
1482 : }
1483 : GEN
1484 4291 : group_to_cc(GEN G)
1485 : {
1486 4291 : GEN elts = checkgroupelts(G), z = cgetg(5,t_VEC);
1487 4270 : long n, flag = 1;
1488 4270 : if (typ(gel(G,1)) == t_POL)
1489 3885 : elts = galois_elts_sorted(G); /* galoisinit */
1490 : else
1491 : {
1492 385 : long i, l = lg(elts);
1493 385 : elts = gen_sort_shallow(elts,(void*)vecsmall_lexcmp,cmp_nodata);
1494 5824 : for (i = 1; i < l; i++)
1495 5586 : if (gel(elts,i)[1] != i) { flag = 0; break; }
1496 : }
1497 4270 : gel(z,1) = elts;
1498 4270 : gel(z,2) = groupelts_conjclasses(elts,&n);
1499 4270 : gel(z,3) = conjclasses_repr(gel(z,2),n);
1500 4270 : gel(z,4) = utoi(flag); return z;
1501 : }
1502 :
1503 : /* S a list of generators */
1504 : GEN
1505 0 : groupelts_abelian_group(GEN S)
1506 : {
1507 0 : pari_sp ltop = avma;
1508 : GEN Qgen, Qord, Qelt;
1509 0 : long i, j, n = lg(gel(S,1))-1, l = lg(S);
1510 0 : Qord = cgetg(l, t_VECSMALL);
1511 0 : Qgen = cgetg(l, t_VEC);
1512 0 : Qelt = mkvec(identity_perm(n));
1513 0 : for (i = 1, j = 1; i < l; ++i)
1514 : {
1515 0 : GEN g = gel(S,i);
1516 0 : long o = perm_relorder(g, groupelts_set(Qelt, n));
1517 0 : gel(Qgen,j) = g;
1518 0 : Qord[j] = o;
1519 0 : if (o != 1) { Qelt = perm_generate(g, Qelt, o); j++; }
1520 : }
1521 0 : setlg(Qgen,j);
1522 0 : setlg(Qord,j);
1523 0 : return gerepilecopy(ltop, mkvec2(Qgen, Qord));
1524 : }
1525 :
1526 : GEN
1527 14 : group_export_GAP(GEN G)
1528 : {
1529 14 : pari_sp av = avma;
1530 14 : GEN s, comma, g = grp_get_gen(G);
1531 14 : long i, k, l = lg(g);
1532 14 : if (l == 1) return strtoGENstr("Group(())");
1533 7 : s = cgetg(2*l, t_VEC);
1534 7 : comma = strtoGENstr(", ");
1535 7 : gel(s,1) = strtoGENstr("Group(");
1536 28 : for (i=1, k=2; i < l; ++i)
1537 : {
1538 21 : if (i > 1) gel(s,k++) = comma;
1539 21 : gel(s,k++) = perm_to_GAP(gel(g,i));
1540 : }
1541 7 : gel(s,k++) = strtoGENstr(")");
1542 7 : return gerepilecopy(av, shallowconcat1(s));
1543 : }
1544 :
1545 : GEN
1546 14 : group_export_MAGMA(GEN G)
1547 : {
1548 14 : pari_sp av = avma;
1549 14 : GEN s, comma, g = grp_get_gen(G);
1550 14 : long i, k, l = lg(g);
1551 14 : if (l == 1) return strtoGENstr("PermutationGroup<1|>");
1552 7 : s = cgetg(2*l, t_VEC);
1553 7 : comma = strtoGENstr(", ");
1554 7 : gel(s,1) = gsprintf("PermutationGroup<%ld|",group_domain(G));
1555 28 : for (i=1, k=2; i < l; ++i)
1556 : {
1557 21 : if (i > 1) gel(s,k++) = comma;
1558 21 : gel(s,k++) = GENtoGENstr( vecsmall_to_vec(gel(g,i)) );
1559 : }
1560 7 : gel(s,k++) = strtoGENstr(">");
1561 7 : return gerepilecopy(av, shallowconcat1(s));
1562 : }
1563 :
1564 : GEN
1565 28 : group_export(GEN G, long format)
1566 : {
1567 28 : switch(format)
1568 : {
1569 14 : case 0: return group_export_GAP(G);
1570 14 : case 1: return group_export_MAGMA(G);
1571 : }
1572 0 : pari_err_FLAG("galoisexport");
1573 0 : return NULL; /*-Wall*/
1574 : }
1575 :
1576 : static GEN
1577 3577 : groupelts_cyclic_subgroups(GEN G)
1578 : {
1579 3577 : pari_sp av = avma;
1580 3577 : long i, j, n = lg(G)-1;
1581 : GEN elts, f, gen, ord;
1582 3577 : if (n==1) return cgetg(1,t_VEC);
1583 3577 : elts = zero_F2v(lg(gel(G,1))-1);
1584 3577 : gen = cgetg(n+1, t_VECSMALL);
1585 3577 : ord = cgetg(n+1, t_VECSMALL);
1586 53109 : for (i=1, j=1; i<=n; i++)
1587 : {
1588 49532 : long k = 1, o, c = 0;
1589 49532 : GEN p = gel(G, i);
1590 49532 : if (F2v_coeff(elts, p[1])) continue;
1591 35903 : o = perm_orderu(p);
1592 35903 : gen[j] = i; ord[j] = o; j++;
1593 : do
1594 : {
1595 96229 : if (cgcd(o, ++c)==1) F2v_set(elts, p[k]);
1596 96229 : k = p[k];
1597 96229 : } while (k!=1);
1598 : }
1599 3577 : setlg(gen, j);
1600 3577 : setlg(ord, j);
1601 3577 : f = vecsmall_indexsort(ord);
1602 3577 : return gerepilecopy(av, mkvec2(vecsmallpermute(gen, f),
1603 : vecsmallpermute(ord, f)));
1604 : }
1605 :
1606 : GEN
1607 3584 : groupelts_to_group(GEN G)
1608 : {
1609 3584 : pari_sp av = avma;
1610 : GEN L, cyc, ord;
1611 3584 : long i, l, n = lg(G)-1;
1612 3584 : if (n==1) return trivialgroup();
1613 3563 : L = groupelts_cyclic_subgroups(G);
1614 3563 : cyc = gel(L,1); ord = gel(L,2);
1615 3563 : l = lg(cyc);
1616 16324 : for (i = l-1; i >= 2; i--)
1617 : {
1618 15645 : GEN p = gel(G,cyc[i]);
1619 15645 : long o = ord[i];
1620 15645 : GEN H = cyclicgroup(p, o);
1621 15645 : if (o == n) return gerepileupto(av, H);
1622 14273 : if (groupelts_subgroup_isnormal(G, H))
1623 : {
1624 1512 : GEN C = groupelts_quotient(G, H);
1625 1512 : GEN Q = quotient_groupelts(C);
1626 1512 : GEN R = groupelts_to_group(Q);
1627 1512 : if (!R) return gc_NULL(av);
1628 1512 : return gerepilecopy(av, quotient_subgroup_lift(C, H, R));
1629 : }
1630 : }
1631 679 : if (n==12 && l==9 && ord[2]==2 && ord[3]==2 && ord[5]==3)
1632 602 : return gerepilecopy(av,
1633 301 : mkvec2(mkvec3(gel(G,cyc[2]), gel(G,cyc[3]), gel(G,cyc[5])), mkvecsmall3(2,2,3)));
1634 378 : if (n==24 && l==18 && ord[11]==3 && ord[15]==4 && ord[16]==4)
1635 : {
1636 350 : GEN t21 = perm_sqr(gel(G,cyc[15]));
1637 350 : GEN t22 = perm_sqr(gel(G,cyc[16]));
1638 350 : GEN s = perm_mul(t22, gel(G,cyc[15]));
1639 700 : return gerepilecopy(av,
1640 350 : mkvec2(mkvec4(t21,t22, gel(G,cyc[11]), s), mkvecsmall4(2,2,3,2)));
1641 : }
1642 28 : if (n==36 && l==24 && ord[11]==3 && ord[15]==4)
1643 : {
1644 7 : GEN t1 = gel(G,cyc[11]), t3 = gel(G,cyc[15]);
1645 7 : return gerepilecopy(av,
1646 : mkvec2(mkvec3(perm_conj(t3, t1), t1, t3), mkvecsmall3(3,3,4)));
1647 : }
1648 21 : return gc_NULL(av);
1649 : }
1650 :
1651 : static GEN
1652 1176 : subg_get_gen(GEN subg) { return gel(subg, 1); }
1653 :
1654 : static GEN
1655 8694 : subg_get_set(GEN subg) { return gel(subg, 2); }
1656 :
1657 : static GEN
1658 812 : groupelt_subg_normalize(GEN elt, GEN subg, GEN cyc)
1659 : {
1660 812 : GEN gen = subg_get_gen(subg), set = subg_get_set(subg);
1661 812 : long i, j, u, n = lg(elt)-1, lgen = lg(gen);
1662 812 : GEN b = F2v_copy(cyc), res = zero_F2v(n);
1663 49532 : for(i = 1; i <= n; i++)
1664 : {
1665 : GEN g;
1666 48720 : if (!F2v_coeff(b, i)) continue;
1667 22386 : g = gel(elt,i);
1668 763532 : for(u=1; u<=n; u++)
1669 763532 : if (g[u]==1) break;
1670 25172 : for(j=1; j<lgen; j++)
1671 : {
1672 23786 : GEN h = gel(elt,gen[j]);
1673 23786 : if (!F2v_coeff(set,g[h[u]])) break;
1674 : }
1675 22386 : if (j < lgen) continue;
1676 1386 : F2v_set(res,i);
1677 84546 : for(j=1; j <= n; j++)
1678 83160 : if (F2v_coeff(set, j))
1679 6720 : F2v_clear(b,g[gel(elt,j)[1]]);
1680 : }
1681 812 : return res;
1682 : }
1683 :
1684 : static GEN
1685 14 : triv_subg(GEN elt)
1686 : {
1687 14 : GEN v = cgetg(3, t_VEC);
1688 14 : gel(v,1) = cgetg(1,t_VECSMALL);
1689 14 : gel(v,2) = zero_F2v(lg(elt)-1);
1690 14 : F2v_set(gel(v,2),1);
1691 14 : return v;
1692 : }
1693 :
1694 : static GEN
1695 364 : subg_extend(GEN U, long e, long o, GEN elt)
1696 : {
1697 364 : long i, j, n = lg(elt)-1;
1698 364 : GEN g = gel(elt, e);
1699 364 : GEN gen = vecsmall_append(subg_get_gen(U), e);
1700 364 : GEN set = subg_get_set(U);
1701 364 : GEN Vset = zv_copy(set);
1702 22204 : for(i = 1; i <= n; i++)
1703 21840 : if (F2v_coeff(set, i))
1704 : {
1705 1260 : long h = gel(elt, i)[1];
1706 2800 : for(j = 1; j < o; j++)
1707 : {
1708 1540 : h = g[h];
1709 1540 : F2v_set(Vset, h);
1710 : }
1711 : }
1712 364 : return mkvec2(gen, Vset);
1713 : }
1714 :
1715 : static GEN
1716 434 : cyclic_subg(long e, long o, GEN elt)
1717 : {
1718 434 : long j, n = lg(elt)-1, h = 1;
1719 434 : GEN g = gel(elt, e);
1720 434 : GEN gen = mkvecsmall(e);
1721 434 : GEN set = zero_F2v(n);
1722 434 : F2v_set(set,1);
1723 1260 : for(j = 1; j < o; j++)
1724 : {
1725 826 : h = g[h];
1726 826 : F2v_set(set, h);
1727 : }
1728 434 : return mkvec2(gen, set);
1729 : }
1730 :
1731 : static GEN
1732 14 : groupelts_to_regular(GEN elt)
1733 : {
1734 14 : long i, j, n = lg(elt)-1;
1735 14 : GEN V = cgetg(n+1,t_VEC);
1736 854 : for (i=1; i<=n; i++)
1737 : {
1738 840 : pari_sp av = avma;
1739 840 : GEN g = gel(elt, i);
1740 840 : GEN W = cgetg(n+1,t_VEC);
1741 51240 : for(j=1; j<=n; j++)
1742 50400 : gel(W,j) = perm_mul(g, gel(elt,j));
1743 840 : gel(V, i) = gerepileuptoleaf(av,vecvecsmall_indexsort(W));
1744 : }
1745 14 : vecvecsmall_sort_inplace(V, NULL);
1746 14 : return V;
1747 : }
1748 :
1749 : static long
1750 434 : groupelts_pow(GEN elt, long j, long n)
1751 : {
1752 434 : GEN g = gel(elt,j);
1753 434 : long i, h = 1;
1754 1694 : for (i=1; i<=n; i++)
1755 1260 : h = g[h];
1756 434 : return h;
1757 : }
1758 :
1759 : static GEN
1760 14 : groupelts_cyclic_primepow(GEN elt, GEN *pt_pr, GEN *pt_po)
1761 : {
1762 14 : GEN R = groupelts_cyclic_subgroups(elt);
1763 14 : GEN gen = gel(R,1), ord = gel(R,2);
1764 14 : long i, n = lg(elt)-1, l = lg(gen);
1765 14 : GEN set = zero_F2v(n);
1766 14 : GEN pr = zero_Flv(n);
1767 14 : GEN po = zero_Flv(n);
1768 462 : for (i = 1; i < l; i++)
1769 : {
1770 448 : long h = gen[i];
1771 : ulong p;
1772 448 : if (uisprimepower(ord[i], &p))
1773 : {
1774 434 : F2v_set(set, h);
1775 434 : uel(pr,h) = p;
1776 434 : po[h] = groupelts_pow(elt, h, p);
1777 : }
1778 : }
1779 14 : *pt_pr = pr; *pt_po = po;
1780 14 : return set;
1781 : }
1782 :
1783 : static GEN
1784 50400 : perm_bracket(GEN p, GEN q)
1785 : {
1786 50400 : return perm_mul(perm_mul(p,q), perm_inv(perm_mul(q,p)));
1787 : }
1788 :
1789 : static GEN
1790 826 : set_groupelts(GEN S, GEN x)
1791 : {
1792 826 : long i, n = F2v_hamming(x), k=1, m = x[1];
1793 826 : GEN v = cgetg(n+1, t_VEC);
1794 50386 : for (i=1; i<=m; i++)
1795 49560 : if (F2v_coeff(x,i))
1796 4914 : gel(v,k++) = gel(S,i);
1797 826 : return v;
1798 : }
1799 :
1800 : static GEN
1801 14 : set_idx(GEN x)
1802 : {
1803 14 : long i, n = F2v_hamming(x), k=1, m = x[1];
1804 14 : GEN v = cgetg(n+1, t_VECSMALL);
1805 854 : for (i=1; i<=m; i++)
1806 840 : if (F2v_coeff(x,i))
1807 840 : uel(v,k++) = i;
1808 14 : return v;
1809 : }
1810 :
1811 : static GEN
1812 14 : set_derived(GEN set, GEN elts)
1813 : {
1814 14 : long i, j, l = lg(elts);
1815 14 : GEN V = zero_F2v(l-1);
1816 854 : for(i = 1; i < l; i++)
1817 840 : if (F2v_coeff(set, i))
1818 51240 : for(j = 1; j < l; j++)
1819 50400 : if (F2v_coeff(set, j))
1820 50400 : F2v_set(V, perm_bracket(gel(elts,i),gel(elts,j))[1]);
1821 14 : return V;
1822 : }
1823 :
1824 : static GEN
1825 14 : groupelts_residuum(GEN elts)
1826 : {
1827 14 : pari_sp av = avma;
1828 14 : long o = lg(elts)-1, oo;
1829 14 : GEN set = const_F2v(o);
1830 : do
1831 : {
1832 14 : oo = o;
1833 14 : set = set_derived(set, elts);
1834 14 : o = F2v_hamming(set);
1835 14 : } while (o > 1 && o < oo);
1836 14 : if (o==1) return NULL;
1837 14 : return gerepilecopy(av,mkvec2(set_idx(set), set));
1838 : }
1839 :
1840 : static GEN
1841 14 : all_cyclic_subg(GEN pr, GEN po, GEN elt)
1842 : {
1843 14 : long i, n = lg(pr)-1, m = 0, k = 1;
1844 : GEN W;
1845 854 : for (i=1; i <= n; i++)
1846 840 : m += po[i]==1;
1847 14 : W = cgetg(m+1, t_VEC);
1848 854 : for (i=1; i <= n; i++)
1849 840 : if (po[i]==1)
1850 434 : gel(W, k++) = cyclic_subg(i, pr[i], elt);
1851 14 : return W;
1852 : }
1853 :
1854 : static GEN
1855 14 : groupelts_subgroups_raw(GEN elts)
1856 : {
1857 14 : pari_sp av = avma;
1858 14 : GEN elt = groupelts_to_regular(elts);
1859 14 : GEN pr, po, cyc = groupelts_cyclic_primepow(elt, &pr, &po);
1860 14 : long n = lg(elt)-1;
1861 14 : long i, j, nS = 1;
1862 14 : GEN S, L, R = NULL;
1863 14 : S = cgetg(1+bigomegau(n)+1, t_VEC);
1864 14 : gel(S, nS++) = mkvec(triv_subg(elt));
1865 14 : gel(S, nS++) = L = all_cyclic_subg(pr, po, elt);
1866 14 : if (DEBUGLEVEL) err_printf("subgroups: level %ld: %ld\n",nS-1,lg(L)-1);
1867 70 : while (lg(L) > 1)
1868 : {
1869 56 : pari_sp av2 = avma;
1870 56 : long nW = 1, lL = lg(L);
1871 56 : long ng = n;
1872 56 : GEN W = cgetg(1+ng, t_VEC);
1873 868 : for (i=1; i<lL; i++)
1874 : {
1875 812 : GEN U = gel(L, i), set = subg_get_set(U);
1876 812 : GEN G = groupelt_subg_normalize(elt, U, cyc);
1877 7154 : for (j=1; j<nW; j++)
1878 : {
1879 6342 : GEN Wj = subg_get_set(gel(W, j));
1880 6342 : if (F2v_subset(set, Wj))
1881 728 : F2v_negimply_inplace(G, Wj);
1882 : }
1883 49532 : for (j=1; j<=n; j++)
1884 48720 : if(F2v_coeff(G,j))
1885 : {
1886 770 : long p = pr[j];
1887 770 : if (F2v_coeff(set, j)) continue;
1888 364 : if (F2v_coeff(set, po[j]))
1889 : {
1890 364 : GEN U2 = subg_extend(U, j, p, elt);
1891 364 : F2v_negimply_inplace(G, subg_get_set(U2));
1892 364 : if (nW > ng) { ng<<=1; W = vec_lengthen(W, ng); }
1893 364 : gel(W, nW++) = U2;
1894 : }
1895 : }
1896 : }
1897 56 : setlg(W, nW);
1898 56 : L = W;
1899 56 : if (nW > 1) gel(S, nS++) = L = gerepilecopy(av2, W);
1900 56 : if (DEBUGLEVEL) err_printf("subgroups: level %ld: %ld\n",nS-1,nW-1);
1901 56 : if (lg(L)==1 && !R)
1902 : {
1903 14 : R = groupelts_residuum(elt);
1904 14 : if (!R) break;
1905 14 : gel(S, nS++) = L = mkvec(R);
1906 : }
1907 : }
1908 14 : setlg(S, nS);
1909 14 : return gerepilecopy(av, shallowconcat1(S));
1910 : }
1911 :
1912 : static GEN
1913 14 : subg_to_elts(GEN S, GEN x)
1914 840 : { pari_APPLY_type(t_VEC, set_groupelts(S, gmael(x,i,2))); }
1915 :
1916 : GEN
1917 14 : groupelts_solvablesubgroups(GEN G)
1918 : {
1919 14 : pari_sp av = avma;
1920 14 : GEN S = vecvecsmall_sort(checkgroupelts(G));
1921 14 : GEN L = groupelts_subgroups_raw(S);
1922 14 : return gerepilecopy(av, subg_to_elts(S, L));
1923 : }
|