Line data Source code
1 : /* Copyright (C) 2000 The PARI group.
2 :
3 : This file is part of the PARI/GP package.
4 :
5 : PARI/GP is free software; you can redistribute it and/or modify it under the
6 : terms of the GNU General Public License as published by the Free Software
7 : Foundation; either version 2 of the License, or (at your option) any later
8 : version. It is distributed in the hope that it will be useful, but WITHOUT
9 : ANY WARRANTY WHATSOEVER.
10 :
11 : Check the License for details. You should have received a copy of it, along
12 : with the package; see the file 'COPYING'. If not, write to the Free Software
13 : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
14 :
15 : /*******************************************************************/
16 : /* */
17 : /* BASIC NF OPERATIONS */
18 : /* (continued) */
19 : /* */
20 : /*******************************************************************/
21 : #include "pari.h"
22 : #include "paripriv.h"
23 :
24 : #define DEBUGLEVEL DEBUGLEVEL_nf
25 :
26 : /*******************************************************************/
27 : /* */
28 : /* IDEAL OPERATIONS */
29 : /* */
30 : /*******************************************************************/
31 :
32 : /* A valid ideal is either principal (valid nf_element), or prime, or a matrix
33 : * on the integer basis in HNF.
34 : * A prime ideal is of the form [p,a,e,f,b], where the ideal is p.Z_K+a.Z_K,
35 : * p is a rational prime, a belongs to Z_K, e=e(P/p), f=f(P/p), and b
36 : * is Lenstra's constant, such that p.P^(-1)= p Z_K + b Z_K.
37 : *
38 : * An extended ideal is a couple [I,F] where I is an ideal and F is either an
39 : * algebraic number, or a factorization matrix attached to an algebraic number.
40 : * All routines work with either extended ideals or ideals (an omitted F is
41 : * assumed to be factor(1)). All ideals are output in HNF form. */
42 :
43 : /* types and conversions */
44 :
45 : long
46 16234071 : idealtyp(GEN *ideal, GEN *arch)
47 : {
48 16234071 : GEN x = *ideal;
49 16234071 : long t,lx,tx = typ(x);
50 :
51 16234071 : if (tx!=t_VEC || lg(x)!=3) { if (arch) *arch = NULL; }
52 : else
53 : {
54 1239710 : GEN a = gel(x,2);
55 1239710 : if (typ(a) == t_MAT && lg(a) != 3)
56 : { /* allow [;] */
57 14 : if (lg(a) != 1) pari_err_TYPE("idealtyp [extended ideal]",x);
58 7 : if (arch) *arch = trivial_fact();
59 : }
60 : else
61 1239696 : if (arch) *arch = a;
62 1239703 : x = gel(x,1); tx = typ(x);
63 : }
64 16234064 : switch(tx)
65 : {
66 12229629 : case t_MAT: lx = lg(x);
67 12229629 : if (lx == 1) { t = id_PRINCIPAL; x = gen_0; break; }
68 12229454 : if (lx != lgcols(x)) pari_err_TYPE("idealtyp [nonsquare t_MAT]",x);
69 12229411 : t = id_MAT;
70 12229411 : break;
71 :
72 3047959 : case t_VEC:
73 3047959 : if (!checkprid_i(x)) pari_err_TYPE("idealtyp [fake prime ideal]",x);
74 3047920 : t = id_PRIME; break;
75 :
76 956762 : case t_POL: case t_POLMOD: case t_COL:
77 : case t_INT: case t_FRAC:
78 956762 : t = id_PRINCIPAL; break;
79 6 : default:
80 6 : pari_err_TYPE("idealtyp",x);
81 : return 0; /*LCOV_EXCL_LINE*/
82 : }
83 16234268 : *ideal = x; return t;
84 : }
85 :
86 : /* true nf; v = [a,x,...], a in Z. Return (a,x) */
87 : GEN
88 802258 : idealhnf_two(GEN nf, GEN v)
89 : {
90 802258 : GEN p = gel(v,1), pi = gel(v,2), m = zk_scalar_or_multable(nf, pi);
91 802230 : if (typ(m) == t_INT) return scalarmat(gcdii(m,p), nf_get_degree(nf));
92 728836 : return ZM_hnfmodid(m, p);
93 : }
94 : /* true nf */
95 : GEN
96 3801278 : pr_hnf(GEN nf, GEN pr)
97 : {
98 3801278 : GEN p = pr_get_p(pr), m;
99 3801243 : if (pr_is_inert(pr)) return scalarmat(p, nf_get_degree(nf));
100 3364399 : m = zk_scalar_or_multable(nf, pr_get_gen(pr));
101 3364068 : return ZM_hnfmodprime(m, p);
102 : }
103 :
104 : GEN
105 1391734 : idealhnf_principal(GEN nf, GEN x)
106 : {
107 : GEN cx;
108 1391734 : x = nf_to_scalar_or_basis(nf, x);
109 1391732 : switch(typ(x))
110 : {
111 1130531 : case t_COL: break;
112 242957 : case t_INT: if (!signe(x)) return cgetg(1,t_MAT);
113 241291 : return scalarmat(absi_shallow(x), nf_get_degree(nf));
114 18245 : case t_FRAC:
115 18245 : return scalarmat(Q_abs_shallow(x), nf_get_degree(nf));
116 0 : default: pari_err_TYPE("idealhnf",x);
117 : }
118 1130531 : x = Q_primitive_part(x, &cx);
119 1130524 : RgV_check_ZV(x, "idealhnf");
120 1130524 : x = zk_multable(nf, x);
121 1130524 : x = ZM_hnfmodid(x, zkmultable_capZ(x));
122 1130529 : return cx? ZM_Q_mul(x,cx): x;
123 : }
124 :
125 : /* true nf; x integral Z_K-module as t_MAT generated by its columns.
126 : * Return square hnf representation */
127 : static GEN
128 91 : vec_mulid(GEN nf, GEN x)
129 : {
130 91 : long i, j, l = lg(x);
131 91 : GEN D = NULL, v;
132 91 : v = cgetg(l, t_VEC);
133 196 : for (i = j = 1; i < l; i++)
134 : {
135 : GEN m, d;
136 161 : if (D && ZV_Z_dvd(gel(x,i), D)) l--;
137 161 : gel(v,j++) = m = zk_multable(nf, gel(x,i));
138 161 : d = zkmultable_capZ(m);
139 161 : D = D? gcdii(D, d): d;
140 161 : if (is_pm1(D)) return matid(lg(m)-1);
141 : }
142 35 : setlg(v, l); if (l == 1) return cgetg(1, t_MAT);
143 35 : return ZM_hnfmodid(shallowconcat1(v), D);
144 : }
145 :
146 : static GEN
147 84 : nfV_idealhnf(GEN nf, GEN v, GEN *pden)
148 : {
149 84 : GEN H = ZM_hnf(Q_remove_denom(matalgtobasis(nf, v), pden));
150 84 : return vec_mulid(nf, H);
151 : }
152 :
153 : /* true nf */
154 : GEN
155 1855929 : idealhnf_shallow(GEN nf, GEN x)
156 : {
157 1855929 : long tx = typ(x), lx = lg(x), N;
158 :
159 : /* cannot use idealtyp because here we allow nonsquare matrices */
160 1855929 : if (tx == t_VEC && lx == 3) { x = gel(x,1); tx = typ(x); lx = lg(x); }
161 1855929 : if (tx == t_VEC && lx == 6)
162 : {
163 538829 : if (!checkprid_i(x)) pari_err_TYPE("idealhnf [fake prime ideal]",x);
164 538818 : return pr_hnf(nf,x); /* PRIME */
165 : }
166 1317100 : switch(tx)
167 : {
168 91845 : case t_MAT:
169 : {
170 : GEN cx;
171 91845 : long nx = lx-1;
172 91845 : N = nf_get_degree(nf);
173 91845 : if (nx == 0) return cgetg(1, t_MAT);
174 91824 : if (nbrows(x) != N) pari_err_TYPE("idealhnf [wrong dimension]",x);
175 91817 : if (nx == 1) return idealhnf_principal(nf, gel(x,1));
176 :
177 74934 : if (nx == N && RgM_is_ZM(x) && ZM_ishnf(x)) return x;
178 46382 : x = Q_primitive_part(x, &cx);
179 46382 : if (nx < N)
180 7 : x = vec_mulid(nf, x); /* build ZK-module generated from cols */
181 : else
182 46375 : x = ZM_hnfmod(x, ZM_detmult(x)); /* assume Z-span cols is ZK-module */
183 46382 : return cx? ZM_Q_mul(x,cx): x;
184 : }
185 14 : case t_QFB:
186 : {
187 14 : pari_sp av = avma;
188 14 : GEN u, D = nf_get_disc(nf), T = nf_get_pol(nf), f = nf_get_index(nf);
189 14 : GEN A = gel(x,1), B = gel(x,2);
190 14 : N = nf_get_degree(nf);
191 14 : if (N != 2)
192 0 : pari_err_TYPE("idealhnf [Qfb for nonquadratic fields]", x);
193 14 : if (!equalii(qfb_disc(x), D))
194 7 : pari_err_DOMAIN("idealhnf [Qfb]", "disc(q)", "!=", D, x);
195 : /* x -> A Z + (-B + sqrt(D)) / 2 Z
196 : K = Q[t]/T(t), t^2 + ut + v = 0, u^2 - 4v = Df^2
197 : => t = (-u + sqrt(D) f)/2
198 : => sqrt(D)/2 = (t + u/2)/f */
199 7 : u = gel(T,3);
200 7 : B = deg1pol_shallow(ginv(f),
201 : gsub(gdiv(u, shifti(f,1)), gdiv(B,gen_2)),
202 7 : varn(T));
203 7 : return gc_upto(av, idealhnf_two(nf, mkvec2(A,B)));
204 : }
205 1225241 : default: return idealhnf_principal(nf, x); /* PRINCIPAL */
206 : }
207 : }
208 : /* true nf */
209 : GEN
210 665 : idealhnf(GEN nf, GEN x)
211 : {
212 665 : pari_sp av = avma;
213 665 : GEN y = idealhnf_shallow(nf, x);
214 651 : return (avma == av)? gcopy(y): gc_upto(av, y);
215 : }
216 :
217 : static GEN
218 84 : nfV_eltembed(GEN nf, GEN x, long prec)
219 518 : { pari_APPLY_type(t_VEC, nfeltembed(nf, gel(x,i), NULL, prec)) }
220 :
221 : /* true nf */
222 : static GEN
223 84 : nfweilheight_i(GEN nf, GEN v, long prec)
224 : {
225 84 : long i, j, r1, r2, u, N, l = lg(v);
226 84 : GEN den, h = gen_1, id = nfV_idealhnf(nf, v, &den);
227 84 : GEN V = nfV_eltembed(nf, v, prec);
228 :
229 84 : nf_get_sign(nf, &r1, &r2); u = r1 + r2; N = u + r2;
230 259 : for (i = 1; i <= r1; i++)
231 1029 : for (j = 1; j < l; j++) gmael(V,j,i) = gabs(gmael(V,j,i), prec);
232 343 : for ( ; i <= u; i++)
233 1771 : for (j = 1; j < l; j++) gmael(V,j,i) = gnorm(gmael(V,j,i));
234 518 : for (i = 1; i <= u; i++)
235 : {
236 434 : long j0 = 1;
237 2366 : for (j = 2; j < l; j++)
238 1932 : if (gcmp(gmael(V,j,i), gmael(V,j0,i)) > 0) j0 = j;
239 434 : h = gmul(h, gmael(V,j0,i));
240 : }
241 84 : if (den) h = gmul(h, powiu(den, N));
242 84 : return divru(glog(gdiv(h, idealnorm(nf, id)), prec), N);
243 : }
244 :
245 : GEN
246 84 : nfweilheight(GEN nf, GEN v, long prec)
247 : {
248 84 : pari_sp av = avma;
249 84 : nf = checknf(nf);
250 84 : if (!is_vec_t(typ(v)) || lg(v) < 2) pari_err_TYPE("nfweilheight",v);
251 84 : return gc_upto(av, nfweilheight_i(nf, v, prec));
252 : }
253 :
254 : /* GP functions */
255 :
256 : GEN
257 2485 : idealtwoelt0(GEN nf, GEN x, GEN a)
258 : {
259 2485 : if (!a) return idealtwoelt(nf,x);
260 42 : return idealtwoelt2(nf,x,a);
261 : }
262 :
263 : GEN
264 2499 : idealpow0(GEN nf, GEN x, GEN n, long flag)
265 : {
266 2499 : if (flag) return idealpowred(nf,x,n);
267 2492 : return idealpow(nf,x,n);
268 : }
269 :
270 : GEN
271 70 : idealmul0(GEN nf, GEN x, GEN y, long flag)
272 : {
273 70 : if (flag) return idealmulred(nf,x,y);
274 63 : return idealmul(nf,x,y);
275 : }
276 :
277 : GEN
278 56 : idealdiv0(GEN nf, GEN x, GEN y, long flag)
279 : {
280 56 : switch(flag)
281 : {
282 28 : case 0: return idealdiv(nf,x,y);
283 28 : case 1: return idealdivexact(nf,x,y);
284 0 : default: pari_err_FLAG("idealdiv");
285 : }
286 : return NULL; /* LCOV_EXCL_LINE */
287 : }
288 :
289 : GEN
290 70 : idealaddtoone0(GEN nf, GEN arg1, GEN arg2)
291 : {
292 70 : if (!arg2) return idealaddmultoone(nf,arg1);
293 35 : return idealaddtoone(nf,arg1,arg2);
294 : }
295 :
296 : /* b not a scalar */
297 : static GEN
298 77 : hnf_Z_ZC(GEN nf, GEN a, GEN b) { return hnfmodid(zk_multable(nf,b), a); }
299 : /* b not a scalar */
300 : static GEN
301 70 : hnf_Z_QC(GEN nf, GEN a, GEN b)
302 : {
303 : GEN db;
304 70 : b = Q_remove_denom(b, &db);
305 70 : if (db) a = mulii(a, db);
306 70 : b = hnf_Z_ZC(nf,a,b);
307 70 : return db? RgM_Rg_div(b, db): b;
308 : }
309 : /* b not a scalar (not point in trying to optimize for this case) */
310 : static GEN
311 77 : hnf_Q_QC(GEN nf, GEN a, GEN b)
312 : {
313 : GEN da, db;
314 77 : if (typ(a) == t_INT) return hnf_Z_QC(nf, a, b);
315 7 : da = gel(a,2);
316 7 : a = gel(a,1);
317 7 : b = Q_remove_denom(b, &db);
318 : /* write da = d*A, db = d*B, gcd(A,B) = 1
319 : * gcd(a/(d A), b/(d B)) = gcd(a B, A b) / A B d = gcd(a B, b) / A B d */
320 7 : if (db)
321 : {
322 7 : GEN d = gcdii(da,db);
323 7 : if (!is_pm1(d)) db = diviiexact(db,d); /* B */
324 7 : if (!is_pm1(db))
325 : {
326 7 : a = mulii(a, db); /* a B */
327 7 : da = mulii(da, db); /* A B d = lcm(denom(a),denom(b)) */
328 : }
329 : }
330 7 : return RgM_Rg_div(hnf_Z_ZC(nf,a,b), da);
331 : }
332 : static GEN
333 7 : hnf_QC_QC(GEN nf, GEN a, GEN b)
334 : {
335 : GEN da, db, d, x;
336 7 : a = Q_remove_denom(a, &da);
337 7 : b = Q_remove_denom(b, &db);
338 7 : if (da) b = ZC_Z_mul(b, da);
339 7 : if (db) a = ZC_Z_mul(a, db);
340 7 : d = mul_denom(da, db);
341 7 : a = zk_multable(nf,a); da = zkmultable_capZ(a);
342 7 : b = zk_multable(nf,b); db = zkmultable_capZ(b);
343 7 : x = ZM_hnfmodid(shallowconcat(a,b), gcdii(da,db));
344 7 : return d? RgM_Rg_div(x, d): x;
345 : }
346 : static GEN
347 21 : hnf_Q_Q(GEN nf, GEN a, GEN b) {return scalarmat(Q_gcd(a,b), nf_get_degree(nf));}
348 : GEN
349 413 : idealhnf0(GEN nf, GEN a, GEN b)
350 : {
351 : long ta, tb;
352 : pari_sp av;
353 : GEN x;
354 413 : nf = checknf(nf);
355 413 : if (!b) return idealhnf(nf,a);
356 :
357 : /* HNF of aZ_K+bZ_K */
358 112 : av = avma;
359 112 : a = nf_to_scalar_or_basis(nf,a); ta = typ(a);
360 112 : b = nf_to_scalar_or_basis(nf,b); tb = typ(b);
361 105 : if (ta == t_COL)
362 14 : x = (tb==t_COL)? hnf_QC_QC(nf, a,b): hnf_Q_QC(nf, b,a);
363 : else
364 91 : x = (tb==t_COL)? hnf_Q_QC(nf, a,b): hnf_Q_Q(nf, a,b);
365 105 : return gc_upto(av, x);
366 : }
367 :
368 : /*******************************************************************/
369 : /* */
370 : /* TWO-ELEMENT FORM */
371 : /* */
372 : /*******************************************************************/
373 : static GEN idealapprfact_i(GEN nf, GEN x, int nored);
374 :
375 : static int
376 225027 : ok_elt(GEN x, GEN xZ, GEN y)
377 : {
378 225027 : pari_sp av = avma;
379 225027 : return gc_bool(av, ZM_equal(x, ZM_hnfmodid(y, xZ)));
380 : }
381 :
382 : /* a + s * b, a and b ZM, s integer */
383 : static GEN
384 66763 : addmul_mat(GEN a, GEN s, GEN b)
385 : {
386 66763 : if (!signe(s)) return a;
387 58126 : if (!equali1(s)) b = ZM_Z_mul(b, s);
388 58126 : return a? ZM_add(a, b): b;
389 : }
390 :
391 : static GEN
392 117658 : get_random_a(GEN nf, GEN x, GEN xZ)
393 : {
394 : pari_sp av;
395 117658 : long i, lm, l = lg(x);
396 : GEN z, beta, mul;
397 :
398 117658 : beta= cgetg(l, t_MAT);
399 117658 : mul = cgetg(l, t_VEC); lm = 1; /* = lg(mul) */
400 : /* look for a in x such that a O/xZ = x O/xZ */
401 249825 : for (i = 2; i < l; i++)
402 : {
403 239729 : GEN xi = gel(x,i);
404 239729 : GEN t = FpM_red(zk_multable(nf,xi), xZ); /* ZM, cannot be a scalar */
405 239724 : if (gequal0(t)) continue;
406 196710 : if (ok_elt(x,xZ, t)) return xi;
407 89151 : gel(beta,lm) = xi;
408 : /* mul[i] = { canonical generators for x[i] O/xZ as Z-module } */
409 89151 : gel(mul,lm) = t; lm++;
410 : }
411 10096 : setlg(mul, lm);
412 10096 : setlg(beta,lm); z = cgetg(lm, t_VEC);
413 30175 : for(av = avma;; set_avma(av))
414 20079 : {
415 30175 : GEN a = NULL;
416 96938 : for (i = 1; i < lm; i++)
417 : {
418 66763 : gel(z,i) = randomi(xZ);
419 66763 : a = addmul_mat(a, gel(z,i), gel(mul,i));
420 : }
421 : /* a = matrix (NOT HNF) of ideal generated by beta.z in O/xZ */
422 30175 : if (a && ok_elt(x,xZ, a)) break;
423 : }
424 10096 : return ZM_ZC_mul(beta, z);
425 : }
426 :
427 : /* x square matrix, assume it is HNF */
428 : static GEN
429 250508 : mat_ideal_two_elt(GEN nf, GEN x)
430 : {
431 : GEN y, a, cx, xZ;
432 250508 : long N = nf_get_degree(nf);
433 : pari_sp av, tetpil;
434 :
435 250508 : if (lg(x)-1 != N) pari_err_DIM("idealtwoelt");
436 250494 : if (N == 2) return mkvec2copy(gcoeff(x,1,1), gel(x,2));
437 :
438 136232 : y = cgetg(3,t_VEC); av = avma;
439 136232 : cx = Q_content(x);
440 136233 : xZ = gcoeff(x,1,1);
441 136233 : if (gequal(xZ, cx)) /* x = (cx) */
442 : {
443 7566 : gel(y,1) = cx;
444 7566 : gel(y,2) = gen_0; return y;
445 : }
446 128667 : if (equali1(cx)) cx = NULL;
447 : else
448 : {
449 1026 : x = Q_div_to_int(x, cx);
450 1026 : xZ = gcoeff(x,1,1);
451 : }
452 128667 : if (N < 6)
453 109186 : a = get_random_a(nf, x, xZ);
454 : else
455 : {
456 19481 : const long FB[] = { _evallg(15+1) | evaltyp(t_VECSMALL),
457 : 2, 3, 5, 7, 11, 13, 17, 19, 23, 29, 31, 37, 41, 43, 47
458 : };
459 19481 : GEN P, E, a1 = Z_lsmoothen(xZ, (GEN)FB, &P, &E);
460 19481 : if (!a1) /* factors completely */
461 11009 : a = idealapprfact_i(nf, idealfactor(nf,x), 1);
462 8472 : else if (lg(P) == 1) /* no small factors */
463 2611 : a = get_random_a(nf, x, xZ);
464 : else /* general case */
465 : {
466 : GEN A0, A1, a0, u0, u1, v0, v1, pi0, pi1, t, u;
467 5861 : a0 = diviiexact(xZ, a1);
468 5861 : A0 = ZM_hnfmodid(x, a0); /* smooth part of x */
469 5861 : A1 = ZM_hnfmodid(x, a1); /* cofactor */
470 5861 : pi0 = idealapprfact_i(nf, idealfactor(nf,A0), 1);
471 5861 : pi1 = get_random_a(nf, A1, a1);
472 5861 : (void)bezout(a0, a1, &v0,&v1);
473 5861 : u0 = mulii(a0, v0);
474 5861 : u1 = mulii(a1, v1);
475 5861 : if (typ(pi0) != t_COL) t = addmulii(u0, pi0, u1);
476 : else
477 5861 : { t = ZC_Z_mul(pi0, u1); gel(t,1) = addii(gel(t,1), u0); }
478 5861 : u = ZC_Z_mul(pi1, u0); gel(u,1) = addii(gel(u,1), u1);
479 5861 : a = nfmuli(nf, centermod(u, xZ), centermod(t, xZ));
480 : }
481 : }
482 128666 : if (cx)
483 : {
484 1026 : a = centermod(a, xZ);
485 1026 : tetpil = avma;
486 1026 : if (typ(cx) == t_INT)
487 : {
488 91 : gel(y,1) = mulii(xZ, cx);
489 91 : gel(y,2) = ZC_Z_mul(a, cx);
490 : }
491 : else
492 : {
493 935 : gel(y,1) = gmul(xZ, cx);
494 935 : gel(y,2) = RgC_Rg_mul(a, cx);
495 : }
496 : }
497 : else
498 : {
499 127640 : tetpil = avma;
500 127640 : gel(y,1) = icopy(xZ);
501 127640 : gel(y,2) = centermod(a, xZ);
502 : }
503 128665 : gc_slice_unsafe(av,tetpil,y+1,2); return y;
504 : }
505 :
506 : /* Given an ideal x, returns [a,alpha] such that a is in Q,
507 : * x = a Z_K + alpha Z_K, alpha in K^*
508 : * a = 0 or alpha = 0 are possible, but do not try to determine whether
509 : * x is principal. */
510 : GEN
511 97578 : idealtwoelt(GEN nf, GEN x)
512 : {
513 : pari_sp av;
514 97578 : long tx = idealtyp(&x, NULL);
515 97571 : nf = checknf(nf);
516 97571 : if (tx == id_MAT) return mat_ideal_two_elt(nf,x);
517 735 : if (tx == id_PRIME) return mkvec2copy(gel(x,1), gel(x,2));
518 : /* id_PRINCIPAL */
519 714 : av = avma; x = nf_to_scalar_or_basis(nf, x);
520 1232 : return gc_GEN(av, typ(x)==t_COL? mkvec2(gen_0,x):
521 609 : mkvec2(Q_abs_shallow(x),gen_0));
522 : }
523 :
524 : /*******************************************************************/
525 : /* */
526 : /* FACTORIZATION */
527 : /* */
528 : /*******************************************************************/
529 : /* x integral ideal in HNF, Zval = v_p(x \cap Z) > 0; return v_p(Nx) */
530 : static long
531 4101110 : idealHNF_norm_pval(GEN x, GEN p, long Zval)
532 : {
533 4101110 : long i, v = Zval, l = lg(x);
534 31884483 : for (i = 2; i < l; i++) v += Z_pval(gcoeff(x,i,i), p);
535 4101132 : return v;
536 : }
537 :
538 : /* x integral in HNF, f0 = partial factorization of a multiple of
539 : * x[1,1] = x\cap Z */
540 : GEN
541 293999 : idealHNF_Z_factor_i(GEN x, GEN f0, GEN *pvN, GEN *pvZ)
542 : {
543 293999 : GEN P, E, vN, vZ, xZ = gcoeff(x,1,1), f = f0? f0: Z_factor(xZ);
544 : long i, l;
545 294007 : P = gel(f,1); l = lg(P);
546 294007 : E = gel(f,2);
547 294007 : *pvN = vN = cgetg(l, t_VECSMALL);
548 294018 : *pvZ = vZ = cgetg(l, t_VECSMALL);
549 734576 : for (i = 1; i < l; i++)
550 : {
551 440556 : GEN p = gel(P,i);
552 440556 : vZ[i] = f0? Z_pval(xZ, p): (long) itou(gel(E,i));
553 440558 : vN[i] = idealHNF_norm_pval(x,p, vZ[i]);
554 : }
555 294020 : return P;
556 : }
557 : /* return P, primes dividing Nx and xZ = x\cap Z, set v_p(Nx), v_p(xZ);
558 : * x integral in HNF */
559 : GEN
560 0 : idealHNF_Z_factor(GEN x, GEN *pvN, GEN *pvZ)
561 0 : { return idealHNF_Z_factor_i(x, NULL, pvN, pvZ); }
562 :
563 : /* v_P(A)*f(P) <= Nval [e.g. Nval = v_p(Norm A)], Zval = v_p(A \cap Z).
564 : * Return v_P(A) */
565 : static long
566 4238704 : idealHNF_val(GEN A, GEN P, long Nval, long Zval)
567 : {
568 4238704 : long f = pr_get_f(P), vmax, v, e, i, j, k, l;
569 : GEN mul, B, a, y, r, p, pk, cx, vals;
570 : pari_sp av;
571 :
572 4238703 : if (Nval < f) return 0;
573 4235533 : p = pr_get_p(P);
574 4235538 : e = pr_get_e(P);
575 : /* v_P(A) <= max [ e * v_p(A \cap Z), floor[v_p(Nix) / f ] */
576 4235539 : vmax = minss(Zval * e, Nval / f);
577 4235534 : mul = pr_get_tau(P);
578 4235532 : l = lg(mul);
579 4235532 : B = cgetg(l,t_MAT);
580 : /* B[1] not needed: v_pr(A[1]) = v_pr(A \cap Z) is known already */
581 4241686 : gel(B,1) = gen_0; /* dummy */
582 23393526 : for (j = 2; j < l; j++)
583 : {
584 21002506 : GEN x = gel(A,j);
585 21002506 : gel(B,j) = y = cgetg(l, t_COL);
586 228181126 : for (i = 1; i < l; i++)
587 : { /* compute a = (x.t0)_i, A in HNF ==> x[j+1..l-1] = 0 */
588 209029286 : a = mulii(gel(x,1), gcoeff(mul,i,1));
589 1464256611 : for (k = 2; k <= j; k++) a = addii(a, mulii(gel(x,k), gcoeff(mul,i,k)));
590 : /* p | a ? */
591 209021840 : gel(y,i) = dvmdii(a,p,&r); if (signe(r)) return 0;
592 : }
593 : }
594 2391020 : vals = cgetg(l, t_VECSMALL);
595 : /* vals[1] not needed */
596 16692136 : for (j = 2; j < l; j++)
597 : {
598 14300945 : gel(B,j) = Q_primitive_part(gel(B,j), &cx);
599 14300975 : vals[j] = cx? 1 + e * Q_pval(cx, p): 1;
600 : }
601 2391191 : pk = powiu(p, ceildivuu(vmax, e));
602 2391151 : av = avma; y = cgetg(l,t_COL);
603 : /* can compute mod p^ceil((vmax-v)/e) */
604 3942135 : for (v = 1; v < vmax; v++)
605 : { /* we know v_pr(Bj) >= v for all j */
606 1581348 : if (e == 1 || (vmax - v) % e == 0) pk = diviiexact(pk, p);
607 8255708 : for (j = 2; j < l; j++)
608 : {
609 6704757 : GEN x = gel(B,j); if (v < vals[j]) continue;
610 44312953 : for (i = 1; i < l; i++)
611 : {
612 39828351 : pari_sp av2 = avma;
613 39828351 : a = mulii(gel(x,1), gcoeff(mul,i,1));
614 515470848 : for (k = 2; k < l; k++) a = addii(a, mulii(gel(x,k), gcoeff(mul,i,k)));
615 : /* a = (x.t_0)_i; p | a ? */
616 39827610 : a = dvmdii(a,p,&r); if (signe(r)) return v;
617 39797083 : if (lgefint(a) > lgefint(pk)) a = remii(a, pk);
618 39797076 : gel(y,i) = gc_INT(av2, a);
619 : }
620 4484602 : gel(B,j) = y; y = x;
621 4484602 : if (gc_needed(av,3))
622 : {
623 0 : if(DEBUGMEM>1) pari_warn(warnmem,"idealval");
624 0 : (void)gc_all(av,3, &y,&B,&pk);
625 : }
626 : }
627 : }
628 2360787 : return v;
629 : }
630 : /* true nf, x != 0 integral ideal in HNF, cx t_INT or NULL,
631 : * FA integer factorization matrix or NULL. Return partial factorization of
632 : * cx * x above primes in FA (complete factorization if !FA)*/
633 : static GEN
634 294001 : idealHNF_factor_i(GEN nf, GEN x, GEN cx, GEN FA)
635 : {
636 294001 : const long N = lg(x)-1;
637 : long i, j, k, l, v;
638 294001 : GEN vN, vZ, vP, vE, vp = idealHNF_Z_factor_i(x, FA, &vN,&vZ);
639 :
640 294020 : l = lg(vp);
641 294020 : i = cx? expi(cx)+1: 1;
642 294022 : vP = cgetg((l+i-2)*N+1, t_COL);
643 294020 : vE = cgetg((l+i-2)*N+1, t_COL);
644 734569 : for (i = k = 1; i < l; i++)
645 : {
646 440555 : GEN L, p = gel(vp,i);
647 440555 : long Nval = vN[i], Zval = vZ[i], vc = cx? Z_pvalrem(cx,p,&cx): 0;
648 440551 : if (vc)
649 : {
650 48283 : L = idealprimedec(nf,p);
651 48283 : if (is_pm1(cx)) cx = NULL;
652 : }
653 : else
654 392268 : L = idealprimedec_limit_f(nf,p,Nval);
655 1018677 : for (j = 1; Nval && j < lg(L); j++) /* !Nval => only cx contributes */
656 : {
657 578125 : GEN P = gel(L,j);
658 578125 : pari_sp av = avma;
659 578125 : v = idealHNF_val(x, P, Nval, Zval);
660 578102 : set_avma(av);
661 578100 : Nval -= v*pr_get_f(P);
662 578104 : v += vc * pr_get_e(P); if (!v) continue;
663 483601 : gel(vP,k) = P;
664 483601 : gel(vE,k) = utoipos(v); k++;
665 : }
666 491193 : if (vc) for (; j<lg(L); j++)
667 : {
668 50646 : GEN P = gel(L,j);
669 50646 : gel(vP,k) = P;
670 50646 : gel(vE,k) = utoipos(vc * pr_get_e(P)); k++;
671 : }
672 : }
673 294014 : if (cx && !FA)
674 : { /* complete factorization */
675 73364 : GEN f = Z_factor(cx), cP = gel(f,1), cE = gel(f,2);
676 73365 : long lc = lg(cP);
677 159421 : for (i=1; i<lc; i++)
678 : {
679 86056 : GEN p = gel(cP,i), L = idealprimedec(nf,p);
680 86056 : long vc = itos(gel(cE,i));
681 188317 : for (j=1; j<lg(L); j++)
682 : {
683 102261 : GEN P = gel(L,j);
684 102261 : gel(vP,k) = P;
685 102261 : gel(vE,k) = utoipos(vc * pr_get_e(P)); k++;
686 : }
687 : }
688 : }
689 294015 : setlg(vP, k);
690 294015 : setlg(vE, k); return mkmat2(vP, vE);
691 : }
692 : /* true nf, x integral ideal */
693 : static GEN
694 240110 : idealHNF_factor(GEN nf, GEN x, ulong lim)
695 : {
696 240110 : GEN cx, F = NULL;
697 240110 : if (lim)
698 : {
699 : GEN P, E;
700 : long i;
701 : /* strict useless because of prime table */
702 77 : F = absZ_factor_limit(gcoeff(x,1,1), lim);
703 77 : P = gel(F,1);
704 77 : E = gel(F,2);
705 : /* filter out entries > lim */
706 126 : for (i = lg(P)-1; i; i--)
707 126 : if (cmpiu(gel(P,i), lim) <= 0) break;
708 77 : setlg(P, i+1);
709 77 : setlg(E, i+1);
710 : }
711 240110 : x = Q_primitive_part(x, &cx);
712 240095 : return idealHNF_factor_i(nf, x, cx, F);
713 : }
714 : /* c * vector(#L,i,L[i].e), assume results fit in ulong */
715 : static GEN
716 28495 : prV_e_muls(GEN L, long c)
717 : {
718 28495 : long j, l = lg(L);
719 28495 : GEN z = cgetg(l, t_COL);
720 58682 : for (j = 1; j < l; j++) gel(z,j) = stoi(c * pr_get_e(gel(L,j)));
721 28478 : return z;
722 : }
723 : /* true nf, y in Q */
724 : static GEN
725 28485 : Q_nffactor(GEN nf, GEN y, ulong lim)
726 : {
727 : GEN f, P, E;
728 : long l, i;
729 28485 : if (typ(y) == t_INT)
730 : {
731 28457 : if (!signe(y)) pari_err_DOMAIN("idealfactor", "ideal", "=",gen_0,y);
732 28443 : if (is_pm1(y)) return trivial_fact();
733 : }
734 17992 : y = Q_abs_shallow(y);
735 17993 : if (!lim) f = Q_factor(y);
736 : else
737 : {
738 189 : f = Q_factor_limit(y, lim);
739 189 : P = gel(f,1);
740 189 : E = gel(f,2);
741 273 : for (i = lg(P)-1; i > 0; i--)
742 238 : if (abscmpiu(gel(P,i), lim) < 0) break;
743 189 : setlg(P,i+1); setlg(E,i+1);
744 : }
745 18005 : P = gel(f,1); l = lg(P); if (l == 1) return f;
746 17970 : E = gel(f,2);
747 46484 : for (i = 1; i < l; i++)
748 : {
749 28534 : gel(P,i) = idealprimedec(nf, gel(P,i));
750 28481 : gel(E,i) = prV_e_muls(gel(P,i), itos(gel(E,i)));
751 : }
752 17950 : P = shallowconcat1(P); gel(f,1) = P; settyp(P, t_COL);
753 17975 : E = shallowconcat1(E); gel(f,2) = E; return f;
754 : }
755 :
756 : GEN
757 25515 : idealfactor_partial(GEN nf, GEN x, GEN L)
758 : {
759 25515 : pari_sp av = avma;
760 : long i, j, l;
761 : GEN P, E;
762 25515 : if (!L) return idealfactor(nf, x);
763 24661 : if (typ(L) == t_INT) return idealfactor_limit(nf, x, itou(L));
764 24633 : l = lg(L); if (l == 1) return trivial_fact();
765 23863 : P = cgetg(l, t_VEC);
766 89964 : for (i = 1; i < l; i++)
767 : {
768 66101 : GEN p = gel(L,i);
769 66101 : gel(P,i) = typ(p) == t_INT? idealprimedec(nf, p): mkvec(p);
770 : }
771 23863 : P = shallowconcat1(P); settyp(P, t_COL);
772 23863 : P = gen_sort_uniq(P, (void*)&cmp_prime_ideal, &cmp_nodata);
773 23863 : E = cgetg_copy(P, &l);
774 114695 : for (i = j = 1; i < l; i++)
775 : {
776 90832 : long v = idealval(nf, x, gel(P,i));
777 90832 : if (v) { gel(P,j) = gel(P,i); gel(E,j) = stoi(v); j++; }
778 : }
779 23863 : setlg(P,j);
780 23863 : setlg(E,j); return gc_GEN(av, mkmat2(P, E));
781 : }
782 : GEN
783 268711 : idealfactor_limit(GEN nf, GEN x, ulong lim)
784 : {
785 268711 : pari_sp av = avma;
786 : GEN fa, y;
787 268711 : long tx = idealtyp(&x, NULL);
788 :
789 268689 : if (tx == id_PRIME)
790 : {
791 119 : if (lim && abscmpiu(pr_get_p(x), lim) >= 0) return trivial_fact();
792 112 : retmkmat2(mkcolcopy(x), mkcol(gen_1));
793 : }
794 268570 : nf = checknf(nf);
795 268567 : if (tx == id_PRINCIPAL)
796 : {
797 29835 : y = nf_to_scalar_or_basis(nf, x);
798 29835 : if (typ(y) != t_COL) return gc_GEN(av, Q_nffactor(nf, y, lim));
799 : }
800 240083 : y = idealnumden(nf, x);
801 240098 : fa = idealHNF_factor(nf, gel(y,1), lim);
802 240090 : if (!isint1(gel(y,2)))
803 14 : fa = famat_div_shallow(fa, idealHNF_factor(nf, gel(y,2), lim));
804 240089 : fa = gc_GEN(av, fa);
805 240096 : return sort_factor(fa, (void*)&cmp_prime_ideal, &cmp_nodata);
806 : }
807 : GEN
808 268284 : idealfactor(GEN nf, GEN x) { return idealfactor_limit(nf, x, 0); }
809 : GEN
810 189 : gpidealfactor(GEN nf, GEN x, GEN lim)
811 : {
812 189 : ulong L = 0;
813 189 : if (lim)
814 : {
815 70 : if (typ(lim) != t_INT || signe(lim) < 0) pari_err_FLAG("idealfactor");
816 70 : L = itou(lim);
817 : }
818 189 : return idealfactor_limit(nf, x, L);
819 : }
820 :
821 : static GEN
822 7735 : ramified_root(GEN nf, GEN R, GEN A, long n)
823 : {
824 7735 : GEN v, P = gel(idealfactor(nf, R), 1);
825 7735 : long i, l = lg(P);
826 7735 : v = cgetg(l, t_VECSMALL);
827 8393 : for (i = 1; i < l; i++)
828 : {
829 665 : long w = idealval(nf, A, gel(P,i));
830 665 : if (w % n) return NULL;
831 658 : v[i] = w / n;
832 : }
833 7728 : return idealfactorback(nf, P, v, 0);
834 : }
835 : static int
836 7 : ramified_root_simple(GEN nf, long n, GEN P, GEN v)
837 : {
838 7 : long i, l = lg(v);
839 21 : for (i = 1; i < l; i++)
840 : {
841 14 : long w = v[i] % n;
842 14 : if (w)
843 : {
844 7 : GEN vpr = idealprimedec(nf, gel(P,i));
845 7 : long lpr = lg(vpr), j;
846 14 : for (j = 1; j < lpr; j++)
847 : {
848 7 : long e = pr_get_e(gel(vpr,j));
849 7 : if ((e * w) % n) return 0;
850 : }
851 : }
852 : }
853 7 : return 1;
854 : }
855 : /* true nf, n > 1, A a non-zero integral ideal; check whether A is the n-th
856 : * power of an ideal and set *pB to its n-th root if so */
857 : static long
858 7742 : idealsqrtn_int(GEN nf, GEN A, long n, GEN *pB)
859 : {
860 : GEN C, root;
861 : long i, l;
862 :
863 7742 : if (typ(A) == t_MAT && ZM_isscalar(A, NULL)) A = gcoeff(A,1,1);
864 7742 : if (typ(A) == t_INT) /* > 0 */
865 : {
866 5544 : GEN P = nf_get_ramified_primes(nf), v, q;
867 5544 : l = lg(P); v = cgetg(l, t_VECSMALL);
868 24920 : for (i = 1; i < l; i++) v[i] = Z_pvalrem(A, gel(P,i), &A);
869 5544 : C = gen_1;
870 5544 : if (!isint1(A) && !Z_ispowerall(A, n, pB? &C: NULL)) return 0;
871 5544 : if (!pB) return ramified_root_simple(nf, n, P, v);
872 5537 : q = factorback2(P, v);
873 5537 : root = ramified_root(nf, q, q, n);
874 5537 : if (!root) return 0;
875 5537 : if (!equali1(C)) root = isint1(root)? C: ZM_Z_mul(root, C);
876 5537 : *pB = root; return 1;
877 : }
878 : /* compute valuations at ramified primes */
879 2198 : root = ramified_root(nf, idealadd(nf, nf_get_diff(nf), A), A, n);
880 2198 : if (!root) return 0;
881 : /* remove ramified primes */
882 2191 : if (isint1(root))
883 1834 : root = matid(nf_get_degree(nf));
884 : else
885 357 : A = idealdivexact(nf, A, idealpows(nf,root,n));
886 2191 : A = Q_primitive_part(A, &C);
887 2191 : if (C)
888 : {
889 7 : if (!Z_ispowerall(C,n,&C)) return 0;
890 0 : if (pB) root = ZM_Z_mul(root, C);
891 : }
892 :
893 : /* compute final n-th root, at most degree(nf)-1 iterations */
894 2184 : for (i = 0;; i++)
895 2044 : {
896 4228 : GEN J, b, a = gcoeff(A,1,1); /* A \cap Z */
897 4228 : if (is_pm1(a)) break;
898 2072 : if (!Z_ispowerall(a,n,&b)) return 0;
899 2044 : J = idealadd(nf, b, A);
900 2044 : A = idealdivexact(nf, idealpows(nf,J,n), A);
901 : /* div and not divexact here */
902 2044 : if (pB) root = odd(i)? idealdiv(nf, root, J): idealmul(nf, root, J);
903 : }
904 2156 : if (pB) *pB = root;
905 2156 : return 1;
906 : }
907 :
908 : /* A is assumed to be the n-th power of an ideal in nf
909 : returns its n-th root. */
910 : long
911 3899 : idealispower(GEN nf, GEN A, long n, GEN *pB)
912 : {
913 3899 : pari_sp av = avma;
914 : GEN v, N, D;
915 3899 : nf = checknf(nf);
916 3899 : if (n <= 0) pari_err_DOMAIN("idealispower", "n", "<=", gen_0, stoi(n));
917 3899 : if (n == 1) { if (pB) *pB = idealhnf(nf,A); return 1; }
918 3892 : v = idealnumden(nf,A);
919 3892 : if (gequal0(gel(v,1))) { set_avma(av); if (pB) *pB = cgetg(1,t_MAT); return 1; }
920 3892 : if (!idealsqrtn_int(nf, gel(v,1), n, pB? &N: NULL)) return 0;
921 3850 : if (!idealsqrtn_int(nf, gel(v,2), n, pB? &D: NULL)) return 0;
922 3850 : if (pB) *pB = gc_upto(av, idealdiv(nf,N,D)); else set_avma(av);
923 3850 : return 1;
924 : }
925 :
926 : /* x t_INT or integral nonzero ideal in HNF */
927 : static GEN
928 118370 : idealredmodpower_i(GEN nf, GEN x, ulong k, ulong B)
929 : {
930 : GEN cx, y, U, N, F, Q;
931 118370 : if (typ(x) == t_INT)
932 : {
933 63644 : if (!signe(x) || is_pm1(x)) return gen_1;
934 3549 : F = Z_factor_limit(x, B);
935 3549 : gel(F,2) = gdiventgs(gel(F,2), k);
936 3549 : return ginv(factorback(F));
937 : }
938 54726 : N = gcoeff(x,1,1); if (is_pm1(N)) return gen_1;
939 53913 : F = absZ_factor_limit_strict(N, B, &U);
940 53913 : if (U)
941 : {
942 147 : GEN M = powii(gel(U,1), gel(U,2));
943 147 : y = hnfmodid(x, M); /* coprime part to B! */
944 147 : if (!idealispower(nf, y, k, &U)) U = NULL;
945 147 : x = hnfmodid(x, diviiexact(N, M));
946 : }
947 : /* x = B-smooth part of initial x */
948 53913 : x = Q_primitive_part(x, &cx);
949 53913 : F = idealHNF_factor_i(nf, x, cx, F);
950 53913 : gel(F,2) = gdiventgs(gel(F,2), k);
951 53913 : Q = idealfactorback(nf, gel(F,1), gel(F,2), 0);
952 53913 : if (U) Q = idealmul(nf,Q,U);
953 53913 : if (typ(Q) == t_INT) return Q;
954 13361 : y = idealred_elt(nf, idealHNF_inv_Z(nf, Q));
955 13361 : return gdiv(y, gcoeff(Q,1,1));
956 : }
957 : GEN
958 59190 : idealredmodpower(GEN nf, GEN x, ulong n, ulong B)
959 : {
960 59190 : pari_sp av = avma;
961 : GEN a, b;
962 59190 : nf = checknf(nf);
963 59190 : if (!n) pari_err_DOMAIN("idealredmodpower","n", "=", gen_0, gen_0);
964 59190 : x = idealnumden(nf, x);
965 59192 : a = gel(x,1);
966 59192 : if (isintzero(a)) { set_avma(av); return gen_1; }
967 59185 : a = idealredmodpower_i(nf, gel(x,1), n, B);
968 59185 : b = idealredmodpower_i(nf, gel(x,2), n, B);
969 59185 : if (!isint1(b)) a = nf_to_scalar_or_basis(nf, nfdiv(nf, a, b));
970 59185 : return gc_GEN(av, a);
971 : }
972 :
973 : /* P prime ideal in idealprimedec format. Return valuation(A) at P */
974 : long
975 9492185 : idealval(GEN nf, GEN A, GEN P)
976 : {
977 9492185 : pari_sp av = avma;
978 : GEN p, cA;
979 9492185 : long vcA, v, Zval, tx = idealtyp(&A, NULL);
980 :
981 9492193 : if (tx == id_PRINCIPAL) return nfval(nf,A,P);
982 9384896 : checkprid(P);
983 9385146 : if (tx == id_PRIME) return pr_equal(P, A)? 1: 0;
984 : /* id_MAT */
985 9385118 : nf = checknf(nf);
986 9385263 : A = Q_primitive_part(A, &cA);
987 9385427 : p = pr_get_p(P);
988 9385442 : vcA = cA? Q_pval(cA,p): 0;
989 9385443 : if (pr_is_inert(P)) return gc_long(av,vcA);
990 8987358 : Zval = Z_pval(gcoeff(A,1,1), p);
991 8987400 : if (!Zval) v = 0;
992 : else
993 : {
994 3660449 : long Nval = idealHNF_norm_pval(A, p, Zval);
995 3660583 : v = idealHNF_val(A, P, Nval, Zval);
996 : }
997 8987368 : return gc_long(av, vcA? v + vcA*pr_get_e(P): v);
998 : }
999 : GEN
1000 7119 : gpidealval(GEN nf, GEN ix, GEN P)
1001 : {
1002 7119 : long v = idealval(nf,ix,P);
1003 7105 : return v == LONG_MAX? mkoo(): stoi(v);
1004 : }
1005 :
1006 : /* gcd and generalized Bezout */
1007 :
1008 : GEN
1009 83536 : idealadd(GEN nf, GEN x, GEN y)
1010 : {
1011 83536 : pari_sp av = avma;
1012 : long tx, ty;
1013 : GEN z, a, dx, dy, dz;
1014 :
1015 83536 : tx = idealtyp(&x, NULL);
1016 83536 : ty = idealtyp(&y, NULL); nf = checknf(nf);
1017 83536 : if (tx != id_MAT) x = idealhnf_shallow(nf,x);
1018 83536 : if (ty != id_MAT) y = idealhnf_shallow(nf,y);
1019 83536 : if (lg(x) == 1) return gc_GEN(av,y);
1020 82661 : if (lg(y) == 1) return gc_GEN(av,x); /* check for 0 ideal */
1021 82066 : dx = Q_denom(x);
1022 82066 : dy = Q_denom(y); dz = lcmii(dx,dy);
1023 82066 : if (is_pm1(dz)) dz = NULL; else {
1024 5975 : x = Q_muli_to_int(x, dz);
1025 5975 : y = Q_muli_to_int(y, dz);
1026 : }
1027 82066 : a = gcdii(gcoeff(x,1,1), gcoeff(y,1,1));
1028 82066 : if (is_pm1(a))
1029 : {
1030 22889 : long N = lg(x)-1;
1031 22889 : if (!dz) { set_avma(av); return matid(N); }
1032 1085 : return gc_upto(av, scalarmat(ginv(dz), N));
1033 : }
1034 59177 : z = ZM_hnfmodid(shallowconcat(x,y), a);
1035 59177 : if (dz) z = RgM_Rg_div(z,dz);
1036 59177 : return gc_upto(av,z);
1037 : }
1038 :
1039 : static GEN
1040 28 : trivial_merge(GEN x)
1041 28 : { return (lg(x) == 1 || !is_pm1(gcoeff(x,1,1)))? NULL: gen_1; }
1042 : /* true nf */
1043 : static GEN
1044 803379 : _idealaddtoone(GEN nf, GEN x, GEN y, long red)
1045 : {
1046 : GEN a;
1047 803379 : long tx = idealtyp(&x, NULL);
1048 803370 : long ty = idealtyp(&y, NULL);
1049 : long ea;
1050 803372 : if (tx != id_MAT) x = idealhnf_shallow(nf, x);
1051 803377 : if (ty != id_MAT) y = idealhnf_shallow(nf, y);
1052 803377 : if (lg(x) == 1)
1053 14 : a = trivial_merge(y);
1054 803363 : else if (lg(y) == 1)
1055 14 : a = trivial_merge(x);
1056 : else
1057 803349 : a = hnfmerge_get_1(x, y);
1058 803330 : if (!a) pari_err_COPRIME("idealaddtoone",x,y);
1059 803315 : if (red && (ea = gexpo(a)) > 10)
1060 : {
1061 4557 : GEN b = (typ(a) == t_COL)? a: scalarcol_shallow(a, nf_get_degree(nf));
1062 4557 : b = ZC_reducemodlll(b, idealHNF_mul(nf,x,y));
1063 4557 : if (gexpo(b) < ea) a = b;
1064 : }
1065 803315 : return a;
1066 : }
1067 : /* true nf */
1068 : GEN
1069 17633 : idealaddtoone_i(GEN nf, GEN x, GEN y)
1070 17633 : { return _idealaddtoone(nf, x, y, 1); }
1071 : /* true nf */
1072 : GEN
1073 785751 : idealaddtoone_raw(GEN nf, GEN x, GEN y)
1074 785751 : { return _idealaddtoone(nf, x, y, 0); }
1075 :
1076 : GEN
1077 98 : idealaddtoone(GEN nf, GEN x, GEN y)
1078 : {
1079 98 : GEN z = cgetg(3,t_VEC), a;
1080 98 : pari_sp av = avma;
1081 98 : nf = checknf(nf);
1082 98 : a = gc_upto(av, idealaddtoone_i(nf,x,y));
1083 84 : gel(z,1) = a;
1084 84 : gel(z,2) = typ(a) == t_COL? Z_ZC_sub(gen_1,a): subui(1,a);
1085 84 : return z;
1086 : }
1087 :
1088 : /* assume elements of list are integral ideals */
1089 : GEN
1090 35 : idealaddmultoone(GEN nf, GEN list)
1091 : {
1092 35 : pari_sp av = avma;
1093 35 : long N, i, l, nz, tx = typ(list);
1094 : GEN H, U, perm, L;
1095 :
1096 35 : nf = checknf(nf); N = nf_get_degree(nf);
1097 35 : if (!is_vec_t(tx)) pari_err_TYPE("idealaddmultoone",list);
1098 35 : l = lg(list);
1099 35 : L = cgetg(l, t_VEC);
1100 35 : if (l == 1)
1101 0 : pari_err_DOMAIN("idealaddmultoone", "sum(ideals)", "!=", gen_1, L);
1102 35 : nz = 0; /* number of nonzero ideals in L */
1103 98 : for (i=1; i<l; i++)
1104 : {
1105 70 : GEN I = gel(list,i);
1106 70 : if (typ(I) != t_MAT) I = idealhnf_shallow(nf,I);
1107 70 : if (lg(I) != 1)
1108 : {
1109 42 : nz++; RgM_check_ZM(I,"idealaddmultoone");
1110 35 : if (lgcols(I) != N+1) pari_err_TYPE("idealaddmultoone [not an ideal]", I);
1111 : }
1112 63 : gel(L,i) = I;
1113 : }
1114 28 : H = ZM_hnfperm(shallowconcat1(L), &U, &perm);
1115 28 : if (lg(H) == 1 || !equali1(gcoeff(H,1,1)))
1116 7 : pari_err_DOMAIN("idealaddmultoone", "sum(ideals)", "!=", gen_1, L);
1117 49 : for (i=1; i<=N; i++)
1118 49 : if (perm[i] == 1) break;
1119 21 : U = gel(U,(nz-1)*N + i); /* (L[1]|...|L[nz]) U = 1 */
1120 21 : nz = 0;
1121 63 : for (i=1; i<l; i++)
1122 : {
1123 42 : GEN c = gel(L,i);
1124 42 : if (lg(c) == 1)
1125 14 : c = gen_0;
1126 : else {
1127 28 : c = ZM_ZC_mul(c, vecslice(U, nz*N + 1, (nz+1)*N));
1128 28 : nz++;
1129 : }
1130 42 : gel(L,i) = c;
1131 : }
1132 21 : return gc_GEN(av, L);
1133 : }
1134 :
1135 : /* multiplication */
1136 :
1137 : /* x integral ideal (without archimedean component) in HNF form
1138 : * y = [a,alpha] corresponds to the integral ideal aZ_K+alpha Z_K, a in Z,
1139 : * alpha a ZV or a ZM (multiplication table). Multiply them */
1140 : static GEN
1141 970611 : idealHNF_mul_two(GEN nf, GEN x, GEN y)
1142 : {
1143 970611 : GEN m, a = gel(y,1), alpha = gel(y,2);
1144 : long i, N;
1145 :
1146 970611 : if (typ(alpha) != t_MAT)
1147 : {
1148 621379 : alpha = zk_scalar_or_multable(nf, alpha);
1149 621402 : if (typ(alpha) == t_INT) /* e.g. y inert ? 0 should not (but may) occur */
1150 14278 : return signe(a)? ZM_Z_mul(x, gcdii(a, alpha)): cgetg(1,t_MAT);
1151 : }
1152 956356 : N = lg(x)-1; m = cgetg((N<<1)+1,t_MAT);
1153 3791983 : for (i=1; i<=N; i++) gel(m,i) = ZM_ZC_mul(alpha,gel(x,i));
1154 3790822 : for (i=1; i<=N; i++) gel(m,i+N) = ZC_Z_mul(gel(x,i), a);
1155 955499 : return ZM_hnfmodid(m, mulii(a, gcoeff(x,1,1)));
1156 : }
1157 :
1158 : /* Assume x and y are integral in HNF form [NOT extended]. Not memory clean.
1159 : * HACK: ideal in y can be of the form [a,b], a in Z, b in Z_K */
1160 : GEN
1161 473300 : idealHNF_mul(GEN nf, GEN x, GEN y)
1162 : {
1163 : GEN z;
1164 473300 : if (typ(y) == t_VEC)
1165 300062 : z = idealHNF_mul_two(nf,x,y);
1166 : else
1167 : { /* reduce one ideal to two-elt form. The smallest */
1168 173238 : GEN xZ = gcoeff(x,1,1), yZ = gcoeff(y,1,1);
1169 173238 : if (cmpii(xZ, yZ) < 0)
1170 : {
1171 35905 : if (is_pm1(xZ)) return gcopy(y);
1172 20303 : z = idealHNF_mul_two(nf, y, mat_ideal_two_elt(nf,x));
1173 : }
1174 : else
1175 : {
1176 137334 : if (is_pm1(yZ)) return gcopy(x);
1177 36270 : z = idealHNF_mul_two(nf, x, mat_ideal_two_elt(nf,y));
1178 : }
1179 : }
1180 356657 : return z;
1181 : }
1182 :
1183 : /* operations on elements in factored form */
1184 :
1185 : GEN
1186 224855 : famat_mul_shallow(GEN f, GEN g)
1187 : {
1188 224855 : if (typ(f) != t_MAT) f = to_famat_shallow(f,gen_1);
1189 224855 : if (typ(g) != t_MAT) g = to_famat_shallow(g,gen_1);
1190 224855 : if (lgcols(f) == 1) return g;
1191 167595 : if (lgcols(g) == 1) return f;
1192 165673 : return mkmat2(shallowconcat(gel(f,1), gel(g,1)),
1193 165673 : shallowconcat(gel(f,2), gel(g,2)));
1194 : }
1195 : GEN
1196 88237 : famat_mulpow_shallow(GEN f, GEN g, GEN e)
1197 : {
1198 88237 : if (!signe(e)) return f;
1199 54531 : return famat_mul_shallow(f, famat_pow_shallow(g, e));
1200 : }
1201 :
1202 : GEN
1203 147572 : famat_mulpows_shallow(GEN f, GEN g, long e)
1204 : {
1205 147572 : if (e==0) return f;
1206 118874 : return famat_mul_shallow(f, famat_pows_shallow(g, e));
1207 : }
1208 :
1209 : GEN
1210 10304 : famat_div_shallow(GEN f, GEN g)
1211 10304 : { return famat_mul_shallow(f, famat_inv_shallow(g)); }
1212 :
1213 : GEN
1214 376306 : Z_to_famat(GEN x)
1215 : {
1216 : long k;
1217 376306 : if (equali1(x)) return trivial_fact();
1218 192485 : k = Z_isanypower(x, &x) ;
1219 192485 : return to_famat_shallow(x, k? utoi(k): gen_1);
1220 : }
1221 : GEN
1222 197112 : Q_to_famat(GEN x)
1223 : {
1224 197112 : if (typ(x) == t_INT) return Z_to_famat(x);
1225 179194 : return famat_div(Z_to_famat(gel(x,1)), Z_to_famat(gel(x,2)));
1226 : }
1227 : GEN
1228 0 : to_famat(GEN x, GEN y) { retmkmat2(mkcolcopy(x), mkcolcopy(y)); }
1229 : GEN
1230 2768043 : to_famat_shallow(GEN x, GEN y) { return mkmat2(mkcol(x), mkcol(y)); }
1231 :
1232 : /* concat the single elt x; not gconcat since x may be a t_COL */
1233 : static GEN
1234 154142 : append(GEN v, GEN x)
1235 : {
1236 154142 : long i, l = lg(v);
1237 154142 : GEN w = cgetg(l+1, typ(v));
1238 658366 : for (i=1; i<l; i++) gel(w,i) = gcopy(gel(v,i));
1239 154142 : gel(w,i) = gcopy(x); return w;
1240 : }
1241 : /* add x^1 to famat f */
1242 : static GEN
1243 160787 : famat_add(GEN f, GEN x)
1244 : {
1245 160787 : GEN h = cgetg(3,t_MAT);
1246 160787 : if (lgcols(f) == 1)
1247 : {
1248 13186 : gel(h,1) = mkcolcopy(x);
1249 13186 : gel(h,2) = mkcol(gen_1);
1250 : }
1251 : else
1252 : {
1253 147601 : gel(h,1) = append(gel(f,1), x);
1254 147601 : gel(h,2) = gconcat(gel(f,2), gen_1);
1255 : }
1256 160787 : return h;
1257 : }
1258 : /* add x^-1 to famat f */
1259 : static GEN
1260 20927 : famat_sub(GEN f, GEN x)
1261 : {
1262 20927 : GEN h = cgetg(3,t_MAT);
1263 20927 : if (lgcols(f) == 1)
1264 : {
1265 14386 : gel(h,1) = mkcolcopy(x);
1266 14386 : gel(h,2) = mkcol(gen_m1);
1267 : }
1268 : else
1269 : {
1270 6541 : gel(h,1) = append(gel(f,1), x);
1271 6541 : gel(h,2) = gconcat(gel(f,2), gen_m1);
1272 : }
1273 20927 : return h;
1274 : }
1275 :
1276 : GEN
1277 449586 : famat_mul(GEN f, GEN g)
1278 : {
1279 : GEN h;
1280 449586 : if (typ(g) != t_MAT) {
1281 30482 : if (typ(f) == t_MAT) return famat_add(f, g);
1282 0 : h = cgetg(3, t_MAT);
1283 0 : gel(h,1) = mkcol2(gcopy(f), gcopy(g));
1284 0 : gel(h,2) = mkcol2(gen_1, gen_1);
1285 0 : return h;
1286 : }
1287 419104 : if (typ(f) != t_MAT) return famat_add(g, f);
1288 288799 : if (lgcols(f) == 1) return gcopy(g);
1289 265481 : if (lgcols(g) == 1) return gcopy(f);
1290 259494 : h = cgetg(3,t_MAT);
1291 259494 : gel(h,1) = gconcat(gel(f,1), gel(g,1));
1292 259494 : gel(h,2) = gconcat(gel(f,2), gel(g,2));
1293 259494 : return h;
1294 : }
1295 :
1296 : GEN
1297 200128 : famat_div(GEN f, GEN g)
1298 : {
1299 : GEN h;
1300 200128 : if (typ(g) != t_MAT) {
1301 20885 : if (typ(f) == t_MAT) return famat_sub(f, g);
1302 0 : h = cgetg(3, t_MAT);
1303 0 : gel(h,1) = mkcol2(gcopy(f), gcopy(g));
1304 0 : gel(h,2) = mkcol2(gen_1, gen_m1);
1305 0 : return h;
1306 : }
1307 179243 : if (typ(f) != t_MAT) return famat_sub(g, f);
1308 179201 : if (lgcols(f) == 1) return famat_inv(g);
1309 260 : if (lgcols(g) == 1) return gcopy(f);
1310 260 : h = cgetg(3,t_MAT);
1311 260 : gel(h,1) = gconcat(gel(f,1), gel(g,1));
1312 260 : gel(h,2) = gconcat(gel(f,2), gneg(gel(g,2)));
1313 260 : return h;
1314 : }
1315 :
1316 : GEN
1317 22923 : famat_sqr(GEN f)
1318 : {
1319 : GEN h;
1320 22923 : if (typ(f) != t_MAT) return to_famat(f,gen_2);
1321 22923 : if (lgcols(f) == 1) return gcopy(f);
1322 13131 : h = cgetg(3,t_MAT);
1323 13131 : gel(h,1) = gcopy(gel(f,1));
1324 13131 : gel(h,2) = gmul2n(gel(f,2),1);
1325 13131 : return h;
1326 : }
1327 :
1328 : GEN
1329 26945 : famat_inv_shallow(GEN f)
1330 : {
1331 26945 : if (typ(f) != t_MAT) return to_famat_shallow(f,gen_m1);
1332 10444 : if (lgcols(f) == 1) return f;
1333 10444 : return mkmat2(gel(f,1), ZC_neg(gel(f,2)));
1334 : }
1335 : GEN
1336 199322 : famat_inv(GEN f)
1337 : {
1338 199322 : if (typ(f) != t_MAT) return to_famat(f,gen_m1);
1339 199322 : if (lgcols(f) == 1) return gcopy(f);
1340 180803 : retmkmat2(gcopy(gel(f,1)), ZC_neg(gel(f,2)));
1341 : }
1342 : GEN
1343 60642 : famat_pow(GEN f, GEN n)
1344 : {
1345 60642 : if (typ(f) != t_MAT) return to_famat(f,n);
1346 60642 : if (lgcols(f) == 1) return gcopy(f);
1347 60642 : retmkmat2(gcopy(gel(f,1)), ZC_Z_mul(gel(f,2),n));
1348 : }
1349 : GEN
1350 61964 : famat_pow_shallow(GEN f, GEN n)
1351 : {
1352 61964 : if (is_pm1(n)) return signe(n) > 0? f: famat_inv_shallow(f);
1353 35637 : if (typ(f) != t_MAT) return to_famat_shallow(f,n);
1354 7988 : if (lgcols(f) == 1) return f;
1355 6075 : return mkmat2(gel(f,1), ZC_Z_mul(gel(f,2),n));
1356 : }
1357 :
1358 : GEN
1359 151998 : famat_pows_shallow(GEN f, long n)
1360 : {
1361 151998 : if (n==1) return f;
1362 35003 : if (n==-1) return famat_inv_shallow(f);
1363 34996 : if (typ(f) != t_MAT) return to_famat_shallow(f, stoi(n));
1364 26712 : if (lgcols(f) == 1) return f;
1365 26712 : return mkmat2(gel(f,1), ZC_z_mul(gel(f,2),n));
1366 : }
1367 :
1368 : GEN
1369 0 : famat_Z_gcd(GEN M, GEN n)
1370 : {
1371 0 : pari_sp av=avma;
1372 0 : long i, j, l=lgcols(M);
1373 0 : GEN F=cgetg(3,t_MAT);
1374 0 : gel(F,1)=cgetg(l,t_COL);
1375 0 : gel(F,2)=cgetg(l,t_COL);
1376 0 : for (i=1, j=1; i<l; i++)
1377 : {
1378 0 : GEN p = gcoeff(M,i,1);
1379 0 : GEN e = gminsg(Z_pval(n,p),gcoeff(M,i,2));
1380 0 : if (signe(e))
1381 : {
1382 0 : gcoeff(F,j,1)=p;
1383 0 : gcoeff(F,j,2)=e;
1384 0 : j++;
1385 : }
1386 : }
1387 0 : setlg(gel(F,1),j); setlg(gel(F,2),j);
1388 0 : return gc_GEN(av,F);
1389 : }
1390 :
1391 : /* x assumed to be a t_MATs (factorization matrix), or compatible with
1392 : * the element_* functions. */
1393 : static GEN
1394 33885 : ext_sqr(GEN nf, GEN x)
1395 33885 : { return (typ(x)==t_MAT)? famat_sqr(x): nfsqr(nf, x); }
1396 : static GEN
1397 60786 : ext_mul(GEN nf, GEN x, GEN y)
1398 60786 : { return (typ(x)==t_MAT)? famat_mul(x,y): nfmul(nf, x, y); }
1399 : static GEN
1400 20381 : ext_inv(GEN nf, GEN x)
1401 20381 : { return (typ(x)==t_MAT)? famat_inv(x): nfinv(nf, x); }
1402 : static GEN
1403 0 : ext_pow(GEN nf, GEN x, GEN n)
1404 0 : { return (typ(x)==t_MAT)? famat_pow(x,n): nfpow(nf, x, n); }
1405 :
1406 : GEN
1407 0 : famat_to_nf(GEN nf, GEN f)
1408 : {
1409 : GEN t, x, e;
1410 : long i;
1411 0 : if (lgcols(f) == 1) return gen_1;
1412 0 : x = gel(f,1);
1413 0 : e = gel(f,2);
1414 0 : t = nfpow(nf, gel(x,1), gel(e,1));
1415 0 : for (i=lg(x)-1; i>1; i--)
1416 0 : t = nfmul(nf, t, nfpow(nf, gel(x,i), gel(e,i)));
1417 0 : return t;
1418 : }
1419 :
1420 : GEN
1421 0 : famat_idealfactor(GEN nf, GEN x)
1422 : {
1423 : long i, l;
1424 0 : GEN g = gel(x,1), e = gel(x,2), h = cgetg_copy(g, &l);
1425 0 : for (i = 1; i < l; i++) gel(h,i) = idealfactor(nf, gel(g,i));
1426 0 : h = famat_reduce(famatV_factorback(h,e));
1427 0 : return sort_factor(h, (void*)&cmp_prime_ideal, &cmp_nodata);
1428 : }
1429 :
1430 : GEN
1431 306316 : famat_reduce(GEN fa)
1432 : {
1433 : GEN E, G, L, g, e;
1434 : long i, k, l;
1435 :
1436 306316 : if (typ(fa) != t_MAT || lgcols(fa) == 1) return fa;
1437 295394 : g = gel(fa,1); l = lg(g);
1438 295394 : e = gel(fa,2);
1439 295394 : L = gen_indexsort(g, (void*)&cmp_universal, &cmp_nodata);
1440 295393 : G = cgetg(l, t_COL);
1441 295393 : E = cgetg(l, t_COL);
1442 : /* merge */
1443 2356836 : for (k=i=1; i<l; i++,k++)
1444 : {
1445 2061445 : gel(G,k) = gel(g,L[i]);
1446 2061445 : gel(E,k) = gel(e,L[i]);
1447 2061445 : if (k > 1 && gidentical(gel(G,k), gel(G,k-1)))
1448 : {
1449 838772 : gel(E,k-1) = addii(gel(E,k), gel(E,k-1));
1450 838771 : k--;
1451 : }
1452 : }
1453 : /* kill 0 exponents */
1454 295391 : l = k;
1455 1518063 : for (k=i=1; i<l; i++)
1456 1222672 : if (!gequal0(gel(E,i)))
1457 : {
1458 1198081 : gel(G,k) = gel(G,i);
1459 1198081 : gel(E,k) = gel(E,i); k++;
1460 : }
1461 295391 : setlg(G, k);
1462 295391 : setlg(E, k); return mkmat2(G,E);
1463 : }
1464 : GEN
1465 63 : matreduce(GEN f)
1466 63 : { pari_sp av = avma;
1467 63 : switch(typ(f))
1468 : {
1469 21 : case t_VEC: case t_COL:
1470 : {
1471 21 : GEN e; f = vec_reduce(f, &e); settyp(f, t_COL);
1472 21 : return gc_GEN(av, mkmat2(f, zc_to_ZC(e)));
1473 : }
1474 35 : case t_MAT:
1475 35 : if (lg(f) == 3) break;
1476 : default:
1477 14 : pari_err_TYPE("matreduce", f);
1478 : }
1479 28 : if (typ(gel(f,1)) == t_VECSMALL)
1480 0 : f = famatsmall_reduce(f);
1481 : else
1482 : {
1483 28 : if (!RgV_is_ZV(gel(f,2))) pari_err_TYPE("matreduce",f);
1484 21 : f = famat_reduce(f);
1485 : }
1486 21 : return gc_GEN(av, f);
1487 : }
1488 :
1489 : GEN
1490 189693 : famatsmall_reduce(GEN fa)
1491 : {
1492 : GEN E, G, L, g, e;
1493 : long i, k, l;
1494 189693 : if (lgcols(fa) == 1) return fa;
1495 189693 : g = gel(fa,1); l = lg(g);
1496 189693 : e = gel(fa,2);
1497 189693 : L = vecsmall_indexsort(g);
1498 189693 : G = cgetg(l, t_VECSMALL);
1499 189693 : E = cgetg(l, t_VECSMALL);
1500 : /* merge */
1501 512985 : for (k=i=1; i<l; i++,k++)
1502 : {
1503 323292 : G[k] = g[L[i]];
1504 323292 : E[k] = e[L[i]];
1505 323292 : if (k > 1 && G[k] == G[k-1])
1506 : {
1507 7956 : E[k-1] += E[k];
1508 7956 : k--;
1509 : }
1510 : }
1511 : /* kill 0 exponents */
1512 189693 : l = k;
1513 505029 : for (k=i=1; i<l; i++)
1514 315336 : if (E[i])
1515 : {
1516 311688 : G[k] = G[i];
1517 311688 : E[k] = E[i]; k++;
1518 : }
1519 189693 : setlg(G, k);
1520 189693 : setlg(E, k); return mkmat2(G,E);
1521 : }
1522 :
1523 : GEN
1524 66958 : famat_remove_trivial(GEN fa)
1525 : {
1526 66958 : GEN P, E, p = gel(fa,1), e = gel(fa,2);
1527 66958 : long j, k, l = lg(p);
1528 66958 : P = cgetg(l, t_COL);
1529 66958 : E = cgetg(l, t_COL);
1530 1821763 : for (j = k = 1; j < l; j++)
1531 1754805 : if (signe(gel(e,j))) { gel(P,k) = gel(p,j); gel(E,k++) = gel(e,j); }
1532 66958 : setlg(P, k); setlg(E, k); return mkmat2(P,E);
1533 : }
1534 :
1535 : GEN
1536 13765 : famatV_factorback(GEN v, GEN e)
1537 : {
1538 13765 : long i, l = lg(e);
1539 : GEN V;
1540 13765 : if (l == 1) return trivial_fact();
1541 13380 : V = signe(gel(e,1))? famat_pow_shallow(gel(v,1), gel(e,1)): trivial_fact();
1542 56655 : for (i = 2; i < l; i++) V = famat_mulpow_shallow(V, gel(v,i), gel(e,i));
1543 13380 : return V;
1544 : }
1545 :
1546 : GEN
1547 56273 : famatV_zv_factorback(GEN v, GEN e)
1548 : {
1549 56273 : long i, l = lg(e);
1550 : GEN V;
1551 56273 : if (l == 1) return trivial_fact();
1552 53816 : V = uel(e,1)? famat_pows_shallow(gel(v,1), uel(e,1)): trivial_fact();
1553 177730 : for (i = 2; i < l; i++) V = famat_mulpows_shallow(V, gel(v,i), uel(e,i));
1554 53816 : return V;
1555 : }
1556 :
1557 : GEN
1558 568161 : ZM_famat_limit(GEN fa, GEN limit)
1559 : {
1560 : pari_sp av;
1561 : GEN E, G, g, e, r;
1562 : long i, k, l, n, lG;
1563 :
1564 568161 : if (lgcols(fa) == 1) return fa;
1565 568154 : g = gel(fa,1); l = lg(g);
1566 568154 : e = gel(fa,2);
1567 1137456 : for(n=0, i=1; i<l; i++)
1568 569302 : if (cmpii(gel(g,i),limit)<=0) n++;
1569 568154 : lG = n<l-1 ? n+2 : n+1;
1570 568154 : G = cgetg(lG, t_COL);
1571 568154 : E = cgetg(lG, t_COL);
1572 568154 : av = avma;
1573 1137456 : for (i=1, k=1, r = gen_1; i<l; i++)
1574 : {
1575 569302 : if (cmpii(gel(g,i),limit)<=0)
1576 : {
1577 569169 : gel(G,k) = gel(g,i);
1578 569169 : gel(E,k) = gel(e,i);
1579 569169 : k++;
1580 133 : } else r = mulii(r, powii(gel(g,i), gel(e,i)));
1581 : }
1582 568154 : if (k<i)
1583 : {
1584 133 : gel(G, k) = gc_INT(av, r);
1585 133 : gel(E, k) = gen_1;
1586 : }
1587 568154 : return mkmat2(G,E);
1588 : }
1589 :
1590 : /* assume pr has degree 1 and coprime to Q_denom(x) */
1591 : static GEN
1592 123032 : to_Fp_coprime(GEN nf, GEN x, GEN modpr)
1593 : {
1594 123032 : GEN d, r, p = modpr_get_p(modpr);
1595 123032 : x = nf_to_scalar_or_basis(nf,x);
1596 123032 : if (typ(x) != t_COL) return Rg_to_Fp(x,p);
1597 121394 : x = Q_remove_denom(x, &d);
1598 121394 : r = zk_to_Fq(x, modpr);
1599 121394 : if (d) r = Fp_div(r, d, p);
1600 121394 : return r;
1601 : }
1602 :
1603 : /* pr coprime to all denominators occurring in x */
1604 : static GEN
1605 693 : famat_to_Fp_coprime(GEN nf, GEN x, GEN modpr)
1606 : {
1607 693 : GEN p = modpr_get_p(modpr);
1608 693 : GEN t = NULL, g = gel(x,1), e = gel(x,2), q = subiu(p,1);
1609 693 : long i, l = lg(g);
1610 4571 : for (i = 1; i < l; i++)
1611 : {
1612 3878 : GEN n = modii(gel(e,i), q);
1613 3878 : if (signe(n))
1614 : {
1615 3871 : GEN h = to_Fp_coprime(nf, gel(g,i), modpr);
1616 3871 : h = Fp_pow(h, n, p);
1617 3871 : t = t? Fp_mul(t, h, p): h;
1618 : }
1619 : }
1620 693 : return t? modii(t, p): gen_1;
1621 : }
1622 :
1623 : /* cf famat_to_nf_modideal_coprime, modpr attached to prime of degree 1 */
1624 : GEN
1625 119854 : nf_to_Fp_coprime(GEN nf, GEN x, GEN modpr)
1626 : {
1627 693 : return typ(x)==t_MAT? famat_to_Fp_coprime(nf, x, modpr)
1628 120547 : : to_Fp_coprime(nf, x, modpr);
1629 : }
1630 :
1631 : static long
1632 4325469 : zk_pvalrem(GEN x, GEN p, GEN *py)
1633 4325469 : { return (typ(x) == t_INT)? Z_pvalrem(x, p, py): ZV_pvalrem(x, p, py); }
1634 : /* x a QC or Q. Return a ZC or Z, whose content is coprime to Z. Set v, dx
1635 : * such that x = p^v (newx / dx); dx = NULL if 1 */
1636 : static GEN
1637 4371263 : nf_remove_denom_p(GEN nf, GEN x, GEN p, GEN *pdx, long *pv)
1638 : {
1639 : long vcx;
1640 : GEN dx;
1641 4371263 : x = nf_to_scalar_or_basis(nf, x);
1642 4371265 : x = Q_remove_denom(x, &dx);
1643 4371270 : if (dx)
1644 : {
1645 61201 : vcx = - Z_pvalrem(dx, p, &dx);
1646 61201 : if (!vcx) vcx = zk_pvalrem(x, p, &x);
1647 61201 : if (isint1(dx)) dx = NULL;
1648 : }
1649 : else
1650 : {
1651 4310069 : vcx = zk_pvalrem(x, p, &x);
1652 4310075 : dx = NULL;
1653 : }
1654 4371276 : *pv = vcx;
1655 4371276 : *pdx = dx; return x;
1656 : }
1657 : /* x = b^e/p^(e-1) in Z_K; x = 0 mod p/pr^e, (x,pr) = 1. Return NULL
1658 : * if p inert (instead of 1) */
1659 : static GEN
1660 94879 : p_makecoprime(GEN pr)
1661 : {
1662 94879 : GEN B = pr_get_tau(pr), b;
1663 : long i, e;
1664 :
1665 94879 : if (typ(B) == t_INT) return NULL;
1666 72122 : b = gel(B,1); /* B = multiplication table by b */
1667 72122 : e = pr_get_e(pr);
1668 72122 : if (e == 1) return b;
1669 : /* one could also divide (exactly) by p in each iteration */
1670 45198 : for (i = 1; i < e; i++) b = ZM_ZC_mul(B, b);
1671 22085 : return ZC_Z_divexact(b, powiu(pr_get_p(pr), e-1));
1672 : }
1673 :
1674 : /* Compute A = prod g[i]^e[i] mod pr^k, assuming (A, pr) = 1.
1675 : * Method: modify each g[i] so that it becomes coprime to pr,
1676 : * g[i] *= (b/p)^v_pr(g[i]), where b/p = pr^(-1) times something integral
1677 : * and prime to p; globally, we multiply by (b/p)^v_pr(A) = 1.
1678 : * Optimizations:
1679 : * 1) remove all powers of p from contents, and consider extra generator p^vp;
1680 : * modified as p * (b/p)^e = b^e / p^(e-1)
1681 : * 2) remove denominators, coprime to p, by multiplying by inverse mod prk\cap Z
1682 : *
1683 : * EX = multiple of exponent of (O_K / pr^k)^* used to reduce the product in
1684 : * case the e[i] are large */
1685 : GEN
1686 2220895 : famat_makecoprime(GEN nf, GEN g, GEN e, GEN pr, GEN prk, GEN EX)
1687 : {
1688 2220895 : GEN G, E, t, vp = NULL, p = pr_get_p(pr), prkZ = gcoeff(prk, 1,1);
1689 2220891 : long i, l = lg(g);
1690 :
1691 2220891 : G = cgetg(l+1, t_VEC);
1692 2220906 : E = cgetg(l+1, t_VEC); /* l+1: room for "modified p" */
1693 6592166 : for (i=1; i < l; i++)
1694 : {
1695 : long vcx;
1696 4371263 : GEN dx, x = nf_remove_denom_p(nf, gel(g,i), p, &dx, &vcx);
1697 4371270 : if (vcx) /* = v_p(content(g[i])) */
1698 : {
1699 142050 : GEN a = mulsi(vcx, gel(e,i));
1700 142064 : vp = vp? addii(vp, a): a;
1701 : }
1702 : /* x integral, content coprime to p; dx coprime to p */
1703 4371284 : if (typ(x) == t_INT)
1704 : { /* x coprime to p, hence to pr */
1705 1133548 : x = modii(x, prkZ);
1706 1133548 : if (dx) x = Fp_div(x, dx, prkZ);
1707 : }
1708 : else
1709 : {
1710 3237736 : (void)ZC_nfvalrem(x, pr, &x); /* x *= (b/p)^v_pr(x) */
1711 3237690 : x = ZC_hnfrem(FpC_red(x,prkZ), prk);
1712 3237728 : if (dx) x = FpC_Fp_mul(x, Fp_inv(dx,prkZ), prkZ);
1713 : }
1714 4371253 : gel(G,i) = x;
1715 4371253 : gel(E,i) = gel(e,i);
1716 : }
1717 :
1718 2220903 : t = vp? p_makecoprime(pr): NULL;
1719 2220909 : if (!t)
1720 : { /* no need for extra generator */
1721 2148863 : setlg(G,l);
1722 2148864 : setlg(E,l);
1723 : }
1724 : else
1725 : {
1726 72046 : gel(G,i) = FpC_red(t, prkZ);
1727 72046 : gel(E,i) = vp;
1728 : }
1729 2220909 : return famat_to_nf_modideal_coprime(nf, G, E, prk, EX);
1730 : }
1731 :
1732 : /* simplified version of famat_makecoprime for X = SUnits[1] */
1733 : GEN
1734 98 : sunits_makecoprime(GEN X, GEN pr, GEN prk)
1735 : {
1736 98 : GEN G, p = pr_get_p(pr), prkZ = gcoeff(prk,1,1);
1737 98 : long i, l = lg(X);
1738 :
1739 98 : G = cgetg(l, t_VEC);
1740 9205 : for (i = 1; i < l; i++)
1741 : {
1742 9107 : GEN x = gel(X,i);
1743 9107 : if (typ(x) == t_INT) /* a prime */
1744 1491 : x = equalii(x,p)? p_makecoprime(pr): modii(x, prkZ);
1745 : else
1746 : {
1747 7616 : (void)ZC_nfvalrem(x, pr, &x); /* x *= (b/p)^v_pr(x) */
1748 7616 : x = ZC_hnfrem(FpC_red(x,prkZ), prk);
1749 : }
1750 9107 : gel(G,i) = x;
1751 : }
1752 98 : return G;
1753 : }
1754 :
1755 : /* prod g[i]^e[i] mod bid, assume (g[i], id) = 1 and 1 < lg(g) <= lg(e) */
1756 : GEN
1757 18606 : famat_to_nf_moddivisor(GEN nf, GEN g, GEN e, GEN bid)
1758 : {
1759 18606 : GEN t, cyc = bid_get_cyc(bid);
1760 18606 : if (lg(cyc) == 1)
1761 0 : t = gen_1;
1762 : else
1763 18606 : t = famat_to_nf_modideal_coprime(nf, g, e, bid_get_ideal(bid),
1764 : cyc_get_expo(cyc));
1765 18606 : return set_sign_mod_divisor(nf, mkmat2(g,e), t, bid_get_sarch(bid));
1766 : }
1767 :
1768 : GEN
1769 16016383 : vecmul(GEN x, GEN y)
1770 : {
1771 16016383 : if (!is_vec_t(typ(x))) return gmul(x,y);
1772 3530528 : pari_APPLY_same(vecmul(gel(x,i), gel(y,i)))
1773 : }
1774 :
1775 : GEN
1776 185983 : vecsqr(GEN x)
1777 : {
1778 185983 : if (!is_vec_t(typ(x))) return gsqr(x);
1779 46606 : pari_APPLY_same(vecsqr(gel(x,i)))
1780 : }
1781 :
1782 : GEN
1783 826 : vecinv(GEN x)
1784 : {
1785 826 : if (!is_vec_t(typ(x))) return ginv(x);
1786 56 : pari_APPLY_same(vecinv(gel(x,i)))
1787 : }
1788 :
1789 : GEN
1790 0 : vecpow(GEN x, GEN n)
1791 : {
1792 0 : if (!is_vec_t(typ(x))) return powgi(x,n);
1793 0 : pari_APPLY_same(vecpow(gel(x,i), n))
1794 : }
1795 :
1796 : GEN
1797 903 : vecdiv(GEN x, GEN y)
1798 : {
1799 903 : if (!is_vec_t(typ(x))) return gdiv(x,y);
1800 903 : pari_APPLY_same(vecdiv(gel(x,i), gel(y,i)))
1801 : }
1802 :
1803 : /* A ideal as a square t_MAT */
1804 : static GEN
1805 266235 : idealmulelt(GEN nf, GEN x, GEN A)
1806 : {
1807 : long i, lx;
1808 : GEN dx, dA, D;
1809 266235 : if (lg(A) == 1) return cgetg(1, t_MAT);
1810 266235 : x = nf_to_scalar_or_basis(nf,x);
1811 266235 : if (typ(x) != t_COL)
1812 71772 : return isintzero(x)? cgetg(1,t_MAT): RgM_Rg_mul(A, Q_abs_shallow(x));
1813 194463 : x = Q_remove_denom(x, &dx);
1814 194463 : A = Q_remove_denom(A, &dA);
1815 194463 : x = zk_multable(nf, x);
1816 194463 : D = mulii(zkmultable_capZ(x), gcoeff(A,1,1));
1817 194463 : x = zkC_multable_mul(A, x);
1818 194463 : settyp(x, t_MAT); lx = lg(x);
1819 : /* x may contain scalars (at most 1 since the ideal is nonzero)*/
1820 751434 : for (i=1; i<lx; i++)
1821 566234 : if (typ(gel(x,i)) == t_INT)
1822 : {
1823 9263 : if (i > 1) swap(gel(x,1), gel(x,i)); /* help HNF */
1824 9263 : gel(x,1) = scalarcol_shallow(gel(x,1), lx-1);
1825 9263 : break;
1826 : }
1827 194463 : x = ZM_hnfmodid(x, D);
1828 194463 : dx = mul_denom(dx,dA);
1829 194463 : return dx? gdiv(x,dx): x;
1830 : }
1831 :
1832 : /* nf a true nf, tx <= ty */
1833 : static GEN
1834 511798 : idealmul_aux(GEN nf, GEN x, GEN y, long tx, long ty)
1835 : {
1836 : GEN z, cx, cy;
1837 511798 : switch(tx)
1838 : {
1839 304269 : case id_PRINCIPAL:
1840 304269 : switch(ty)
1841 : {
1842 37586 : case id_PRINCIPAL:
1843 37586 : return idealhnf_principal(nf, nfmul(nf,x,y));
1844 448 : case id_PRIME:
1845 : {
1846 448 : GEN p = pr_get_p(y), pi = pr_get_gen(y), cx;
1847 448 : if (pr_is_inert(y)) return RgM_Rg_mul(idealhnf_principal(nf,x),p);
1848 :
1849 217 : x = nf_to_scalar_or_basis(nf, x);
1850 217 : switch(typ(x))
1851 : {
1852 203 : case t_INT:
1853 203 : if (!signe(x)) return cgetg(1,t_MAT);
1854 203 : return ZM_Z_mul(pr_hnf(nf,y), absi_shallow(x));
1855 7 : case t_FRAC:
1856 7 : return RgM_Rg_mul(pr_hnf(nf,y), Q_abs_shallow(x));
1857 : }
1858 : /* t_COL */
1859 7 : x = Q_primitive_part(x, &cx);
1860 7 : x = zk_multable(nf, x);
1861 7 : z = shallowconcat(ZM_Z_mul(x,p), ZM_ZC_mul(x,pi));
1862 7 : z = ZM_hnfmodid(z, mulii(p, zkmultable_capZ(x)));
1863 7 : return cx? ZM_Q_mul(z, cx): z;
1864 : }
1865 266235 : default: /* id_MAT */
1866 266235 : return idealmulelt(nf, x,y);
1867 : }
1868 42635 : case id_PRIME:
1869 42635 : if (ty==id_PRIME)
1870 4347 : { y = pr_hnf(nf,y); cy = NULL; }
1871 : else
1872 38288 : y = Q_primitive_part(y, &cy);
1873 42635 : y = idealHNF_mul_two(nf,y,x);
1874 42635 : return cy? ZM_Q_mul(y,cy): y;
1875 :
1876 164894 : default: /* id_MAT */
1877 : {
1878 164894 : long N = nf_get_degree(nf);
1879 164894 : if (lg(x)-1 != N || lg(y)-1 != N) pari_err_DIM("idealmul");
1880 164880 : x = Q_primitive_part(x, &cx);
1881 164880 : y = Q_primitive_part(y, &cy); cx = mul_content(cx,cy);
1882 164880 : y = idealHNF_mul(nf,x,y);
1883 164880 : return cx? ZM_Q_mul(y,cx): y;
1884 : }
1885 : }
1886 : }
1887 :
1888 : /* output the ideal product x.y */
1889 : GEN
1890 511798 : idealmul(GEN nf, GEN x, GEN y)
1891 : {
1892 : pari_sp av;
1893 : GEN res, ax, ay, z;
1894 511798 : long tx = idealtyp(&x,&ax);
1895 511798 : long ty = idealtyp(&y,&ay), f;
1896 511798 : if (tx>ty) { swap(ax,ay); swap(x,y); lswap(tx,ty); }
1897 511798 : f = (ax||ay); res = f? cgetg(3,t_VEC): NULL; /*product is an extended ideal*/
1898 511798 : av = avma;
1899 511798 : z = gc_upto(av, idealmul_aux(checknf(nf), x,y, tx,ty));
1900 511784 : if (!f) return z;
1901 28024 : if (ax && ay)
1902 26526 : ax = ext_mul(nf, ax, ay);
1903 : else
1904 1498 : ax = gcopy(ax? ax: ay);
1905 28024 : gel(res,1) = z; gel(res,2) = ax; return res;
1906 : }
1907 :
1908 : /* Return x, integral in 2-elt form, such that pr^2 = c * x. cf idealpowprime
1909 : * nf = true nf */
1910 : static GEN
1911 317317 : idealsqrprime(GEN nf, GEN pr, GEN *pc)
1912 : {
1913 317317 : GEN p = pr_get_p(pr), q, gen;
1914 317316 : long e = pr_get_e(pr), f = pr_get_f(pr);
1915 :
1916 317317 : q = (e == 1)? sqri(p): p;
1917 317307 : if (e <= 2 && e * f == nf_get_degree(nf))
1918 : { /* pr^e = (p) */
1919 45682 : *pc = q;
1920 45682 : return mkvec2(gen_1,gen_0);
1921 : }
1922 271625 : gen = nfsqr(nf, pr_get_gen(pr));
1923 271629 : gen = FpC_red(gen, q);
1924 271623 : *pc = NULL;
1925 271623 : return mkvec2(q, gen);
1926 : }
1927 : /* cf idealpow_aux */
1928 : static GEN
1929 39149 : idealsqr_aux(GEN nf, GEN x, long tx)
1930 : {
1931 39149 : GEN T = nf_get_pol(nf), m, cx, a, alpha;
1932 39149 : long N = degpol(T);
1933 39149 : switch(tx)
1934 : {
1935 385 : case id_PRINCIPAL:
1936 385 : return idealhnf_principal(nf, nfsqr(nf,x));
1937 10773 : case id_PRIME:
1938 10773 : if (pr_is_inert(x)) return scalarmat(sqri(gel(x,1)), N);
1939 10605 : x = idealsqrprime(nf, x, &cx);
1940 10605 : x = idealhnf_two(nf,x);
1941 10605 : return cx? ZM_Z_mul(x, cx): x;
1942 27991 : default:
1943 27991 : x = Q_primitive_part(x, &cx);
1944 27991 : a = mat_ideal_two_elt(nf,x); alpha = gel(a,2); a = gel(a,1);
1945 27991 : alpha = nfsqr(nf,alpha);
1946 27991 : m = zk_scalar_or_multable(nf, alpha);
1947 27991 : if (typ(m) == t_INT) {
1948 1642 : x = gcdii(sqri(a), m);
1949 1642 : if (cx) x = gmul(x, gsqr(cx));
1950 1642 : x = scalarmat(x, N);
1951 : }
1952 : else
1953 : { /* could use gcdii(sqri(a), zkmultable_capZ(m)), but costly */
1954 26349 : x = ZM_hnfmodid(m, sqri(a));
1955 26349 : if (cx) cx = gsqr(cx);
1956 26349 : if (cx) x = ZM_Q_mul(x, cx);
1957 : }
1958 27991 : return x;
1959 : }
1960 : }
1961 : GEN
1962 39149 : idealsqr(GEN nf, GEN x)
1963 : {
1964 : pari_sp av;
1965 : GEN res, ax, z;
1966 39149 : long tx = idealtyp(&x,&ax);
1967 39149 : res = ax? cgetg(3,t_VEC): NULL; /*product is an extended ideal*/
1968 39149 : av = avma;
1969 39149 : z = gc_upto(av, idealsqr_aux(checknf(nf), x, tx));
1970 39149 : if (!ax) return z;
1971 33885 : gel(res,1) = z;
1972 33885 : gel(res,2) = ext_sqr(nf, ax); return res;
1973 : }
1974 :
1975 : /* norm of an ideal */
1976 : GEN
1977 105700 : idealnorm(GEN nf, GEN x)
1978 : {
1979 : pari_sp av;
1980 : long tx;
1981 :
1982 105700 : switch(idealtyp(&x, NULL))
1983 : {
1984 4914 : case id_PRIME: return pr_norm(x);
1985 11088 : case id_MAT: return RgM_det_triangular(x);
1986 : }
1987 : /* id_PRINCIPAL */
1988 89698 : nf = checknf(nf); av = avma;
1989 89698 : x = nfnorm(nf, x);
1990 89698 : tx = typ(x);
1991 89698 : if (tx == t_INT) return gc_INT(av, absi(x));
1992 420 : if (tx != t_FRAC) pari_err_TYPE("idealnorm",x);
1993 420 : return gc_upto(av, Q_abs(x));
1994 : }
1995 :
1996 : /* x \cap Z */
1997 : GEN
1998 3031 : idealdown(GEN nf, GEN x)
1999 : {
2000 3031 : pari_sp av = avma;
2001 : GEN y, c;
2002 3031 : switch(idealtyp(&x, NULL))
2003 : {
2004 7 : case id_PRIME: return icopy(pr_get_p(x));
2005 2121 : case id_MAT: return gcopy(gcoeff(x,1,1));
2006 : }
2007 : /* id_PRINCIPAL */
2008 903 : nf = checknf(nf); av = avma;
2009 903 : x = nf_to_scalar_or_basis(nf, x);
2010 903 : if (is_rational_t(typ(x))) return Q_abs(x);
2011 14 : x = Q_primitive_part(x, &c);
2012 14 : y = zkmultable_capZ(zk_multable(nf, x));
2013 14 : return gc_GEN(av, mul_content(c, y));
2014 : }
2015 :
2016 : /* true nf */
2017 : static GEN
2018 42 : idealismaximal_int(GEN nf, GEN p)
2019 : {
2020 : GEN L;
2021 42 : if (!BPSW_psp(p)) return NULL;
2022 77 : if (!dvdii(nf_get_index(nf), p) &&
2023 49 : !FpX_is_irred(FpX_red(nf_get_pol(nf),p), p)) return NULL;
2024 28 : L = idealprimedec(nf, p);
2025 28 : return (lg(L) == 2 && pr_get_e(gel(L,1)) == 1)? gel(L,1): NULL;
2026 : }
2027 : /* true nf */
2028 : static GEN
2029 21 : idealismaximal_mat(GEN nf, GEN x)
2030 : {
2031 : GEN p, c, L;
2032 : long i, l, f;
2033 21 : x = Q_primitive_part(x, &c);
2034 21 : p = gcoeff(x,1,1);
2035 21 : if (c)
2036 : {
2037 7 : if (typ(c) == t_FRAC || !equali1(p)) return NULL;
2038 7 : return idealismaximal_int(nf, c);
2039 : }
2040 14 : if (!BPSW_psp(p)) return NULL;
2041 14 : l = lg(x); f = 1;
2042 35 : for (i = 2; i < l; i++)
2043 : {
2044 21 : c = gcoeff(x,i,i);
2045 21 : if (equalii(c, p)) f++; else if (!equali1(c)) return NULL;
2046 : }
2047 14 : L = idealprimedec_limit_f(nf, p, f);
2048 28 : for (i = lg(L)-1; i; i--)
2049 : {
2050 28 : GEN pr = gel(L,i);
2051 28 : if (pr_get_f(pr) != f) break;
2052 28 : if (idealval(nf, x, pr) == 1) return pr;
2053 : }
2054 0 : return NULL;
2055 : }
2056 : /* true nf */
2057 : static GEN
2058 77 : idealismaximal_i(GEN nf, GEN x)
2059 : {
2060 : GEN L, p, pr, c;
2061 : long i, l;
2062 77 : switch(idealtyp(&x, NULL))
2063 : {
2064 7 : case id_PRIME: return x;
2065 21 : case id_MAT: return idealismaximal_mat(nf, x);
2066 : }
2067 : /* id_PRINCIPAL */
2068 49 : x = nf_to_scalar_or_basis(nf, x);
2069 49 : switch(typ(x))
2070 : {
2071 35 : case t_INT: return idealismaximal_int(nf, absi_shallow(x));
2072 0 : case t_FRAC: return NULL;
2073 : }
2074 14 : x = Q_primitive_part(x, &c);
2075 14 : if (c) return NULL;
2076 14 : p = zkmultable_capZ(zk_multable(nf, x));
2077 14 : if (!BPSW_psp(p)) return NULL;
2078 7 : L = idealprimedec(nf, p); l = lg(L); pr = NULL;
2079 21 : for (i = 1; i < l; i++)
2080 : {
2081 14 : long v = ZC_nfval(x, gel(L,i));
2082 14 : if (v > 1 || (v && pr)) return NULL;
2083 14 : pr = gel(L,i);
2084 : }
2085 7 : return pr;
2086 : }
2087 : GEN
2088 77 : idealismaximal(GEN nf, GEN x)
2089 : {
2090 77 : pari_sp av = avma;
2091 77 : x = idealismaximal_i(checknf(nf), x);
2092 77 : if (!x) { set_avma(av); return gen_0; }
2093 49 : return gc_GEN(av, x);
2094 : }
2095 :
2096 : /* I^(-1) = { x \in K, Tr(x D^(-1) I) \in Z }, D different of K/Q
2097 : *
2098 : * nf[5][6] = pp( D^(-1) ) = pp( HNF( T^(-1) ) ), T = (Tr(wi wj))
2099 : * nf[5][7] = same in 2-elt form.
2100 : * Assume I integral. Return the integral ideal (I\cap Z) I^(-1) */
2101 : GEN
2102 211905 : idealHNF_inv_Z(GEN nf, GEN I)
2103 : {
2104 211905 : GEN J, dual, IZ = gcoeff(I,1,1); /* I \cap Z */
2105 211905 : if (isint1(IZ)) return matid(lg(I)-1);
2106 197442 : J = idealHNF_mul(nf,I, gmael(nf,5,7));
2107 : /* I in HNF, hence easily inverted; multiply by IZ to get integer coeffs
2108 : * missing content cancels while solving the linear equation */
2109 197443 : dual = shallowtrans( hnf_divscale(J, gmael(nf,5,6), IZ) );
2110 197442 : return ZM_hnfmodid(dual, IZ);
2111 : }
2112 : /* I HNF with rational coefficients (denominator d). */
2113 : GEN
2114 69848 : idealHNF_inv(GEN nf, GEN I)
2115 : {
2116 69848 : GEN J, IQ = gcoeff(I,1,1); /* I \cap Q; d IQ = dI \cap Z */
2117 69848 : J = idealHNF_inv_Z(nf, Q_remove_denom(I, NULL)); /* = (dI)^(-1) * (d IQ) */
2118 69848 : return equali1(IQ)? J: RgM_Rg_div(J, IQ);
2119 : }
2120 :
2121 : /* return p * P^(-1) [integral] */
2122 : GEN
2123 38485 : pr_inv_p(GEN pr)
2124 : {
2125 38485 : if (pr_is_inert(pr)) return matid(pr_get_f(pr));
2126 37813 : return ZM_hnfmodid(pr_get_tau(pr), pr_get_p(pr));
2127 : }
2128 : GEN
2129 17836 : pr_inv(GEN pr)
2130 : {
2131 17836 : GEN p = pr_get_p(pr);
2132 17836 : if (pr_is_inert(pr)) return scalarmat(ginv(p), pr_get_f(pr));
2133 17563 : return RgM_Rg_div(ZM_hnfmodid(pr_get_tau(pr),p), p);
2134 : }
2135 :
2136 : GEN
2137 115785 : idealinv(GEN nf, GEN x)
2138 : {
2139 : GEN res, ax;
2140 : pari_sp av;
2141 115785 : long tx = idealtyp(&x,&ax), N;
2142 :
2143 115785 : res = ax? cgetg(3,t_VEC): NULL;
2144 115785 : nf = checknf(nf); av = avma;
2145 115785 : N = nf_get_degree(nf);
2146 115785 : switch (tx)
2147 : {
2148 62948 : case id_MAT:
2149 62948 : if (lg(x)-1 != N) pari_err_DIM("idealinv");
2150 62948 : x = idealHNF_inv(nf,x); break;
2151 35833 : case id_PRINCIPAL:
2152 35833 : x = nf_to_scalar_or_basis(nf, x);
2153 35833 : if (typ(x) != t_COL)
2154 35784 : x = idealhnf_principal(nf,ginv(x));
2155 : else
2156 : { /* nfinv + idealhnf where we already know (x) \cap Z */
2157 : GEN c, d;
2158 49 : x = Q_remove_denom(x, &c);
2159 49 : x = zk_inv(nf, x);
2160 49 : x = Q_remove_denom(x, &d); /* true inverse is c/d * x */
2161 49 : if (!d) /* x and x^(-1) integral => x a unit */
2162 14 : x = c? scalarmat(c, N): matid(N);
2163 : else
2164 : {
2165 35 : c = c? gdiv(c,d): ginv(d);
2166 35 : x = zk_multable(nf, x);
2167 35 : x = ZM_Q_mul(ZM_hnfmodid(x,d), c);
2168 : }
2169 : }
2170 35833 : break;
2171 17004 : case id_PRIME:
2172 17004 : x = pr_inv(x); break;
2173 : }
2174 115785 : x = gc_upto(av,x); if (!ax) return x;
2175 20381 : gel(res,1) = x;
2176 20381 : gel(res,2) = ext_inv(nf, ax); return res;
2177 : }
2178 :
2179 : /* write x = A/B, A,B coprime integral ideals */
2180 : GEN
2181 389730 : idealnumden(GEN nf, GEN x)
2182 : {
2183 389730 : pari_sp av = avma;
2184 : GEN x0, c, d, A, B, J;
2185 389730 : long tx = idealtyp(&x, NULL);
2186 389729 : nf = checknf(nf);
2187 389735 : switch (tx)
2188 : {
2189 7 : case id_PRIME:
2190 7 : retmkvec2(idealhnf(nf, x), gen_1);
2191 148175 : case id_PRINCIPAL:
2192 : {
2193 : GEN xZ, mx;
2194 148175 : x = nf_to_scalar_or_basis(nf, x);
2195 148175 : switch(typ(x))
2196 : {
2197 88060 : case t_INT: return gc_GEN(av, mkvec2(absi_shallow(x),gen_1));
2198 2639 : case t_FRAC:return gc_GEN(av, mkvec2(absi_shallow(gel(x,1)), gel(x,2)));
2199 : }
2200 : /* t_COL */
2201 57476 : x = Q_remove_denom(x, &d);
2202 57476 : if (!d) return gc_GEN(av, mkvec2(idealhnf_shallow(nf, x), gen_1));
2203 105 : mx = zk_multable(nf, x);
2204 105 : xZ = zkmultable_capZ(mx);
2205 105 : x = ZM_hnfmodid(mx, xZ); /* principal ideal (x) */
2206 105 : x0 = mkvec2(xZ, mx); /* same, for fast multiplication */
2207 105 : break;
2208 : }
2209 241553 : default: /* id_MAT */
2210 : {
2211 241553 : long n = lg(x)-1;
2212 241553 : if (n == 0) return mkvec2(gen_0, gen_1);
2213 241553 : if (n != nf_get_degree(nf)) pari_err_DIM("idealnumden");
2214 241552 : x0 = x = Q_remove_denom(x, &d);
2215 241551 : if (!d) return gc_GEN(av, mkvec2(x, gen_1));
2216 21 : break;
2217 : }
2218 : }
2219 126 : J = hnfmodid(x, d); /* = d/B */
2220 126 : c = gcoeff(J,1,1); /* (d/B) \cap Z, divides d */
2221 126 : B = idealHNF_inv_Z(nf, J); /* (d/B \cap Z) B/d */
2222 126 : if (!equalii(c,d)) B = ZM_Z_mul(B, diviiexact(d,c)); /* = B ! */
2223 126 : A = idealHNF_mul(nf, B, x0); /* d * (original x) * B = d A */
2224 126 : A = ZM_Z_divexact(A, d); /* = A ! */
2225 126 : return gc_GEN(av, mkvec2(A, B));
2226 : }
2227 :
2228 : /* Return x, integral in 2-elt form, such that pr^n = c * x. Assume n != 0.
2229 : * nf = true nf */
2230 : static GEN
2231 1305574 : idealpowprime(GEN nf, GEN pr, GEN n, GEN *pc)
2232 : {
2233 1305574 : GEN p = pr_get_p(pr), q, gen;
2234 :
2235 1305554 : *pc = NULL;
2236 1305554 : if (is_pm1(n)) /* n = 1 special cased for efficiency */
2237 : {
2238 634935 : q = p;
2239 634935 : if (typ(pr_get_tau(pr)) == t_INT) /* inert */
2240 : {
2241 0 : *pc = (signe(n) >= 0)? p: ginv(p);
2242 0 : return mkvec2(gen_1,gen_0);
2243 : }
2244 634920 : if (signe(n) >= 0) gen = pr_get_gen(pr);
2245 : else
2246 : {
2247 170800 : gen = pr_get_tau(pr); /* possibly t_MAT */
2248 170814 : *pc = ginv(p);
2249 : }
2250 : }
2251 670698 : else if (equalis(n,2)) return idealsqrprime(nf, pr, pc);
2252 : else
2253 : {
2254 363984 : long e = pr_get_e(pr), f = pr_get_f(pr);
2255 363998 : GEN r, m = truedvmdis(n, e, &r);
2256 363980 : if (e * f == nf_get_degree(nf))
2257 : { /* pr^e = (p) */
2258 76624 : if (signe(m)) *pc = powii(p,m);
2259 76628 : if (!signe(r)) return mkvec2(gen_1,gen_0);
2260 36507 : q = p;
2261 36507 : gen = nfpow(nf, pr_get_gen(pr), r);
2262 : }
2263 : else
2264 : {
2265 287359 : m = absi_shallow(m);
2266 287362 : if (signe(r)) m = addiu(m,1);
2267 287362 : q = powii(p,m); /* m = ceil(|n|/e) */
2268 287370 : if (signe(n) >= 0) gen = nfpow(nf, pr_get_gen(pr), n);
2269 : else
2270 : {
2271 43151 : gen = pr_get_tau(pr);
2272 43151 : if (typ(gen) == t_MAT) gen = gel(gen,1);
2273 43151 : n = negi(n);
2274 43151 : gen = ZC_Z_divexact(nfpow(nf, gen, n), powii(p, subii(n,m)));
2275 43151 : *pc = ginv(q);
2276 : }
2277 : }
2278 323883 : gen = FpC_red(gen, q);
2279 : }
2280 958779 : return mkvec2(q, gen);
2281 : }
2282 :
2283 : /* True nf. x * pr^n. Assume x in HNF or scalar (possibly nonintegral) */
2284 : GEN
2285 825031 : idealmulpowprime(GEN nf, GEN x, GEN pr, GEN n)
2286 : {
2287 : GEN c, cx, y;
2288 825031 : long N = nf_get_degree(nf);
2289 :
2290 825034 : if (!signe(n)) return typ(x) == t_MAT? x: scalarmat_shallow(x, N);
2291 :
2292 : /* inert, special cased for efficiency */
2293 825027 : if (pr_is_inert(pr))
2294 : {
2295 76988 : GEN q = powii(pr_get_p(pr), n);
2296 74924 : return typ(x) == t_MAT? RgM_Rg_mul(x,q)
2297 151914 : : scalarmat_shallow(gmul(Q_abs(x),q), N);
2298 : }
2299 :
2300 748051 : y = idealpowprime(nf, pr, n, &c);
2301 748015 : if (typ(x) == t_MAT)
2302 745672 : { x = Q_primitive_part(x, &cx); if (is_pm1(gcoeff(x,1,1))) x = NULL; }
2303 : else
2304 2343 : { cx = x; x = NULL; }
2305 747835 : cx = mul_content(c,cx);
2306 747831 : if (x)
2307 571333 : x = idealHNF_mul_two(nf,x,y);
2308 : else
2309 176498 : x = idealhnf_two(nf,y);
2310 748258 : if (cx) x = ZM_Q_mul(x,cx);
2311 747910 : return x;
2312 : }
2313 : GEN
2314 13007 : idealdivpowprime(GEN nf, GEN x, GEN pr, GEN n)
2315 : {
2316 13007 : return idealmulpowprime(nf,x,pr, negi(n));
2317 : }
2318 :
2319 : /* nf = true nf */
2320 : static GEN
2321 941584 : idealpow_aux(GEN nf, GEN x, long tx, GEN n)
2322 : {
2323 941584 : GEN T = nf_get_pol(nf), m, cx, n1, a, alpha;
2324 941584 : long N = degpol(T), s = signe(n);
2325 941586 : if (!s) return matid(N);
2326 926479 : switch(tx)
2327 : {
2328 75528 : case id_PRINCIPAL:
2329 75528 : return idealhnf_principal(nf, nfpow(nf,x,n));
2330 656018 : case id_PRIME:
2331 656018 : if (pr_is_inert(x)) return scalarmat(powii(gel(x,1), n), N);
2332 557496 : x = idealpowprime(nf, x, n, &cx);
2333 557480 : x = idealhnf_two(nf,x);
2334 557501 : return cx? ZM_Q_mul(x, cx): x;
2335 194933 : default:
2336 194933 : if (is_pm1(n)) return (s < 0)? idealinv(nf, x): gcopy(x);
2337 69109 : n1 = (s < 0)? negi(n): n;
2338 :
2339 69109 : x = Q_primitive_part(x, &cx);
2340 69109 : a = mat_ideal_two_elt(nf,x); alpha = gel(a,2); a = gel(a,1);
2341 69109 : alpha = nfpow(nf,alpha,n1);
2342 69109 : m = zk_scalar_or_multable(nf, alpha);
2343 69109 : if (typ(m) == t_INT) {
2344 553 : x = gcdii(powii(a,n1), m);
2345 553 : if (s<0) x = ginv(x);
2346 553 : if (cx) x = gmul(x, powgi(cx,n));
2347 553 : x = scalarmat(x, N);
2348 : }
2349 : else
2350 : { /* could use gcdii(powii(a,n1), zkmultable_capZ(m)), but costly */
2351 68556 : x = ZM_hnfmodid(m, powii(a,n1));
2352 68556 : if (cx) cx = powgi(cx,n);
2353 68556 : if (s<0) {
2354 7 : GEN xZ = gcoeff(x,1,1);
2355 7 : cx = cx ? gdiv(cx, xZ): ginv(xZ);
2356 7 : x = idealHNF_inv_Z(nf,x);
2357 : }
2358 68556 : if (cx) x = ZM_Q_mul(x, cx);
2359 : }
2360 69109 : return x;
2361 : }
2362 : }
2363 :
2364 : /* raise the ideal x to the power n (in Z) */
2365 : GEN
2366 941584 : idealpow(GEN nf, GEN x, GEN n)
2367 : {
2368 : pari_sp av;
2369 : long tx;
2370 : GEN res, ax;
2371 :
2372 941584 : if (typ(n) != t_INT) pari_err_TYPE("idealpow",n);
2373 941584 : tx = idealtyp(&x,&ax);
2374 941585 : res = ax? cgetg(3,t_VEC): NULL;
2375 941585 : av = avma;
2376 941585 : x = gc_upto(av, idealpow_aux(checknf(nf), x, tx, n));
2377 941590 : if (!ax) return x;
2378 0 : gel(res,1) = x;
2379 0 : gel(res,2) = ext_pow(nf, ax, n);
2380 0 : return res;
2381 : }
2382 :
2383 : /* Return ideal^e in number field nf. e is a C integer. */
2384 : GEN
2385 313276 : idealpows(GEN nf, GEN ideal, long e)
2386 : {
2387 313276 : long court[] = {evaltyp(t_INT) | _evallg(3),0,0};
2388 313276 : affsi(e,court); return idealpow(nf,ideal,court);
2389 : }
2390 :
2391 : static GEN
2392 28591 : _idealmulred(GEN nf, GEN x, GEN y)
2393 28591 : { return idealred(nf,idealmul(nf,x,y)); }
2394 : static GEN
2395 35663 : _idealsqrred(GEN nf, GEN x)
2396 35663 : { return idealred(nf,idealsqr(nf,x)); }
2397 : static GEN
2398 11370 : _mul(void *data, GEN x, GEN y) { return _idealmulred((GEN)data,x,y); }
2399 : static GEN
2400 35663 : _sqr(void *data, GEN x) { return _idealsqrred((GEN)data, x); }
2401 :
2402 : /* compute x^n (x ideal, n integer), reducing along the way */
2403 : GEN
2404 80206 : idealpowred(GEN nf, GEN x, GEN n)
2405 : {
2406 80206 : pari_sp av = avma, av2;
2407 : long s;
2408 : GEN y;
2409 :
2410 80206 : if (typ(n) != t_INT) pari_err_TYPE("idealpowred",n);
2411 80206 : s = signe(n); if (s == 0) return idealpow(nf,x,n);
2412 80206 : y = gen_pow_i(x, n, (void*)nf, &_sqr, &_mul);
2413 80206 : av2 = avma;
2414 80206 : if (s < 0) y = idealinv(nf,y);
2415 80206 : if (s < 0 || is_pm1(n)) y = idealred(nf,y);
2416 80206 : return avma == av2? gc_GEN(av,y): gc_upto(av,y);
2417 : }
2418 :
2419 : GEN
2420 17221 : idealmulred(GEN nf, GEN x, GEN y)
2421 : {
2422 17221 : pari_sp av = avma;
2423 17221 : return gc_upto(av, _idealmulred(nf,x,y));
2424 : }
2425 :
2426 : long
2427 91 : isideal(GEN nf,GEN x)
2428 : {
2429 91 : long N, i, j, lx, tx = typ(x);
2430 : pari_sp av;
2431 : GEN T, xZ;
2432 :
2433 91 : nf = checknf(nf); T = nf_get_pol(nf); lx = lg(x);
2434 91 : if (tx==t_VEC && lx==3) { x = gel(x,1); tx = typ(x); lx = lg(x); }
2435 91 : switch(tx)
2436 : {
2437 14 : case t_INT: case t_FRAC: return 1;
2438 7 : case t_POL: return varn(x) == varn(T);
2439 7 : case t_POLMOD: return RgX_equal_var(T, gel(x,1));
2440 14 : case t_VEC: return get_prid(x)? 1 : 0;
2441 42 : case t_MAT: break;
2442 7 : default: return 0;
2443 : }
2444 42 : N = degpol(T);
2445 42 : if (lx-1 != N) return (lx == 1);
2446 28 : if (nbrows(x) != N) return 0;
2447 :
2448 28 : av = avma; x = Q_primpart(x);
2449 28 : if (!ZM_ishnf(x)) return 0;
2450 14 : xZ = gcoeff(x,1,1);
2451 21 : for (j=2; j<=N; j++)
2452 14 : if (!dvdii(xZ, gcoeff(x,j,j))) return gc_long(av,0);
2453 14 : for (i=2; i<=N; i++)
2454 14 : for (j=2; j<=N; j++)
2455 7 : if (! hnf_invimage(x, zk_ei_mul(nf,gel(x,i),j))) return gc_long(av,0);
2456 7 : return gc_long(av,1);
2457 : }
2458 :
2459 : GEN
2460 39716 : idealdiv(GEN nf, GEN x, GEN y)
2461 : {
2462 39716 : pari_sp av = avma;
2463 39716 : return gc_upto(av, idealmul(nf, x, idealinv(nf,y)));
2464 : }
2465 :
2466 : /* This routine computes the quotient x/y of two ideals in the number field nf.
2467 : * It assumes that the quotient is an integral ideal. The idea is to find an
2468 : * ideal z dividing y such that gcd(Nx/Nz, Nz) = 1. Then
2469 : *
2470 : * x + (Nx/Nz) x
2471 : * ----------- = ---
2472 : * y + (Ny/Nz) y
2473 : *
2474 : * Proof: we can assume x and y are integral. Let p be any prime ideal
2475 : *
2476 : * If p | Nz, then it divides neither Nx/Nz nor Ny/Nz (since Nx/Nz is the
2477 : * product of the integers N(x/y) and N(y/z)). Both the numerator and the
2478 : * denominator on the left will be coprime to p. So will x/y, since x/y is
2479 : * assumed integral and its norm N(x/y) is coprime to p.
2480 : *
2481 : * If instead p does not divide Nz, then v_p (Nx/Nz) = v_p (Nx) >= v_p(x).
2482 : * Hence v_p (x + Nx/Nz) = v_p(x). Likewise for the denominators. QED.
2483 : *
2484 : * Peter Montgomery. July, 1994. */
2485 : static void
2486 7 : err_divexact(GEN x, GEN y)
2487 7 : { pari_err_DOMAIN("idealdivexact","denominator(x/y)", "!=",
2488 0 : gen_1,mkvec2(x,y)); }
2489 : GEN
2490 5221 : idealdivexact(GEN nf, GEN x0, GEN y0)
2491 : {
2492 5221 : pari_sp av = avma;
2493 : GEN x, y, xZ, yZ, Nx, Ny, Nz, cy, q, r;
2494 :
2495 5221 : nf = checknf(nf);
2496 5221 : x = idealhnf_shallow(nf, x0);
2497 5221 : y = idealhnf_shallow(nf, y0);
2498 5221 : if (lg(y) == 1) pari_err_INV("idealdivexact", y0);
2499 5214 : if (lg(x) == 1) retgc_const(av, cgetg(1, t_MAT)); /* numerator is zero */
2500 5214 : y = Q_primitive_part(y, &cy);
2501 5214 : if (cy) x = RgM_Rg_div(x,cy);
2502 5214 : xZ = gcoeff(x,1,1); if (typ(xZ) != t_INT) err_divexact(x,y);
2503 5207 : yZ = gcoeff(y,1,1); if (isint1(yZ)) return gc_GEN(av, x);
2504 2828 : Nx = idealnorm(nf,x);
2505 2828 : Ny = idealnorm(nf,y);
2506 2828 : if (typ(Nx) != t_INT) err_divexact(x,y);
2507 2828 : q = dvmdii(Nx,Ny, &r);
2508 2828 : if (signe(r)) err_divexact(x,y);
2509 2828 : if (is_pm1(q)) { set_avma(av); return matid(nf_get_degree(nf)); }
2510 : /* Find a norm Nz | Ny such that gcd(Nx/Nz, Nz) = 1 */
2511 581 : for (Nz = Ny;;) /* q = Nx/Nz */
2512 470 : {
2513 1051 : GEN p1 = gcdii(Nz, q);
2514 1051 : if (is_pm1(p1)) break;
2515 470 : Nz = diviiexact(Nz,p1);
2516 470 : q = mulii(q,p1);
2517 : }
2518 581 : xZ = gcoeff(x,1,1); q = gcdii(q, xZ);
2519 581 : if (!equalii(xZ,q))
2520 : { /* Replace x/y by x+(Nx/Nz) / y+(Ny/Nz) */
2521 433 : x = ZM_hnfmodid(x, q);
2522 : /* y reduced to unit ideal ? */
2523 433 : if (Nz == Ny) return gc_upto(av, x);
2524 :
2525 125 : yZ = gcoeff(y,1,1); q = gcdii(diviiexact(Ny,Nz), yZ);
2526 125 : y = ZM_hnfmodid(y, q);
2527 : }
2528 273 : yZ = gcoeff(y,1,1);
2529 273 : y = idealHNF_mul(nf,x, idealHNF_inv_Z(nf,y));
2530 273 : return gc_upto(av, ZM_Z_divexact(y, yZ));
2531 : }
2532 :
2533 : GEN
2534 21 : idealintersect(GEN nf, GEN x, GEN y)
2535 : {
2536 21 : pari_sp av = avma;
2537 : long lz, lx, i;
2538 : GEN z, dx, dy, xZ, yZ;;
2539 :
2540 21 : nf = checknf(nf);
2541 21 : x = idealhnf_shallow(nf,x);
2542 21 : y = idealhnf_shallow(nf,y);
2543 21 : if (lg(x) == 1 || lg(y) == 1) retgc_const(av, cgetg(1, t_MAT));
2544 14 : x = Q_remove_denom(x, &dx);
2545 14 : y = Q_remove_denom(y, &dy);
2546 14 : if (dx) y = ZM_Z_mul(y, dx);
2547 14 : if (dy) x = ZM_Z_mul(x, dy);
2548 14 : xZ = gcoeff(x,1,1);
2549 14 : yZ = gcoeff(y,1,1);
2550 14 : dx = mul_denom(dx,dy);
2551 14 : z = ZM_lll(shallowconcat(x,y), 0.99, LLL_KER); lz = lg(z);
2552 14 : lx = lg(x);
2553 63 : for (i=1; i<lz; i++) setlg(z[i], lx);
2554 14 : z = ZM_hnfmodid(ZM_mul(x,z), lcmii(xZ, yZ));
2555 14 : if (dx) z = RgM_Rg_div(z,dx);
2556 14 : return gc_upto(av,z);
2557 : }
2558 :
2559 : /*******************************************************************/
2560 : /* */
2561 : /* T2-IDEAL REDUCTION */
2562 : /* */
2563 : /*******************************************************************/
2564 :
2565 : static GEN
2566 21 : chk_vdir(GEN nf, GEN vdir)
2567 : {
2568 21 : long i, l = lg(vdir);
2569 : GEN v;
2570 21 : if (l != lg(nf_get_roots(nf))) pari_err_DIM("idealred");
2571 14 : switch(typ(vdir))
2572 : {
2573 0 : case t_VECSMALL: return vdir;
2574 14 : case t_VEC: break;
2575 0 : default: pari_err_TYPE("idealred",vdir);
2576 : }
2577 14 : v = cgetg(l, t_VECSMALL);
2578 56 : for (i = 1; i < l; i++) v[i] = itos(gceil(gel(vdir,i)));
2579 14 : return v;
2580 : }
2581 :
2582 : static void
2583 12639 : twistG(GEN G, long r1, long i, long v)
2584 : {
2585 12639 : long j, lG = lg(G);
2586 12639 : if (i <= r1) {
2587 37275 : for (j=1; j<lG; j++) gcoeff(G,i,j) = gmul2n(gcoeff(G,i,j), v);
2588 : } else {
2589 578 : long k = (i<<1) - r1;
2590 4402 : for (j=1; j<lG; j++)
2591 : {
2592 3824 : gcoeff(G,k-1,j) = gmul2n(gcoeff(G,k-1,j), v);
2593 3824 : gcoeff(G,k ,j) = gmul2n(gcoeff(G,k ,j), v);
2594 : }
2595 : }
2596 12639 : }
2597 :
2598 : GEN
2599 138582 : nf_get_Gtwist(GEN nf, GEN vdir)
2600 : {
2601 : long i, l, v, r1;
2602 : GEN G;
2603 :
2604 138582 : if (!vdir) return nf_get_roundG(nf);
2605 21 : if (typ(vdir) == t_MAT)
2606 : {
2607 0 : long N = nf_get_degree(nf);
2608 0 : if (lg(vdir) != N+1 || lgcols(vdir) != N+1) pari_err_DIM("idealred");
2609 0 : return vdir;
2610 : }
2611 21 : vdir = chk_vdir(nf, vdir);
2612 14 : G = RgM_shallowcopy(nf_get_G(nf));
2613 14 : r1 = nf_get_r1(nf);
2614 14 : l = lg(vdir);
2615 56 : for (i=1; i<l; i++)
2616 : {
2617 42 : v = vdir[i]; if (!v) continue;
2618 42 : twistG(G, r1, i, v);
2619 : }
2620 14 : return RM_round_maxrank(G);
2621 : }
2622 : GEN
2623 12597 : nf_get_Gtwist1(GEN nf, long i)
2624 : {
2625 12597 : GEN G = RgM_shallowcopy( nf_get_G(nf) );
2626 12597 : long r1 = nf_get_r1(nf);
2627 12597 : twistG(G, r1, i, 10);
2628 12597 : return RM_round_maxrank(G);
2629 : }
2630 :
2631 : GEN
2632 98091 : RM_round_maxrank(GEN G0)
2633 : {
2634 98091 : long e, r = lg(G0)-1;
2635 98091 : pari_sp av = avma;
2636 98091 : for (e = 4; ; e <<= 1, set_avma(av))
2637 0 : {
2638 98091 : GEN G = gmul2n(G0, e), H = ground(G);
2639 98090 : if (ZM_rank(H) == r) return H; /* maximal rank ? */
2640 : }
2641 : }
2642 :
2643 : GEN
2644 138575 : idealred0(GEN nf, GEN I, GEN vdir)
2645 : {
2646 138575 : pari_sp av = avma;
2647 138575 : GEN G, aI, IZ, J, y, my, dyi, yi, c1 = NULL;
2648 : long N;
2649 :
2650 138575 : nf = checknf(nf);
2651 138575 : N = nf_get_degree(nf);
2652 : /* put first for sanity checks, unused when I obviously principal */
2653 138575 : G = nf_get_Gtwist(nf, vdir);
2654 138568 : switch (idealtyp(&I,&aI))
2655 : {
2656 37209 : case id_PRIME:
2657 37209 : if (pr_is_inert(I)) {
2658 655 : if (!aI) { set_avma(av); return matid(N); }
2659 655 : c1 = gel(I,1); I = matid(N);
2660 655 : goto END;
2661 : }
2662 36554 : IZ = pr_get_p(I);
2663 36554 : J = pr_inv_p(I);
2664 36554 : I = idealhnf_two(nf,I);
2665 36554 : break;
2666 101331 : case id_MAT:
2667 101331 : if (lg(I)-1 != N) pari_err_DIM("idealred");
2668 101324 : I = Q_primitive_part(I, &c1);
2669 101324 : IZ = gcoeff(I,1,1);
2670 101324 : if (is_pm1(IZ))
2671 : {
2672 8838 : if (!aI) { set_avma(av); return matid(N); }
2673 8754 : goto END;
2674 : }
2675 92486 : J = idealHNF_inv_Z(nf, I);
2676 92486 : break;
2677 21 : default: /* id_PRINCIPAL, silly case */
2678 21 : if (gequal0(I)) I = cgetg(1,t_MAT); else { c1 = I; I = matid(N); }
2679 21 : if (!aI) return I;
2680 14 : goto END;
2681 : }
2682 : /* now I integral, HNF; and J = (I\cap Z) I^(-1), integral */
2683 129040 : y = idealpseudomin(J, G); /* small elt in (I\cap Z)I^(-1), integral */
2684 129039 : if (ZV_isscalar(y))
2685 : { /* already reduced */
2686 71008 : if (!aI) return gc_GEN(av, I);
2687 67822 : goto END;
2688 : }
2689 :
2690 58031 : my = zk_multable(nf, y);
2691 58032 : I = ZM_Z_divexact(ZM_mul(my, I), IZ); /* y I / (I\cap Z), integral */
2692 58031 : c1 = mul_content(c1, IZ);
2693 58031 : if (equali1(c1)) c1 = NULL; /* can be simplified with IZ */
2694 58031 : yi = ZM_gauss(my, col_ei(N,1)); /* y^-1 */
2695 58032 : dyi = Q_denom(yi); /* generates (y) \cap Z */
2696 58032 : I = hnfmodid(I, dyi);
2697 58032 : if (!aI) return gc_upto(av, I);
2698 56080 : if (typ(aI) == t_MAT)
2699 : {
2700 39308 : GEN nyi = Q_muli_to_int(yi, dyi);
2701 39308 : if (gexpo(nyi) >= gexpo(y))
2702 20843 : aI = famat_div(aI, y); /* yi "larger" than y, keep the latter */
2703 : else
2704 : { /* use yi */
2705 18465 : aI = famat_mul(aI, nyi);
2706 18465 : c1 = div_content(c1, dyi);
2707 : }
2708 39308 : if (c1) { aI = famat_mul(aI, Q_to_famat(c1)); c1 = NULL; }
2709 : }
2710 : else
2711 16772 : c1 = c1? RgC_Rg_mul(yi, c1): yi;
2712 133325 : END:
2713 133325 : if (c1) aI = ext_mul(nf, aI,c1);
2714 133325 : return gc_GEN(av, mkvec2(I, aI));
2715 : }
2716 :
2717 : /* I integral ZM (not HNF), G ZM, rounded Cholesky form of a weighted
2718 : * T2 matrix. Reduce I wrt G */
2719 : GEN
2720 1344417 : idealpseudored(GEN I, GEN G)
2721 1344417 : { return ZM_mul(I, ZM_lll(ZM_mul(G, I), 0.99, LLL_IM)); }
2722 :
2723 : /* Same I, G; m in I with T2(m) small */
2724 : GEN
2725 142443 : idealpseudomin(GEN I, GEN G)
2726 : {
2727 142443 : GEN u = ZM_lll(ZM_mul(G, I), 0.99, LLL_IM);
2728 142442 : return ZM_ZC_mul(I, gel(u,1));
2729 : }
2730 : /* Same I,G; irrational m in I with T2(m) small */
2731 : GEN
2732 0 : idealpseudomin_nonscalar(GEN I, GEN G)
2733 : {
2734 0 : GEN u = ZM_lll(ZM_mul(G, I), 0.99, LLL_IM);
2735 0 : GEN m = ZM_ZC_mul(I, gel(u,1));
2736 0 : if (ZV_isscalar(m) && lg(u) > 2) m = ZM_ZC_mul(I, gel(u,2));
2737 0 : return m;
2738 : }
2739 : /* Same I,G; t_VEC of irrational m in I with T2(m) small */
2740 : GEN
2741 1252510 : idealpseudominvec(GEN I, GEN G)
2742 : {
2743 1252510 : long i, j, k, n = lg(I)-1;
2744 1252510 : GEN x, L, b = idealpseudored(I, G);
2745 1252509 : L = cgetg(1 + (n*(n+1))/2, t_VEC);
2746 4413538 : for (i = k = 1; i <= n; i++)
2747 : {
2748 3161026 : x = gel(b,i);
2749 3161026 : if (!ZV_isscalar(x)) gel(L,k++) = x;
2750 : }
2751 3161033 : for (i = 2; i <= n; i++)
2752 : {
2753 1908520 : long J = minss(i, 4);
2754 4729898 : for (j = 1; j < J; j++)
2755 : {
2756 2821377 : x = ZC_add(gel(b,i),gel(b,j));
2757 2821378 : if (!ZV_isscalar(x)) gel(L,k++) = x;
2758 : }
2759 : }
2760 1252513 : setlg(L,k); return L;
2761 : }
2762 :
2763 : GEN
2764 13396 : idealred_elt(GEN nf, GEN I)
2765 : {
2766 13396 : pari_sp av = avma;
2767 13396 : GEN u = idealpseudomin(I, nf_get_roundG(nf));
2768 13396 : return gc_upto(av, u);
2769 : }
2770 :
2771 : GEN
2772 7 : idealmin(GEN nf, GEN x, GEN vdir)
2773 : {
2774 7 : pari_sp av = avma;
2775 : GEN y, dx;
2776 7 : nf = checknf(nf);
2777 7 : switch( idealtyp(&x, NULL) )
2778 : {
2779 0 : case id_PRINCIPAL: return gcopy(x);
2780 0 : case id_PRIME: x = pr_hnf(nf,x); break;
2781 7 : case id_MAT: if (lg(x) == 1) return gen_0;
2782 : }
2783 7 : x = Q_remove_denom(x, &dx);
2784 7 : y = idealpseudomin(x, nf_get_Gtwist(nf,vdir));
2785 7 : if (dx) y = RgC_Rg_div(y, dx);
2786 7 : return gc_upto(av, y);
2787 : }
2788 :
2789 : /*******************************************************************/
2790 : /* */
2791 : /* APPROXIMATION THEOREM */
2792 : /* */
2793 : /*******************************************************************/
2794 : /* a = ppi(a,b) ppo(a,b), where ppi regroups primes common to a and b
2795 : * and ppo(a,b) = Z_ppo(a,b) */
2796 : /* return gcd(a,b),ppi(a,b),ppo(a,b) */
2797 : GEN
2798 986062 : Z_ppio(GEN a, GEN b)
2799 : {
2800 986062 : GEN x, y, d = gcdii(a,b);
2801 986062 : if (is_pm1(d)) return mkvec3(gen_1, gen_1, a);
2802 757428 : x = d; y = diviiexact(a,d);
2803 : for(;;)
2804 131075 : {
2805 888503 : GEN g = gcdii(x,y);
2806 888503 : if (is_pm1(g)) return mkvec3(d, x, y);
2807 131075 : x = mulii(x,g); y = diviiexact(y,g);
2808 : }
2809 : }
2810 : /* a = ppg(a,b)pple(a,b), where ppg regroups primes such that v(a) > v(b)
2811 : * and pple all others */
2812 : /* return gcd(a,b),ppg(a,b),pple(a,b) */
2813 : GEN
2814 0 : Z_ppgle(GEN a, GEN b)
2815 : {
2816 0 : GEN x, y, g, d = gcdii(a,b);
2817 0 : if (equalii(a, d)) return mkvec3(a, gen_1, a);
2818 0 : x = diviiexact(a,d); y = d;
2819 : for(;;)
2820 : {
2821 0 : g = gcdii(x,y);
2822 0 : if (is_pm1(g)) return mkvec3(d, x, y);
2823 0 : x = mulii(x,g); y = diviiexact(y,g);
2824 : }
2825 : }
2826 : static void
2827 0 : Z_dcba_rec(GEN L, GEN a, GEN b)
2828 : {
2829 : GEN x, r, v, g, h, c, c0;
2830 : long n;
2831 0 : if (is_pm1(b)) {
2832 0 : if (!is_pm1(a)) vectrunc_append(L, a);
2833 0 : return;
2834 : }
2835 0 : v = Z_ppio(a,b);
2836 0 : a = gel(v,2);
2837 0 : r = gel(v,3);
2838 0 : if (!is_pm1(r)) vectrunc_append(L, r);
2839 0 : v = Z_ppgle(a,b);
2840 0 : g = gel(v,1);
2841 0 : h = gel(v,2);
2842 0 : x = c0 = gel(v,3);
2843 0 : for (n = 1; !is_pm1(h); n++)
2844 : {
2845 : GEN d, y;
2846 : long i;
2847 0 : v = Z_ppgle(h,sqri(g));
2848 0 : g = gel(v,1);
2849 0 : h = gel(v,2);
2850 0 : c = gel(v,3); if (is_pm1(c)) continue;
2851 0 : d = gcdii(c,b);
2852 0 : x = mulii(x,d);
2853 0 : y = d; for (i=1; i < n; i++) y = sqri(y);
2854 0 : Z_dcba_rec(L, diviiexact(c,y), d);
2855 : }
2856 0 : Z_dcba_rec(L,diviiexact(b,x), c0);
2857 : }
2858 : static GEN
2859 6827597 : Z_cba_rec(GEN L, GEN a, GEN b)
2860 : {
2861 : GEN g;
2862 : /* a few naive steps before switching to dcba */
2863 6827597 : if (lg(L) > 10) { Z_dcba_rec(L, a, b); return veclast(L); }
2864 6827597 : if (is_pm1(a)) return b;
2865 4061792 : g = gcdii(a,b);
2866 4061792 : if (is_pm1(g)) { vectrunc_append(L, a); return b; }
2867 3035165 : a = diviiexact(a,g);
2868 3035165 : b = diviiexact(b,g);
2869 3035165 : return Z_cba_rec(L, Z_cba_rec(L, a, g), b);
2870 : }
2871 : GEN
2872 757267 : Z_cba(GEN a, GEN b)
2873 : {
2874 757267 : GEN L = vectrunc_init(expi(a) + expi(b) + 2);
2875 757267 : GEN t = Z_cba_rec(L, a, b);
2876 757267 : if (!is_pm1(t)) vectrunc_append(L, t);
2877 757267 : return L;
2878 : }
2879 : /* P = coprime base, extend it by b; TODO: quadratic for now */
2880 : GEN
2881 49 : ZV_cba_extend(GEN P, GEN b)
2882 : {
2883 49 : long i, l = lg(P);
2884 49 : GEN w = cgetg(l+1, t_VEC);
2885 175 : for (i = 1; i < l; i++)
2886 : {
2887 126 : GEN v = Z_cba(gel(P,i), b);
2888 126 : long nv = lg(v)-1;
2889 126 : gel(w,i) = vecslice(v, 1, nv-1); /* those divide P[i] but not b */
2890 126 : b = gel(v,nv);
2891 : }
2892 49 : gel(w,l) = b; return shallowconcat1(w);
2893 : }
2894 : GEN
2895 28 : ZV_cba(GEN v)
2896 : {
2897 28 : long i, l = lg(v);
2898 : GEN P;
2899 28 : if (l <= 2) return v;
2900 14 : P = Z_cba(gel(v,1), gel(v,2));
2901 42 : for (i = 3; i < l; i++) P = ZV_cba_extend(P, gel(v,i));
2902 14 : return P;
2903 : }
2904 :
2905 : /* write x = x1 x2, x2 maximal s.t. (x2,f) = 1, return x2 */
2906 : GEN
2907 18702624 : Z_ppo(GEN x, GEN f)
2908 : {
2909 18702624 : (void)Z_pvalrem(x, f, &x);
2910 : for (;;)
2911 : {
2912 56849945 : f = gcdii(x, f); if (is_pm1(f)) break;
2913 38147728 : x = diviiexact(x, f);
2914 : }
2915 18702392 : return x;
2916 : }
2917 : /* write x = x1 x2, x2 maximal s.t. (x2,f) = 1, return x2 */
2918 : ulong
2919 70509713 : u_ppo(ulong x, ulong f)
2920 : {
2921 : for (;;)
2922 : {
2923 70509713 : f = ugcd(x, f); if (f == 1) break;
2924 16213819 : x /= f;
2925 : }
2926 54295852 : return x;
2927 : }
2928 :
2929 : /* result known to be representable as an ulong */
2930 : static ulong
2931 1636664 : lcmuu(ulong a, ulong b) { ulong d = ugcd(a,b); return (a/d) * b; }
2932 :
2933 : /* assume 0 < x < N; return u in (Z/NZ)^* such that u x = gcd(x,N) (mod N);
2934 : * set *pd = gcd(x,N) */
2935 : ulong
2936 5897717 : Fl_invgen(ulong x, ulong N, ulong *pd)
2937 : {
2938 : ulong d, d0, e, v, v1;
2939 : long s;
2940 5897717 : *pd = d = xgcduu(N, x, 0, &v, &v1, &s);
2941 5898536 : if (s > 0) v = N - v;
2942 5898536 : if (d == 1) return v;
2943 : /* vx = gcd(x,N) (mod N), v coprime to N/d but need not be coprime to N */
2944 2749474 : e = N / d;
2945 2749474 : d0 = u_ppo(d, e); /* d = d0 d1, d0 coprime to N/d, rad(d1) | N/d */
2946 2749630 : if (d0 == 1) return v;
2947 1636618 : e = lcmuu(e, d / d0);
2948 1636665 : return u_chinese_coprime(v, 1, e, d0, e*d0);
2949 : }
2950 :
2951 : /* x t_INT, f ideal. Write x = x1 x2, sqf(x1) | f, (x2,f) = 1. Return x2 */
2952 : static GEN
2953 126 : nf_coprime_part(GEN nf, GEN x, GEN listpr)
2954 : {
2955 126 : long v, j, lp = lg(listpr), N = nf_get_degree(nf);
2956 : GEN x1, x2, ex;
2957 :
2958 : #if 0 /*1) via many gcds. Expensive ! */
2959 : GEN f = idealprodprime(nf, listpr);
2960 : f = ZM_hnfmodid(f, x); /* first gcd is less expensive since x in Z */
2961 : x = scalarmat(x, N);
2962 : for (;;)
2963 : {
2964 : if (gequal1(gcoeff(f,1,1))) break;
2965 : x = idealdivexact(nf, x, f);
2966 : f = ZM_hnfmodid(shallowconcat(f,x), gcoeff(x,1,1)); /* gcd(f,x) */
2967 : }
2968 : x2 = x;
2969 : #else /*2) from prime decomposition */
2970 126 : x1 = NULL;
2971 350 : for (j=1; j<lp; j++)
2972 : {
2973 224 : GEN pr = gel(listpr,j);
2974 224 : v = Z_pval(x, pr_get_p(pr)); if (!v) continue;
2975 :
2976 126 : ex = muluu(v, pr_get_e(pr)); /* = v_pr(x) > 0 */
2977 126 : x1 = x1? idealmulpowprime(nf, x1, pr, ex)
2978 126 : : idealpow(nf, pr, ex);
2979 : }
2980 126 : x = scalarmat(x, N);
2981 126 : x2 = x1? idealdivexact(nf, x, x1): x;
2982 : #endif
2983 126 : return x2;
2984 : }
2985 :
2986 : /* L0 in K^*, assume (L0,f) = 1. Return L integral, L0 = L mod f */
2987 : GEN
2988 10920 : make_integral(GEN nf, GEN L0, GEN f, GEN listpr)
2989 : {
2990 : GEN fZ, t, L, D2, d1, d2, d;
2991 :
2992 10920 : L = Q_remove_denom(L0, &d);
2993 10920 : if (!d) return L0;
2994 :
2995 : /* L0 = L / d, L integral */
2996 518 : fZ = gcoeff(f,1,1);
2997 518 : if (typ(L) == t_INT) return Fp_mul(L, Fp_inv(d, fZ), fZ);
2998 : /* Kill denom part coprime to fZ */
2999 126 : d2 = Z_ppo(d, fZ);
3000 126 : t = Fp_inv(d2, fZ); if (!is_pm1(t)) L = ZC_Z_mul(L,t);
3001 126 : if (equalii(d, d2)) return L;
3002 :
3003 126 : d1 = diviiexact(d, d2);
3004 : /* L0 = (L / d1) mod f. d1 not coprime to f
3005 : * write (d1) = D1 D2, D2 minimal, (D2,f) = 1. */
3006 126 : D2 = nf_coprime_part(nf, d1, listpr);
3007 126 : t = idealaddtoone_i(nf, D2, f); /* in D2, 1 mod f */
3008 126 : L = nfmuli(nf,t,L);
3009 :
3010 : /* if (L0, f) = 1, then L in D1 ==> in D1 D2 = (d1) */
3011 126 : return Q_div_to_int(L, d1); /* exact division */
3012 : }
3013 :
3014 : /* assume L is a list of prime ideals. Return the product */
3015 : GEN
3016 666 : idealprodprime(GEN nf, GEN L)
3017 : {
3018 666 : long l = lg(L), i;
3019 : GEN z;
3020 666 : if (l == 1) return matid(nf_get_degree(nf));
3021 666 : z = pr_hnf(nf, gel(L,1));
3022 694 : for (i=2; i<l; i++) z = idealHNF_mul_two(nf,z, gel(L,i));
3023 666 : return z;
3024 : }
3025 :
3026 : /* optimize for the frequent case I = nfhnf()[2]: lots of them are 1 */
3027 : GEN
3028 1064 : idealprod(GEN nf, GEN I)
3029 : {
3030 1064 : long i, l = lg(I);
3031 : GEN z;
3032 2450 : for (i = 1; i < l; i++)
3033 2443 : if (!equali1(gel(I,i))) break;
3034 1064 : if (i == l) return gen_1;
3035 1057 : z = gel(I,i);
3036 1855 : for (i++; i<l; i++) z = idealmul(nf, z, gel(I,i));
3037 1057 : return z;
3038 : }
3039 :
3040 : /* v_pr(idealprod(nf,I)) */
3041 : long
3042 1917 : idealprodval(GEN nf, GEN I, GEN pr)
3043 : {
3044 1917 : long i, l = lg(I), v = 0;
3045 10371 : for (i = 1; i < l; i++)
3046 8454 : if (!equali1(gel(I,i))) v += idealval(nf, gel(I,i), pr);
3047 1917 : return v;
3048 : }
3049 :
3050 : /* assume L is a list of prime ideals. Return prod L[i]^e[i] */
3051 : GEN
3052 75596 : factorbackprime(GEN nf, GEN L, GEN e)
3053 : {
3054 75596 : long l = lg(L), i;
3055 : GEN z;
3056 :
3057 75596 : if (l == 1) return matid(nf_get_degree(nf));
3058 60973 : z = idealpow(nf, gel(L,1), gel(e,1));
3059 113897 : for (i=2; i<l; i++)
3060 52923 : if (signe(gel(e,i))) z = idealmulpowprime(nf,z, gel(L,i),gel(e,i));
3061 60974 : return z;
3062 : }
3063 :
3064 : /* F in Z, divisible exactly by pr.p. Return F-uniformizer for pr, i.e.
3065 : * a t in Z_K such that v_pr(t) = 1 and (t, F/pr) = 1 */
3066 : GEN
3067 58448 : pr_uniformizer(GEN pr, GEN F)
3068 : {
3069 58448 : GEN p = pr_get_p(pr), t = pr_get_gen(pr);
3070 58448 : if (!equalii(F, p))
3071 : {
3072 36948 : long e = pr_get_e(pr);
3073 36948 : GEN u, v, q = (e == 1)? sqri(p): p;
3074 36948 : u = mulii(q, Fp_inv(q, diviiexact(F,p))); /* 1 mod F/p, 0 mod q */
3075 36948 : v = subui(1UL, u); /* 0 mod F/p, 1 mod q */
3076 36948 : if (pr_is_inert(pr))
3077 28 : t = addii(mulii(p, v), u);
3078 : else
3079 : {
3080 36920 : t = ZC_Z_mul(t, v);
3081 36920 : gel(t,1) = addii(gel(t,1), u); /* return u + vt */
3082 : }
3083 : }
3084 58448 : return t;
3085 : }
3086 : /* L = list of prime ideals, return lcm_i (L[i] \cap \ZM) */
3087 : GEN
3088 80962 : prV_lcm_capZ(GEN L)
3089 : {
3090 80962 : long i, r = lg(L);
3091 : GEN F;
3092 80962 : if (r == 1) return gen_1;
3093 68467 : F = pr_get_p(gel(L,1));
3094 121887 : for (i = 2; i < r; i++)
3095 : {
3096 53420 : GEN pr = gel(L,i), p = pr_get_p(pr);
3097 53420 : if (!dvdii(F, p)) F = mulii(F,p);
3098 : }
3099 68467 : return F;
3100 : }
3101 : /* v vector of prid. Return underlying list of rational primes */
3102 : GEN
3103 66241 : prV_primes(GEN v)
3104 : {
3105 66241 : long i, l = lg(v);
3106 66241 : GEN w = cgetg(l,t_VEC);
3107 219775 : for (i=1; i<l; i++) gel(w,i) = pr_get_p(gel(v,i));
3108 66241 : return ZV_sort_uniq(w);
3109 : }
3110 :
3111 : /* Given a prime ideal factorization with possibly zero or negative
3112 : * exponents, gives b such that v_p(b) = v_p(x) for all prime ideals pr | x
3113 : * and v_pr(b) >= 0 for all other pr.
3114 : * For optimal performance, all [anti-]uniformizers should be precomputed,
3115 : * but no support for this yet. If nored, do not reduce result. */
3116 : static GEN
3117 54095 : idealapprfact_i(GEN nf, GEN x, int nored)
3118 : {
3119 54095 : GEN d = NULL, z, L, e, e2, F;
3120 : long i, r;
3121 54095 : int hasden = 0;
3122 :
3123 54095 : nf = checknf(nf);
3124 54096 : L = gel(x,1);
3125 54096 : e = gel(x,2);
3126 54096 : F = prV_lcm_capZ(L);
3127 54096 : z = NULL; r = lg(e);
3128 136965 : for (i = 1; i < r; i++)
3129 : {
3130 82869 : long s = signe(gel(e,i));
3131 : GEN pi, q;
3132 82869 : if (!s) continue;
3133 54080 : if (s < 0) hasden = 1;
3134 54080 : pi = pr_uniformizer(gel(L,i), F);
3135 54080 : q = nfpow(nf, pi, gel(e,i));
3136 54080 : z = z? nfmul(nf, z, q): q;
3137 : }
3138 54096 : if (!z) return gen_1;
3139 26973 : if (hasden) /* denominator */
3140 : {
3141 10075 : z = Q_remove_denom(z, &d);
3142 10075 : d = diviiexact(d, Z_ppo(d, F));
3143 : }
3144 26973 : if (nored || typ(z) != t_COL) return d? gdiv(z, d): z;
3145 10075 : e2 = cgetg(r, t_VEC);
3146 28611 : for (i = 1; i < r; i++) gel(e2,i) = addiu(gel(e,i), 1);
3147 10075 : x = factorbackprime(nf, L, e2);
3148 10075 : if (d) x = RgM_Rg_mul(x, d);
3149 10075 : z = ZC_reducemodlll(z, x);
3150 10075 : return d? RgC_Rg_div(z,d): z;
3151 : }
3152 :
3153 : GEN
3154 0 : idealapprfact(GEN nf, GEN x) {
3155 0 : pari_sp av = avma;
3156 0 : return gc_upto(av, idealapprfact_i(nf, x, 0));
3157 : }
3158 : GEN
3159 14 : idealappr(GEN nf, GEN x) {
3160 14 : pari_sp av = avma;
3161 14 : if (!is_nf_extfactor(x)) x = idealfactor(nf, x);
3162 14 : return gc_upto(av, idealapprfact_i(nf, x, 0));
3163 : }
3164 :
3165 : /* OBSOLETE */
3166 : GEN
3167 14 : idealappr0(GEN nf, GEN x, long fl) { (void)fl; return idealappr(nf, x); }
3168 :
3169 : static GEN
3170 21 : mat_ideal_two_elt2(GEN nf, GEN x, GEN a)
3171 : {
3172 21 : GEN F = idealfactor(nf,a), P = gel(F,1), E = gel(F,2);
3173 21 : long i, r = lg(E);
3174 84 : for (i=1; i<r; i++) gel(E,i) = stoi( idealval(nf,x,gel(P,i)) );
3175 21 : return idealapprfact_i(nf,F,1);
3176 : }
3177 :
3178 : static void
3179 14 : not_in_ideal(GEN a) {
3180 14 : pari_err_DOMAIN("idealtwoelt2","element mod ideal", "!=", gen_0, a);
3181 0 : }
3182 : /* x integral in HNF, a an 'nf' */
3183 : static int
3184 28 : in_ideal(GEN x, GEN a)
3185 : {
3186 28 : switch(typ(a))
3187 : {
3188 14 : case t_INT: return dvdii(a, gcoeff(x,1,1));
3189 7 : case t_COL: return RgV_is_ZV(a) && !!hnf_invimage(x, a);
3190 7 : default: return 0;
3191 : }
3192 : }
3193 :
3194 : /* Given an integral ideal x and a in x, gives a b such that
3195 : * x = aZ_K + bZ_K using the approximation theorem */
3196 : GEN
3197 42 : idealtwoelt2(GEN nf, GEN x, GEN a)
3198 : {
3199 42 : pari_sp av = avma;
3200 : GEN cx, b;
3201 :
3202 42 : nf = checknf(nf);
3203 42 : a = nf_to_scalar_or_basis(nf, a);
3204 42 : x = idealhnf_shallow(nf,x);
3205 42 : if (lg(x) == 1)
3206 : {
3207 14 : if (!isintzero(a)) not_in_ideal(a);
3208 7 : set_avma(av); return gen_0;
3209 : }
3210 28 : x = Q_primitive_part(x, &cx);
3211 28 : if (cx) a = gdiv(a, cx);
3212 28 : if (!in_ideal(x, a)) not_in_ideal(a);
3213 21 : b = mat_ideal_two_elt2(nf, x, a);
3214 21 : if (typ(b) == t_COL)
3215 : {
3216 14 : GEN mod = idealhnf_principal(nf,a);
3217 14 : b = ZC_hnfrem(b,mod);
3218 14 : if (ZV_isscalar(b)) b = gel(b,1);
3219 : }
3220 : else
3221 : {
3222 7 : GEN aZ = typ(a) == t_COL? Q_denom(zk_inv(nf,a)): a; /* (a) \cap Z */
3223 7 : b = centermodii(b, aZ, shifti(aZ,-1));
3224 : }
3225 21 : b = cx? gmul(b,cx): gcopy(b);
3226 21 : return gc_upto(av, b);
3227 : }
3228 :
3229 : /* Given 2 integral ideals x and y in nf, returns a beta in nf such that
3230 : * beta * x is an integral ideal coprime to y */
3231 : GEN
3232 37191 : idealcoprimefact(GEN nf, GEN x, GEN fy)
3233 : {
3234 37191 : GEN L = gel(fy,1), e;
3235 37191 : long i, r = lg(L);
3236 :
3237 37191 : e = cgetg(r, t_COL);
3238 76053 : for (i=1; i<r; i++) gel(e,i) = stoi( -idealval(nf,x,gel(L,i)) );
3239 37190 : return idealapprfact_i(nf, mkmat2(L,e), 0);
3240 : }
3241 : GEN
3242 84 : idealcoprime(GEN nf, GEN x, GEN y)
3243 : {
3244 84 : pari_sp av = avma;
3245 84 : return gc_upto(av, idealcoprimefact(nf, x, idealfactor(nf,y)));
3246 : }
3247 :
3248 : GEN
3249 7 : nfmulmodpr(GEN nf, GEN x, GEN y, GEN modpr)
3250 : {
3251 7 : pari_sp av = avma;
3252 7 : GEN z, p, pr = modpr, T;
3253 :
3254 7 : nf = checknf(nf); modpr = nf_to_Fq_init(nf,&pr,&T,&p);
3255 0 : x = nf_to_Fq(nf,x,modpr);
3256 0 : y = nf_to_Fq(nf,y,modpr);
3257 0 : z = Fq_mul(x,y,T,p);
3258 0 : return gc_upto(av, algtobasis(nf, Fq_to_nf(z,modpr)));
3259 : }
3260 :
3261 : GEN
3262 0 : nfdivmodpr(GEN nf, GEN x, GEN y, GEN modpr)
3263 : {
3264 0 : pari_sp av = avma;
3265 0 : nf = checknf(nf);
3266 0 : return gc_upto(av, nfreducemodpr(nf, nfdiv(nf,x,y), modpr));
3267 : }
3268 :
3269 : GEN
3270 0 : nfpowmodpr(GEN nf, GEN x, GEN k, GEN modpr)
3271 : {
3272 0 : pari_sp av=avma;
3273 0 : GEN z, T, p, pr = modpr;
3274 :
3275 0 : nf = checknf(nf); modpr = nf_to_Fq_init(nf,&pr,&T,&p);
3276 0 : z = nf_to_Fq(nf,x,modpr);
3277 0 : z = Fq_pow(z,k,T,p);
3278 0 : return gc_upto(av, algtobasis(nf, Fq_to_nf(z,modpr)));
3279 : }
3280 :
3281 : GEN
3282 0 : nfkermodpr(GEN nf, GEN x, GEN modpr)
3283 : {
3284 0 : pari_sp av = avma;
3285 0 : GEN T, p, pr = modpr;
3286 :
3287 0 : nf = checknf(nf); modpr = nf_to_Fq_init(nf, &pr,&T,&p);
3288 0 : if (typ(x)!=t_MAT) pari_err_TYPE("nfkermodpr",x);
3289 0 : x = nfM_to_FqM(x, nf, modpr);
3290 0 : return gc_GEN(av, FqM_to_nfM(FqM_ker(x,T,p), modpr));
3291 : }
3292 :
3293 : GEN
3294 0 : nfsolvemodpr(GEN nf, GEN a, GEN b, GEN pr)
3295 : {
3296 0 : const char *f = "nfsolvemodpr";
3297 0 : pari_sp av = avma;
3298 : GEN T, p, modpr;
3299 :
3300 0 : nf = checknf(nf);
3301 0 : modpr = nf_to_Fq_init(nf, &pr,&T,&p);
3302 0 : if (typ(a)!=t_MAT) pari_err_TYPE(f,a);
3303 0 : a = nfM_to_FqM(a, nf, modpr);
3304 0 : switch(typ(b))
3305 : {
3306 0 : case t_MAT:
3307 0 : b = nfM_to_FqM(b, nf, modpr);
3308 0 : b = FqM_gauss(a,b,T,p);
3309 0 : if (!b) pari_err_INV(f,a);
3310 0 : a = FqM_to_nfM(b, modpr);
3311 0 : break;
3312 0 : case t_COL:
3313 0 : b = nfV_to_FqV(b, nf, modpr);
3314 0 : b = FqM_FqC_gauss(a,b,T,p);
3315 0 : if (!b) pari_err_INV(f,a);
3316 0 : a = FqV_to_nfV(b, modpr);
3317 0 : break;
3318 0 : default: pari_err_TYPE(f,b);
3319 : }
3320 0 : return gc_GEN(av, a);
3321 : }
|