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 - base4.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 22296-eb24562) Lines: 1400 1552 90.2 %
Date: 2018-04-19 06:16:12 Functions: 140 154 90.9 %
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. 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             : /*******************************************************************/
      15             : /*                                                                 */
      16             : /*                       BASIC NF OPERATIONS                       */
      17             : /*                           (continued)                           */
      18             : /*                                                                 */
      19             : /*******************************************************************/
      20             : #include "pari.h"
      21             : #include "paripriv.h"
      22             : 
      23             : /*******************************************************************/
      24             : /*                                                                 */
      25             : /*                     IDEAL OPERATIONS                            */
      26             : /*                                                                 */
      27             : /*******************************************************************/
      28             : 
      29             : /* A valid ideal is either principal (valid nf_element), or prime, or a matrix
      30             :  * on the integer basis in HNF.
      31             :  * A prime ideal is of the form [p,a,e,f,b], where the ideal is p.Z_K+a.Z_K,
      32             :  * p is a rational prime, a belongs to Z_K, e=e(P/p), f=f(P/p), and b
      33             :  * is Lenstra's constant, such that p.P^(-1)= p Z_K + b Z_K.
      34             :  *
      35             :  * An extended ideal is a couple [I,F] where I is a valid ideal and F is
      36             :  * either an algebraic number, or a factorization matrix attached to an
      37             :  * algebraic number. All routines work with either extended ideals or ideals
      38             :  * (an omitted F is assumed to be [;] <-> 1).
      39             :  * All ideals are output in HNF form. */
      40             : 
      41             : /* types and conversions */
      42             : 
      43             : long
      44     3135734 : idealtyp(GEN *ideal, GEN *arch)
      45             : {
      46     3135734 :   GEN x = *ideal;
      47     3135734 :   long t,lx,tx = typ(x);
      48             : 
      49     3135734 :   if (tx==t_VEC && lg(x)==3)
      50      324897 :   { *arch = gel(x,2); x = gel(x,1); tx = typ(x); }
      51             :   else
      52     2810837 :     *arch = NULL;
      53     3135734 :   switch(tx)
      54             :   {
      55     1577172 :     case t_MAT: lx = lg(x);
      56     1577172 :       if (lx == 1) { t = id_PRINCIPAL; x = gen_0; break; }
      57     1577095 :       if (lx != lgcols(x)) pari_err_TYPE("idealtyp [non-square t_MAT]",x);
      58     1577088 :       t = id_MAT;
      59     1577088 :       break;
      60             : 
      61     1157454 :     case t_VEC: if (lg(x)!=6) pari_err_TYPE("idealtyp",x);
      62     1157440 :       t = id_PRIME; break;
      63             : 
      64             :     case t_POL: case t_POLMOD: case t_COL:
      65             :     case t_INT: case t_FRAC:
      66      401108 :       t = id_PRINCIPAL; break;
      67             :     default:
      68           0 :       pari_err_TYPE("idealtyp",x);
      69             :       return 0; /*LCOV_EXCL_LINE*/
      70             :   }
      71     3135713 :   *ideal = x; return t;
      72             : }
      73             : 
      74             : /* true nf; v = [a,x,...], a in Z. Return (a,x) */
      75             : GEN
      76      108490 : idealhnf_two(GEN nf, GEN v)
      77             : {
      78      108490 :   GEN p = gel(v,1), pi = gel(v,2), m = zk_scalar_or_multable(nf, pi);
      79      108490 :   if (typ(m) == t_INT) return scalarmat(gcdii(m,p), nf_get_degree(nf));
      80       94287 :   return ZM_hnfmodid(m, p);
      81             : }
      82             : /* true nf */
      83             : GEN
      84     1167695 : pr_hnf(GEN nf, GEN pr)
      85             : {
      86     1167695 :   GEN p = pr_get_p(pr), m;
      87     1167695 :   if (pr_is_inert(pr)) return scalarmat(p, nf_get_degree(nf));
      88      910249 :   m = zk_scalar_or_multable(nf, pr_get_gen(pr));
      89      910249 :   return ZM_hnfmodprime(m, p);
      90             : }
      91             : 
      92             : GEN
      93      270427 : idealhnf_principal(GEN nf, GEN x)
      94             : {
      95             :   GEN cx;
      96      270427 :   x = nf_to_scalar_or_basis(nf, x);
      97      270427 :   switch(typ(x))
      98             :   {
      99      154380 :     case t_COL: break;
     100       92066 :     case t_INT:  if (!signe(x)) return cgetg(1,t_MAT);
     101       91653 :       return scalarmat(absi(x), nf_get_degree(nf));
     102             :     case t_FRAC:
     103       23981 :       return scalarmat(Q_abs_shallow(x), nf_get_degree(nf));
     104           0 :     default: pari_err_TYPE("idealhnf",x);
     105             :   }
     106      154380 :   x = Q_primitive_part(x, &cx);
     107      154380 :   RgV_check_ZV(x, "idealhnf");
     108      154380 :   x = zk_multable(nf, x);
     109      154380 :   x = ZM_hnfmodid(x, zkmultable_capZ(x));
     110      154380 :   return cx? ZM_Q_mul(x,cx): x;
     111             : }
     112             : 
     113             : /* x integral ideal in t_MAT form, nx columns */
     114             : static GEN
     115           7 : vec_mulid(GEN nf, GEN x, long nx, long N)
     116             : {
     117           7 :   GEN m = cgetg(nx*N + 1, t_MAT);
     118             :   long i, j, k;
     119          21 :   for (i=k=1; i<=nx; i++)
     120          14 :     for (j=1; j<=N; j++) gel(m, k++) = zk_ei_mul(nf, gel(x,i),j);
     121           7 :   return m;
     122             : }
     123             : /* true nf */
     124             : GEN
     125      323206 : idealhnf_shallow(GEN nf, GEN x)
     126             : {
     127      323206 :   long tx = typ(x), lx = lg(x), N;
     128             : 
     129             :   /* cannot use idealtyp because here we allow non-square matrices */
     130      323206 :   if (tx == t_VEC && lx == 3) { x = gel(x,1); tx = typ(x); lx = lg(x); }
     131      323206 :   if (tx == t_VEC && lx == 6) return pr_hnf(nf,x); /* PRIME */
     132      223704 :   switch(tx)
     133             :   {
     134             :     case t_MAT:
     135             :     {
     136             :       GEN cx;
     137       47334 :       long nx = lx-1;
     138       47334 :       N = nf_get_degree(nf);
     139       47334 :       if (nx == 0) return cgetg(1, t_MAT);
     140       47313 :       if (nbrows(x) != N) pari_err_TYPE("idealhnf [wrong dimension]",x);
     141       47306 :       if (nx == 1) return idealhnf_principal(nf, gel(x,1));
     142             : 
     143       46109 :       if (nx == N && RgM_is_ZM(x) && ZM_ishnf(x)) return x;
     144       22435 :       x = Q_primitive_part(x, &cx);
     145       22435 :       if (nx < N) x = vec_mulid(nf, x, nx, N);
     146       22435 :       x = ZM_hnfmod(x, ZM_detmult(x));
     147       22435 :       return cx? ZM_Q_mul(x,cx): x;
     148             :     }
     149             :     case t_QFI:
     150             :     case t_QFR:
     151             :     {
     152          14 :       pari_sp av = avma;
     153          14 :       GEN u, D = nf_get_disc(nf), T = nf_get_pol(nf), f = nf_get_index(nf);
     154          14 :       GEN A = gel(x,1), B = gel(x,2);
     155          14 :       N = nf_get_degree(nf);
     156          14 :       if (N != 2)
     157           0 :         pari_err_TYPE("idealhnf [Qfb for non-quadratic fields]", x);
     158          14 :       if (!equalii(qfb_disc(x), D))
     159           7 :         pari_err_DOMAIN("idealhnf [Qfb]", "disc(q)", "!=", D, x);
     160             :       /* x -> A Z + (-B + sqrt(D)) / 2 Z
     161             :          K = Q[t]/T(t), t^2 + ut + v = 0,  u^2 - 4v = Df^2
     162             :          => t = (-u + sqrt(D) f)/2
     163             :          => sqrt(D)/2 = (t + u/2)/f */
     164           7 :       u = gel(T,3);
     165           7 :       B = deg1pol_shallow(ginv(f),
     166             :                           gsub(gdiv(u, shifti(f,1)), gdiv(B,gen_2)),
     167           7 :                           varn(T));
     168           7 :       return gerepileupto(av, idealhnf_two(nf, mkvec2(A,B)));
     169             :     }
     170      176356 :     default: return idealhnf_principal(nf, x); /* PRINCIPAL */
     171             :   }
     172             : }
     173             : GEN
     174        3157 : idealhnf(GEN nf, GEN x)
     175             : {
     176        3157 :   pari_sp av = avma;
     177        3157 :   GEN y = idealhnf_shallow(checknf(nf), x);
     178        3143 :   return (avma == av)? gcopy(y): gerepileupto(av, y);
     179             : }
     180             : 
     181             : /* GP functions */
     182             : 
     183             : GEN
     184          63 : idealtwoelt0(GEN nf, GEN x, GEN a)
     185             : {
     186          63 :   if (!a) return idealtwoelt(nf,x);
     187          42 :   return idealtwoelt2(nf,x,a);
     188             : }
     189             : 
     190             : GEN
     191          42 : idealpow0(GEN nf, GEN x, GEN n, long flag)
     192             : {
     193          42 :   if (flag) return idealpowred(nf,x,n);
     194          35 :   return idealpow(nf,x,n);
     195             : }
     196             : 
     197             : GEN
     198          56 : idealmul0(GEN nf, GEN x, GEN y, long flag)
     199             : {
     200          56 :   if (flag) return idealmulred(nf,x,y);
     201          49 :   return idealmul(nf,x,y);
     202             : }
     203             : 
     204             : GEN
     205          49 : idealdiv0(GEN nf, GEN x, GEN y, long flag)
     206             : {
     207          49 :   switch(flag)
     208             :   {
     209          21 :     case 0: return idealdiv(nf,x,y);
     210          28 :     case 1: return idealdivexact(nf,x,y);
     211           0 :     default: pari_err_FLAG("idealdiv");
     212             :   }
     213             :   return NULL; /* LCOV_EXCL_LINE */
     214             : }
     215             : 
     216             : GEN
     217          70 : idealaddtoone0(GEN nf, GEN arg1, GEN arg2)
     218             : {
     219          70 :   if (!arg2) return idealaddmultoone(nf,arg1);
     220          35 :   return idealaddtoone(nf,arg1,arg2);
     221             : }
     222             : 
     223             : /* b not a scalar */
     224             : static GEN
     225          28 : hnf_Z_ZC(GEN nf, GEN a, GEN b) { return hnfmodid(zk_multable(nf,b), a); }
     226             : /* b not a scalar */
     227             : static GEN
     228          21 : hnf_Z_QC(GEN nf, GEN a, GEN b)
     229             : {
     230             :   GEN db;
     231          21 :   b = Q_remove_denom(b, &db);
     232          21 :   if (db) a = mulii(a, db);
     233          21 :   b = hnf_Z_ZC(nf,a,b);
     234          21 :   return db? RgM_Rg_div(b, db): b;
     235             : }
     236             : /* b not a scalar (not point in trying to optimize for this case) */
     237             : static GEN
     238          28 : hnf_Q_QC(GEN nf, GEN a, GEN b)
     239             : {
     240             :   GEN da, db;
     241          28 :   if (typ(a) == t_INT) return hnf_Z_QC(nf, a, b);
     242           7 :   da = gel(a,2);
     243           7 :   a = gel(a,1);
     244           7 :   b = Q_remove_denom(b, &db);
     245             :   /* write da = d*A, db = d*B, gcd(A,B) = 1
     246             :    * gcd(a/(d A), b/(d B)) = gcd(a B, A b) / A B d = gcd(a B, b) / A B d */
     247           7 :   if (db)
     248             :   {
     249           7 :     GEN d = gcdii(da,db);
     250           7 :     if (!is_pm1(d)) db = diviiexact(db,d); /* B */
     251           7 :     if (!is_pm1(db))
     252             :     {
     253           7 :       a = mulii(a, db); /* a B */
     254           7 :       da = mulii(da, db); /* A B d = lcm(denom(a),denom(b)) */
     255             :     }
     256             :   }
     257           7 :   return RgM_Rg_div(hnf_Z_ZC(nf,a,b), da);
     258             : }
     259             : static GEN
     260           7 : hnf_QC_QC(GEN nf, GEN a, GEN b)
     261             : {
     262             :   GEN da, db, d, x;
     263           7 :   a = Q_remove_denom(a, &da);
     264           7 :   b = Q_remove_denom(b, &db);
     265           7 :   if (da) b = ZC_Z_mul(b, da);
     266           7 :   if (db) a = ZC_Z_mul(a, db);
     267           7 :   d = mul_denom(da, db);
     268           7 :   a = zk_multable(nf,a); da = zkmultable_capZ(a);
     269           7 :   b = zk_multable(nf,b); db = zkmultable_capZ(b);
     270           7 :   x = ZM_hnfmodid(shallowconcat(a,b), gcdii(da,db));
     271           7 :   return d? RgM_Rg_div(x, d): x;
     272             : }
     273             : static GEN
     274          21 : hnf_Q_Q(GEN nf, GEN a, GEN b) {return scalarmat(Q_gcd(a,b), nf_get_degree(nf));}
     275             : GEN
     276         119 : idealhnf0(GEN nf, GEN a, GEN b)
     277             : {
     278             :   long ta, tb;
     279             :   pari_sp av;
     280             :   GEN x;
     281         119 :   if (!b) return idealhnf(nf,a);
     282             : 
     283             :   /* HNF of aZ_K+bZ_K */
     284          56 :   av = avma; nf = checknf(nf);
     285          56 :   a = nf_to_scalar_or_basis(nf,a); ta = typ(a);
     286          56 :   b = nf_to_scalar_or_basis(nf,b); tb = typ(b);
     287          56 :   if (ta == t_COL)
     288          14 :     x = (tb==t_COL)? hnf_QC_QC(nf, a,b): hnf_Q_QC(nf, b,a);
     289             :   else
     290          42 :     x = (tb==t_COL)? hnf_Q_QC(nf, a,b): hnf_Q_Q(nf, a,b);
     291          56 :   return gerepileupto(av, x);
     292             : }
     293             : 
     294             : /*******************************************************************/
     295             : /*                                                                 */
     296             : /*                       TWO-ELEMENT FORM                          */
     297             : /*                                                                 */
     298             : /*******************************************************************/
     299             : static GEN idealapprfact_i(GEN nf, GEN x, int nored);
     300             : 
     301             : static int
     302      154173 : ok_elt(GEN x, GEN xZ, GEN y)
     303             : {
     304      154173 :   pari_sp av = avma;
     305      154173 :   int r = ZM_equal(x, ZM_hnfmodid(y, xZ));
     306      154173 :   avma = av; return r;
     307             : }
     308             : 
     309             : static GEN
     310       54446 : addmul_col(GEN a, long s, GEN b)
     311             : {
     312             :   long i,l;
     313       54446 :   if (!s) return a? leafcopy(a): a;
     314       54279 :   if (!a) return gmulsg(s,b);
     315       51135 :   l = lg(a);
     316      267096 :   for (i=1; i<l; i++)
     317      215961 :     if (signe(gel(b,i))) gel(a,i) = addii(gel(a,i), mulsi(s, gel(b,i)));
     318       51135 :   return a;
     319             : }
     320             : 
     321             : /* a <-- a + s * b, all coeffs integers */
     322             : static GEN
     323       23999 : addmul_mat(GEN a, long s, GEN b)
     324             : {
     325             :   long j,l;
     326             :   /* copy otherwise next call corrupts a */
     327       23999 :   if (!s) return a? RgM_shallowcopy(a): a;
     328       22516 :   if (!a) return gmulsg(s,b);
     329       12265 :   l = lg(a);
     330       58695 :   for (j=1; j<l; j++)
     331       46430 :     (void)addmul_col(gel(a,j), s, gel(b,j));
     332       12265 :   return a;
     333             : }
     334             : 
     335             : static GEN
     336       79192 : get_random_a(GEN nf, GEN x, GEN xZ)
     337             : {
     338             :   pari_sp av;
     339       79192 :   long i, lm, l = lg(x);
     340             :   GEN a, z, beta, mul;
     341             : 
     342       79192 :   beta= cgetg(l, t_VEC);
     343       79192 :   mul = cgetg(l, t_VEC); lm = 1; /* = lg(mul) */
     344             :   /* look for a in x such that a O/xZ = x O/xZ */
     345      158165 :   for (i = 2; i < l; i++)
     346             :   {
     347      155021 :     GEN xi = gel(x,i);
     348      155021 :     GEN t = FpM_red(zk_multable(nf,xi), xZ); /* ZM, cannot be a scalar */
     349      155021 :     if (gequal0(t)) continue;
     350      143922 :     if (ok_elt(x,xZ, t)) return xi;
     351       67874 :     gel(beta,lm) = xi;
     352             :     /* mul[i] = { canonical generators for x[i] O/xZ as Z-module } */
     353       67874 :     gel(mul,lm) = t; lm++;
     354             :   }
     355        3144 :   setlg(mul, lm);
     356        3144 :   setlg(beta,lm);
     357        3144 :   z = cgetg(lm, t_VECSMALL);
     358       10272 :   for(av = avma;; avma = av)
     359             :   {
     360       34271 :     for (a=NULL,i=1; i<lm; i++)
     361             :     {
     362       23999 :       long t = random_bits(4) - 7; /* in [-7,8] */
     363       23999 :       z[i] = t;
     364       23999 :       a = addmul_mat(a, t, gel(mul,i));
     365             :     }
     366             :     /* a = matrix (NOT HNF) of ideal generated by beta.z in O/xZ */
     367       10272 :     if (a && ok_elt(x,xZ, a)) break;
     368        7128 :   }
     369       11160 :   for (a=NULL,i=1; i<lm; i++)
     370        8016 :     a = addmul_col(a, z[i], gel(beta,i));
     371        3144 :   return a;
     372             : }
     373             : 
     374             : /* x square matrix, assume it is HNF */
     375             : static GEN
     376      196565 : mat_ideal_two_elt(GEN nf, GEN x)
     377             : {
     378             :   GEN y, a, cx, xZ;
     379      196565 :   long N = nf_get_degree(nf);
     380             :   pari_sp av, tetpil;
     381             : 
     382      196565 :   if (lg(x)-1 != N) pari_err_DIM("idealtwoelt");
     383      196551 :   if (N == 2) return mkvec2copy(gcoeff(x,1,1), gel(x,2));
     384             : 
     385       89556 :   y = cgetg(3,t_VEC); av = avma;
     386       89556 :   cx = Q_content(x);
     387       89556 :   xZ = gcoeff(x,1,1);
     388       89556 :   if (gequal(xZ, cx)) /* x = (cx) */
     389             :   {
     390        3423 :     gel(y,1) = cx;
     391        3423 :     gel(y,2) = gen_0; return y;
     392             :   }
     393       86133 :   if (equali1(cx)) cx = NULL;
     394             :   else
     395             :   {
     396        1680 :     x = Q_div_to_int(x, cx);
     397        1680 :     xZ = gcoeff(x,1,1);
     398             :   }
     399       86133 :   if (N < 6)
     400       74027 :     a = get_random_a(nf, x, xZ);
     401             :   else
     402             :   {
     403       12106 :     const long FB[] = { _evallg(15+1) | evaltyp(t_VECSMALL),
     404             :       2, 3, 5, 7, 11, 13, 17, 19, 23, 29, 31, 37, 41, 43, 47
     405             :     };
     406       12106 :     GEN P, E, a1 = Z_smoothen(xZ, (GEN)FB, &P, &E);
     407       12106 :     if (!a1) /* factors completely */
     408        6941 :       a = idealapprfact_i(nf, idealfactor(nf,x), 1);
     409        5165 :     else if (lg(P) == 1) /* no small factors */
     410        3820 :       a = get_random_a(nf, x, xZ);
     411             :     else /* general case */
     412             :     {
     413             :       GEN A0, A1, a0, u0, u1, v0, v1, pi0, pi1, t, u;
     414        1345 :       a0 = diviiexact(xZ, a1);
     415        1345 :       A0 = ZM_hnfmodid(x, a0); /* smooth part of x */
     416        1345 :       A1 = ZM_hnfmodid(x, a1); /* cofactor */
     417        1345 :       pi0 = idealapprfact_i(nf, idealfactor(nf,A0), 1);
     418        1345 :       pi1 = get_random_a(nf, A1, a1);
     419        1345 :       (void)bezout(a0, a1, &v0,&v1);
     420        1345 :       u0 = mulii(a0, v0);
     421        1345 :       u1 = mulii(a1, v1);
     422        1345 :       if (typ(pi0) != t_COL) t = addmulii(u0, pi0, u1);
     423             :       else
     424        1345 :       { t = ZC_Z_mul(pi0, u1); gel(t,1) = addii(gel(t,1), u0); }
     425        1345 :       u = ZC_Z_mul(pi1, u0); gel(u,1) = addii(gel(u,1), u1);
     426        1345 :       a = nfmuli(nf, centermod(u, xZ), centermod(t, xZ));
     427             :     }
     428             :   }
     429       86133 :   if (cx)
     430             :   {
     431        1680 :     a = centermod(a, xZ);
     432        1680 :     tetpil = avma;
     433        1680 :     if (typ(cx) == t_INT)
     434             :     {
     435         469 :       gel(y,1) = mulii(xZ, cx);
     436         469 :       gel(y,2) = ZC_Z_mul(a, cx);
     437             :     }
     438             :     else
     439             :     {
     440        1211 :       gel(y,1) = gmul(xZ, cx);
     441        1211 :       gel(y,2) = RgC_Rg_mul(a, cx);
     442             :     }
     443             :   }
     444             :   else
     445             :   {
     446       84453 :     tetpil = avma;
     447       84453 :     gel(y,1) = icopy(xZ);
     448       84453 :     gel(y,2) = centermod(a, xZ);
     449             :   }
     450       86133 :   gerepilecoeffssp(av,tetpil,y+1,2); return y;
     451             : }
     452             : 
     453             : /* Given an ideal x, returns [a,alpha] such that a is in Q,
     454             :  * x = a Z_K + alpha Z_K, alpha in K^*
     455             :  * a = 0 or alpha = 0 are possible, but do not try to determine whether
     456             :  * x is principal. */
     457             : GEN
     458       40588 : idealtwoelt(GEN nf, GEN x)
     459             : {
     460             :   pari_sp av;
     461             :   GEN z;
     462       40588 :   long tx = idealtyp(&x,&z);
     463       40581 :   nf = checknf(nf);
     464       40581 :   if (tx == id_MAT) return mat_ideal_two_elt(nf,x);
     465        1715 :   if (tx == id_PRIME) return mkvec2copy(gel(x,1), gel(x,2));
     466             :   /* id_PRINCIPAL */
     467         875 :   av = avma; x = nf_to_scalar_or_basis(nf, x);
     468        1554 :   return gerepilecopy(av, typ(x)==t_COL? mkvec2(gen_0,x):
     469         770 :                                          mkvec2(Q_abs_shallow(x),gen_0));
     470             : }
     471             : 
     472             : /*******************************************************************/
     473             : /*                                                                 */
     474             : /*                         FACTORIZATION                           */
     475             : /*                                                                 */
     476             : /*******************************************************************/
     477             : /* x integral ideal in HNF, Zval = v_p(x \cap Z) > 0; return v_p(Nx) */
     478             : static long
     479      200221 : idealHNF_norm_pval(GEN x, GEN p, long Zval)
     480             : {
     481      200221 :   long i, v = Zval, l = lg(x);
     482      200221 :   for (i = 2; i < l; i++) v += Z_pval(gcoeff(x,i,i), p);
     483      200221 :   return v;
     484             : }
     485             : 
     486             : /* x integral in HNF, f = partial factorization of x[1,1] = x\cap Z */
     487             : GEN
     488       37469 : idealHNF_Z_factor_i(GEN x, GEN f, GEN *pvN, GEN *pvZ)
     489             : {
     490             :   GEN P, E, vN, vZ;
     491             :   long i, l;
     492       37469 :   if (!f) f = Z_factor(gcoeff(x,1,1));
     493       37469 :   P = gel(f,1); l = lg(P);
     494       37469 :   E = gel(f,2);
     495       37469 :   *pvN = vN = cgetg(l, t_VECSMALL);
     496       37469 :   *pvZ = vZ = cgetg(l, t_VECSMALL);
     497       70864 :   for (i = 1; i < l; i++)
     498             :   {
     499       33395 :     vZ[i] = itou(gel(E,i));
     500       33395 :     vN[i] = idealHNF_norm_pval(x,gel(P,i), vZ[i]);
     501             :   }
     502       37469 :   return P;
     503             : }
     504             : /* return P, primes dividing Nx and xZ = x\cap Z, set v_p(Nx), v_p(xZ);
     505             :  * x integral in HNF */
     506             : GEN
     507           0 : idealHNF_Z_factor(GEN x, GEN *pvN, GEN *pvZ)
     508           0 : { return idealHNF_Z_factor_i(x, NULL, pvN, pvZ); }
     509             : 
     510             : /* v_P(A)*f(P) <= Nval [e.g. Nval = v_p(Norm A)], Zval = v_p(A \cap Z).
     511             :  * Return v_P(A) */
     512             : static long
     513      217970 : idealHNF_val(GEN A, GEN P, long Nval, long Zval)
     514             : {
     515      217970 :   long f = pr_get_f(P), vmax, v, e, i, j, k, l;
     516             :   GEN mul, B, a, y, r, p, pk, cx, vals;
     517             :   pari_sp av;
     518             : 
     519      217970 :   if (Nval < f) return 0;
     520      217906 :   p = pr_get_p(P);
     521      217906 :   e = pr_get_e(P);
     522             :   /* v_P(A) <= max [ e * v_p(A \cap Z), floor[v_p(Nix) / f ] */
     523      217906 :   vmax = minss(Zval * e, Nval / f);
     524      217906 :   mul = pr_get_tau(P);
     525      217906 :   l = lg(mul);
     526      217906 :   B = cgetg(l,t_MAT);
     527             :   /* B[1] not needed: v_pr(A[1]) = v_pr(A \cap Z) is known already */
     528      217906 :   gel(B,1) = gen_0; /* dummy */
     529      655158 :   for (j = 2; j < l; j++)
     530             :   {
     531      508505 :     GEN x = gel(A,j);
     532      508505 :     gel(B,j) = y = cgetg(l, t_COL);
     533     4074452 :     for (i = 1; i < l; i++)
     534             :     { /* compute a = (x.t0)_i, A in HNF ==> x[j+1..l-1] = 0 */
     535     3637200 :       a = mulii(gel(x,1), gcoeff(mul,i,1));
     536     3637200 :       for (k = 2; k <= j; k++) a = addii(a, mulii(gel(x,k), gcoeff(mul,i,k)));
     537             :       /* p | a ? */
     538     3637200 :       gel(y,i) = dvmdii(a,p,&r); if (signe(r)) return 0;
     539             :     }
     540             :   }
     541      146653 :   vals = cgetg(l, t_VECSMALL);
     542             :   /* vals[1] not needed */
     543      525028 :   for (j = 2; j < l; j++)
     544             :   {
     545      378375 :     gel(B,j) = Q_primitive_part(gel(B,j), &cx);
     546      378375 :     vals[j] = cx? 1 + e * Q_pval(cx, p): 1;
     547             :   }
     548      146653 :   pk = powiu(p, ceildivuu(vmax, e));
     549      146653 :   av = avma; y = cgetg(l,t_COL);
     550             :   /* can compute mod p^ceil((vmax-v)/e) */
     551      209384 :   for (v = 1; v < vmax; v++)
     552             :   { /* we know v_pr(Bj) >= v for all j */
     553       65987 :     if (e == 1 || (vmax - v) % e == 0) pk = diviiexact(pk, p);
     554      510200 :     for (j = 2; j < l; j++)
     555             :     {
     556      447469 :       GEN x = gel(B,j); if (v < vals[j]) continue;
     557     4469952 :       for (i = 1; i < l; i++)
     558             :       {
     559     4145678 :         pari_sp av2 = avma;
     560     4145678 :         a = mulii(gel(x,1), gcoeff(mul,i,1));
     561     4145678 :         for (k = 2; k < l; k++) a = addii(a, mulii(gel(x,k), gcoeff(mul,i,k)));
     562             :         /* a = (x.t_0)_i; p | a ? */
     563     4145678 :         a = dvmdii(a,p,&r); if (signe(r)) return v;
     564     4142422 :         if (lgefint(a) > lgefint(pk)) a = remii(a, pk);
     565     4142422 :         gel(y,i) = gerepileuptoint(av2, a);
     566             :       }
     567      324274 :       gel(B,j) = y; y = x;
     568      324274 :       if (gc_needed(av,3))
     569             :       {
     570           0 :         if(DEBUGMEM>1) pari_warn(warnmem,"idealval");
     571           0 :         gerepileall(av,3, &y,&B,&pk);
     572             :       }
     573             :     }
     574             :   }
     575      143397 :   return v;
     576             : }
     577             : /* true nf, x != 0 integral ideal in HNF, cx t_INT or NULL,
     578             :  * FA integer factorization matrix or NULL. Return partial factorization of
     579             :  * cx * x above primes in FA (complete factorization if !FA)*/
     580             : static GEN
     581       37469 : idealHNF_factor_i(GEN nf, GEN x, GEN cx, GEN FA)
     582             : {
     583       37469 :   const long N = lg(x)-1;
     584             :   long i, j, k, l, v;
     585       37469 :   GEN vN, vZ, vP, vE, vp = idealHNF_Z_factor_i(x, FA, &vN,&vZ);
     586             : 
     587       37469 :   l = lg(vp);
     588       37469 :   i = cx? expi(cx)+1: 1;
     589       37469 :   vP = cgetg((l+i-2)*N+1, t_COL);
     590       37469 :   vE = cgetg((l+i-2)*N+1, t_COL);
     591       70864 :   for (i = k = 1; i < l; i++)
     592             :   {
     593       33395 :     GEN L, p = gel(vp,i);
     594       33395 :     long Nval = vN[i], Zval = vZ[i], vc = cx? Z_pvalrem(cx,p,&cx): 0;
     595       33395 :     if (vc)
     596             :     {
     597        1813 :       L = idealprimedec(nf,p);
     598        1813 :       if (is_pm1(cx)) cx = NULL;
     599             :     }
     600             :     else
     601       31582 :       L = idealprimedec_limit_f(nf,p,Nval);
     602       51207 :     for (j = 1; j < lg(L); j++)
     603             :     {
     604       51144 :       GEN P = gel(L,j);
     605       51144 :       pari_sp av = avma;
     606       51144 :       v = idealHNF_val(x, P, Nval, Zval);
     607       51144 :       avma = av;
     608       51144 :       Nval -= v*pr_get_f(P);
     609       51144 :       v += vc * pr_get_e(P); if (!v) continue;
     610       37715 :       gel(vP,k) = P;
     611       37715 :       gel(vE,k) = utoipos(v); k++;
     612       37715 :       if (!Nval) break; /* now only the content contributes */
     613             :     }
     614       34110 :     if (vc) for (j++; j<lg(L); j++)
     615             :     {
     616         715 :       GEN P = gel(L,j);
     617         715 :       gel(vP,k) = P;
     618         715 :       gel(vE,k) = utoipos(vc * pr_get_e(P)); k++;
     619             :     }
     620             :   }
     621       37469 :   if (cx && !FA)
     622             :   { /* complete factorization */
     623        7504 :     GEN f = Z_factor(cx), cP = gel(f,1), cE = gel(f,2);
     624        7504 :     long lc = lg(cP);
     625       15757 :     for (i=1; i<lc; i++)
     626             :     {
     627        8253 :       GEN p = gel(cP,i), L = idealprimedec(nf,p);
     628        8253 :       long vc = itos(gel(cE,i));
     629       18284 :       for (j=1; j<lg(L); j++)
     630             :       {
     631       10031 :         GEN P = gel(L,j);
     632       10031 :         gel(vP,k) = P;
     633       10031 :         gel(vE,k) = utoipos(vc * pr_get_e(P)); k++;
     634             :       }
     635             :     }
     636             :   }
     637       37469 :   setlg(vP, k);
     638       37469 :   setlg(vE, k); return mkmat2(vP, vE);
     639             : }
     640             : /* true nf, x integral ideal */
     641             : static GEN
     642       37413 : idealHNF_factor(GEN nf, GEN x, ulong lim)
     643             : {
     644       37413 :   GEN cx, F = NULL;
     645       37413 :   if (lim)
     646             :   {
     647             :     GEN P, E;
     648             :     long l;
     649          21 :     F = Z_factor_limit(gcoeff(x,1,1), lim);
     650          21 :     P = gel(F,1); l = lg(P);
     651          21 :     E = gel(F,2);
     652          21 :     if (l > 1 && cmpiu(gel(P,l-1), lim) >= 0) { setlg(P,l-1); setlg(E,l-1); }
     653             :   }
     654       37413 :   x = Q_primitive_part(x, &cx);
     655       37413 :   return idealHNF_factor_i(nf, x, cx, F);
     656             : }
     657             : /* c * vector(#L,i,L[i].e), assume results fit in ulong */
     658             : static GEN
     659        3122 : prV_e_muls(GEN L, long c)
     660             : {
     661        3122 :   long j, l = lg(L);
     662        3122 :   GEN z = cgetg(l, t_COL);
     663        3122 :   for (j = 1; j < l; j++) gel(z,j) = stoi(c * pr_get_e(gel(L,j)));
     664        3122 :   return z;
     665             : }
     666             : /* true nf, y in Q */
     667             : static GEN
     668        3164 : Q_nffactor(GEN nf, GEN y, ulong lim)
     669             : {
     670             :   GEN f, P, E;
     671             :   long lfa, i;
     672        3164 :   if (typ(y) == t_INT)
     673             :   {
     674        3150 :     if (!signe(y)) pari_err_DOMAIN("idealfactor", "ideal", "=",gen_0,y);
     675        3129 :     if (is_pm1(y)) return trivial_fact();
     676             :   }
     677        2359 :   y = Q_abs_shallow(y);
     678        2359 :   f = lim? Q_factor_limit(y, lim): Q_factor(y);
     679        2359 :   P = gel(f,1); lfa = lg(P);
     680        2359 :   E = gel(f,2);
     681        5481 :   for (i = 1; i < lfa; i++)
     682             :   {
     683        3122 :     gel(P,i) = idealprimedec(nf, gel(P,i));
     684        3122 :     gel(E,i) = prV_e_muls(gel(P,i), itos(gel(E,i)));
     685             :   }
     686        2359 :   settyp(P,t_VEC); P = shallowconcat1(P);
     687        2359 :   settyp(E,t_VEC); E = shallowconcat1(E);
     688        2359 :   gel(f,1) = P; settyp(P, t_COL);
     689        2359 :   gel(f,2) = E; return f;
     690             : }
     691             : 
     692             : GEN
     693       40619 : idealfactor_limit(GEN nf, GEN x, ulong lim)
     694             : {
     695       40619 :   pari_sp av = avma;
     696             :   GEN fa, y;
     697       40619 :   long tx = idealtyp(&x,&y);
     698             : 
     699       40619 :   nf = checknf(nf);
     700       40619 :   if (tx == id_PRIME)
     701             :   {
     702          49 :     if (lim && cmpiu(pr_get_p(x), lim) >= 0) return trivial_fact();
     703          42 :     retmkmat2(mkcolcopy(x), mkcol(gen_1));
     704             :   }
     705       40570 :   if (tx == id_PRINCIPAL)
     706             :   {
     707        5117 :     y = nf_to_scalar_or_basis(nf, x);
     708        5117 :     if (typ(y) != t_COL) return gerepilecopy(av, Q_nffactor(nf, y, lim));
     709             :   }
     710       37406 :   y = idealnumden(nf, x);
     711       37406 :   fa = idealHNF_factor(nf, gel(y,1), lim);
     712       37406 :   if (!isint1(gel(y,2)))
     713           7 :     fa = famat_div_shallow(fa, idealHNF_factor(nf, gel(y,2), lim));
     714       37406 :   fa = gerepilecopy(av, fa);
     715       37406 :   return sort_factor(fa, (void*)&cmp_prime_ideal, &cmp_nodata);
     716             : }
     717             : GEN
     718       40521 : idealfactor(GEN nf, GEN x) { return idealfactor_limit(nf, x, 0); }
     719             : GEN
     720          98 : gpidealfactor(GEN nf, GEN x, GEN lim)
     721             : {
     722          98 :   ulong L = 0;
     723          98 :   if (lim)
     724             :   {
     725          35 :     if (typ(lim) != t_INT || signe(lim) < 0) pari_err_FLAG("idealfactor");
     726          35 :     L = itou(lim);
     727             :   }
     728          98 :   return idealfactor_limit(nf, x, L);
     729             : }
     730             : 
     731             : /* true nf; A is assumed to be the n-th power of an integral ideal,
     732             :  * return its n-th root; n > 1 */
     733             : static long
     734          91 : idealsqrtn_int(GEN nf, GEN A, long n, GEN *pB)
     735             : {
     736             :   GEN C, ram, vram, root;
     737             :   long i, l;
     738             : 
     739          91 :   if (typ(A) == t_INT) return Z_ispowerall(A, n, pB);
     740             :   /* compute valuations at ramified primes */
     741          42 :   ram = gel(idealfactor(nf, idealadd(nf, nf_get_diff(nf),A)), 1);
     742          42 :   l = lg(ram); vram = cgetg(l, t_VECSMALL);
     743          56 :   for (i = 1; i < l; i++)
     744             :   {
     745          14 :     long v = idealval(nf,A,gel(ram,i));
     746          14 :     if (v % n) return 0;
     747          14 :     vram[i] = v / n;
     748             :   }
     749          42 :   root = idealfactorback(nf, ram, vram, 0);
     750             :   /* remove ramified primes */
     751          42 :   if (isint1(root))
     752          28 :     root = matid(nf_get_degree(nf));
     753             :   else
     754          14 :     A = idealdivexact(nf, A, idealpows(nf,root,n));
     755          42 :   A = Q_primitive_part(A, &C);
     756          42 :   if (C)
     757             :   {
     758           0 :     if (!Z_ispowerall(C,n,&C)) return 0;
     759           0 :     if (pB) root = ZM_Z_mul(root, C);
     760             :   }
     761             : 
     762             :   /* compute final n-th root, at most degree(nf)-1 iterations */
     763          77 :   for (i = 0;; i++)
     764             :   {
     765          77 :     GEN J, b, a = gcoeff(A,1,1); /* A \cap Z */
     766          77 :     if (is_pm1(a)) break;
     767          42 :     if (!Z_ispowerall(a,n,&b)) return 0;
     768          35 :     J = idealadd(nf, b, A);
     769          35 :     A = idealdivexact(nf, idealpows(nf,J,n), A);
     770          35 :     if (pB) root = odd(i)? idealdivexact(nf, root, J): idealmul(nf, root, J);
     771          35 :   }
     772          70 :   if (pB) *pB = root;
     773          35 :   return 1;
     774             : }
     775             : 
     776             : /* A is assumed to be the n-th power of an ideal in nf
     777             :  returns its n-th root. */
     778             : long
     779          56 : idealispower(GEN nf, GEN A, long n, GEN *pB)
     780             : {
     781          56 :   pari_sp av = avma;
     782             :   GEN v, N, D;
     783          56 :   nf = checknf(nf);
     784          56 :   if (n <= 0) pari_err_DOMAIN("idealispower", "n", "<=", gen_0, stoi(n));
     785          56 :   if (n == 1) { if (pB) *pB = idealhnf(nf,A); return 1; }
     786          49 :   v = idealnumden(nf,A);
     787          49 :   if (gequal0(gel(v,1))) { avma = av; if (pB) *pB = cgetg(1,t_MAT); return 1; }
     788          49 :   if (!idealsqrtn_int(nf, gel(v,1), n, pB? &N: NULL)) return 0;
     789          42 :   if (!idealsqrtn_int(nf, gel(v,2), n, pB? &D: NULL)) return 0;
     790          42 :   if (pB) *pB = gerepileupto(av, idealdiv(nf,N,D)); else avma = av;
     791          42 :   return 1;
     792             : }
     793             : 
     794             : /* x t_INT or integral non-0 ideal in HNF */
     795             : static GEN
     796         112 : idealredmodpower_i(GEN nf, GEN x, ulong k, ulong B)
     797             : {
     798             :   GEN cx, y, U, N, F, Q;
     799             :   long nF;
     800         112 :   if (typ(x) == t_INT)
     801             :   {
     802          49 :     if (!signe(x) || is_pm1(x)) return gen_1;
     803           0 :     F = Z_factor_limit(x, B);
     804           0 :     gel(F,2) = gdiventgs(gel(F,2), k);
     805           0 :     return factorback(F);
     806             :   }
     807          63 :   N = gcoeff(x,1,1); if (is_pm1(N)) return gen_1;
     808          56 :   F = Z_factor_limit(N, B); nF=lg(gel(F,1))-1;
     809          56 :   if (BPSW_psp(gcoeff(F,nF,1))) U = NULL;
     810             :   else
     811             :   {
     812          28 :     GEN M = powii(gcoeff(F,nF,1), gcoeff(F,nF,2));
     813          28 :     y = hnfmodid(x, M); /* coprime part to B! */
     814          28 :     if (!idealispower(nf, y, k, &U)) U = NULL;
     815          28 :     x = hnfmodid(x, diviiexact(N, M));
     816          28 :     setlg(gel(F,1), nF); /* remove last entry (unfactored part) */
     817          28 :     setlg(gel(F,2), nF);
     818             :   }
     819             :   /* x = B-smooth part of initial x */
     820          56 :   x = Q_primitive_part(x, &cx);
     821          56 :   F = idealHNF_factor_i(nf, x, cx, F);
     822          56 :   gel(F,2) = gdiventgs(gel(F,2), k);
     823          56 :   Q = idealfactorback(nf, gel(F,1), gel(F,2), 0);
     824          56 :   if (U) Q = idealmul(nf,Q,U);
     825          56 :   if (typ(Q) == t_INT) return Q;
     826          56 :   y = idealred_elt(nf, idealHNF_inv_Z(nf, Q));
     827          56 :   return gdiv(y, gcoeff(Q,1,1));
     828             : }
     829             : GEN
     830          63 : idealredmodpower(GEN nf, GEN x, ulong n, ulong B)
     831             : {
     832          63 :   pari_sp av = avma;
     833             :   GEN a, b;
     834          63 :   nf = checknf(nf);
     835          63 :   if (!n) pari_err_DOMAIN("idealredmodpower","n", "=", gen_0, gen_0);
     836          63 :   x = idealnumden(nf, x);
     837          63 :   a = gel(x,1);
     838          63 :   if (isintzero(a)) { avma = av; return gen_1; }
     839          56 :   a = idealredmodpower_i(nf, gel(x,1), n, B);
     840          56 :   b = idealredmodpower_i(nf, gel(x,2), n, B);
     841          56 :   if (!isint1(b)) a = nf_to_scalar_or_basis(nf, nfdiv(nf, a, b));
     842          56 :   return gerepilecopy(av, a);
     843             : }
     844             : 
     845             : /* P prime ideal in idealprimedec format. Return valuation(A) at P */
     846             : long
     847      492362 : idealval(GEN nf, GEN A, GEN P)
     848             : {
     849      492362 :   pari_sp av = avma;
     850             :   GEN a, p, cA;
     851      492362 :   long vcA, v, Zval, tx = idealtyp(&A,&a);
     852             : 
     853      492362 :   if (tx == id_PRINCIPAL) return nfval(nf,A,P);
     854      487931 :   checkprid(P);
     855      487931 :   if (tx == id_PRIME) return pr_equal(P, A)? 1: 0;
     856             :   /* id_MAT */
     857      487903 :   nf = checknf(nf);
     858      487903 :   A = Q_primitive_part(A, &cA);
     859      487903 :   p = pr_get_p(P);
     860      487903 :   vcA = cA? Q_pval(cA,p): 0;
     861      487903 :   if (pr_is_inert(P)) { avma = av; return vcA; }
     862      480294 :   Zval = Z_pval(gcoeff(A,1,1), p);
     863      480294 :   if (!Zval) v = 0;
     864             :   else
     865             :   {
     866      166826 :     long Nval = idealHNF_norm_pval(A, p, Zval);
     867      166826 :     v = idealHNF_val(A, P, Nval, Zval);
     868             :   }
     869      480294 :   avma = av; return vcA? v + vcA*pr_get_e(P): v;
     870             : }
     871             : GEN
     872        6573 : gpidealval(GEN nf, GEN ix, GEN P)
     873             : {
     874        6573 :   long v = idealval(nf,ix,P);
     875        6573 :   return v == LONG_MAX? mkoo(): stoi(v);
     876             : }
     877             : 
     878             : /* gcd and generalized Bezout */
     879             : 
     880             : GEN
     881       60368 : idealadd(GEN nf, GEN x, GEN y)
     882             : {
     883       60368 :   pari_sp av = avma;
     884             :   long tx, ty;
     885             :   GEN z, a, dx, dy, dz;
     886             : 
     887       60368 :   tx = idealtyp(&x,&z);
     888       60368 :   ty = idealtyp(&y,&z); nf = checknf(nf);
     889       60368 :   if (tx != id_MAT) x = idealhnf_shallow(nf,x);
     890       60368 :   if (ty != id_MAT) y = idealhnf_shallow(nf,y);
     891       60368 :   if (lg(x) == 1) return gerepilecopy(av,y);
     892       60354 :   if (lg(y) == 1) return gerepilecopy(av,x); /* check for 0 ideal */
     893       60053 :   dx = Q_denom(x);
     894       60053 :   dy = Q_denom(y); dz = lcmii(dx,dy);
     895       60053 :   if (is_pm1(dz)) dz = NULL; else {
     896       12621 :     x = Q_muli_to_int(x, dz);
     897       12621 :     y = Q_muli_to_int(y, dz);
     898             :   }
     899       60053 :   a = gcdii(gcoeff(x,1,1), gcoeff(y,1,1));
     900       60053 :   if (is_pm1(a))
     901             :   {
     902       27971 :     long N = lg(x)-1;
     903       27971 :     if (!dz) { avma = av; return matid(N); }
     904        3611 :     return gerepileupto(av, scalarmat(ginv(dz), N));
     905             :   }
     906       32082 :   z = ZM_hnfmodid(shallowconcat(x,y), a);
     907       32082 :   if (dz) z = RgM_Rg_div(z,dz);
     908       32082 :   return gerepileupto(av,z);
     909             : }
     910             : 
     911             : static GEN
     912          28 : trivial_merge(GEN x)
     913          28 : { return (lg(x) == 1 || !is_pm1(gcoeff(x,1,1)))? NULL: gen_1; }
     914             : /* true nf */
     915             : static GEN
     916      120715 : _idealaddtoone(GEN nf, GEN x, GEN y, long red)
     917             : {
     918             :   GEN a;
     919      120715 :   long tx = idealtyp(&x, &a/*junk*/);
     920      120715 :   long ty = idealtyp(&y, &a/*junk*/);
     921             :   long ea;
     922      120715 :   if (tx != id_MAT) x = idealhnf_shallow(nf, x);
     923      120715 :   if (ty != id_MAT) y = idealhnf_shallow(nf, y);
     924      120715 :   if (lg(x) == 1)
     925          14 :     a = trivial_merge(y);
     926      120701 :   else if (lg(y) == 1)
     927          14 :     a = trivial_merge(x);
     928             :   else
     929      120687 :     a = hnfmerge_get_1(x, y);
     930      120715 :   if (!a) pari_err_COPRIME("idealaddtoone",x,y);
     931      120701 :   if (red && (ea = gexpo(a)) > 10)
     932             :   {
     933        6687 :     GEN b = (typ(a) == t_COL)? a: scalarcol_shallow(a, nf_get_degree(nf));
     934        6687 :     b = ZC_reducemodlll(b, idealHNF_mul(nf,x,y));
     935        6687 :     if (gexpo(b) < ea) a = b;
     936             :   }
     937      120701 :   return a;
     938             : }
     939             : /* true nf */
     940             : GEN
     941       12334 : idealaddtoone_i(GEN nf, GEN x, GEN y)
     942       12334 : { return _idealaddtoone(nf, x, y, 1); }
     943             : /* true nf */
     944             : GEN
     945      108381 : idealaddtoone_raw(GEN nf, GEN x, GEN y)
     946      108381 : { return _idealaddtoone(nf, x, y, 0); }
     947             : 
     948             : GEN
     949          98 : idealaddtoone(GEN nf, GEN x, GEN y)
     950             : {
     951          98 :   GEN z = cgetg(3,t_VEC), a;
     952          98 :   pari_sp av = avma;
     953          98 :   nf = checknf(nf);
     954          98 :   a = gerepileupto(av, idealaddtoone_i(nf,x,y));
     955          84 :   gel(z,1) = a;
     956          84 :   gel(z,2) = typ(a) == t_COL? Z_ZC_sub(gen_1,a): subui(1,a);
     957          84 :   return z;
     958             : }
     959             : 
     960             : /* assume elements of list are integral ideals */
     961             : GEN
     962          35 : idealaddmultoone(GEN nf, GEN list)
     963             : {
     964          35 :   pari_sp av = avma;
     965          35 :   long N, i, l, nz, tx = typ(list);
     966             :   GEN H, U, perm, L;
     967             : 
     968          35 :   nf = checknf(nf); N = nf_get_degree(nf);
     969          35 :   if (!is_vec_t(tx)) pari_err_TYPE("idealaddmultoone",list);
     970          35 :   l = lg(list);
     971          35 :   L = cgetg(l, t_VEC);
     972          35 :   if (l == 1)
     973           0 :     pari_err_DOMAIN("idealaddmultoone", "sum(ideals)", "!=", gen_1, L);
     974          35 :   nz = 0; /* number of non-zero ideals in L */
     975          98 :   for (i=1; i<l; i++)
     976             :   {
     977          70 :     GEN I = gel(list,i);
     978          70 :     if (typ(I) != t_MAT) I = idealhnf_shallow(nf,I);
     979          70 :     if (lg(I) != 1)
     980             :     {
     981          42 :       nz++; RgM_check_ZM(I,"idealaddmultoone");
     982          35 :       if (lgcols(I) != N+1) pari_err_TYPE("idealaddmultoone [not an ideal]", I);
     983             :     }
     984          63 :     gel(L,i) = I;
     985             :   }
     986          28 :   H = ZM_hnfperm(shallowconcat1(L), &U, &perm);
     987          28 :   if (lg(H) == 1 || !equali1(gcoeff(H,1,1)))
     988           7 :     pari_err_DOMAIN("idealaddmultoone", "sum(ideals)", "!=", gen_1, L);
     989          49 :   for (i=1; i<=N; i++)
     990          49 :     if (perm[i] == 1) break;
     991          21 :   U = gel(U,(nz-1)*N + i); /* (L[1]|...|L[nz]) U = 1 */
     992          21 :   nz = 0;
     993          63 :   for (i=1; i<l; i++)
     994             :   {
     995          42 :     GEN c = gel(L,i);
     996          42 :     if (lg(c) == 1)
     997          14 :       c = gen_0;
     998             :     else {
     999          28 :       c = ZM_ZC_mul(c, vecslice(U, nz*N + 1, (nz+1)*N));
    1000          28 :       nz++;
    1001             :     }
    1002          42 :     gel(L,i) = c;
    1003             :   }
    1004          21 :   return gerepilecopy(av, L);
    1005             : }
    1006             : 
    1007             : /* multiplication */
    1008             : 
    1009             : /* x integral ideal (without archimedean component) in HNF form
    1010             :  * y = [a,alpha] corresponds to the integral ideal aZ_K+alpha Z_K, a in Z,
    1011             :  * alpha a ZV or a ZM (multiplication table). Multiply them */
    1012             : static GEN
    1013      633032 : idealHNF_mul_two(GEN nf, GEN x, GEN y)
    1014             : {
    1015      633032 :   GEN m, a = gel(y,1), alpha = gel(y,2);
    1016             :   long i, N;
    1017             : 
    1018      633032 :   if (typ(alpha) != t_MAT)
    1019             :   {
    1020      429074 :     alpha = zk_scalar_or_multable(nf, alpha);
    1021      429074 :     if (typ(alpha) == t_INT) /* e.g. y inert ? 0 should not (but may) occur */
    1022        3206 :       return signe(a)? ZM_Z_mul(x, gcdii(a, alpha)): cgetg(1,t_MAT);
    1023             :   }
    1024      629826 :   N = lg(x)-1; m = cgetg((N<<1)+1,t_MAT);
    1025      629826 :   for (i=1; i<=N; i++) gel(m,i)   = ZM_ZC_mul(alpha,gel(x,i));
    1026      629826 :   for (i=1; i<=N; i++) gel(m,i+N) = ZC_Z_mul(gel(x,i), a);
    1027      629826 :   return ZM_hnfmodid(m, mulii(a, gcoeff(x,1,1)));
    1028             : }
    1029             : 
    1030             : /* Assume ix and iy are integral in HNF form [NOT extended]. Not memory clean.
    1031             :  * HACK: ideal in iy can be of the form [a,b], a in Z, b in Z_K */
    1032             : GEN
    1033      310928 : idealHNF_mul(GEN nf, GEN x, GEN y)
    1034             : {
    1035             :   GEN z;
    1036      310928 :   if (typ(y) == t_VEC)
    1037      212602 :     z = idealHNF_mul_two(nf,x,y);
    1038             :   else
    1039             :   { /* reduce one ideal to two-elt form. The smallest */
    1040       98326 :     GEN xZ = gcoeff(x,1,1), yZ = gcoeff(y,1,1);
    1041       98326 :     if (cmpii(xZ, yZ) < 0)
    1042             :     {
    1043       35026 :       if (is_pm1(xZ)) return gcopy(y);
    1044       24597 :       z = idealHNF_mul_two(nf, y, mat_ideal_two_elt(nf,x));
    1045             :     }
    1046             :     else
    1047             :     {
    1048       63300 :       if (is_pm1(yZ)) return gcopy(x);
    1049       39003 :       z = idealHNF_mul_two(nf, x, mat_ideal_two_elt(nf,y));
    1050             :     }
    1051             :   }
    1052      276202 :   return z;
    1053             : }
    1054             : 
    1055             : /* operations on elements in factored form */
    1056             : 
    1057             : GEN
    1058       93385 : famat_mul_shallow(GEN f, GEN g)
    1059             : {
    1060       93385 :   if (typ(f) != t_MAT) f = to_famat_shallow(f,gen_1);
    1061       93385 :   if (typ(g) != t_MAT) g = to_famat_shallow(g,gen_1);
    1062       93385 :   if (lg(f) == 1) return g;
    1063       74206 :   if (lg(g) == 1) return f;
    1064      145890 :   return mkmat2(shallowconcat(gel(f,1), gel(g,1)),
    1065      145890 :                 shallowconcat(gel(f,2), gel(g,2)));
    1066             : }
    1067             : GEN
    1068       59234 : famat_mulpow_shallow(GEN f, GEN g, GEN e)
    1069             : {
    1070       59234 :   if (!signe(e)) return f;
    1071       59108 :   return famat_mul_shallow(f, famat_pow_shallow(g, e));
    1072             : }
    1073             : 
    1074             : GEN
    1075        4466 : famat_mulpows_shallow(GEN f, GEN g, long e)
    1076             : {
    1077        4466 :   if (e==0) return f;
    1078        2820 :   return famat_mul_shallow(f, famat_pows_shallow(g, e));
    1079             : }
    1080             : 
    1081             : GEN
    1082           7 : famat_div_shallow(GEN f, GEN g)
    1083           7 : { return famat_mul_shallow(f, famat_inv_shallow(g)); }
    1084             : 
    1085             : GEN
    1086           0 : to_famat(GEN x, GEN y) { retmkmat2(mkcolcopy(x), mkcolcopy(y)); }
    1087             : GEN
    1088      879563 : to_famat_shallow(GEN x, GEN y) { return mkmat2(mkcol(x), mkcol(y)); }
    1089             : 
    1090             : /* concat the single elt x; not gconcat since x may be a t_COL */
    1091             : static GEN
    1092       58824 : append(GEN v, GEN x)
    1093             : {
    1094       58824 :   long i, l = lg(v);
    1095       58824 :   GEN w = cgetg(l+1, typ(v));
    1096       58824 :   for (i=1; i<l; i++) gel(w,i) = gcopy(gel(v,i));
    1097       58824 :   gel(w,i) = gcopy(x); return w;
    1098             : }
    1099             : /* add x^1 to famat f */
    1100             : static GEN
    1101       85251 : famat_add(GEN f, GEN x)
    1102             : {
    1103       85251 :   GEN h = cgetg(3,t_MAT);
    1104       85251 :   if (lg(f) == 1)
    1105             :   {
    1106       26427 :     gel(h,1) = mkcolcopy(x);
    1107       26427 :     gel(h,2) = mkcol(gen_1);
    1108             :   }
    1109             :   else
    1110             :   {
    1111       58824 :     gel(h,1) = append(gel(f,1), x);
    1112       58824 :     gel(h,2) = gconcat(gel(f,2), gen_1);
    1113             :   }
    1114       85251 :   return h;
    1115             : }
    1116             : 
    1117             : GEN
    1118      108975 : famat_mul(GEN f, GEN g)
    1119             : {
    1120             :   GEN h;
    1121      108975 :   if (typ(g) != t_MAT) {
    1122       85251 :     if (typ(f) == t_MAT) return famat_add(f, g);
    1123           0 :     h = cgetg(3, t_MAT);
    1124           0 :     gel(h,1) = mkcol2(gcopy(f), gcopy(g));
    1125           0 :     gel(h,2) = mkcol2(gen_1, gen_1);
    1126             :   }
    1127       23724 :   if (typ(f) != t_MAT) return famat_add(g, f);
    1128       23724 :   if (lg(f) == 1) return gcopy(g);
    1129        4447 :   if (lg(g) == 1) return gcopy(f);
    1130        1955 :   h = cgetg(3,t_MAT);
    1131        1955 :   gel(h,1) = gconcat(gel(f,1), gel(g,1));
    1132        1955 :   gel(h,2) = gconcat(gel(f,2), gel(g,2));
    1133        1955 :   return h;
    1134             : }
    1135             : 
    1136             : GEN
    1137       50136 : famat_sqr(GEN f)
    1138             : {
    1139             :   GEN h;
    1140       50136 :   if (lg(f) == 1) return cgetg(1,t_MAT);
    1141       25155 :   if (typ(f) != t_MAT) return to_famat(f,gen_2);
    1142       25155 :   h = cgetg(3,t_MAT);
    1143       25155 :   gel(h,1) = gcopy(gel(f,1));
    1144       25155 :   gel(h,2) = gmul2n(gel(f,2),1);
    1145       25155 :   return h;
    1146             : }
    1147             : 
    1148             : GEN
    1149       27027 : famat_inv_shallow(GEN f)
    1150             : {
    1151       27027 :   if (lg(f) == 1) return f;
    1152       27027 :   if (typ(f) != t_MAT) return to_famat_shallow(f,gen_m1);
    1153          21 :   return mkmat2(gel(f,1), ZC_neg(gel(f,2)));
    1154             : }
    1155             : GEN
    1156       11060 : famat_inv(GEN f)
    1157             : {
    1158       11060 :   if (lg(f) == 1) return cgetg(1,t_MAT);
    1159        4159 :   if (typ(f) != t_MAT) return to_famat(f,gen_m1);
    1160        4159 :   retmkmat2(gcopy(gel(f,1)), ZC_neg(gel(f,2)));
    1161             : }
    1162             : GEN
    1163        1167 : famat_pow(GEN f, GEN n)
    1164             : {
    1165        1167 :   if (lg(f) == 1) return cgetg(1,t_MAT);
    1166           0 :   if (typ(f) != t_MAT) return to_famat(f,n);
    1167           0 :   retmkmat2(gcopy(gel(f,1)), ZC_Z_mul(gel(f,2),n));
    1168             : }
    1169             : GEN
    1170       59108 : famat_pow_shallow(GEN f, GEN n)
    1171             : {
    1172       59108 :   if (is_pm1(n)) return signe(n) > 0? f: famat_inv_shallow(f);
    1173       30296 :   if (lg(f) == 1) return f;
    1174       30296 :   if (typ(f) != t_MAT) return to_famat_shallow(f,n);
    1175         168 :   return mkmat2(gel(f,1), ZC_Z_mul(gel(f,2),n));
    1176             : }
    1177             : 
    1178             : GEN
    1179        2820 : famat_pows_shallow(GEN f, long n)
    1180             : {
    1181        2820 :   if (n==1) return f;
    1182        1330 :   if (n==-1) return famat_inv_shallow(f);
    1183        1197 :   if (lg(f) == 1) return f;
    1184        1197 :   if (typ(f) != t_MAT) return to_famat_shallow(f, stoi(n));
    1185        1112 :   return mkmat2(gel(f,1), ZC_z_mul(gel(f,2),n));
    1186             : }
    1187             : 
    1188             : GEN
    1189           0 : famat_Z_gcd(GEN M, GEN n)
    1190             : {
    1191           0 :   pari_sp av=avma;
    1192           0 :   long i, j, l=lgcols(M);
    1193           0 :   GEN F=cgetg(3,t_MAT);
    1194           0 :   gel(F,1)=cgetg(l,t_COL);
    1195           0 :   gel(F,2)=cgetg(l,t_COL);
    1196           0 :   for (i=1, j=1; i<l; i++)
    1197             :   {
    1198           0 :     GEN p = gcoeff(M,i,1);
    1199           0 :     GEN e = gminsg(Z_pval(n,p),gcoeff(M,i,2));
    1200           0 :     if (signe(e))
    1201             :     {
    1202           0 :       gcoeff(F,j,1)=p;
    1203           0 :       gcoeff(F,j,2)=e;
    1204           0 :       j++;
    1205             :     }
    1206             :   }
    1207           0 :   setlg(gel(F,1),j); setlg(gel(F,2),j);
    1208           0 :   return gerepilecopy(av,F);
    1209             : }
    1210             : 
    1211             : /* x assumed to be a t_MATs (factorization matrix), or compatible with
    1212             :  * the element_* functions. */
    1213             : static GEN
    1214       60713 : ext_sqr(GEN nf, GEN x)
    1215       60713 : { return (typ(x)==t_MAT)? famat_sqr(x): nfsqr(nf, x); }
    1216             : static GEN
    1217      143681 : ext_mul(GEN nf, GEN x, GEN y)
    1218      143681 : { return (typ(x)==t_MAT)? famat_mul(x,y): nfmul(nf, x, y); }
    1219             : static GEN
    1220       10920 : ext_inv(GEN nf, GEN x)
    1221       10920 : { return (typ(x)==t_MAT)? famat_inv(x): nfinv(nf, x); }
    1222             : static GEN
    1223        1167 : ext_pow(GEN nf, GEN x, GEN n)
    1224        1167 : { return (typ(x)==t_MAT)? famat_pow(x,n): nfpow(nf, x, n); }
    1225             : 
    1226             : GEN
    1227           0 : famat_to_nf(GEN nf, GEN f)
    1228             : {
    1229             :   GEN t, x, e;
    1230             :   long i;
    1231           0 :   if (lg(f) == 1) return gen_1;
    1232             : 
    1233           0 :   x = gel(f,1);
    1234           0 :   e = gel(f,2);
    1235           0 :   t = nfpow(nf, gel(x,1), gel(e,1));
    1236           0 :   for (i=lg(x)-1; i>1; i--)
    1237           0 :     t = nfmul(nf, t, nfpow(nf, gel(x,i), gel(e,i)));
    1238           0 :   return t;
    1239             : }
    1240             : 
    1241             : GEN
    1242       18032 : famat_reduce(GEN fa)
    1243             : {
    1244             :   GEN E, G, L, g, e;
    1245             :   long i, k, l;
    1246             : 
    1247       18032 :   if (lg(fa) == 1) return fa;
    1248       15477 :   g = gel(fa,1); l = lg(g);
    1249       15477 :   e = gel(fa,2);
    1250       15477 :   L = gen_indexsort(g, (void*)&cmp_universal, &cmp_nodata);
    1251       15477 :   G = cgetg(l, t_COL);
    1252       15477 :   E = cgetg(l, t_COL);
    1253             :   /* merge */
    1254       37932 :   for (k=i=1; i<l; i++,k++)
    1255             :   {
    1256       22455 :     gel(G,k) = gel(g,L[i]);
    1257       22455 :     gel(E,k) = gel(e,L[i]);
    1258       22455 :     if (k > 1 && gidentical(gel(G,k), gel(G,k-1)))
    1259             :     {
    1260         735 :       gel(E,k-1) = addii(gel(E,k), gel(E,k-1));
    1261         735 :       k--;
    1262             :     }
    1263             :   }
    1264             :   /* kill 0 exponents */
    1265       15477 :   l = k;
    1266       37197 :   for (k=i=1; i<l; i++)
    1267       21720 :     if (!gequal0(gel(E,i)))
    1268             :     {
    1269       20985 :       gel(G,k) = gel(G,i);
    1270       20985 :       gel(E,k) = gel(E,i); k++;
    1271             :     }
    1272       15477 :   setlg(G, k);
    1273       15477 :   setlg(E, k); return mkmat2(G,E);
    1274             : }
    1275             : 
    1276             : GEN
    1277       12628 : famatsmall_reduce(GEN fa)
    1278             : {
    1279             :   GEN E, G, L, g, e;
    1280             :   long i, k, l;
    1281       12628 :   if (lg(fa) == 1) return fa;
    1282       12628 :   g = gel(fa,1); l = lg(g);
    1283       12628 :   e = gel(fa,2);
    1284       12628 :   L = vecsmall_indexsort(g);
    1285       12628 :   G = cgetg(l, t_VECSMALL);
    1286       12628 :   E = cgetg(l, t_VECSMALL);
    1287             :   /* merge */
    1288      112948 :   for (k=i=1; i<l; i++,k++)
    1289             :   {
    1290      100320 :     G[k] = g[L[i]];
    1291      100320 :     E[k] = e[L[i]];
    1292      100320 :     if (k > 1 && G[k] == G[k-1])
    1293             :     {
    1294        5848 :       E[k-1] += E[k];
    1295        5848 :       k--;
    1296             :     }
    1297             :   }
    1298             :   /* kill 0 exponents */
    1299       12628 :   l = k;
    1300      107098 :   for (k=i=1; i<l; i++)
    1301       94470 :     if (E[i])
    1302             :     {
    1303       91576 :       G[k] = G[i];
    1304       91576 :       E[k] = E[i]; k++;
    1305             :     }
    1306       12628 :   setlg(G, k);
    1307       12628 :   setlg(E, k); return mkmat2(G,E);
    1308             : }
    1309             : 
    1310             : GEN
    1311       55930 : ZM_famat_limit(GEN fa, GEN limit)
    1312             : {
    1313             :   pari_sp av;
    1314             :   GEN E, G, g, e, r;
    1315             :   long i, k, l, n, lG;
    1316             : 
    1317       55930 :   if (lg(fa) == 1) return fa;
    1318       55930 :   g = gel(fa,1); l = lg(g);
    1319       55930 :   e = gel(fa,2);
    1320      124152 :   for(n=0, i=1; i<l; i++)
    1321       68222 :     if (cmpii(gel(g,i),limit)<=0) n++;
    1322       55930 :   lG = n<l-1 ? n+2 : n+1;
    1323       55930 :   G = cgetg(lG, t_COL);
    1324       55930 :   E = cgetg(lG, t_COL);
    1325       55930 :   av = avma;
    1326      124152 :   for (i=1, k=1, r = gen_1; i<l; i++)
    1327             :   {
    1328       68222 :     if (cmpii(gel(g,i),limit)<=0)
    1329             :     {
    1330       68138 :       gel(G,k) = gel(g,i);
    1331       68138 :       gel(E,k) = gel(e,i);
    1332       68138 :       k++;
    1333          84 :     } else r = mulii(r, powii(gel(g,i), gel(e,i)));
    1334             :   }
    1335       55930 :   if (k<i)
    1336             :   {
    1337          84 :     gel(G, k) = gerepileuptoint(av, r);
    1338          84 :     gel(E, k) = gen_1;
    1339             :   }
    1340       55930 :   return mkmat2(G,E);
    1341             : }
    1342             : 
    1343             : /* assume pr has degree 1 and coprime to Q_denom(x) */
    1344             : static GEN
    1345        4755 : to_Fp_coprime(GEN nf, GEN x, GEN modpr)
    1346             : {
    1347        4755 :   GEN d, r, p = modpr_get_p(modpr);
    1348        4755 :   x = nf_to_scalar_or_basis(nf,x);
    1349        4755 :   if (typ(x) != t_COL) return Rg_to_Fp(x,p);
    1350        4447 :   x = Q_remove_denom(x, &d);
    1351        4447 :   r = zk_to_Fq(x, modpr);
    1352        4447 :   if (d) r = Fp_div(r, d, p);
    1353        4447 :   return r;
    1354             : }
    1355             : 
    1356             : /* pr coprime to all denominators occurring in x */
    1357             : static GEN
    1358         621 : famat_to_Fp_coprime(GEN nf, GEN x, GEN modpr)
    1359             : {
    1360         621 :   GEN p = modpr_get_p(modpr);
    1361         621 :   GEN t = NULL, g = gel(x,1), e = gel(x,2), q = subiu(p,1);
    1362         621 :   long i, l = lg(g);
    1363        1967 :   for (i = 1; i < l; i++)
    1364             :   {
    1365        1346 :     GEN n = modii(gel(e,i), q);
    1366        1346 :     if (signe(n))
    1367             :     {
    1368        1346 :       GEN h = to_Fp_coprime(nf, gel(g,i), modpr);
    1369        1346 :       h = Fp_pow(h, n, p);
    1370        1346 :       t = t? Fp_mul(t, h, p): h;
    1371             :     }
    1372             :   }
    1373         621 :   return t? modii(t, p): gen_1;
    1374             : }
    1375             : 
    1376             : /* cf famat_to_nf_modideal_coprime, modpr attached to prime of degree 1 */
    1377             : GEN
    1378        4030 : nf_to_Fp_coprime(GEN nf, GEN x, GEN modpr)
    1379             : {
    1380        8060 :   return typ(x)==t_MAT? famat_to_Fp_coprime(nf, x, modpr)
    1381        4030 :                       : to_Fp_coprime(nf, x, modpr);
    1382             : }
    1383             : 
    1384             : static long
    1385      135815 : zk_pvalrem(GEN x, GEN p, GEN *py)
    1386      135815 : { return (typ(x) == t_INT)? Z_pvalrem(x, p, py): ZV_pvalrem(x, p, py); }
    1387             : /* x a QC or Q. Return a ZC or Z, whose content is coprime to Z. Set v, dx
    1388             :  * such that x = p^v (newx / dx); dx = NULL if 1 */
    1389             : static GEN
    1390      263418 : nf_remove_denom_p(GEN nf, GEN x, GEN p, GEN *pdx, long *pv)
    1391             : {
    1392             :   long vcx;
    1393             :   GEN dx;
    1394      263418 :   x = nf_to_scalar_or_basis(nf, x);
    1395      263418 :   x = Q_remove_denom(x, &dx);
    1396      263418 :   if (dx)
    1397             :   {
    1398      170940 :     vcx = - Z_pvalrem(dx, p, &dx);
    1399      170940 :     if (!vcx) vcx = zk_pvalrem(x, p, &x);
    1400      170940 :     if (isint1(dx)) dx = NULL;
    1401             :   }
    1402             :   else
    1403             :   {
    1404       92478 :     vcx = zk_pvalrem(x, p, &x);
    1405       92478 :     dx = NULL;
    1406             :   }
    1407      263418 :   *pv = vcx;
    1408      263418 :   *pdx = dx; return x;
    1409             : }
    1410             : /* x = b^e/p^(e-1) in Z_K; x = 0 mod p/pr^e, (x,pr) = 1. Return NULL
    1411             :  * if p inert (instead of 1) */
    1412             : static GEN
    1413       62069 : p_makecoprime(GEN pr)
    1414             : {
    1415       62069 :   GEN B = pr_get_tau(pr), b;
    1416             :   long i, e;
    1417             : 
    1418       62069 :   if (typ(B) == t_INT) return NULL;
    1419       61929 :   b = gel(B,1); /* B = multiplication table by b */
    1420       61929 :   e = pr_get_e(pr);
    1421       61929 :   if (e == 1) return b;
    1422             :   /* one could also divide (exactly) by p in each iteration */
    1423       17150 :   for (i = 1; i < e; i++) b = ZM_ZC_mul(B, b);
    1424       17150 :   return ZC_Z_divexact(b, powiu(pr_get_p(pr), e-1));
    1425             : }
    1426             : 
    1427             : /* Compute A = prod g[i]^e[i] mod pr^k, assuming (A, pr) = 1.
    1428             :  * Method: modify each g[i] so that it becomes coprime to pr,
    1429             :  * g[i] *= (b/p)^v_pr(g[i]), where b/p = pr^(-1) times something integral
    1430             :  * and prime to p; globally, we multiply by (b/p)^v_pr(A) = 1.
    1431             :  * Optimizations:
    1432             :  * 1) remove all powers of p from contents, and consider extra generator p^vp;
    1433             :  * modified as p * (b/p)^e = b^e / p^(e-1)
    1434             :  * 2) remove denominators, coprime to p, by multiplying by inverse mod prk\cap Z
    1435             :  *
    1436             :  * EX = multiple of exponent of (O_K / pr^k)^* used to reduce the product in
    1437             :  * case the e[i] are large */
    1438             : GEN
    1439      111633 : famat_makecoprime(GEN nf, GEN g, GEN e, GEN pr, GEN prk, GEN EX)
    1440             : {
    1441      111633 :   GEN G, E, t, vp = NULL, p = pr_get_p(pr), prkZ = gcoeff(prk, 1,1);
    1442      111633 :   long i, l = lg(g);
    1443             : 
    1444      111633 :   G = cgetg(l+1, t_VEC);
    1445      111633 :   E = cgetg(l+1, t_VEC); /* l+1: room for "modified p" */
    1446      375051 :   for (i=1; i < l; i++)
    1447             :   {
    1448             :     long vcx;
    1449      263418 :     GEN dx, x = nf_remove_denom_p(nf, gel(g,i), p, &dx, &vcx);
    1450      263418 :     if (vcx) /* = v_p(content(g[i])) */
    1451             :     {
    1452      129262 :       GEN a = mulsi(vcx, gel(e,i));
    1453      129262 :       vp = vp? addii(vp, a): a;
    1454             :     }
    1455             :     /* x integral, content coprime to p; dx coprime to p */
    1456      263418 :     if (typ(x) == t_INT)
    1457             :     { /* x coprime to p, hence to pr */
    1458       38595 :       x = modii(x, prkZ);
    1459       38595 :       if (dx) x = Fp_div(x, dx, prkZ);
    1460             :     }
    1461             :     else
    1462             :     {
    1463      224823 :       (void)ZC_nfvalrem(x, pr, &x); /* x *= (b/p)^v_pr(x) */
    1464      224823 :       x = ZC_hnfrem(FpC_red(x,prkZ), prk);
    1465      224823 :       if (dx) x = FpC_Fp_mul(x, Fp_inv(dx,prkZ), prkZ);
    1466             :     }
    1467      263418 :     gel(G,i) = x;
    1468      263418 :     gel(E,i) = gel(e,i);
    1469             :   }
    1470             : 
    1471      111633 :   t = vp? p_makecoprime(pr): NULL;
    1472      111633 :   if (!t)
    1473             :   { /* no need for extra generator */
    1474       49704 :     setlg(G,l);
    1475       49704 :     setlg(E,l);
    1476             :   }
    1477             :   else
    1478             :   {
    1479       61929 :     gel(G,i) = FpC_red(t, prkZ);
    1480       61929 :     gel(E,i) = vp;
    1481             :   }
    1482      111633 :   return famat_to_nf_modideal_coprime(nf, G, E, prk, EX);
    1483             : }
    1484             : 
    1485             : /* prod g[i]^e[i] mod bid, assume (g[i], id) = 1 */
    1486             : GEN
    1487       11011 : famat_to_nf_moddivisor(GEN nf, GEN g, GEN e, GEN bid)
    1488             : {
    1489             :   GEN t, cyc;
    1490       11011 :   if (lg(g) == 1) return gen_1;
    1491       11011 :   cyc = bid_get_cyc(bid);
    1492       11011 :   if (lg(cyc) == 1)
    1493           0 :     t = gen_1;
    1494             :   else
    1495       11011 :     t = famat_to_nf_modideal_coprime(nf, g, e, bid_get_ideal(bid), gel(cyc,1));
    1496       11011 :   return set_sign_mod_divisor(nf, mkmat2(g,e), t, bid_get_sarch(bid));
    1497             : }
    1498             : 
    1499             : GEN
    1500      196280 : vecmul(GEN x, GEN y)
    1501             : {
    1502      196280 :   if (is_scalar_t(typ(x))) return gmul(x,y);
    1503       17619 :   pari_APPLY_same(vecmul(gel(x,i), gel(y,i)))
    1504             : }
    1505             : 
    1506             : GEN
    1507           0 : vecinv(GEN x)
    1508             : {
    1509           0 :   if (is_scalar_t(typ(x))) return ginv(x);
    1510           0 :   pari_APPLY_same(vecinv(gel(x,i)))
    1511             : }
    1512             : 
    1513             : GEN
    1514       15729 : vecpow(GEN x, GEN n)
    1515             : {
    1516       15729 :   if (is_scalar_t(typ(x))) return powgi(x,n);
    1517        4270 :   pari_APPLY_same(vecpow(gel(x,i), n))
    1518             : }
    1519             : 
    1520             : GEN
    1521         903 : vecdiv(GEN x, GEN y)
    1522             : {
    1523         903 :   if (is_scalar_t(typ(x))) return gdiv(x,y);
    1524         301 :   pari_APPLY_same(vecdiv(gel(x,i), gel(y,i)))
    1525             : }
    1526             : 
    1527             : /* A ideal as a square t_MAT */
    1528             : static GEN
    1529      197404 : idealmulelt(GEN nf, GEN x, GEN A)
    1530             : {
    1531             :   long i, lx;
    1532             :   GEN dx, dA, D;
    1533      197404 :   if (lg(A) == 1) return cgetg(1, t_MAT);
    1534      197404 :   x = nf_to_scalar_or_basis(nf,x);
    1535      197404 :   if (typ(x) != t_COL)
    1536       68081 :     return isintzero(x)? cgetg(1,t_MAT): RgM_Rg_mul(A, Q_abs_shallow(x));
    1537      129323 :   x = Q_remove_denom(x, &dx);
    1538      129323 :   A = Q_remove_denom(A, &dA);
    1539      129323 :   x = zk_multable(nf, x);
    1540      129323 :   D = mulii(zkmultable_capZ(x), gcoeff(A,1,1));
    1541      129323 :   x = zkC_multable_mul(A, x);
    1542      129323 :   settyp(x, t_MAT); lx = lg(x);
    1543             :   /* x may contain scalars (at most 1 since the ideal is non-0)*/
    1544      447365 :   for (i=1; i<lx; i++)
    1545      326553 :     if (typ(gel(x,i)) == t_INT)
    1546             :     {
    1547        8511 :       if (i > 1) swap(gel(x,1), gel(x,i)); /* help HNF */
    1548        8511 :       gel(x,1) = scalarcol_shallow(gel(x,1), lx-1);
    1549        8511 :       break;
    1550             :     }
    1551      129323 :   x = ZM_hnfmodid(x, D);
    1552      129323 :   dx = mul_denom(dx,dA);
    1553      129323 :   return dx? gdiv(x,dx): x;
    1554             : }
    1555             : 
    1556             : /* nf a true nf, tx <= ty */
    1557             : static GEN
    1558      644675 : idealmul_aux(GEN nf, GEN x, GEN y, long tx, long ty)
    1559             : {
    1560             :   GEN z, cx, cy;
    1561      644675 :   switch(tx)
    1562             :   {
    1563             :     case id_PRINCIPAL:
    1564      247138 :       switch(ty)
    1565             :       {
    1566             :         case id_PRINCIPAL:
    1567       49608 :           return idealhnf_principal(nf, nfmul(nf,x,y));
    1568             :         case id_PRIME:
    1569             :         {
    1570         126 :           GEN p = pr_get_p(y), pi = pr_get_gen(y), cx;
    1571         126 :           if (pr_is_inert(y)) return RgM_Rg_mul(idealhnf_principal(nf,x),p);
    1572             : 
    1573          42 :           x = nf_to_scalar_or_basis(nf, x);
    1574          42 :           switch(typ(x))
    1575             :           {
    1576             :             case t_INT:
    1577          28 :               if (!signe(x)) return cgetg(1,t_MAT);
    1578          28 :               return ZM_Z_mul(pr_hnf(nf,y), absi_shallow(x));
    1579             :             case t_FRAC:
    1580           7 :               return RgM_Rg_mul(pr_hnf(nf,y), Q_abs_shallow(x));
    1581             :           }
    1582             :           /* t_COL */
    1583           7 :           x = Q_primitive_part(x, &cx);
    1584           7 :           x = zk_multable(nf, x);
    1585           7 :           z = shallowconcat(ZM_Z_mul(x,p), ZM_ZC_mul(x,pi));
    1586           7 :           z = ZM_hnfmodid(z, mulii(p, zkmultable_capZ(x)));
    1587           7 :           return cx? ZM_Q_mul(z, cx): z;
    1588             :         }
    1589             :         default: /* id_MAT */
    1590      197404 :           return idealmulelt(nf, x,y);
    1591             :       }
    1592             :     case id_PRIME:
    1593      322733 :       if (ty==id_PRIME)
    1594      301112 :       { y = pr_hnf(nf,y); cy = NULL; }
    1595             :       else
    1596       21621 :         y = Q_primitive_part(y, &cy);
    1597      322733 :       y = idealHNF_mul_two(nf,y,x);
    1598      322733 :       return cy? ZM_Q_mul(y,cy): y;
    1599             : 
    1600             :     default: /* id_MAT */
    1601             :     {
    1602       74804 :       long N = nf_get_degree(nf);
    1603       74804 :       if (lg(x)-1 != N || lg(y)-1 != N) pari_err_DIM("idealmul");
    1604       74790 :       x = Q_primitive_part(x, &cx);
    1605       74790 :       y = Q_primitive_part(y, &cy); cx = mul_content(cx,cy);
    1606       74790 :       y = idealHNF_mul(nf,x,y);
    1607       74790 :       return cx? ZM_Q_mul(y,cx): y;
    1608             :     }
    1609             :   }
    1610             : }
    1611             : 
    1612             : /* output the ideal product ix.iy */
    1613             : GEN
    1614      644675 : idealmul(GEN nf, GEN x, GEN y)
    1615             : {
    1616             :   pari_sp av;
    1617             :   GEN res, ax, ay, z;
    1618      644675 :   long tx = idealtyp(&x,&ax);
    1619      644675 :   long ty = idealtyp(&y,&ay), f;
    1620      644675 :   if (tx>ty) { swap(ax,ay); swap(x,y); lswap(tx,ty); }
    1621      644675 :   f = (ax||ay); res = f? cgetg(3,t_VEC): NULL; /*product is an extended ideal*/
    1622      644675 :   av = avma;
    1623      644675 :   z = gerepileupto(av, idealmul_aux(checknf(nf), x,y, tx,ty));
    1624      644661 :   if (!f) return z;
    1625       42371 :   if (ax && ay)
    1626       40622 :     ax = ext_mul(nf, ax, ay);
    1627             :   else
    1628        1749 :     ax = gcopy(ax? ax: ay);
    1629       42371 :   gel(res,1) = z; gel(res,2) = ax; return res;
    1630             : }
    1631             : 
    1632             : /* Return x, integral in 2-elt form, such that pr^2 = c * x. cf idealpowprime
    1633             :  * nf = true nf */
    1634             : static GEN
    1635       44260 : idealsqrprime(GEN nf, GEN pr, GEN *pc)
    1636             : {
    1637       44260 :   GEN p = pr_get_p(pr), q, gen;
    1638       44260 :   long e = pr_get_e(pr), f = pr_get_f(pr);
    1639             : 
    1640       44260 :   q = (e == 1)? sqri(p): p;
    1641       44260 :   if (e <= 2 && e * f == nf_get_degree(nf))
    1642             :   { /* pr^e = (p) */
    1643        9947 :     *pc = q;
    1644        9947 :     return mkvec2(gen_1,gen_0);
    1645             :   }
    1646       34313 :   gen = nfsqr(nf, pr_get_gen(pr));
    1647       34313 :   gen = FpC_red(gen, q);
    1648       34313 :   *pc = NULL;
    1649       34313 :   return mkvec2(q, gen);
    1650             : }
    1651             : /* cf idealpow_aux */
    1652             : static GEN
    1653       61448 : idealsqr_aux(GEN nf, GEN x, long tx)
    1654             : {
    1655       61448 :   GEN T = nf_get_pol(nf), m, cx, a, alpha;
    1656       61448 :   long N = degpol(T);
    1657       61448 :   switch(tx)
    1658             :   {
    1659             :     case id_PRINCIPAL:
    1660          77 :       return idealhnf_principal(nf, nfsqr(nf,x));
    1661             :     case id_PRIME:
    1662       23029 :       if (pr_is_inert(x)) return scalarmat(sqri(gel(x,1)), N);
    1663       22861 :       x = idealsqrprime(nf, x, &cx);
    1664       22861 :       x = idealhnf_two(nf,x);
    1665       22861 :       return cx? ZM_Z_mul(x, cx): x;
    1666             :     default:
    1667       38342 :       x = Q_primitive_part(x, &cx);
    1668       38342 :       a = mat_ideal_two_elt(nf,x); alpha = gel(a,2); a = gel(a,1);
    1669       38342 :       alpha = nfsqr(nf,alpha);
    1670       38342 :       m = zk_scalar_or_multable(nf, alpha);
    1671       38342 :       if (typ(m) == t_INT) {
    1672        1435 :         x = gcdii(sqri(a), m);
    1673        1435 :         if (cx) x = gmul(x, gsqr(cx));
    1674        1435 :         x = scalarmat(x, N);
    1675             :       }
    1676             :       else
    1677             :       {
    1678       36907 :         x = ZM_hnfmodid(m, gcdii(sqri(a), zkmultable_capZ(m)));
    1679       36907 :         if (cx) cx = gsqr(cx);
    1680       36907 :         if (cx) x = ZM_Q_mul(x, cx);
    1681             :       }
    1682       38342 :       return x;
    1683             :   }
    1684             : }
    1685             : GEN
    1686       61448 : idealsqr(GEN nf, GEN x)
    1687             : {
    1688             :   pari_sp av;
    1689             :   GEN res, ax, z;
    1690       61448 :   long tx = idealtyp(&x,&ax);
    1691       61448 :   res = ax? cgetg(3,t_VEC): NULL; /*product is an extended ideal*/
    1692       61448 :   av = avma;
    1693       61448 :   z = gerepileupto(av, idealsqr_aux(checknf(nf), x, tx));
    1694       61448 :   if (!ax) return z;
    1695       60713 :   gel(res,1) = z;
    1696       60713 :   gel(res,2) = ext_sqr(nf, ax); return res;
    1697             : }
    1698             : 
    1699             : /* norm of an ideal */
    1700             : GEN
    1701        5271 : idealnorm(GEN nf, GEN x)
    1702             : {
    1703             :   pari_sp av;
    1704             :   GEN y, T;
    1705             :   long tx;
    1706             : 
    1707        5271 :   switch(idealtyp(&x,&y))
    1708             :   {
    1709         175 :     case id_PRIME: return pr_norm(x);
    1710        3038 :     case id_MAT: return RgM_det_triangular(x);
    1711             :   }
    1712             :   /* id_PRINCIPAL */
    1713        2058 :   nf = checknf(nf); T = nf_get_pol(nf); av = avma;
    1714        2058 :   x = nf_to_scalar_or_alg(nf, x);
    1715        2058 :   x = (typ(x) == t_POL)? RgXQ_norm(x, T): gpowgs(x, degpol(T));
    1716        2058 :   tx = typ(x);
    1717        2058 :   if (tx == t_INT) return gerepileuptoint(av, absi(x));
    1718         532 :   if (tx != t_FRAC) pari_err_TYPE("idealnorm",x);
    1719         532 :   return gerepileupto(av, Q_abs(x));
    1720             : }
    1721             : 
    1722             : /* I^(-1) = { x \in K, Tr(x D^(-1) I) \in Z }, D different of K/Q
    1723             :  *
    1724             :  * nf[5][6] = pp( D^(-1) ) = pp( HNF( T^(-1) ) ), T = (Tr(wi wj))
    1725             :  * nf[5][7] = same in 2-elt form.
    1726             :  * Assume I integral. Return the integral ideal (I\cap Z) I^(-1) */
    1727             : GEN
    1728      190829 : idealHNF_inv_Z(GEN nf, GEN I)
    1729             : {
    1730      190829 :   GEN J, dual, IZ = gcoeff(I,1,1); /* I \cap Z */
    1731      190829 :   if (isint1(IZ)) return matid(lg(I)-1);
    1732      179755 :   J = idealHNF_mul(nf,I, gmael(nf,5,7));
    1733             :  /* I in HNF, hence easily inverted; multiply by IZ to get integer coeffs
    1734             :   * missing content cancels while solving the linear equation */
    1735      179755 :   dual = shallowtrans( hnf_divscale(J, gmael(nf,5,6), IZ) );
    1736      179755 :   return ZM_hnfmodid(dual, IZ);
    1737             : }
    1738             : /* I HNF with rational coefficients (denominator d). */
    1739             : GEN
    1740       56111 : idealHNF_inv(GEN nf, GEN I)
    1741             : {
    1742       56111 :   GEN J, IQ = gcoeff(I,1,1); /* I \cap Q; d IQ = dI \cap Z */
    1743       56111 :   J = idealHNF_inv_Z(nf, Q_remove_denom(I, NULL)); /* = (dI)^(-1) * (d IQ) */
    1744       56111 :   return equali1(IQ)? J: RgM_Rg_div(J, IQ);
    1745             : }
    1746             : 
    1747             : /* return p * P^(-1)  [integral] */
    1748             : GEN
    1749       24502 : pr_inv_p(GEN pr)
    1750             : {
    1751       24502 :   if (pr_is_inert(pr)) return matid(pr_get_f(pr));
    1752       24012 :   return ZM_hnfmodid(pr_get_tau(pr), pr_get_p(pr));
    1753             : }
    1754             : GEN
    1755        3572 : pr_inv(GEN pr)
    1756             : {
    1757        3572 :   GEN p = pr_get_p(pr);
    1758        3572 :   if (pr_is_inert(pr)) return scalarmat(ginv(p), pr_get_f(pr));
    1759        3222 :   return RgM_Rg_div(ZM_hnfmodid(pr_get_tau(pr),p), p);
    1760             : }
    1761             : 
    1762             : GEN
    1763       96341 : idealinv(GEN nf, GEN x)
    1764             : {
    1765             :   GEN res, ax;
    1766             :   pari_sp av;
    1767       96341 :   long tx = idealtyp(&x,&ax), N;
    1768             : 
    1769       96341 :   res = ax? cgetg(3,t_VEC): NULL;
    1770       96341 :   nf = checknf(nf); av = avma;
    1771       96341 :   N = nf_get_degree(nf);
    1772       96341 :   switch (tx)
    1773             :   {
    1774             :     case id_MAT:
    1775       50574 :       if (lg(x)-1 != N) pari_err_DIM("idealinv");
    1776       50574 :       x = idealHNF_inv(nf,x); break;
    1777             :     case id_PRINCIPAL:
    1778       43063 :       x = nf_to_scalar_or_basis(nf, x);
    1779       43063 :       if (typ(x) != t_COL)
    1780       43021 :         x = idealhnf_principal(nf,ginv(x));
    1781             :       else
    1782             :       { /* nfinv + idealhnf where we already know (x) \cap Z */
    1783             :         GEN c, d;
    1784          42 :         x = Q_remove_denom(x, &c);
    1785          42 :         x = zk_inv(nf, x);
    1786          42 :         x = Q_remove_denom(x, &d); /* true inverse is c/d * x */
    1787          42 :         if (!d) /* x and x^(-1) integral => x a unit */
    1788           7 :           x = scalarmat_shallow(c? c: gen_1, N);
    1789             :         else
    1790             :         {
    1791          35 :           c = c? gdiv(c,d): ginv(d);
    1792          35 :           x = zk_multable(nf, x);
    1793          35 :           x = ZM_Q_mul(ZM_hnfmodid(x,d), c);
    1794             :         }
    1795             :       }
    1796       43063 :       break;
    1797             :     case id_PRIME:
    1798        2704 :       x = pr_inv(x); break;
    1799             :   }
    1800       96341 :   x = gerepileupto(av,x); if (!ax) return x;
    1801       10920 :   gel(res,1) = x;
    1802       10920 :   gel(res,2) = ext_inv(nf, ax); return res;
    1803             : }
    1804             : 
    1805             : /* write x = A/B, A,B coprime integral ideals */
    1806             : GEN
    1807       37742 : idealnumden(GEN nf, GEN x)
    1808             : {
    1809       37742 :   pari_sp av = avma;
    1810             :   GEN x0, ax, c, d, A, B, J;
    1811       37742 :   long tx = idealtyp(&x,&ax);
    1812       37742 :   nf = checknf(nf);
    1813       37742 :   switch (tx)
    1814             :   {
    1815             :     case id_PRIME:
    1816           7 :       retmkvec2(idealhnf(nf, x), gen_1);
    1817             :     case id_PRINCIPAL:
    1818             :     {
    1819             :       GEN xZ, mx;
    1820        2233 :       x = nf_to_scalar_or_basis(nf, x);
    1821        2233 :       switch(typ(x))
    1822             :       {
    1823          98 :         case t_INT: return gerepilecopy(av, mkvec2(absi(x),gen_1));
    1824          14 :         case t_FRAC:return gerepilecopy(av, mkvec2(absi(gel(x,1)), gel(x,2)));
    1825             :       }
    1826             :       /* t_COL */
    1827        2121 :       x = Q_remove_denom(x, &d);
    1828        2121 :       if (!d) return gerepilecopy(av, mkvec2(idealhnf(nf, x), gen_1));
    1829          35 :       mx = zk_multable(nf, x);
    1830          35 :       xZ = zkmultable_capZ(mx);
    1831          35 :       x = ZM_hnfmodid(mx, xZ); /* principal ideal (x) */
    1832          35 :       x0 = mkvec2(xZ, mx); /* same, for fast multiplication */
    1833          35 :       break;
    1834             :     }
    1835             :     default: /* id_MAT */
    1836             :     {
    1837       35502 :       long n = lg(x)-1;
    1838       35502 :       if (n == 0) return mkvec2(gen_0, gen_1);
    1839       35502 :       if (n != nf_get_degree(nf)) pari_err_DIM("idealnumden");
    1840       35502 :       x0 = x = Q_remove_denom(x, &d);
    1841       35502 :       if (!d) return gerepilecopy(av, mkvec2(x, gen_1));
    1842          14 :       break;
    1843             :     }
    1844             :   }
    1845          49 :   J = hnfmodid(x, d); /* = d/B */
    1846          49 :   c = gcoeff(J,1,1); /* (d/B) \cap Z, divides d */
    1847          49 :   B = idealHNF_inv_Z(nf, J); /* (d/B \cap Z) B/d */
    1848          49 :   if (!equalii(c,d)) B = ZM_Z_mul(B, diviiexact(d,c)); /* = B ! */
    1849          49 :   A = idealHNF_mul(nf, B, x0); /* d * (original x) * B = d A */
    1850          49 :   A = ZM_Z_divexact(A, d); /* = A ! */
    1851          49 :   return gerepilecopy(av, mkvec2(A, B));
    1852             : }
    1853             : 
    1854             : /* Return x, integral in 2-elt form, such that pr^n = c * x. Assume n != 0.
    1855             :  * nf = true nf */
    1856             : static GEN
    1857       89284 : idealpowprime(GEN nf, GEN pr, GEN n, GEN *pc)
    1858             : {
    1859       89284 :   GEN p = pr_get_p(pr), q, gen;
    1860             : 
    1861       89284 :   *pc = NULL;
    1862       89284 :   if (is_pm1(n)) /* n = 1 special cased for efficiency */
    1863             :   {
    1864       50888 :     q = p;
    1865       50888 :     if (typ(pr_get_tau(pr)) == t_INT) /* inert */
    1866             :     {
    1867           0 :       *pc = (signe(n) >= 0)? p: ginv(p);
    1868           0 :       return mkvec2(gen_1,gen_0);
    1869             :     }
    1870       50888 :     if (signe(n) >= 0) gen = pr_get_gen(pr);
    1871             :     else
    1872             :     {
    1873        8176 :       gen = pr_get_tau(pr); /* possibly t_MAT */
    1874        8176 :       *pc = ginv(p);
    1875             :     }
    1876             :   }
    1877       38396 :   else if (equalis(n,2)) return idealsqrprime(nf, pr, pc);
    1878             :   else
    1879             :   {
    1880       16997 :     long e = pr_get_e(pr), f = pr_get_f(pr);
    1881       16997 :     GEN r, m = truedvmdis(n, e, &r);
    1882       16997 :     if (e * f == nf_get_degree(nf))
    1883             :     { /* pr^e = (p) */
    1884        7812 :       if (signe(m)) *pc = powii(p,m);
    1885        7812 :       if (!signe(r)) return mkvec2(gen_1,gen_0);
    1886        3206 :       q = p;
    1887        3206 :       gen = nfpow(nf, pr_get_gen(pr), r);
    1888             :     }
    1889             :     else
    1890             :     {
    1891        9185 :       m = absi(m);
    1892        9185 :       if (signe(r)) m = addiu(m,1);
    1893        9185 :       q = powii(p,m); /* m = ceil(|n|/e) */
    1894        9185 :       if (signe(n) >= 0) gen = nfpow(nf, pr_get_gen(pr), n);
    1895             :       else
    1896             :       {
    1897        2261 :         gen = pr_get_tau(pr);
    1898        2261 :         if (typ(gen) == t_MAT) gen = gel(gen,1);
    1899        2261 :         n = negi(n);
    1900        2261 :         gen = ZC_Z_divexact(nfpow(nf, gen, n), powii(p, subii(n,m)));
    1901        2261 :         *pc = ginv(q);
    1902             :       }
    1903             :     }
    1904       12391 :     gen = FpC_red(gen, q);
    1905             :   }
    1906       63279 :   return mkvec2(q, gen);
    1907             : }
    1908             : 
    1909             : /* x * pr^n. Assume x in HNF or scalar (possibly non-integral) */
    1910             : GEN
    1911       67753 : idealmulpowprime(GEN nf, GEN x, GEN pr, GEN n)
    1912             : {
    1913             :   GEN c, cx, y;
    1914             :   long N;
    1915             : 
    1916       67753 :   nf = checknf(nf);
    1917       67753 :   N = nf_get_degree(nf);
    1918       67753 :   if (!signe(n)) return typ(x) == t_MAT? x: scalarmat_shallow(x, N);
    1919             : 
    1920             :   /* inert, special cased for efficiency */
    1921       67641 :   if (pr_is_inert(pr))
    1922             :   {
    1923        5565 :     GEN q = powii(pr_get_p(pr), n);
    1924        5565 :     return typ(x) == t_MAT? RgM_Rg_mul(x,q): scalarmat_shallow(gmul(x,q), N);
    1925             :   }
    1926             : 
    1927       62076 :   y = idealpowprime(nf, pr, n, &c);
    1928       62076 :   if (typ(x) == t_MAT)
    1929       60704 :   { x = Q_primitive_part(x, &cx); if (is_pm1(gcoeff(x,1,1))) x = NULL; }
    1930             :   else
    1931        1372 :   { cx = x; x = NULL; }
    1932       62076 :   cx = mul_content(c,cx);
    1933       62076 :   if (x)
    1934       34076 :     x = idealHNF_mul_two(nf,x,y);
    1935             :   else
    1936       28000 :     x = idealhnf_two(nf,y);
    1937       62076 :   if (cx) x = ZM_Q_mul(x,cx);
    1938       62076 :   return x;
    1939             : }
    1940             : GEN
    1941       13545 : idealdivpowprime(GEN nf, GEN x, GEN pr, GEN n)
    1942             : {
    1943       13545 :   return idealmulpowprime(nf,x,pr, negi(n));
    1944             : }
    1945             : 
    1946             : /* nf = true nf */
    1947             : static GEN
    1948      178746 : idealpow_aux(GEN nf, GEN x, long tx, GEN n)
    1949             : {
    1950      178746 :   GEN T = nf_get_pol(nf), m, cx, n1, a, alpha;
    1951      178746 :   long N = degpol(T), s = signe(n);
    1952      178746 :   if (!s) return matid(N);
    1953      173176 :   switch(tx)
    1954             :   {
    1955             :     case id_PRINCIPAL:
    1956           0 :       return idealhnf_principal(nf, nfpow(nf,x,n));
    1957             :     case id_PRIME:
    1958       71252 :       if (pr_is_inert(x)) return scalarmat(powii(gel(x,1), n), N);
    1959       27208 :       x = idealpowprime(nf, x, n, &cx);
    1960       27208 :       x = idealhnf_two(nf,x);
    1961       27208 :       return cx? ZM_Q_mul(x, cx): x;
    1962             :     default:
    1963      101924 :       if (is_pm1(n)) return (s < 0)? idealinv(nf, x): gcopy(x);
    1964       55757 :       n1 = (s < 0)? negi(n): n;
    1965             : 
    1966       55757 :       x = Q_primitive_part(x, &cx);
    1967       55757 :       a = mat_ideal_two_elt(nf,x); alpha = gel(a,2); a = gel(a,1);
    1968       55757 :       alpha = nfpow(nf,alpha,n1);
    1969       55757 :       m = zk_scalar_or_multable(nf, alpha);
    1970       55757 :       if (typ(m) == t_INT) {
    1971         189 :         x = gcdii(powii(a,n1), m);
    1972         189 :         if (s<0) x = ginv(x);
    1973         189 :         if (cx) x = gmul(x, powgi(cx,n));
    1974         189 :         x = scalarmat(x, N);
    1975             :       }
    1976             :       else
    1977             :       {
    1978       55568 :         x = ZM_hnfmodid(m, gcdii(powii(a,n1), zkmultable_capZ(m)));
    1979       55568 :         if (cx) cx = powgi(cx,n);
    1980       55568 :         if (s<0) {
    1981           7 :           GEN xZ = gcoeff(x,1,1);
    1982           7 :           cx = cx ? gdiv(cx, xZ): ginv(xZ);
    1983           7 :           x = idealHNF_inv_Z(nf,x);
    1984             :         }
    1985       55568 :         if (cx) x = ZM_Q_mul(x, cx);
    1986             :       }
    1987       55757 :       return x;
    1988             :   }
    1989             : }
    1990             : 
    1991             : /* raise the ideal x to the power n (in Z) */
    1992             : GEN
    1993      178746 : idealpow(GEN nf, GEN x, GEN n)
    1994             : {
    1995             :   pari_sp av;
    1996             :   long tx;
    1997             :   GEN res, ax;
    1998             : 
    1999      178746 :   if (typ(n) != t_INT) pari_err_TYPE("idealpow",n);
    2000      178746 :   tx = idealtyp(&x,&ax);
    2001      178746 :   res = ax? cgetg(3,t_VEC): NULL;
    2002      178746 :   av = avma;
    2003      178746 :   x = gerepileupto(av, idealpow_aux(checknf(nf), x, tx, n));
    2004      178746 :   if (!ax) return x;
    2005        1167 :   ax = ext_pow(nf, ax, n);
    2006        1167 :   gel(res,1) = x;
    2007        1167 :   gel(res,2) = ax;
    2008        1167 :   return res;
    2009             : }
    2010             : 
    2011             : /* Return ideal^e in number field nf. e is a C integer. */
    2012             : GEN
    2013       21245 : idealpows(GEN nf, GEN ideal, long e)
    2014             : {
    2015       21245 :   long court[] = {evaltyp(t_INT) | _evallg(3),0,0};
    2016       21245 :   affsi(e,court); return idealpow(nf,ideal,court);
    2017             : }
    2018             : 
    2019             : static GEN
    2020       42392 : _idealmulred(GEN nf, GEN x, GEN y)
    2021       42392 : { return idealred(nf,idealmul(nf,x,y)); }
    2022             : static GEN
    2023       60727 : _idealsqrred(GEN nf, GEN x)
    2024       60727 : { return idealred(nf,idealsqr(nf,x)); }
    2025             : static GEN
    2026       25832 : _mul(void *data, GEN x, GEN y) { return _idealmulred((GEN)data,x,y); }
    2027             : static GEN
    2028       60727 : _sqr(void *data, GEN x) { return _idealsqrred((GEN)data, x); }
    2029             : 
    2030             : /* compute x^n (x ideal, n integer), reducing along the way */
    2031             : GEN
    2032       58591 : idealpowred(GEN nf, GEN x, GEN n)
    2033             : {
    2034       58591 :   pari_sp av = avma;
    2035             :   long s;
    2036             :   GEN y;
    2037             : 
    2038       58591 :   if (typ(n) != t_INT) pari_err_TYPE("idealpowred",n);
    2039       58591 :   s = signe(n); if (s == 0) return idealpow(nf,x,n);
    2040       57424 :   y = gen_pow(x, n, (void*)nf, &_sqr, &_mul);
    2041             : 
    2042       57424 :   if (s < 0) y = idealinv(nf,y);
    2043       57424 :   if (s < 0 || is_pm1(n)) y = idealred(nf,y);
    2044       57424 :   return gerepileupto(av,y);
    2045             : }
    2046             : 
    2047             : GEN
    2048       16560 : idealmulred(GEN nf, GEN x, GEN y)
    2049             : {
    2050       16560 :   pari_sp av = avma;
    2051       16560 :   return gerepileupto(av, _idealmulred(nf,x,y));
    2052             : }
    2053             : 
    2054             : long
    2055          91 : isideal(GEN nf,GEN x)
    2056             : {
    2057          91 :   long N, i, j, lx, tx = typ(x);
    2058             :   pari_sp av;
    2059             :   GEN T, xZ;
    2060             : 
    2061          91 :   nf = checknf(nf); T = nf_get_pol(nf); lx = lg(x);
    2062          91 :   if (tx==t_VEC && lx==3) { x = gel(x,1); tx = typ(x); lx = lg(x); }
    2063          91 :   switch(tx)
    2064             :   {
    2065          14 :     case t_INT: case t_FRAC: return 1;
    2066           7 :     case t_POL: return varn(x) == varn(T);
    2067           7 :     case t_POLMOD: return RgX_equal_var(T, gel(x,1));
    2068          14 :     case t_VEC: return get_prid(x)? 1 : 0;
    2069          42 :     case t_MAT: break;
    2070           7 :     default: return 0;
    2071             :   }
    2072          42 :   N = degpol(T);
    2073          42 :   if (lx-1 != N) return (lx == 1);
    2074          28 :   if (nbrows(x) != N) return 0;
    2075             : 
    2076          28 :   av = avma; x = Q_primpart(x);
    2077          28 :   if (!ZM_ishnf(x)) return 0;
    2078          14 :   xZ = gcoeff(x,1,1);
    2079          21 :   for (j=2; j<=N; j++)
    2080          14 :     if (!dvdii(xZ, gcoeff(x,j,j))) { avma = av; return 0; }
    2081          14 :   for (i=2; i<=N; i++)
    2082          14 :     for (j=2; j<=N; j++)
    2083           7 :       if (! hnf_invimage(x, zk_ei_mul(nf,gel(x,i),j))) { avma = av; return 0; }
    2084           7 :   avma=av; return 1;
    2085             : }
    2086             : 
    2087             : GEN
    2088       21056 : idealdiv(GEN nf, GEN x, GEN y)
    2089             : {
    2090       21056 :   pari_sp av = avma, tetpil;
    2091       21056 :   GEN z = idealinv(nf,y);
    2092       21056 :   tetpil = avma; return gerepile(av,tetpil, idealmul(nf,x,z));
    2093             : }
    2094             : 
    2095             : /* This routine computes the quotient x/y of two ideals in the number field nf.
    2096             :  * It assumes that the quotient is an integral ideal.  The idea is to find an
    2097             :  * ideal z dividing y such that gcd(Nx/Nz, Nz) = 1.  Then
    2098             :  *
    2099             :  *   x + (Nx/Nz)    x
    2100             :  *   ----------- = ---
    2101             :  *   y + (Ny/Nz)    y
    2102             :  *
    2103             :  * Proof: we can assume x and y are integral. Let p be any prime ideal
    2104             :  *
    2105             :  * If p | Nz, then it divides neither Nx/Nz nor Ny/Nz (since Nx/Nz is the
    2106             :  * product of the integers N(x/y) and N(y/z)).  Both the numerator and the
    2107             :  * denominator on the left will be coprime to p.  So will x/y, since x/y is
    2108             :  * assumed integral and its norm N(x/y) is coprime to p.
    2109             :  *
    2110             :  * If instead p does not divide Nz, then v_p (Nx/Nz) = v_p (Nx) >= v_p(x).
    2111             :  * Hence v_p (x + Nx/Nz) = v_p(x).  Likewise for the denominators.  QED.
    2112             :  *
    2113             :  *                Peter Montgomery.  July, 1994. */
    2114             : static void
    2115           7 : err_divexact(GEN x, GEN y)
    2116           7 : { pari_err_DOMAIN("idealdivexact","denominator(x/y)", "!=",
    2117           0 :                   gen_1,mkvec2(x,y)); }
    2118             : GEN
    2119         854 : idealdivexact(GEN nf, GEN x0, GEN y0)
    2120             : {
    2121         854 :   pari_sp av = avma;
    2122             :   GEN x, y, xZ, yZ, Nx, Ny, Nz, cy, q, r;
    2123             : 
    2124         854 :   nf = checknf(nf);
    2125         854 :   x = idealhnf_shallow(nf, x0);
    2126         854 :   y = idealhnf_shallow(nf, y0);
    2127         854 :   if (lg(y) == 1) pari_err_INV("idealdivexact", y0);
    2128         847 :   if (lg(x) == 1) { avma = av; return cgetg(1, t_MAT); } /* numerator is zero */
    2129         847 :   y = Q_primitive_part(y, &cy);
    2130         847 :   if (cy) x = RgM_Rg_div(x,cy);
    2131         847 :   xZ = gcoeff(x,1,1); if (typ(xZ) != t_INT) err_divexact(x,y);
    2132         840 :   yZ = gcoeff(y,1,1); if (isint1(yZ)) return gerepilecopy(av, x);
    2133         231 :   Nx = idealnorm(nf,x);
    2134         231 :   Ny = idealnorm(nf,y);
    2135         231 :   if (typ(Nx) != t_INT) err_divexact(x,y);
    2136         231 :   q = dvmdii(Nx,Ny, &r);
    2137         231 :   if (signe(r)) err_divexact(x,y);
    2138         231 :   if (is_pm1(q)) { avma = av; return matid(nf_get_degree(nf)); }
    2139             :   /* Find a norm Nz | Ny such that gcd(Nx/Nz, Nz) = 1 */
    2140         196 :   for (Nz = Ny;;) /* q = Nx/Nz */
    2141             :   {
    2142         357 :     GEN p1 = gcdii(Nz, q);
    2143         357 :     if (is_pm1(p1)) break;
    2144         161 :     Nz = diviiexact(Nz,p1);
    2145         161 :     q = mulii(q,p1);
    2146         161 :   }
    2147         196 :   xZ = gcoeff(x,1,1); q = gcdii(q, xZ);
    2148         196 :   if (!equalii(xZ,q))
    2149             :   { /* Replace x/y  by  x+(Nx/Nz) / y+(Ny/Nz) */
    2150          42 :     x = ZM_hnfmodid(x, q);
    2151             :     /* y reduced to unit ideal ? */
    2152          42 :     if (Nz == Ny) return gerepileupto(av, x);
    2153             : 
    2154           7 :     yZ = gcoeff(y,1,1); q = gcdii(diviiexact(Ny,Nz), yZ);
    2155           7 :     y = ZM_hnfmodid(y, q);
    2156             :   }
    2157         161 :   yZ = gcoeff(y,1,1);
    2158         161 :   y = idealHNF_mul(nf,x, idealHNF_inv_Z(nf,y));
    2159         161 :   return gerepileupto(av, ZM_Z_divexact(y, yZ));
    2160             : }
    2161             : 
    2162             : GEN
    2163          21 : idealintersect(GEN nf, GEN x, GEN y)
    2164             : {
    2165          21 :   pari_sp av = avma;
    2166             :   long lz, lx, i;
    2167             :   GEN z, dx, dy, xZ, yZ;;
    2168             : 
    2169          21 :   nf = checknf(nf);
    2170          21 :   x = idealhnf_shallow(nf,x);
    2171          21 :   y = idealhnf_shallow(nf,y);
    2172          21 :   if (lg(x) == 1 || lg(y) == 1) { avma = av; return cgetg(1,t_MAT); }
    2173          14 :   x = Q_remove_denom(x, &dx);
    2174          14 :   y = Q_remove_denom(y, &dy);
    2175          14 :   if (dx) y = ZM_Z_mul(y, dx);
    2176          14 :   if (dy) x = ZM_Z_mul(x, dy);
    2177          14 :   xZ = gcoeff(x,1,1);
    2178          14 :   yZ = gcoeff(y,1,1);
    2179          14 :   dx = mul_denom(dx,dy);
    2180          14 :   z = ZM_lll(shallowconcat(x,y), 0.99, LLL_KER); lz = lg(z);
    2181          14 :   lx = lg(x);
    2182          14 :   for (i=1; i<lz; i++) setlg(z[i], lx);
    2183          14 :   z = ZM_hnfmodid(ZM_mul(x,z), lcmii(xZ, yZ));
    2184          14 :   if (dx) z = RgM_Rg_div(z,dx);
    2185          14 :   return gerepileupto(av,z);
    2186             : }
    2187             : 
    2188             : /*******************************************************************/
    2189             : /*                                                                 */
    2190             : /*                      T2-IDEAL REDUCTION                         */
    2191             : /*                                                                 */
    2192             : /*******************************************************************/
    2193             : 
    2194             : static GEN
    2195          21 : chk_vdir(GEN nf, GEN vdir)
    2196             : {
    2197          21 :   long i, l = lg(vdir);
    2198             :   GEN v;
    2199          21 :   if (l != lg(nf_get_roots(nf))) pari_err_DIM("idealred");
    2200          14 :   switch(typ(vdir))
    2201             :   {
    2202           0 :     case t_VECSMALL: return vdir;
    2203          14 :     case t_VEC: break;
    2204           0 :     default: pari_err_TYPE("idealred",vdir);
    2205             :   }
    2206          14 :   v = cgetg(l, t_VECSMALL);
    2207          14 :   for (i = 1; i < l; i++) v[i] = itos(gceil(gel(vdir,i)));
    2208          14 :   return v;
    2209             : }
    2210             : 
    2211             : static void
    2212       26987 : twistG(GEN G, long r1, long i, long v)
    2213             : {
    2214       26987 :   long j, lG = lg(G);
    2215       26987 :   if (i <= r1) {
    2216       23560 :     for (j=1; j<lG; j++) gcoeff(G,i,j) = gmul2n(gcoeff(G,i,j), v);
    2217             :   } else {
    2218        3427 :     long k = (i<<1) - r1;
    2219       18285 :     for (j=1; j<lG; j++)
    2220             :     {
    2221       14858 :       gcoeff(G,k-1,j) = gmul2n(gcoeff(G,k-1,j), v);
    2222       14858 :       gcoeff(G,k  ,j) = gmul2n(gcoeff(G,k  ,j), v);
    2223             :     }
    2224             :   }
    2225       26987 : }
    2226             : 
    2227             : GEN
    2228      162972 : nf_get_Gtwist(GEN nf, GEN vdir)
    2229             : {
    2230             :   long i, l, v, r1;
    2231             :   GEN G;
    2232             : 
    2233      162972 :   if (!vdir) return nf_get_roundG(nf);
    2234        2990 :   if (typ(vdir) == t_MAT)
    2235             :   {
    2236        2969 :     long N = nf_get_degree(nf);
    2237        2969 :     if (lg(vdir) != N+1 || lgcols(vdir) != N+1) pari_err_DIM("idealred");
    2238        2969 :     return vdir;
    2239             :   }
    2240          21 :   vdir = chk_vdir(nf, vdir);
    2241          14 :   G = RgM_shallowcopy(nf_get_G(nf));
    2242          14 :   r1 = nf_get_r1(nf);
    2243          14 :   l = lg(vdir);
    2244          56 :   for (i=1; i<l; i++)
    2245             :   {
    2246          42 :     v = vdir[i]; if (!v) continue;
    2247          42 :     twistG(G, r1, i, v);
    2248             :   }
    2249          14 :   return RM_round_maxrank(G);
    2250             : }
    2251             : GEN
    2252       26945 : nf_get_Gtwist1(GEN nf, long i)
    2253             : {
    2254       26945 :   GEN G = RgM_shallowcopy( nf_get_G(nf) );
    2255       26945 :   long r1 = nf_get_r1(nf);
    2256       26945 :   twistG(G, r1, i, 10);
    2257       26945 :   return RM_round_maxrank(G);
    2258             : }
    2259             : 
    2260             : GEN
    2261       40231 : RM_round_maxrank(GEN G0)
    2262             : {
    2263       40231 :   long e, r = lg(G0)-1;
    2264       40231 :   pari_sp av = avma;
    2265       40231 :   GEN G = G0;
    2266       40231 :   for (e = 4; ; e <<= 1)
    2267             :   {
    2268       40231 :     GEN H = ground(G);
    2269       80462 :     if (ZM_rank(H) == r) return H; /* maximal rank ? */
    2270           0 :     avma = av;
    2271           0 :     G = gmul2n(G0, e);
    2272           0 :   }
    2273             : }
    2274             : 
    2275             : GEN
    2276      162965 : idealred0(GEN nf, GEN I, GEN vdir)
    2277             : {
    2278      162965 :   pari_sp av = avma;
    2279      162965 :   GEN G, aI, IZ, J, y, yZ, my, c1 = NULL;
    2280             :   long N;
    2281             : 
    2282      162965 :   nf = checknf(nf);
    2283      162965 :   N = nf_get_degree(nf);
    2284             :   /* put first for sanity checks, unused when I obviously principal */
    2285      162965 :   G = nf_get_Gtwist(nf, vdir);
    2286      162958 :   switch (idealtyp(&I,&aI))
    2287             :   {
    2288             :     case id_PRIME:
    2289       23221 :       if (pr_is_inert(I)) {
    2290         581 :         if (!aI) { avma = av; return matid(N); }
    2291         581 :         c1 = gel(I,1); I = matid(N);
    2292         581 :         goto END;
    2293             :       }
    2294       22640 :       IZ = pr_get_p(I);
    2295       22640 :       J = pr_inv_p(I);
    2296       22640 :       I = idealhnf_two(nf,I);
    2297       22640 :       break;
    2298             :     case id_MAT:
    2299      139723 :       I = Q_primitive_part(I, &c1);
    2300      139723 :       IZ = gcoeff(I,1,1);
    2301      139723 :       if (is_pm1(IZ))
    2302             :       {
    2303        7987 :         if (!aI) { avma = av; return matid(N); }
    2304        7931 :         goto END;
    2305             :       }
    2306      131736 :       J = idealHNF_inv_Z(nf, I);
    2307      131736 :       break;
    2308             :     default: /* id_PRINCIPAL, silly case */
    2309          14 :       if (gequal0(I)) I = cgetg(1,t_MAT); else { c1 = I; I = matid(N); }
    2310          14 :       if (!aI) return I;
    2311           7 :       goto END;
    2312             :   }
    2313             :   /* now I integral, HNF; and J = (I\cap Z) I^(-1), integral */
    2314      154376 :   y = idealpseudomin(J, G); /* small elt in (I\cap Z)I^(-1), integral */
    2315      154376 :   if (ZV_isscalar(y))
    2316             :   { /* already reduced */
    2317       65080 :     if (!aI) return gerepilecopy(av, I);
    2318       64681 :     goto END;
    2319             :   }
    2320             : 
    2321       89296 :   my = zk_multable(nf, y);
    2322       89296 :   I = ZM_Z_divexact(ZM_mul(my, I), IZ); /* y I / (I\cap Z), integral */
    2323       89296 :   c1 = mul_content(c1, IZ);
    2324       89296 :   my = ZM_gauss(my, col_ei(N,1)); /* y^-1 */
    2325       89296 :   yZ = Q_denom(my); /* (y) \cap Z */
    2326       89296 :   I = hnfmodid(I, yZ);
    2327       89296 :   if (!aI) return gerepileupto(av, I);
    2328       89100 :   c1 = RgC_Rg_mul(my, c1);
    2329             : END:
    2330      162300 :   if (c1) aI = ext_mul(nf, aI,c1);
    2331      162300 :   return gerepilecopy(av, mkvec2(I, aI));
    2332             : }
    2333             : 
    2334             : GEN
    2335           7 : idealmin(GEN nf, GEN x, GEN vdir)
    2336             : {
    2337           7 :   pari_sp av = avma;
    2338             :   GEN y, dx;
    2339           7 :   nf = checknf(nf);
    2340           7 :   switch( idealtyp(&x,&y) )
    2341             :   {
    2342           0 :     case id_PRINCIPAL: return gcopy(x);
    2343           0 :     case id_PRIME: x = pr_hnf(nf,x); break;
    2344           7 :     case id_MAT: if (lg(x) == 1) return gen_0;
    2345             :   }
    2346           7 :   x = Q_remove_denom(x, &dx);
    2347           7 :   y = idealpseudomin(x, nf_get_Gtwist(nf,vdir));
    2348           7 :   if (dx) y = RgC_Rg_div(y, dx);
    2349           7 :   return gerepileupto(av, y);
    2350             : }
    2351             : 
    2352             : /*******************************************************************/
    2353             : /*                                                                 */
    2354             : /*                   APPROXIMATION THEOREM                         */
    2355             : /*                                                                 */
    2356             : /*******************************************************************/
    2357             : /* a = ppi(a,b) ppo(a,b), where ppi regroups primes common to a and b
    2358             :  * and ppo(a,b) = Z_ppo(a,b) */
    2359             : /* return gcd(a,b),ppi(a,b),ppo(a,b) */
    2360             : GEN
    2361      452865 : Z_ppio(GEN a, GEN b)
    2362             : {
    2363      452865 :   GEN x, y, d = gcdii(a,b);
    2364      452865 :   if (is_pm1(d)) return mkvec3(gen_1, gen_1, a);
    2365      344568 :   x = d; y = diviiexact(a,d);
    2366             :   for(;;)
    2367             :   {
    2368      407253 :     GEN g = gcdii(x,y);
    2369      407253 :     if (is_pm1(g)) return mkvec3(d, x, y);
    2370       62685 :     x = mulii(x,g); y = diviiexact(y,g);
    2371       62685 :   }
    2372             : }
    2373             : /* a = ppg(a,b)pple(a,b), where ppg regroups primes such that v(a) > v(b)
    2374             :  * and pple all others */
    2375             : /* return gcd(a,b),ppg(a,b),pple(a,b) */
    2376             : GEN
    2377           0 : Z_ppgle(GEN a, GEN b)
    2378             : {
    2379           0 :   GEN x, y, g, d = gcdii(a,b);
    2380           0 :   if (equalii(a, d)) return mkvec3(a, gen_1, a);
    2381           0 :   x = diviiexact(a,d); y = d;
    2382             :   for(;;)
    2383             :   {
    2384           0 :     g = gcdii(x,y);
    2385           0 :     if (is_pm1(g)) return mkvec3(d, x, y);
    2386           0 :     x = mulii(x,g); y = diviiexact(y,g);
    2387           0 :   }
    2388             : }
    2389             : static void
    2390           0 : Z_dcba_rec(GEN L, GEN a, GEN b)
    2391             : {
    2392             :   GEN x, r, v, g, h, c, c0;
    2393             :   long n;
    2394           0 :   if (is_pm1(b)) {
    2395           0 :     if (!is_pm1(a)) vectrunc_append(L, a);
    2396           0 :     return;
    2397             :   }
    2398           0 :   v = Z_ppio(a,b);
    2399           0 :   a = gel(v,2);
    2400           0 :   r = gel(v,3);
    2401           0 :   if (!is_pm1(r)) vectrunc_append(L, r);
    2402           0 :   v = Z_ppgle(a,b);
    2403           0 :   g = gel(v,1);
    2404           0 :   h = gel(v,2);
    2405           0 :   x = c0 = gel(v,3);
    2406           0 :   for (n = 1; !is_pm1(h); n++)
    2407             :   {
    2408             :     GEN d, y;
    2409             :     long i;
    2410           0 :     v = Z_ppgle(h,sqri(g));
    2411           0 :     g = gel(v,1);
    2412           0 :     h = gel(v,2);
    2413           0 :     c = gel(v,3); if (is_pm1(c)) continue;
    2414           0 :     d = gcdii(c,b);
    2415           0 :     x = mulii(x,d);
    2416           0 :     y = d; for (i=1; i < n; i++) y = sqri(y);
    2417           0 :     Z_dcba_rec(L, diviiexact(c,y), d);
    2418             :   }
    2419           0 :   Z_dcba_rec(L,diviiexact(b,x), c0);
    2420             : }
    2421             : static GEN
    2422     3068653 : Z_cba_rec(GEN L, GEN a, GEN b)
    2423             : {
    2424             :   GEN g;
    2425     3068653 :   if (lg(L) > 10)
    2426             :   { /* a few naive steps before switching to dcba */
    2427           0 :     Z_dcba_rec(L, a, b);
    2428           0 :     return gel(L, lg(L)-1);
    2429             :   }
    2430     3068653 :   if (is_pm1(a)) return b;
    2431     1823311 :   g = gcdii(a,b);
    2432     1823311 :   if (is_pm1(g)) { vectrunc_append(L, a); return b; }
    2433     1362018 :   a = diviiexact(a,g);
    2434     1362018 :   b = diviiexact(b,g);
    2435     1362018 :   return Z_cba_rec(L, Z_cba_rec(L, a, g), b);
    2436             : }
    2437             : GEN
    2438      344617 : Z_cba(GEN a, GEN b)
    2439             : {
    2440      344617 :   GEN L = vectrunc_init(expi(a) + expi(b) + 2);
    2441      344617 :   GEN t = Z_cba_rec(L, a, b);
    2442      344617 :   if (!is_pm1(t)) vectrunc_append(L, t);
    2443      344617 :   return L;
    2444             : }
    2445             : /* P = coprime base, extend it by b; TODO: quadratic for now */
    2446             : GEN
    2447           0 : ZV_cba_extend(GEN P, GEN b)
    2448             : {
    2449           0 :   long i, l = lg(P);
    2450           0 :   GEN w = cgetg(l+1, t_VEC);
    2451           0 :   for (i = 1; i < l; i++)
    2452             :   {
    2453           0 :     GEN v = Z_cba(gel(P,i), b);
    2454           0 :     long nv = lg(v)-1;
    2455           0 :     gel(w,i) = vecslice(v, 1, nv-1); /* those divide P[i] but not b */
    2456           0 :     b = gel(v,nv);
    2457             :   }
    2458           0 :   gel(w,l) = b; return shallowconcat1(w);
    2459             : }
    2460             : GEN
    2461           0 : ZV_cba(GEN v)
    2462             : {
    2463           0 :   long i, l = lg(v);
    2464             :   GEN P;
    2465           0 :   if (l <= 2) return v;
    2466           0 :   P = Z_cba(gel(v,1), gel(v,2));
    2467           0 :   for (i = 3; i < l; i++) P = ZV_cba_extend(P, gel(v,i));
    2468           0 :   return P;
    2469             : }
    2470             : 
    2471             : /* write x = x1 x2, x2 maximal s.t. (x2,f) = 1, return x2 */
    2472             : GEN
    2473     1541075 : Z_ppo(GEN x, GEN f)
    2474             : {
    2475             :   for (;;)
    2476             :   {
    2477     1541075 :     f = gcdii(x, f); if (is_pm1(f)) break;
    2478      953442 :     x = diviiexact(x, f);
    2479      953442 :   }
    2480      587633 :   return x;
    2481             : }
    2482             : /* write x = x1 x2, x2 maximal s.t. (x2,f) = 1, return x2 */
    2483             : ulong
    2484    40321635 : u_ppo(ulong x, ulong f)
    2485             : {
    2486             :   for (;;)
    2487             :   {
    2488    40321635 :     f = ugcd(x, f); if (f == 1) break;
    2489     8105269 :     x /= f;
    2490     8105269 :   }
    2491    32216366 :   return x;
    2492             : }
    2493             : 
    2494             : /* x t_INT, f ideal. Write x = x1 x2, sqf(x1) | f, (x2,f) = 1. Return x2 */
    2495             : static GEN
    2496         140 : nf_coprime_part(GEN nf, GEN x, GEN listpr)
    2497             : {
    2498         140 :   long v, j, lp = lg(listpr), N = nf_get_degree(nf);
    2499             :   GEN x1, x2, ex;
    2500             : 
    2501             : #if 0 /*1) via many gcds. Expensive ! */
    2502             :   GEN f = idealprodprime(nf, listpr);
    2503             :   f = ZM_hnfmodid(f, x); /* first gcd is less expensive since x in Z */
    2504             :   x = scalarmat(x, N);
    2505             :   for (;;)
    2506             :   {
    2507             :     if (gequal1(gcoeff(f,1,1))) break;
    2508             :     x = idealdivexact(nf, x, f);
    2509             :     f = ZM_hnfmodid(shallowconcat(f,x), gcoeff(x,1,1)); /* gcd(f,x) */
    2510             :   }
    2511             :   x2 = x;
    2512             : #else /*2) from prime decomposition */
    2513         140 :   x1 = NULL;
    2514         399 :   for (j=1; j<lp; j++)
    2515             :   {
    2516         259 :     GEN pr = gel(listpr,j);
    2517         259 :     v = Z_pval(x, pr_get_p(pr)); if (!v) continue;
    2518             : 
    2519         140 :     ex = muluu(v, pr_get_e(pr)); /* = v_pr(x) > 0 */
    2520         140 :     x1 = x1? idealmulpowprime(nf, x1, pr, ex)
    2521         140 :            : idealpow(nf, pr, ex);
    2522             :   }
    2523         140 :   x = scalarmat(x, N);
    2524         140 :   x2 = x1? idealdivexact(nf, x, x1): x;
    2525             : #endif
    2526         140 :   return x2;
    2527             : }
    2528             : 
    2529             : /* L0 in K^*, assume (L0,f) = 1. Return L integral, L0 = L mod f  */
    2530             : GEN
    2531        5600 : make_integral(GEN nf, GEN L0, GEN f, GEN listpr)
    2532             : {
    2533             :   GEN fZ, t, L, D2, d1, d2, d;
    2534             : 
    2535        5600 :   L = Q_remove_denom(L0, &d);
    2536        5600 :   if (!d) return L0;
    2537             : 
    2538             :   /* L0 = L / d, L integral */
    2539        2149 :   fZ = gcoeff(f,1,1);
    2540        2149 :   if (typ(L) == t_INT) return Fp_mul(L, Fp_inv(d, fZ), fZ);
    2541             :   /* Kill denom part coprime to fZ */
    2542        1918 :   d2 = Z_ppo(d, fZ);
    2543        1918 :   t = Fp_inv(d2, fZ); if (!is_pm1(t)) L = ZC_Z_mul(L,t);
    2544        1918 :   if (equalii(d, d2)) return L;
    2545             : 
    2546         140 :   d1 = diviiexact(d, d2);
    2547             :   /* L0 = (L / d1) mod f. d1 not coprime to f
    2548             :    * write (d1) = D1 D2, D2 minimal, (D2,f) = 1. */
    2549         140 :   D2 = nf_coprime_part(nf, d1, listpr);
    2550         140 :   t = idealaddtoone_i(nf, D2, f); /* in D2, 1 mod f */
    2551         140 :   L = nfmuli(nf,t,L);
    2552             : 
    2553             :   /* if (L0, f) = 1, then L in D1 ==> in D1 D2 = (d1) */
    2554         140 :   return Q_div_to_int(L, d1); /* exact division */
    2555             : }
    2556             : 
    2557             : /* assume L is a list of prime ideals. Return the product */
    2558             : GEN
    2559         126 : idealprodprime(GEN nf, GEN L)
    2560             : {
    2561         126 :   long l = lg(L), i;
    2562             :   GEN z;
    2563         126 :   if (l == 1) return matid(nf_get_degree(nf));
    2564         126 :   z = pr_hnf(nf, gel(L,1));
    2565         126 :   for (i=2; i<l; i++) z = idealHNF_mul_two(nf,z, gel(L,i));
    2566         126 :   return z;
    2567             : }
    2568             : 
    2569             : /* optimize for the frequent case I = nfhnf()[2]: lots of them are 1 */
    2570             : GEN
    2571        1484 : idealprod(GEN nf, GEN I)
    2572             : {
    2573        1484 :   long i, l = lg(I);
    2574             :   GEN z;
    2575        2576 :   for (i = 1; i < l; i++)
    2576        2555 :     if (!equali1(gel(I,i))) break;
    2577        1484 :   if (i == l) return gen_1;
    2578        1463 :   z = gel(I,i);
    2579        1463 :   for (i++; i<l; i++) z = idealmul(nf, z, gel(I,i));
    2580        1463 :   return z;
    2581             : }
    2582             : 
    2583             : /* assume L is a list of prime ideals. Return prod L[i]^e[i] */
    2584             : GEN
    2585        7266 : factorbackprime(GEN nf, GEN L, GEN e)
    2586             : {
    2587        7266 :   long l = lg(L), i;
    2588             :   GEN z;
    2589             : 
    2590        7266 :   if (l == 1) return matid(nf_get_degree(nf));
    2591        7252 :   z = idealpow(nf, gel(L,1), gel(e,1));
    2592       11130 :   for (i=2; i<l; i++)
    2593        3878 :     if (signe(gel(e,i))) z = idealmulpowprime(nf,z, gel(L,i),gel(e,i));
    2594        7252 :   return z;
    2595             : }
    2596             : 
    2597             : /* F in Z, divisible exactly by pr.p. Return F-uniformizer for pr, i.e.
    2598             :  * a t in Z_K such that v_pr(t) = 1 and (t, F/pr) = 1 */
    2599             : GEN
    2600       18102 : pr_uniformizer(GEN pr, GEN F)
    2601             : {
    2602       18102 :   GEN p = pr_get_p(pr), t = pr_get_gen(pr);
    2603       18102 :   if (!equalii(F, p))
    2604             :   {
    2605        7630 :     long e = pr_get_e(pr);
    2606        7630 :     GEN u, v, q = (e == 1)? sqri(p): p;
    2607        7630 :     u = mulii(q, Fp_inv(q, diviiexact(F,p))); /* 1 mod F/p, 0 mod q */
    2608        7630 :     v = subui(1UL, u); /* 0 mod F/p, 1 mod q */
    2609        7630 :     if (pr_is_inert(pr))
    2610           0 :       t = addii(mulii(p, v), u);
    2611             :     else
    2612             :     {
    2613        7630 :       t = ZC_Z_mul(t, v);
    2614        7630 :       gel(t,1) = addii(gel(t,1), u); /* return u + vt */
    2615             :     }
    2616             :   }
    2617       18102 :   return t;
    2618             : }
    2619             : /* L = list of prime ideals, return lcm_i (L[i] \cap \ZM) */
    2620             : GEN
    2621       35320 : prV_lcm_capZ(GEN L)
    2622             : {
    2623       35320 :   long i, r = lg(L);
    2624             :   GEN F;
    2625       35320 :   if (r == 1) return gen_1;
    2626       29888 :   F = pr_get_p(gel(L,1));
    2627       44471 :   for (i = 2; i < r; i++)
    2628             :   {
    2629       14583 :     GEN pr = gel(L,i), p = pr_get_p(pr);
    2630       14583 :     if (!dvdii(F, p)) F = mulii(F,p);
    2631             :   }
    2632       29888 :   return F;
    2633             : }
    2634             : 
    2635             : /* Given a prime ideal factorization with possibly zero or negative
    2636             :  * exponents, gives b such that v_p(b) = v_p(x) for all prime ideals pr | x
    2637             :  * and v_pr(b) >= 0 for all other pr.
    2638             :  * For optimal performance, all [anti-]uniformizers should be precomputed,
    2639             :  * but no support for this yet.
    2640             :  *
    2641             :  * If nored, do not reduce result.
    2642             :  * No garbage collecting */
    2643             : static GEN
    2644       20893 : idealapprfact_i(GEN nf, GEN x, int nored)
    2645             : {
    2646             :   GEN z, d, L, e, e2, F;
    2647             :   long i, r;
    2648             :   int flagden;
    2649             : 
    2650       20893 :   nf = checknf(nf);
    2651       20893 :   L = gel(x,1);
    2652       20893 :   e = gel(x,2);
    2653       20893 :   F = prV_lcm_capZ(L);
    2654       20893 :   flagden = 0;
    2655       20893 :   z = NULL; r = lg(e);
    2656       44315 :   for (i = 1; i < r; i++)
    2657             :   {
    2658       23422 :     long s = signe(gel(e,i));
    2659             :     GEN pi, q;
    2660       23422 :     if (!s) continue;
    2661       15946 :     if (s < 0) flagden = 1;
    2662       15946 :     pi = pr_uniformizer(gel(L,i), F);
    2663       15946 :     q = nfpow(nf, pi, gel(e,i));
    2664       15946 :     z = z? nfmul(nf, z, q): q;
    2665             :   }
    2666       20893 :   if (!z) return gen_1;
    2667       11016 :   if (nored || typ(z) != t_COL) return z;
    2668        2716 :   e2 = cgetg(r, t_VEC);
    2669        2716 :   for (i=1; i<r; i++) gel(e2,i) = addiu(gel(e,i), 1);
    2670        2716 :   x = factorbackprime(nf, L,e2);
    2671        2716 :   if (flagden) /* denominator */
    2672             :   {
    2673        2702 :     z = Q_remove_denom(z, &d);
    2674        2702 :     d = diviiexact(d, Z_ppo(d, F));
    2675        2702 :     x = RgM_Rg_mul(x, d);
    2676             :   }
    2677             :   else
    2678          14 :     d = NULL;
    2679        2716 :   z = ZC_reducemodlll(z, x);
    2680        2716 :   return d? RgC_Rg_div(z,d): z;
    2681             : }
    2682             : 
    2683             : GEN
    2684           0 : idealapprfact(GEN nf, GEN x) {
    2685           0 :   pari_sp av = avma;
    2686           0 :   return gerepileupto(av, idealapprfact_i(nf, x, 0));
    2687             : }
    2688             : GEN
    2689          14 : idealappr(GEN nf, GEN x) {
    2690          14 :   pari_sp av = avma;
    2691          14 :   if (!is_nf_extfactor(x)) x = idealfactor(nf, x);
    2692          14 :   return gerepileupto(av, idealapprfact_i(nf, x, 0));
    2693             : }
    2694             : 
    2695             : /* OBSOLETE */
    2696             : GEN
    2697          14 : idealappr0(GEN nf, GEN x, long fl) { (void)fl; return idealappr(nf, x); }
    2698             : 
    2699             : static GEN
    2700          21 : mat_ideal_two_elt2(GEN nf, GEN x, GEN a)
    2701             : {
    2702          21 :   GEN F = idealfactor(nf,a), P = gel(F,1), E = gel(F,2);
    2703          21 :   long i, r = lg(E);
    2704          21 :   for (i=1; i<r; i++) gel(E,i) = stoi( idealval(nf,x,gel(P,i)) );
    2705          21 :   return idealapprfact_i(nf,F,1);
    2706             : }
    2707             : 
    2708             : static void
    2709          14 : not_in_ideal(GEN a) {
    2710          14 :   pari_err_DOMAIN("idealtwoelt2","element mod ideal", "!=", gen_0, a);
    2711           0 : }
    2712             : /* x integral in HNF, a an 'nf' */
    2713             : static int
    2714          28 : in_ideal(GEN x, GEN a)
    2715             : {
    2716          28 :   switch(typ(a))
    2717             :   {
    2718          14 :     case t_INT: return dvdii(a, gcoeff(x,1,1));
    2719           7 :     case t_COL: return RgV_is_ZV(a) && !!hnf_invimage(x, a);
    2720           7 :     default: return 0;
    2721             :   }
    2722             : }
    2723             : 
    2724             : /* Given an integral ideal x and a in x, gives a b such that
    2725             :  * x = aZ_K + bZ_K using the approximation theorem */
    2726             : GEN
    2727          42 : idealtwoelt2(GEN nf, GEN x, GEN a)
    2728             : {
    2729          42 :   pari_sp av = avma;
    2730             :   GEN cx, b;
    2731             : 
    2732          42 :   nf = checknf(nf);
    2733          42 :   a = nf_to_scalar_or_basis(nf, a);
    2734          42 :   x = idealhnf_shallow(nf,x);
    2735          42 :   if (lg(x) == 1)
    2736             :   {
    2737          14 :     if (!isintzero(a)) not_in_ideal(a);
    2738           7 :     avma = av; return gen_0;
    2739             :   }
    2740          28 :   x = Q_primitive_part(x, &cx);
    2741          28 :   if (cx) a = gdiv(a, cx);
    2742          28 :   if (!in_ideal(x, a)) not_in_ideal(a);
    2743          21 :   b = mat_ideal_two_elt2(nf, x, a);
    2744          21 :   if (typ(b) == t_COL)
    2745             :   {
    2746          14 :     GEN mod = idealhnf_principal(nf,a);
    2747          14 :     b = ZC_hnfrem(b,mod);
    2748          14 :     if (ZV_isscalar(b)) b = gel(b,1);
    2749             :   }
    2750             :   else
    2751             :   {
    2752           7 :     GEN aZ = typ(a) == t_COL? Q_denom(zk_inv(nf,a)): a; /* (a) \cap Z */
    2753           7 :     b = centermodii(b, aZ, shifti(aZ,-1));
    2754             :   }
    2755          21 :   b = cx? gmul(b,cx): gcopy(b);
    2756          21 :   return gerepileupto(av, b);
    2757             : }
    2758             : 
    2759             : /* Given 2 integral ideals x and y in nf, returns a beta in nf such that
    2760             :  * beta * x is an integral ideal coprime to y */
    2761             : GEN
    2762       12572 : idealcoprimefact(GEN nf, GEN x, GEN fy)
    2763             : {
    2764       12572 :   GEN L = gel(fy,1), e;
    2765       12572 :   long i, r = lg(L);
    2766             : 
    2767       12572 :   e = cgetg(r, t_COL);
    2768       12572 :   for (i=1; i<r; i++) gel(e,i) = stoi( -idealval(nf,x,gel(L,i)) );
    2769       12572 :   return idealapprfact_i(nf, mkmat2(L,e), 0);
    2770             : }
    2771             : GEN
    2772          70 : idealcoprime(GEN nf, GEN x, GEN y)
    2773             : {
    2774          70 :   pari_sp av = avma;
    2775          70 :   return gerepileupto(av, idealcoprimefact(nf, x, idealfactor(nf,y)));
    2776             : }
    2777             : 
    2778             : GEN
    2779           7 : nfmulmodpr(GEN nf, GEN x, GEN y, GEN modpr)
    2780             : {
    2781           7 :   pari_sp av = avma;
    2782           7 :   GEN z, p, pr = modpr, T;
    2783             : 
    2784           7 :   nf = checknf(nf); modpr = nf_to_Fq_init(nf,&pr,&T,&p);
    2785           0 :   x = nf_to_Fq(nf,x,modpr);
    2786           0 :   y = nf_to_Fq(nf,y,modpr);
    2787           0 :   z = Fq_mul(x,y,T,p);
    2788           0 :   return gerepileupto(av, algtobasis(nf, Fq_to_nf(z,modpr)));
    2789             : }
    2790             : 
    2791             : GEN
    2792           0 : nfdivmodpr(GEN nf, GEN x, GEN y, GEN modpr)
    2793             : {
    2794           0 :   pari_sp av = avma;
    2795           0 :   nf = checknf(nf);
    2796           0 :   return gerepileupto(av, nfreducemodpr(nf, nfdiv(nf,x,y), modpr));
    2797             : }
    2798             : 
    2799             : GEN
    2800           0 : nfpowmodpr(GEN nf, GEN x, GEN k, GEN modpr)
    2801             : {
    2802           0 :   pari_sp av=avma;
    2803           0 :   GEN z, T, p, pr = modpr;
    2804             : 
    2805           0 :   nf = checknf(nf); modpr = nf_to_Fq_init(nf,&pr,&T,&p);
    2806           0 :   z = nf_to_Fq(nf,x,modpr);
    2807           0 :   z = Fq_pow(z,k,T,p);
    2808           0 :   return gerepileupto(av, algtobasis(nf, Fq_to_nf(z,modpr)));
    2809             : }
    2810             : 
    2811             : GEN
    2812           0 : nfkermodpr(GEN nf, GEN x, GEN modpr)
    2813             : {
    2814           0 :   pari_sp av = avma;
    2815           0 :   GEN T, p, pr = modpr;
    2816             : 
    2817           0 :   nf = checknf(nf); modpr = nf_to_Fq_init(nf, &pr,&T,&p);
    2818           0 :   if (typ(x)!=t_MAT) pari_err_TYPE("nfkermodpr",x);
    2819           0 :   x = nfM_to_FqM(x, nf, modpr);
    2820           0 :   return gerepilecopy(av, FqM_to_nfM(FqM_ker(x,T,p), modpr));
    2821             : }
    2822             : 
    2823             : GEN
    2824           0 : nfsolvemodpr(GEN nf, GEN a, GEN b, GEN pr)
    2825             : {
    2826           0 :   const char *f = "nfsolvemodpr";
    2827           0 :   pari_sp av = avma;
    2828             :   GEN T, p, modpr;
    2829             : 
    2830           0 :   nf = checknf(nf);
    2831           0 :   modpr = nf_to_Fq_init(nf, &pr,&T,&p);
    2832           0 :   if (typ(a)!=t_MAT) pari_err_TYPE(f,a);
    2833           0 :   a = nfM_to_FqM(a, nf, modpr);
    2834           0 :   switch(typ(b))
    2835             :   {
    2836             :     case t_MAT:
    2837           0 :       b = nfM_to_FqM(b, nf, modpr);
    2838           0 :       b = FqM_gauss(a,b,T,p);
    2839           0 :       if (!b) pari_err_INV(f,a);
    2840           0 :       a = FqM_to_nfM(b, modpr);
    2841           0 :       break;
    2842             :     case t_COL:
    2843           0 :       b = nfV_to_FqV(b, nf, modpr);
    2844           0 :       b = FqM_FqC_gauss(a,b,T,p);
    2845           0 :       if (!b) pari_err_INV(f,a);
    2846           0 :       a = FqV_to_nfV(b, modpr);
    2847           0 :       break;
    2848           0 :     default: pari_err_TYPE(f,b);
    2849             :   }
    2850           0 :   return gerepilecopy(av, a);
    2851             : }

Generated by: LCOV version 1.11