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 to exceed 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 - char.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.18.1 lcov report (development 30288-703288f808) Lines: 827 864 95.7 %
Date: 2025-05-19 09:23:07 Functions: 64 64 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000  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             : #include "pari.h"
      15             : #include "paripriv.h"
      16             : 
      17             : /*********************************************************************/
      18             : /**                                                                 **/
      19             : /**                  GENERIC ABELIAN CHARACTERS                     **/
      20             : /**                                                                 **/
      21             : /*********************************************************************/
      22             : /* check whether G is a znstar */
      23             : int
      24     3870833 : checkznstar_i(GEN G)
      25             : {
      26     3870616 :   return (typ(G) == t_VEC && lg(G) == 6
      27     3870238 :       && typ(znstar_get_N(G)) == t_INT
      28     3870231 :       && typ(znstar_get_faN(G)) == t_VEC
      29     7741449 :       && typ(gel(G,1)) == t_VEC && lg(gel(G,1)) == 3);
      30             : }
      31             : 
      32             : int
      33      121170 : char_check(GEN cyc, GEN chi)
      34      121170 : { return typ(chi) == t_VEC && lg(chi) == lg(cyc) && RgV_is_ZV(chi); }
      35             : 
      36             : /* Shallow; return [ d[1],  d[1]/d[2],...,d[1]/d[n] ] */
      37             : GEN
      38      254737 : cyc_normalize(GEN d)
      39             : {
      40      254737 :   long i, l = lg(d);
      41             :   GEN C, D;
      42      254737 :   if (l == 1) return mkvec(gen_1);
      43      254716 :   D = cgetg(l, t_VEC); gel(D,1) = C = gel(d,1);
      44      603099 :   for (i = 2; i < l; i++) gel(D,i) = diviiexact(C, gel(d,i));
      45      254716 :   return D;
      46             : }
      47             : 
      48             : /* chi character [D,C] given by chi(g_i) = \zeta_D^C[i] for all i, return
      49             :  * [d,c] such that chi(g_i) = \zeta_d^c[i] for all i and d minimal */
      50             : GEN
      51      296884 : char_simplify(GEN D, GEN C)
      52             : {
      53      296884 :   GEN d = D;
      54      296884 :   if (lg(C) == 1) d = gen_1;
      55             :   else
      56             :   {
      57      296324 :     GEN t = gcdii(d, ZV_content(C));
      58      296324 :     if (!equali1(t))
      59             :     {
      60      204729 :       long tc = typ(C);
      61      204729 :       C = ZC_Z_divexact(C, t); settyp(C, tc);
      62      204729 :       d = diviiexact(d, t);
      63             :     }
      64             :   }
      65      296884 :   return mkvec2(d,C);
      66             : }
      67             : 
      68             : /* Shallow; ncyc from cyc_normalize(): ncyc[1] = cyc[1],
      69             :  * ncyc[i] = cyc[i]/cyc[1] for i > 1; chi character on G ~ cyc.
      70             :  * Return [d,c] such that: chi( g_i ) = e(chi[i] / cyc[i]) = e(c[i]/ d) */
      71             : GEN
      72      295435 : char_normalize(GEN chi, GEN ncyc)
      73             : {
      74      295435 :   long i, l = lg(chi);
      75      295435 :   GEN c = cgetg(l, t_VEC);
      76      295435 :   if (l > 1) {
      77      295414 :     gel(c,1) = gel(chi,1);
      78      745661 :     for (i = 2; i < l; i++) gel(c,i) = mulii(gel(chi,i), gel(ncyc,i));
      79             :   }
      80      295435 :   return char_simplify(gel(ncyc,1), c);
      81             : }
      82             : 
      83             : /* \chi(gen[i]) = zeta_D^chic[i])
      84             :  * denormalize: express chi(gen[i]) in terms of zeta_{cyc[i]} */
      85             : GEN
      86      215754 : char_denormalize(GEN cyc, GEN D, GEN chic)
      87             : {
      88      215754 :   long i, l = lg(chic);
      89      215754 :   GEN chi = cgetg(l, t_VEC);
      90             :   /* \chi(gen[i]) = e(chic[i] / D) = e(chi[i] / cyc[i])
      91             :    * hence chi[i] = chic[i]cyc[i]/ D  mod cyc[i] */
      92      828709 :   for (i = 1; i < l; ++i)
      93             :   {
      94      612955 :     GEN di = gel(cyc, i), t = diviiexact(mulii(di, gel(chic,i)), D);
      95      612955 :     gel(chi, i) = modii(t, di);
      96             :   }
      97      215754 :   return chi;
      98             : }
      99             : 
     100             : /* Called by function 's'. x is a group object affording ".cyc" method, and
     101             :  * chi an abelian character. Return NULL if the group is (Z/nZ)^* [special
     102             :  * case more character types allowed] and x.cyc otherwise */
     103             : static GEN
     104         532 : get_cyc(GEN x, GEN chi, const char *s)
     105             : {
     106         532 :   switch(nftyp(x))
     107             :   {
     108         448 :     case typ_BIDZ:
     109         448 :       if (!zncharcheck(x, chi)) pari_err_TYPE(s, chi);
     110         448 :       return NULL;
     111           0 :     case typ_GCHAR:
     112           0 :       x = gchar_get_cyc(x);
     113           0 :       if (!is_vec_t(typ(chi)) || lg(chi) != lg(x) || !RgV_is_ZV(chi))
     114           0 :         pari_err_TYPE(s, chi); /* FIXME: handle norm component */
     115           0 :       return x;
     116          84 :     default:
     117          84 :       if (typ(x) != t_VEC || !RgV_is_ZV(x)) x = member_cyc(x);
     118          84 :       if (!char_check(x, chi)) pari_err_TYPE(s, chi);
     119          84 :       return x;
     120             :   }
     121             : }
     122             : 
     123             : /* conjugate character [ZV/ZC] */
     124             : GEN
     125        6664 : charconj(GEN cyc, GEN chi)
     126             : {
     127             :   long i, l;
     128        6664 :   GEN z = cgetg_copy(chi, &l);
     129       11998 :   for (i = 1; i < l; i++)
     130             :   {
     131        5334 :     GEN c = gel(chi,i);
     132        5334 :     gel(z,i) = signe(c)? subii(gel(cyc,i), c): gen_0;
     133             :   }
     134        6664 :   return z;
     135             : }
     136             : GEN
     137          28 : charconj0(GEN x, GEN chi)
     138             : {
     139          28 :   GEN cyc = get_cyc(x, chi, "charconj");
     140          28 :   return cyc? charconj(cyc, chi): zncharconj(x, chi);
     141             : }
     142             : 
     143             : GEN
     144     1383382 : charorder(GEN cyc, GEN x)
     145             : {
     146     1383382 :   pari_sp av = avma;
     147     1383382 :   long i, l = lg(cyc);
     148     1383382 :   GEN f = gen_1;
     149     3616928 :   for (i = 1; i < l; i++)
     150     2233546 :     if (signe(gel(x,i)))
     151             :     {
     152     1276240 :       GEN c, o = gel(cyc,i);
     153     1276240 :       if (!signe(o))
     154           0 :         return gc_upto(av,mkoo());
     155     1276240 :       c = gcdii(o, gel(x,i));
     156     1276240 :       if (!is_pm1(c)) o = diviiexact(o,c);
     157     1276240 :       f = lcmii(f, o);
     158             :     }
     159     1383382 :   return gc_INT(av, f);
     160             : }
     161             : GEN
     162         210 : charorder0(GEN x, GEN chi)
     163             : {
     164         210 :   GEN cyc = get_cyc(x, chi, "charorder");
     165         210 :   return cyc? charorder(cyc, chi): zncharorder(x, chi);
     166             : }
     167             : 
     168             : /* chi character of abelian G: chi[i] = chi(z_i), where G = \oplus Z/cyc[i] z_i.
     169             :  * Return Ker chi */
     170             : GEN
     171       99211 : charker(GEN cyc, GEN chi)
     172             : {
     173       99211 :   long i, l = lg(cyc);
     174             :   GEN nchi, ncyc, m, U;
     175             : 
     176       99211 :   if (l == 1) return cgetg(1,t_MAT); /* trivial subgroup */
     177       99162 :   ncyc = cyc_normalize(cyc);
     178       99162 :   nchi = char_normalize(chi, ncyc);
     179       99162 :   m = shallowconcat(gel(nchi,2), gel(nchi,1));
     180       99162 :   U = gel(ZV_extgcd(m), 2); setlg(U,l);
     181      272503 :   for (i = 1; i < l; i++) setlg(U[i], l);
     182       99162 :   return hnfmodid(U, cyc);
     183             : }
     184             : GEN
     185          35 : charker0(GEN x, GEN chi)
     186             : {
     187          35 :   GEN cyc = get_cyc(x, chi, "charker");
     188          35 :   return cyc? charker(cyc, chi): zncharker(x, chi);
     189             : }
     190             : 
     191             : GEN
     192         189 : charpow(GEN cyc, GEN a, GEN N)
     193             : {
     194             :   long i, l;
     195         189 :   GEN v = cgetg_copy(a, &l);
     196         350 :   for (i = 1; i < l; i++) gel(v,i) = Fp_mul(gel(a,i), N, gel(cyc,i));
     197         189 :   return v;
     198             : }
     199             : GEN
     200      302008 : charmul(GEN cyc, GEN a, GEN b)
     201             : {
     202             :   long i, l;
     203      302008 :   GEN v = cgetg_copy(a, &l);
     204     1350328 :   for (i = 1; i < l; i++) gel(v,i) = Fp_add(gel(a,i), gel(b,i), gel(cyc,i));
     205      302008 :   return v;
     206             : }
     207             : GEN
     208        5544 : chardiv(GEN cyc, GEN a, GEN b)
     209             : {
     210             :   long i, l;
     211        5544 :   GEN v = cgetg_copy(a, &l);
     212       11991 :   for (i = 1; i < l; i++) gel(v,i) = Fp_sub(gel(a,i), gel(b,i), gel(cyc,i));
     213        5544 :   return v;
     214             : }
     215             : GEN
     216          63 : charpow0(GEN x, GEN a, GEN N)
     217             : {
     218          63 :   GEN cyc = get_cyc(x, a, "charpow");
     219          63 :   return cyc? charpow(cyc, a, N): zncharpow(x, a, N);
     220             : }
     221             : GEN
     222         154 : charmul0(GEN x, GEN a, GEN b)
     223             : {
     224         154 :   const char *s = "charmul";
     225         154 :   GEN cyc = get_cyc(x, a, s);
     226         154 :   if (!cyc)
     227             :   {
     228         154 :     if (!zncharcheck(x, b)) pari_err_TYPE(s, b);
     229         154 :     return zncharmul(x, a, b);
     230             :   }
     231             :   else
     232             :   {
     233           0 :     if (!char_check(cyc, b)) pari_err_TYPE(s, b);
     234           0 :     return charmul(cyc, a, b);
     235             :   }
     236             : }
     237             : GEN
     238          42 : chardiv0(GEN x, GEN a, GEN b)
     239             : {
     240          42 :   const char *s = "chardiv";
     241          42 :   GEN cyc = get_cyc(x, a, s);
     242          42 :   if (!cyc)
     243             :   {
     244          42 :     if (!zncharcheck(x, b)) pari_err_TYPE(s, b);
     245          42 :     return znchardiv(x, a, b);
     246             :   }
     247             :   else
     248             :   {
     249           0 :     if (!char_check(cyc, b)) pari_err_TYPE(s, b);
     250           0 :     return chardiv(cyc, a, b);
     251             :   }
     252             : }
     253             : 
     254             : static GEN
     255     2520917 : chareval_i(GEN nchi, GEN dlog, GEN z)
     256             : {
     257     2520917 :   GEN o, q, r, b = gel(nchi,1);
     258     2520917 :   GEN a = FpV_dotproduct(gel(nchi,2), dlog, b);
     259             :   /* image is a/b in Q/Z */
     260     2520917 :   if (!z) return gdiv(a,b);
     261     2520539 :   if (typ(z) == t_INT)
     262             :   {
     263     2509143 :     q = dvmdii(z, b, &r);
     264     2509143 :     if (signe(r)) pari_err_TYPE("chareval", z);
     265     2509143 :     return mulii(a, q);
     266             :   }
     267             :   /* return z^(a*o/b), assuming z^o = 1 and b | o */
     268       11396 :   if (typ(z) != t_VEC || lg(z) != 3) pari_err_TYPE("chareval", z);
     269       11396 :   o = gel(z,2); if (typ(o) != t_INT) pari_err_TYPE("chareval", z);
     270       11396 :   q = dvmdii(o, b, &r); if (signe(r)) pari_err_TYPE("chareval", z);
     271       11396 :   q = mulii(a, q); /* in [0, o[ since a is reduced mod b */
     272       11396 :   z = gel(z,1);
     273       11396 :   if (typ(z) == t_VEC)
     274             :   {
     275       10080 :     if (itos_or_0(o) != lg(z)-1) pari_err_TYPE("chareval", z);
     276       10080 :     return gcopy(gel(z, itos(q)+1));
     277             :   }
     278             :   else
     279        1316 :     return gpow(z, q, DEFAULTPREC);
     280             : }
     281             : 
     282             : static GEN
     283        1855 : not_coprime(GEN z)
     284        1855 : { return (!z || typ(z) == t_INT)? gen_m1: gen_0; }
     285             : 
     286             : static GEN
     287          35 : get_chi(GEN cyc, GEN chi)
     288             : {
     289          35 :   if (!char_check(cyc,chi)) pari_err_TYPE("chareval", chi);
     290          35 :   return char_normalize(chi, cyc_normalize(cyc));
     291             : }
     292             : /* G a bnr.  FIXME: horribly inefficient to check that (x,N)=1, what to do ? */
     293             : static int
     294          42 : bnr_coprime(GEN G, GEN x)
     295             : {
     296          42 :   GEN t, N = gel(bnr_get_mod(G), 1);
     297          42 :   if (typ(x) == t_INT) /* shortcut */
     298             :   {
     299          14 :     t = gcdii(gcoeff(N,1,1), x);
     300          14 :     if (equali1(t)) return 1;
     301           0 :     t = idealadd(G, N, x);
     302           0 :     return equali1(gcoeff(t,1,1));
     303             :   }
     304          28 :   x = idealnumden(G, x);
     305          28 :   t = idealadd(G, N, gel(x,1));
     306          28 :   if (!equali1(gcoeff(t,1,1))) return 0;
     307          21 :   t = idealadd(G, N, gel(x,2));
     308          21 :   return equali1(gcoeff(t,1,1));
     309             : }
     310             : GEN
     311        3605 : chareval(GEN G, GEN chi, GEN x, GEN z)
     312             : {
     313        3605 :   pari_sp av = avma;
     314             :   GEN nchi, L;
     315             : 
     316        3605 :   switch(nftyp(G))
     317             :   {
     318          42 :     case typ_BNR:
     319          42 :       if (!bnr_coprime(G, x)) return not_coprime(z);
     320          28 :       L = isprincipalray(G, x);
     321          28 :       nchi = get_chi(bnr_get_cyc(G), chi);
     322          28 :       break;
     323           7 :     case typ_BNF:
     324           7 :       L = isprincipal(G, x);
     325           7 :       nchi = get_chi(bnf_get_cyc(G), chi);
     326           7 :       break;
     327        3542 :     case typ_BIDZ:
     328        3542 :       if (checkznstar_i(G)) return gc_upto(av, znchareval(G, chi, x, z));
     329             :       /* don't implement chars on general bid: need an nf... */
     330             :     case typ_GCHAR:
     331           7 :         pari_err_TYPE("chareval [use gchareval to evaluate a grossencharacter]", G);
     332           7 :     default:
     333           7 :       pari_err_TYPE("chareval", G);
     334             :       return NULL;/* LCOV_EXCL_LINE */
     335             :   }
     336          35 :   return gc_upto(av, chareval_i(nchi, L, z));
     337             : }
     338             : 
     339             : /* nchi = [ord,D] a quasi-normalized character (ord may be a multiple of
     340             :  * the character order); return v such that v[n] = -1 if (n,N) > 1 else
     341             :  * chi(n) = e(v[n]/ord), 1 <= n <= N */
     342             : GEN
     343       37331 : ncharvecexpo(GEN G, GEN nchi)
     344             : {
     345       37331 :   long N = itou(znstar_get_N(G)), ord = itou(gel(nchi,1)), i, j, l;
     346             :   GEN cyc, gen, d, t, t1, t2, t3, e, u, u1, u2, u3;
     347       37331 :   GEN D = gel(nchi,2), v = const_vecsmall(N,-1);
     348       37331 :   pari_sp av = avma;
     349       37331 :   if (typ(D) == t_COL) {
     350       37331 :     cyc = znstar_get_conreycyc(G);
     351       37331 :     gen = znstar_get_conreygen(G);
     352             :   } else {
     353           0 :     cyc = znstar_get_cyc(G);
     354           0 :     gen = znstar_get_gen(G);
     355             :   }
     356       37331 :   l = lg(cyc);
     357       37331 :   e = u = cgetg(N+1,t_VECSMALL);
     358       37331 :   d = t = cgetg(N+1,t_VECSMALL);
     359       37331 :   *++d = 1;
     360       37331 :   *++e = 0; v[*d] = *e;
     361       78715 :   for (i = 1; i < l; i++)
     362             :   {
     363       41384 :     ulong g = itou(gel(gen,i)), c = itou(gel(cyc,i)), x = itou(gel(D,i));
     364      551908 :     for (t1=t,u1=u,j=c-1; j; j--,t1=t2,u1=u2)
     365     1132880 :       for (t2=d,u2=e, t3=t1,u3=u1; t3<t2; )
     366             :       {
     367      622356 :         *++d = Fl_mul(*++t3, g, N);
     368      622356 :         *++e = Fl_add(*++u3, x, ord); v[*d] = *e;
     369             :       }
     370             :   }
     371       37331 :   set_avma(av); return v;
     372             : }
     373             : 
     374             : /*****************************************************************************/
     375             : 
     376             : static ulong
     377      257075 : lcmuu(ulong a, ulong b) { return (a/ugcd(a,b)) * b; }
     378             : static ulong
     379      101248 : zv_charorder(GEN cyc, GEN x)
     380             : {
     381      101248 :   long i, l = lg(cyc);
     382      101248 :   ulong f = 1;
     383      332115 :   for (i = 1; i < l; i++)
     384      230867 :     if (x[i])
     385             :     {
     386      155827 :       ulong o = cyc[i];
     387      155827 :       f = lcmuu(f, o / ugcd(o, x[i]));
     388             :     }
     389      101248 :   return f;
     390             : }
     391             : 
     392             : /* N > 0 */
     393             : GEN
     394      549969 : coprimes_zv(ulong N)
     395             : {
     396      549969 :   GEN v = const_vecsmall(N,1);
     397      549969 :   pari_sp av = avma;
     398      549969 :   GEN P = gel(factoru(N),1);
     399      549969 :   long i, l = lg(P);
     400     1282309 :   for (i = 1; i < l; i++)
     401             :   {
     402      732340 :     ulong p = P[i], j;
     403     3414453 :     for (j = p; j <= N; j += p) v[j] = 0;
     404             :   }
     405      549969 :   set_avma(av); return v;
     406             : }
     407             : /* cf zv_cyc_minimal: return k such that g*k is minimal (wrt lex) */
     408             : long
     409       46627 : zv_cyc_minimize(GEN cyc, GEN g, GEN coprime)
     410             : {
     411       46627 :   pari_sp av = avma;
     412       46627 :   long d, k, e, i, maxi, k0, bestk, l = lg(g), o = lg(coprime)-1;
     413             :   GEN best, gk, gd;
     414             :   ulong t;
     415       46627 :   if (o == 1) return 1;
     416       54187 :   for (i = 1; i < l; i++)
     417       54187 :     if (g[i]) break;
     418       46627 :   if (g[i] == 1) return 1;
     419       38423 :   k0 = Fl_invgen(g[i], cyc[i], &t);
     420       38423 :   d = cyc[i] / (long)t;
     421       38423 :   if (k0 > 1) g = vecmoduu(Flv_Fl_mul(g, k0, cyc[i]), cyc);
     422       50785 :   for (i++; i < l; i++)
     423       44751 :     if (g[i]) break;
     424       38423 :   if (i == l) return k0;
     425       32389 :   cyc = vecslice(cyc,i,l-1);
     426       32389 :   g   = vecslice(g,  i,l-1);
     427       32389 :   e = cyc[1];
     428       32389 :   gd = Flv_Fl_mul(g, d, e);
     429       32389 :   bestk = 1; best = g; maxi = e/ugcd(d,e);
     430       48643 :   for (gk = g, k = d+1, i = 1; i < maxi; k += d, i++)
     431             :   {
     432       16254 :     long ko = k % o;
     433       16254 :     gk = Flv_add(gk, gd, e); if (!ko || !coprime[ko]) continue;
     434        7357 :     gk = vecmoduu(gk, cyc);
     435        7357 :     if (vecsmall_lexcmp(gk, best) < 0) { best = gk; bestk = k; }
     436             :   }
     437       32389 :   return gc_long(av, bestk == 1? k0: (long) Fl_mul(k0, bestk, o));
     438             : }
     439             : /* g of order o in abelian group G attached to cyc. Is g a minimal generator
     440             :  * [wrt lex order] of the cyclic subgroup it generates;
     441             :  * coprime = coprimes_zv(o) */
     442             : long
     443      100800 : zv_cyc_minimal(GEN cyc, GEN g, GEN coprime)
     444             : {
     445      100800 :   pari_sp av = avma;
     446      100800 :   long i, maxi, d, k, e, l = lg(g), o = lg(coprime)-1; /* elt order */
     447             :   GEN gd, gk;
     448      100800 :   if (o == 1) return 1;
     449      106407 :   for (k = 1; k < l; k++)
     450      106407 :     if (g[k]) break;
     451      100800 :   if (g[k] == 1) return 1;
     452       80332 :   if (cyc[k] % g[k]) return 0;
     453       80332 :   d = cyc[k] / g[k]; /* > 1 */
     454       97069 :   for (k++; k < l; k++) /* skip following 0s */
     455       97069 :     if (g[k]) break;
     456       80332 :   if (k == l) return 1;
     457       80332 :   cyc = vecslice(cyc,k,l-1);
     458       80332 :   g   = vecslice(g,  k,l-1);
     459       80332 :   e = cyc[1];
     460             :   /* find k in (Z/e)^* such that g*k mod cyc is lexicographically minimal,
     461             :    * k = 1 mod d to fix the first nonzero entry */
     462       80332 :   gd = Flv_Fl_mul(g, d, e); maxi = e/ugcd(d,e);
     463      138159 :   for (gk = g, k = d+1, i = 1; i < maxi; i++, k += d)
     464             :   {
     465       67823 :     long ko = k % o;
     466       67823 :     gk = Flv_add(gk, gd, e); if (!coprime[ko]) continue;
     467       36141 :     gk = vecmoduu(gk, cyc);
     468       36141 :     if (vecsmall_lexcmp(gk, g) < 0) return gc_long(av,0);
     469             :   }
     470       70336 :   return gc_long(av,1);
     471             : }
     472             : 
     473             : static GEN
     474        8008 : coprime_tables(long N)
     475             : {
     476        8008 :   GEN D = divisorsu(N), v = const_vec(N, NULL);
     477        8008 :   long i, l = lg(D);
     478       47481 :   for (i = 1; i < l; i++) gel(v, D[i]) = coprimes_zv(D[i]);
     479        8008 :   return v;
     480             : }
     481             : /* enumerate all group elements, modulo (Z/cyc[1])^* */
     482             : static GEN
     483        8036 : cyc2elts_normal(GEN cyc, long maxord, GEN ORD)
     484             : {
     485        8036 :   long i, n, o, N, j = 1;
     486             :   GEN z, vcoprime;
     487             : 
     488        8036 :   if (typ(cyc) != t_VECSMALL) cyc = vec_to_vecsmall(cyc);
     489        8036 :   n = lg(cyc)-1;
     490        8036 :   N = zv_prod(cyc);
     491        8036 :   z = cgetg(N+1, t_VEC);
     492        8036 :   if (1 <= maxord && (!ORD|| zv_search(ORD,1)))
     493        7651 :     gel(z,j++) = zero_zv(n);
     494        8036 :   if (n == 0) { setlg(z, j); return z; }
     495        8008 :   vcoprime = coprime_tables(cyc[1]);
     496       23184 :   for (i = n; i > 0; i--)
     497             :   {
     498       15176 :     GEN cyc0 = vecslice(cyc,i+1,n), pre = zero_zv(i);
     499       15176 :     GEN D = divisorsu(cyc[i]), C = cyc2elts(cyc0);
     500       15176 :     long s, t, lD = lg(D), nC = lg(C)-1; /* remove last element */
     501       54810 :     for (s = 1; s < lD-1; s++)
     502             :     {
     503       39634 :       long o0 = D[lD-s]; /* cyc[i] / D[s] */
     504       39634 :       if (o0 > maxord) continue;
     505       38143 :       pre[i] = D[s];
     506       38143 :       if (!ORD || zv_search(ORD,o0))
     507             :       {
     508       37842 :         GEN c = vecsmall_concat(pre, zero_zv(n-i));
     509       37842 :         gel(z,j++) = c;
     510             :       }
     511      139391 :       for (t = 1; t < nC; t++)
     512             :       {
     513      101248 :         GEN chi0 = gel(C,t);
     514      101248 :         o = lcmuu(o0, zv_charorder(cyc0,chi0));
     515      101248 :         if (o <= maxord && (!ORD || zv_search(ORD,o)))
     516             :         {
     517      100800 :           GEN c = vecsmall_concat(pre, chi0);
     518      100800 :           if (zv_cyc_minimal(cyc, c, gel(vcoprime,o))) gel(z,j++) = c;
     519             :         }
     520             :       }
     521             :     }
     522             :   }
     523        8008 :   setlg(z,j); return z;
     524             : }
     525             : 
     526             : GEN
     527       10479 : chargalois(GEN G, GEN ORD)
     528             : {
     529       10479 :   pari_sp av = avma;
     530             :   long maxord, i, l;
     531       10479 :   GEN v, cyc = (typ(G) == t_VEC && RgV_is_ZVpos(G))? G: member_cyc(G);
     532       10479 :   if (lg(cyc) == 1 && !ORD) retmkvec(cgetg(1,t_VEC));
     533        8043 :   maxord = itou(cyc_get_expo(cyc));
     534        8043 :   if (ORD)
     535         441 :     switch(typ(ORD))
     536             :     {
     537             :       long l;
     538          42 :       case t_VEC:
     539          42 :         ORD = ZV_to_zv(ORD);
     540         406 :       case t_VECSMALL:
     541         406 :         l = lg(ORD);
     542         406 :         if (l > 2)
     543             :         {
     544         371 :           ORD = leafcopy(ORD);
     545         371 :           vecsmall_sort(ORD);
     546             :         }
     547         406 :         if (l == 1) retgc_const(av, cgetg(1, t_VEC));
     548         399 :         maxord = minss(maxord, ORD[l-1]);
     549         399 :         break;
     550          35 :       case t_INT:
     551          35 :         maxord = minss(maxord, itos(ORD));
     552          35 :         ORD = NULL;
     553          35 :         break;
     554           0 :       default: pari_err_TYPE("chargalois", ORD);
     555             :     }
     556        8036 :   v = cyc2elts_normal(cyc, maxord, ORD); l = lg(v);
     557      144333 :   for(i = 1; i < l; i++) gel(v,i) = zv_to_ZV(gel(v,i));
     558        8036 :   return gc_upto(av, v);
     559             : }
     560             : 
     561             : /*********************************************************************/
     562             : /**                                                                 **/
     563             : /**                  (Z/NZ)^* AND DIRICHLET CHARACTERS              **/
     564             : /**                                                                 **/
     565             : /*********************************************************************/
     566             : 
     567             : GEN
     568      131416 : znstar0(GEN N, long flag)
     569             : {
     570      131416 :   GEN F = NULL, P, E, cyc, gen, mod, G;
     571             :   long i, i0, l, nbprimes;
     572      131416 :   pari_sp av = avma;
     573             : 
     574      131416 :   if (flag && flag != 1) pari_err_FLAG("znstar");
     575      131416 :   if ((F = check_arith_all(N,"znstar")))
     576             :   {
     577       10073 :     F = clean_Z_factor(F);
     578       10073 :     N = typ(N) == t_VEC? gel(N,1): factorback(F);
     579             :   }
     580      131416 :   if (!signe(N))
     581             :   {
     582          21 :     if (flag) pari_err_IMPL("znstar(0,1)");
     583          14 :     set_avma(av);
     584          14 :     retmkvec3(gen_2, mkvec(gen_2), mkvec(gen_m1));
     585             :   }
     586      131395 :   N = absi_shallow(N);
     587      131395 :   if (abscmpiu(N,2) <= 0)
     588             :   {
     589       10626 :     G = mkvec3(gen_1, cgetg(1,t_VEC), cgetg(1,t_VEC));
     590       10626 :     if (flag)
     591             :     {
     592       10556 :       GEN v = const_vec(6,cgetg(1,t_VEC));
     593       10556 :       gel(v,3) = cgetg(1,t_MAT);
     594       21077 :       F = equali1(N)? mkvec2(cgetg(1,t_COL),cgetg(1,t_VECSMALL))
     595       10556 :                     : mkvec2(mkcol(gen_2), mkvecsmall(1));
     596       10556 :       G = mkvec5(mkvec2(N,mkvec(gen_0)), G, F, v, cgetg(1,t_MAT));
     597             :     }
     598       10626 :     return gc_GEN(av,G);
     599             :   }
     600      120769 :   if (!F) F = Z_factor(N);
     601      120769 :   P = gel(F,1); nbprimes = lg(P)-1;
     602      120769 :   E = ZV_to_nv( gel(F,2) );
     603      120769 :   switch(mod8(N))
     604             :   {
     605       20368 :     case 0:
     606       20368 :       P = shallowconcat(gen_2,P);
     607       20368 :       E = vecsmall_prepend(E, E[1]); /* add a copy of p=2 row */
     608       20368 :       i = 2; /* 2 generators at 2 */
     609       20368 :       break;
     610       19235 :     case 4:
     611       19235 :       i = 1; /* 1 generator at 2 */
     612       19235 :       break;
     613        7021 :     case 2: case 6:
     614        7021 :       P = vecsplice(P,1);
     615        7021 :       E = vecsplice(E,1); /* remove 2 */
     616        7021 :       i = 0; /* no generator at 2 */
     617        7021 :       break;
     618       74145 :     default:
     619       74145 :       i = 0; /* no generator at 2 */
     620       74145 :       break;
     621             :   }
     622      120769 :   l = lg(P);
     623      120769 :   cyc = cgetg(l,t_VEC);
     624      120769 :   gen = cgetg(l,t_VEC);
     625      120769 :   mod = cgetg(l,t_VEC);
     626             :   /* treat p=2 first */
     627      120769 :   if (i == 2)
     628             :   {
     629       20368 :     long v2 = E[1];
     630       20368 :     GEN q = int2n(v2);
     631       20368 :     gel(cyc,1) = gen_2;
     632       20368 :     gel(gen,1) = subiu(q,1); /* -1 */
     633       20368 :     gel(mod,1) = q;
     634       20368 :     gel(cyc,2) = int2n(v2-2);
     635       20368 :     gel(gen,2) = utoipos(5); /* Conrey normalization */
     636       20368 :     gel(mod,2) = q;
     637       20368 :     i0 = 3;
     638             :   }
     639      100401 :   else if (i == 1)
     640             :   {
     641       19235 :     gel(cyc,1) = gen_2;
     642       19235 :     gel(gen,1) = utoipos(3);
     643       19235 :     gel(mod,1) = utoipos(4);
     644       19235 :     i0 = 2;
     645             :   }
     646             :   else
     647       81166 :     i0 = 1;
     648             :   /* odd primes, fill remaining entries */
     649      291716 :   for (i = i0; i < l; i++)
     650             :   {
     651      170947 :     long e = E[i];
     652      170947 :     GEN p = gel(P,i), q = powiu(p, e-1), Q = mulii(p, q);
     653      170947 :     gel(cyc,i) = subii(Q, q); /* phi(p^e) */
     654      170947 :     gel(gen,i) = pgener_Zp(p);/* Conrey normalization, for e = 1 also */
     655      170947 :     gel(mod,i) = Q;
     656             :   }
     657             :   /* gen[i] has order cyc[i] and generates (Z/mod[i]Z)^* */
     658      120769 :   if (nbprimes > 1) /* lift generators to (Z/NZ)^*, = 1 mod N/mod[i] */
     659      243977 :     for (i=1; i<l; i++)
     660             :     {
     661      175663 :       GEN Q = gel(mod,i), g = gel(gen,i), qinv = Fp_inv(Q, diviiexact(N,Q));
     662      175663 :       g = addii(g, mulii(mulii(subsi(1,g),qinv),Q));
     663      175663 :       gel(gen,i) = modii(g, N);
     664             :     }
     665             : 
     666             :   /* cyc[i] > 1 and remain so in the loop, gen[i] = 1 mod (N/mod[i]) */
     667      120769 :   if (!flag)
     668             :   { /* update generators in place; about twice faster */
     669       36145 :     G = gen;
     670       36390 :     for (i=l-1; i>=2; i--)
     671             :     {
     672         245 :       GEN ci = gel(cyc,i), gi = gel(G,i);
     673             :       long j;
     674         588 :       for (j=i-1; j>=1; j--) /* we want cyc[i] | cyc[j] */
     675             :       {
     676         343 :         GEN cj = gel(cyc,j), gj, qj, v, d;
     677             : 
     678         343 :         d = bezout(ci,cj,NULL,&v); /* > 1 */
     679         518 :         if (absequalii(ci, d)) continue; /* ci | cj */
     680         210 :         if (absequalii(cj, d)) { /* cj | ci */
     681         175 :           swap(gel(G,j),gel(G,i));
     682         175 :           gi = gel(G,i);
     683         175 :           swap(gel(cyc,j),gel(cyc,i));
     684         175 :           ci = gel(cyc,i); continue;
     685             :         }
     686             : 
     687          35 :         qj = diviiexact(cj,d);
     688          35 :         gel(cyc,j) = mulii(ci,qj);
     689          35 :         gel(cyc,i) = d;
     690             : 
     691             :         /* [1,v*cj/d; 0,1]*[1,0;-1,1]*diag(cj,ci)*[ci/d,-v; cj/d,u]
     692             :          * = diag(lcm,gcd), with u ci + v cj = d */
     693          35 :         gj = gel(G,j);
     694             :         /* (gj, gi) *= [1,0; -1,1]^-1 */
     695          35 :         gj = Fp_mul(gj, gi, N); /* order ci*qj = lcm(ci,cj) */
     696             :         /* (gj,gi) *= [1,v*qj; 0,1]^-1 */
     697          35 :         togglesign_safe(&v);
     698          35 :         if (signe(v) < 0) v = modii(v,ci); /* >= 0 to avoid inversions */
     699          35 :         gel(G,i) = gi = Fp_mul(gi, Fp_pow(gj, mulii(qj, v), N), N);
     700          35 :         gel(G,j) = gj;
     701          35 :         ci = d; if (absequaliu(ci, 2)) break;
     702             :       }
     703             :     }
     704       36145 :     G = mkvec3(ZV_prod(cyc), cyc, FpV_to_mod(G,N));
     705             :   }
     706             :   else
     707             :   { /* keep matrices between generators, return an 'init' structure */
     708       84624 :     GEN D, U, Ui, fao = cgetg(l, t_VEC), lo = cgetg(l, t_VEC);
     709       84624 :     F = mkvec2(P, E);
     710       84624 :     D = ZV_snf_group(cyc,&U,&Ui);
     711      279152 :     for (i = 1; i < l; i++)
     712             :     {
     713      194528 :       GEN t = gen_0, p = gel(P,i), p_1 = subiu(p,1);
     714      194528 :       long e = E[i];
     715      194528 :       gel(fao,i) = get_arith_ZZM(p_1);
     716      194528 :       if (e >= 2 && !absequaliu(p,2))
     717             :       {
     718       12944 :         GEN q = gel(mod,i), g = Fp_pow(gel(gen,i),p_1,q);
     719       12944 :         if (e == 2)
     720       10564 :           t = Fp_inv(diviiexact(subiu(g,1), p), p);
     721             :         else
     722        2380 :           t = ginv(Qp_log(cvtop(g,p,e)));
     723             :       }
     724      194528 :       gel(lo,i) = t;
     725             :     }
     726       84624 :     G = cgetg(l, t_VEC);
     727      279152 :     for (i = 1; i < l; i++) gel(G,i) = FpV_factorback(gen, gel(Ui,i), N);
     728       84624 :     G = mkvec3(ZV_prod(D), D, G);
     729       84624 :     G = mkvec5(mkvec2(N,mkvec(gen_0)), G, F,
     730             :                mkvecn(6,mod,fao,Ui,gen,cyc,lo), U);
     731             :   }
     732      120769 :   return gc_GEN(av, G);
     733             : }
     734             : GEN
     735       36138 : znstar(GEN N) { return znstar0(N, 0); }
     736             : 
     737             : /* g has order 2^(e-2), g,h = 1 (mod 4); return x s.t. g^x = h (mod 2^e) */
     738             : static GEN
     739      949576 : Zideallog_2k(GEN h, GEN g, long e, GEN pe)
     740             : {
     741      949576 :   GEN a = Fp_log(h, g, int2n(e-2), pe);
     742      949576 :   if (typ(a) != t_INT) return NULL;
     743      949576 :   return a;
     744             : }
     745             : 
     746             : /* ord = get_arith_ZZM(p-1), simplified form of znlog_rec: g is known
     747             :  * to be a primitive root mod p^e; lo = 1/log_p(g^(p-1)) */
     748             : static GEN
     749     2000680 : Zideallog_pk(GEN h, GEN g, GEN p, long e, GEN pe, GEN ord, GEN lo)
     750             : {
     751     2000680 :   GEN gp = (e == 1)? g: modii(g, p);
     752     2000680 :   GEN hp = (e == 1)? h: modii(h, p);
     753     2000680 :   GEN a = Fp_log(hp, gp, ord, p);
     754     2000680 :   if (typ(a) != t_INT) return NULL;
     755     2000673 :   if (e > 1)
     756             :   { /* find a s.t. g^a = h (mod p^e), p odd prime, e > 0, (h,p) = 1 */
     757             :     /* use p-adic log: O(log p + e) mul*/
     758       53572 :     GEN b, p_1 = gel(ord,1);
     759       53572 :     h = Fp_mul(h, Fp_pow(g, negi(a), pe), pe);
     760             :     /* g,h = 1 mod p; compute b s.t. h = g^b */
     761       53572 :     if (e == 2) /* simpler */
     762       46145 :       b = Fp_mul(diviiexact(subiu(h,1), p), lo, p);
     763             :     else
     764        7427 :       b = padic_to_Q(gmul(Qp_log(cvtop(h, p, e)), lo));
     765       53572 :     a = addii(a, mulii(p_1, b));
     766             :   }
     767     2000673 :   return a;
     768             : }
     769             : 
     770             : int
     771      903777 : znconrey_check(GEN cyc, GEN chi)
     772      903777 : { return typ(chi) == t_COL && lg(chi) == lg(cyc) && RgV_is_ZV(chi); }
     773             : 
     774             : int
     775      207207 : zncharcheck(GEN G, GEN chi)
     776             : {
     777      207207 :   switch(typ(chi))
     778             :   {
     779         875 :     case t_INT: return 1;
     780      197834 :     case t_COL: return znconrey_check(znstar_get_conreycyc(G), chi);
     781        8498 :     case t_VEC: return char_check(znstar_get_cyc(G), chi);
     782             :   }
     783           0 :   return 0;
     784             : }
     785             : 
     786             : GEN
     787      153055 : znconreyfromchar_normalized(GEN bid, GEN chi)
     788             : {
     789      153055 :   GEN nchi, U = znstar_get_U(bid);
     790      153055 :   long l = lg(chi);
     791      153055 :   if (l == 1) retmkvec2(gen_1,cgetg(1,t_VEC));
     792      150612 :   if (!RgV_is_ZV(chi) || lgcols(U) != l) pari_err_TYPE("lfunchiZ", chi);
     793      150605 :   nchi = char_normalize(chi, cyc_normalize(znstar_get_cyc(bid)));
     794      150605 :   gel(nchi,2) = ZV_ZM_mul(gel(nchi,2),U); return nchi;
     795             : }
     796             : 
     797             : GEN
     798      142919 : znconreyfromchar(GEN bid, GEN chi)
     799             : {
     800      142919 :   GEN nchi = znconreyfromchar_normalized(bid, chi);
     801      142912 :   GEN v = char_denormalize(znstar_get_conreycyc(bid), gel(nchi,1), gel(nchi,2));
     802      142912 :   settyp(v, t_COL); return v;
     803             : }
     804             : 
     805             : /* discrete log on canonical "primitive root" generators
     806             :  * Allow log(x) instead of x [usual discrete log on bid's generators] */
     807             : GEN
     808     2606528 : znconreylog(GEN bid, GEN x)
     809             : {
     810     2606528 :   pari_sp av = avma;
     811             :   GEN N, L, F, P,E, y, pe, fao, gen, lo, cycg;
     812             :   long i, l;
     813     2606528 :   if (!checkznstar_i(bid)) pari_err_TYPE("znconreylog", bid);
     814     2606521 :   N = znstar_get_N(bid);
     815     2606521 :   if (abscmpiu(N, 2) <= 0)
     816             :   {
     817       12390 :     switch(typ(x))
     818             :     {
     819        2688 :     case t_INT: break;
     820           7 :     case t_INTMOD:
     821           7 :       if (!equalii(N, gel(x,1))) pari_err_TYPE("znconreylog", x);
     822           0 :       x = gel(x,2); break;
     823        9688 :     case t_COL:
     824             :     case t_VEC:
     825        9688 :       if (lg(x) != 1) pari_err_TYPE("znconreylog", x);
     826        9674 :       break;
     827           7 :     default: pari_err_TYPE("znconreylog", x);
     828             :     }
     829       12362 :     return cgetg(1, t_COL);
     830             :   }
     831     2594131 :   cycg = znstar_get_conreycyc(bid);
     832     2594131 :   switch(typ(x))
     833             :   {
     834             :     GEN Ui;
     835           0 :     case t_INTMOD:
     836           0 :       if (!equalii(N, gel(x,1))) pari_err_TYPE("znconreylog", x);
     837           0 :       x = gel(x,2); /* fall through */
     838     2560013 :     case t_INT:
     839     2560013 :       if (!signe(x)) pari_err_COPRIME("znconreylog", x, N);
     840     2560006 :       break;
     841          35 :     case t_COL: /* log_bid(x) */
     842          35 :       Ui = znstar_get_Ui(bid);
     843          35 :       if (!RgV_is_ZV(x) || lg(x) != lg(Ui)) pari_err_TYPE("znconreylog", x);
     844          35 :       return gc_upto(av, ZV_ZV_mod(ZM_ZC_mul(Ui,x), cycg));
     845       34083 :     case t_VEC:
     846       34083 :       return gc_GEN(av, znconreyfromchar(bid, x));
     847           0 :     default: pari_err_TYPE("znconreylog", x);
     848             :   }
     849     2560006 :   F = znstar_get_faN(bid); /* factor(N) */
     850     2560006 :   P = gel(F, 1); /* prime divisors of N */
     851     2560006 :   E = gel(F, 2); /* exponents */
     852     2560006 :   L = gel(bid,4);
     853     2560006 :   pe = znstar_get_pe(bid);
     854     2560006 :   fao = gel(L,2);
     855     2560006 :   gen = znstar_get_conreygen(bid); /* local generators of (Z/p^k)^* */
     856     2560006 :   lo = gel(L,6); /* 1/log_p((g_i)^(p_i-1)) */
     857             : 
     858     2560006 :   l = lg(gen); i = 1;
     859     2560006 :   y = cgetg(l, t_COL);
     860     2560006 :   if (!mod2(N) && !mod2(x)) pari_err_COPRIME("znconreylog", x, N);
     861     2559992 :   if (absequaliu(gel(P,1), 2) && E[1] >= 2)
     862             :   {
     863     1259759 :     if (E[1] == 2)
     864      310183 :       gel(y,i++) = mod4(x) == 1? gen_0: gen_1;
     865             :     else
     866             :     {
     867      949576 :       GEN a, x2, q2 = gel(pe,1);
     868      949576 :       x2 = modii(x, q2);
     869      949576 :       if (mod4(x) == 1) /* 1 or 5 mod 8*/
     870      592362 :         gel(y,i++) = gen_0;
     871             :       else /* 3 or 7 */
     872      357214 :       { gel(y,i++) = gen_1; x2 = subii(q2, x2); }
     873             :       /* x2 = 5^x mod q */
     874      949576 :       a = Zideallog_2k(x2, gel(gen,i), E[1], q2);
     875      949576 :       if (!a) pari_err_COPRIME("znconreylog", x, N);
     876      949576 :       gel(y, i++) = a;
     877             :     }
     878             :   }
     879     4560665 :   while (i < l)
     880             :   {
     881     2000680 :     GEN p = gel(P,i), q = gel(pe,i), xpe = modii(x, q);
     882     2000680 :     GEN a = Zideallog_pk(xpe, gel(gen,i), p, E[i], q, gel(fao,i), gel(lo,i));
     883     2000680 :     if (!a) pari_err_COPRIME("znconreylog", x, N);
     884     2000673 :     gel(y, i++) = a;
     885             :   }
     886     2559985 :   return gc_GEN(av, y);
     887             : }
     888             : GEN
     889       23024 : Zideallog(GEN bid, GEN x)
     890             : {
     891       23024 :   pari_sp av = avma;
     892       23024 :   GEN y = znconreylog(bid, x), U = znstar_get_U(bid);
     893       22996 :   return gc_upto(av, ZM_ZC_mul(U, y));
     894             : }
     895             : GEN
     896         294 : znlog0(GEN h, GEN g, GEN o)
     897             : {
     898         294 :   if (typ(g) == t_VEC)
     899             :   {
     900             :     GEN N;
     901          56 :     if (o) pari_err_TYPE("znlog [with znstar]", o);
     902          56 :     if (!checkznstar_i(g)) pari_err_TYPE("znlog", g);
     903          56 :     N = znstar_get_N(g);
     904          56 :     h = Rg_to_Fp(h,N);
     905          49 :     return Zideallog(g, h);
     906             :   }
     907         238 :   return znlog(h, g, o);
     908             : }
     909             : 
     910             : GEN
     911      260792 : znconreyexp(GEN bid, GEN x)
     912             : {
     913      260792 :   pari_sp av = avma;
     914             :   long i, l;
     915             :   GEN N, pe, gen, cycg, v, vmod;
     916             :   int e2;
     917      260792 :   if (!checkznstar_i(bid)) pari_err_TYPE("znconreyexp", bid);
     918      260792 :   cycg = znstar_get_conreycyc(bid);
     919      260792 :   switch(typ(x))
     920             :   {
     921          21 :     case t_VEC:
     922          21 :       x = znconreylog(bid, x);
     923          21 :       break;
     924      260771 :     case t_COL:
     925      260771 :       if (RgV_is_ZV(x) && lg(x) == lg(cycg)) break;
     926           7 :     default: pari_err_TYPE("znconreyexp",x);
     927             :   }
     928      260785 :   pe = znstar_get_pe(bid);
     929      260785 :   gen = znstar_get_conreygen(bid); /* local generators of (Z/p^k)^* */
     930      260785 :   cycg = znstar_get_conreycyc(bid);
     931      260785 :   l = lg(x); v = cgetg(l, t_VEC);
     932      260785 :   N = znstar_get_N(bid);
     933      260785 :   e2 = !mod8(N); /* 2 generators at p = 2 */
     934      805707 :   for (i = 1; i < l; i++)
     935             :   {
     936             :     GEN q, g, m;
     937      544922 :     if (i == 1 && e2) { gel(v,1) = NULL; continue; }
     938      512638 :     q = gel(pe,i);
     939      512638 :     g = gel(gen,i);
     940      512638 :     m = modii(gel(x,i), gel(cycg,i));
     941      512638 :     m = Fp_pow(g, m, q);
     942      512638 :     if (i == 2 && e2 && signe(gel(x,1))) m = Fp_neg(m, q);
     943      512638 :     gel(v,i) = mkintmod(m, q);
     944             :   }
     945      260785 :   if (e2) v = vecsplice(v, 1);
     946      260785 :   v = chinese1_coprime_Z(v);
     947      260785 :   vmod = gel(v,1);
     948      260785 :   v = gel(v,2);
     949      260785 :   if (mpodd(v) || mpodd(N)) return gc_GEN(av, v);
     950             :   /* handle N = 2 mod 4 */
     951         231 :   return gc_INT(av, addii(v, vmod));
     952             : }
     953             : 
     954             : /* Return Dirichlet character \chi_q(m,.), where bid = znstar(q);
     955             :  * m is either a t_INT, or a t_COL [Conrey logarithm] */
     956             : GEN
     957       46767 : znconreychar(GEN bid, GEN m)
     958             : {
     959       46767 :   pari_sp av = avma;
     960             :   GEN c, d, nchi;
     961             : 
     962       46767 :   if (!checkznstar_i(bid)) pari_err_TYPE("znconreychar", bid);
     963       46760 :   switch(typ(m))
     964             :   {
     965           0 :     case t_INTMOD:
     966           0 :       if (!equalii(gel(m,1), znstar_get_N(bid)))
     967           0 :         pari_err_TYPE("znconreychar",m);
     968           0 :       m = gel(m,2); /* fall through */
     969       46760 :     case t_INT:
     970             :     case t_COL:
     971       46760 :       nchi = znconrey_normalized(bid,m); /* images of primroot gens */
     972       46753 :       break;
     973           0 :     default:
     974           0 :       pari_err_TYPE("znconreychar",m);
     975             :       return NULL;/*LCOV_EXCL_LINE*/
     976             :   }
     977       46753 :   d = gel(nchi,1);
     978       46753 :   c = ZV_ZM_mul(gel(nchi,2), znstar_get_Ui(bid)); /* images of bid gens */
     979       46753 :   return gc_GEN(av, char_denormalize(znstar_get_cyc(bid),d,c));
     980             : }
     981             : 
     982             : /* chi a t_INT or Conrey log describing a character. Return conductor, as an
     983             :  * integer if primitive; as a t_VEC [N,factor(N)] if not. Set *pm=m to the
     984             :  * attached primitive character: chi(g_i) = m[i]/ord(g_i)
     985             :  * Caller should use znconreylog_normalize(BID, m), once BID(conductor) is
     986             :  * computed (wasteful to do it here since BID is shared by many characters) */
     987             : GEN
     988      736484 : znconreyconductor(GEN bid, GEN chi, GEN *pm)
     989             : {
     990      736484 :   pari_sp av = avma;
     991             :   GEN q, m, F, P, E;
     992             :   long i, j, l;
     993      736484 :   int e2, primitive = 1;
     994             : 
     995      736484 :   if (!checkznstar_i(bid)) pari_err_TYPE("znconreyconductor", bid);
     996      736484 :   if (typ(chi) == t_COL)
     997             :   {
     998      705943 :     if (!znconrey_check(znstar_get_conreycyc(bid), chi))
     999           0 :       pari_err_TYPE("znconreyconductor",chi);
    1000             :   }
    1001             :   else
    1002       30541 :     chi = znconreylog(bid, chi);
    1003      736477 :   l = lg(chi);
    1004      736477 :   F = znstar_get_faN(bid);
    1005      736477 :   P = gel(F,1);
    1006      736477 :   E = gel(F,2);
    1007      736477 :   if (l == 1)
    1008             :   {
    1009      105266 :     set_avma(av);
    1010      105266 :     if (pm) *pm = cgetg(1,t_COL);
    1011      105266 :     if (lg(P) == 1) return gen_1;
    1012          14 :     retmkvec2(gen_1, trivial_fact());
    1013             :   }
    1014      631211 :   P = leafcopy(P);
    1015      631211 :   E = leafcopy(E);
    1016      631211 :   m = cgetg(l, t_COL);
    1017      631211 :   e2 = (E[1] >= 3 && absequaliu(gel(P,1),2));
    1018      631211 :   i = j = 1;
    1019      631211 :   if (e2)
    1020             :   { /* two generators at p=2 */
    1021      286153 :     GEN a1 = gel(chi,1), a = gel(chi,2);
    1022      286153 :     i = 3;
    1023      286153 :     if (!signe(a))
    1024             :     {
    1025       98014 :       e2 =  primitive = 0;
    1026       98014 :       if (signe(a1))
    1027             :       { /* lose one generator */
    1028       45507 :         E[1] = 2;
    1029       45507 :         gel(m,1) = a1;
    1030       45507 :         j = 2;
    1031             :       }
    1032             :       /* else lose both */
    1033             :     }
    1034             :     else
    1035             :     {
    1036      188139 :       long v = Z_pvalrem(a, gen_2, &a);
    1037      188139 :       if (v) { E[1] -= v; E[2] = E[1]; primitive = 0; }
    1038      188139 :       gel(m,1) = a1;
    1039      188139 :       gel(m,2) = a;
    1040      188139 :       j = 3;
    1041             :     }
    1042             :   }
    1043      631211 :   l = lg(P);
    1044     1829590 :   for (; i < l; i++)
    1045             :   {
    1046     1198379 :     GEN p = gel(P,i), a = gel(chi,i);
    1047             :     /* image of g_i in Q/Z is a/cycg[i], cycg[i] = order(g_i) */
    1048     1198379 :     if (!signe(a)) primitive = 0;
    1049             :     else
    1050             :     {
    1051      937818 :       long v = Z_pvalrem(a, p, &a);
    1052      937818 :       E[j] = E[i]; if (v) { E[j] -= v; primitive = 0; }
    1053      937818 :       gel(P,j) = gel(P,i);
    1054      937818 :       gel(m,j) = a; j++;
    1055             :     }
    1056             :   }
    1057      631211 :   setlg(m,j);
    1058      631211 :   setlg(P,j);
    1059      631211 :   setlg(E,j);
    1060      631211 :   if (pm) *pm = m; /* attached primitive  character */
    1061      631211 :   if (primitive)
    1062             :   {
    1063      295099 :     q = znstar_get_N(bid);
    1064      295099 :     if (mod4(q) == 2) primitive = 0;
    1065             :   }
    1066      631211 :   if (!primitive)
    1067             :   {
    1068      336889 :     if (e2)
    1069             :     { /* remove duplicate p=2 row from factorization */
    1070      114744 :       P = vecsplice(P,1);
    1071      114744 :       E = vecsplice(E,1);
    1072             :     }
    1073      336889 :     E = zc_to_ZC(E);
    1074      336889 :     q = mkvec2(factorback2(P,E), mkmat2(P,E));
    1075             :   }
    1076      631211 :   return gc_all(av, pm? 2: 1, &q, pm);
    1077             : }
    1078             : 
    1079             : GEN
    1080        8477 : zncharinduce(GEN G, GEN chi, GEN N)
    1081             : {
    1082        8477 :   pari_sp av = avma;
    1083             :   GEN q, faq, P, E, Pq, Eq, CHI;
    1084             :   long i, j, l;
    1085             :   int e2;
    1086             : 
    1087        8477 :   if (!checkznstar_i(G)) pari_err_TYPE("zncharinduce", G);
    1088        8477 :   if (!zncharcheck(G, chi)) pari_err_TYPE("zncharinduce", chi);
    1089        8477 :   q = znstar_get_N(G);
    1090        8477 :   if (typ(chi) != t_COL) chi = znconreylog(G, chi);
    1091        8477 :   if (checkznstar_i(N))
    1092             :   {
    1093        8288 :     GEN faN = znstar_get_faN(N);
    1094        8288 :     P = gel(faN,1); l = lg(P);
    1095        8288 :     E = gel(faN,2);
    1096        8288 :     N = znstar_get_N(N);
    1097        8288 :     if (l > 2 && equalii(gel(P,1),gel(P,2)))
    1098             :     { /* remove duplicate 2 */
    1099        2583 :       l--;
    1100        2583 :       P = vecsplice(P,1);
    1101        2583 :       E = vecsplice(E,1);
    1102             :     }
    1103             :   }
    1104             :   else
    1105             :   {
    1106         189 :     GEN faN = check_arith_pos(N, "zncharinduce");
    1107         189 :     if (!faN) faN = Z_factor(N);
    1108             :     else
    1109           0 :       N = (typ(N) == t_VEC)? gel(N,1): factorback(faN);
    1110         189 :     P = gel(faN,1);
    1111         189 :     E = gel(faN,2);
    1112             :   }
    1113        8477 :   if (!dvdii(N,q)) pari_err_DOMAIN("zncharinduce", "N % q", "!=", gen_0, N);
    1114        8470 :   if (mod4(N) == 2)
    1115             :   { /* remove 2 */
    1116          77 :     if (lg(P) > 1 && absequaliu(gel(P,1), 2))
    1117             :     {
    1118          35 :       P = vecsplice(P,1);
    1119          35 :       E = vecsplice(E,1);
    1120             :     }
    1121          77 :     N = shifti(N,-1);
    1122             :   }
    1123        8470 :   l = lg(P);
    1124             :   /* q = N or q = 2N, N odd */
    1125        8470 :   if (cmpii(N,q) <= 0) return gc_GEN(av, chi);
    1126             :   /* N > 1 => l > 1*/
    1127        8344 :   if (typ(E) != t_VECSMALL) E = ZV_to_zv(E);
    1128        8344 :   e2 = (E[1] >= 3 && absequaliu(gel(P,1),2)); /* 2 generators at 2 mod N */
    1129        8344 :   if (ZV_equal0(chi))
    1130             :   {
    1131        5446 :     set_avma(av);
    1132        5446 :     return equali1(N)? cgetg(1, t_COL): zerocol(l+e2 - 1);
    1133             :   }
    1134             : 
    1135        2898 :   faq = znstar_get_faN(G);
    1136        2898 :   Pq = gel(faq,1);
    1137        2898 :   Eq = gel(faq,2);
    1138        2898 :   CHI = cgetg(l+e2, t_COL);
    1139        2898 :   i = j = 1;
    1140        2898 :   if (e2)
    1141             :   {
    1142        1183 :     i = 2; j = 3;
    1143        1183 :     if (absequaliu(gel(Pq,1), 2))
    1144             :     {
    1145        1015 :       if (Eq[1] >= 3)
    1146             :       { /* 2 generators at 2 mod q */
    1147         560 :         gel(CHI,1) = gel(chi,1);
    1148         560 :         gel(CHI,2) = shifti(gel(chi,2), E[1]-Eq[1]);
    1149             :       }
    1150         455 :       else if (Eq[1] == 2)
    1151             :       { /* 1 generator at 2 mod q */
    1152         455 :         gel(CHI,1) = gel(chi,1);
    1153         455 :         gel(CHI,2) = gen_0;
    1154             :       }
    1155             :       else
    1156           0 :         gel(CHI,1) = gel(CHI,2) = gen_0;
    1157             :     }
    1158             :     else
    1159         168 :       gel(CHI,1) = gel(CHI,2) = gen_0;
    1160             :   }
    1161        7315 :   for (; i < l; i++,j++)
    1162             :   {
    1163        4417 :     GEN p = gel(P,i);
    1164        4417 :     long k = ZV_search(Pq, p);
    1165        4417 :     gel(CHI,j) = k? mulii(gel(chi,k), powiu(p, E[i]-Eq[k])): gen_0;
    1166             :   }
    1167        2898 :   return gc_GEN(av, CHI);
    1168             : }
    1169             : 
    1170             : /* m a Conrey log [on the canonical primitive roots], cycg the primitive
    1171             :  * roots orders */
    1172             : GEN
    1173     2594340 : znconreylog_normalize(GEN G, GEN m)
    1174             : {
    1175     2594340 :   GEN cycg = znstar_get_conreycyc(G);
    1176             :   long i, l;
    1177     2594340 :   GEN d, M = cgetg_copy(m, &l);
    1178     2594340 :   if (typ(cycg) != t_VEC || lg(cycg) != l)
    1179           0 :     pari_err_TYPE("znconreylog_normalize",mkvec2(m,cycg));
    1180     6881798 :   for (i = 1; i < l; i++) gel(M,i) = gdiv(gel(m,i), gel(cycg,i));
    1181             :   /* m[i]: image of primroot generators g_i in Q/Z */
    1182     2594340 :   M = Q_remove_denom(M, &d);
    1183     2594340 :   return mkvec2(d? d: gen_1, M);
    1184             : }
    1185             : 
    1186             : /* return normalized character on Conrey generators attached to chi: Conrey
    1187             :  * label (t_INT), char on (SNF) G.gen* (t_VEC), or Conrey log (t_COL) */
    1188             : GEN
    1189     2587515 : znconrey_normalized(GEN G, GEN chi)
    1190             : {
    1191     2587515 :   switch(typ(chi))
    1192             :   {
    1193         420 :     case t_INT: /* Conrey label */
    1194         420 :       return znconreylog_normalize(G, znconreylog(G, chi));
    1195     2576959 :     case t_COL: /* Conrey log */
    1196     2576959 :       if (!RgV_is_ZV(chi)) break;
    1197     2576959 :       return znconreylog_normalize(G, chi);
    1198       10136 :     case t_VEC: /* char on G.gen */
    1199       10136 :       if (!RgV_is_ZV(chi)) break;
    1200       10136 :       return znconreyfromchar_normalized(G, chi);
    1201             :   }
    1202           0 :   pari_err_TYPE("znchareval",chi);
    1203             :   return NULL;/* LCOV_EXCL_LINE */
    1204             : }
    1205             : 
    1206             : /* return 1 iff chi(-1) = -1, and 0 otherwise */
    1207             : long
    1208      197498 : zncharisodd(GEN G, GEN chi)
    1209             : {
    1210             :   long i, l, s;
    1211             :   GEN N;
    1212      197498 :   if (!checkznstar_i(G)) pari_err_TYPE("zncharisodd", G);
    1213      197498 :   if (!zncharcheck(G, chi)) pari_err_TYPE("zncharisodd", chi);
    1214      197498 :   if (typ(chi) != t_COL) chi = znconreylog(G, chi);
    1215      197498 :   N = znstar_get_N(G);
    1216      197498 :   l = lg(chi);
    1217      197498 :   s = 0;
    1218      197498 :   if (!mod8(N))
    1219             :   {
    1220       91371 :     s = mpodd(gel(chi,1));
    1221       91371 :     i = 3;
    1222             :   }
    1223             :   else
    1224      106127 :     i = 1;
    1225      537698 :   for (; i < l; i++) s += mpodd(gel(chi,i));
    1226      197498 :   return odd(s);
    1227             : }
    1228             : 
    1229             : GEN
    1230         847 : znchartokronecker(GEN G, GEN chi, long flag)
    1231             : {
    1232         847 :   pari_sp av = avma;
    1233             :   long s;
    1234             :   GEN F, o;
    1235             : 
    1236         847 :   if (flag && flag != 1) pari_err_FLAG("znchartokronecker");
    1237         847 :   s = zncharisodd(G, chi)? -1: 1;
    1238         847 :   if (typ(chi) != t_COL) chi = znconreylog(G, chi);
    1239         847 :   o = zncharorder(G, chi);
    1240         847 :   if (abscmpiu(o,2) > 0) { set_avma(av); return gen_0; }
    1241         581 :   F = znconreyconductor(G, chi, NULL);
    1242         581 :   if (typ(F) == t_INT)
    1243             :   {
    1244         469 :     if (s < 0) F = negi(F);
    1245         469 :     return gc_INT(av, F);
    1246             :   }
    1247         112 :   F = gel(F,1);
    1248         112 :   F = (s < 0)? negi(F): icopy(F);
    1249         112 :   if (!flag)
    1250             :   {
    1251          49 :     GEN MF = znstar_get_faN(G), P = gel(MF,1);
    1252          49 :     long i, l = lg(P);
    1253         140 :     for (i = 1; i < l; i++)
    1254             :     {
    1255          91 :       GEN p = gel(P,i);
    1256          91 :       if (!dvdii(F,p)) F = mulii(F,sqri(p));
    1257             :     }
    1258             :   }
    1259         112 :   return gc_INT(av, F);
    1260             : }
    1261             : 
    1262             : /* (D/.) as a character mod N; assume |D| divides N and D = 0,1 mod 4*/
    1263             : GEN
    1264      303331 : znchar_quad(GEN G, GEN D)
    1265             : {
    1266      303331 :   GEN cyc = znstar_get_conreycyc(G);
    1267      303331 :   GEN gen = znstar_get_conreygen(G);
    1268      303331 :   long i, l = lg(cyc);
    1269      303331 :   GEN chi = cgetg(l, t_COL);
    1270     1351686 :   for (i = 1; i < l; i++)
    1271             :   {
    1272     1048355 :     long k = kronecker(D, gel(gen,i));
    1273     1048355 :     gel(chi,i) = (k==1)? gen_0: shifti(gel(cyc,i), -1);
    1274             :   }
    1275      303331 :   return chi;
    1276             : }
    1277             : 
    1278             : GEN
    1279        3304 : znchar(GEN D)
    1280             : {
    1281        3304 :   pari_sp av = avma;
    1282             :   GEN G, chi;
    1283        3304 :   switch(typ(D))
    1284             :   {
    1285        2513 :     case t_INT:
    1286        2513 :       if (!signe(D) || Mod4(D) > 1) pari_err_TYPE("znchar", D);
    1287        2492 :       G = znstar0(D,1);
    1288        2492 :       chi = mkvec2(G, znchar_quad(G,D));
    1289        2492 :       break;
    1290         721 :     case t_INTMOD:
    1291         721 :       G = znstar0(gel(D,1), 1);
    1292         721 :       chi = mkvec2(G, znconreylog(G, gel(D,2)));
    1293         721 :       break;
    1294          56 :     case t_VEC:
    1295          56 :       if (checkMF_i(D)) { chi = vecslice(MF_get_CHI(D),1,2); break; }
    1296          49 :       else if (checkmf_i(D)) { chi = vecslice(mf_get_CHI(D),1,2); break; }
    1297          42 :       if (lg(D) != 3) pari_err_TYPE("znchar", D);
    1298          35 :       G = gel(D,1);
    1299          35 :       if (!checkznstar_i(G)) pari_err_TYPE("znchar", D);
    1300          28 :       chi = gel(D,2);
    1301          28 :       if (typ(chi) == t_VEC && lg(chi) == 3 && is_vec_t(typ(gel(chi,2))))
    1302             :       { /* normalized character */
    1303           7 :         GEN n = gel(chi,1), chic = gel(chi,2);
    1304           7 :         GEN cyc = typ(chic)==t_VEC? znstar_get_cyc(G): znstar_get_conreycyc(G);
    1305           7 :         if (!char_check(cyc, chic)) pari_err_TYPE("znchar",D);
    1306           7 :         chi = char_denormalize(cyc, n, chic);
    1307             :       }
    1308          28 :       if (!zncharcheck(G, chi)) pari_err_TYPE("znchar", D);
    1309          21 :       chi = mkvec2(G,chi); break;
    1310          14 :     default:
    1311          14 :       pari_err_TYPE("znchar", D);
    1312             :       return NULL; /*LCOV_EXCL_LINE*/
    1313             :   }
    1314        3248 :   return gc_GEN(av, chi);
    1315             : }
    1316             : 
    1317             : /* G a znstar, not stack clean */
    1318             : GEN
    1319     2522723 : znchareval(GEN G, GEN chi, GEN n, GEN z)
    1320             : {
    1321     2522723 :   GEN nchi, N = znstar_get_N(G);
    1322             :   /* avoid division by 0 */
    1323     2522723 :   if (typ(n) == t_FRAC && !equali1(gcdii(gel(n,2), N))) return not_coprime(z);
    1324     2522716 :   n = Rg_to_Fp(n, N);
    1325     2522716 :   if (!equali1(gcdii(n, N))) return not_coprime(z);
    1326             :   /* nchi: normalized character on Conrey generators */
    1327     2520882 :   nchi = znconrey_normalized(G, chi);
    1328     2520882 :   return chareval_i(nchi, znconreylog(G,n), z);
    1329             : }
    1330             : 
    1331             : /* G is a znstar, chi a Dirichlet character */
    1332             : GEN
    1333        5810 : zncharconj(GEN G, GEN chi)
    1334             : {
    1335        5810 :   switch(typ(chi))
    1336             :   {
    1337           7 :     case t_INT: chi = znconreylog(G, chi); /* fall through */
    1338         665 :     case t_COL: return charconj(znstar_get_conreycyc(G), chi);
    1339        5145 :     case t_VEC: return charconj(znstar_get_cyc(G), chi);
    1340             :   }
    1341           0 :   pari_err_TYPE("zncharconj",chi);
    1342             :   return NULL; /*LCOV_EXCL_LINE*/
    1343             : }
    1344             : 
    1345             : /* G is a znstar, chi a Dirichlet character */
    1346             : GEN
    1347      388661 : zncharorder(GEN G,  GEN chi)
    1348             : {
    1349      388661 :   switch(typ(chi))
    1350             :   {
    1351          21 :     case t_INT: chi = znconreylog(G, chi); /*fall through*/
    1352      383201 :     case t_COL: return charorder(znstar_get_conreycyc(G), chi);
    1353        5460 :     case t_VEC: return charorder(znstar_get_cyc(G), chi);
    1354           0 :     default: pari_err_TYPE("zncharorder",chi);
    1355             :              return NULL; /* LCOV_EXCL_LINE */
    1356             :   }
    1357             : }
    1358             : 
    1359             : /* G is a znstar, chi a Dirichlet character */
    1360             : GEN
    1361          21 : zncharker(GEN G, GEN chi)
    1362             : {
    1363          21 :   if (typ(chi) != t_VEC) chi = znconreychar(G, chi);
    1364          21 :   return charker(znstar_get_cyc(G), chi);
    1365             : }
    1366             : 
    1367             : /* G is a znstar, 'a' is a Dirichlet character */
    1368             : GEN
    1369         210 : zncharpow(GEN G, GEN a, GEN n)
    1370             : {
    1371         210 :   switch(typ(a))
    1372             :   {
    1373          21 :     case t_INT: return Fp_pow(a, n, znstar_get_N(G));
    1374          21 :     case t_VEC: return charpow(znstar_get_cyc(G), a, n);
    1375         168 :     case t_COL: return charpow(znstar_get_conreycyc(G), a, n);
    1376           0 :     default: pari_err_TYPE("znchapow",a);
    1377             :              return NULL; /* LCOV_EXCL_LINE */
    1378             :   }
    1379             : }
    1380             : /* G is a znstar, 'a' and 'b' are Dirichlet character */
    1381             : GEN
    1382      302015 : zncharmul(GEN G, GEN a, GEN b)
    1383             : {
    1384      302015 :   long ta = typ(a), tb = typ(b);
    1385      302015 :   if (ta == tb) switch(ta)
    1386             :   {
    1387           7 :     case t_INT: return Fp_mul(a, b, znstar_get_N(G));
    1388           7 :     case t_VEC: return charmul(znstar_get_cyc(G), a, b);
    1389      301980 :     case t_COL: return charmul(znstar_get_conreycyc(G), a, b);
    1390           0 :     default: pari_err_TYPE("zncharmul",a);
    1391             :              return NULL; /* LCOV_EXCL_LINE */
    1392             :   }
    1393          21 :   if (ta != t_COL) a = znconreylog(G, a);
    1394          21 :   if (tb != t_COL) b = znconreylog(G, b);
    1395          21 :   return charmul(znstar_get_conreycyc(G), a, b);
    1396             : }
    1397             : 
    1398             : /* G is a znstar, 'a' and 'b' are Dirichlet character */
    1399             : GEN
    1400        5551 : znchardiv(GEN G, GEN a, GEN b)
    1401             : {
    1402        5551 :   long ta = typ(a), tb = typ(b);
    1403        5551 :   if (ta == tb) switch(ta)
    1404             :   {
    1405           7 :     case t_INT: return Fp_div(a, b, znstar_get_N(G));
    1406           7 :     case t_VEC: return chardiv(znstar_get_cyc(G), a, b);
    1407        5516 :     case t_COL: return chardiv(znstar_get_conreycyc(G), a, b);
    1408           0 :     default: pari_err_TYPE("znchardiv",a);
    1409             :              return NULL; /* LCOV_EXCL_LINE */
    1410             :   }
    1411          21 :   if (ta != t_COL) a = znconreylog(G, a);
    1412          21 :   if (tb != t_COL) b = znconreylog(G, b);
    1413          21 :   return chardiv(znstar_get_conreycyc(G), a, b);
    1414             : }
    1415             : 
    1416             : /* CHI mod N = \prod_p p^e; let CHI = \prod CHI_p, CHI_p mod p^e
    1417             :  * return \prod_{p | (Q,N)} CHI_p. E.g if Q = p, return chi_p */
    1418             : GEN
    1419         791 : znchardecompose(GEN G, GEN chi, GEN Q)
    1420             : {
    1421             :   GEN c, P, E, F;
    1422             :   long l, lP, i;
    1423             : 
    1424         791 :   if (!checkznstar_i(G)) pari_err_TYPE("znchardecompose", G);
    1425         791 :   if (typ(Q) != t_INT) pari_err_TYPE("znchardecompose", Q);
    1426         791 :   if (typ(chi) == t_COL)
    1427         560 :   { if (!zncharcheck(G, chi)) pari_err_TYPE("znchardecompose", chi); }
    1428             :   else
    1429         231 :     chi = znconreylog(G, chi);
    1430         791 :   l = lg(chi); if (l == 1) return cgetg(1, t_VEC);
    1431         784 :   F = znstar_get_faN(G);
    1432         784 :   c = zerocol(l-1);
    1433         784 :   P = gel(F,1); /* prime divisors of N */
    1434         784 :   lP = lg(P);
    1435         784 :   E = gel(F,2); /* exponents */
    1436        2471 :   for (i = 1; i < lP; i++)
    1437             :   {
    1438        1687 :     GEN p = gel(P,i);
    1439        1687 :     if (i == 1 && equaliu(p,2) && E[1] >= 3)
    1440             :     {
    1441         567 :       if (!mpodd(Q))
    1442             :       {
    1443         203 :         gel(c,1) = icopy(gel(chi,1));
    1444         203 :         gel(c,2) = icopy(gel(chi,2));
    1445             :       }
    1446         567 :       i = 2; /* skip P[2] = P[1] = 2 */
    1447             :     }
    1448             :     else
    1449        1120 :       if (dvdii(Q, p)) gel(c,i) = icopy(gel(chi,i));
    1450             :   }
    1451         784 :   return c;
    1452             : }
    1453             : 
    1454             : GEN
    1455       24325 : zncharconductor(GEN G, GEN chi)
    1456             : {
    1457       24325 :   pari_sp av = avma;
    1458       24325 :   GEN F = znconreyconductor(G, chi, NULL);
    1459       24325 :   if (typ(F) == t_INT) return F;
    1460       12873 :   return gc_GEN(av, gel(F,1));
    1461             : }
    1462             : GEN
    1463        5558 : znchartoprimitive(GEN G, GEN chi)
    1464             : {
    1465        5558 :   pari_sp av = avma;
    1466        5558 :   GEN chi0, F = znconreyconductor(G, chi, &chi0);
    1467        5558 :   if (typ(F) == t_INT)
    1468        5299 :     chi = mkvec2(G,chi);
    1469             :   else
    1470         259 :     chi = mkvec2(znstar0(F,1), chi0);
    1471        5558 :   return gc_GEN(av, chi);
    1472             : }

Generated by: LCOV version 1.16