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 16096248 : idealtyp(GEN *ideal, GEN *arch)
47 : {
48 16096248 : GEN x = *ideal;
49 16096248 : long t,lx,tx = typ(x);
50 :
51 16096248 : if (tx!=t_VEC || lg(x)!=3) { if (arch) *arch = NULL; }
52 : else
53 : {
54 1239719 : GEN a = gel(x,2);
55 1239719 : 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 1239705 : if (arch) *arch = a;
62 1239712 : x = gel(x,1); tx = typ(x);
63 : }
64 16096241 : switch(tx)
65 : {
66 12163327 : case t_MAT: lx = lg(x);
67 12163327 : if (lx == 1) { t = id_PRINCIPAL; x = gen_0; break; }
68 12163166 : if (lx != lgcols(x)) pari_err_TYPE("idealtyp [nonsquare t_MAT]",x);
69 12163135 : t = id_MAT;
70 12163135 : break;
71 :
72 2884764 : case t_VEC:
73 2884764 : if (!checkprid_i(x)) pari_err_TYPE("idealtyp [fake prime ideal]",x);
74 2884718 : t = id_PRIME; break;
75 :
76 1048438 : case t_POL: case t_POLMOD: case t_COL:
77 : case t_INT: case t_FRAC:
78 1048438 : t = id_PRINCIPAL; break;
79 6 : default:
80 6 : pari_err_TYPE("idealtyp",x);
81 : return 0; /*LCOV_EXCL_LINE*/
82 : }
83 16096452 : *ideal = x; return t;
84 : }
85 :
86 : /* true nf; v = [a,x,...], a in Z. Return (a,x) */
87 : GEN
88 739306 : idealhnf_two(GEN nf, GEN v)
89 : {
90 739306 : GEN p = gel(v,1), pi = gel(v,2), m = zk_scalar_or_multable(nf, pi);
91 739277 : if (typ(m) == t_INT) return scalarmat(gcdii(m,p), nf_get_degree(nf));
92 666070 : return ZM_hnfmodid(m, p);
93 : }
94 : /* true nf */
95 : GEN
96 3581372 : pr_hnf(GEN nf, GEN pr)
97 : {
98 3581372 : GEN p = pr_get_p(pr), m;
99 3581338 : if (pr_is_inert(pr)) return scalarmat(p, nf_get_degree(nf));
100 3145026 : m = zk_scalar_or_multable(nf, pr_get_gen(pr));
101 3144646 : return ZM_hnfmodprime(m, p);
102 : }
103 :
104 : GEN
105 1370798 : idealhnf_principal(GEN nf, GEN x)
106 : {
107 : GEN cx;
108 1370798 : x = nf_to_scalar_or_basis(nf, x);
109 1370793 : switch(typ(x))
110 : {
111 1052846 : case t_COL: break;
112 284729 : case t_INT: if (!signe(x)) return cgetg(1,t_MAT);
113 283098 : return scalarmat(absi_shallow(x), nf_get_degree(nf));
114 33218 : case t_FRAC:
115 33218 : return scalarmat(Q_abs_shallow(x), nf_get_degree(nf));
116 0 : default: pari_err_TYPE("idealhnf",x);
117 : }
118 1052846 : x = Q_primitive_part(x, &cx);
119 1052838 : RgV_check_ZV(x, "idealhnf");
120 1052838 : x = zk_multable(nf, x);
121 1052841 : x = ZM_hnfmodid(x, zkmultable_capZ(x));
122 1052845 : 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 1747621 : idealhnf_shallow(GEN nf, GEN x)
156 : {
157 1747621 : long tx = typ(x), lx = lg(x), N;
158 :
159 : /* cannot use idealtyp because here we allow nonsquare matrices */
160 1747621 : if (tx == t_VEC && lx == 3) { x = gel(x,1); tx = typ(x); lx = lg(x); }
161 1747621 : if (tx == t_VEC && lx == 6)
162 : {
163 484985 : if (!checkprid_i(x)) pari_err_TYPE("idealhnf [fake prime ideal]",x);
164 484974 : return pr_hnf(nf,x); /* PRIME */
165 : }
166 1262636 : switch(tx)
167 : {
168 103098 : case t_MAT:
169 : {
170 : GEN cx;
171 103098 : long nx = lx-1;
172 103098 : N = nf_get_degree(nf);
173 103098 : if (nx == 0) return cgetg(1, t_MAT);
174 103077 : if (nbrows(x) != N) pari_err_TYPE("idealhnf [wrong dimension]",x);
175 103070 : if (nx == 1) return idealhnf_principal(nf, gel(x,1));
176 :
177 86334 : if (nx == N && RgM_is_ZM(x) && ZM_ishnf(x)) return x;
178 49609 : x = Q_primitive_part(x, &cx);
179 49609 : if (nx < N)
180 7 : x = vec_mulid(nf, x); /* build ZK-module generated from cols */
181 : else
182 49602 : x = ZM_hnfmod(x, ZM_detmult(x)); /* assume Z-span cols is ZK-module */
183 49609 : 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 gerepileupto(av, idealhnf_two(nf, mkvec2(A,B)));
204 : }
205 1159524 : default: return idealhnf_principal(nf, x); /* PRINCIPAL */
206 : }
207 : }
208 : /* true nf */
209 : GEN
210 574 : idealhnf(GEN nf, GEN x)
211 : {
212 574 : pari_sp av = avma;
213 574 : GEN y = idealhnf_shallow(nf, x);
214 560 : return (avma == av)? gcopy(y): gerepileupto(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 gerepileupto(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 gerepileupto(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 229093 : ok_elt(GEN x, GEN xZ, GEN y)
377 : {
378 229093 : pari_sp av = avma;
379 229093 : 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 65503 : addmul_mat(GEN a, GEN s, GEN b)
385 : {
386 65503 : if (!signe(s)) return a;
387 56666 : if (!equali1(s)) b = ZM_Z_mul(b, s);
388 56665 : return a? ZM_add(a, b): b;
389 : }
390 :
391 : static GEN
392 122186 : get_random_a(GEN nf, GEN x, GEN xZ)
393 : {
394 : pari_sp av;
395 122186 : long i, lm, l = lg(x);
396 : GEN z, beta, mul;
397 :
398 122186 : beta= cgetg(l, t_MAT);
399 122186 : mul = cgetg(l, t_VEC); lm = 1; /* = lg(mul) */
400 : /* look for a in x such that a O/xZ = x O/xZ */
401 256347 : for (i = 2; i < l; i++)
402 : {
403 246434 : GEN xi = gel(x,i);
404 246434 : GEN t = FpM_red(zk_multable(nf,xi), xZ); /* ZM, cannot be a scalar */
405 246426 : if (gequal0(t)) continue;
406 201189 : if (ok_elt(x,xZ, t)) return xi;
407 88918 : gel(beta,lm) = xi;
408 : /* mul[i] = { canonical generators for x[i] O/xZ as Z-module } */
409 88918 : gel(mul,lm) = t; lm++;
410 : }
411 9913 : setlg(mul, lm);
412 9913 : setlg(beta,lm); z = cgetg(lm, t_VEC);
413 29810 : for(av = avma;; set_avma(av))
414 19897 : {
415 29810 : GEN a = NULL;
416 95313 : for (i = 1; i < lm; i++)
417 : {
418 65503 : gel(z,i) = randomi(xZ);
419 65503 : 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 29810 : if (a && ok_elt(x,xZ, a)) break;
423 : }
424 9913 : return ZM_ZC_mul(beta, z);
425 : }
426 :
427 : /* x square matrix, assume it is HNF */
428 : static GEN
429 263938 : mat_ideal_two_elt(GEN nf, GEN x)
430 : {
431 : GEN y, a, cx, xZ;
432 263938 : long N = nf_get_degree(nf);
433 : pari_sp av, tetpil;
434 :
435 263939 : if (lg(x)-1 != N) pari_err_DIM("idealtwoelt");
436 263925 : if (N == 2) return mkvec2copy(gcoeff(x,1,1), gel(x,2));
437 :
438 141052 : y = cgetg(3,t_VEC); av = avma;
439 141052 : cx = Q_content(x);
440 141053 : xZ = gcoeff(x,1,1);
441 141053 : if (gequal(xZ, cx)) /* x = (cx) */
442 : {
443 7973 : gel(y,1) = cx;
444 7973 : gel(y,2) = gen_0; return y;
445 : }
446 133080 : if (equali1(cx)) cx = NULL;
447 : else
448 : {
449 3744 : x = Q_div_to_int(x, cx);
450 3744 : xZ = gcoeff(x,1,1);
451 : }
452 133080 : if (N < 6)
453 113778 : a = get_random_a(nf, x, xZ);
454 : else
455 : {
456 19302 : 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 19302 : GEN P, E, a1 = Z_lsmoothen(xZ, (GEN)FB, &P, &E);
460 19302 : if (!a1) /* factors completely */
461 10894 : a = idealapprfact_i(nf, idealfactor(nf,x), 1);
462 8408 : else if (lg(P) == 1) /* no small factors */
463 2603 : 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 5805 : a0 = diviiexact(xZ, a1);
468 5805 : A0 = ZM_hnfmodid(x, a0); /* smooth part of x */
469 5805 : A1 = ZM_hnfmodid(x, a1); /* cofactor */
470 5805 : pi0 = idealapprfact_i(nf, idealfactor(nf,A0), 1);
471 5805 : pi1 = get_random_a(nf, A1, a1);
472 5805 : (void)bezout(a0, a1, &v0,&v1);
473 5805 : u0 = mulii(a0, v0);
474 5805 : u1 = mulii(a1, v1);
475 5805 : if (typ(pi0) != t_COL) t = addmulii(u0, pi0, u1);
476 : else
477 5805 : { t = ZC_Z_mul(pi0, u1); gel(t,1) = addii(gel(t,1), u0); }
478 5805 : u = ZC_Z_mul(pi1, u0); gel(u,1) = addii(gel(u,1), u1);
479 5805 : a = nfmuli(nf, centermod(u, xZ), centermod(t, xZ));
480 : }
481 : }
482 133079 : if (cx)
483 : {
484 3744 : a = centermod(a, xZ);
485 3744 : tetpil = avma;
486 3744 : if (typ(cx) == t_INT)
487 : {
488 98 : gel(y,1) = mulii(xZ, cx);
489 98 : gel(y,2) = ZC_Z_mul(a, cx);
490 : }
491 : else
492 : {
493 3646 : gel(y,1) = gmul(xZ, cx);
494 3646 : gel(y,2) = RgC_Rg_mul(a, cx);
495 : }
496 : }
497 : else
498 : {
499 129335 : tetpil = avma;
500 129335 : gel(y,1) = icopy(xZ);
501 129335 : gel(y,2) = centermod(a, xZ);
502 : }
503 133079 : gerepilecoeffssp(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 108949 : idealtwoelt(GEN nf, GEN x)
512 : {
513 : pari_sp av;
514 108949 : long tx = idealtyp(&x, NULL);
515 108942 : nf = checknf(nf);
516 108941 : if (tx == id_MAT) return mat_ideal_two_elt(nf,x);
517 1099 : if (tx == id_PRIME) return mkvec2copy(gel(x,1), gel(x,2));
518 : /* id_PRINCIPAL */
519 1078 : av = avma; x = nf_to_scalar_or_basis(nf, x);
520 1960 : return gerepilecopy(av, typ(x)==t_COL? mkvec2(gen_0,x):
521 973 : 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 4045851 : idealHNF_norm_pval(GEN x, GEN p, long Zval)
532 : {
533 4045851 : long i, v = Zval, l = lg(x);
534 31744836 : for (i = 2; i < l; i++) v += Z_pval(gcoeff(x,i,i), p);
535 4045886 : 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 285639 : idealHNF_Z_factor_i(GEN x, GEN f0, GEN *pvN, GEN *pvZ)
542 : {
543 285639 : GEN P, E, vN, vZ, xZ = gcoeff(x,1,1), f = f0? f0: Z_factor(xZ);
544 : long i, l;
545 285658 : P = gel(f,1); l = lg(P);
546 285658 : E = gel(f,2);
547 285658 : *pvN = vN = cgetg(l, t_VECSMALL);
548 285663 : *pvZ = vZ = cgetg(l, t_VECSMALL);
549 714387 : for (i = 1; i < l; i++)
550 : {
551 428725 : GEN p = gel(P,i);
552 428725 : vZ[i] = f0? Z_pval(xZ, p): (long) itou(gel(E,i));
553 428724 : vN[i] = idealHNF_norm_pval(x,p, vZ[i]);
554 : }
555 285662 : 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 4181553 : idealHNF_val(GEN A, GEN P, long Nval, long Zval)
567 : {
568 4181553 : 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 4181552 : if (Nval < f) return 0;
573 4178375 : p = pr_get_p(P);
574 4178372 : e = pr_get_e(P);
575 : /* v_P(A) <= max [ e * v_p(A \cap Z), floor[v_p(Nix) / f ] */
576 4178375 : vmax = minss(Zval * e, Nval / f);
577 4178374 : mul = pr_get_tau(P);
578 4178369 : l = lg(mul);
579 4178369 : B = cgetg(l,t_MAT);
580 : /* B[1] not needed: v_pr(A[1]) = v_pr(A \cap Z) is known already */
581 4185474 : gel(B,1) = gen_0; /* dummy */
582 23274999 : for (j = 2; j < l; j++)
583 : {
584 20921251 : GEN x = gel(A,j);
585 20921251 : gel(B,j) = y = cgetg(l, t_COL);
586 227807051 : for (i = 1; i < l; i++)
587 : { /* compute a = (x.t0)_i, A in HNF ==> x[j+1..l-1] = 0 */
588 208717526 : a = mulii(gel(x,1), gcoeff(mul,i,1));
589 1461356623 : for (k = 2; k <= j; k++) a = addii(a, mulii(gel(x,k), gcoeff(mul,i,k)));
590 : /* p | a ? */
591 208710279 : gel(y,i) = dvmdii(a,p,&r); if (signe(r)) return 0;
592 : }
593 : }
594 2353748 : vals = cgetg(l, t_VECSMALL);
595 : /* vals[1] not needed */
596 16599377 : for (j = 2; j < l; j++)
597 : {
598 14245444 : gel(B,j) = Q_primitive_part(gel(B,j), &cx);
599 14245465 : vals[j] = cx? 1 + e * Q_pval(cx, p): 1;
600 : }
601 2353933 : pk = powiu(p, ceildivuu(vmax, e));
602 2353921 : av = avma; y = cgetg(l,t_COL);
603 : /* can compute mod p^ceil((vmax-v)/e) */
604 3914562 : for (v = 1; v < vmax; v++)
605 : { /* we know v_pr(Bj) >= v for all j */
606 1590640 : if (e == 1 || (vmax - v) % e == 0) pk = diviiexact(pk, p);
607 8307606 : for (j = 2; j < l; j++)
608 : {
609 6747000 : GEN x = gel(B,j); if (v < vals[j]) continue;
610 44574726 : for (i = 1; i < l; i++)
611 : {
612 40069205 : pari_sp av2 = avma;
613 40069205 : a = mulii(gel(x,1), gcoeff(mul,i,1));
614 518287239 : 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 40068440 : a = dvmdii(a,p,&r); if (signe(r)) return v;
617 40038257 : if (lgefint(a) > lgefint(pk)) a = remii(a, pk);
618 40038257 : gel(y,i) = gerepileuptoint(av2, a);
619 : }
620 4505521 : gel(B,j) = y; y = x;
621 4505521 : if (gc_needed(av,3))
622 : {
623 0 : if(DEBUGMEM>1) pari_warn(warnmem,"idealval");
624 0 : gerepileall(av,3, &y,&B,&pk);
625 : }
626 : }
627 : }
628 2323922 : 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 285638 : idealHNF_factor_i(GEN nf, GEN x, GEN cx, GEN FA)
635 : {
636 285638 : const long N = lg(x)-1;
637 : long i, j, k, l, v;
638 285638 : GEN vN, vZ, vP, vE, vp = idealHNF_Z_factor_i(x, FA, &vN,&vZ);
639 :
640 285662 : l = lg(vp);
641 285662 : i = cx? expi(cx)+1: 1;
642 285663 : vP = cgetg((l+i-2)*N+1, t_COL);
643 285664 : vE = cgetg((l+i-2)*N+1, t_COL);
644 714379 : for (i = k = 1; i < l; i++)
645 : {
646 428725 : GEN L, p = gel(vp,i);
647 428725 : long Nval = vN[i], Zval = vZ[i], vc = cx? Z_pvalrem(cx,p,&cx): 0;
648 428719 : if (vc)
649 : {
650 45908 : L = idealprimedec(nf,p);
651 45908 : if (is_pm1(cx)) cx = NULL;
652 : }
653 : else
654 382811 : L = idealprimedec_limit_f(nf,p,Nval);
655 993105 : for (j = 1; Nval && j < lg(L); j++) /* !Nval => only cx contributes */
656 : {
657 564397 : GEN P = gel(L,j);
658 564397 : pari_sp av = avma;
659 564397 : v = idealHNF_val(x, P, Nval, Zval);
660 564370 : set_avma(av);
661 564368 : Nval -= v*pr_get_f(P);
662 564371 : v += vc * pr_get_e(P); if (!v) continue;
663 474159 : gel(vP,k) = P;
664 474159 : gel(vE,k) = utoipos(v); k++;
665 : }
666 473474 : if (vc) for (; j<lg(L); j++)
667 : {
668 44768 : GEN P = gel(L,j);
669 44768 : gel(vP,k) = P;
670 44768 : gel(vE,k) = utoipos(vc * pr_get_e(P)); k++;
671 : }
672 : }
673 285654 : if (cx && !FA)
674 : { /* complete factorization */
675 73800 : GEN f = Z_factor(cx), cP = gel(f,1), cE = gel(f,2);
676 73800 : long lc = lg(cP);
677 160599 : for (i=1; i<lc; i++)
678 : {
679 86799 : GEN p = gel(cP,i), L = idealprimedec(nf,p);
680 86799 : long vc = itos(gel(cE,i));
681 190069 : for (j=1; j<lg(L); j++)
682 : {
683 103270 : GEN P = gel(L,j);
684 103270 : gel(vP,k) = P;
685 103270 : gel(vE,k) = utoipos(vc * pr_get_e(P)); k++;
686 : }
687 : }
688 : }
689 285654 : setlg(vP, k);
690 285659 : setlg(vE, k); return mkmat2(vP, vE);
691 : }
692 : /* true nf, x integral ideal */
693 : static GEN
694 239943 : idealHNF_factor(GEN nf, GEN x, ulong lim)
695 : {
696 239943 : GEN cx, F = NULL;
697 239943 : if (lim)
698 : {
699 : GEN P, E;
700 : long i;
701 : /* strict useless because of prime table */
702 70 : F = absZ_factor_limit(gcoeff(x,1,1), lim);
703 70 : P = gel(F,1);
704 70 : E = gel(F,2);
705 : /* filter out entries > lim */
706 119 : for (i = lg(P)-1; i; i--)
707 119 : if (cmpiu(gel(P,i), lim) <= 0) break;
708 70 : setlg(P, i+1);
709 70 : setlg(E, i+1);
710 : }
711 239943 : x = Q_primitive_part(x, &cx);
712 239927 : 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 28077 : prV_e_muls(GEN L, long c)
717 : {
718 28077 : long j, l = lg(L);
719 28077 : GEN z = cgetg(l, t_COL);
720 57785 : for (j = 1; j < l; j++) gel(z,j) = stoi(c * pr_get_e(gel(L,j)));
721 28071 : return z;
722 : }
723 : /* true nf, y in Q */
724 : static GEN
725 28059 : Q_nffactor(GEN nf, GEN y, ulong lim)
726 : {
727 : GEN f, P, E;
728 : long l, i;
729 28059 : if (typ(y) == t_INT)
730 : {
731 28031 : if (!signe(y)) pari_err_DOMAIN("idealfactor", "ideal", "=",gen_0,y);
732 28017 : if (is_pm1(y)) return trivial_fact();
733 : }
734 17642 : y = Q_abs_shallow(y);
735 17643 : if (!lim) f = Q_factor(y);
736 : else
737 : {
738 154 : f = Q_factor_limit(y, lim);
739 154 : P = gel(f,1);
740 154 : E = gel(f,2);
741 224 : for (i = lg(P)-1; i > 0; i--)
742 203 : if (abscmpiu(gel(P,i), lim) < 0) break;
743 154 : setlg(P,i+1); setlg(E,i+1);
744 : }
745 17645 : P = gel(f,1); l = lg(P); if (l == 1) return f;
746 17624 : E = gel(f,2);
747 45732 : for (i = 1; i < l; i++)
748 : {
749 28124 : gel(P,i) = idealprimedec(nf, gel(P,i));
750 28075 : gel(E,i) = prV_e_muls(gel(P,i), itos(gel(E,i)));
751 : }
752 17608 : P = shallowconcat1(P); gel(f,1) = P; settyp(P, t_COL);
753 17628 : E = shallowconcat1(E); gel(f,2) = E; return f;
754 : }
755 :
756 : GEN
757 25501 : idealfactor_partial(GEN nf, GEN x, GEN L)
758 : {
759 25501 : pari_sp av = avma;
760 : long i, j, l;
761 : GEN P, E;
762 25501 : 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 23842 : P = cgetg(l, t_VEC);
766 89915 : for (i = 1; i < l; i++)
767 : {
768 66073 : GEN p = gel(L,i);
769 66073 : gel(P,i) = typ(p) == t_INT? idealprimedec(nf, p): mkvec(p);
770 : }
771 23842 : P = shallowconcat1(P); settyp(P, t_COL);
772 23842 : P = gen_sort_uniq(P, (void*)&cmp_prime_ideal, &cmp_nodata);
773 23842 : E = cgetg_copy(P, &l);
774 114604 : for (i = j = 1; i < l; i++)
775 : {
776 90762 : long v = idealval(nf, x, gel(P,i));
777 90762 : if (v) { gel(P,j) = gel(P,i); gel(E,j) = stoi(v); j++; }
778 : }
779 23842 : setlg(P,j);
780 23842 : setlg(E,j); return gerepilecopy(av, mkmat2(P, E));
781 : }
782 : GEN
783 268124 : idealfactor_limit(GEN nf, GEN x, ulong lim)
784 : {
785 268124 : pari_sp av = avma;
786 : GEN fa, y;
787 268124 : long tx = idealtyp(&x, NULL);
788 :
789 268101 : 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 267982 : nf = checknf(nf);
795 267979 : if (tx == id_PRINCIPAL)
796 : {
797 29537 : y = nf_to_scalar_or_basis(nf, x);
798 29536 : if (typ(y) != t_COL) return gerepilecopy(av, Q_nffactor(nf, y, lim));
799 : }
800 239919 : y = idealnumden(nf, x);
801 239930 : fa = idealHNF_factor(nf, gel(y,1), lim);
802 239929 : if (!isint1(gel(y,2)))
803 14 : fa = famat_div_shallow(fa, idealHNF_factor(nf, gel(y,2), lim));
804 239929 : fa = gerepilecopy(av, fa);
805 239935 : return sort_factor(fa, (void*)&cmp_prime_ideal, &cmp_nodata);
806 : }
807 : GEN
808 267746 : idealfactor(GEN nf, GEN x) { return idealfactor_limit(nf, x, 0); }
809 : GEN
810 182 : gpidealfactor(GEN nf, GEN x, GEN lim)
811 : {
812 182 : ulong L = 0;
813 182 : if (lim)
814 : {
815 70 : if (typ(lim) != t_INT || signe(lim) < 0) pari_err_FLAG("idealfactor");
816 70 : L = itou(lim);
817 : }
818 182 : return idealfactor_limit(nf, x, L);
819 : }
820 :
821 : static GEN
822 7775 : ramified_root(GEN nf, GEN R, GEN A, long n)
823 : {
824 7775 : GEN v, P = gel(idealfactor(nf, R), 1);
825 7775 : long i, l = lg(P);
826 7775 : v = cgetg(l, t_VECSMALL);
827 8426 : for (i = 1; i < l; i++)
828 : {
829 658 : long w = idealval(nf, A, gel(P,i));
830 658 : if (w % n) return NULL;
831 651 : v[i] = w / n;
832 : }
833 7768 : 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 7782 : idealsqrtn_int(GEN nf, GEN A, long n, GEN *pB)
859 : {
860 : GEN C, root;
861 : long i, l;
862 :
863 7782 : if (typ(A) == t_MAT && ZM_isscalar(A, NULL)) A = gcoeff(A,1,1);
864 7782 : if (typ(A) == t_INT) /* > 0 */
865 : {
866 5564 : GEN P = nf_get_ramified_primes(nf), v, q;
867 5564 : l = lg(P); v = cgetg(l, t_VECSMALL);
868 24981 : for (i = 1; i < l; i++) v[i] = Z_pvalrem(A, gel(P,i), &A);
869 5564 : C = gen_1;
870 5564 : if (!isint1(A) && !Z_ispowerall(A, n, pB? &C: NULL)) return 0;
871 5564 : if (!pB) return ramified_root_simple(nf, n, P, v);
872 5557 : q = factorback2(P, v);
873 5557 : root = ramified_root(nf, q, q, n);
874 5557 : if (!root) return 0;
875 5557 : if (!equali1(C)) root = isint1(root)? C: ZM_Z_mul(root, C);
876 5557 : *pB = root; return 1;
877 : }
878 : /* compute valuations at ramified primes */
879 2218 : root = ramified_root(nf, idealadd(nf, nf_get_diff(nf), A), A, n);
880 2218 : if (!root) return 0;
881 : /* remove ramified primes */
882 2211 : if (isint1(root))
883 1861 : root = matid(nf_get_degree(nf));
884 : else
885 350 : A = idealdivexact(nf, A, idealpows(nf,root,n));
886 2211 : A = Q_primitive_part(A, &C);
887 2211 : 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 2204 : for (i = 0;; i++)
895 2071 : {
896 4275 : GEN J, b, a = gcoeff(A,1,1); /* A \cap Z */
897 4275 : if (is_pm1(a)) break;
898 2099 : if (!Z_ispowerall(a,n,&b)) return 0;
899 2071 : J = idealadd(nf, b, A);
900 2071 : A = idealdivexact(nf, idealpows(nf,J,n), A);
901 : /* div and not divexact here */
902 2071 : if (pB) root = odd(i)? idealdiv(nf, root, J): idealmul(nf, root, J);
903 : }
904 2176 : if (pB) *pB = root;
905 2176 : 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 3919 : idealispower(GEN nf, GEN A, long n, GEN *pB)
912 : {
913 3919 : pari_sp av = avma;
914 : GEN v, N, D;
915 3919 : nf = checknf(nf);
916 3919 : if (n <= 0) pari_err_DOMAIN("idealispower", "n", "<=", gen_0, stoi(n));
917 3919 : if (n == 1) { if (pB) *pB = idealhnf(nf,A); return 1; }
918 3912 : v = idealnumden(nf,A);
919 3912 : if (gequal0(gel(v,1))) { set_avma(av); if (pB) *pB = cgetg(1,t_MAT); return 1; }
920 3912 : if (!idealsqrtn_int(nf, gel(v,1), n, pB? &N: NULL)) return 0;
921 3870 : if (!idealsqrtn_int(nf, gel(v,2), n, pB? &D: NULL)) return 0;
922 3870 : if (pB) *pB = gerepileupto(av, idealdiv(nf,N,D)); else set_avma(av);
923 3870 : return 1;
924 : }
925 :
926 : /* x t_INT or integral nonzero ideal in HNF */
927 : static GEN
928 98112 : idealredmodpower_i(GEN nf, GEN x, ulong k, ulong B)
929 : {
930 : GEN cx, y, U, N, F, Q;
931 98112 : if (typ(x) == t_INT)
932 : {
933 51912 : if (!signe(x) || is_pm1(x)) return gen_1;
934 2142 : F = Z_factor_limit(x, B);
935 2142 : gel(F,2) = gdiventgs(gel(F,2), k);
936 2142 : return ginv(factorback(F));
937 : }
938 46200 : N = gcoeff(x,1,1); if (is_pm1(N)) return gen_1;
939 45717 : F = absZ_factor_limit_strict(N, B, &U);
940 45716 : if (U)
941 : {
942 146 : GEN M = powii(gel(U,1), gel(U,2));
943 146 : y = hnfmodid(x, M); /* coprime part to B! */
944 146 : if (!idealispower(nf, y, k, &U)) U = NULL;
945 146 : x = hnfmodid(x, diviiexact(N, M));
946 : }
947 : /* x = B-smooth part of initial x */
948 45716 : x = Q_primitive_part(x, &cx);
949 45717 : F = idealHNF_factor_i(nf, x, cx, F);
950 45717 : gel(F,2) = gdiventgs(gel(F,2), k);
951 45717 : Q = idealfactorback(nf, gel(F,1), gel(F,2), 0);
952 45716 : if (U) Q = idealmul(nf,Q,U);
953 45717 : if (typ(Q) == t_INT) return Q;
954 13886 : y = idealred_elt(nf, idealHNF_inv_Z(nf, Q));
955 13886 : return gdiv(y, gcoeff(Q,1,1));
956 : }
957 : GEN
958 49063 : idealredmodpower(GEN nf, GEN x, ulong n, ulong B)
959 : {
960 49063 : pari_sp av = avma;
961 : GEN a, b;
962 49063 : nf = checknf(nf);
963 49063 : if (!n) pari_err_DOMAIN("idealredmodpower","n", "=", gen_0, gen_0);
964 49063 : x = idealnumden(nf, x);
965 49063 : a = gel(x,1);
966 49063 : if (isintzero(a)) { set_avma(av); return gen_1; }
967 49056 : a = idealredmodpower_i(nf, gel(x,1), n, B);
968 49056 : b = idealredmodpower_i(nf, gel(x,2), n, B);
969 49056 : if (!isint1(b)) a = nf_to_scalar_or_basis(nf, nfdiv(nf, a, b));
970 49056 : return gerepilecopy(av, a);
971 : }
972 :
973 : /* P prime ideal in idealprimedec format. Return valuation(A) at P */
974 : long
975 9399369 : idealval(GEN nf, GEN A, GEN P)
976 : {
977 9399369 : pari_sp av = avma;
978 : GEN p, cA;
979 9399369 : long vcA, v, Zval, tx = idealtyp(&A, NULL);
980 :
981 9399364 : if (tx == id_PRINCIPAL) return nfval(nf,A,P);
982 9292115 : checkprid(P);
983 9292358 : if (tx == id_PRIME) return pr_equal(P, A)? 1: 0;
984 : /* id_MAT */
985 9292330 : nf = checknf(nf);
986 9292517 : A = Q_primitive_part(A, &cA);
987 9292676 : p = pr_get_p(P);
988 9292696 : vcA = cA? Q_pval(cA,p): 0;
989 9292697 : if (pr_is_inert(P)) return gc_long(av,vcA);
990 8894255 : Zval = Z_pval(gcoeff(A,1,1), p);
991 8894283 : if (!Zval) v = 0;
992 : else
993 : {
994 3617050 : long Nval = idealHNF_norm_pval(A, p, Zval);
995 3617158 : v = idealHNF_val(A, P, Nval, Zval);
996 : }
997 8894190 : 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 109602 : idealadd(GEN nf, GEN x, GEN y)
1010 : {
1011 109602 : pari_sp av = avma;
1012 : long tx, ty;
1013 : GEN z, a, dx, dy, dz;
1014 :
1015 109602 : tx = idealtyp(&x, NULL);
1016 109602 : ty = idealtyp(&y, NULL); nf = checknf(nf);
1017 109602 : if (tx != id_MAT) x = idealhnf_shallow(nf,x);
1018 109602 : if (ty != id_MAT) y = idealhnf_shallow(nf,y);
1019 109602 : if (lg(x) == 1) return gerepilecopy(av,y);
1020 108608 : if (lg(y) == 1) return gerepilecopy(av,x); /* check for 0 ideal */
1021 108153 : dx = Q_denom(x);
1022 108153 : dy = Q_denom(y); dz = lcmii(dx,dy);
1023 108153 : if (is_pm1(dz)) dz = NULL; else {
1024 15995 : x = Q_muli_to_int(x, dz);
1025 15995 : y = Q_muli_to_int(y, dz);
1026 : }
1027 108153 : a = gcdii(gcoeff(x,1,1), gcoeff(y,1,1));
1028 108153 : if (is_pm1(a))
1029 : {
1030 38639 : long N = lg(x)-1;
1031 38639 : if (!dz) { set_avma(av); return matid(N); }
1032 3997 : return gerepileupto(av, scalarmat(ginv(dz), N));
1033 : }
1034 69514 : z = ZM_hnfmodid(shallowconcat(x,y), a);
1035 69514 : if (dz) z = RgM_Rg_div(z,dz);
1036 69514 : return gerepileupto(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 732639 : _idealaddtoone(GEN nf, GEN x, GEN y, long red)
1045 : {
1046 : GEN a;
1047 732639 : long tx = idealtyp(&x, NULL);
1048 732609 : long ty = idealtyp(&y, NULL);
1049 : long ea;
1050 732612 : if (tx != id_MAT) x = idealhnf_shallow(nf, x);
1051 732615 : if (ty != id_MAT) y = idealhnf_shallow(nf, y);
1052 732615 : if (lg(x) == 1)
1053 14 : a = trivial_merge(y);
1054 732601 : else if (lg(y) == 1)
1055 14 : a = trivial_merge(x);
1056 : else
1057 732587 : a = hnfmerge_get_1(x, y);
1058 732555 : if (!a) pari_err_COPRIME("idealaddtoone",x,y);
1059 732542 : if (red && (ea = gexpo(a)) > 10)
1060 : {
1061 5229 : GEN b = (typ(a) == t_COL)? a: scalarcol_shallow(a, nf_get_degree(nf));
1062 5229 : b = ZC_reducemodlll(b, idealHNF_mul(nf,x,y));
1063 5229 : if (gexpo(b) < ea) a = b;
1064 : }
1065 732542 : return a;
1066 : }
1067 : /* true nf */
1068 : GEN
1069 19754 : idealaddtoone_i(GEN nf, GEN x, GEN y)
1070 19754 : { return _idealaddtoone(nf, x, y, 1); }
1071 : /* true nf */
1072 : GEN
1073 712882 : idealaddtoone_raw(GEN nf, GEN x, GEN y)
1074 712882 : { 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 = gerepileupto(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 gerepilecopy(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 931509 : idealHNF_mul_two(GEN nf, GEN x, GEN y)
1142 : {
1143 931509 : GEN m, a = gel(y,1), alpha = gel(y,2);
1144 : long i, N;
1145 :
1146 931509 : if (typ(alpha) != t_MAT)
1147 : {
1148 589379 : alpha = zk_scalar_or_multable(nf, alpha);
1149 589403 : if (typ(alpha) == t_INT) /* e.g. y inert ? 0 should not (but may) occur */
1150 14153 : return signe(a)? ZM_Z_mul(x, gcdii(a, alpha)): cgetg(1,t_MAT);
1151 : }
1152 917380 : N = lg(x)-1; m = cgetg((N<<1)+1,t_MAT);
1153 3668132 : for (i=1; i<=N; i++) gel(m,i) = ZM_ZC_mul(alpha,gel(x,i));
1154 3666830 : for (i=1; i<=N; i++) gel(m,i+N) = ZC_Z_mul(gel(x,i), a);
1155 916501 : 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 485685 : idealHNF_mul(GEN nf, GEN x, GEN y)
1162 : {
1163 : GEN z;
1164 485685 : if (typ(y) == t_VEC)
1165 303404 : z = idealHNF_mul_two(nf,x,y);
1166 : else
1167 : { /* reduce one ideal to two-elt form. The smallest */
1168 182281 : GEN xZ = gcoeff(x,1,1), yZ = gcoeff(y,1,1);
1169 182281 : if (cmpii(xZ, yZ) < 0)
1170 : {
1171 39356 : if (is_pm1(xZ)) return gcopy(y);
1172 23210 : z = idealHNF_mul_two(nf, y, mat_ideal_two_elt(nf,x));
1173 : }
1174 : else
1175 : {
1176 142926 : if (is_pm1(yZ)) return gcopy(x);
1177 36046 : z = idealHNF_mul_two(nf, x, mat_ideal_two_elt(nf,y));
1178 : }
1179 : }
1180 362683 : return z;
1181 : }
1182 :
1183 : /* operations on elements in factored form */
1184 :
1185 : GEN
1186 200420 : famat_mul_shallow(GEN f, GEN g)
1187 : {
1188 200420 : if (typ(f) != t_MAT) f = to_famat_shallow(f,gen_1);
1189 200420 : if (typ(g) != t_MAT) g = to_famat_shallow(g,gen_1);
1190 200420 : if (lgcols(f) == 1) return g;
1191 148098 : if (lgcols(g) == 1) return f;
1192 146173 : return mkmat2(shallowconcat(gel(f,1), gel(g,1)),
1193 146173 : shallowconcat(gel(f,2), gel(g,2)));
1194 : }
1195 : GEN
1196 88283 : famat_mulpow_shallow(GEN f, GEN g, GEN e)
1197 : {
1198 88283 : if (!signe(e)) return f;
1199 54577 : return famat_mul_shallow(f, famat_pow_shallow(g, e));
1200 : }
1201 :
1202 : GEN
1203 118508 : famat_mulpows_shallow(GEN f, GEN g, long e)
1204 : {
1205 118508 : if (e==0) return f;
1206 94757 : 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 376396 : Z_to_famat(GEN x)
1215 : {
1216 : long k;
1217 376396 : if (equali1(x)) return trivial_fact();
1218 192586 : k = Z_isanypower(x, &x) ;
1219 192586 : return to_famat_shallow(x, k? utoi(k): gen_1);
1220 : }
1221 : GEN
1222 197168 : Q_to_famat(GEN x)
1223 : {
1224 197168 : if (typ(x) == t_INT) return Z_to_famat(x);
1225 179228 : 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 2688575 : 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 151933 : append(GEN v, GEN x)
1235 : {
1236 151933 : long i, l = lg(v);
1237 151933 : GEN w = cgetg(l+1, typ(v));
1238 653861 : for (i=1; i<l; i++) gel(w,i) = gcopy(gel(v,i));
1239 151933 : gel(w,i) = gcopy(x); return w;
1240 : }
1241 : /* add x^1 to famat f */
1242 : static GEN
1243 158586 : famat_add(GEN f, GEN x)
1244 : {
1245 158586 : GEN h = cgetg(3,t_MAT);
1246 158586 : if (lgcols(f) == 1)
1247 : {
1248 13194 : gel(h,1) = mkcolcopy(x);
1249 13194 : gel(h,2) = mkcol(gen_1);
1250 : }
1251 : else
1252 : {
1253 145392 : gel(h,1) = append(gel(f,1), x);
1254 145392 : gel(h,2) = gconcat(gel(f,2), gen_1);
1255 : }
1256 158586 : return h;
1257 : }
1258 : /* add x^-1 to famat f */
1259 : static GEN
1260 20957 : famat_sub(GEN f, GEN x)
1261 : {
1262 20957 : GEN h = cgetg(3,t_MAT);
1263 20957 : if (lgcols(f) == 1)
1264 : {
1265 14416 : gel(h,1) = mkcolcopy(x);
1266 14416 : 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 20957 : return h;
1274 : }
1275 :
1276 : GEN
1277 447445 : famat_mul(GEN f, GEN g)
1278 : {
1279 : GEN h;
1280 447445 : if (typ(g) != t_MAT) {
1281 30486 : 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 : }
1286 416959 : if (typ(f) != t_MAT) return famat_add(g, f);
1287 288859 : if (lgcols(f) == 1) return gcopy(g);
1288 265545 : if (lgcols(g) == 1) return gcopy(f);
1289 259561 : h = cgetg(3,t_MAT);
1290 259561 : gel(h,1) = gconcat(gel(f,1), gel(g,1));
1291 259561 : gel(h,2) = gconcat(gel(f,2), gel(g,2));
1292 259561 : return h;
1293 : }
1294 :
1295 : GEN
1296 200192 : famat_div(GEN f, GEN g)
1297 : {
1298 : GEN h;
1299 200192 : if (typ(g) != t_MAT) {
1300 20915 : if (typ(f) == t_MAT) return famat_sub(f, g);
1301 0 : h = cgetg(3, t_MAT);
1302 0 : gel(h,1) = mkcol2(gcopy(f), gcopy(g));
1303 0 : gel(h,2) = mkcol2(gen_1, gen_m1);
1304 : }
1305 179277 : if (typ(f) != t_MAT) return famat_sub(g, f);
1306 179235 : if (lgcols(f) == 1) return famat_inv(g);
1307 261 : if (lgcols(g) == 1) return gcopy(f);
1308 261 : h = cgetg(3,t_MAT);
1309 261 : gel(h,1) = gconcat(gel(f,1), gel(g,1));
1310 261 : gel(h,2) = gconcat(gel(f,2), gneg(gel(g,2)));
1311 261 : return h;
1312 : }
1313 :
1314 : GEN
1315 22917 : famat_sqr(GEN f)
1316 : {
1317 : GEN h;
1318 22917 : if (typ(f) != t_MAT) return to_famat(f,gen_2);
1319 22917 : if (lgcols(f) == 1) return gcopy(f);
1320 13127 : h = cgetg(3,t_MAT);
1321 13127 : gel(h,1) = gcopy(gel(f,1));
1322 13127 : gel(h,2) = gmul2n(gel(f,2),1);
1323 13127 : return h;
1324 : }
1325 :
1326 : GEN
1327 26974 : famat_inv_shallow(GEN f)
1328 : {
1329 26974 : if (typ(f) != t_MAT) return to_famat_shallow(f,gen_m1);
1330 10444 : if (lgcols(f) == 1) return f;
1331 10444 : return mkmat2(gel(f,1), ZC_neg(gel(f,2)));
1332 : }
1333 : GEN
1334 199362 : famat_inv(GEN f)
1335 : {
1336 199362 : if (typ(f) != t_MAT) return to_famat(f,gen_m1);
1337 199362 : if (lgcols(f) == 1) return gcopy(f);
1338 180835 : retmkmat2(gcopy(gel(f,1)), ZC_neg(gel(f,2)));
1339 : }
1340 : GEN
1341 60642 : famat_pow(GEN f, GEN n)
1342 : {
1343 60642 : if (typ(f) != t_MAT) return to_famat(f,n);
1344 60642 : if (lgcols(f) == 1) return gcopy(f);
1345 60642 : retmkmat2(gcopy(gel(f,1)), ZC_Z_mul(gel(f,2),n));
1346 : }
1347 : GEN
1348 62010 : famat_pow_shallow(GEN f, GEN n)
1349 : {
1350 62010 : if (is_pm1(n)) return signe(n) > 0? f: famat_inv_shallow(f);
1351 35647 : if (typ(f) != t_MAT) return to_famat_shallow(f,n);
1352 7988 : if (lgcols(f) == 1) return f;
1353 6073 : return mkmat2(gel(f,1), ZC_Z_mul(gel(f,2),n));
1354 : }
1355 :
1356 : GEN
1357 123127 : famat_pows_shallow(GEN f, long n)
1358 : {
1359 123127 : if (n==1) return f;
1360 34956 : if (n==-1) return famat_inv_shallow(f);
1361 34949 : if (typ(f) != t_MAT) return to_famat_shallow(f, stoi(n));
1362 26702 : if (lgcols(f) == 1) return f;
1363 26702 : return mkmat2(gel(f,1), ZC_z_mul(gel(f,2),n));
1364 : }
1365 :
1366 : GEN
1367 0 : famat_Z_gcd(GEN M, GEN n)
1368 : {
1369 0 : pari_sp av=avma;
1370 0 : long i, j, l=lgcols(M);
1371 0 : GEN F=cgetg(3,t_MAT);
1372 0 : gel(F,1)=cgetg(l,t_COL);
1373 0 : gel(F,2)=cgetg(l,t_COL);
1374 0 : for (i=1, j=1; i<l; i++)
1375 : {
1376 0 : GEN p = gcoeff(M,i,1);
1377 0 : GEN e = gminsg(Z_pval(n,p),gcoeff(M,i,2));
1378 0 : if (signe(e))
1379 : {
1380 0 : gcoeff(F,j,1)=p;
1381 0 : gcoeff(F,j,2)=e;
1382 0 : j++;
1383 : }
1384 : }
1385 0 : setlg(gel(F,1),j); setlg(gel(F,2),j);
1386 0 : return gerepilecopy(av,F);
1387 : }
1388 :
1389 : /* x assumed to be a t_MATs (factorization matrix), or compatible with
1390 : * the element_* functions. */
1391 : static GEN
1392 33865 : ext_sqr(GEN nf, GEN x)
1393 33865 : { return (typ(x)==t_MAT)? famat_sqr(x): nfsqr(nf, x); }
1394 : static GEN
1395 60731 : ext_mul(GEN nf, GEN x, GEN y)
1396 60731 : { return (typ(x)==t_MAT)? famat_mul(x,y): nfmul(nf, x, y); }
1397 : static GEN
1398 20388 : ext_inv(GEN nf, GEN x)
1399 20388 : { return (typ(x)==t_MAT)? famat_inv(x): nfinv(nf, x); }
1400 : static GEN
1401 0 : ext_pow(GEN nf, GEN x, GEN n)
1402 0 : { return (typ(x)==t_MAT)? famat_pow(x,n): nfpow(nf, x, n); }
1403 :
1404 : GEN
1405 0 : famat_to_nf(GEN nf, GEN f)
1406 : {
1407 : GEN t, x, e;
1408 : long i;
1409 0 : if (lgcols(f) == 1) return gen_1;
1410 0 : x = gel(f,1);
1411 0 : e = gel(f,2);
1412 0 : t = nfpow(nf, gel(x,1), gel(e,1));
1413 0 : for (i=lg(x)-1; i>1; i--)
1414 0 : t = nfmul(nf, t, nfpow(nf, gel(x,i), gel(e,i)));
1415 0 : return t;
1416 : }
1417 :
1418 : GEN
1419 0 : famat_idealfactor(GEN nf, GEN x)
1420 : {
1421 : long i, l;
1422 0 : GEN g = gel(x,1), e = gel(x,2), h = cgetg_copy(g, &l);
1423 0 : for (i = 1; i < l; i++) gel(h,i) = idealfactor(nf, gel(g,i));
1424 0 : h = famat_reduce(famatV_factorback(h,e));
1425 0 : return sort_factor(h, (void*)&cmp_prime_ideal, &cmp_nodata);
1426 : }
1427 :
1428 : GEN
1429 295800 : famat_reduce(GEN fa)
1430 : {
1431 : GEN E, G, L, g, e;
1432 : long i, k, l;
1433 :
1434 295800 : if (typ(fa) != t_MAT || lgcols(fa) == 1) return fa;
1435 285057 : g = gel(fa,1); l = lg(g);
1436 285057 : e = gel(fa,2);
1437 285057 : L = gen_indexsort(g, (void*)&cmp_universal, &cmp_nodata);
1438 285058 : G = cgetg(l, t_COL);
1439 285058 : E = cgetg(l, t_COL);
1440 : /* merge */
1441 2316023 : for (k=i=1; i<l; i++,k++)
1442 : {
1443 2030966 : gel(G,k) = gel(g,L[i]);
1444 2030966 : gel(E,k) = gel(e,L[i]);
1445 2030966 : if (k > 1 && gidentical(gel(G,k), gel(G,k-1)))
1446 : {
1447 838794 : gel(E,k-1) = addii(gel(E,k), gel(E,k-1));
1448 838796 : k--;
1449 : }
1450 : }
1451 : /* kill 0 exponents */
1452 285057 : l = k;
1453 1477230 : for (k=i=1; i<l; i++)
1454 1192173 : if (!gequal0(gel(E,i)))
1455 : {
1456 1168080 : gel(G,k) = gel(G,i);
1457 1168080 : gel(E,k) = gel(E,i); k++;
1458 : }
1459 285057 : setlg(G, k);
1460 285059 : setlg(E, k); return mkmat2(G,E);
1461 : }
1462 : GEN
1463 63 : matreduce(GEN f)
1464 63 : { pari_sp av = avma;
1465 63 : switch(typ(f))
1466 : {
1467 21 : case t_VEC: case t_COL:
1468 : {
1469 21 : GEN e; f = vec_reduce(f, &e); settyp(f, t_COL);
1470 21 : return gerepilecopy(av, mkmat2(f, zc_to_ZC(e)));
1471 : }
1472 35 : case t_MAT:
1473 35 : if (lg(f) == 3) break;
1474 : default:
1475 14 : pari_err_TYPE("matreduce", f);
1476 : }
1477 28 : if (typ(gel(f,1)) == t_VECSMALL)
1478 0 : f = famatsmall_reduce(f);
1479 : else
1480 : {
1481 28 : if (!RgV_is_ZV(gel(f,2))) pari_err_TYPE("matreduce",f);
1482 21 : f = famat_reduce(f);
1483 : }
1484 21 : return gerepilecopy(av, f);
1485 : }
1486 :
1487 : GEN
1488 176580 : famatsmall_reduce(GEN fa)
1489 : {
1490 : GEN E, G, L, g, e;
1491 : long i, k, l;
1492 176580 : if (lgcols(fa) == 1) return fa;
1493 176580 : g = gel(fa,1); l = lg(g);
1494 176580 : e = gel(fa,2);
1495 176580 : L = vecsmall_indexsort(g);
1496 176581 : G = cgetg(l, t_VECSMALL);
1497 176581 : E = cgetg(l, t_VECSMALL);
1498 : /* merge */
1499 484407 : for (k=i=1; i<l; i++,k++)
1500 : {
1501 307826 : G[k] = g[L[i]];
1502 307826 : E[k] = e[L[i]];
1503 307826 : if (k > 1 && G[k] == G[k-1])
1504 : {
1505 7988 : E[k-1] += E[k];
1506 7988 : k--;
1507 : }
1508 : }
1509 : /* kill 0 exponents */
1510 176581 : l = k;
1511 476419 : for (k=i=1; i<l; i++)
1512 299838 : if (E[i])
1513 : {
1514 296170 : G[k] = G[i];
1515 296170 : E[k] = E[i]; k++;
1516 : }
1517 176581 : setlg(G, k);
1518 176581 : setlg(E, k); return mkmat2(G,E);
1519 : }
1520 :
1521 : GEN
1522 53784 : famat_remove_trivial(GEN fa)
1523 : {
1524 53784 : GEN P, E, p = gel(fa,1), e = gel(fa,2);
1525 53784 : long j, k, l = lg(p);
1526 53784 : P = cgetg(l, t_COL);
1527 53784 : E = cgetg(l, t_COL);
1528 1740245 : for (j = k = 1; j < l; j++)
1529 1686461 : if (signe(gel(e,j))) { gel(P,k) = gel(p,j); gel(E,k++) = gel(e,j); }
1530 53784 : setlg(P, k); setlg(E, k); return mkmat2(P,E);
1531 : }
1532 :
1533 : GEN
1534 13765 : famatV_factorback(GEN v, GEN e)
1535 : {
1536 13765 : long i, l = lg(e);
1537 : GEN V;
1538 13765 : if (l == 1) return trivial_fact();
1539 13380 : V = signe(gel(e,1))? famat_pow_shallow(gel(v,1), gel(e,1)): trivial_fact();
1540 56655 : for (i = 2; i < l; i++) V = famat_mulpow_shallow(V, gel(v,i), gel(e,i));
1541 13380 : return V;
1542 : }
1543 :
1544 : GEN
1545 46724 : famatV_zv_factorback(GEN v, GEN e)
1546 : {
1547 46724 : long i, l = lg(e);
1548 : GEN V;
1549 46724 : if (l == 1) return trivial_fact();
1550 44295 : V = uel(e,1)? famat_pows_shallow(gel(v,1), uel(e,1)): trivial_fact();
1551 139151 : for (i = 2; i < l; i++) V = famat_mulpows_shallow(V, gel(v,i), uel(e,i));
1552 44295 : return V;
1553 : }
1554 :
1555 : GEN
1556 568185 : ZM_famat_limit(GEN fa, GEN limit)
1557 : {
1558 : pari_sp av;
1559 : GEN E, G, g, e, r;
1560 : long i, k, l, n, lG;
1561 :
1562 568185 : if (lgcols(fa) == 1) return fa;
1563 568178 : g = gel(fa,1); l = lg(g);
1564 568178 : e = gel(fa,2);
1565 1137504 : for(n=0, i=1; i<l; i++)
1566 569326 : if (cmpii(gel(g,i),limit)<=0) n++;
1567 568178 : lG = n<l-1 ? n+2 : n+1;
1568 568178 : G = cgetg(lG, t_COL);
1569 568178 : E = cgetg(lG, t_COL);
1570 568179 : av = avma;
1571 1137506 : for (i=1, k=1, r = gen_1; i<l; i++)
1572 : {
1573 569327 : if (cmpii(gel(g,i),limit)<=0)
1574 : {
1575 569194 : gel(G,k) = gel(g,i);
1576 569194 : gel(E,k) = gel(e,i);
1577 569194 : k++;
1578 133 : } else r = mulii(r, powii(gel(g,i), gel(e,i)));
1579 : }
1580 568179 : if (k<i)
1581 : {
1582 133 : gel(G, k) = gerepileuptoint(av, r);
1583 133 : gel(E, k) = gen_1;
1584 : }
1585 568179 : return mkmat2(G,E);
1586 : }
1587 :
1588 : /* assume pr has degree 1 and coprime to Q_denom(x) */
1589 : static GEN
1590 123270 : to_Fp_coprime(GEN nf, GEN x, GEN modpr)
1591 : {
1592 123270 : GEN d, r, p = modpr_get_p(modpr);
1593 123270 : x = nf_to_scalar_or_basis(nf,x);
1594 123270 : if (typ(x) != t_COL) return Rg_to_Fp(x,p);
1595 121632 : x = Q_remove_denom(x, &d);
1596 121632 : r = zk_to_Fq(x, modpr);
1597 121632 : if (d) r = Fp_div(r, d, p);
1598 121632 : return r;
1599 : }
1600 :
1601 : /* pr coprime to all denominators occurring in x */
1602 : static GEN
1603 693 : famat_to_Fp_coprime(GEN nf, GEN x, GEN modpr)
1604 : {
1605 693 : GEN p = modpr_get_p(modpr);
1606 693 : GEN t = NULL, g = gel(x,1), e = gel(x,2), q = subiu(p,1);
1607 693 : long i, l = lg(g);
1608 4571 : for (i = 1; i < l; i++)
1609 : {
1610 3878 : GEN n = modii(gel(e,i), q);
1611 3878 : if (signe(n))
1612 : {
1613 3871 : GEN h = to_Fp_coprime(nf, gel(g,i), modpr);
1614 3871 : h = Fp_pow(h, n, p);
1615 3871 : t = t? Fp_mul(t, h, p): h;
1616 : }
1617 : }
1618 693 : return t? modii(t, p): gen_1;
1619 : }
1620 :
1621 : /* cf famat_to_nf_modideal_coprime, modpr attached to prime of degree 1 */
1622 : GEN
1623 120092 : nf_to_Fp_coprime(GEN nf, GEN x, GEN modpr)
1624 : {
1625 693 : return typ(x)==t_MAT? famat_to_Fp_coprime(nf, x, modpr)
1626 120785 : : to_Fp_coprime(nf, x, modpr);
1627 : }
1628 :
1629 : static long
1630 4098985 : zk_pvalrem(GEN x, GEN p, GEN *py)
1631 4098985 : { return (typ(x) == t_INT)? Z_pvalrem(x, p, py): ZV_pvalrem(x, p, py); }
1632 : /* x a QC or Q. Return a ZC or Z, whose content is coprime to Z. Set v, dx
1633 : * such that x = p^v (newx / dx); dx = NULL if 1 */
1634 : static GEN
1635 4144791 : nf_remove_denom_p(GEN nf, GEN x, GEN p, GEN *pdx, long *pv)
1636 : {
1637 : long vcx;
1638 : GEN dx;
1639 4144791 : x = nf_to_scalar_or_basis(nf, x);
1640 4144794 : x = Q_remove_denom(x, &dx);
1641 4144785 : if (dx)
1642 : {
1643 61197 : vcx = - Z_pvalrem(dx, p, &dx);
1644 61197 : if (!vcx) vcx = zk_pvalrem(x, p, &x);
1645 61197 : if (isint1(dx)) dx = NULL;
1646 : }
1647 : else
1648 : {
1649 4083588 : vcx = zk_pvalrem(x, p, &x);
1650 4083591 : dx = NULL;
1651 : }
1652 4144788 : *pv = vcx;
1653 4144788 : *pdx = dx; return x;
1654 : }
1655 : /* x = b^e/p^(e-1) in Z_K; x = 0 mod p/pr^e, (x,pr) = 1. Return NULL
1656 : * if p inert (instead of 1) */
1657 : static GEN
1658 94947 : p_makecoprime(GEN pr)
1659 : {
1660 94947 : GEN B = pr_get_tau(pr), b;
1661 : long i, e;
1662 :
1663 94947 : if (typ(B) == t_INT) return NULL;
1664 72190 : b = gel(B,1); /* B = multiplication table by b */
1665 72190 : e = pr_get_e(pr);
1666 72190 : if (e == 1) return b;
1667 : /* one could also divide (exactly) by p in each iteration */
1668 45305 : for (i = 1; i < e; i++) b = ZM_ZC_mul(B, b);
1669 22139 : return ZC_Z_divexact(b, powiu(pr_get_p(pr), e-1));
1670 : }
1671 :
1672 : /* Compute A = prod g[i]^e[i] mod pr^k, assuming (A, pr) = 1.
1673 : * Method: modify each g[i] so that it becomes coprime to pr,
1674 : * g[i] *= (b/p)^v_pr(g[i]), where b/p = pr^(-1) times something integral
1675 : * and prime to p; globally, we multiply by (b/p)^v_pr(A) = 1.
1676 : * Optimizations:
1677 : * 1) remove all powers of p from contents, and consider extra generator p^vp;
1678 : * modified as p * (b/p)^e = b^e / p^(e-1)
1679 : * 2) remove denominators, coprime to p, by multiplying by inverse mod prk\cap Z
1680 : *
1681 : * EX = multiple of exponent of (O_K / pr^k)^* used to reduce the product in
1682 : * case the e[i] are large */
1683 : GEN
1684 2003774 : famat_makecoprime(GEN nf, GEN g, GEN e, GEN pr, GEN prk, GEN EX)
1685 : {
1686 2003774 : GEN G, E, t, vp = NULL, p = pr_get_p(pr), prkZ = gcoeff(prk, 1,1);
1687 2003776 : long i, l = lg(g);
1688 :
1689 2003776 : G = cgetg(l+1, t_VEC);
1690 2003802 : E = cgetg(l+1, t_VEC); /* l+1: room for "modified p" */
1691 6148561 : for (i=1; i < l; i++)
1692 : {
1693 : long vcx;
1694 4144785 : GEN dx, x = nf_remove_denom_p(nf, gel(g,i), p, &dx, &vcx);
1695 4144782 : if (vcx) /* = v_p(content(g[i])) */
1696 : {
1697 142122 : GEN a = mulsi(vcx, gel(e,i));
1698 142133 : vp = vp? addii(vp, a): a;
1699 : }
1700 : /* x integral, content coprime to p; dx coprime to p */
1701 4144794 : if (typ(x) == t_INT)
1702 : { /* x coprime to p, hence to pr */
1703 1129145 : x = modii(x, prkZ);
1704 1129145 : if (dx) x = Fp_div(x, dx, prkZ);
1705 : }
1706 : else
1707 : {
1708 3015649 : (void)ZC_nfvalrem(x, pr, &x); /* x *= (b/p)^v_pr(x) */
1709 3015623 : x = ZC_hnfrem(FpC_red(x,prkZ), prk);
1710 3015612 : if (dx) x = FpC_Fp_mul(x, Fp_inv(dx,prkZ), prkZ);
1711 : }
1712 4144742 : gel(G,i) = x;
1713 4144742 : gel(E,i) = gel(e,i);
1714 : }
1715 :
1716 2003776 : t = vp? p_makecoprime(pr): NULL;
1717 2003788 : if (!t)
1718 : { /* no need for extra generator */
1719 1931675 : setlg(G,l);
1720 1931677 : setlg(E,l);
1721 : }
1722 : else
1723 : {
1724 72113 : gel(G,i) = FpC_red(t, prkZ);
1725 72113 : gel(E,i) = vp;
1726 : }
1727 2003794 : return famat_to_nf_modideal_coprime(nf, G, E, prk, EX);
1728 : }
1729 :
1730 : /* simplified version of famat_makecoprime for X = SUnits[1] */
1731 : GEN
1732 98 : sunits_makecoprime(GEN X, GEN pr, GEN prk)
1733 : {
1734 98 : GEN G, p = pr_get_p(pr), prkZ = gcoeff(prk,1,1);
1735 98 : long i, l = lg(X);
1736 :
1737 98 : G = cgetg(l, t_VEC);
1738 9205 : for (i = 1; i < l; i++)
1739 : {
1740 9107 : GEN x = gel(X,i);
1741 9107 : if (typ(x) == t_INT) /* a prime */
1742 1491 : x = equalii(x,p)? p_makecoprime(pr): modii(x, prkZ);
1743 : else
1744 : {
1745 7616 : (void)ZC_nfvalrem(x, pr, &x); /* x *= (b/p)^v_pr(x) */
1746 7616 : x = ZC_hnfrem(FpC_red(x,prkZ), prk);
1747 : }
1748 9107 : gel(G,i) = x;
1749 : }
1750 98 : return G;
1751 : }
1752 :
1753 : /* prod g[i]^e[i] mod bid, assume (g[i], id) = 1 and 1 < lg(g) <= lg(e) */
1754 : GEN
1755 20118 : famat_to_nf_moddivisor(GEN nf, GEN g, GEN e, GEN bid)
1756 : {
1757 20118 : GEN t, cyc = bid_get_cyc(bid);
1758 20118 : if (lg(cyc) == 1)
1759 0 : t = gen_1;
1760 : else
1761 20118 : t = famat_to_nf_modideal_coprime(nf, g, e, bid_get_ideal(bid),
1762 : cyc_get_expo(cyc));
1763 20118 : return set_sign_mod_divisor(nf, mkmat2(g,e), t, bid_get_sarch(bid));
1764 : }
1765 :
1766 : GEN
1767 15996217 : vecmul(GEN x, GEN y)
1768 : {
1769 15996217 : if (!is_vec_t(typ(x))) return gmul(x,y);
1770 3524274 : pari_APPLY_same(vecmul(gel(x,i), gel(y,i)))
1771 : }
1772 :
1773 : GEN
1774 185983 : vecsqr(GEN x)
1775 : {
1776 185983 : if (!is_vec_t(typ(x))) return gsqr(x);
1777 46606 : pari_APPLY_same(vecsqr(gel(x,i)))
1778 : }
1779 :
1780 : GEN
1781 826 : vecinv(GEN x)
1782 : {
1783 826 : if (!is_vec_t(typ(x))) return ginv(x);
1784 56 : pari_APPLY_same(vecinv(gel(x,i)))
1785 : }
1786 :
1787 : GEN
1788 0 : vecpow(GEN x, GEN n)
1789 : {
1790 0 : if (!is_vec_t(typ(x))) return powgi(x,n);
1791 0 : pari_APPLY_same(vecpow(gel(x,i), n))
1792 : }
1793 :
1794 : GEN
1795 903 : vecdiv(GEN x, GEN y)
1796 : {
1797 903 : if (!is_vec_t(typ(x))) return gdiv(x,y);
1798 903 : pari_APPLY_same(vecdiv(gel(x,i), gel(y,i)))
1799 : }
1800 :
1801 : /* A ideal as a square t_MAT */
1802 : static GEN
1803 306076 : idealmulelt(GEN nf, GEN x, GEN A)
1804 : {
1805 : long i, lx;
1806 : GEN dx, dA, D;
1807 306076 : if (lg(A) == 1) return cgetg(1, t_MAT);
1808 306076 : x = nf_to_scalar_or_basis(nf,x);
1809 306076 : if (typ(x) != t_COL)
1810 100355 : return isintzero(x)? cgetg(1,t_MAT): RgM_Rg_mul(A, Q_abs_shallow(x));
1811 205721 : x = Q_remove_denom(x, &dx);
1812 205721 : A = Q_remove_denom(A, &dA);
1813 205721 : x = zk_multable(nf, x);
1814 205721 : D = mulii(zkmultable_capZ(x), gcoeff(A,1,1));
1815 205721 : x = zkC_multable_mul(A, x);
1816 205721 : settyp(x, t_MAT); lx = lg(x);
1817 : /* x may contain scalars (at most 1 since the ideal is nonzero)*/
1818 783422 : for (i=1; i<lx; i++)
1819 592616 : if (typ(gel(x,i)) == t_INT)
1820 : {
1821 14915 : if (i > 1) swap(gel(x,1), gel(x,i)); /* help HNF */
1822 14915 : gel(x,1) = scalarcol_shallow(gel(x,1), lx-1);
1823 14915 : break;
1824 : }
1825 205721 : x = ZM_hnfmodid(x, D);
1826 205721 : dx = mul_denom(dx,dA);
1827 205721 : return dx? gdiv(x,dx): x;
1828 : }
1829 :
1830 : /* nf a true nf, tx <= ty */
1831 : static GEN
1832 579271 : idealmul_aux(GEN nf, GEN x, GEN y, long tx, long ty)
1833 : {
1834 : GEN z, cx, cy;
1835 579271 : switch(tx)
1836 : {
1837 363891 : case id_PRINCIPAL:
1838 363891 : switch(ty)
1839 : {
1840 57388 : case id_PRINCIPAL:
1841 57388 : return idealhnf_principal(nf, nfmul(nf,x,y));
1842 427 : case id_PRIME:
1843 : {
1844 427 : GEN p = pr_get_p(y), pi = pr_get_gen(y), cx;
1845 427 : if (pr_is_inert(y)) return RgM_Rg_mul(idealhnf_principal(nf,x),p);
1846 :
1847 217 : x = nf_to_scalar_or_basis(nf, x);
1848 217 : switch(typ(x))
1849 : {
1850 203 : case t_INT:
1851 203 : if (!signe(x)) return cgetg(1,t_MAT);
1852 203 : return ZM_Z_mul(pr_hnf(nf,y), absi_shallow(x));
1853 7 : case t_FRAC:
1854 7 : return RgM_Rg_mul(pr_hnf(nf,y), Q_abs_shallow(x));
1855 : }
1856 : /* t_COL */
1857 7 : x = Q_primitive_part(x, &cx);
1858 7 : x = zk_multable(nf, x);
1859 7 : z = shallowconcat(ZM_Z_mul(x,p), ZM_ZC_mul(x,pi));
1860 7 : z = ZM_hnfmodid(z, mulii(p, zkmultable_capZ(x)));
1861 7 : return cx? ZM_Q_mul(z, cx): z;
1862 : }
1863 306076 : default: /* id_MAT */
1864 306076 : return idealmulelt(nf, x,y);
1865 : }
1866 42119 : case id_PRIME:
1867 42119 : if (ty==id_PRIME)
1868 4347 : { y = pr_hnf(nf,y); cy = NULL; }
1869 : else
1870 37772 : y = Q_primitive_part(y, &cy);
1871 42119 : y = idealHNF_mul_two(nf,y,x);
1872 42119 : return cy? ZM_Q_mul(y,cy): y;
1873 :
1874 173261 : default: /* id_MAT */
1875 : {
1876 173261 : long N = nf_get_degree(nf);
1877 173261 : if (lg(x)-1 != N || lg(y)-1 != N) pari_err_DIM("idealmul");
1878 173247 : x = Q_primitive_part(x, &cx);
1879 173247 : y = Q_primitive_part(y, &cy); cx = mul_content(cx,cy);
1880 173247 : y = idealHNF_mul(nf,x,y);
1881 173247 : return cx? ZM_Q_mul(y,cx): y;
1882 : }
1883 : }
1884 : }
1885 :
1886 : /* output the ideal product x.y */
1887 : GEN
1888 579270 : idealmul(GEN nf, GEN x, GEN y)
1889 : {
1890 : pari_sp av;
1891 : GEN res, ax, ay, z;
1892 579270 : long tx = idealtyp(&x,&ax);
1893 579270 : long ty = idealtyp(&y,&ay), f;
1894 579270 : if (tx>ty) { swap(ax,ay); swap(x,y); lswap(tx,ty); }
1895 579270 : f = (ax||ay); res = f? cgetg(3,t_VEC): NULL; /*product is an extended ideal*/
1896 579270 : av = avma;
1897 579270 : z = gerepileupto(av, idealmul_aux(checknf(nf), x,y, tx,ty));
1898 579257 : if (!f) return z;
1899 27997 : if (ax && ay)
1900 26513 : ax = ext_mul(nf, ax, ay);
1901 : else
1902 1484 : ax = gcopy(ax? ax: ay);
1903 27997 : gel(res,1) = z; gel(res,2) = ax; return res;
1904 : }
1905 :
1906 : /* Return x, integral in 2-elt form, such that pr^2 = c * x. cf idealpowprime
1907 : * nf = true nf */
1908 : static GEN
1909 285450 : idealsqrprime(GEN nf, GEN pr, GEN *pc)
1910 : {
1911 285450 : GEN p = pr_get_p(pr), q, gen;
1912 285450 : long e = pr_get_e(pr), f = pr_get_f(pr);
1913 :
1914 285453 : q = (e == 1)? sqri(p): p;
1915 285446 : if (e <= 2 && e * f == nf_get_degree(nf))
1916 : { /* pr^e = (p) */
1917 45507 : *pc = q;
1918 45507 : return mkvec2(gen_1,gen_0);
1919 : }
1920 239937 : gen = nfsqr(nf, pr_get_gen(pr));
1921 239947 : gen = FpC_red(gen, q);
1922 239933 : *pc = NULL;
1923 239933 : return mkvec2(q, gen);
1924 : }
1925 : /* cf idealpow_aux */
1926 : static GEN
1927 38534 : idealsqr_aux(GEN nf, GEN x, long tx)
1928 : {
1929 38534 : GEN T = nf_get_pol(nf), m, cx, a, alpha;
1930 38534 : long N = degpol(T);
1931 38534 : switch(tx)
1932 : {
1933 84 : case id_PRINCIPAL:
1934 84 : return idealhnf_principal(nf, nfsqr(nf,x));
1935 10758 : case id_PRIME:
1936 10758 : if (pr_is_inert(x)) return scalarmat(sqri(gel(x,1)), N);
1937 10590 : x = idealsqrprime(nf, x, &cx);
1938 10590 : x = idealhnf_two(nf,x);
1939 10590 : return cx? ZM_Z_mul(x, cx): x;
1940 27692 : default:
1941 27692 : x = Q_primitive_part(x, &cx);
1942 27692 : a = mat_ideal_two_elt(nf,x); alpha = gel(a,2); a = gel(a,1);
1943 27692 : alpha = nfsqr(nf,alpha);
1944 27692 : m = zk_scalar_or_multable(nf, alpha);
1945 27692 : if (typ(m) == t_INT) {
1946 1635 : x = gcdii(sqri(a), m);
1947 1635 : if (cx) x = gmul(x, gsqr(cx));
1948 1635 : x = scalarmat(x, N);
1949 : }
1950 : else
1951 : { /* could use gcdii(sqri(a), zkmultable_capZ(m)), but costly */
1952 26057 : x = ZM_hnfmodid(m, sqri(a));
1953 26057 : if (cx) cx = gsqr(cx);
1954 26057 : if (cx) x = ZM_Q_mul(x, cx);
1955 : }
1956 27692 : return x;
1957 : }
1958 : }
1959 : GEN
1960 38534 : idealsqr(GEN nf, GEN x)
1961 : {
1962 : pari_sp av;
1963 : GEN res, ax, z;
1964 38534 : long tx = idealtyp(&x,&ax);
1965 38534 : res = ax? cgetg(3,t_VEC): NULL; /*product is an extended ideal*/
1966 38534 : av = avma;
1967 38534 : z = gerepileupto(av, idealsqr_aux(checknf(nf), x, tx));
1968 38533 : if (!ax) return z;
1969 33865 : gel(res,1) = z;
1970 33865 : gel(res,2) = ext_sqr(nf, ax); return res;
1971 : }
1972 :
1973 : /* norm of an ideal */
1974 : GEN
1975 104323 : idealnorm(GEN nf, GEN x)
1976 : {
1977 : pari_sp av;
1978 : long tx;
1979 :
1980 104323 : switch(idealtyp(&x, NULL))
1981 : {
1982 4865 : case id_PRIME: return pr_norm(x);
1983 10306 : case id_MAT: return RgM_det_triangular(x);
1984 : }
1985 : /* id_PRINCIPAL */
1986 89152 : nf = checknf(nf); av = avma;
1987 89152 : x = nfnorm(nf, x);
1988 89152 : tx = typ(x);
1989 89152 : if (tx == t_INT) return gerepileuptoint(av, absi(x));
1990 406 : if (tx != t_FRAC) pari_err_TYPE("idealnorm",x);
1991 406 : return gerepileupto(av, Q_abs(x));
1992 : }
1993 :
1994 : /* x \cap Z */
1995 : GEN
1996 3031 : idealdown(GEN nf, GEN x)
1997 : {
1998 3031 : pari_sp av = avma;
1999 : GEN y, c;
2000 3031 : switch(idealtyp(&x, NULL))
2001 : {
2002 7 : case id_PRIME: return icopy(pr_get_p(x));
2003 2135 : case id_MAT: return gcopy(gcoeff(x,1,1));
2004 : }
2005 : /* id_PRINCIPAL */
2006 889 : nf = checknf(nf); av = avma;
2007 889 : x = nf_to_scalar_or_basis(nf, x);
2008 889 : if (is_rational_t(typ(x))) return Q_abs(x);
2009 14 : x = Q_primitive_part(x, &c);
2010 14 : y = zkmultable_capZ(zk_multable(nf, x));
2011 14 : return gerepilecopy(av, mul_content(c, y));
2012 : }
2013 :
2014 : /* true nf */
2015 : static GEN
2016 42 : idealismaximal_int(GEN nf, GEN p)
2017 : {
2018 : GEN L;
2019 42 : if (!BPSW_psp(p)) return NULL;
2020 77 : if (!dvdii(nf_get_index(nf), p) &&
2021 49 : !FpX_is_irred(FpX_red(nf_get_pol(nf),p), p)) return NULL;
2022 28 : L = idealprimedec(nf, p);
2023 28 : return (lg(L) == 2 && pr_get_e(gel(L,1)) == 1)? gel(L,1): NULL;
2024 : }
2025 : /* true nf */
2026 : static GEN
2027 21 : idealismaximal_mat(GEN nf, GEN x)
2028 : {
2029 : GEN p, c, L;
2030 : long i, l, f;
2031 21 : x = Q_primitive_part(x, &c);
2032 21 : p = gcoeff(x,1,1);
2033 21 : if (c)
2034 : {
2035 7 : if (typ(c) == t_FRAC || !equali1(p)) return NULL;
2036 7 : return idealismaximal_int(nf, c);
2037 : }
2038 14 : if (!BPSW_psp(p)) return NULL;
2039 14 : l = lg(x); f = 1;
2040 35 : for (i = 2; i < l; i++)
2041 : {
2042 21 : c = gcoeff(x,i,i);
2043 21 : if (equalii(c, p)) f++; else if (!equali1(c)) return NULL;
2044 : }
2045 14 : L = idealprimedec_limit_f(nf, p, f);
2046 28 : for (i = lg(L)-1; i; i--)
2047 : {
2048 28 : GEN pr = gel(L,i);
2049 28 : if (pr_get_f(pr) != f) break;
2050 28 : if (idealval(nf, x, pr) == 1) return pr;
2051 : }
2052 0 : return NULL;
2053 : }
2054 : /* true nf */
2055 : static GEN
2056 77 : idealismaximal_i(GEN nf, GEN x)
2057 : {
2058 : GEN L, p, pr, c;
2059 : long i, l;
2060 77 : switch(idealtyp(&x, NULL))
2061 : {
2062 7 : case id_PRIME: return x;
2063 21 : case id_MAT: return idealismaximal_mat(nf, x);
2064 : }
2065 : /* id_PRINCIPAL */
2066 49 : x = nf_to_scalar_or_basis(nf, x);
2067 49 : switch(typ(x))
2068 : {
2069 35 : case t_INT: return idealismaximal_int(nf, absi_shallow(x));
2070 0 : case t_FRAC: return NULL;
2071 : }
2072 14 : x = Q_primitive_part(x, &c);
2073 14 : if (c) return NULL;
2074 14 : p = zkmultable_capZ(zk_multable(nf, x));
2075 14 : if (!BPSW_psp(p)) return NULL;
2076 7 : L = idealprimedec(nf, p); l = lg(L); pr = NULL;
2077 21 : for (i = 1; i < l; i++)
2078 : {
2079 14 : long v = ZC_nfval(x, gel(L,i));
2080 14 : if (v > 1 || (v && pr)) return NULL;
2081 14 : pr = gel(L,i);
2082 : }
2083 7 : return pr;
2084 : }
2085 : GEN
2086 77 : idealismaximal(GEN nf, GEN x)
2087 : {
2088 77 : pari_sp av = avma;
2089 77 : x = idealismaximal_i(checknf(nf), x);
2090 77 : if (!x) { set_avma(av); return gen_0; }
2091 49 : return gerepilecopy(av, x);
2092 : }
2093 :
2094 : /* I^(-1) = { x \in K, Tr(x D^(-1) I) \in Z }, D different of K/Q
2095 : *
2096 : * nf[5][6] = pp( D^(-1) ) = pp( HNF( T^(-1) ) ), T = (Tr(wi wj))
2097 : * nf[5][7] = same in 2-elt form.
2098 : * Assume I integral. Return the integral ideal (I\cap Z) I^(-1) */
2099 : GEN
2100 223258 : idealHNF_inv_Z(GEN nf, GEN I)
2101 : {
2102 223258 : GEN J, dual, IZ = gcoeff(I,1,1); /* I \cap Z */
2103 223258 : if (isint1(IZ)) return matid(lg(I)-1);
2104 201727 : J = idealHNF_mul(nf,I, gmael(nf,5,7));
2105 : /* I in HNF, hence easily inverted; multiply by IZ to get integer coeffs
2106 : * missing content cancels while solving the linear equation */
2107 201728 : dual = shallowtrans( hnf_divscale(J, gmael(nf,5,6), IZ) );
2108 201728 : return ZM_hnfmodid(dual, IZ);
2109 : }
2110 : /* I HNF with rational coefficients (denominator d). */
2111 : GEN
2112 81175 : idealHNF_inv(GEN nf, GEN I)
2113 : {
2114 81175 : GEN J, IQ = gcoeff(I,1,1); /* I \cap Q; d IQ = dI \cap Z */
2115 81175 : J = idealHNF_inv_Z(nf, Q_remove_denom(I, NULL)); /* = (dI)^(-1) * (d IQ) */
2116 81175 : return equali1(IQ)? J: RgM_Rg_div(J, IQ);
2117 : }
2118 :
2119 : /* return p * P^(-1) [integral] */
2120 : GEN
2121 39501 : pr_inv_p(GEN pr)
2122 : {
2123 39501 : if (pr_is_inert(pr)) return matid(pr_get_f(pr));
2124 38689 : return ZM_hnfmodid(pr_get_tau(pr), pr_get_p(pr));
2125 : }
2126 : GEN
2127 18376 : pr_inv(GEN pr)
2128 : {
2129 18376 : GEN p = pr_get_p(pr);
2130 18376 : if (pr_is_inert(pr)) return scalarmat(ginv(p), pr_get_f(pr));
2131 17963 : return RgM_Rg_div(ZM_hnfmodid(pr_get_tau(pr),p), p);
2132 : }
2133 :
2134 : GEN
2135 152321 : idealinv(GEN nf, GEN x)
2136 : {
2137 : GEN res, ax;
2138 : pari_sp av;
2139 152321 : long tx = idealtyp(&x,&ax), N;
2140 :
2141 152321 : res = ax? cgetg(3,t_VEC): NULL;
2142 152321 : nf = checknf(nf); av = avma;
2143 152321 : N = nf_get_degree(nf);
2144 152321 : switch (tx)
2145 : {
2146 74030 : case id_MAT:
2147 74030 : if (lg(x)-1 != N) pari_err_DIM("idealinv");
2148 74030 : x = idealHNF_inv(nf,x); break;
2149 61279 : case id_PRINCIPAL:
2150 61279 : x = nf_to_scalar_or_basis(nf, x);
2151 61279 : if (typ(x) != t_COL)
2152 61230 : x = idealhnf_principal(nf,ginv(x));
2153 : else
2154 : { /* nfinv + idealhnf where we already know (x) \cap Z */
2155 : GEN c, d;
2156 49 : x = Q_remove_denom(x, &c);
2157 49 : x = zk_inv(nf, x);
2158 49 : x = Q_remove_denom(x, &d); /* true inverse is c/d * x */
2159 49 : if (!d) /* x and x^(-1) integral => x a unit */
2160 14 : x = c? scalarmat(c, N): matid(N);
2161 : else
2162 : {
2163 35 : c = c? gdiv(c,d): ginv(d);
2164 35 : x = zk_multable(nf, x);
2165 35 : x = ZM_Q_mul(ZM_hnfmodid(x,d), c);
2166 : }
2167 : }
2168 61279 : break;
2169 17012 : case id_PRIME:
2170 17012 : x = pr_inv(x); break;
2171 : }
2172 152321 : x = gerepileupto(av,x); if (!ax) return x;
2173 20388 : gel(res,1) = x;
2174 20388 : gel(res,2) = ext_inv(nf, ax); return res;
2175 : }
2176 :
2177 : /* write x = A/B, A,B coprime integral ideals */
2178 : GEN
2179 379460 : idealnumden(GEN nf, GEN x)
2180 : {
2181 379460 : pari_sp av = avma;
2182 : GEN x0, c, d, A, B, J;
2183 379460 : long tx = idealtyp(&x, NULL);
2184 379459 : nf = checknf(nf);
2185 379463 : switch (tx)
2186 : {
2187 7 : case id_PRIME:
2188 7 : retmkvec2(idealhnf(nf, x), gen_1);
2189 138194 : case id_PRINCIPAL:
2190 : {
2191 : GEN xZ, mx;
2192 138194 : x = nf_to_scalar_or_basis(nf, x);
2193 138194 : switch(typ(x))
2194 : {
2195 86457 : case t_INT: return gerepilecopy(av, mkvec2(absi_shallow(x),gen_1));
2196 2639 : case t_FRAC:return gerepilecopy(av, mkvec2(absi_shallow(gel(x,1)), gel(x,2)));
2197 : }
2198 : /* t_COL */
2199 49098 : x = Q_remove_denom(x, &d);
2200 49098 : if (!d) return gerepilecopy(av, mkvec2(idealhnf_shallow(nf, x), gen_1));
2201 105 : mx = zk_multable(nf, x);
2202 105 : xZ = zkmultable_capZ(mx);
2203 105 : x = ZM_hnfmodid(mx, xZ); /* principal ideal (x) */
2204 105 : x0 = mkvec2(xZ, mx); /* same, for fast multiplication */
2205 105 : break;
2206 : }
2207 241262 : default: /* id_MAT */
2208 : {
2209 241262 : long n = lg(x)-1;
2210 241262 : if (n == 0) return mkvec2(gen_0, gen_1);
2211 241262 : if (n != nf_get_degree(nf)) pari_err_DIM("idealnumden");
2212 241261 : x0 = x = Q_remove_denom(x, &d);
2213 241261 : if (!d) return gerepilecopy(av, mkvec2(x, gen_1));
2214 21 : break;
2215 : }
2216 : }
2217 126 : J = hnfmodid(x, d); /* = d/B */
2218 126 : c = gcoeff(J,1,1); /* (d/B) \cap Z, divides d */
2219 126 : B = idealHNF_inv_Z(nf, J); /* (d/B \cap Z) B/d */
2220 126 : if (!equalii(c,d)) B = ZM_Z_mul(B, diviiexact(d,c)); /* = B ! */
2221 126 : A = idealHNF_mul(nf, B, x0); /* d * (original x) * B = d A */
2222 126 : A = ZM_Z_divexact(A, d); /* = A ! */
2223 126 : return gerepilecopy(av, mkvec2(A, B));
2224 : }
2225 :
2226 : /* Return x, integral in 2-elt form, such that pr^n = c * x. Assume n != 0.
2227 : * nf = true nf */
2228 : static GEN
2229 1198065 : idealpowprime(GEN nf, GEN pr, GEN n, GEN *pc)
2230 : {
2231 1198065 : GEN p = pr_get_p(pr), q, gen;
2232 :
2233 1198051 : *pc = NULL;
2234 1198051 : if (is_pm1(n)) /* n = 1 special cased for efficiency */
2235 : {
2236 604835 : q = p;
2237 604835 : if (typ(pr_get_tau(pr)) == t_INT) /* inert */
2238 : {
2239 0 : *pc = (signe(n) >= 0)? p: ginv(p);
2240 0 : return mkvec2(gen_1,gen_0);
2241 : }
2242 604829 : if (signe(n) >= 0) gen = pr_get_gen(pr);
2243 : else
2244 : {
2245 156432 : gen = pr_get_tau(pr); /* possibly t_MAT */
2246 156447 : *pc = ginv(p);
2247 : }
2248 : }
2249 593297 : else if (equalis(n,2)) return idealsqrprime(nf, pr, pc);
2250 : else
2251 : {
2252 318436 : long e = pr_get_e(pr), f = pr_get_f(pr);
2253 318440 : GEN r, m = truedvmdis(n, e, &r);
2254 318424 : if (e * f == nf_get_degree(nf))
2255 : { /* pr^e = (p) */
2256 76528 : if (signe(m)) *pc = powii(p,m);
2257 76527 : if (!signe(r)) return mkvec2(gen_1,gen_0);
2258 36469 : q = p;
2259 36469 : gen = nfpow(nf, pr_get_gen(pr), r);
2260 : }
2261 : else
2262 : {
2263 241909 : m = absi_shallow(m);
2264 241913 : if (signe(r)) m = addiu(m,1);
2265 241913 : q = powii(p,m); /* m = ceil(|n|/e) */
2266 241914 : if (signe(n) >= 0) gen = nfpow(nf, pr_get_gen(pr), n);
2267 : else
2268 : {
2269 24418 : gen = pr_get_tau(pr);
2270 24418 : if (typ(gen) == t_MAT) gen = gel(gen,1);
2271 24418 : n = negi(n);
2272 24418 : gen = ZC_Z_divexact(nfpow(nf, gen, n), powii(p, subii(n,m)));
2273 24418 : *pc = ginv(q);
2274 : }
2275 : }
2276 278387 : gen = FpC_red(gen, q);
2277 : }
2278 883196 : return mkvec2(q, gen);
2279 : }
2280 :
2281 : /* True nf. x * pr^n. Assume x in HNF or scalar (possibly nonintegral) */
2282 : GEN
2283 770138 : idealmulpowprime(GEN nf, GEN x, GEN pr, GEN n)
2284 : {
2285 : GEN c, cx, y;
2286 770138 : long N = nf_get_degree(nf);
2287 :
2288 770143 : if (!signe(n)) return typ(x) == t_MAT? x: scalarmat_shallow(x, N);
2289 :
2290 : /* inert, special cased for efficiency */
2291 770136 : if (pr_is_inert(pr))
2292 : {
2293 76472 : GEN q = powii(pr_get_p(pr), n);
2294 74129 : return typ(x) == t_MAT? RgM_Rg_mul(x,q)
2295 150605 : : scalarmat_shallow(gmul(Q_abs(x),q), N);
2296 : }
2297 :
2298 693690 : y = idealpowprime(nf, pr, n, &c);
2299 693716 : if (typ(x) == t_MAT)
2300 690344 : { x = Q_primitive_part(x, &cx); if (is_pm1(gcoeff(x,1,1))) x = NULL; }
2301 : else
2302 3372 : { cx = x; x = NULL; }
2303 693442 : cx = mul_content(c,cx);
2304 693464 : if (x)
2305 526700 : x = idealHNF_mul_two(nf,x,y);
2306 : else
2307 166764 : x = idealhnf_two(nf,y);
2308 693912 : if (cx) x = ZM_Q_mul(x,cx);
2309 693513 : return x;
2310 : }
2311 : GEN
2312 11537 : idealdivpowprime(GEN nf, GEN x, GEN pr, GEN n)
2313 : {
2314 11537 : return idealmulpowprime(nf,x,pr, negi(n));
2315 : }
2316 :
2317 : /* nf = true nf */
2318 : static GEN
2319 887706 : idealpow_aux(GEN nf, GEN x, long tx, GEN n)
2320 : {
2321 887706 : GEN T = nf_get_pol(nf), m, cx, n1, a, alpha;
2322 887704 : long N = degpol(T), s = signe(n);
2323 887705 : if (!s) return matid(N);
2324 871498 : switch(tx)
2325 : {
2326 75528 : case id_PRINCIPAL:
2327 75528 : return idealhnf_principal(nf, nfpow(nf,x,n));
2328 600975 : case id_PRIME:
2329 600975 : if (pr_is_inert(x)) return scalarmat(powii(gel(x,1), n), N);
2330 504326 : x = idealpowprime(nf, x, n, &cx);
2331 504320 : x = idealhnf_two(nf,x);
2332 504334 : return cx? ZM_Q_mul(x, cx): x;
2333 194995 : default:
2334 194995 : if (is_pm1(n)) return (s < 0)? idealinv(nf, x): gcopy(x);
2335 69148 : n1 = (s < 0)? negi(n): n;
2336 :
2337 69148 : x = Q_primitive_part(x, &cx);
2338 69148 : a = mat_ideal_two_elt(nf,x); alpha = gel(a,2); a = gel(a,1);
2339 69148 : alpha = nfpow(nf,alpha,n1);
2340 69148 : m = zk_scalar_or_multable(nf, alpha);
2341 69148 : if (typ(m) == t_INT) {
2342 553 : x = gcdii(powii(a,n1), m);
2343 553 : if (s<0) x = ginv(x);
2344 553 : if (cx) x = gmul(x, powgi(cx,n));
2345 553 : x = scalarmat(x, N);
2346 : }
2347 : else
2348 : { /* could use gcdii(powii(a,n1), zkmultable_capZ(m)), but costly */
2349 68595 : x = ZM_hnfmodid(m, powii(a,n1));
2350 68595 : if (cx) cx = powgi(cx,n);
2351 68595 : if (s<0) {
2352 7 : GEN xZ = gcoeff(x,1,1);
2353 7 : cx = cx ? gdiv(cx, xZ): ginv(xZ);
2354 7 : x = idealHNF_inv_Z(nf,x);
2355 : }
2356 68595 : if (cx) x = ZM_Q_mul(x, cx);
2357 : }
2358 69148 : return x;
2359 : }
2360 : }
2361 :
2362 : /* raise the ideal x to the power n (in Z) */
2363 : GEN
2364 887708 : idealpow(GEN nf, GEN x, GEN n)
2365 : {
2366 : pari_sp av;
2367 : long tx;
2368 : GEN res, ax;
2369 :
2370 887708 : if (typ(n) != t_INT) pari_err_TYPE("idealpow",n);
2371 887708 : tx = idealtyp(&x,&ax);
2372 887708 : res = ax? cgetg(3,t_VEC): NULL;
2373 887708 : av = avma;
2374 887708 : x = gerepileupto(av, idealpow_aux(checknf(nf), x, tx, n));
2375 887706 : if (!ax) return x;
2376 0 : gel(res,1) = x;
2377 0 : gel(res,2) = ext_pow(nf, ax, n);
2378 0 : return res;
2379 : }
2380 :
2381 : /* Return ideal^e in number field nf. e is a C integer. */
2382 : GEN
2383 274635 : idealpows(GEN nf, GEN ideal, long e)
2384 : {
2385 274635 : long court[] = {evaltyp(t_INT) | _evallg(3),0,0};
2386 274635 : affsi(e,court); return idealpow(nf,ideal,court);
2387 : }
2388 :
2389 : static GEN
2390 28571 : _idealmulred(GEN nf, GEN x, GEN y)
2391 28571 : { return idealred(nf,idealmul(nf,x,y)); }
2392 : static GEN
2393 35650 : _idealsqrred(GEN nf, GEN x)
2394 35650 : { return idealred(nf,idealsqr(nf,x)); }
2395 : static GEN
2396 11364 : _mul(void *data, GEN x, GEN y) { return _idealmulred((GEN)data,x,y); }
2397 : static GEN
2398 35650 : _sqr(void *data, GEN x) { return _idealsqrred((GEN)data, x); }
2399 :
2400 : /* compute x^n (x ideal, n integer), reducing along the way */
2401 : GEN
2402 80263 : idealpowred(GEN nf, GEN x, GEN n)
2403 : {
2404 80263 : pari_sp av = avma, av2;
2405 : long s;
2406 : GEN y;
2407 :
2408 80263 : if (typ(n) != t_INT) pari_err_TYPE("idealpowred",n);
2409 80263 : s = signe(n); if (s == 0) return idealpow(nf,x,n);
2410 80263 : y = gen_pow_i(x, n, (void*)nf, &_sqr, &_mul);
2411 80263 : av2 = avma;
2412 80263 : if (s < 0) y = idealinv(nf,y);
2413 80263 : if (s < 0 || is_pm1(n)) y = idealred(nf,y);
2414 80264 : return avma == av2? gerepilecopy(av,y): gerepileupto(av,y);
2415 : }
2416 :
2417 : GEN
2418 17207 : idealmulred(GEN nf, GEN x, GEN y)
2419 : {
2420 17207 : pari_sp av = avma;
2421 17207 : return gerepileupto(av, _idealmulred(nf,x,y));
2422 : }
2423 :
2424 : long
2425 91 : isideal(GEN nf,GEN x)
2426 : {
2427 91 : long N, i, j, lx, tx = typ(x);
2428 : pari_sp av;
2429 : GEN T, xZ;
2430 :
2431 91 : nf = checknf(nf); T = nf_get_pol(nf); lx = lg(x);
2432 91 : if (tx==t_VEC && lx==3) { x = gel(x,1); tx = typ(x); lx = lg(x); }
2433 91 : switch(tx)
2434 : {
2435 14 : case t_INT: case t_FRAC: return 1;
2436 7 : case t_POL: return varn(x) == varn(T);
2437 7 : case t_POLMOD: return RgX_equal_var(T, gel(x,1));
2438 14 : case t_VEC: return get_prid(x)? 1 : 0;
2439 42 : case t_MAT: break;
2440 7 : default: return 0;
2441 : }
2442 42 : N = degpol(T);
2443 42 : if (lx-1 != N) return (lx == 1);
2444 28 : if (nbrows(x) != N) return 0;
2445 :
2446 28 : av = avma; x = Q_primpart(x);
2447 28 : if (!ZM_ishnf(x)) return 0;
2448 14 : xZ = gcoeff(x,1,1);
2449 21 : for (j=2; j<=N; j++)
2450 14 : if (!dvdii(xZ, gcoeff(x,j,j))) return gc_long(av,0);
2451 14 : for (i=2; i<=N; i++)
2452 14 : for (j=2; j<=N; j++)
2453 7 : if (! hnf_invimage(x, zk_ei_mul(nf,gel(x,i),j))) return gc_long(av,0);
2454 7 : return gc_long(av,1);
2455 : }
2456 :
2457 : GEN
2458 40555 : idealdiv(GEN nf, GEN x, GEN y)
2459 : {
2460 40555 : pari_sp av = avma, tetpil;
2461 40555 : GEN z = idealinv(nf,y);
2462 40555 : tetpil = avma; return gerepile(av,tetpil, idealmul(nf,x,z));
2463 : }
2464 :
2465 : /* This routine computes the quotient x/y of two ideals in the number field nf.
2466 : * It assumes that the quotient is an integral ideal. The idea is to find an
2467 : * ideal z dividing y such that gcd(Nx/Nz, Nz) = 1. Then
2468 : *
2469 : * x + (Nx/Nz) x
2470 : * ----------- = ---
2471 : * y + (Ny/Nz) y
2472 : *
2473 : * Proof: we can assume x and y are integral. Let p be any prime ideal
2474 : *
2475 : * If p | Nz, then it divides neither Nx/Nz nor Ny/Nz (since Nx/Nz is the
2476 : * product of the integers N(x/y) and N(y/z)). Both the numerator and the
2477 : * denominator on the left will be coprime to p. So will x/y, since x/y is
2478 : * assumed integral and its norm N(x/y) is coprime to p.
2479 : *
2480 : * If instead p does not divide Nz, then v_p (Nx/Nz) = v_p (Nx) >= v_p(x).
2481 : * Hence v_p (x + Nx/Nz) = v_p(x). Likewise for the denominators. QED.
2482 : *
2483 : * Peter Montgomery. July, 1994. */
2484 : static void
2485 7 : err_divexact(GEN x, GEN y)
2486 7 : { pari_err_DOMAIN("idealdivexact","denominator(x/y)", "!=",
2487 0 : gen_1,mkvec2(x,y)); }
2488 : GEN
2489 4880 : idealdivexact(GEN nf, GEN x0, GEN y0)
2490 : {
2491 4880 : pari_sp av = avma;
2492 : GEN x, y, xZ, yZ, Nx, Ny, Nz, cy, q, r;
2493 :
2494 4880 : nf = checknf(nf);
2495 4880 : x = idealhnf_shallow(nf, x0);
2496 4880 : y = idealhnf_shallow(nf, y0);
2497 4880 : if (lg(y) == 1) pari_err_INV("idealdivexact", y0);
2498 4873 : if (lg(x) == 1) { set_avma(av); return cgetg(1, t_MAT); } /* numerator is zero */
2499 4873 : y = Q_primitive_part(y, &cy);
2500 4873 : if (cy) x = RgM_Rg_div(x,cy);
2501 4873 : xZ = gcoeff(x,1,1); if (typ(xZ) != t_INT) err_divexact(x,y);
2502 4866 : yZ = gcoeff(y,1,1); if (isint1(yZ)) return gerepilecopy(av, x);
2503 2731 : Nx = idealnorm(nf,x);
2504 2731 : Ny = idealnorm(nf,y);
2505 2731 : if (typ(Nx) != t_INT) err_divexact(x,y);
2506 2731 : q = dvmdii(Nx,Ny, &r);
2507 2731 : if (signe(r)) err_divexact(x,y);
2508 2731 : if (is_pm1(q)) { set_avma(av); return matid(nf_get_degree(nf)); }
2509 : /* Find a norm Nz | Ny such that gcd(Nx/Nz, Nz) = 1 */
2510 527 : for (Nz = Ny;;) /* q = Nx/Nz */
2511 465 : {
2512 992 : GEN p1 = gcdii(Nz, q);
2513 992 : if (is_pm1(p1)) break;
2514 465 : Nz = diviiexact(Nz,p1);
2515 465 : q = mulii(q,p1);
2516 : }
2517 527 : xZ = gcoeff(x,1,1); q = gcdii(q, xZ);
2518 527 : if (!equalii(xZ,q))
2519 : { /* Replace x/y by x+(Nx/Nz) / y+(Ny/Nz) */
2520 377 : x = ZM_hnfmodid(x, q);
2521 : /* y reduced to unit ideal ? */
2522 377 : if (Nz == Ny) return gerepileupto(av, x);
2523 :
2524 118 : yZ = gcoeff(y,1,1); q = gcdii(diviiexact(Ny,Nz), yZ);
2525 118 : y = ZM_hnfmodid(y, q);
2526 : }
2527 268 : yZ = gcoeff(y,1,1);
2528 268 : y = idealHNF_mul(nf,x, idealHNF_inv_Z(nf,y));
2529 268 : return gerepileupto(av, ZM_Z_divexact(y, yZ));
2530 : }
2531 :
2532 : GEN
2533 21 : idealintersect(GEN nf, GEN x, GEN y)
2534 : {
2535 21 : pari_sp av = avma;
2536 : long lz, lx, i;
2537 : GEN z, dx, dy, xZ, yZ;;
2538 :
2539 21 : nf = checknf(nf);
2540 21 : x = idealhnf_shallow(nf,x);
2541 21 : y = idealhnf_shallow(nf,y);
2542 21 : if (lg(x) == 1 || lg(y) == 1) { set_avma(av); return cgetg(1,t_MAT); }
2543 14 : x = Q_remove_denom(x, &dx);
2544 14 : y = Q_remove_denom(y, &dy);
2545 14 : if (dx) y = ZM_Z_mul(y, dx);
2546 14 : if (dy) x = ZM_Z_mul(x, dy);
2547 14 : xZ = gcoeff(x,1,1);
2548 14 : yZ = gcoeff(y,1,1);
2549 14 : dx = mul_denom(dx,dy);
2550 14 : z = ZM_lll(shallowconcat(x,y), 0.99, LLL_KER); lz = lg(z);
2551 14 : lx = lg(x);
2552 63 : for (i=1; i<lz; i++) setlg(z[i], lx);
2553 14 : z = ZM_hnfmodid(ZM_mul(x,z), lcmii(xZ, yZ));
2554 14 : if (dx) z = RgM_Rg_div(z,dx);
2555 14 : return gerepileupto(av,z);
2556 : }
2557 :
2558 : /*******************************************************************/
2559 : /* */
2560 : /* T2-IDEAL REDUCTION */
2561 : /* */
2562 : /*******************************************************************/
2563 :
2564 : static GEN
2565 21 : chk_vdir(GEN nf, GEN vdir)
2566 : {
2567 21 : long i, l = lg(vdir);
2568 : GEN v;
2569 21 : if (l != lg(nf_get_roots(nf))) pari_err_DIM("idealred");
2570 14 : switch(typ(vdir))
2571 : {
2572 0 : case t_VECSMALL: return vdir;
2573 14 : case t_VEC: break;
2574 0 : default: pari_err_TYPE("idealred",vdir);
2575 : }
2576 14 : v = cgetg(l, t_VECSMALL);
2577 56 : for (i = 1; i < l; i++) v[i] = itos(gceil(gel(vdir,i)));
2578 14 : return v;
2579 : }
2580 :
2581 : static void
2582 12639 : twistG(GEN G, long r1, long i, long v)
2583 : {
2584 12639 : long j, lG = lg(G);
2585 12639 : if (i <= r1) {
2586 37275 : for (j=1; j<lG; j++) gcoeff(G,i,j) = gmul2n(gcoeff(G,i,j), v);
2587 : } else {
2588 578 : long k = (i<<1) - r1;
2589 4402 : for (j=1; j<lG; j++)
2590 : {
2591 3824 : gcoeff(G,k-1,j) = gmul2n(gcoeff(G,k-1,j), v);
2592 3824 : gcoeff(G,k ,j) = gmul2n(gcoeff(G,k ,j), v);
2593 : }
2594 : }
2595 12639 : }
2596 :
2597 : GEN
2598 138657 : nf_get_Gtwist(GEN nf, GEN vdir)
2599 : {
2600 : long i, l, v, r1;
2601 : GEN G;
2602 :
2603 138657 : if (!vdir) return nf_get_roundG(nf);
2604 21 : if (typ(vdir) == t_MAT)
2605 : {
2606 0 : long N = nf_get_degree(nf);
2607 0 : if (lg(vdir) != N+1 || lgcols(vdir) != N+1) pari_err_DIM("idealred");
2608 0 : return vdir;
2609 : }
2610 21 : vdir = chk_vdir(nf, vdir);
2611 14 : G = RgM_shallowcopy(nf_get_G(nf));
2612 14 : r1 = nf_get_r1(nf);
2613 14 : l = lg(vdir);
2614 56 : for (i=1; i<l; i++)
2615 : {
2616 42 : v = vdir[i]; if (!v) continue;
2617 42 : twistG(G, r1, i, v);
2618 : }
2619 14 : return RM_round_maxrank(G);
2620 : }
2621 : GEN
2622 12597 : nf_get_Gtwist1(GEN nf, long i)
2623 : {
2624 12597 : GEN G = RgM_shallowcopy( nf_get_G(nf) );
2625 12597 : long r1 = nf_get_r1(nf);
2626 12597 : twistG(G, r1, i, 10);
2627 12597 : return RM_round_maxrank(G);
2628 : }
2629 :
2630 : GEN
2631 97272 : RM_round_maxrank(GEN G0)
2632 : {
2633 97272 : long e, r = lg(G0)-1;
2634 97272 : pari_sp av = avma;
2635 97272 : for (e = 4; ; e <<= 1, set_avma(av))
2636 0 : {
2637 97272 : GEN G = gmul2n(G0, e), H = ground(G);
2638 97270 : if (ZM_rank(H) == r) return H; /* maximal rank ? */
2639 : }
2640 : }
2641 :
2642 : GEN
2643 138650 : idealred0(GEN nf, GEN I, GEN vdir)
2644 : {
2645 138650 : pari_sp av = avma;
2646 138650 : GEN G, aI, IZ, J, y, my, dyi, yi, c1 = NULL;
2647 : long N;
2648 :
2649 138650 : nf = checknf(nf);
2650 138650 : N = nf_get_degree(nf);
2651 : /* put first for sanity checks, unused when I obviously principal */
2652 138650 : G = nf_get_Gtwist(nf, vdir);
2653 138643 : switch (idealtyp(&I,&aI))
2654 : {
2655 37238 : case id_PRIME:
2656 37238 : if (pr_is_inert(I)) {
2657 655 : if (!aI) { set_avma(av); return matid(N); }
2658 655 : c1 = gel(I,1); I = matid(N);
2659 655 : goto END;
2660 : }
2661 36583 : IZ = pr_get_p(I);
2662 36583 : J = pr_inv_p(I);
2663 36583 : I = idealhnf_two(nf,I);
2664 36584 : break;
2665 101377 : case id_MAT:
2666 101377 : if (lg(I)-1 != N) pari_err_DIM("idealred");
2667 101370 : I = Q_primitive_part(I, &c1);
2668 101369 : IZ = gcoeff(I,1,1);
2669 101369 : if (is_pm1(IZ))
2670 : {
2671 8831 : if (!aI) { set_avma(av); return matid(N); }
2672 8747 : goto END;
2673 : }
2674 92538 : J = idealHNF_inv_Z(nf, I);
2675 92539 : break;
2676 21 : default: /* id_PRINCIPAL, silly case */
2677 21 : if (gequal0(I)) I = cgetg(1,t_MAT); else { c1 = I; I = matid(N); }
2678 21 : if (!aI) return I;
2679 14 : goto END;
2680 : }
2681 : /* now I integral, HNF; and J = (I\cap Z) I^(-1), integral */
2682 129123 : y = idealpseudomin(J, G); /* small elt in (I\cap Z)I^(-1), integral */
2683 129120 : if (ZV_isscalar(y))
2684 : { /* already reduced */
2685 71070 : if (!aI) return gerepilecopy(av, I);
2686 67891 : goto END;
2687 : }
2688 :
2689 58051 : my = zk_multable(nf, y);
2690 58052 : I = ZM_Z_divexact(ZM_mul(my, I), IZ); /* y I / (I\cap Z), integral */
2691 58050 : c1 = mul_content(c1, IZ);
2692 58050 : if (equali1(c1)) c1 = NULL; /* can be simplified with IZ */
2693 58050 : yi = ZM_gauss(my, col_ei(N,1)); /* y^-1 */
2694 58052 : dyi = Q_denom(yi); /* generates (y) \cap Z */
2695 58052 : I = hnfmodid(I, dyi);
2696 58052 : if (!aI) return gerepileupto(av, I);
2697 56079 : if (typ(aI) == t_MAT)
2698 : {
2699 39342 : GEN nyi = Q_muli_to_int(yi, dyi);
2700 39342 : if (gexpo(nyi) >= gexpo(y))
2701 20873 : aI = famat_div(aI, y); /* yi "larger" than y, keep the latter */
2702 : else
2703 : { /* use yi */
2704 18469 : aI = famat_mul(aI, nyi);
2705 18469 : c1 = div_content(c1, dyi);
2706 : }
2707 39342 : if (c1) { aI = famat_mul(aI, Q_to_famat(c1)); c1 = NULL; }
2708 : }
2709 : else
2710 16737 : c1 = c1? RgC_Rg_mul(yi, c1): yi;
2711 133386 : END:
2712 133386 : if (c1) aI = ext_mul(nf, aI,c1);
2713 133386 : return gerepilecopy(av, mkvec2(I, aI));
2714 : }
2715 :
2716 : /* I integral ZM (not HNF), G ZM, rounded Cholesky form of a weighted
2717 : * T2 matrix. Reduce I wrt G */
2718 : GEN
2719 1285180 : idealpseudored(GEN I, GEN G)
2720 1285180 : { return ZM_mul(I, ZM_lll(ZM_mul(G, I), 0.99, LLL_IM)); }
2721 :
2722 : /* Same I, G; m in I with T2(m) small */
2723 : GEN
2724 143037 : idealpseudomin(GEN I, GEN G)
2725 : {
2726 143037 : GEN u = ZM_lll(ZM_mul(G, I), 0.99, LLL_IM);
2727 143037 : return ZM_ZC_mul(I, gel(u,1));
2728 : }
2729 : /* Same I,G; irrational m in I with T2(m) small */
2730 : GEN
2731 0 : idealpseudomin_nonscalar(GEN I, GEN G)
2732 : {
2733 0 : GEN u = ZM_lll(ZM_mul(G, I), 0.99, LLL_IM);
2734 0 : GEN m = ZM_ZC_mul(I, gel(u,1));
2735 0 : if (ZV_isscalar(m) && lg(u) > 2) m = ZM_ZC_mul(I, gel(u,2));
2736 0 : return m;
2737 : }
2738 : /* Same I,G; t_VEC of irrational m in I with T2(m) small */
2739 : GEN
2740 1204348 : idealpseudominvec(GEN I, GEN G)
2741 : {
2742 1204348 : long i, j, k, n = lg(I)-1;
2743 1204348 : GEN x, L, b = idealpseudored(I, G);
2744 1204348 : L = cgetg(1 + (n*(n+1))/2, t_VEC);
2745 4257241 : for (i = k = 1; i <= n; i++)
2746 : {
2747 3052890 : x = gel(b,i);
2748 3052890 : if (!ZV_isscalar(x)) gel(L,k++) = x;
2749 : }
2750 3052892 : for (i = 2; i <= n; i++)
2751 : {
2752 1848545 : long J = minss(i, 4);
2753 4592198 : for (j = 1; j < J; j++)
2754 : {
2755 2743657 : x = ZC_add(gel(b,i),gel(b,j));
2756 2743653 : if (!ZV_isscalar(x)) gel(L,k++) = x;
2757 : }
2758 : }
2759 1204347 : setlg(L,k); return L;
2760 : }
2761 :
2762 : GEN
2763 13907 : idealred_elt(GEN nf, GEN I)
2764 : {
2765 13907 : pari_sp av = avma;
2766 13907 : GEN u = idealpseudomin(I, nf_get_roundG(nf));
2767 13907 : return gerepileupto(av, u);
2768 : }
2769 :
2770 : GEN
2771 7 : idealmin(GEN nf, GEN x, GEN vdir)
2772 : {
2773 7 : pari_sp av = avma;
2774 : GEN y, dx;
2775 7 : nf = checknf(nf);
2776 7 : switch( idealtyp(&x, NULL) )
2777 : {
2778 0 : case id_PRINCIPAL: return gcopy(x);
2779 0 : case id_PRIME: x = pr_hnf(nf,x); break;
2780 7 : case id_MAT: if (lg(x) == 1) return gen_0;
2781 : }
2782 7 : x = Q_remove_denom(x, &dx);
2783 7 : y = idealpseudomin(x, nf_get_Gtwist(nf,vdir));
2784 7 : if (dx) y = RgC_Rg_div(y, dx);
2785 7 : return gerepileupto(av, y);
2786 : }
2787 :
2788 : /*******************************************************************/
2789 : /* */
2790 : /* APPROXIMATION THEOREM */
2791 : /* */
2792 : /*******************************************************************/
2793 : /* a = ppi(a,b) ppo(a,b), where ppi regroups primes common to a and b
2794 : * and ppo(a,b) = Z_ppo(a,b) */
2795 : /* return gcd(a,b),ppi(a,b),ppo(a,b) */
2796 : GEN
2797 532735 : Z_ppio(GEN a, GEN b)
2798 : {
2799 532735 : GEN x, y, d = gcdii(a,b);
2800 532735 : if (is_pm1(d)) return mkvec3(gen_1, gen_1, a);
2801 413028 : x = d; y = diviiexact(a,d);
2802 : for(;;)
2803 68425 : {
2804 481453 : GEN g = gcdii(x,y);
2805 481453 : if (is_pm1(g)) return mkvec3(d, x, y);
2806 68425 : x = mulii(x,g); y = diviiexact(y,g);
2807 : }
2808 : }
2809 : /* a = ppg(a,b)pple(a,b), where ppg regroups primes such that v(a) > v(b)
2810 : * and pple all others */
2811 : /* return gcd(a,b),ppg(a,b),pple(a,b) */
2812 : GEN
2813 0 : Z_ppgle(GEN a, GEN b)
2814 : {
2815 0 : GEN x, y, g, d = gcdii(a,b);
2816 0 : if (equalii(a, d)) return mkvec3(a, gen_1, a);
2817 0 : x = diviiexact(a,d); y = d;
2818 : for(;;)
2819 : {
2820 0 : g = gcdii(x,y);
2821 0 : if (is_pm1(g)) return mkvec3(d, x, y);
2822 0 : x = mulii(x,g); y = diviiexact(y,g);
2823 : }
2824 : }
2825 : static void
2826 0 : Z_dcba_rec(GEN L, GEN a, GEN b)
2827 : {
2828 : GEN x, r, v, g, h, c, c0;
2829 : long n;
2830 0 : if (is_pm1(b)) {
2831 0 : if (!is_pm1(a)) vectrunc_append(L, a);
2832 0 : return;
2833 : }
2834 0 : v = Z_ppio(a,b);
2835 0 : a = gel(v,2);
2836 0 : r = gel(v,3);
2837 0 : if (!is_pm1(r)) vectrunc_append(L, r);
2838 0 : v = Z_ppgle(a,b);
2839 0 : g = gel(v,1);
2840 0 : h = gel(v,2);
2841 0 : x = c0 = gel(v,3);
2842 0 : for (n = 1; !is_pm1(h); n++)
2843 : {
2844 : GEN d, y;
2845 : long i;
2846 0 : v = Z_ppgle(h,sqri(g));
2847 0 : g = gel(v,1);
2848 0 : h = gel(v,2);
2849 0 : c = gel(v,3); if (is_pm1(c)) continue;
2850 0 : d = gcdii(c,b);
2851 0 : x = mulii(x,d);
2852 0 : y = d; for (i=1; i < n; i++) y = sqri(y);
2853 0 : Z_dcba_rec(L, diviiexact(c,y), d);
2854 : }
2855 0 : Z_dcba_rec(L,diviiexact(b,x), c0);
2856 : }
2857 : static GEN
2858 3765118 : Z_cba_rec(GEN L, GEN a, GEN b)
2859 : {
2860 : GEN g;
2861 : /* a few naive steps before switching to dcba */
2862 3765118 : if (lg(L) > 10) { Z_dcba_rec(L, a, b); return veclast(L); }
2863 3765118 : if (is_pm1(a)) return b;
2864 2242975 : g = gcdii(a,b);
2865 2242975 : if (is_pm1(g)) { vectrunc_append(L, a); return b; }
2866 1676143 : a = diviiexact(a,g);
2867 1676143 : b = diviiexact(b,g);
2868 1676143 : return Z_cba_rec(L, Z_cba_rec(L, a, g), b);
2869 : }
2870 : GEN
2871 412832 : Z_cba(GEN a, GEN b)
2872 : {
2873 412832 : GEN L = vectrunc_init(expi(a) + expi(b) + 2);
2874 412832 : GEN t = Z_cba_rec(L, a, b);
2875 412832 : if (!is_pm1(t)) vectrunc_append(L, t);
2876 412832 : return L;
2877 : }
2878 : /* P = coprime base, extend it by b; TODO: quadratic for now */
2879 : GEN
2880 35 : ZV_cba_extend(GEN P, GEN b)
2881 : {
2882 35 : long i, l = lg(P);
2883 35 : GEN w = cgetg(l+1, t_VEC);
2884 133 : for (i = 1; i < l; i++)
2885 : {
2886 98 : GEN v = Z_cba(gel(P,i), b);
2887 98 : long nv = lg(v)-1;
2888 98 : gel(w,i) = vecslice(v, 1, nv-1); /* those divide P[i] but not b */
2889 98 : b = gel(v,nv);
2890 : }
2891 35 : gel(w,l) = b; return shallowconcat1(w);
2892 : }
2893 : GEN
2894 28 : ZV_cba(GEN v)
2895 : {
2896 28 : long i, l = lg(v);
2897 : GEN P;
2898 28 : if (l <= 2) return v;
2899 14 : P = Z_cba(gel(v,1), gel(v,2));
2900 42 : for (i = 3; i < l; i++) P = ZV_cba_extend(P, gel(v,i));
2901 14 : return P;
2902 : }
2903 :
2904 : /* write x = x1 x2, x2 maximal s.t. (x2,f) = 1, return x2 */
2905 : GEN
2906 56728214 : Z_ppo(GEN x, GEN f)
2907 : {
2908 : for (;;)
2909 : {
2910 56728214 : f = gcdii(x, f); if (is_pm1(f)) break;
2911 38583798 : x = diviiexact(x, f);
2912 : }
2913 18144322 : return x;
2914 : }
2915 : /* write x = x1 x2, x2 maximal s.t. (x2,f) = 1, return x2 */
2916 : ulong
2917 70412940 : u_ppo(ulong x, ulong f)
2918 : {
2919 : for (;;)
2920 : {
2921 70412940 : f = ugcd(x, f); if (f == 1) break;
2922 16165489 : x /= f;
2923 : }
2924 54247409 : return x;
2925 : }
2926 :
2927 : /* result known to be representable as an ulong */
2928 : static ulong
2929 1598653 : lcmuu(ulong a, ulong b) { ulong d = ugcd(a,b); return (a/d) * b; }
2930 :
2931 : /* assume 0 < x < N; return u in (Z/NZ)^* such that u x = gcd(x,N) (mod N);
2932 : * set *pd = gcd(x,N) */
2933 : ulong
2934 5744640 : Fl_invgen(ulong x, ulong N, ulong *pd)
2935 : {
2936 : ulong d, d0, e, v, v1;
2937 : long s;
2938 5744640 : *pd = d = xgcduu(N, x, 0, &v, &v1, &s);
2939 5745327 : if (s > 0) v = N - v;
2940 5745327 : if (d == 1) return v;
2941 : /* vx = gcd(x,N) (mod N), v coprime to N/d but need not be coprime to N */
2942 2709478 : e = N / d;
2943 2709478 : d0 = u_ppo(d, e); /* d = d0 d1, d0 coprime to N/d, rad(d1) | N/d */
2944 2709632 : if (d0 == 1) return v;
2945 1598609 : e = lcmuu(e, d / d0);
2946 1598639 : return u_chinese_coprime(v, 1, e, d0, e*d0);
2947 : }
2948 :
2949 : /* x t_INT, f ideal. Write x = x1 x2, sqf(x1) | f, (x2,f) = 1. Return x2 */
2950 : static GEN
2951 126 : nf_coprime_part(GEN nf, GEN x, GEN listpr)
2952 : {
2953 126 : long v, j, lp = lg(listpr), N = nf_get_degree(nf);
2954 : GEN x1, x2, ex;
2955 :
2956 : #if 0 /*1) via many gcds. Expensive ! */
2957 : GEN f = idealprodprime(nf, listpr);
2958 : f = ZM_hnfmodid(f, x); /* first gcd is less expensive since x in Z */
2959 : x = scalarmat(x, N);
2960 : for (;;)
2961 : {
2962 : if (gequal1(gcoeff(f,1,1))) break;
2963 : x = idealdivexact(nf, x, f);
2964 : f = ZM_hnfmodid(shallowconcat(f,x), gcoeff(x,1,1)); /* gcd(f,x) */
2965 : }
2966 : x2 = x;
2967 : #else /*2) from prime decomposition */
2968 126 : x1 = NULL;
2969 350 : for (j=1; j<lp; j++)
2970 : {
2971 224 : GEN pr = gel(listpr,j);
2972 224 : v = Z_pval(x, pr_get_p(pr)); if (!v) continue;
2973 :
2974 126 : ex = muluu(v, pr_get_e(pr)); /* = v_pr(x) > 0 */
2975 126 : x1 = x1? idealmulpowprime(nf, x1, pr, ex)
2976 126 : : idealpow(nf, pr, ex);
2977 : }
2978 126 : x = scalarmat(x, N);
2979 126 : x2 = x1? idealdivexact(nf, x, x1): x;
2980 : #endif
2981 126 : return x2;
2982 : }
2983 :
2984 : /* L0 in K^*, assume (L0,f) = 1. Return L integral, L0 = L mod f */
2985 : GEN
2986 10920 : make_integral(GEN nf, GEN L0, GEN f, GEN listpr)
2987 : {
2988 : GEN fZ, t, L, D2, d1, d2, d;
2989 :
2990 10920 : L = Q_remove_denom(L0, &d);
2991 10920 : if (!d) return L0;
2992 :
2993 : /* L0 = L / d, L integral */
2994 518 : fZ = gcoeff(f,1,1);
2995 518 : if (typ(L) == t_INT) return Fp_mul(L, Fp_inv(d, fZ), fZ);
2996 : /* Kill denom part coprime to fZ */
2997 126 : d2 = Z_ppo(d, fZ);
2998 126 : t = Fp_inv(d2, fZ); if (!is_pm1(t)) L = ZC_Z_mul(L,t);
2999 126 : if (equalii(d, d2)) return L;
3000 :
3001 126 : d1 = diviiexact(d, d2);
3002 : /* L0 = (L / d1) mod f. d1 not coprime to f
3003 : * write (d1) = D1 D2, D2 minimal, (D2,f) = 1. */
3004 126 : D2 = nf_coprime_part(nf, d1, listpr);
3005 126 : t = idealaddtoone_i(nf, D2, f); /* in D2, 1 mod f */
3006 126 : L = nfmuli(nf,t,L);
3007 :
3008 : /* if (L0, f) = 1, then L in D1 ==> in D1 D2 = (d1) */
3009 126 : return Q_div_to_int(L, d1); /* exact division */
3010 : }
3011 :
3012 : /* assume L is a list of prime ideals. Return the product */
3013 : GEN
3014 483 : idealprodprime(GEN nf, GEN L)
3015 : {
3016 483 : long l = lg(L), i;
3017 : GEN z;
3018 483 : if (l == 1) return matid(nf_get_degree(nf));
3019 483 : z = pr_hnf(nf, gel(L,1));
3020 511 : for (i=2; i<l; i++) z = idealHNF_mul_two(nf,z, gel(L,i));
3021 483 : return z;
3022 : }
3023 :
3024 : /* optimize for the frequent case I = nfhnf()[2]: lots of them are 1 */
3025 : GEN
3026 462 : idealprod(GEN nf, GEN I)
3027 : {
3028 462 : long i, l = lg(I);
3029 : GEN z;
3030 1134 : for (i = 1; i < l; i++)
3031 1127 : if (!equali1(gel(I,i))) break;
3032 462 : if (i == l) return gen_1;
3033 455 : z = gel(I,i);
3034 763 : for (i++; i<l; i++) z = idealmul(nf, z, gel(I,i));
3035 455 : return z;
3036 : }
3037 :
3038 : /* v_pr(idealprod(nf,I)) */
3039 : long
3040 2904 : idealprodval(GEN nf, GEN I, GEN pr)
3041 : {
3042 2904 : long i, l = lg(I), v = 0;
3043 16363 : for (i = 1; i < l; i++)
3044 13459 : if (!equali1(gel(I,i))) v += idealval(nf, gel(I,i), pr);
3045 2904 : return v;
3046 : }
3047 :
3048 : /* assume L is a list of prime ideals. Return prod L[i]^e[i] */
3049 : GEN
3050 61152 : factorbackprime(GEN nf, GEN L, GEN e)
3051 : {
3052 61152 : long l = lg(L), i;
3053 : GEN z;
3054 :
3055 61152 : if (l == 1) return matid(nf_get_degree(nf));
3056 47719 : z = idealpow(nf, gel(L,1), gel(e,1));
3057 78449 : for (i=2; i<l; i++)
3058 30730 : if (signe(gel(e,i))) z = idealmulpowprime(nf,z, gel(L,i),gel(e,i));
3059 47719 : return z;
3060 : }
3061 :
3062 : /* F in Z, divisible exactly by pr.p. Return F-uniformizer for pr, i.e.
3063 : * a t in Z_K such that v_pr(t) = 1 and (t, F/pr) = 1 */
3064 : GEN
3065 58002 : pr_uniformizer(GEN pr, GEN F)
3066 : {
3067 58002 : GEN p = pr_get_p(pr), t = pr_get_gen(pr);
3068 58002 : if (!equalii(F, p))
3069 : {
3070 36645 : long e = pr_get_e(pr);
3071 36645 : GEN u, v, q = (e == 1)? sqri(p): p;
3072 36645 : u = mulii(q, Fp_inv(q, diviiexact(F,p))); /* 1 mod F/p, 0 mod q */
3073 36645 : v = subui(1UL, u); /* 0 mod F/p, 1 mod q */
3074 36645 : if (pr_is_inert(pr))
3075 28 : t = addii(mulii(p, v), u);
3076 : else
3077 : {
3078 36616 : t = ZC_Z_mul(t, v);
3079 36616 : gel(t,1) = addii(gel(t,1), u); /* return u + vt */
3080 : }
3081 : }
3082 58004 : return t;
3083 : }
3084 : /* L = list of prime ideals, return lcm_i (L[i] \cap \ZM) */
3085 : GEN
3086 81903 : prV_lcm_capZ(GEN L)
3087 : {
3088 81903 : long i, r = lg(L);
3089 : GEN F;
3090 81903 : if (r == 1) return gen_1;
3091 69345 : F = pr_get_p(gel(L,1));
3092 122825 : for (i = 2; i < r; i++)
3093 : {
3094 53481 : GEN pr = gel(L,i), p = pr_get_p(pr);
3095 53481 : if (!dvdii(F, p)) F = mulii(F,p);
3096 : }
3097 69344 : return F;
3098 : }
3099 : /* v vector of prid. Return underlying list of rational primes */
3100 : GEN
3101 61285 : prV_primes(GEN v)
3102 : {
3103 61285 : long i, l = lg(v);
3104 61285 : GEN w = cgetg(l,t_VEC);
3105 198247 : for (i=1; i<l; i++) gel(w,i) = pr_get_p(gel(v,i));
3106 61285 : return ZV_sort_uniq(w);
3107 : }
3108 :
3109 : /* Given a prime ideal factorization with possibly zero or negative
3110 : * exponents, gives b such that v_p(b) = v_p(x) for all prime ideals pr | x
3111 : * and v_pr(b) >= 0 for all other pr.
3112 : * For optimal performance, all [anti-]uniformizers should be precomputed,
3113 : * but no support for this yet. If nored, do not reduce result. */
3114 : static GEN
3115 53945 : idealapprfact_i(GEN nf, GEN x, int nored)
3116 : {
3117 53945 : GEN d = NULL, z, L, e, e2, F;
3118 : long i, r;
3119 53945 : int hasden = 0;
3120 :
3121 53945 : nf = checknf(nf);
3122 53945 : L = gel(x,1);
3123 53945 : e = gel(x,2);
3124 53945 : F = prV_lcm_capZ(L);
3125 53944 : z = NULL; r = lg(e);
3126 136370 : for (i = 1; i < r; i++)
3127 : {
3128 82424 : long s = signe(gel(e,i));
3129 : GEN pi, q;
3130 82424 : if (!s) continue;
3131 53634 : if (s < 0) hasden = 1;
3132 53634 : pi = pr_uniformizer(gel(L,i), F);
3133 53636 : q = nfpow(nf, pi, gel(e,i));
3134 53636 : z = z? nfmul(nf, z, q): q;
3135 : }
3136 53946 : if (!z) return gen_1;
3137 26801 : if (hasden) /* denominator */
3138 : {
3139 10074 : z = Q_remove_denom(z, &d);
3140 10074 : d = diviiexact(d, Z_ppo(d, F));
3141 : }
3142 26801 : if (nored || typ(z) != t_COL) return d? gdiv(z, d): z;
3143 10074 : e2 = cgetg(r, t_VEC);
3144 28608 : for (i = 1; i < r; i++) gel(e2,i) = addiu(gel(e,i), 1);
3145 10074 : x = factorbackprime(nf, L, e2);
3146 10074 : if (d) x = RgM_Rg_mul(x, d);
3147 10074 : z = ZC_reducemodlll(z, x);
3148 10074 : return d? RgC_Rg_div(z,d): z;
3149 : }
3150 :
3151 : GEN
3152 0 : idealapprfact(GEN nf, GEN x) {
3153 0 : pari_sp av = avma;
3154 0 : return gerepileupto(av, idealapprfact_i(nf, x, 0));
3155 : }
3156 : GEN
3157 14 : idealappr(GEN nf, GEN x) {
3158 14 : pari_sp av = avma;
3159 14 : if (!is_nf_extfactor(x)) x = idealfactor(nf, x);
3160 14 : return gerepileupto(av, idealapprfact_i(nf, x, 0));
3161 : }
3162 :
3163 : /* OBSOLETE */
3164 : GEN
3165 14 : idealappr0(GEN nf, GEN x, long fl) { (void)fl; return idealappr(nf, x); }
3166 :
3167 : static GEN
3168 21 : mat_ideal_two_elt2(GEN nf, GEN x, GEN a)
3169 : {
3170 21 : GEN F = idealfactor(nf,a), P = gel(F,1), E = gel(F,2);
3171 21 : long i, r = lg(E);
3172 84 : for (i=1; i<r; i++) gel(E,i) = stoi( idealval(nf,x,gel(P,i)) );
3173 21 : return idealapprfact_i(nf,F,1);
3174 : }
3175 :
3176 : static void
3177 14 : not_in_ideal(GEN a) {
3178 14 : pari_err_DOMAIN("idealtwoelt2","element mod ideal", "!=", gen_0, a);
3179 0 : }
3180 : /* x integral in HNF, a an 'nf' */
3181 : static int
3182 28 : in_ideal(GEN x, GEN a)
3183 : {
3184 28 : switch(typ(a))
3185 : {
3186 14 : case t_INT: return dvdii(a, gcoeff(x,1,1));
3187 7 : case t_COL: return RgV_is_ZV(a) && !!hnf_invimage(x, a);
3188 7 : default: return 0;
3189 : }
3190 : }
3191 :
3192 : /* Given an integral ideal x and a in x, gives a b such that
3193 : * x = aZ_K + bZ_K using the approximation theorem */
3194 : GEN
3195 42 : idealtwoelt2(GEN nf, GEN x, GEN a)
3196 : {
3197 42 : pari_sp av = avma;
3198 : GEN cx, b;
3199 :
3200 42 : nf = checknf(nf);
3201 42 : a = nf_to_scalar_or_basis(nf, a);
3202 42 : x = idealhnf_shallow(nf,x);
3203 42 : if (lg(x) == 1)
3204 : {
3205 14 : if (!isintzero(a)) not_in_ideal(a);
3206 7 : set_avma(av); return gen_0;
3207 : }
3208 28 : x = Q_primitive_part(x, &cx);
3209 28 : if (cx) a = gdiv(a, cx);
3210 28 : if (!in_ideal(x, a)) not_in_ideal(a);
3211 21 : b = mat_ideal_two_elt2(nf, x, a);
3212 21 : if (typ(b) == t_COL)
3213 : {
3214 14 : GEN mod = idealhnf_principal(nf,a);
3215 14 : b = ZC_hnfrem(b,mod);
3216 14 : if (ZV_isscalar(b)) b = gel(b,1);
3217 : }
3218 : else
3219 : {
3220 7 : GEN aZ = typ(a) == t_COL? Q_denom(zk_inv(nf,a)): a; /* (a) \cap Z */
3221 7 : b = centermodii(b, aZ, shifti(aZ,-1));
3222 : }
3223 21 : b = cx? gmul(b,cx): gcopy(b);
3224 21 : return gerepileupto(av, b);
3225 : }
3226 :
3227 : /* Given 2 integral ideals x and y in nf, returns a beta in nf such that
3228 : * beta * x is an integral ideal coprime to y */
3229 : GEN
3230 37209 : idealcoprimefact(GEN nf, GEN x, GEN fy)
3231 : {
3232 37209 : GEN L = gel(fy,1), e;
3233 37209 : long i, r = lg(L);
3234 :
3235 37209 : e = cgetg(r, t_COL);
3236 76072 : for (i=1; i<r; i++) gel(e,i) = stoi( -idealval(nf,x,gel(L,i)) );
3237 37211 : return idealapprfact_i(nf, mkmat2(L,e), 0);
3238 : }
3239 : GEN
3240 84 : idealcoprime(GEN nf, GEN x, GEN y)
3241 : {
3242 84 : pari_sp av = avma;
3243 84 : return gerepileupto(av, idealcoprimefact(nf, x, idealfactor(nf,y)));
3244 : }
3245 :
3246 : GEN
3247 7 : nfmulmodpr(GEN nf, GEN x, GEN y, GEN modpr)
3248 : {
3249 7 : pari_sp av = avma;
3250 7 : GEN z, p, pr = modpr, T;
3251 :
3252 7 : nf = checknf(nf); modpr = nf_to_Fq_init(nf,&pr,&T,&p);
3253 0 : x = nf_to_Fq(nf,x,modpr);
3254 0 : y = nf_to_Fq(nf,y,modpr);
3255 0 : z = Fq_mul(x,y,T,p);
3256 0 : return gerepileupto(av, algtobasis(nf, Fq_to_nf(z,modpr)));
3257 : }
3258 :
3259 : GEN
3260 0 : nfdivmodpr(GEN nf, GEN x, GEN y, GEN modpr)
3261 : {
3262 0 : pari_sp av = avma;
3263 0 : nf = checknf(nf);
3264 0 : return gerepileupto(av, nfreducemodpr(nf, nfdiv(nf,x,y), modpr));
3265 : }
3266 :
3267 : GEN
3268 0 : nfpowmodpr(GEN nf, GEN x, GEN k, GEN modpr)
3269 : {
3270 0 : pari_sp av=avma;
3271 0 : GEN z, T, p, pr = modpr;
3272 :
3273 0 : nf = checknf(nf); modpr = nf_to_Fq_init(nf,&pr,&T,&p);
3274 0 : z = nf_to_Fq(nf,x,modpr);
3275 0 : z = Fq_pow(z,k,T,p);
3276 0 : return gerepileupto(av, algtobasis(nf, Fq_to_nf(z,modpr)));
3277 : }
3278 :
3279 : GEN
3280 0 : nfkermodpr(GEN nf, GEN x, GEN modpr)
3281 : {
3282 0 : pari_sp av = avma;
3283 0 : GEN T, p, pr = modpr;
3284 :
3285 0 : nf = checknf(nf); modpr = nf_to_Fq_init(nf, &pr,&T,&p);
3286 0 : if (typ(x)!=t_MAT) pari_err_TYPE("nfkermodpr",x);
3287 0 : x = nfM_to_FqM(x, nf, modpr);
3288 0 : return gerepilecopy(av, FqM_to_nfM(FqM_ker(x,T,p), modpr));
3289 : }
3290 :
3291 : GEN
3292 0 : nfsolvemodpr(GEN nf, GEN a, GEN b, GEN pr)
3293 : {
3294 0 : const char *f = "nfsolvemodpr";
3295 0 : pari_sp av = avma;
3296 : GEN T, p, modpr;
3297 :
3298 0 : nf = checknf(nf);
3299 0 : modpr = nf_to_Fq_init(nf, &pr,&T,&p);
3300 0 : if (typ(a)!=t_MAT) pari_err_TYPE(f,a);
3301 0 : a = nfM_to_FqM(a, nf, modpr);
3302 0 : switch(typ(b))
3303 : {
3304 0 : case t_MAT:
3305 0 : b = nfM_to_FqM(b, nf, modpr);
3306 0 : b = FqM_gauss(a,b,T,p);
3307 0 : if (!b) pari_err_INV(f,a);
3308 0 : a = FqM_to_nfM(b, modpr);
3309 0 : break;
3310 0 : case t_COL:
3311 0 : b = nfV_to_FqV(b, nf, modpr);
3312 0 : b = FqM_FqC_gauss(a,b,T,p);
3313 0 : if (!b) pari_err_INV(f,a);
3314 0 : a = FqV_to_nfV(b, modpr);
3315 0 : break;
3316 0 : default: pari_err_TYPE(f,b);
3317 : }
3318 0 : return gerepilecopy(av, a);
3319 : }
|