Code coverage tests

This page documents the degree to which the PARI/GP source code is tested by our public test suite, distributed with the source distribution in directory src/test/. This is measured by the gcov utility; we then process gcov output using the lcov frond-end.

We test a few variants depending on Configure flags on the pari.math.u-bordeaux.fr machine (x86_64 architecture), and agregate them in the final report:

The target is 90% coverage for all mathematical modules (given that branches depending on DEBUGLEVEL or DEBUGMEM are not covered). This script is run to produce the results below.

LCOV - code coverage report
Current view: top level - basemath - perm.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.11.0 lcov report (development 22860-5579deb0b) Lines: 653 700 93.3 %
Date: 2018-07-18 05:36:42 Functions: 75 79 94.9 %
Legend: Lines: hit not hit

          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. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : 
      14             : #include "pari.h"
      15             : #include "paripriv.h"
      16             : 
      17             : /*************************************************************************/
      18             : /**                                                                     **/
      19             : /**                   Routines for handling VEC/COL                     **/
      20             : /**                                                                     **/
      21             : /*************************************************************************/
      22             : int
      23        5784 : vec_isconst(GEN v)
      24             : {
      25        5784 :   long i, l = lg(v);
      26             :   GEN w;
      27        5784 :   if (l==1) return 1;
      28        5784 :   w = gel(v,1);
      29       13030 :   for(i=2;i<l;i++)
      30       11083 :     if (!gequal(gel(v,i), w)) return 0;
      31        1947 :   return 1;
      32             : }
      33             : 
      34             : /* Check if all the elements of v are different.
      35             :  * Use a quadratic algorithm. Could be done in n*log(n) by sorting. */
      36             : int
      37        2731 : vec_is1to1(GEN v)
      38             : {
      39        2731 :   long i, j, l = lg(v);
      40       16922 :   for (i=1; i<l; i++)
      41             :   {
      42       14331 :     GEN w = gel(v,i);
      43       91801 :     for(j=i+1; j<l; j++)
      44       77610 :       if (gequal(gel(v,j), w)) return 0;
      45             :   }
      46        2591 :   return 1;
      47             : }
      48             : 
      49             : GEN
      50       97538 : vec_insert(GEN v, long n, GEN x)
      51             : {
      52       97538 :   long i, l=lg(v);
      53       97538 :   GEN V = cgetg(l+1,t_VEC);
      54       97538 :   for(i=1; i<n; i++) gel(V,i) = gel(v,i);
      55       97538 :   gel(V,n) = x;
      56       97538 :   for(i=n+1; i<=l; i++) gel(V,i) = gel(v,i-1);
      57       97538 :   return V;
      58             : }
      59             : /*************************************************************************/
      60             : /**                                                                     **/
      61             : /**                   Routines for handling VECSMALL                    **/
      62             : /**                                                                     **/
      63             : /*************************************************************************/
      64             : /* Sort v[0]...v[n-1] and put result in w[0]...w[n-1].
      65             :  * We accept v==w. w must be allocated. */
      66             : static void
      67    97165495 : vecsmall_sortspec(GEN v, long n, GEN w)
      68             : {
      69    97165495 :   pari_sp ltop=avma;
      70    97165495 :   long nx=n>>1, ny=n-nx;
      71             :   long m, ix, iy;
      72             :   GEN x, y;
      73    97165495 :   if (n<=2)
      74             :   {
      75    54680033 :     if (n==1)
      76    10680281 :       w[0]=v[0];
      77    43999752 :     else if (n==2)
      78             :     {
      79    43999756 :       long v0=v[0], v1=v[1];
      80    43999756 :       if (v0<=v1) { w[0]=v0; w[1]=v1; }
      81     1871281 :       else        { w[0]=v1; w[1]=v0; }
      82             :     }
      83    54680033 :     return;
      84             :   }
      85    42485462 :   x=new_chunk(nx); y=new_chunk(ny);
      86    42485463 :   vecsmall_sortspec(v,nx,x);
      87    42485463 :   vecsmall_sortspec(v+nx,ny,y);
      88   228589798 :   for (m=0, ix=0, iy=0; ix<nx && iy<ny; )
      89   143618868 :     if (x[ix]<=y[iy])
      90   120442243 :       w[m++]=x[ix++];
      91             :     else
      92    23176625 :       w[m++]=y[iy++];
      93    42485465 :   for(;ix<nx;) w[m++]=x[ix++];
      94    42485465 :   for(;iy<ny;) w[m++]=y[iy++];
      95    42485465 :   avma=ltop;
      96             : }
      97             : 
      98             : /*in place sort.*/
      99             : void
     100    18140866 : vecsmall_sort(GEN V)
     101             : {
     102    18140866 :   long l = lg(V)-1;
     103    18140866 :   if (l<=1) return;
     104    12194571 :   vecsmall_sortspec(V+1,l,V+1);
     105             : }
     106             : 
     107             : /* cf gen_sortspec */
     108             : static GEN
     109    21084789 : vecsmall_indexsortspec(GEN v, long n)
     110             : {
     111             :   long nx, ny, m, ix, iy;
     112             :   GEN x, y, w;
     113    21084789 :   switch(n)
     114             :   {
     115       50057 :     case 1: return mkvecsmall(1);
     116     5574249 :     case 2: return (v[1] <= v[2])? mkvecsmall2(1,2): mkvecsmall2(2,1);
     117             :     case 3:
     118     6369027 :       if (v[1] <= v[2]) {
     119     5458631 :         if (v[2] <= v[3]) return mkvecsmall3(1,2,3);
     120     1116465 :         return (v[1] <= v[3])? mkvecsmall3(1,3,2)
     121     1116465 :                              : mkvecsmall3(3,1,2);
     122             :       } else {
     123      910396 :         if (v[1] <= v[3]) return mkvecsmall3(2,1,3);
     124      643850 :         return (v[2] <= v[3])? mkvecsmall3(2,3,1)
     125      643850 :                              : mkvecsmall3(3,2,1);
     126             :       }
     127             :   }
     128     9091456 :   nx = n>>1; ny = n-nx;
     129     9091456 :   w = cgetg(n+1,t_VECSMALL);
     130     9091456 :   x = vecsmall_indexsortspec(v,nx);
     131     9091457 :   y = vecsmall_indexsortspec(v+nx,ny);
     132   156631014 :   for (m=1, ix=1, iy=1; ix<=nx && iy<=ny; )
     133   138448098 :     if (v[x[ix]] <= v[y[iy]+nx])
     134    89913681 :       w[m++] = x[ix++];
     135             :     else
     136    48534417 :       w[m++] = y[iy++]+nx;
     137     9091458 :   for(;ix<=nx;) w[m++] = x[ix++];
     138     9091458 :   for(;iy<=ny;) w[m++] = y[iy++]+nx;
     139     9091458 :   avma = (pari_sp)w; return w;
     140             : }
     141             : 
     142             : /*indirect sort.*/
     143             : GEN
     144     2901937 : vecsmall_indexsort(GEN V)
     145             : {
     146     2901937 :   long l=lg(V)-1;
     147     2901937 :   if (l==0) return cgetg(1, t_VECSMALL);
     148     2901874 :   return vecsmall_indexsortspec(V,l);
     149             : }
     150             : 
     151             : /* assume V sorted */
     152             : GEN
     153        1449 : vecsmall_uniq_sorted(GEN V)
     154             : {
     155             :   GEN W;
     156        1449 :   long i,j, l = lg(V);
     157        1449 :   if (l == 1) return vecsmall_copy(V);
     158        1428 :   W = cgetg(l,t_VECSMALL);
     159        1428 :   W[1] = V[1];
     160        2737 :   for(i=j=2; i<l; i++)
     161        1309 :     if (V[i] != W[j-1]) W[j++] = V[i];
     162        1428 :   stackdummy((pari_sp)(W + l), (pari_sp)(W + j));
     163        1428 :   setlg(W, j); return W;
     164             : }
     165             : 
     166             : GEN
     167        1372 : vecsmall_uniq(GEN V)
     168             : {
     169        1372 :   pari_sp av = avma;
     170        1372 :   V = zv_copy(V); vecsmall_sort(V);
     171        1372 :   return gerepileuptoleaf(av, vecsmall_uniq_sorted(V));
     172             : }
     173             : 
     174             : /* assume x sorted */
     175             : long
     176           0 : vecsmall_duplicate_sorted(GEN x)
     177             : {
     178           0 :   long i,k,l=lg(x);
     179           0 :   if (l==1) return 0;
     180           0 :   for (k=x[1],i=2; i<l; k=x[i++])
     181           0 :     if (x[i] == k) return i;
     182           0 :   return 0;
     183             : }
     184             : 
     185             : long
     186       14362 : vecsmall_duplicate(GEN x)
     187             : {
     188       14362 :   pari_sp av=avma;
     189       14362 :   GEN p=vecsmall_indexsort(x);
     190       14362 :   long k,i,r=0,l=lg(x);
     191       14362 :   if (l==1) return 0;
     192       19519 :   for (k=x[p[1]],i=2; i<l; k=x[p[i++]])
     193        5157 :     if (x[p[i]] == k) { r=p[i]; break; }
     194       14362 :   avma=av;
     195       14362 :   return r;
     196             : }
     197             : 
     198             : /*************************************************************************/
     199             : /**                                                                     **/
     200             : /**             Routines for handling vectors of VECSMALL               **/
     201             : /**                                                                     **/
     202             : /*************************************************************************/
     203             : 
     204             : GEN
     205      121072 : vecvecsmall_sort(GEN x)
     206             : {
     207      121072 :   return gen_sort(x, (void*)&vecsmall_lexcmp, cmp_nodata);
     208             : }
     209             : 
     210             : GEN
     211         392 : vecvecsmall_sort_uniq(GEN x)
     212             : {
     213         392 :   return gen_sort_uniq(x, (void*)&vecsmall_lexcmp, cmp_nodata);
     214             : }
     215             : 
     216             : GEN
     217          21 : vecvecsmall_indexsort(GEN x)
     218             : {
     219          21 :   return gen_indexsort(x, (void*)&vecsmall_lexcmp, cmp_nodata);
     220             : }
     221             : 
     222             : long
     223    19730179 : vecvecsmall_search(GEN x, GEN y, long flag)
     224             : {
     225    19730179 :   return gen_search(x,y,flag,(void*)vecsmall_prefixcmp, cmp_nodata);
     226             : }
     227             : 
     228             : /* assume x non empty */
     229             : long
     230         133 : vecvecsmall_max(GEN x)
     231             : {
     232         133 :   long i, l = lg(x), m = vecsmall_max(gel(x,1));
     233        1099 :   for (i = 2; i < l; i++)
     234             :   {
     235         966 :     long t = vecsmall_max(gel(x,i));
     236         966 :     if (t > m) m = t;
     237             :   }
     238         133 :   return m;
     239             : }
     240             : 
     241             : /*************************************************************************/
     242             : /**                                                                     **/
     243             : /**                  Routines for handling permutations                 **/
     244             : /**                                                                     **/
     245             : /*************************************************************************/
     246             : 
     247             : /* Permutations may be given by
     248             :  * perm (VECSMALL): a bijection from 1...n to 1...n i-->perm[i]
     249             :  * cyc (VEC of VECSMALL): a product of disjoint cycles. */
     250             : 
     251             : /* Multiply (compose) two permutations, putting the result in the second one. */
     252             : static void
     253           7 : perm_mul_inplace2(GEN s, GEN t)
     254             : {
     255           7 :   long i, l = lg(s);
     256           7 :   for (i = 1; i < l; i++) t[i] = s[t[i]];
     257           7 : }
     258             : 
     259             : /* Orbits of the subgroup generated by v on {1,..,n} */
     260             : static GEN
     261      571430 : vecperm_orbits_i(GEN v, long n)
     262             : {
     263      571430 :   long mj = 1, lv = lg(v), k, l;
     264      571430 :   GEN cycle = cgetg(n+1, t_VEC), bit = const_vecsmall(n, 0);
     265     4653242 :   for (k = 1, l = 1; k <= n;)
     266             :   {
     267     3510375 :     long m = 1;
     268     3510375 :     GEN cy = cgetg(n+1, t_VECSMALL);
     269     3510373 :     for (  ; bit[mj]; mj++) /*empty*/;
     270     3510373 :     k++; cy[m++] = mj;
     271     3510373 :     bit[mj++] = 1;
     272             :     for(;;)
     273     1627126 :     {
     274     5137499 :       long o, mold = m;
     275    10278005 :       for (o = 1; o < lv; o++)
     276             :       {
     277     5140506 :         GEN vo = gel(v,o);
     278             :         long p;
     279    19028265 :         for (p = 1; p < m; p++) /* m increases! */
     280             :         {
     281    13887759 :           long j = vo[ cy[p] ];
     282    13887759 :           if (!bit[j]) cy[m++] = j;
     283    13887759 :           bit[j] = 1;
     284             :         }
     285             :       }
     286     5137499 :       if (m == mold) break;
     287     1627126 :       k += m - mold;
     288             :     }
     289     3510373 :     setlg(cy, m); gel(cycle,l++) = cy;
     290             :   }
     291      571432 :   setlg(cycle, l); return cycle;
     292             : }
     293             : /* memory clean version */
     294             : GEN
     295         826 : vecperm_orbits(GEN v, long n)
     296             : {
     297         826 :   pari_sp av = avma;
     298         826 :   return gerepilecopy(av, vecperm_orbits_i(v, n));
     299             : }
     300             : 
     301             : /* Compute the cyclic decomposition of a permutation */
     302             : GEN
     303        4055 : perm_cycles(GEN v)
     304             : {
     305        4055 :   pari_sp av = avma;
     306        4055 :   return gerepilecopy(av, vecperm_orbits_i(mkvec(v), lg(v)-1));
     307             : }
     308             : 
     309             : static long
     310        1890 : isperm(GEN v)
     311             : {
     312        1890 :   pari_sp av = avma;
     313        1890 :   long i, n = lg(v)-1;
     314             :   GEN w;
     315        1890 :   if (typ(v) != t_VECSMALL) return 0;
     316        1890 :   w = zero_zv(n);
     317       12404 :   for (i=1; i<=n; i++)
     318             :   {
     319       10542 :     long d = v[i];
     320       10542 :     if (d < 1 || d > n || w[d]) { avma = av; return 0; }
     321       10514 :     w[d] = 1;
     322             :   }
     323        1862 :   avma = av; return 1;
     324             : }
     325             : 
     326             : /* Output the order of p */
     327             : long
     328      388920 : perm_order(GEN v)
     329             : {
     330      388920 :   pari_sp ltop = avma;
     331      388920 :   GEN c = vecperm_orbits_i(mkvec(v), lg(v)-1);
     332             :   long i, d;
     333      388920 :   for(i=1, d=1; i<lg(c); i++) d = ulcm(d, lg(gel(c,i))-1);
     334      388920 :   avma = ltop; return d;
     335             : }
     336             : 
     337             : long
     338          91 : permorder(GEN v)
     339             : {
     340          91 :   if (!isperm(v)) pari_err_TYPE("permorder",v);
     341          84 :   return perm_order(v);
     342             : }
     343             : 
     344             : /* sign of a permutation */
     345             : long
     346      177629 : perm_sign(GEN v)
     347             : {
     348      177629 :   pari_sp av = avma;
     349      177629 :   GEN c = vecperm_orbits_i(mkvec(v), lg(v)-1);
     350      177630 :   long i, l = lg(c), s = 1;
     351     1512714 :   for (i = 1; i < l; i++)
     352     1335085 :     if (odd(lg(gel(c, i)))) s = -s;
     353      177629 :   avma = av; return s;
     354             : }
     355             : 
     356             : long
     357          98 : permsign(GEN v)
     358             : {
     359          98 :   if (!isperm(v)) pari_err_TYPE("permsign",v);
     360          84 :   return perm_sign(v);
     361             : }
     362             : 
     363             : GEN
     364        5915 : Z_to_perm(long n, GEN x)
     365             : {
     366             :   pari_sp av;
     367             :   ulong i, r;
     368        5915 :   GEN v = cgetg(n+1, t_VECSMALL);
     369        5915 :   if (n==0) return v;
     370        5908 :   uel(v,n) = 1; av = avma;
     371        5908 :   if (signe(x) <= 0) x = modii(x, mpfact(n));
     372       27146 :   for (r=n-1; r>=1; r--)
     373             :   {
     374             :     ulong a;
     375       21238 :     x = absdiviu_rem(x, n+1-r, &a);
     376       71687 :     for (i=r+1; i<=(ulong)n; i++)
     377       50449 :       if (uel(v,i) > a) uel(v,i)++;
     378       21238 :     uel(v,r) = a+1;
     379             :   }
     380        5908 :   avma = av; return v;
     381             : }
     382             : GEN
     383        5915 : numtoperm(long n, GEN x)
     384             : {
     385        5915 :   if (n < 0) pari_err_DOMAIN("numtoperm", "n", "<", gen_0, stoi(n));
     386        5915 :   if (typ(x) != t_INT) pari_err_TYPE("numtoperm",x);
     387        5915 :   return Z_to_perm(n, x);
     388             : }
     389             : 
     390             : /* destroys v */
     391             : static GEN
     392        1701 : perm_to_Z_inplace(GEN v)
     393             : {
     394        1701 :   long l = lg(v), i, r;
     395        1701 :   GEN x = gen_0;
     396        1701 :   if (!isperm(v)) pari_err_TYPE("permsign",v);
     397       10143 :   for (i = 1; i < l; i++)
     398             :   {
     399        8449 :     long vi = v[i];
     400        8449 :     if (vi <= 0) return NULL;
     401        8449 :     x = i==1 ? utoi(vi-1): addiu(muliu(x,l-i), vi-1);
     402       25396 :     for (r = i+1; r < l; r++)
     403       16947 :       if (v[r] > vi) v[r]--;
     404             :   }
     405        1694 :   return x;
     406             : }
     407             : GEN
     408        1680 : perm_to_Z(GEN v)
     409             : {
     410        1680 :   pari_sp av = avma;
     411        1680 :   GEN x = perm_to_Z_inplace(leafcopy(v));
     412        1680 :   if (!x) pari_err_TYPE("permtonum",v);
     413        1680 :   return gerepileuptoint(av, x);
     414             : }
     415             : GEN
     416        1708 : permtonum(GEN p)
     417             : {
     418        1708 :   pari_sp av = avma;
     419             :   GEN v, x;
     420        1708 :   switch(typ(p))
     421             :   {
     422        1680 :     case t_VECSMALL: return perm_to_Z(p);
     423             :     case t_VEC: case t_COL:
     424          21 :       if (RgV_is_ZV(p)) { v = ZV_to_zv(p); break; }
     425           7 :     default: pari_err_TYPE("permtonum",p); return NULL;
     426             :   }
     427          21 :   x = perm_to_Z_inplace(v);
     428          14 :   if (!x) pari_err_TYPE("permtonum",p);
     429          14 :   return gerepileuptoint(av, x);
     430             : }
     431             : 
     432             : 
     433             : GEN
     434        1891 : cyc_pow(GEN cyc, long exp)
     435             : {
     436             :   long i, j, k, l, r;
     437             :   GEN c;
     438        6344 :   for (r = j = 1; j < lg(cyc); j++)
     439             :   {
     440        4453 :     long n = lg(gel(cyc,j)) - 1;
     441        4453 :     r += cgcd(n, exp);
     442             :   }
     443        1891 :   c = cgetg(r, t_VEC);
     444        6344 :   for (r = j = 1; j < lg(cyc); j++)
     445             :   {
     446        4453 :     GEN v = gel(cyc,j);
     447        4453 :     long n = lg(v) - 1, e = smodss(exp,n), g = (long)ugcd(n, e), m = n / g;
     448        9725 :     for (i = 0; i < g; i++)
     449             :     {
     450        5272 :       GEN p = cgetg(m+1, t_VECSMALL);
     451        5272 :       gel(c,r++) = p;
     452       17790 :       for (k = 1, l = i; k <= m; k++)
     453             :       {
     454       12518 :         p[k] = v[l+1];
     455       12518 :         l += e; if (l >= n) l -= n;
     456             :       }
     457             :     }
     458             :   }
     459        1891 :   return c;
     460             : }
     461             : 
     462             : /* Compute the power of a permutation given by product of cycles
     463             :  * Ouput a perm, not a cyc */
     464             : GEN
     465           0 : cyc_pow_perm(GEN cyc, long exp)
     466             : {
     467             :   long e, j, k, l, n;
     468             :   GEN p;
     469           0 :   for (n = 0, j = 1; j < lg(cyc); j++) n += lg(gel(cyc,j))-1;
     470           0 :   p = cgetg(n + 1, t_VECSMALL);
     471           0 :   for (j = 1; j < lg(cyc); j++)
     472             :   {
     473           0 :     GEN v = gel(cyc,j);
     474           0 :     n = lg(v) - 1; e = smodss(exp, n);
     475           0 :     for (k = 1, l = e; k <= n; k++)
     476             :     {
     477           0 :       p[v[k]] = v[l+1];
     478           0 :       if (++l == n) l = 0;
     479             :     }
     480             :   }
     481           0 :   return p;
     482             : }
     483             : 
     484             : GEN
     485        8359 : perm_pow(GEN perm, long exp)
     486             : {
     487        8359 :   long i, r = lg(perm)-1;
     488        8359 :   GEN p = zero_zv(r);
     489        8359 :   pari_sp av = avma;
     490        8359 :   GEN v = cgetg(r+1, t_VECSMALL);
     491      142018 :   for (i=1; i<=r; i++)
     492             :   {
     493             :     long e, n, k, l;
     494      133659 :     if (p[i]) continue;
     495       50894 :     v[1] = i;
     496      133659 :     for (n=1, k=perm[i]; k!=i; k=perm[k], n++)
     497       82765 :       v[n+1] = k;
     498       50894 :     e = smodss(exp, n);
     499      184553 :     for (k = 1, l = e; k <= n; k++)
     500             :     {
     501      133659 :       p[v[k]] = v[l+1];
     502      133659 :       if (++l == n) l = 0;
     503             :     }
     504             :   }
     505        8359 :   avma = av; return p;
     506             : }
     507             : 
     508             : static GEN
     509          21 : perm_to_GAP(GEN p)
     510             : {
     511          21 :   pari_sp ltop=avma;
     512             :   GEN gap;
     513             :   GEN x;
     514             :   long i;
     515          21 :   long nb, c=0;
     516             :   char *s;
     517             :   long sz;
     518          21 :   long lp=lg(p)-1;
     519          21 :   if (typ(p) != t_VECSMALL)  pari_err_TYPE("perm_to_GAP",p);
     520          21 :   x = perm_cycles(p);
     521          21 :   sz = (long) ((bfffo(lp)+1) * LOG10_2 + 1);
     522             :   /*Dry run*/
     523         133 :   for (i = 1, nb = 1; i < lg(x); ++i)
     524             :   {
     525         112 :     GEN z = gel(x,i);
     526         112 :     long lz = lg(z)-1;
     527         112 :     nb += 1+lz*(sz+2);
     528             :   }
     529          21 :   nb++;
     530             :   /*Real run*/
     531          21 :   gap = cgetg(nchar2nlong(nb) + 1, t_STR);
     532          21 :   s = GSTR(gap);
     533         133 :   for (i = 1; i < lg(x); ++i)
     534             :   {
     535             :     long j;
     536         112 :     GEN z = gel(x,i);
     537         112 :     if (lg(z) > 2)
     538             :     {
     539         112 :       s[c++] = '(';
     540         364 :       for (j = 1; j < lg(z); ++j)
     541             :       {
     542         252 :         if (j > 1)
     543             :         {
     544         140 :           s[c++] = ','; s[c++] = ' ';
     545             :         }
     546         252 :         sprintf(s+c,"%ld",z[j]);
     547         252 :         while(s[c++]) /* empty */;
     548         252 :         c--;
     549             :       }
     550         112 :       s[c++] = ')';
     551             :     }
     552             :   }
     553          21 :   if (!c) { s[c++]='('; s[c++]=')'; }
     554          21 :   s[c] = '\0';
     555          21 :   return gerepileupto(ltop,gap);
     556             : }
     557             : 
     558             : int
     559      536711 : perm_commute(GEN s, GEN t)
     560             : {
     561      536711 :   long i, l = lg(t);
     562    38374455 :   for (i = 1; i < l; i++)
     563    37851499 :     if (t[ s[i] ] != s[ t[i] ]) return 0;
     564      522956 :   return 1;
     565             : }
     566             : 
     567             : /*************************************************************************/
     568             : /**                                                                     **/
     569             : /**                  Routines for handling groups                       **/
     570             : /**                                                                     **/
     571             : /*************************************************************************/
     572             : /* A Group is a t_VEC [gen,orders]
     573             :  * gen (vecvecsmall): list of generators given by permutations
     574             :  * orders (vecsmall): relatives orders of generators. */
     575      432699 : INLINE GEN grp_get_gen(GEN G) { return gel(G,1); }
     576      738551 : INLINE GEN grp_get_ord(GEN G) { return gel(G,2); }
     577             : 
     578             : /* A Quotient Group is a t_VEC [gen,coset]
     579             :  * gen (vecvecsmall): coset generators
     580             :  * coset (vecsmall): gen[coset[p[1]]] generate the p-coset.
     581             :  */
     582       76328 : INLINE GEN quo_get_gen(GEN C) { return gel(C,1); }
     583       12215 : INLINE GEN quo_get_coset(GEN C) { return gel(C,2); }
     584             : 
     585             : static GEN
     586       30744 : trivialsubgroups(void)
     587       30744 : { GEN L = cgetg(2, t_VEC); gel(L,1) = trivialgroup(); return L; }
     588             : 
     589             : /* Compute the order of p modulo the group given by a set */
     590             : long
     591      122906 : perm_relorder(GEN p, GEN set)
     592             : {
     593      122906 :   pari_sp ltop = avma;
     594      122906 :   long n = 1;
     595      122906 :   long q = p[1];
     596      122906 :   while (!F2v_coeff(set,q)) { q = p[q]; n++; }
     597      122906 :   avma = ltop; return n;
     598             : }
     599             : 
     600             : GEN
     601        8071 : perm_generate(GEN S, GEN H, long o)
     602             : {
     603        8071 :   long i, n = lg(H)-1;
     604        8071 :   GEN L = cgetg(n*o + 1, t_VEC);
     605        8071 :   for(i=1; i<=n;     i++) gel(L,i) = vecsmall_copy(gel(H,i));
     606        8071 :   for(   ; i <= n*o; i++) gel(L,i) = perm_mul(gel(L,i-n), S);
     607        8071 :   return L;
     608             : }
     609             : 
     610             : /*Return the order (cardinality) of a group */
     611             : long
     612      318550 : group_order(GEN G)
     613             : {
     614      318550 :   return zv_prod(grp_get_ord(G));
     615             : }
     616             : 
     617             : /* G being a subgroup of S_n, output n */
     618             : long
     619        6587 : group_domain(GEN G)
     620             : {
     621        6587 :   GEN gen = grp_get_gen(G);
     622        6587 :   if (lg(gen) < 2) pari_err_DOMAIN("group_domain", "#G", "=", gen_1,G);
     623        6587 :   return lg(gel(gen,1)) - 1;
     624             : }
     625             : 
     626             : /*Left coset of g mod G: gG*/
     627             : GEN
     628      136451 : group_leftcoset(GEN G, GEN g)
     629             : {
     630      136451 :   GEN gen = grp_get_gen(G), ord = grp_get_ord(G);
     631      136451 :   GEN res = cgetg(group_order(G)+1, t_VEC);
     632             :   long i, j, k;
     633      136451 :   gel(res,1) = vecsmall_copy(g);
     634      136451 :   k = 1;
     635      253652 :   for (i = 1; i < lg(gen); i++)
     636             :   {
     637      117201 :     long c = k * (ord[i] - 1);
     638      117201 :     for (j = 1; j <= c; j++) gel(res,++k) = perm_mul(gel(res,j), gel(gen,i));
     639             :   }
     640      136451 :   return res;
     641             : }
     642             : /*Right coset of g mod G: Gg*/
     643             : GEN
     644       60480 : group_rightcoset(GEN G, GEN g)
     645             : {
     646       60480 :   GEN gen = grp_get_gen(G), ord = grp_get_ord(G);
     647       60480 :   GEN res = cgetg(group_order(G)+1, t_VEC);
     648             :   long i, j, k;
     649       60480 :   gel(res,1) = vecsmall_copy(g);
     650       60480 :   k = 1;
     651      101598 :   for (i = 1; i < lg(gen); i++)
     652             :   {
     653       41118 :     long c = k * (ord[i] - 1);
     654       41118 :     for (j = 1; j <= c; j++) gel(res,++k) = perm_mul(gel(gen,i), gel(res,j));
     655             :   }
     656       60480 :   return res;
     657             : }
     658             : /*Elements of a group from the generators, cf group_leftcoset*/
     659             : GEN
     660       70995 : group_elts(GEN G, long n)
     661             : {
     662       70995 :   GEN gen = grp_get_gen(G), ord = grp_get_ord(G);
     663       70995 :   GEN res = cgetg(group_order(G)+1, t_VEC);
     664             :   long i, j, k;
     665       70995 :   gel(res,1) = identity_perm(n);
     666       70995 :   k = 1;
     667      142992 :   for (i = 1; i < lg(gen); i++)
     668             :   {
     669       71997 :     long c = k * (ord[i] - 1);
     670             :     /* j = 1, use res[1] = identity */
     671       71997 :     gel(res,++k) = vecsmall_copy(gel(gen,i));
     672       71997 :     for (j = 2; j <= c; j++) gel(res,++k) = perm_mul(gel(res,j), gel(gen,i));
     673             :   }
     674       70995 :   return res;
     675             : }
     676             : 
     677             : GEN
     678       13188 : groupelts_set(GEN elts, long n)
     679             : {
     680       13188 :   GEN res = zero_F2v(n);
     681       13188 :   long i, l = lg(elts);
     682       65625 :   for(i=1; i<l; i++)
     683       52437 :     F2v_set(res,mael(elts,i,1));
     684       13188 :   return res;
     685             : }
     686             : 
     687             : /*Elements of a group from the generators, returned as a set (bitmap)*/
     688             : GEN
     689       58338 : group_set(GEN G, long n)
     690             : {
     691       58338 :   GEN res = zero_F2v(n);
     692       58338 :   pari_sp av = avma;
     693       58338 :   GEN elts = group_elts(G, n);
     694       58338 :   long i, l = lg(elts);
     695      189203 :   for(i=1; i<l; i++)
     696      130865 :     F2v_set(res,mael(elts,i,1));
     697       58338 :   avma = av;
     698       58338 :   return res;
     699             : }
     700             : 
     701             : static int
     702       17206 : sgcmp(GEN a, GEN b) { return vecsmall_lexcmp(gel(a,1),gel(b,1)); }
     703             : 
     704             : GEN
     705         490 : subgroups_tableset(GEN S, long n)
     706             : {
     707         490 :   long i, l = lg(S);
     708         490 :   GEN  v = cgetg(l, t_VEC);
     709        5390 :   for(i=1; i<l; i++)
     710        4900 :     gel(v,i) = mkvec2(group_set(gel(S,i), n), mkvecsmall(i));
     711         490 :   gen_sort_inplace(v,(void*)sgcmp,cmp_nodata, NULL);
     712         490 :   return v;
     713             : }
     714             : 
     715             : long
     716        1981 : tableset_find_index(GEN tbl, GEN set)
     717             : {
     718        1981 :   long i = tablesearch(tbl,mkvec2(set,mkvecsmall(0)),sgcmp);
     719        1981 :   if (!i) return 0;
     720        1981 :   return mael3(tbl,i,2,1);
     721             : }
     722             : 
     723             : GEN
     724       30744 : trivialgroup(void) { retmkvec2(cgetg(1,t_VEC), cgetg(1,t_VECSMALL)); }
     725             : /*Cyclic group generated by g of order s*/
     726             : GEN
     727        7049 : cyclicgroup(GEN g, long s)
     728        7049 : { retmkvec2(mkvec( vecsmall_copy(g) ),
     729             :             mkvecsmall(s)); }
     730             : /*Return the group generated by g1,g2 of relative orders s1,s2*/
     731             : GEN
     732         973 : dicyclicgroup(GEN g1, GEN g2, long s1, long s2)
     733         973 : { retmkvec2( mkvec2(vecsmall_copy(g1), vecsmall_copy(g2)),
     734             :              mkvecsmall2(s1, s2) ); }
     735             : 
     736             : /* return the quotient map G --> G/H */
     737             : /*The ouput is [gen,hash]*/
     738             : /* gen (vecvecsmall): coset generators
     739             :  * coset (vecsmall): vecsmall of coset number) */
     740             : GEN
     741        5124 : group_quotient(GEN G, GEN H)
     742             : {
     743        5124 :   pari_sp ltop = avma;
     744             :   GEN  p2, p3;
     745        5124 :   long i, j, a = 1;
     746        5124 :   long n = group_domain(G), o = group_order(H);
     747        5124 :   GEN  elt = group_elts(G,n), el;
     748        5124 :   long le = lg(elt)-1;
     749        5124 :   GEN used = zero_F2v(le+1);
     750        5124 :   long l = le/o;
     751        5124 :   p2 = cgetg(l+1, t_VEC);
     752        5124 :   p3 = zero_zv(n);
     753        5124 :   el = zero_zv(n);
     754       69160 :   for (i = 1; i<=le; i++)
     755       64036 :     el[mael(elt,i,1)]=i;
     756       36204 :   for (i = 1; i <= l; ++i)
     757             :   {
     758             :     GEN V;
     759       31087 :     while(F2v_coeff(used,a)) a++;
     760       31087 :     V = group_leftcoset(H,gel(elt,a));
     761       31087 :     gel(p2,i) = gel(V,1);
     762       95018 :     for(j=1;j<lg(V);j++)
     763             :     {
     764       63938 :       long b = el[mael(V,j,1)];
     765       63938 :       if (b==0) pari_err_IMPL("group_quotient for a non-WSS group");
     766       63931 :       F2v_set(used,b);
     767             :     }
     768       95004 :     for (j = 1; j <= o; j++)
     769       63924 :       p3[mael(V, j, 1)] = i;
     770             :   }
     771        5117 :   return gerepilecopy(ltop,mkvec2(p2,p3));
     772             : }
     773             : 
     774             : /*Compute the image of a permutation by a quotient map.*/
     775             : GEN
     776       12215 : quotient_perm(GEN C, GEN p)
     777             : {
     778       12215 :   GEN gen = quo_get_gen(C);
     779       12215 :   GEN coset = quo_get_coset(C);
     780       12215 :   long j, l = lg(gen);
     781       12215 :   GEN p3 = cgetg(l, t_VECSMALL);
     782      115850 :   for (j = 1; j < l; ++j)
     783             :   {
     784      103635 :     p3[j] = coset[p[mael(gen,j,1)]];
     785      103635 :     if (p3[j]==0) pari_err_IMPL("quotient_perm for a non-WSS group");
     786             :   }
     787       12215 :   return p3;
     788             : }
     789             : 
     790             : /* H is a subgroup of G, C is the quotient map G --> G/H
     791             :  *
     792             :  * Lift a subgroup S of G/H to a subgroup of G containing H */
     793             : GEN
     794       29498 : quotient_subgroup_lift(GEN C, GEN H, GEN S)
     795             : {
     796       29498 :   GEN genH = grp_get_gen(H);
     797       29498 :   GEN genS = grp_get_gen(S);
     798       29498 :   GEN genC = quo_get_gen(C);
     799       29498 :   long l1 = lg(genH)-1;
     800       29498 :   long l2 = lg(genS)-1, j;
     801       29498 :   GEN p1 = cgetg(3, t_VEC), L = cgetg(l1+l2+1, t_VEC);
     802       29498 :   for (j = 1; j <= l1; ++j) gel(L,j) = gel(genH,j);
     803       29498 :   for (j = 1; j <= l2; ++j) gel(L,l1+j) = gel(genC, mael(genS,j,1));
     804       29498 :   gel(p1,1) = L;
     805       29498 :   gel(p1,2) = vecsmall_concat(grp_get_ord(H), grp_get_ord(S));
     806       29498 :   return p1;
     807             : }
     808             : 
     809             : /* Let G a group and C a quotient map G --> G/H
     810             :  * Assume H is normal, return the group G/H */
     811             : GEN
     812        5117 : quotient_group(GEN C, GEN G)
     813             : {
     814        5117 :   pari_sp ltop = avma;
     815             :   GEN Qgen, Qord, Qelt, Qset, Q;
     816        5117 :   GEN Cgen = quo_get_gen(C);
     817        5117 :   GEN Ggen = grp_get_gen(G);
     818        5117 :   long i,j, n = lg(Cgen)-1, l = lg(Ggen);
     819        5117 :   Qord = cgetg(l, t_VECSMALL);
     820        5117 :   Qgen = cgetg(l, t_VEC);
     821        5117 :   Qelt = mkvec(identity_perm(n));
     822        5117 :   Qset = groupelts_set(Qelt, n);
     823       17332 :   for (i = 1, j = 1; i < l; ++i)
     824             :   {
     825       12215 :     GEN  g = quotient_perm(C, gel(Ggen,i));
     826       12215 :     long o = perm_relorder(g, Qset);
     827       12215 :     gel(Qgen,j) = g;
     828       12215 :     Qord[j] = o;
     829       12215 :     if (o != 1)
     830             :     {
     831        8071 :       Qelt = perm_generate(g, Qelt, o);
     832        8071 :       Qset = groupelts_set(Qelt, n);
     833        8071 :       j++;
     834             :     }
     835             :   }
     836        5117 :   setlg(Qgen,j);
     837        5117 :   setlg(Qord,j); Q = mkvec2(Qgen, Qord);
     838        5117 :   return gerepilecopy(ltop,Q);
     839             : }
     840             : 
     841             : /* Return 1 if g normalizes N, 0 otherwise */
     842             : long
     843       60480 : group_perm_normalize(GEN N, GEN g)
     844             : {
     845       60480 :   pari_sp ltop = avma;
     846       60480 :   long r = gequal(vecvecsmall_sort(group_leftcoset(N, g)),
     847             :                   vecvecsmall_sort(group_rightcoset(N, g)));
     848       60480 :   avma = ltop; return r;
     849             : }
     850             : 
     851             : /* L is a list of subgroups, C is a coset and r a relative order.*/
     852             : static GEN
     853       44884 : liftlistsubgroups(GEN L, GEN C, long r)
     854             : {
     855       44884 :   pari_sp ltop = avma;
     856       44884 :   long c = lg(C)-1, l = lg(L)-1, n = lg(gel(C,1))-1, i, k;
     857             :   GEN R;
     858       44884 :   if (!l) return cgetg(1,t_VEC);
     859       37443 :   R = cgetg(l*c+1, t_VEC);
     860       89565 :   for (i = 1, k = 1; i <= l; ++i)
     861             :   {
     862       52122 :     GEN S = gel(L,i), Selt = group_set(S,n);
     863       52122 :     GEN gen = grp_get_gen(S);
     864       52122 :     GEN ord = grp_get_ord(S);
     865             :     long j;
     866      158270 :     for (j = 1; j <= c; ++j)
     867             :     {
     868      106148 :       GEN p = gel(C,j);
     869      106148 :       if (perm_relorder(p, Selt) == r && group_perm_normalize(S, p))
     870       58275 :         gel(R,k++) = mkvec2(vec_append(gen, p),
     871             :                             vecsmall_append(ord, r));
     872             :     }
     873             :   }
     874       37443 :   setlg(R, k);
     875       37443 :   return gerepilecopy(ltop, R);
     876             : }
     877             : 
     878             : /* H is a normal subgroup, C is the quotient map G -->G/H,
     879             :  * S is a subgroup of G/H, and G is embedded in Sym(l)
     880             :  * Return all the subgroups K of G such that
     881             :  * S= K mod H and K inter H={1} */
     882             : static GEN
     883       29498 : liftsubgroup(GEN C, GEN H, GEN S)
     884             : {
     885       29498 :   pari_sp ltop = avma;
     886       29498 :   GEN V = trivialsubgroups();
     887       29498 :   GEN Sgen = grp_get_gen(S);
     888       29498 :   GEN Sord = grp_get_ord(S);
     889       29498 :   GEN Cgen = quo_get_gen(C);
     890       29498 :   long n = lg(Sgen), i;
     891       74382 :   for (i = 1; i < n; ++i)
     892             :   { /*loop over generators of S*/
     893       44884 :     GEN W = group_leftcoset(H, gel(Cgen, mael(Sgen, i, 1)));
     894       44884 :     V = liftlistsubgroups(V, W, Sord[i]);
     895             :   }
     896       29498 :   return gerepilecopy(ltop,V);
     897             : }
     898             : 
     899             : /* 1:A4 2:S4 0: other */
     900             : long
     901        4942 : group_isA4S4(GEN G)
     902             : {
     903        4942 :   GEN elt = grp_get_gen(G);
     904        4942 :   GEN ord = grp_get_ord(G);
     905        4942 :   long n = lg(ord);
     906        4942 :   if (n != 4 && n != 5) return 0;
     907        1330 :   if (ord[1]!=2 || ord[2]!=2 || ord[3]!=3) return 0;
     908          14 :   if (perm_commute(gel(elt,1),gel(elt,3))) return 0;
     909          14 :   if (n==4) return 1;
     910           7 :   if (ord[4]!=2) return 0;
     911           7 :   if (perm_commute(gel(elt,3),gel(elt,4))) return 0;
     912           7 :   return 2;
     913             : }
     914             : /* compute all the subgroups of a group G */
     915             : GEN
     916        6188 : group_subgroups(GEN G)
     917             : {
     918        6188 :   pari_sp ltop = avma;
     919             :   GEN p1, H, C, Q, M, sg1, sg2, sg3;
     920        6188 :   GEN gen = grp_get_gen(G);
     921        6188 :   GEN ord = grp_get_ord(G);
     922        6188 :   long lM, i, j, n = lg(gen);
     923        6188 :   if (n == 1) return trivialsubgroups();
     924        4942 :   if (group_isA4S4(G))
     925             :   {
     926          14 :     GEN s = gel(gen,1);       /*s = (1,2)(3,4) */
     927          14 :     GEN t = gel(gen,2);       /*t = (1,3)(2,4) */
     928          14 :     GEN st = perm_mul(s, t); /*st = (1,4)(2,3) */
     929          14 :     H = dicyclicgroup(s, t, 2, 2);
     930             :     /* sg3 is the list of subgroups intersecting only partially with H*/
     931          14 :     sg3 = cgetg((n==4)?4: 10, t_VEC);
     932          14 :     gel(sg3,1) = cyclicgroup(s, 2);
     933          14 :     gel(sg3,2) = cyclicgroup(t, 2);
     934          14 :     gel(sg3,3) = cyclicgroup(st, 2);
     935          14 :     if (n==5)
     936             :     {
     937           7 :       GEN u = gel(gen,3);
     938           7 :       GEN v = gel(gen,4), w, u2;
     939           7 :       if (zv_equal(perm_conj(u,s), t)) /*u=(2,3,4)*/
     940           7 :         u2 = perm_mul(u,u);
     941             :       else
     942             :       {
     943           0 :         u2 = u;
     944           0 :         u = perm_mul(u,u);
     945             :       }
     946           7 :       if (perm_order(v)==2)
     947             :       {
     948           7 :         if (!perm_commute(s,v)) /*v=(1,2)*/
     949             :         {
     950           0 :           v = perm_conj(u,v);
     951           0 :           if (!perm_commute(s,v)) v = perm_conj(u,v);
     952             :         }
     953           7 :         w = perm_mul(v,t); /*w=(1,4,2,3)*/
     954             :       }
     955             :       else
     956             :       {
     957           0 :         w = v;
     958           0 :         if (!zv_equal(perm_mul(w,w), s)) /*w=(1,4,2,3)*/
     959             :         {
     960           0 :           w = perm_conj(u,w);
     961           0 :           if (!zv_equal(perm_mul(w,w), s)) w = perm_conj(u,w);
     962             :         }
     963           0 :         v = perm_mul(w,t); /*v=(1,2)*/
     964             :       }
     965           7 :       gel(sg3,4) = dicyclicgroup(s,v,2,2);
     966           7 :       gel(sg3,5) = dicyclicgroup(t,perm_conj(u,v),2,2);
     967           7 :       gel(sg3,6) = dicyclicgroup(st,perm_conj(u2,v),2,2);
     968           7 :       gel(sg3,7) = dicyclicgroup(s,w,2,2);
     969           7 :       gel(sg3,8) = dicyclicgroup(t,perm_conj(u,w),2,2);
     970           7 :       gel(sg3,9) = dicyclicgroup(st,perm_conj(u2,w),2,2);
     971             :     }
     972             :   }
     973             :   else
     974             :   {
     975        4928 :     long osig = mael(factoru(ord[1]), 1, 1);
     976        4928 :     GEN sig = perm_pow(gel(gen,1), ord[1]/osig);
     977        4928 :     H = cyclicgroup(sig,osig);
     978        4928 :     sg3 = NULL;
     979             :   }
     980        4942 :   C = group_quotient(G,H);
     981        4935 :   Q = quotient_group(C,G);
     982        4935 :   M = group_subgroups(Q); lM = lg(M);
     983             :   /* sg1 is the list of subgroups containing H*/
     984        4928 :   sg1 = cgetg(lM, t_VEC);
     985        4928 :   for (i = 1; i < lM; ++i) gel(sg1,i) = quotient_subgroup_lift(C,H,gel(M,i));
     986             :   /*sg2 is a list of lists of subgroups not intersecting with H*/
     987        4928 :   sg2 = cgetg(lM, t_VEC);
     988             :   /* Loop over all subgroups of G/H */
     989        4928 :   for (j = 1; j < lM; ++j) gel(sg2,j) = liftsubgroup(C, H, gel(M,j));
     990        4928 :   p1 = gconcat(sg1, shallowconcat1(sg2));
     991        4928 :   if (sg3)
     992             :   {
     993          14 :     p1 = gconcat(p1, sg3);
     994          14 :     if (n==5) /*ensure that the D4 subgroups of S4 are in supersolvable format*/
     995          28 :       for(j = 3; j <= 5; j++)
     996             :       {
     997          21 :         GEN c = gmael(p1,j,1);
     998          21 :         if (!perm_commute(gel(c,1),gel(c,3)))
     999             :         {
    1000          14 :           if (perm_commute(gel(c,2),gel(c,3))) { swap(gel(c,1), gel(c,2)); }
    1001             :           else
    1002           7 :             perm_mul_inplace2(gel(c,2), gel(c,1));
    1003             :         }
    1004             :       }
    1005             :   }
    1006        4928 :   return gerepileupto(ltop,p1);
    1007             : }
    1008             : 
    1009             : /*return 1 if G is abelian, else 0*/
    1010             : long
    1011         854 : group_isabelian(GEN G)
    1012             : {
    1013         854 :   GEN g = grp_get_gen(G);
    1014         854 :   long i, j, n = lg(g);
    1015        1407 :   for(i=2; i<n; i++)
    1016        1477 :     for(j=1; j<i; j++)
    1017         924 :       if (!perm_commute(gel(g,i), gel(g,j))) return 0;
    1018         728 :   return 1;
    1019             : }
    1020             : 
    1021             : /*If G is abelian, return its HNF matrix*/
    1022             : GEN
    1023         329 : group_abelianHNF(GEN G, GEN S)
    1024             : {
    1025         329 :   GEN M, g = grp_get_gen(G), o = grp_get_ord(G);
    1026         329 :   long i, j, k, n = lg(g);
    1027         329 :   if (!group_isabelian(G)) return NULL;
    1028         259 :   if (n==1) return cgetg(1,t_MAT);
    1029         245 :   if (!S) S = group_elts(G, group_domain(G));
    1030         245 :   M = cgetg(n,t_MAT);
    1031         868 :   for(i=1; i<n; i++)
    1032             :   {
    1033         623 :     GEN P, C = cgetg(n,t_COL);
    1034         623 :     pari_sp av = avma;
    1035         623 :     gel(M,i) = C;
    1036         623 :     P = perm_pow(gel(g,i), o[i]);
    1037         903 :     for(j=1; j<lg(S); j++)
    1038         903 :       if (zv_equal(P, gel(S,j))) break;
    1039         623 :     avma = av;
    1040         623 :     if (j==lg(S)) pari_err_BUG("galoisisabelian [inconsistent group]");
    1041         623 :     j--;
    1042        1162 :     for(k=1; k<i; k++)
    1043             :     {
    1044         539 :       long q = j / o[k];
    1045         539 :       gel(C,k) = stoi(j - q*o[k]);
    1046         539 :       j = q;
    1047             :     }
    1048         623 :     gel(C,k) = stoi(o[i]);
    1049         623 :     for (k++; k<n; k++) gel(C,k) = gen_0;
    1050             :   }
    1051         245 :   return M;
    1052             : }
    1053             : 
    1054             : /*If G is abelian, return its abstract SNF matrix*/
    1055             : GEN
    1056         280 : group_abelianSNF(GEN G, GEN L)
    1057             : {
    1058         280 :   pari_sp ltop = avma;
    1059         280 :   GEN H = group_abelianHNF(G,L);
    1060         280 :   if (!H) return NULL;
    1061         210 :   return gerepileupto(ltop, smithclean( ZM_snf(H) ));
    1062             : }
    1063             : 
    1064             : GEN
    1065         189 : abelian_group(GEN v)
    1066             : {
    1067         189 :   long card = zv_prod(v), i, d = 1, l = lg(v);
    1068         189 :   GEN G = cgetg(3,t_VEC), gen = cgetg(l,t_VEC);
    1069         189 :   gel(G,1) = gen;
    1070         189 :   gel(G,2) = vecsmall_copy(v);
    1071         420 :   for(i=1; i<l; i++)
    1072             :   {
    1073         231 :     GEN p = cgetg(card+1, t_VECSMALL);
    1074         231 :     long o = v[i], u = d*(o-1), j, k, l;
    1075         231 :     gel(gen, i) = p;
    1076             :     /* The following loop is over-optimized. Remember that I wrote it for
    1077             :      * testpermutation. Something has survived... BA */
    1078         826 :     for(j=1;j<=card;)
    1079             :     {
    1080        1582 :       for(k=1;k<o;k++)
    1081        1218 :         for(l=1;l<=d; l++,j++) p[j] = j+d;
    1082         364 :       for (l=1; l<=d; l++,j++) p[j] = j-u;
    1083             :     }
    1084         231 :     d += u;
    1085             :   }
    1086         189 :   return G;
    1087             : }
    1088             : 
    1089             : /*return 1 if H is a normal subgroup of G*/
    1090             : long
    1091          56 : group_subgroup_isnormal(GEN G, GEN H)
    1092             : {
    1093          56 :   GEN g = grp_get_gen(G);
    1094          56 :   long i, n = lg(g);
    1095          56 :   if (lg(grp_get_gen(H)) > 1 && group_domain(G) != group_domain(H))
    1096           0 :     pari_err_DOMAIN("group_subgroup_isnormal","domain(H)","!=",
    1097             :                     strtoGENstr("domain(G)"), H);
    1098         126 :   for(i=1; i<n; i++)
    1099          91 :     if (!group_perm_normalize(H, gel(g,i))) return 0;
    1100          35 :   return 1;
    1101             : }
    1102             : 
    1103             : long
    1104           0 : groupelts_exponent(GEN elts)
    1105             : {
    1106           0 :   long i, n = lg(elts)-1, expo = 1;
    1107           0 :   for(i=1; i<=n; i++) expo = ulcm(expo, perm_order(gel(elts,i)));
    1108           0 :   return expo;
    1109             : }
    1110             : 
    1111             : GEN
    1112         693 : groupelts_center(GEN S)
    1113             : {
    1114         693 :   pari_sp ltop = avma;
    1115         693 :   long i, j, n = lg(S)-1, l = n;
    1116         693 :   GEN V, elts = zero_F2v(n+1);
    1117       24969 :   for(i=1; i<=n; i++)
    1118             :   {
    1119       24276 :     if (F2v_coeff(elts,i)) { l--;  continue; }
    1120      545895 :     for(j=1; j<=n; j++)
    1121      535724 :       if (!perm_commute(gel(S,i),gel(S,j)))
    1122             :       {
    1123       13587 :         F2v_set(elts,i);
    1124       13587 :         F2v_set(elts,j); l--; break;
    1125             :       }
    1126             :   }
    1127         693 :   V = cgetg(l+1,t_VEC);
    1128       24969 :   for (i=1, j=1; i<=n ;i++)
    1129       24276 :     if (!F2v_coeff(elts,i)) gel(V,j++) = vecsmall_copy(gel(S,i));
    1130         693 :   return gerepileupto(ltop,V);
    1131             : }
    1132             : 
    1133             : GEN
    1134        4249 : groupelts_conjclasses(GEN elts, long *pnbcl)
    1135             : {
    1136        4249 :   long i, j, cl = 0, n = lg(elts)-1;
    1137        4249 :   GEN c = const_vecsmall(n,0);
    1138        4249 :   pari_sp av = avma;
    1139       52787 :   for (i=1; i<=n; i++)
    1140             :   {
    1141       48538 :     GEN g = gel(elts,i);
    1142       48538 :     if (c[i]) continue;
    1143       34923 :     c[i] = ++cl;
    1144      486745 :     for(j=1; j<=n; j++)
    1145      451822 :       if (j != i)
    1146             :       {
    1147      416899 :         GEN h = perm_conj(gel(elts,j), g);
    1148      416899 :         long i2 = gen_search(elts,h,0,(void*)&vecsmall_lexcmp,&cmp_nodata);
    1149      416899 :         c[i2] = cl;
    1150      416899 :         avma = av;
    1151             :       }
    1152             :   }
    1153        4249 :   if (pnbcl) *pnbcl = cl;
    1154        4249 :   return c;
    1155             : }
    1156             : 
    1157             : GEN
    1158        4249 : conjclasses_repr(GEN conj, long nb)
    1159             : {
    1160        4249 :   long i, l = lg(conj);
    1161        4249 :   GEN e = const_vecsmall(nb, 0);
    1162       52787 :   for(i=1; i<l; i++)
    1163             :   {
    1164       48538 :     long ci = conj[i];
    1165       48538 :     if (!e[ci]) e[ci] = i;
    1166             :   }
    1167        4249 :   return e;
    1168             : }
    1169             : 
    1170             : /* elts of G sorted wrt vecsmall_lexcmp order: g in G is determined by g[1]
    1171             :  * so sort by increasing g[1] */
    1172             : static GEN
    1173        3864 : galois_elts_sorted(GEN gal)
    1174             : {
    1175             :   long i, l;
    1176        3864 :   GEN elts = gal_get_group(gal), v = cgetg_copy(elts, &l);
    1177        3864 :   for (i = 1; i < l; i++) { GEN g = gel(elts,i); gel(v, g[1]) = g; }
    1178        3864 :   return v;
    1179             : }
    1180             : GEN
    1181        4263 : group_to_cc(GEN G)
    1182             : {
    1183        4263 :   GEN elts = checkgroupelts(G), z = cgetg(5,t_VEC);
    1184        4249 :   long n, flag = 1;
    1185        4249 :   if (typ(gel(G,1)) == t_POL)
    1186        3864 :     elts = galois_elts_sorted(G); /* galoisinit */
    1187             :   else
    1188             :   {
    1189         385 :     long i, l = lg(elts);
    1190         385 :     elts = gen_sort(elts,(void*)vecsmall_lexcmp,cmp_nodata); /* general case */
    1191        5824 :     for (i = 1; i < l; i++)
    1192        5586 :       if (gel(elts,i)[1] != i) { flag = 0; break; }
    1193             :   }
    1194        4249 :   gel(z,1) = elts;
    1195        4249 :   gel(z,2) = groupelts_conjclasses(elts,&n);
    1196        4249 :   gel(z,3) = conjclasses_repr(gel(z,2),n);
    1197        4249 :   gel(z,4) = utoi(flag); return z;
    1198             : }
    1199             : 
    1200             : /* S a list of generators */
    1201             : GEN
    1202           0 : groupelts_abelian_group(GEN S)
    1203             : {
    1204           0 :   pari_sp ltop = avma;
    1205             :   GEN Qgen, Qord, Qelt;
    1206           0 :   long i, j, n = lg(gel(S,1))-1, l = lg(S);
    1207           0 :   Qord = cgetg(l, t_VECSMALL);
    1208           0 :   Qgen = cgetg(l, t_VEC);
    1209           0 :   Qelt = mkvec(identity_perm(n));
    1210           0 :   for (i = 1, j = 1; i < l; ++i)
    1211             :   {
    1212           0 :     GEN  g = gel(S,i);
    1213           0 :     long o = perm_relorder(g, groupelts_set(Qelt, n));
    1214           0 :     gel(Qgen,j) = g;
    1215           0 :     Qord[j] = o;
    1216           0 :     if (o != 1) { Qelt = perm_generate(g, Qelt, o); j++; }
    1217             :   }
    1218           0 :   setlg(Qgen,j);
    1219           0 :   setlg(Qord,j);
    1220           0 :   return gerepilecopy(ltop, mkvec2(Qgen, Qord));
    1221             : }
    1222             : 
    1223             : GEN
    1224          14 : group_export_GAP(GEN G)
    1225             : {
    1226          14 :   pari_sp av = avma;
    1227          14 :   GEN s, comma, g = grp_get_gen(G);
    1228          14 :   long i, k, l = lg(g);
    1229          14 :   if (l == 1) return strtoGENstr("Group(())");
    1230           7 :   s = cgetg(2*l, t_VEC);
    1231           7 :   comma = strtoGENstr(", ");
    1232           7 :   gel(s,1) = strtoGENstr("Group(");
    1233          28 :   for (i=1, k=2; i < l; ++i)
    1234             :   {
    1235          21 :     if (i > 1) gel(s,k++) = comma;
    1236          21 :     gel(s,k++) = perm_to_GAP(gel(g,i));
    1237             :   }
    1238           7 :   gel(s,k++) = strtoGENstr(")");
    1239           7 :   return gerepilecopy(av, shallowconcat1(s));
    1240             : }
    1241             : 
    1242             : GEN
    1243          14 : group_export_MAGMA(GEN G)
    1244             : {
    1245          14 :   pari_sp av = avma;
    1246          14 :   GEN s, comma, g = grp_get_gen(G);
    1247          14 :   long i, k, l = lg(g);
    1248          14 :   if (l == 1) return strtoGENstr("PermutationGroup<1|>");
    1249           7 :   s = cgetg(2*l, t_VEC);
    1250           7 :   comma = strtoGENstr(", ");
    1251           7 :   gel(s,1) = gsprintf("PermutationGroup<%ld|",group_domain(G));
    1252          28 :   for (i=1, k=2; i < l; ++i)
    1253             :   {
    1254          21 :     if (i > 1) gel(s,k++) = comma;
    1255          21 :     gel(s,k++) = GENtoGENstr( vecsmall_to_vec(gel(g,i)) );
    1256             :   }
    1257           7 :   gel(s,k++) = strtoGENstr(">");
    1258           7 :   return gerepilecopy(av, shallowconcat1(s));
    1259             : }
    1260             : 
    1261             : GEN
    1262          28 : group_export(GEN G, long format)
    1263             : {
    1264          28 :   switch(format)
    1265             :   {
    1266          14 :   case 0: return group_export_GAP(G);
    1267          14 :   case 1: return group_export_MAGMA(G);
    1268             :   }
    1269           0 :   pari_err_FLAG("galoisexport");
    1270           0 :   return NULL; /*-Wall*/
    1271             : }

Generated by: LCOV version 1.13