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., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. */
14 : #include "pari.h"
15 : #include "paripriv.h"
16 :
17 : #define DEBUGLEVEL DEBUGLEVEL_genus2red
18 :
19 : /* extract coefficients of a polynomial a0 X^6 + ... + a6, of degree <= 6 */
20 : static void
21 1750 : RgX_to_06(GEN q, GEN *a0, GEN *a1, GEN *a2, GEN *a3, GEN *a4, GEN *a5, GEN *a6)
22 : {
23 1750 : *a0 = gen_0;
24 1750 : *a1 = gen_0;
25 1750 : *a2 = gen_0;
26 1750 : *a3 = gen_0;
27 1750 : *a4 = gen_0;
28 1750 : *a5 = gen_0;
29 1750 : *a6 = gen_0;
30 1750 : switch(degpol(q))
31 : {
32 1274 : case 6: *a0 = gel(q,8); /*fall through*/
33 1750 : case 5: *a1 = gel(q,7); /*fall through*/
34 1750 : case 4: *a2 = gel(q,6); /*fall through*/
35 1750 : case 3: *a3 = gel(q,5); /*fall through*/
36 1750 : case 2: *a4 = gel(q,4); /*fall through*/
37 1750 : case 1: *a5 = gel(q,3); /*fall through*/
38 1750 : case 0: *a6 = gel(q,2); /*fall through*/
39 : }
40 1750 : }
41 :
42 : /********************************************************************/
43 : /** **/
44 : /** IGUSA INVARIANTS **/
45 : /** (GP2C-generated) **/
46 : /** **/
47 : /********************************************************************/
48 : /*
49 : j2(a0,a1,a2,a3,a4,a5,a6) = (-120*a0*a6+20*a1*a5-8*a2*a4+3*a3^2) / 4;
50 : */
51 : static GEN
52 1463 : igusaj2(GEN a0, GEN a1, GEN a2, GEN a3, GEN a4, GEN a5, GEN a6)
53 : {
54 1463 : pari_sp av = avma;
55 1463 : return gerepileupto(av, gmul2n(gadd(gsub(gadd(gmul(gmulsg(-120, a0), a6), gmul(gmulsg(20, a1), a5)), gmul(gmulsg(8, a2), a4)), gmulsg(3, gsqr(a3))), -2));
56 : }
57 :
58 : /*
59 : j4(a0,a1,a2,a3,a4,a5,a6) = (240*(a0*a3*a4*a5+a1*a2*a3*a6)-400*(a0*a2*a5^2+a1^2*a4*a6)-64*(a0*a4^3+a2^3*a6)+16*(a1*a3*a4^2+a2^2*a3*a5)-672*a0*a3^2*a6+240*a1^2*a5^2-112*a1*a2*a4*a5-8*a1*a3^2*a5+16*a2^2*a4^2-16*a2*a3^2*a4+3*a3^4+2640*a0^2*a6^2-880*a0*a1*a5*a6+1312*a0*a2*a4*a6) / 2^7
60 : */
61 : static GEN
62 1463 : igusaj4(GEN a0, GEN a1, GEN a2, GEN a3, GEN a4, GEN a5, GEN a6)
63 : {
64 1463 : pari_sp av = avma;
65 1463 : return gerepileupto(av,
66 : gmul2n(gadd(gsub(gadd(gadd(gsub(gadd(gsub(gsub(gadd(gsub(gadd(gsub(gsub(gmulsg(240,
67 : gadd(gmul(gmul(gmul(a0, a3), a4), a5), gmul(gmul(gmul(a1, a2), a3), a6))),
68 : gmulsg(400, gadd(gmul(gmul(a0, a2), gsqr(a5)), gmul(gmul(gsqr(a1), a4), a6)))),
69 : gmulsg(64, gadd(gmul(a0, gpowgs(a4, 3)), gmul(gpowgs(a2, 3), a6)))), gmulsg(16,
70 : gadd(gmul(gmul(a1, a3), gsqr(a4)), gmul(gmul(gsqr(a2), a3), a5)))),
71 : gmul(gmul(gmulsg(672, a0), gsqr(a3)), a6)), gmul(gmulsg(240, gsqr(a1)),
72 : gsqr(a5))), gmul(gmul(gmul(gmulsg(112, a1), a2), a4), a5)), gmul(gmul(gmulsg(8,
73 : a1), gsqr(a3)), a5)), gmul(gmulsg(16, gsqr(a2)), gsqr(a4))),
74 : gmul(gmul(gmulsg(16, a2), gsqr(a3)), a4)), gmulsg(3, gpowgs(a3, 4))),
75 : gmul(gmulsg(2640, gsqr(a0)), gsqr(a6))), gmul(gmul(gmul(gmulsg(880, a0), a1),
76 : a5), a6)), gmul(gmul(gmul(gmulsg(1312, a0), a2), a4), a6)), -7));
77 : }
78 :
79 : /*
80 : j6(a0,a1,a2,a3,a4,a5,a6) = (1600*(a0^2*a4^2*a5^2+a1^2*a2^2*a6^2)+1600*(a0*a1*a2*a5^3+a1^3*a4*a5*a6)+640*(a0*a1*a3*a4*a5^2+a1^2*a2*a3*a5*a6)-4000*(a0^2*a3*a5^3+a1^3*a3*a6^2)-384*(a0*a1*a4^3*a5+a1*a2^3*a5*a6)-640*(a0*a2^2*a4*a5^2+a1^2*a2*a4^2*a6)+80*(a0*a2*a3^2*a5^2+a1^2*a3^2*a4*a6)+192*(a0*a2*a3*a4^2*a5+a1*a2^2*a3*a4*a6)-48*(a0*a3^3*a4*a5+a1*a2*a3^3*a6)-224*(a1^2*a3*a4^2*a5+a1*a2^2*a3*a5^2)+64*(a1^2*a4^4+a2^4*a5^2)-64*(a1*a2*a3*a4^3+a2^3*a3*a4*a5)+16*(a1*a3^3*a4^2+a2^2*a3^3*a5)-4096*(a0^2*a4^3*a6+a0*a2^3*a6^2)+6400*(a0^2*a2*a5^2*a6+a0*a1^2*a4*a6^2)+10560*(a0^2*a3*a4*a5*a6+a0*a1*a2*a3*a6^2)+2624*(a0*a1*a3*a4^2*a6+a0*a2^2*a3*a5*a6)-4432*a0*a1*a3^2*a5*a6-8*a2*a3^4*a4+a3^6-320*a1^3*a5^3+64*a1^2*a2*a4*a5^2+176*a1^2*a3^2*a5^2+128*a1*a2^2*a4^2*a5+112*a1*a2*a3^2*a4*a5-28*a1*a3^4*a5+16*a2^2*a3^2*a4^2+5120*a0^3*a6^3-2544*a0^2*a3^2*a6^2+312*a0*a3^4*a6-14336*a0^2*a2*a4*a6^2+1024*a0*a2^2*a4^2*a6-2560*a0^2*a1*a5*a6^2-2240*a0*a1^2*a5^2*a6-6528*a0*a1*a2*a4*a5*a6-1568*a0*a2*a3^2*a4*a6) / 2^10
81 : */
82 : static GEN
83 1463 : igusaj6(GEN a0, GEN a1, GEN a2, GEN a3, GEN a4, GEN a5, GEN a6)
84 : {
85 1463 : pari_sp av = avma;
86 1463 : return gerepileupto(av,
87 : gmul2n(gsub(gsub(gsub(gsub(gadd(gsub(gadd(gsub(gadd(gadd(gsub(gadd(gadd(gadd(gadd(gsub(gadd(gsub(gsub(gadd(gadd(gadd(gsub(gadd(gsub(gadd(gsub(gsub(gadd(gadd(gsub(gsub(gsub(gadd(gadd(gmulsg(1600,
88 : gadd(gmul(gmul(gsqr(a0), gsqr(a4)), gsqr(a5)), gmul(gmul(gsqr(a1), gsqr(a2)),
89 : gsqr(a6)))), gmulsg(1600, gadd(gmul(gmul(gmul(a0, a1), a2), gpowgs(a5, 3)),
90 : gmul(gmul(gmul(gpowgs(a1, 3), a4), a5), a6)))), gmulsg(640,
91 : gadd(gmul(gmul(gmul(gmul(a0, a1), a3), a4), gsqr(a5)),
92 : gmul(gmul(gmul(gmul(gsqr(a1), a2), a3), a5), a6)))), gmulsg(4000,
93 : gadd(gmul(gmul(gsqr(a0), a3), gpowgs(a5, 3)), gmul(gmul(gpowgs(a1, 3), a3),
94 : gsqr(a6))))), gmulsg(384, gadd(gmul(gmul(gmul(a0, a1), gpowgs(a4, 3)), a5),
95 : gmul(gmul(gmul(a1, gpowgs(a2, 3)), a5), a6)))), gmulsg(640,
96 : gadd(gmul(gmul(gmul(a0, gsqr(a2)), a4), gsqr(a5)), gmul(gmul(gmul(gsqr(a1),
97 : a2), gsqr(a4)), a6)))), gmulsg(80, gadd(gmul(gmul(gmul(a0, a2), gsqr(a3)),
98 : gsqr(a5)), gmul(gmul(gmul(gsqr(a1), gsqr(a3)), a4), a6)))), gmulsg(192,
99 : gadd(gmul(gmul(gmul(gmul(a0, a2), a3), gsqr(a4)), a5), gmul(gmul(gmul(gmul(a1,
100 : gsqr(a2)), a3), a4), a6)))), gmulsg(48, gadd(gmul(gmul(gmul(a0, gpowgs(a3, 3)),
101 : a4), a5), gmul(gmul(gmul(a1, a2), gpowgs(a3, 3)), a6)))), gmulsg(224,
102 : gadd(gmul(gmul(gmul(gsqr(a1), a3), gsqr(a4)), a5), gmul(gmul(gmul(a1,
103 : gsqr(a2)), a3), gsqr(a5))))), gmulsg(64, gadd(gmul(gsqr(a1), gpowgs(a4, 4)),
104 : gmul(gpowgs(a2, 4), gsqr(a5))))), gmulsg(64, gadd(gmul(gmul(gmul(a1, a2), a3),
105 : gpowgs(a4, 3)), gmul(gmul(gmul(gpowgs(a2, 3), a3), a4), a5)))), gmulsg(16,
106 : gadd(gmul(gmul(a1, gpowgs(a3, 3)), gsqr(a4)), gmul(gmul(gsqr(a2), gpowgs(a3,
107 : 3)), a5)))), gmulsg(4096, gadd(gmul(gmul(gsqr(a0), gpowgs(a4, 3)), a6),
108 : gmul(gmul(a0, gpowgs(a2, 3)), gsqr(a6))))), gmulsg(6400,
109 : gadd(gmul(gmul(gmul(gsqr(a0), a2), gsqr(a5)), a6), gmul(gmul(gmul(a0,
110 : gsqr(a1)), a4), gsqr(a6))))), gmulsg(10560, gadd(gmul(gmul(gmul(gmul(gsqr(a0),
111 : a3), a4), a5), a6), gmul(gmul(gmul(gmul(a0, a1), a2), a3), gsqr(a6))))),
112 : gmulsg(2624, gadd(gmul(gmul(gmul(gmul(a0, a1), a3), gsqr(a4)), a6),
113 : gmul(gmul(gmul(gmul(a0, gsqr(a2)), a3), a5), a6)))),
114 : gmul(gmul(gmul(gmul(gmulsg(4432, a0), a1), gsqr(a3)), a5), a6)),
115 : gmul(gmul(gmulsg(8, a2), gpowgs(a3, 4)), a4)), gpowgs(a3, 6)), gmul(gmulsg(320,
116 : gpowgs(a1, 3)), gpowgs(a5, 3))), gmul(gmul(gmul(gmulsg(64, gsqr(a1)), a2), a4),
117 : gsqr(a5))), gmul(gmul(gmulsg(176, gsqr(a1)), gsqr(a3)), gsqr(a5))),
118 : gmul(gmul(gmul(gmulsg(128, a1), gsqr(a2)), gsqr(a4)), a5)),
119 : gmul(gmul(gmul(gmul(gmulsg(112, a1), a2), gsqr(a3)), a4), a5)),
120 : gmul(gmul(gmulsg(28, a1), gpowgs(a3, 4)), a5)), gmul(gmul(gmulsg(16, gsqr(a2)),
121 : gsqr(a3)), gsqr(a4))), gmul(gmulsg(5120, gpowgs(a0, 3)), gpowgs(a6, 3))),
122 : gmul(gmul(gmulsg(2544, gsqr(a0)), gsqr(a3)), gsqr(a6))), gmul(gmul(gmulsg(312,
123 : a0), gpowgs(a3, 4)), a6)), gmul(gmul(gmul(gmulsg(14336, gsqr(a0)), a2), a4),
124 : gsqr(a6))), gmul(gmul(gmul(gmulsg(1024, a0), gsqr(a2)), gsqr(a4)), a6)),
125 : gmul(gmul(gmul(gmulsg(2560, gsqr(a0)), a1), a5), gsqr(a6))),
126 : gmul(gmul(gmul(gmulsg(2240, a0), gsqr(a1)), gsqr(a5)), a6)),
127 : gmul(gmul(gmul(gmul(gmul(gmulsg(6528, a0), a1), a2), a4), a5), a6)),
128 : gmul(gmul(gmul(gmul(gmulsg(1568, a0), a2), gsqr(a3)), a4), a6)), -10));
129 : }
130 :
131 : static GEN
132 42 : igusaj8_fromj246(GEN j2, GEN j4, GEN j6)
133 42 : { return gmul2n(gsub(gmul(j2,j6), gsqr(j4)), -2); }
134 :
135 : static GEN
136 21 : igusaj8(GEN a0, GEN a1, GEN a2, GEN a3, GEN a4, GEN a5, GEN a6)
137 : {
138 21 : pari_sp av = avma;
139 21 : GEN j2 = igusaj2(a0,a1,a2,a3,a4,a5,a6);
140 21 : GEN j4 = igusaj4(a0,a1,a2,a3,a4,a5,a6);
141 21 : GEN j6 = igusaj6(a0,a1,a2,a3,a4,a5,a6);
142 21 : return gerepileupto(av, igusaj8_fromj246(j2,j4,j6));
143 : }
144 :
145 : static GEN
146 42 : igusaj10(GEN a0, GEN a1, GEN a2, GEN a3, GEN a4, GEN a5, GEN a6)
147 : {
148 42 : pari_sp av = avma;
149 42 : GEN polr = mkpoln(7, a0, a1, a2, a3, a4, a5, a6);
150 42 : GEN disc = RgX_disc(polr);
151 42 : GEN j10 = degpol(polr) < 6? gmul(gsqr(a1), disc): disc;
152 42 : return gerepileupto(av, gmul2n(j10, -12));
153 : }
154 :
155 : static GEN
156 21 : igusaall(GEN a0, GEN a1, GEN a2, GEN a3, GEN a4, GEN a5, GEN a6)
157 : {
158 21 : GEN j2, j4, j6, V = cgetg(6,t_VEC);
159 : pari_sp av;
160 21 : gel(V,1) = j2 = igusaj2(a0,a1,a2,a3,a4,a5,a6);
161 21 : gel(V,2) = j4 = igusaj4(a0,a1,a2,a3,a4,a5,a6);
162 21 : gel(V,3) = j6 = igusaj6(a0,a1,a2,a3,a4,a5,a6);
163 21 : av = avma;
164 21 : gel(V,4) = gerepileupto(av, igusaj8_fromj246(j2, j4, j6));
165 21 : gel(V,5) = igusaj10(a0,a1,a2,a3,a4,a5,a6);
166 21 : return V;
167 : }
168 :
169 : GEN
170 126 : genus2igusa(GEN P, long n)
171 : {
172 126 : pari_sp av = avma;
173 : GEN a0, a1, a2, a3, a4, a5, a6, r;
174 126 : if (typ(P) == t_VEC && lg(P) == 3)
175 42 : P = gadd(gmul2n(gel(P,1), 2), gsqr(gel(P,2)));
176 : else
177 84 : P = gmul2n(P, 2);
178 126 : if (typ(P)!=t_POL || degpol(P)> 6) pari_err_TYPE("genus2igusa",P);
179 126 : RgX_to_06(P, &a0,&a1,&a2,&a3,&a4,&a5,&a6);
180 126 : switch(n)
181 : {
182 21 : case 0: r = igusaall(a0,a1,a2,a3,a4,a5,a6); break;
183 21 : case 2: r = igusaj2(a0,a1,a2,a3,a4,a5,a6); break;
184 21 : case 4: r = igusaj4(a0,a1,a2,a3,a4,a5,a6); break;
185 21 : case 6: r = igusaj6(a0,a1,a2,a3,a4,a5,a6); break;
186 21 : case 8: r = igusaj8(a0,a1,a2,a3,a4,a5,a6); break;
187 21 : case 10:r = igusaj10(a0,a1,a2,a3,a4,a5,a6); break;
188 0 : default:
189 0 : pari_err_FLAG("genus2igusa");
190 : return NULL; /* LCOV_EXCL_LINE */
191 : }
192 126 : return gerepileupto(av, r);
193 : }
194 :
195 : /********************************************************************/
196 : /** **/
197 : /** A REDUCTION ALGORITHM "A LA TATE" FOR CURVES OF GENUS 2 **/
198 : /** **/
199 : /********************************************************************/
200 : /* Based on genus2reduction-0.3, http://www.math.u-bordeaux.fr/~liu/G2R/
201 : * by Qing Liu <liu@math.u-bordeaux.fr>
202 : * and Henri Cohen <cohen@math.u-bordeaux.fr>
203 :
204 : * Qing Liu: Modeles minimaux des courbes de genre deux
205 : * J. fuer die Reine und Angew. Math., 453 (1994), 137-164.
206 : * http://www.math.u-bordeaux.fr/~liu/articles/modregE.ps */
207 :
208 : /* some auxiliary polynomials, gp2c-generated */
209 :
210 : /*
211 : apol2(a0,a1,a2) = -5*a1^2+12*a0*a2;
212 : */
213 : static GEN
214 1400 : apol2(GEN a0, GEN a1, GEN a2)
215 : {
216 1400 : return gadd(gmulsg(-5, gsqr(a1)), gmul(gmulsg(12, a0), a2));
217 : }
218 :
219 : /*
220 : apol3(a0,a1,a2,a3) = 5*a1^3+9*a0*(-2*a1*a2+3*a0*a3);
221 : */
222 : static GEN
223 1400 : apol3(GEN a0, GEN a1, GEN a2, GEN a3)
224 : {
225 1400 : return gadd(gmulsg(5, gpowgs(a1, 3)), gmul(gmulsg(9, a0), gadd(gmul(gmulsg(-2, a1), a2), gmul(gmulsg(3, a0), a3))));
226 : }
227 :
228 : /*
229 : apol5(a0,a1,a2,a3,a4,a5) = a1^5+3*a0*(-2*a1^3*a2+9*a0*a1^2*a3-36*a0^2*a1*a4+108*a0^3*a5);
230 : */
231 : static GEN
232 1400 : apol5(GEN a0, GEN a1, GEN a2, GEN a3, GEN a4, GEN a5)
233 : {
234 1400 : return gadd(gpowgs(a1, 5), gmul(gmulsg(3, a0), gadd(gsub(gadd(gmul(gmulsg(-2, gpowgs(a1, 3)), a2), gmul(gmul(gmulsg(9, a0), gsqr(a1)), a3)), gmul(gmul(gmulsg(36, gsqr(a0)), a1), a4)), gmul(gmulsg(108, gpowgs(a0, 3)), a5))));
235 : }
236 :
237 : /*
238 : bpol2(a0,a1,a2,a3,a4) = 2*a2^2-5*a1*a3+10*a0*a4;
239 : */
240 : static GEN
241 1400 : bpol2(GEN a0, GEN a1, GEN a2, GEN a3, GEN a4)
242 : {
243 1400 : return gadd(gsub(gmulsg(2, gsqr(a2)), gmul(gmulsg(5, a1), a3)), gmul(gmulsg(10, a0), a4));
244 : }
245 :
246 : static const long VERYBIG = (1L<<20);
247 : static long
248 23695 : myval(GEN x, GEN p) { return signe(x)? Z_pval(x,p): VERYBIG; }
249 : static long
250 2982 : my3val(GEN x) { return signe(x)? Z_lval(x,3): VERYBIG; }
251 : /* b in Z[i], return v_3(b) */
252 : static long
253 1491 : myval_zi(GEN b) { return minss(my3val(real_i(b)), my3val(imag_i(b))); }
254 : /* b in Z[i, Y]/(Y^2-3), return v_Y(b) */
255 : static long
256 672 : myval_zi2(GEN b)
257 : {
258 : long v0, v1;
259 672 : b = lift_shallow(b);
260 672 : v0 = myval_zi(RgX_coeff(b,0));
261 672 : v1 = myval_zi(RgX_coeff(b,1));
262 672 : return minss(2*v0, 2*v1+1);
263 : }
264 :
265 : /* min(a,b,c) */
266 : static long
267 1932 : min3(long a, long b, long c)
268 : {
269 1932 : long m = a;
270 1932 : if (b < m) m = b;
271 1932 : if (c < m) m = c;
272 1932 : return m;
273 : }
274 :
275 : /* Vector of p-adic factors (over Q_p) to accuracy r of pol. */
276 : static GEN
277 119 : padicfactors(GEN pol, GEN p, long r) { return gel(factorpadic(pol,p,r),1); }
278 :
279 : /* x(1/t)*t^6, deg x <= 6 */
280 : static GEN
281 322 : RgX_recip6(GEN x)
282 : {
283 322 : long lx = lg(x), i, j;
284 322 : GEN y = cgetg(9, t_POL);
285 322 : y[1] = x[1];
286 2429 : for (i=8,j=2; j < lx; i--,j++) gel(y,i) = gel(x,j);
287 469 : for ( ; j < 9; i--,j++) gel(y,i) = gen_0;
288 322 : return normalizepol_lg(y, 9);
289 : }
290 : /* extract coefficients a0,...a3 of a polynomial a0 X^6 + ... + a6 */
291 : static void
292 1400 : RgX_to_03(GEN q, GEN *a0, GEN *a1, GEN *a2, GEN *a3)
293 : {
294 1400 : *a0 = gen_0;
295 1400 : *a1 = gen_0;
296 1400 : *a2 = gen_0;
297 1400 : *a3 = gen_0;
298 1400 : switch(degpol(q))
299 : {
300 938 : case 6: *a0 = gel(q,8); /*fall through*/
301 1400 : case 5: *a1 = gel(q,7); /*fall through*/
302 1400 : case 4: *a2 = gel(q,6); /*fall through*/
303 1400 : case 3: *a3 = gel(q,5); /*fall through*/
304 : }
305 1400 : }
306 :
307 : /* deg(H mod p) = 3, return v_p( disc(corresponding p-adic factor) ) */
308 : static long
309 14 : discpart(GEN H, GEN p, long prec)
310 : {
311 : GEN list, prod, dis;
312 : long i, j;
313 :
314 14 : if (degpol(FpX_red(H,p)) != 3)
315 : pari_err_BUG("discpart [must not reach]"); /* LCOV_EXCL_LINE */
316 14 : list = padicfactors(H,p,prec);
317 14 : prod = pol_1(varn(H));
318 56 : for(i = 1; i < lg(list); i++)
319 : {
320 42 : GEN t = gel(list,i);
321 84 : for(j = 3; j < lg(t); j++) /* include if nonconstant mod p */
322 70 : if (!valp(gel(t,j))) { prod = RgX_mul(prod,t); break; }
323 : }
324 14 : if (degpol(prod) != 3) pari_err_BUG("discpart [prod degree]");
325 14 : dis = RgX_disc(prod);
326 14 : return gequal0(dis)? prec+1: valp(dis);
327 : }
328 :
329 : /* B = b0 X^6 + ... + b6 a ZX, 0 <= j <= 3.
330 : * Let theta_j(H) := min { v_p(b_i) / (i - j), j < i <= 6 } >= 0.
331 : * Return 60 theta \in Z */
332 : static long
333 1932 : theta_j(GEN B, GEN p, long j)
334 : {
335 1932 : long i, t = VERYBIG;
336 9310 : for(i = 1+j; i <= 6; i++)
337 7378 : t = minss(t, myval(RgX_coeff(B,6-i), p) * (60 / (i-j)));
338 1932 : return t;
339 : }
340 : /* compute 6 * theta_3 for B in Z[i][X], p = 3 */
341 : static long
342 28 : theta_3_zi(GEN B)
343 : {
344 28 : long v2 = myval_zi(RgX_coeff(B,2));
345 28 : long v1 = myval_zi(RgX_coeff(B,1));
346 28 : long v0 = myval_zi(RgX_coeff(B,0));
347 28 : return min3(6*v2, 3*v1, 2*v0);
348 : }
349 : /* compute 6 * theta_3 for B in (Z[i,Y]/(Y^2-3))[X], p = 3 */
350 : static long
351 84 : theta_3_zi2(GEN B)
352 : {
353 84 : long v2 = myval_zi2(RgX_coeff(B,2));
354 84 : long v1 = myval_zi2(RgX_coeff(B,1));
355 84 : long v0 = myval_zi2(RgX_coeff(B,0));
356 84 : return min3(6*v2, 3*v1, 2*v0);
357 : }
358 :
359 : /* Set maxord to the maximal multiplicity of a factor. If there is at least
360 : * a triple root (=> maxord >= 3) return it, else return NULL */
361 : static GEN
362 812 : factmz(GEN Q, GEN p, long *maxord)
363 : {
364 812 : GEN z = FpX_factor_squarefree(Q, p);
365 812 : long m = lg(z)-1; /* maximal multiplicity */
366 812 : *maxord = m;
367 812 : return (m >= 3)? FpX_oneroot(gel(z,m), p): NULL;
368 : }
369 : static long
370 1743 : get_lambda(GEN H, GEN p)
371 : {
372 1743 : if (!dvdii(RgX_coeff(H,3), p)) return 3;
373 735 : if (!dvdii(RgX_coeff(H,4), p)) return 2;
374 560 : if (!dvdii(RgX_coeff(H,5), p)) return 1;
375 441 : if (!dvdii(RgX_coeff(H,6), p)) return 0;
376 63 : return -1;
377 : }
378 :
379 : /* H integral ZX of degree 5 or 6, p > 2. Modify until
380 : * y^2 = p^alpha H is minimal over Z_p, alpha = 0,1
381 : * Return [H,lambda,60*theta,alpha,quad,beta], where
382 : * - quad = 1 if H has a root of order 3 in F_p^2 \ F_p, 0 otherwise
383 : * - 0 <= lambda <= 3, index of a coefficient with valuation 0
384 : * - theta = theta_j(H(x + r), p, lambda), 60*theta in Z, where r is a root
385 : * of H mod p
386 : * - beta >= -1 s.t. H = p^n H0(r + p^beta * X) for some n, r in Z, where
387 : * H0 is the initial H or polrecip(H) */
388 : static GEN
389 1680 : polymini(GEN H, GEN p)
390 : {
391 1680 : long t60, alpha, lambda, quad = 0, beta = 0;
392 : GEN Hp;
393 :
394 1680 : alpha = ZX_pvalrem(H, p, &H) & 1;
395 1680 : lambda = get_lambda(H, p);
396 1680 : if (lambda < 0) { H = RgX_recip6(H); lambda = get_lambda(H, p); }
397 0 : for(;;)
398 : {
399 : for(;;)
400 252 : { /* lambda <= 3, t60 = 60*theta */
401 : GEN rac;
402 : long e, maxord;
403 1932 : t60 = theta_j(H,p,lambda); e = t60 / 60;
404 1932 : if (e)
405 : {
406 903 : GEN pe = powiu(p,e);
407 : /* H <- H(p^e X) / p^(e(6-lambda)) */
408 903 : H = ZX_unscale_divpow(H, pe, 6-lambda);
409 903 : alpha = (alpha + lambda*e)&1;
410 903 : beta += e;
411 903 : t60 -= 60*e;
412 : }
413 : /* 0 <= t < 60 */
414 1932 : Hp = FpX_red(H, p); if (t60) break;
415 :
416 805 : rac = factmz(Hp,p, &maxord);
417 805 : if (maxord <= 2)
418 : {
419 532 : if (degpol(Hp) <= 3) break;
420 119 : goto end;
421 : }
422 : /* maxord >= 3 */
423 273 : if (!rac) { quad = 1; goto end; }
424 252 : if (signe(rac)) H = ZX_translate(H, rac);
425 252 : lambda = 6 - maxord;
426 : }
427 1561 : if (lambda <= 2)
428 : {
429 630 : if (myval(RgX_coeff(H,2),p) > 1-alpha &&
430 518 : myval(RgX_coeff(H,1),p) > 2-alpha &&
431 413 : myval(RgX_coeff(H,0),p) > 3-alpha)
432 : {
433 0 : H = ZX_unscale(H, p);
434 0 : if (alpha) H = ZX_Z_mul(H, p);
435 0 : return polymini(H, p);
436 : }
437 630 : break;
438 : }
439 : /* lambda = 3 */
440 931 : if (alpha == 0) break;
441 399 : if (degpol(Hp) == 3)
442 : {
443 721 : if (myval(RgX_coeff(H,6),p) >= 3 &&
444 357 : myval(RgX_coeff(H,5),p) >= 2)
445 : { /* too close to root [Kodaira symbol for y^2 = p^alpha*H not
446 : implemented when alpha = 1]: go back one step */
447 357 : H = ZX_rescale(H, p); /* H(x/p)p^(deg H) */
448 357 : H = ZX_Z_divexact(H, powiu(p, degpol(H)-3)); /* H(x/p)p^3 */
449 357 : t60 += 60; alpha = 0; beta--;
450 : }
451 364 : break;
452 : }
453 35 : if (degpol(Hp) != 6) break;
454 7 : if (t60)
455 : {
456 : long m, maxord;
457 7 : GEN v, T, rac = factmz(RgX_mulXn(Hp, -3), p, &maxord);
458 14 : if (maxord <= 2) break;
459 7 : T = ZX_affine(H, p, rac); /* H(rac + px) */
460 7 : if (ZX_pval(T,p) < 3) break;
461 :
462 0 : H = ZX_Z_divexact(T, powiu(p,3));
463 0 : alpha = 0; beta++;
464 0 : v = FpX_factor_squarefree(FpX_red(H,p), p);
465 0 : m = lg(v)-1; /* maximal multiplicity */
466 0 : if (m > 1)
467 : {
468 0 : rac = FpX_oneroot(gel(v,m), p); /* v[m] is linear */
469 0 : H = ZX_translate(H,rac);
470 0 : t60 = theta_j(H,p,3);
471 0 : if (t60 < 60) break;
472 : }
473 : }
474 : }
475 1680 : end:
476 1680 : return mkvec2(H, mkvecsmall5(lambda,t60,alpha,quad,beta));
477 : }
478 :
479 : /* a in Q[i], return a^3 mod 3 */
480 : static GEN
481 14 : zi_pow3mod(GEN a)
482 : {
483 : GEN x, y;
484 14 : if (typ(a) != t_COMPLEX) return gmodgs(a,3);
485 7 : x = gmodgs(gel(a,1), 3);
486 7 : y = gmodgs(gel(a,2), 3);
487 7 : return mkcomplex(x, negi(y));
488 : }
489 : static GEN
490 21 : polymini_zi(GEN pol) /* polynome minimal dans Z[i] */
491 : {
492 21 : GEN polh, rac, a0, a1, a2, a3, a4, a5, a6, p = utoipos(3);
493 21 : long alpha, beta = 0, t6;
494 :
495 21 : alpha = ZX_pval(pol,p) & 1;
496 21 : polh = alpha? RgX_Rg_div(pol, p): pol;
497 21 : rac = mkcomplex(Fp_div(RgX_coeff(polh,3), RgX_coeff(polh,6), p), gen_1);
498 : for(;;)
499 7 : {
500 : long e;
501 28 : polh = RgX_translate(polh, rac);
502 28 : t6 = theta_3_zi(polh); e = t6 / 6;
503 28 : if (e)
504 : {
505 14 : GEN pe = powiu(p,e);
506 14 : polh = RgX_Rg_div(RgX_unscale(polh,pe), powiu(pe,3));
507 14 : alpha = (alpha+e)&1;
508 14 : t6 -= e * 6; beta += e;
509 : }
510 28 : RgX_to_06(polh, &a0,&a1,&a2,&a3,&a4,&a5,&a6);
511 28 : if (t6 || !myval_zi(a4) || !myval_zi(a5)) break;
512 7 : rac = zi_pow3mod(gdiv(a6, gneg(a3)));
513 : }
514 21 : if (alpha && myval_zi(a0) >= 3 && myval_zi(a1) >= 2 && myval_zi(a2) >= 1)
515 : {
516 14 : t6 += 6; beta--; alpha = 0;
517 : }
518 21 : if (alpha && beta >= 1) pari_err_BUG("quadratic");
519 21 : return mkvecsmall3(t6, alpha, beta);
520 : }
521 :
522 : /* pol is a ZX, minimal polynomial over Z_3[i,Y]/(Y^2-3) */
523 : static GEN
524 84 : polymini_zi2(GEN pol)
525 : {
526 : long alpha, beta, t6;
527 : GEN a0, a1, a2, a3, a4, a5, a6;
528 84 : GEN polh, rac, y = pol_x(fetch_var()), p = utoipos(3);
529 :
530 84 : if (ZX_pval(pol,p)) pari_err_BUG("polymini_zi2 [polynomial not minimal]");
531 84 : y = mkpolmod(y, gsubgs(gsqr(y), 3)); /* mod(y,y^2-3) */
532 84 : polh = gdivgs(RgX_unscale(pol, y),27); /* H(y*x) / 27 */
533 161 : if (myval_zi2(RgX_coeff(polh,4)) <= 0 ||
534 77 : myval_zi2(RgX_coeff(polh,2)) <= 0)
535 : {
536 7 : (void)delete_var();
537 7 : return mkvecsmall2(0,0);
538 : }
539 :
540 77 : if (myval_zi2(gsub(RgX_coeff(polh,6), RgX_coeff(polh,0))) > 0)
541 7 : rac = gen_I();
542 : else
543 70 : rac = gen_1;
544 77 : alpha = 0;
545 77 : beta = 0;
546 : for(;;)
547 7 : {
548 : long e;
549 84 : polh = RgX_translate(polh, rac);
550 84 : t6 = theta_3_zi2(polh); e = t6 / 6;
551 84 : if (e)
552 : {
553 77 : GEN pent = gpowgs(y, e);
554 77 : polh = RgX_Rg_div(RgX_unscale(polh, pent), gpowgs(pent,3));
555 77 : alpha = (alpha+e)&1;
556 77 : t6 -= 6*e; beta += e;
557 : }
558 84 : RgX_to_06(polh, &a0,&a1,&a2,&a3,&a4,&a5,&a6);
559 84 : if (t6 || !myval_zi2(a4) || !myval_zi2(a5)) break;
560 7 : a3 = liftpol_shallow(a3); if (typ(a3)==t_POL) a3 = RgX_coeff(a3,0);
561 7 : a6 = liftpol_shallow(a6); if (typ(a6)==t_POL) a6 = RgX_coeff(a6,0);
562 7 : rac = zi_pow3mod(gdiv(a6,gneg(a3)));
563 : }
564 77 : if (alpha)
565 : {
566 42 : if (myval_zi2(a0) < 3 || myval_zi2(a1) < 2 || myval_zi2(a2) < 1)
567 0 : pari_err_BUG("polymini_zi2 [alpha]");
568 42 : t6 += 6; beta--;
569 : }
570 77 : (void)delete_var();
571 77 : if (odd(beta)) pari_err_BUG("quartic [type over Z[i] must be [K-K-(2*m)]]");
572 77 : return mkvecsmall2(t6, beta);
573 : }
574 :
575 : struct igusa {
576 : GEN j2, i4, j4, j6, j8, j10, i12;
577 : GEN a0, A2, A3, A5, B2;
578 : };
579 : struct igusa_p {
580 : long eps, tt, r1, r2, tame;
581 : GEN p, stable, val, neron;
582 : const char *type;
583 : };
584 :
585 : /* initialize Ip */
586 : static void
587 1442 : stable_reduction(struct igusa *I, struct igusa_p *Ip, GEN p)
588 : {
589 : static const long d[9] = { 0,60,30,30,20,15,12,10 }; /* 120 / deg(X) */
590 1442 : GEN j2 = I->j2, i4 = I->i4, j4 = I->j4, j6 = I->j6, j8 = I->j8;
591 1442 : GEN val, J, v, Ieps, j10 = I->j10, i12 = I->i12;
592 : long s, r1, r2, r3, r4, i, eps;
593 :
594 1442 : Ip->tame = 0;
595 1442 : Ip->neron = NULL;
596 1442 : Ip->type = NULL;
597 1442 : Ip->p = p;
598 1442 : Ip->val = val = cgetg(9, t_VECSMALL);
599 1442 : val[1] = myval(j2,p);
600 1442 : val[2] = myval(j4,p);
601 1442 : val[3] = myval(i4,p);
602 1442 : val[4] = myval(j6,p);
603 1442 : val[5] = myval(j8,p);
604 1442 : val[6] = myval(j10,p);
605 1442 : val[7] = myval(i12,p);
606 1442 : switch(itos_or_0(p))
607 : {
608 21 : case 2: eps = 4; val[8] = val[5]; Ieps = j8; break;
609 476 : case 3: eps = 3; val[8] = val[4]; Ieps = j6; break;
610 945 : default: eps = 1; val[8] = val[1]; Ieps = gdivgs(j2,12); break;
611 : }
612 :
613 1442 : v = cgetg(8,t_VECSMALL);
614 11536 : for(i = 1; i <= 7; i++) v[i] = val[i] * d[i];
615 1442 : s = vecsmall_min(v);
616 1442 : Ip->eps = eps;
617 :
618 1442 : r1 = 3*eps*val[3];
619 1442 : r3 = eps*val[6] + val[8];
620 1442 : r2 = eps*val[7];
621 1442 : r4 = min3(r1, r2, r3);
622 :
623 : /* s = max(v_p(X) / deg(X)) */
624 1442 : J = cgetg(1, t_VEC);
625 1442 : if (s == v[6])
626 161 : Ip->tt = 1;
627 1281 : else if (s == v[7])
628 : {
629 119 : J = mkvec( Fp_to_mod(gmod(gdiv(gpowgs(i4,3),i12), p), p) );
630 119 : Ip->tt = 2;
631 : }
632 1162 : else if (s == v[3])
633 210 : Ip->tt = (val[2] == val[3] || 2*val[4] == 3*val[3])? 3: 4;
634 952 : else if (r3 == r4)
635 : {
636 560 : GEN a,b, P, sj, pj, t = gmul(gpowgs(j10,eps),Ieps);
637 560 : sj = gaddsg(1728, gdiv(gpowgs(i12,eps), t));
638 560 : pj = gdiv(gpowgs(i4,3*eps), t);
639 560 : a = gmod(sj, p);
640 560 : b = gmod(pj, p);
641 560 : P = mkpoln(3, gen_1, Fp_neg(a,p), b, 0); /* X^2 - SX + P: roots j1,j2 */
642 560 : J = FpX_roots(P, p);
643 560 : switch(lg(J)-1)
644 : {
645 0 : case 0:
646 0 : P = FpX_to_mod(P, p);
647 0 : a = FpX_to_mod(pol_x(0), p);
648 0 : b = FpX_to_mod(deg1pol_shallow(b, gen_m1,0), p);
649 0 : J = mkvec2(mkpolmod(a,P), mkpolmod(b,P)); break;
650 378 : case 1:
651 378 : a = Fp_to_mod(gel(J,1), p);
652 378 : J = mkvec2(a, a); break;
653 182 : case 2:
654 182 : settyp(J, t_VEC);
655 182 : J = FpV_to_mod(J, p); break;
656 : }
657 560 : Ip->tt = 5;
658 : }
659 392 : else if (r2 == r4)
660 : {
661 280 : J = mkvec( Fp_to_mod(gmod(gdiv(gpowgs(i4,3),i12), p), p) );
662 280 : Ip->tt = 6;
663 : }
664 : else
665 112 : Ip->tt = 7; /* r1 == r4 */
666 1442 : Ip->stable = mkvec2(stoi(Ip->tt), J);
667 1442 : }
668 :
669 : struct red {
670 : const char *t, *pages;
671 : double tnum;
672 : GEN g;
673 : };
674 :
675 : /* destroy v */
676 : static GEN
677 1421 : zv_snf(GEN v)
678 : {
679 1421 : long i, l = lg(v);
680 3143 : for (i = 1; i < l; i++)
681 : {
682 1722 : long j, a = v[i];
683 2485 : for (j = i+1; j < l; j++)
684 : {
685 763 : long b = v[j], d = ugcd(a,b);
686 763 : v[i] = a = a*(b/d);
687 763 : v[j] = d;
688 : }
689 : }
690 1505 : for (i = l-1; i > 0; i--)
691 1204 : if (v[i] != 1) { setlg(v, i+1); break; }
692 1421 : return zv_to_ZV(v);
693 : }
694 :
695 : static GEN
696 1344 : cyclic(long n)
697 1344 : { return (n <= 1)? cgetg(1, t_VECSMALL): mkvecsmall(n); }
698 : static GEN
699 336 : dicyclic(long a, long b)
700 : {
701 : long d;
702 336 : if (!a) a = 1;
703 336 : if (!b) b = 1;
704 336 : if (a < b) lswap(a,b);
705 336 : d = ugcd(a,b);
706 336 : if (d == 1) return cyclic(a*b);
707 280 : return mkvecsmall2(a*b/d, d);
708 : }
709 : /* Z/2xZ/2, resp Z/4 for n even, resp. odd */
710 : static GEN
711 280 : groupH(long n) { return odd(n)? cyclic(4): dicyclic(2,2); }
712 :
713 : static long
714 252 : get_red(struct red *S, struct igusa_p *Ip, GEN polh, GEN p, long alpha, long r)
715 : {
716 252 : GEN val = Ip->val;
717 : long indice;
718 252 : switch(r)
719 : {
720 42 : case 0:
721 42 : indice = FpX_is_squarefree(FpX_red(polh,p), p)
722 : ? 0
723 42 : : val[6] - val[7] + val[8]/Ip->eps;
724 42 : S->t = stack_sprintf("I{%ld}", indice);
725 42 : S->tnum = 1;
726 42 : S->pages = "159-177";
727 42 : S->g = cyclic(indice);
728 42 : return indice ? indice: 1;
729 35 : case 6:
730 35 : if (alpha == 0) polh = ZX_unscale_divpow(polh, p, 3); /* H(px) /p^3 */
731 35 : indice = FpX_is_squarefree(FpX_red(polh,p), p)
732 : ? 0
733 35 : : val[6] - val[7] + val[8]/Ip->eps;
734 35 : S->t = stack_sprintf("I*{%ld}", indice);
735 35 : S->tnum = 1.5;
736 35 : S->pages = "159-177";
737 35 : S->g = groupH(indice);
738 35 : return indice + 5;
739 21 : case 3:
740 21 : S->t = "III";
741 21 : S->tnum = 3;
742 21 : S->pages = "161-177";
743 21 : S->g = cyclic(2);
744 21 : return 2;
745 35 : case 9:
746 35 : S->t = "III*";
747 35 : S->tnum = 3.5;
748 35 : S->pages = "162-177";
749 35 : S->g = cyclic(2);
750 35 : return 8;
751 28 : case 2:
752 28 : S->t = "II";
753 28 : S->tnum = 2;
754 28 : S->pages = "159-174";
755 28 : S->g = cyclic(1);
756 28 : return 1;
757 56 : case 8:
758 56 : S->t = "IV*";
759 56 : S->tnum = 4.5;
760 56 : S->pages = "160-175";
761 56 : S->g = cyclic(3);
762 56 : return 7;
763 21 : case 4:
764 21 : S->t = "IV";
765 21 : S->tnum = 4;
766 21 : S->pages = "160-174";
767 21 : S->g = cyclic(3);
768 21 : return 3;
769 14 : case 10:
770 14 : S->t = "II*";
771 14 : S->tnum = 2.5;
772 14 : S->pages = "160-174";
773 14 : S->g = cyclic(1);
774 14 : return 9;
775 0 : default: pari_err_BUG("get_red [type]");
776 0 : S->t = "";
777 0 : S->tnum = 0;
778 0 : S->pages = ""; /* gcc -Wall */
779 0 : S->g = NULL;
780 : return -1; /*LCOV_EXCL_LINE*/
781 : }
782 : }
783 :
784 : /* reduce a/b; assume b > 0 */
785 : static void
786 1330 : ssQ_red(long a, long b, long *n, long *d)
787 : {
788 1330 : long g = ugcd(labs(a), b);
789 1330 : if (g > 1) { a /= g; b /= g; }
790 1330 : *n = a; *d = b;
791 1330 : }
792 : /* denom(a/b); assume b > 0 */
793 : static long
794 28 : ssQ_denom(long a, long b)
795 : {
796 28 : long g = ugcd(labs(a), b);
797 28 : return g == 1? b: b / g;
798 : }
799 : /* n = lcm(d, denom(a/b)); r = (a/b * n mod n); assume b > 0 and d > 0 */
800 : static void
801 455 : get_nr(long d, long a, long b, long *n, long *r)
802 : {
803 : long c, A, B;
804 455 : ssQ_red(a, b, &A,&B);
805 455 : c = d / ugcd(d, B);
806 455 : *n = B * c;
807 455 : *r = umodsu(A * c, *n);
808 455 : }
809 : /* n = lcm(denom(a/b), denom(c/d)); r = (a/b * n mod n); q = (c/d * n mod n);
810 : * assume b > 0 and d > 0 */
811 : static void
812 154 : get_nrq(long a, long b, long c, long d, long *n, long *r, long *q)
813 : {
814 : long g, A, B, C, D;
815 154 : ssQ_red(a, b, &A,&B);
816 154 : ssQ_red(c, d, &C,&D);
817 154 : g = ugcd(B,D);
818 154 : *n = B * (D/g);
819 154 : *r = umodsu(A * (D/g), *n);
820 154 : *q = umodsu(C * (B/g), *n);
821 154 : }
822 :
823 : /* Ip->tt = 1 */
824 : static long
825 28 : tame_1(struct igusa *I, struct igusa_p *Ip)
826 : {
827 28 : GEN p = Ip->p, val = Ip->val;
828 28 : long condp = -1, va0, va5, r, n;
829 28 : va0 = myval(I->a0,p);
830 28 : va5 = myval(I->A5,p);
831 28 : if (!gequal0(I->A5) && 20*va0+val[6] > 6*va5)
832 21 : get_nr(ssQ_denom(5*val[6]-6*va5, 40), val[6]-2*va5, 20, &n,&r);
833 : else
834 7 : get_nr(ssQ_denom(5*va0-val[6], 10), 10*va0-val[6], 30, &n,&r);
835 28 : switch(n)
836 : {
837 0 : case 1:
838 0 : condp = 0;
839 0 : Ip->type = "[I{0-0-0}] page 155";
840 0 : Ip->neron = cyclic(1); break;
841 21 : case 2:
842 21 : switch(r)
843 : {
844 14 : case 0:
845 14 : condp = 4;
846 14 : Ip->type = "[I*{0-0-0}] page 155";
847 14 : Ip->neron = mkvecsmall4(2,2,2,2); break;
848 7 : case 1:
849 7 : condp = 2;
850 7 : Ip->type = "[II] page 155";
851 7 : Ip->neron = cyclic(1); break;
852 0 : default: pari_err_BUG("tame_1 [bug1]");
853 : }
854 21 : break;
855 7 : case 4:
856 7 : condp = 4;
857 7 : Ip->type = "[VI] page 156";
858 7 : Ip->neron = dicyclic(2,2); break;
859 0 : default: pari_err_BUG("tame_1 [bug8]");
860 : }
861 28 : return condp;
862 : }
863 :
864 : /* (4.2) */
865 : static long
866 203 : tame_234_init(struct igusa *I, struct igusa_p *Ip, long *n, long *q, long *r)
867 : {
868 203 : long va0, va5, vb2, v12 = -1, flc = 1;
869 203 : GEN p = Ip->p;
870 203 : switch(Ip->tt)
871 : {
872 91 : case 2: v12 = myval(I->i12, Ip->p); break;
873 56 : case 3: v12 = 3*myval(I->i4, Ip->p); break;
874 56 : case 4: v12 = 6*myval(I->j2, Ip->p); break;
875 : }
876 203 : va0 = myval(I->a0,p);
877 203 : va5 = myval(I->A5,p);
878 203 : vb2 = myval(I->B2,p);
879 203 : if (9*vb2 >= 6*va0+v12 && 36*va5 >= 120*va0+5*v12)
880 : {
881 42 : get_nrq(12*va0-v12,36, 6*va0-v12,12, n, r, q);
882 : }
883 161 : else if (120*va0+5*v12 > 36*va5 && 60*vb2 >= 12*va5+5*v12)
884 : {
885 49 : ssQ_red(36*va5-25*v12,240, q,n);
886 49 : *r = umodsu(-2* *q, *n);
887 : }
888 : else /* 6*va0+v12 > 9*vb2 && 12*va5+5*v12 > 60*vb2 */
889 : {
890 112 : get_nrq(v12-6*vb2,12, v12-9*vb2,12, n,r,q);
891 112 : flc = 0;
892 : }
893 203 : return flc;
894 : }
895 :
896 : /* Ip->tt = 2 */
897 : static long
898 91 : tame_2(struct igusa *I, struct igusa_p *Ip)
899 : {
900 91 : long condp = -1, d, n, q, r;
901 91 : GEN val = Ip->val;
902 91 : (void)tame_234_init(I, Ip, &n, &q, &r);
903 91 : d = n * (6*val[6]-5*val[7]) / 6;
904 91 : switch(n)
905 : {
906 7 : case 1: condp = 1;
907 7 : Ip->type = stack_sprintf("[I{%ld-0-0}] page 170", d);
908 7 : Ip->neron = cyclic(d); break;
909 21 : case 2:
910 21 : switch(r)
911 : {
912 7 : case 0: condp = 4;
913 7 : Ip->type = stack_sprintf("[I*{%ld-0-0}] page 171",d/2);
914 7 : Ip->neron = shallowconcat(dicyclic(2,2),groupH(d/2)); break;
915 14 : case 1:
916 14 : switch(q)
917 : {
918 7 : case 0: condp = 2;
919 7 : Ip->type = stack_sprintf("[II*{%ld-0}] page 172",d/2);
920 7 : Ip->neron = cyclic(1); break;
921 7 : case 1: condp = 3;
922 7 : Ip->type = stack_sprintf("[II{%ld-0}] page 171",d/2);
923 7 : Ip->neron = cyclic(2*d); break;
924 0 : default: pari_err_BUG("tame2 [bug10]");
925 : }
926 14 : break;
927 0 : default: pari_err_BUG("tame2 [bug11]");
928 : }
929 21 : break;
930 14 : case 3: condp = 3;
931 14 : Ip->neron = cyclic(d);
932 14 : switch(r)
933 : {
934 7 : case 1:
935 7 : Ip->type = stack_sprintf("[II{%ld}-IV] page 175", (d-2)/3);
936 7 : break;
937 7 : case 2:
938 7 : Ip->type = stack_sprintf("[II{%ld}-IV*] page 175", (d-1)/3);
939 7 : break;
940 0 : default: pari_err_BUG("tame2 [bug12]");
941 : }
942 14 : break;
943 42 : case 4:
944 42 : switch(r)
945 : {
946 21 : case 1:
947 21 : switch(q)
948 : {
949 14 : case 1: condp = 3;
950 14 : Ip->type = stack_sprintf("[II{%ld}-III] page 177",(d-2)/4);
951 14 : Ip->neron = cyclic(d/2); break;
952 7 : case 3: condp = 4;
953 7 : Ip->type = stack_sprintf("[II*{%ld}-III*] page 178",(d-2)/4);
954 7 : Ip->neron = cyclic(8); break;
955 0 : default: pari_err_BUG("tame2 [bug13]");
956 : }
957 21 : break;
958 21 : case 3:
959 21 : switch(q)
960 : {
961 7 : case 1: condp = 4;
962 7 : Ip->type = stack_sprintf("[II*{%ld}-III] page 178",(d-2)/4);
963 7 : Ip->neron = cyclic(8); break;
964 14 : case 3: condp = 3;
965 14 : Ip->type = stack_sprintf("[II{%ld}-III*] page 178",(d-2)/4);
966 14 : Ip->neron = cyclic(d/2); break;
967 0 : default: pari_err_BUG("tame2 [bug14]");
968 : }
969 21 : break;
970 0 : default: pari_err_BUG("tame2 [bug15]");
971 : }
972 42 : break;
973 7 : case 6:
974 7 : switch(r)
975 : {
976 7 : case 2: condp = 4;
977 7 : Ip->type = stack_sprintf("[II*-II*{%ld}] page 176", (d-4)/6);
978 7 : Ip->neron = groupH((d+2)/6); break;
979 0 : case 4: condp = 4;
980 0 : Ip->type = stack_sprintf("[II-II*{%ld}] page 176", (d-2)/6);
981 0 : Ip->neron = groupH((d+4)/6); break;
982 0 : default: pari_err_BUG("tame2 [bug16]");
983 : }
984 7 : break;
985 0 : default: pari_err_BUG("tame2 [bug17]");
986 : }
987 91 : return condp;
988 : }
989 :
990 : /* Ip->tt = 3 */
991 : static long
992 56 : tame_3(struct igusa *I, struct igusa_p *Ip)
993 : {
994 56 : long condp = -1, n, q, r, va5, d1, d2;
995 56 : long flc = tame_234_init(I, Ip, &n, &q, &r);
996 56 : GEN val = Ip->val;
997 :
998 56 : va5 = 2*val[6]-5*val[3];
999 56 : d1 = minss(n * (val[7]-3*val[3]), n * va5 / 4);
1000 56 : d2 = n * va5 / 2 - d1;
1001 56 : switch(n)
1002 : {
1003 14 : case 1: condp = 2;
1004 14 : Ip->type = stack_sprintf("[I{%ld-%ld-0}] page 179", d1,d2);
1005 14 : Ip->neron = dicyclic(d1,d2); break;
1006 28 : case 2:
1007 28 : switch(r)
1008 : {
1009 14 : case 0: condp = 4;
1010 14 : Ip->type = stack_sprintf("[I*{%ld-%ld-0}] page 180", d1/2,d2/2);
1011 14 : Ip->neron = shallowconcat(groupH(d1/2),groupH(d2/2)); break;
1012 14 : case 1: condp = 3;
1013 14 : if (flc)
1014 : {
1015 14 : Ip->type = stack_sprintf("[2I{%ld}-0] page 181", d1);
1016 14 : Ip->neron = cyclic(d1);
1017 : }
1018 : else
1019 : { /* FIXME: "or" same with d1<->d2 */
1020 0 : Ip->type = stack_sprintf("[II{%ld-%ld}] page 182",d1/2,d2/2);
1021 0 : Ip->neron = ((d1*d2-4)&7)? cyclic(2*d1): dicyclic(d1,2);
1022 : }
1023 14 : break;
1024 0 : default: pari_err_BUG("tame3 [bug20]");
1025 : }
1026 28 : break;
1027 14 : case 4: condp = 4;
1028 14 : Ip->type = stack_sprintf("[III{%ld}] page 182", d1/2);
1029 14 : Ip->neron = groupH(d1/2); break;
1030 0 : default: pari_err_BUG("tame3 [bug21]");
1031 : }
1032 56 : return condp;
1033 : }
1034 :
1035 : /* Ip->tt = 4 */
1036 : static long
1037 56 : tame_4(struct igusa *I, struct igusa_p *Ip)
1038 : {
1039 56 : long condp = -1, d1,d2,d3, f1,f2, g, h, n, q, r, vl,vn,vm, e1,e2,e3;
1040 56 : GEN val = Ip->val;
1041 56 : (void)tame_234_init(I, Ip, &n, &q, &r);
1042 56 : vl = val[6]-5*val[1];
1043 56 : vn = val[7]-6*val[1];
1044 56 : vm = val[2]-2*val[1]; /* all >= 0 */
1045 56 : e1 = min3(2*vl, 3*vn, 6*vm);
1046 56 : e2 = minss(6*vl - e1, 12*vn - 2*e1); /* >= 0 */
1047 56 : e3 = 12*vl - (2*e1+e2); /* >= 0 */
1048 56 : d1 = e1*n / 6;
1049 56 : d2 = e2*n / 12;
1050 56 : d3 = e3*n / 12;
1051 56 : g = d1*d2 + d1*d3 + d2*d3;
1052 56 : h = ugcd(ugcd(d1,d2),d3);
1053 56 : switch(n)
1054 : {
1055 7 : case 1: condp = 2;
1056 7 : Ip->type = stack_sprintf("[I{%ld-%ld-%ld}] page 182",d1,d2,d3);
1057 7 : Ip->neron = dicyclic(h,g/h); break;
1058 49 : case 2:
1059 49 : switch(r)
1060 : {
1061 7 : case 0: condp = 4;
1062 7 : Ip->type = stack_sprintf("[I*{%ld-%ld-%ld}] page 183",d1/2,d2/2,d3/2);
1063 7 : Ip->neron = shallowconcat(groupH(g/4), groupH(2-((h&2)>>1))); break;
1064 42 : case 1:
1065 42 : if (d1 == d2 || d1 == d3) f2 = d1;
1066 0 : else if (d2 == d3) f2 = d2;
1067 : else {
1068 0 : pari_err_BUG("tame4 [bug23]");
1069 : return -1; /*LCOV_EXCL_LINE*/
1070 : }
1071 42 : f1 = d1+d2+d3-2*f2;
1072 42 : switch(q)
1073 : {
1074 14 : case 0: condp = 3;
1075 14 : Ip->type = stack_sprintf("[II*{%ld-%ld}] page 184", f1/2,f2);
1076 14 : Ip->neron = cyclic(f2); break;
1077 28 : case 1: condp = 3;
1078 28 : Ip->type = stack_sprintf("[II{%ld-%ld}] page 183", f1/2,f2);
1079 28 : Ip->neron = cyclic(2*f1+f2); break;
1080 0 : default: pari_err_BUG("tame4 [bug24]");
1081 : }
1082 42 : break;
1083 0 : default: pari_err_BUG("tame4 [bug25]");
1084 : }
1085 49 : break;
1086 0 : case 3: condp = 4;
1087 0 : Ip->type = stack_sprintf("[III{%ld}] page 184",d1);
1088 0 : Ip->neron = (d1%3)? cyclic(9): dicyclic(3,3); break;
1089 0 : case 6: condp = 4;
1090 0 : Ip->type = stack_sprintf("[III*{%ld}] page 184",d1/2);
1091 0 : Ip->neron = cyclic(1); break;
1092 0 : default: pari_err_BUG("tame4 [bug26]");
1093 : }
1094 56 : return condp;
1095 : }
1096 :
1097 : /* p = 3 */
1098 : static void
1099 91 : tame_567_init_3(struct igusa_p *Ip, long dk,
1100 : long *pd, long *pn, long *pdm, long *pr)
1101 : {
1102 91 : long n = 1 + Ip->r1/6;
1103 91 : *pd = n * dk / 36; /* / (12*Ip->eps) */
1104 91 : *pn = n;
1105 91 : *pr = -1; /* unused */
1106 91 : *pdm = 0;
1107 91 : }
1108 :
1109 : /* (4.3) */
1110 : static void
1111 609 : tame_567_init(struct igusa *I, struct igusa_p *Ip, long dk,
1112 : long *pd, long *pn, long *pdm, long *pr)
1113 : {
1114 : long ndk, ddk;
1115 609 : GEN p = Ip->p, val = Ip->val;
1116 :
1117 609 : if (equaliu(p,3)) { tame_567_init_3(Ip, dk, pd, pn, pdm, pr); return; }
1118 : /* assume p > 3, Ip->eps = 1 */
1119 518 : ssQ_red(dk, 12, &ndk, &ddk);
1120 518 : if (! odd(val[8]))
1121 : {
1122 427 : long va0 = myval(I->a0,p), va2 = myval(I->A2,p), va3 = myval(I->A3,p);
1123 427 : long va5 = myval(I->A5,p), vb2 = myval(I->B2,p);
1124 427 : long v1 = 2*va3-4*va0-val[1], v2 = 6*va5-20*va0-5*val[1];
1125 427 : long v3 = 3*vb2-2*va0-2*val[1], v4 = 10*vb2-2*va5-5*val[1];
1126 427 : if (v3 >= 0 && v2 >= 0 && v1 >= 0)
1127 : {
1128 245 : if (v1==0 || v2==0) get_nr(ddk, va0+val[1], 6,pn,pr); /* Prop 4.3.1 (a) */
1129 : else
1130 : { /* Prop 4.3.1 (d) */
1131 231 : long v5 = myval(subii(mulii(I->A2,I->A3),mului(3,I->A5)),p);
1132 231 : if (gequal0(I->A2)) pari_err_BUG("tame567 [bug27]");
1133 231 : get_nr(ddk, 12*va0 + min3(dk, 6*va3-9*va2, 4*v5 - 10*va2), 24, pn,pr);
1134 : }
1135 : }
1136 182 : else if (v2 < 0 && v4 >= 0)
1137 182 : get_nr(ddk, 2*va5+val[1], 8, pn,pr); /* Prop 4.3.1 (b) */
1138 : else /* (v3 < 0 && v4 < 0) */
1139 0 : get_nr(ddk, vb2, 4, pn,pr); /* Prop 4.3.1 (c) */
1140 427 : *pd = (*pn/ddk) * ndk;
1141 : }
1142 : else
1143 : {
1144 91 : *pr = ndk;
1145 91 : *pn = 2*ddk;
1146 91 : *pd = 2*ndk;
1147 : }
1148 518 : *pdm = umodsu(*pd, *pn);
1149 : }
1150 :
1151 : static long
1152 329 : tame_5(struct igusa *I, struct igusa_p *Ip)
1153 : {
1154 329 : long condp = -1, d, n, dm, r, dk;
1155 329 : GEN val = Ip->val;
1156 :
1157 329 : dk = Ip->eps*val[6]-5*val[8];
1158 329 : tame_567_init(I, Ip, dk, &d, &n, &dm, &r);
1159 329 : if (! odd(val[8]))
1160 : {
1161 266 : switch(n)
1162 : {
1163 7 : case 1: condp = 0;
1164 7 : Ip->type = stack_sprintf("[I{0}-I{0}-%ld] page 158", d);
1165 7 : Ip->neron = cyclic(1); break;
1166 14 : case 2:
1167 14 : switch(dm)
1168 : {
1169 7 : case 0: condp = 4;
1170 7 : Ip->type = stack_sprintf("[I*{0}-I*{0}-%ld] page 158",(d-2)/2);
1171 7 : Ip->neron = mkvecsmall4(2,2,2,2); break;
1172 7 : case 1: condp = 2;
1173 7 : Ip->type = stack_sprintf("[I{0}-I*{0}-%ld] page 159",(d-1)/2);
1174 7 : Ip->neron = dicyclic(2,2); break;
1175 : }
1176 14 : break;
1177 35 : case 3:
1178 35 : switch(dm)
1179 : {
1180 7 : case 0: condp = 4;
1181 7 : Ip->type = stack_sprintf("[IV-IV*-%ld] page 165",(d-3)/3);
1182 7 : Ip->neron = dicyclic(3,3); break;
1183 14 : case 1:
1184 14 : switch(r)
1185 : {
1186 7 : case 0: case 1: condp = 2;
1187 7 : Ip->type = stack_sprintf("[I{0}-IV-%ld] page 160",(d-1)/3);
1188 7 : Ip->neron = cyclic(3); break;
1189 7 : case 2: condp = 4;
1190 7 : Ip->type = stack_sprintf("[IV*-IV*-%ld] page 166",(d-4)/3);
1191 7 : Ip->neron = dicyclic(3,3); break;
1192 : }
1193 14 : break;
1194 14 : case 2:
1195 14 : switch(r)
1196 : {
1197 7 : case 0: case 2: condp = 2;
1198 7 : Ip->type = stack_sprintf("[I{0}-IV*-%ld] page 160",(d-2)/3);
1199 7 : Ip->neron = cyclic(3); break;
1200 7 : case 1: condp = 4;
1201 7 : Ip->type = stack_sprintf("[IV-IV-%ld] page 165",(d-2)/3);
1202 7 : Ip->neron = dicyclic(3,3); break;
1203 : }
1204 14 : break;
1205 : }
1206 35 : break;
1207 49 : case 4:
1208 49 : switch(dm)
1209 : {
1210 7 : case 0: condp = 4;
1211 7 : Ip->type = stack_sprintf("[III-III*-%ld] page 169",(d-4)/4);
1212 7 : Ip->neron = dicyclic(2,2); break;
1213 14 : case 1:
1214 14 : switch(r)
1215 : {
1216 7 : case 0: case 1: condp = 2;
1217 7 : Ip->type = stack_sprintf("[I{0}-III-%ld] page 161",(d-1)/4);
1218 7 : Ip->neron = cyclic(2); break;
1219 7 : case 2: case 3: condp = 4;
1220 7 : Ip->type = stack_sprintf("[I*{0}-III*-%ld] page 162",(d-5)/4);
1221 7 : Ip->neron = mkvecsmall3(2,2,2); break;
1222 : }
1223 14 : break;
1224 14 : case 2: condp = 4;
1225 14 : Ip->neron = dicyclic(2,2);
1226 14 : switch(r)
1227 : {
1228 7 : case 1:
1229 7 : Ip->type = stack_sprintf("[III-III-%ld] page 169",(d-2)/4);
1230 7 : break;
1231 7 : case 3:
1232 7 : Ip->type = stack_sprintf("[III*-III*-%ld] page 169",(d-6)/4);
1233 7 : break;
1234 0 : default: pari_err_BUG("tame5 [bug29]");
1235 : }
1236 14 : break;
1237 14 : case 3:
1238 14 : switch(r)
1239 : {
1240 7 : case 0: case 3: condp = 2;
1241 7 : Ip->type = stack_sprintf("[I{0}-III*-%ld] page 162",(d-3)/4);
1242 7 : Ip->neron = cyclic(2); break;
1243 7 : case 1: case 2: condp = 4;
1244 7 : Ip->type = stack_sprintf("[I*{0}-III-%ld] page 162",(d-3)/4);
1245 7 : Ip->neron = mkvecsmall3(2,2,2); break;
1246 : }
1247 14 : break;
1248 : }
1249 49 : break;
1250 105 : case 6:
1251 105 : switch(dm)
1252 : {
1253 7 : case 0: condp = 4;
1254 7 : Ip->type = stack_sprintf("[II-II*-%ld] page 163",(d-6)/6);
1255 7 : Ip->neron = cyclic(1); break;
1256 21 : case 1:
1257 21 : switch(r)
1258 : {
1259 7 : case 0: case 1: condp = 2;
1260 7 : Ip->type = stack_sprintf("[I{0}-II-%ld] page 159",(d-1)/6);
1261 7 : Ip->neron = cyclic(1); break;
1262 7 : case 2: case 5: condp = 4;
1263 7 : Ip->type = stack_sprintf("[II*-IV-%ld] page 164",(d-7)/6);
1264 7 : Ip->neron = cyclic(3); break;
1265 7 : case 3: case 4: condp = 4;
1266 7 : Ip->type = stack_sprintf("[I*{0}-IV*-%ld] page 161",(d-7)/6);
1267 7 : Ip->neron = mkvecsmall2(6,2); break;
1268 : }
1269 21 : break;
1270 21 : case 2:
1271 21 : switch(r)
1272 : {
1273 14 : case 1: condp = 4;
1274 14 : Ip->type = stack_sprintf("[II-II-%ld] page 163",(d-2)/6);
1275 14 : Ip->neron = cyclic(1); break;
1276 7 : case 3: case 5: condp = 4;
1277 7 : Ip->type = stack_sprintf("[I*{0}-II*-%ld] page 160",(d-8)/6);
1278 7 : Ip->neron = dicyclic(2,2); break;
1279 0 : default: pari_err_BUG("tame5 [bug30]");
1280 : }
1281 21 : break;
1282 14 : case 3:
1283 14 : Ip->neron = cyclic(3);
1284 14 : switch(r)
1285 : {
1286 7 : case 1: case 2: condp = 4;
1287 7 : Ip->type = stack_sprintf("[II-IV-%ld] page 164",(d-3)/6);
1288 7 : break;
1289 7 : case 4: case 5: condp = 4;
1290 7 : Ip->type = stack_sprintf("[II*-IV*-%ld] page 164",(d-9)/6);
1291 7 : break;
1292 0 : default: pari_err_BUG("tame5 [bug31]");
1293 : }
1294 14 : break;
1295 21 : case 4:
1296 21 : switch(r)
1297 : {
1298 7 : case 1: case 3: condp = 4;
1299 7 : Ip->type = stack_sprintf("[I*{0}-II-%ld] page 160",(d-4)/6);
1300 7 : Ip->neron = dicyclic(2,2); break;
1301 14 : case 5: condp = 4;
1302 14 : Ip->type = stack_sprintf("[II*-II*-%ld] page 163",(d-10)/6);
1303 14 : Ip->neron = cyclic(1); break;
1304 0 : default: pari_err_BUG("tame5 [bug32]");
1305 : }
1306 21 : break;
1307 21 : case 5:
1308 21 : switch(r)
1309 : {
1310 7 : case 0: case 5: condp = 2;
1311 7 : Ip->type = stack_sprintf("[I{0}-II*-%ld] page 160",(d-5)/6);
1312 7 : Ip->neron = cyclic(1); break;
1313 7 : case 1: case 4: condp = 4;
1314 7 : Ip->type = stack_sprintf("[II-IV*-%ld] page 164",(d-5)/6);
1315 7 : Ip->neron = cyclic(3); break;
1316 7 : case 2: case 3: condp = 4;
1317 7 : Ip->type = stack_sprintf("[I*{0}-IV-%ld] page 161",(d-5)/6);
1318 7 : Ip->neron = mkvecsmall2(6,2); break;
1319 : }
1320 21 : break;
1321 0 : default: pari_err_BUG("tame5 [bug33]");
1322 : }
1323 105 : break;
1324 56 : case 12:
1325 56 : condp = 4;
1326 56 : switch(dm)
1327 : {
1328 14 : case 1:
1329 14 : switch(r)
1330 : {
1331 7 : case 3: case 10:
1332 7 : Ip->type = stack_sprintf("[II*-III-%ld] page 166",(d-13)/12);
1333 7 : Ip->neron = cyclic(2); break;
1334 7 : case 4: case 9:
1335 7 : Ip->type = stack_sprintf("[III*-IV-%ld] page 167",(d-13)/12);
1336 7 : Ip->neron = cyclic(6); break;
1337 0 : default: pari_err_BUG("tame5 [bug34]");
1338 : }
1339 14 : break;
1340 14 : case 5:
1341 14 : switch(r)
1342 : {
1343 7 : case 2: case 3:
1344 7 : Ip->type = stack_sprintf("[II-III-%ld] page 166",(d-5)/12);
1345 7 : Ip->neron = cyclic(2); break;
1346 7 : case 8: case 9:
1347 7 : Ip->type = stack_sprintf("[III*-IV*-%ld] page 168",(d-17)/12);
1348 7 : Ip->neron = cyclic(6); break;
1349 0 : default: pari_err_BUG("tame5 [bug35]");
1350 : }
1351 14 : break;
1352 14 : case 7:
1353 14 : switch(r)
1354 : {
1355 7 : case 3: case 4:
1356 7 : Ip->type = stack_sprintf("[III-IV-%ld] page 167",(d-7)/12);
1357 7 : Ip->neron = cyclic(6); break;
1358 7 : case 9: case 10:
1359 7 : Ip->type = stack_sprintf("[II*-III*-%ld] page 167",(d-19)/12);
1360 7 : Ip->neron = cyclic(2); break;
1361 0 : default: pari_err_BUG("tame5 [bug36]");
1362 : }
1363 14 : break;
1364 14 : case 11:
1365 14 : switch(r)
1366 : {
1367 7 : case 3: case 8:
1368 7 : Ip->type = stack_sprintf("[III-IV*-%ld] page 168",(d-11)/12);
1369 7 : Ip->neron = cyclic(6); break;
1370 7 : case 2: case 9:
1371 7 : Ip->type = stack_sprintf("[II-III*-%ld] page 166",(d-11)/12);
1372 7 : Ip->neron = cyclic(2); break;
1373 0 : default: pari_err_BUG("tame5 [bug37]");
1374 : }
1375 14 : break;
1376 0 : default: pari_err_BUG("tame5 [bug38]");
1377 : }
1378 56 : break;
1379 0 : default: pari_err_BUG("tame5 [bug39]");
1380 : }
1381 : }
1382 : else
1383 : {
1384 63 : r %= (n >> 1);
1385 63 : switch(n)
1386 : {
1387 7 : case 2: condp = 2;
1388 7 : Ip->type = stack_sprintf("[2I{0}-%ld] page 159",(d/2));
1389 7 : Ip->neron = cyclic(1); break;
1390 14 : case 4: condp = 4;
1391 14 : Ip->type = stack_sprintf("[2I*{0}-%ld] page 159",(d/2-1)/2);
1392 14 : Ip->neron = dicyclic(2,2); break;
1393 14 : case 6: condp = 4;
1394 14 : Ip->neron = cyclic(3);
1395 14 : switch(r)
1396 : {
1397 7 : case 1:
1398 7 : Ip->type = stack_sprintf("[2IV-%ld] page 165",(d/2-1)/3);
1399 7 : break;
1400 7 : case 2:
1401 7 : Ip->type = stack_sprintf("[2IV*-%ld] page 165",(d/2-2)/3);
1402 7 : break;
1403 0 : default: pari_err_BUG("tame5 [bug40]");
1404 : }
1405 14 : break;
1406 14 : case 8: condp = 4;
1407 14 : Ip->neron = cyclic(2);
1408 14 : switch(r)
1409 : {
1410 7 : case 1:
1411 7 : Ip->type = stack_sprintf("[2III-%ld] page 168",(d/2-1)/4);
1412 7 : break;
1413 7 : case 3:
1414 7 : Ip->type = stack_sprintf("[2III*-%ld] page 168",(d/2-3)/4);
1415 7 : break;
1416 0 : default: pari_err_BUG("tame5 [bug41]");
1417 : }
1418 14 : break;
1419 14 : case 12: condp = 4;
1420 14 : Ip->neron = cyclic(1);
1421 14 : switch(r)
1422 : {
1423 7 : case 1:
1424 7 : Ip->type = stack_sprintf("[2II-%ld] page 162",(d/2-1)/6);
1425 7 : break;
1426 7 : case 5:
1427 7 : Ip->type = stack_sprintf("[2II*-%ld] page 163",(d/2-5)/6);
1428 7 : break;
1429 0 : default: pari_err_BUG("tame5 [bug42]");
1430 : }
1431 14 : break;
1432 0 : default: pari_err_BUG("tame5 [bug43]");
1433 : }
1434 : }
1435 329 : return condp;
1436 : }
1437 :
1438 : static long
1439 189 : tame_6(struct igusa *I, struct igusa_p *Ip)
1440 : {
1441 189 : long condp = -1, d, d1, n, dm, r, dk;
1442 189 : GEN val = Ip->val;
1443 :
1444 189 : dk = Ip->eps*val[7]-6*val[8];
1445 189 : tame_567_init(I, Ip, dk, &d, &n, &dm, &r);
1446 189 : d1 = n * (Ip->eps*(val[6]-val[7])+val[8]) / Ip->eps;
1447 189 : switch(n)
1448 : {
1449 56 : case 1: condp = 1;
1450 56 : Ip->type = stack_sprintf("[I{0}-I{%ld}-%ld] page 170",d1,d);
1451 56 : Ip->neron = cyclic(d1); break;
1452 28 : case 2:
1453 28 : switch(dm)
1454 : {
1455 7 : case 0: condp = 4;
1456 7 : Ip->type=stack_sprintf("[I*{0}-I*{%ld}-%ld] page 171", d1/2,(d-2)/2);
1457 7 : Ip->neron = shallowconcat(groupH(d1/2), dicyclic(2,2)); break;
1458 21 : case 1: return -1;
1459 0 : default: pari_err_BUG("tame6 [bug44]");
1460 : }
1461 7 : break;
1462 14 : case 3: condp = 3;
1463 14 : Ip->neron = dicyclic(3,d1/3);
1464 14 : switch(dm)
1465 : {
1466 7 : case 1:
1467 7 : Ip->type = stack_sprintf("[I{%ld}-IV-%ld] page 173",d1/3,(d-1)/3);
1468 7 : break;
1469 7 : case 2:
1470 7 : Ip->type = stack_sprintf("[I{%ld}-IV*-%ld] page 173",d1/3,(d-2)/3);
1471 7 : break;
1472 0 : default: pari_err_BUG("tame6 [bug45]");
1473 : }
1474 14 : break;
1475 35 : case 4:
1476 35 : switch(dm)
1477 : {
1478 21 : case 1:
1479 21 : switch(r)
1480 : {
1481 7 : case 0: case 1: condp = 3;
1482 7 : Ip->type=stack_sprintf("[I{%ld}-III-%ld] page 176",d1/4,(d-1)/4);
1483 7 : Ip->neron = dicyclic(2,d1/4); break;
1484 14 : case 2: case 3: condp = 4;
1485 14 : Ip->type=stack_sprintf("[I*{%ld}-III*-%ld] page 177",d1/4,(d-5)/4);
1486 14 : Ip->neron = shallowconcat(groupH(d1/4), cyclic(2)); break;
1487 0 : default: pari_err_BUG("tame6 [bug46]");
1488 : }
1489 21 : break;
1490 14 : case 3:
1491 14 : switch(r)
1492 : {
1493 7 : case 0: case 3: condp = 3;
1494 7 : Ip->type=stack_sprintf("[I{%ld}-III*-%ld] page 176",d1/4,(d-3)/4);
1495 7 : Ip->neron = dicyclic(2,d1/4); break;
1496 7 : case 1: case 2: condp = 4;
1497 7 : Ip->type=stack_sprintf("[I*{%ld}-III-%ld] page 177",d1/4,(d-3)/4);
1498 7 : Ip->neron = shallowconcat(groupH(d1/4), cyclic(2)); break;
1499 0 : default: pari_err_BUG("tame6 [bug47]");
1500 : }
1501 14 : break;
1502 0 : default: pari_err_BUG("tame6 [bug48]");
1503 : }
1504 35 : break;
1505 56 : case 6:
1506 56 : switch(dm)
1507 : {
1508 21 : case 1:
1509 21 : switch(r)
1510 : {
1511 7 : case 0: case 1: condp = 3;
1512 7 : Ip->type = stack_sprintf("[I{%ld}-II-%ld] page 172",d1/6,(d-1)/6);
1513 7 : Ip->neron = cyclic(d1/6); break;
1514 14 : case 3: case 4: condp = 4;
1515 14 : Ip->type=stack_sprintf("[I*{%ld}-IV*-%ld] page 174",d1/6,(d-7)/6);
1516 14 : Ip->neron = shallowconcat(groupH(d1/6), cyclic(3)); break;
1517 0 : default: pari_err_BUG("tame6 [bug49]");
1518 : }
1519 21 : break;
1520 14 : case 2: condp = 4;
1521 14 : Ip->type = stack_sprintf("[I*{%ld}-II*-%ld] page 174",d1/6,(d-8)/6);
1522 14 : Ip->neron = groupH(d1/6); break;
1523 7 : case 4: condp = 4;
1524 7 : Ip->type = stack_sprintf("[I*{%ld}-II-%ld] page 173",d1/6,(d-4)/6);
1525 7 : Ip->neron = groupH(d1/6); break;
1526 14 : case 5:
1527 14 : switch(r)
1528 : {
1529 7 : case 0: case 5: condp = 3;
1530 7 : Ip->type=stack_sprintf("[I{%ld}-II*-%ld] page 172",d1/6,(d-5)/6);
1531 7 : Ip->neron = cyclic(d1/6); break;
1532 7 : case 2: case 3: condp = 4;
1533 7 : Ip->type=stack_sprintf("[I*{%ld}-IV-%ld] page 174",d1/6,(d-5)/6);
1534 7 : Ip->neron = shallowconcat(groupH(d1/6), cyclic(3)); break;
1535 0 : default: pari_err_BUG("tame6 [bug50]");
1536 : }
1537 14 : break;
1538 0 : default: pari_err_BUG("tame6 [bug51]");
1539 : }
1540 56 : break;
1541 0 : default: pari_err_BUG("tame6 [bug52]");
1542 : }
1543 168 : return condp;
1544 : }
1545 :
1546 : static long
1547 91 : tame_7(struct igusa *I, struct igusa_p *Ip)
1548 : {
1549 91 : long condp = -1, d, D, d1, d2, n, dm, r, dk;
1550 91 : GEN val = Ip->val;
1551 :
1552 91 : dk = 3*(Ip->eps*val[3]-2*val[8]);
1553 91 : tame_567_init(I, Ip, dk, &d, &n, &dm, &r);
1554 91 : D = n * (Ip->eps*(val[6]-3*val[3])+val[8]) / Ip->eps;
1555 91 : d1 = minss(n * (val[7]-3*val[3]), D/2);
1556 91 : d2 = D - d1;
1557 : /* d1 <= d2 */
1558 91 : switch(n)
1559 : {
1560 42 : case 1: condp = 2;
1561 42 : Ip->type = stack_sprintf("[I{%ld}-I{%ld}-%ld] page 179",d1,d2,d);
1562 42 : Ip->neron = dicyclic(d1,d2); break;
1563 35 : case 2:
1564 35 : if (odd(val[8]))
1565 : {
1566 14 : condp = 3;
1567 14 : Ip->type = stack_sprintf("[2I{%ld}-%ld] page 181",d1,d/2);
1568 14 : Ip->neron = cyclic(d1);
1569 : }
1570 21 : else if (dm == 0)
1571 : {
1572 14 : condp = 4;
1573 14 : Ip->type = stack_sprintf("[I*{%ld}-I*{%ld}-%ld] page 180", d1/2,d2/2,(d-2)/2);
1574 14 : Ip->neron = shallowconcat(groupH(d1/2),groupH(d2/2));
1575 : }
1576 : else
1577 : {
1578 : GEN H;
1579 7 : if (d1 != d2) return -1;
1580 0 : condp = 3; H = groupH(d1/2);
1581 0 : Ip->type = stack_sprintf("[I{%ld}-I*{%ld}-%ld] page 180", d1/2,d1/2,(d-1)/2);
1582 0 : Ip->neron = shallowconcat(H, H);
1583 : }
1584 28 : break;
1585 14 : case 4: condp = 4;
1586 14 : Ip->type = stack_sprintf("[2I*{%ld}-%ld] page 181",d1/2,(d-2)/4);
1587 14 : Ip->neron = groupH(d1/2); break;
1588 0 : default: pari_err_BUG("tame7 [bug55]");
1589 : }
1590 84 : return condp;
1591 : }
1592 :
1593 : static long labelm3(GEN polh, long t60, long alpha, long Dmin, struct igusa *I, struct igusa_p *Ip);
1594 : static long
1595 840 : tame(GEN polh, long t60, long alpha, long Dmin, struct igusa *I, struct igusa_p *Ip)
1596 : {
1597 : long d;
1598 840 : Ip->tame = 1;
1599 840 : switch(Ip->tt)
1600 : {
1601 28 : case 1: return tame_1(I,Ip);
1602 91 : case 2: return tame_2(I,Ip);
1603 56 : case 3: return tame_3(I,Ip);
1604 56 : case 4: return tame_4(I,Ip);
1605 329 : case 5: return tame_5(I,Ip);
1606 189 : case 6: d = tame_6(I,Ip); break;
1607 91 : default:d = tame_7(I,Ip); break;
1608 : }
1609 280 : if (d < 0) d = labelm3(polh,t60,alpha,Dmin,I,Ip); /* => tt=6 or 7 */
1610 280 : return d;
1611 : }
1612 :
1613 : /* maxc = maximum conductor valuation at p */
1614 : static long
1615 504 : get_maxc(GEN p)
1616 : {
1617 504 : switch (itos_or_0(p))
1618 : {
1619 0 : case 2: return 20; break;
1620 301 : case 3: return 10; break;
1621 14 : case 5: return 9; break;
1622 189 : default: return 4; break; /* p > 5 */
1623 : }
1624 : }
1625 :
1626 : /* p = 3 */
1627 : static long
1628 84 : quartic(GEN polh, long alpha, long Dmin, struct igusa_p *Ip)
1629 : {
1630 84 : GEN val = Ip->val, p = Ip->p;
1631 84 : GEN polf = polymini_zi2(ZX_Z_mul(polh, powiu(p, alpha)));
1632 84 : long condp = -1, d, R, r1, beta;
1633 84 : r1 = polf[1];
1634 84 : beta = polf[2];
1635 84 : R = beta/2;
1636 84 : switch(Ip->tt)
1637 : {
1638 70 : case 1: case 5: d = 0;break;
1639 0 : case 3: d = val[6] - 5*val[3]/2;break;
1640 14 : case 7: d = val[6] - 3*val[3] + val[8]/Ip->eps;break;
1641 0 : default: pari_err_BUG("quartic [type choices]");
1642 : d = 0; /*LCOV_EXCL_LINE*/
1643 : }
1644 84 : switch(r1)
1645 : {
1646 21 : case 0:
1647 21 : if (d)
1648 : {
1649 7 : condp = 3;
1650 7 : Ip->type = stack_sprintf("[2I{%ld}-%ld] page 181",d,R);
1651 7 : Ip->neron = cyclic(d);
1652 : }
1653 : else
1654 : {
1655 14 : condp = 2;
1656 14 : Ip->neron = cyclic(1);
1657 14 : if (R) Ip->type = stack_sprintf("[2I{0}-%ld] page 159",R);
1658 7 : else Ip->type = "[II] page 155";
1659 : }
1660 21 : break;
1661 14 : case 6: condp = 4;
1662 14 : Ip->type = stack_sprintf("[2I*{%ld}-%ld] pages 159, 181",d,R);
1663 14 : Ip->neron = dicyclic(2,2); break;
1664 7 : case 3: condp = 4;
1665 7 : Ip->type = stack_sprintf("[2III-%ld] page 168",R);
1666 7 : Ip->neron = cyclic(2); break;
1667 7 : case 9: condp = 4;
1668 7 : Ip->type = stack_sprintf("[2III*-%ld] page 168",R);
1669 7 : Ip->neron = cyclic(2); break;
1670 7 : case 2: condp = Dmin-12*R-13;
1671 7 : Ip->type = stack_sprintf("[2II-%ld] page 162",R);
1672 7 : Ip->neron = cyclic(1); break;
1673 14 : case 8: condp = Dmin-12*R-19;
1674 14 : Ip->type = stack_sprintf("[2IV*-%ld] page 165",R);
1675 14 : Ip->neron = cyclic(3); break;
1676 7 : case 4: condp = Dmin-12*R-15;
1677 7 : Ip->type = stack_sprintf("[2IV-%ld] page 165",R);
1678 7 : Ip->neron = cyclic(3); break;
1679 7 : case 10: condp = Dmin-12*R-21;
1680 7 : Ip->type = stack_sprintf("[2II*-%ld] page 163",R);
1681 7 : Ip->neron = cyclic(1); break;
1682 0 : default: pari_err_BUG("quartic [type1]");
1683 : }
1684 84 : if (condp > get_maxc(p) || condp < 0) pari_err_BUG("quartic [conductor]");
1685 84 : return condp;
1686 : }
1687 :
1688 : static long
1689 280 : litredtp(long alpha, long alpha1, long t60, long t60_1, GEN polh, GEN polh1,
1690 : long Dmin, long R, struct igusa *I, struct igusa_p *Ip)
1691 : {
1692 280 : GEN val = Ip->val, p = Ip->p;
1693 280 : long condp = -1, indice, d;
1694 :
1695 280 : if ((Ip->r1 == 0||Ip->r1 == 6) && (Ip->r2 == 0||Ip->r2 == 6))
1696 : { /* (r1,r2) = (0,0), (0,6), (6,0) or (6,6) */
1697 154 : if (Ip->tt == 5)
1698 : {
1699 21 : switch(Ip->r1 + Ip->r2)
1700 : {
1701 7 : case 0: /* (0,0) */
1702 7 : condp = 0;
1703 7 : Ip->type = stack_sprintf("[I{0}-I{0}-%ld] page 158",R);
1704 7 : Ip->neron = cyclic(1); break;
1705 7 : case 6: /* (0,6) or (6,0) */
1706 7 : condp = 2;
1707 7 : Ip->type = stack_sprintf("[I{0}-I*{0}-%ld] page 159",R);
1708 7 : Ip->neron = dicyclic(2,2); break;
1709 7 : case 12: /* (6,6) */
1710 7 : condp = 4;
1711 7 : Ip->type = stack_sprintf("[I*{0}-I*{0}-%ld] page 158",R);
1712 7 : Ip->neron = mkvecsmall4(2,2,2,2); break;
1713 : }
1714 21 : return condp;
1715 : }
1716 133 : if (Ip->r1 == Ip->r2) return tame(polh, t60, alpha, Dmin, I, Ip);
1717 42 : if (Ip->tt == 6)
1718 : {
1719 28 : d = val[6] - val[7] + val[8]/Ip->eps;
1720 28 : if (Ip->r1 && alpha1 == 0) polh1 = ZX_unscale_divpow(polh1, p, 3);
1721 28 : if (FpX_is_squarefree(FpX_red(polh1,p),p))
1722 7 : { indice = 0; condp = 3-Ip->r2/6; }
1723 : else
1724 21 : { indice = d; condp = 3-Ip->r1/6; }
1725 : }
1726 : else
1727 : { /* Ip->tt == 7 */
1728 : long d1;
1729 14 : d = val[6] - 3*val[3] + val[8]/Ip->eps;
1730 14 : if (t60_1 == 60) polh1 = ZX_unscale_divpow(polh1, p, 3);
1731 14 : d1 = minss(val[7]-3*val[3],d/2);
1732 14 : if (d == 2*d1) indice = d1;
1733 : else
1734 : {
1735 14 : indice = discpart(polh1,p,d1+1);
1736 14 : if (indice>= d1+1) indice = d-d1; else indice = d1;
1737 : }
1738 14 : condp = 3;
1739 : }
1740 42 : if (Ip->r1) indice = d - indice; /* (r1,r2) = (6,0) */
1741 42 : Ip->neron = shallowconcat(cyclic(indice),groupH(d-indice));
1742 42 : Ip->type = stack_sprintf("[I{%ld}-I*{%ld}-%ld] page %ld",
1743 42 : indice,d-indice,R, (Ip->tt==6)? 170L: 180L);
1744 42 : return condp;
1745 : }
1746 126 : if (Ip->tt == 7) pari_err_BUG("litredtp [switch ri]");
1747 : {
1748 126 : struct red __S1, __S2, *S1 = &__S1, *S2 = &__S2;
1749 126 : long f1 = get_red(S1, Ip, polh1, p, alpha1, Ip->r1);
1750 126 : long f2 = get_red(S2, Ip, polh, p, alpha, Ip->r2);
1751 : /* reorder to normalize representation */
1752 126 : if (S1->tnum > S2->tnum || (S1->tnum == S2->tnum && f1 > f2))
1753 56 : { struct red *S = S1; S1 = S2; S2 = S; }
1754 126 : Ip->type = stack_sprintf("[%s-%s-%ld] pages %s", S1->t,S2->t, R, S1->pages);
1755 126 : Ip->neron = shallowconcat(S1->g, S2->g);
1756 126 : condp = Dmin - (f1 + f2) + ((R >= 0)? 2-12*R: 4);
1757 : }
1758 126 : if (condp > get_maxc(p)) pari_err_BUG("litredtp [conductor]");
1759 126 : return condp;
1760 : }
1761 :
1762 : static long
1763 259 : labelm3(GEN h1, long t60_1, long alpha1, long Dmin, struct igusa *I, struct igusa_p *Ip)
1764 : {
1765 259 : GEN h, pm, vs, val = Ip->val, p = Ip->p;
1766 : long alpha, t60, lambda, beta, R;
1767 :
1768 259 : pm = polymini(ZX_Z_mul(RgX_recip6(h1), powiu(p,alpha1)), p);
1769 259 : h = gel(pm,1); vs = gel(pm,2);
1770 259 : lambda= vs[1];
1771 259 : t60 = vs[2];
1772 259 : alpha = vs[3];
1773 259 : beta = vs[5];
1774 259 : if (lambda != 3) pari_err_BUG("labelm3 [lambda != 3]");
1775 259 : R = beta-(alpha1+alpha);
1776 259 : if (odd(R)) pari_err_BUG("labelm3 [R odd]");
1777 259 : R /= 2;
1778 259 : if (R <= -2) pari_err_BUG("labelm3 [R <= -2]");
1779 259 : if (val[8] % (2*Ip->eps)) pari_err_BUG("labelm3 [val(eps2)]");
1780 259 : if (R >= 0 && (alpha+alpha1) >= 1) pari_err_BUG("labelm3 [minimal equation]");
1781 259 : Ip->r1 = t60_1 / 10 + 6*alpha1;
1782 259 : Ip->r2 = t60 / 10 + 6*alpha;
1783 259 : return litredtp(alpha, alpha1, t60, t60_1, h, h1, Dmin, R, I, Ip);
1784 : }
1785 :
1786 : /* p = 3 */
1787 : static long
1788 21 : quadratic(GEN polh, long alpha, long Dmin, struct igusa *I, struct igusa_p *Ip)
1789 : {
1790 21 : long alpha1 = alpha, beta, t6, R;
1791 21 : GEN vs = polymini_zi(ZX_Z_mul(polh, powiu(Ip->p,alpha)));
1792 21 : t6 = vs[1];
1793 21 : alpha = vs[2];
1794 21 : beta = vs[3];
1795 21 : R = beta-alpha;
1796 21 : if (R >= 0 && alpha1)
1797 : {
1798 0 : Dmin -= 10;
1799 0 : if (DEBUGLEVEL)
1800 0 : err_printf("(Care: minimal discriminant over Z[i] smaller than over Z)\n");
1801 : }
1802 21 : Ip->r2 = Ip->r1 = t6 + 6*alpha;
1803 21 : return litredtp(alpha, alpha, t6*10, t6*10, polh, polh, Dmin, R, I, Ip);
1804 : }
1805 :
1806 : static long
1807 1442 : genus2localred(struct igusa *I, struct igusa_p *Ip, GEN p, GEN polmini)
1808 : {
1809 : GEN val, vs, polh, list, c1, c2, c3, c4, c5, c6, prod;
1810 : long i, vb5, vb6, d, Dmin, alpha, lambda, t60;
1811 1442 : long condp = -1, indice, vc6, mm, nb, dism;
1812 :
1813 1442 : stable_reduction(I, Ip, p);
1814 1442 : val = Ip->val; Dmin = val[6];
1815 1442 : if (Dmin == 0)
1816 : {
1817 7 : Ip->tame = 1;
1818 7 : Ip->type = "[I{0-0-0}] page 155";
1819 7 : Ip->neron = cyclic(1); return 0;
1820 : }
1821 1435 : if (Dmin == 1)
1822 : {
1823 14 : Ip->type = "[I{1-0-0}] page 170";
1824 14 : Ip->neron = cyclic(1); return 1;
1825 : }
1826 1421 : if (Dmin == 2) switch(Ip->tt)
1827 : {
1828 0 : case 2:
1829 0 : Ip->type = "[I{2-0-0}] page 170";
1830 0 : Ip->neron = cyclic(2); return 1;
1831 0 : case 3:
1832 0 : Ip->type = "[I{1-1-0}] page 179";
1833 0 : Ip->neron = cyclic(1); return 2;
1834 14 : case 5:
1835 14 : if (cmpis(p,3) <= 0) pari_err_BUG("genus2localred [tt 1]");
1836 14 : Ip->type = "[I{0}-II-0] page 159";
1837 14 : Ip->neron = cyclic(1); return 2;
1838 0 : default: pari_err_BUG("genus2localred [tt 2]");
1839 : }
1840 1407 : if (absequaliu(p,2)) return -1;
1841 1386 : polh = gel(polmini,1); vs = gel(polmini,2);
1842 1386 : lambda = vs[1];
1843 1386 : t60 = vs[2];
1844 1386 : alpha = vs[3];
1845 1386 : if (vs[4]) return equaliu(p,3)? quadratic(polh, alpha, Dmin, I, Ip):
1846 0 : tame(polh, t60, alpha, Dmin, I, Ip);
1847 1365 : if (!t60 && lambda<= 2)
1848 : {
1849 7 : if (Ip->tt >= 5) pari_err_BUG("genus2localred [tt 3]");
1850 7 : return tame(polh, t60, alpha, Dmin, I, Ip);
1851 : }
1852 1358 : if (Dmin == 3)
1853 : {
1854 7 : switch(Ip->tt)
1855 : {
1856 0 : case 2: return tame(polh, t60, alpha, Dmin, I, Ip);
1857 0 : case 3: Ip->type = "[I{2-1-0}] page 179"; Ip->neron = cyclic(2); return 2;
1858 7 : case 4: Ip->type = "[I{1-1-1}] page 182"; Ip->neron = cyclic(3); return 2;
1859 0 : case 5:
1860 0 : if (equaliu(p,3) && t60 != 30)
1861 0 : return labelm3(polh,t60,alpha,Dmin,I,Ip);
1862 0 : Ip->type = "[I{0}-III-0] page 161"; Ip->neron = cyclic(2); return 2;
1863 0 : case 6:
1864 0 : if (equaliu(p,3)) pari_err_BUG("genus2localred [conductor]");
1865 0 : Ip->type = "[I{1}-II-0] page 172"; Ip->neron = cyclic(1); return 3;
1866 : }
1867 0 : pari_err_BUG("genus2localred [switch tt 4]");
1868 : return -1; /* LCOV_EXCL_LINE */
1869 : }
1870 1351 : switch(lambda)
1871 : {
1872 364 : case 0:
1873 364 : switch(t60+alpha)
1874 : {
1875 7 : case 10:
1876 7 : condp = Dmin-1;
1877 7 : Ip->type = "[V] page 156";
1878 7 : Ip->neron = cyclic(3); break;
1879 7 : case 11:
1880 7 : condp = Dmin-11;
1881 7 : Ip->type = "[V*] page 156";
1882 7 : Ip->neron = cyclic(3); break;
1883 7 : case 12:
1884 7 : condp = Dmin-2;
1885 7 : Ip->type = "[IX-2] page 157";
1886 7 : Ip->neron = cyclic(5); break;
1887 14 : case 13:
1888 14 : condp = Dmin-12;
1889 14 : Ip->type = "[VIII-4] page 157";
1890 14 : Ip->neron = cyclic(1); break;
1891 7 : case 24:
1892 7 : condp = Dmin-8;
1893 7 : Ip->type = "[IX-4] page 158";
1894 7 : Ip->neron = cyclic(5);
1895 7 : break;
1896 14 : case 15: case 16:
1897 14 : if (Ip->tt>= 5) pari_err_BUG("genus2localred [tt 6]");
1898 14 : return tame(polh, t60, alpha, Dmin, I, Ip);
1899 112 : case 20: case 21:
1900 : {
1901 : GEN b0, b1, b2, b3, b4, b5, b6, b02, b03, b04, b05;
1902 112 : RgX_to_06(polh, &b0,&b1,&b2,&b3,&b4,&b5,&b6);
1903 112 : vb5 = myval(b5,p);
1904 112 : vb6 = myval(b6,p);
1905 112 : if (vb6 >= 3)
1906 : {
1907 14 : if (vb5 < 2) pari_err_BUG("genus2localred [red1]");
1908 14 : if (vb5 >= 3)
1909 : {
1910 7 : condp = Dmin-8;
1911 7 : Ip->type = "[II*-IV-(-1)] page 164";
1912 7 : Ip->neron = cyclic(3);
1913 : }
1914 : else
1915 : {
1916 7 : condp = Dmin-7;
1917 7 : Ip->type = "[IV-III*-(-1)] page 167";
1918 7 : Ip->neron = cyclic(6);
1919 : }
1920 14 : break;
1921 : }
1922 98 : if (dvdii(b0,p)) pari_err_BUG("genus2localred [b0]");
1923 98 : b02 = gsqr(b0);
1924 98 : b03 = gmul(b02, b0);
1925 98 : b04 = gmul(b03, b0);
1926 98 : b05 = gmul(b04, b0);
1927 98 : c1 = gmul2n(b1,-1);
1928 98 : c2 = gmul2n(gsub(gmul(b0,b2), gsqr(c1)),-1);
1929 98 : c3 = gmul2n(gsub(gmul(b02,b3), gmul2n(gmul(c1,c2),1)),-1);
1930 98 : c4 = gsub(gmul(b03,b4), gadd(gmul2n(gmul(c1,c3),1),gsqr(c2)));
1931 98 : c5 = gsub(gmul(b04,b5), gmul2n(gmul(c2,c3),1));
1932 98 : c6 = gsub(gmul(b05,b6), gsqr(c3));
1933 : /* b0^5*H(x/b0) = (x^3+c1*x^2+c2*x+c3)^2+c4*x^2+c5*x+c6 */
1934 98 : vc6 = myval(c6,p);
1935 98 : if (vc6 == 2)
1936 : {
1937 7 : if (alpha)
1938 : {
1939 0 : condp = Dmin-16;
1940 0 : Ip->type = "[IV] page 155";
1941 0 : Ip->neron = cyclic(1);
1942 : }
1943 : else
1944 : {
1945 7 : condp = Dmin-6;
1946 7 : Ip->type = "[III] page 155";
1947 7 : Ip->neron = dicyclic(3,3);
1948 : }
1949 : }
1950 : else
1951 : {
1952 91 : if (myval(c3,p) > 1) pari_err_BUG("genus2localred [c3]");
1953 91 : mm = min3(3*myval(c4,p)-4, 3*myval(c5,p)-5, 3*vc6-6);
1954 91 : if (alpha)
1955 : {
1956 35 : condp = Dmin-mm-16;
1957 35 : Ip->type = stack_sprintf("[III*{%ld}] page 184", mm);
1958 35 : Ip->neron = cyclic(1);
1959 : }
1960 : else
1961 : {
1962 56 : condp = Dmin-mm-6;
1963 56 : Ip->type = stack_sprintf("[III{%ld}] page 184", mm);
1964 56 : Ip->neron = (mm%3)? cyclic(9): dicyclic(3,3);
1965 : }
1966 : }
1967 : }
1968 98 : break;
1969 196 : case 30:
1970 280 : return equaliu(p,3)? quartic(polh, alpha, Dmin, Ip)
1971 280 : : tame(polh, t60, alpha, Dmin, I, Ip);
1972 0 : default: pari_err_BUG("genus2localred [red2]");
1973 : }
1974 154 : break;
1975 119 : case 1:
1976 119 : switch(t60+alpha)
1977 : {
1978 7 : case 12:
1979 7 : condp = Dmin;
1980 7 : Ip->type = "[VIII-1] page 156";
1981 7 : Ip->neron = cyclic(1); break;
1982 7 : case 13:
1983 7 : condp = Dmin-10;
1984 7 : Ip->type = "[IX-3] page 157";
1985 7 : Ip->neron = cyclic(5); break;
1986 7 : case 24:
1987 7 : condp = Dmin-4;
1988 7 : Ip->type = "[IX-1] page 157";
1989 7 : Ip->neron = cyclic(5); break;
1990 7 : case 25:
1991 7 : condp = Dmin-14;
1992 7 : Ip->type = "[VIII-3] page 157";
1993 7 : Ip->neron = cyclic(1); break;
1994 7 : case 36:
1995 7 : condp = Dmin-8;
1996 7 : Ip->type = "[VIII-2] page 157";
1997 7 : Ip->neron = cyclic(1); break;
1998 21 : case 15:
1999 21 : condp = Dmin-1;
2000 21 : Ip->type = "[VII] page 156";
2001 21 : Ip->neron = cyclic(2); break;
2002 7 : case 16:
2003 7 : condp = Dmin-11;
2004 7 : Ip->type = "[VII*] page 156";
2005 7 : Ip->neron = cyclic(2); break;
2006 14 : case 20:
2007 14 : if (cmpis(p,3))
2008 : {
2009 7 : d = 6*val[6]-5*val[7]-2;
2010 7 : if (d%6) pari_err_BUG("genus2localred [index]");
2011 7 : dism = (d/6);
2012 : }
2013 : else
2014 : {
2015 7 : list = padicfactors(polh,p,Dmin-5);
2016 7 : nb = lg(list);
2017 7 : prod = pol_1(varn(polh));
2018 21 : for(i = 1;i<nb;i++)
2019 : {
2020 14 : GEN c = gel(list,i);
2021 14 : if (valp(gel(c,2)) && degpol(c)<= 2) prod = RgX_mul(prod,c);
2022 : }
2023 7 : if (degpol(prod) > 2) pari_err_BUG("genus2localred [padicfactors]");
2024 7 : dism = valp(RgX_disc(prod)) - 1;
2025 : }
2026 14 : condp = Dmin-dism-3;
2027 14 : Ip->type = stack_sprintf("[II-II*{%ld}] page 176", dism);
2028 14 : Ip->neron = groupH(dism+1); break;
2029 14 : case 21:
2030 14 : vb6 = myval(RgX_coeff(polh,0),p);
2031 14 : if (vb6<2) pari_err_BUG("genus2localred [red3]");
2032 14 : condp = Dmin-14;
2033 14 : Ip->type = "[IV*-II{0}] page 175";
2034 14 : Ip->neron = cyclic(1); break;
2035 28 : case 30:
2036 28 : vb5 = myval(RgX_coeff(polh,1),p);
2037 28 : if (vb5 == 2)
2038 : {
2039 21 : if (Ip->tt >= 5) pari_err_BUG("genus2localred [tt 6]");
2040 21 : return tame(polh, t60, alpha, Dmin, I, Ip);
2041 : }
2042 7 : condp = Dmin-7;
2043 7 : Ip->type = "[II*-III-(-1)] page 167";
2044 7 : Ip->neron = cyclic(2); break;
2045 : }
2046 98 : break;
2047 147 : case 2:
2048 147 : if (ugcd(t60, 60) == 15) /* denom(theta) = 4 */
2049 : {
2050 28 : if (Ip->tt>4) pari_err_BUG("genus2localred [tt 5]");
2051 28 : return tame(polh, t60, alpha, Dmin, I, Ip);
2052 : }
2053 119 : if (!equaliu(p,3) && ugcd(t60, 60) == 20) /* denom(theta) = 3 */
2054 21 : return tame(polh, t60, alpha, Dmin, I, Ip);
2055 98 : list = padicfactors(polh,p,Dmin-10*alpha);
2056 98 : nb = lg(list); prod = pol_1(varn(polh));
2057 336 : for(i = 1;i<nb;i++)
2058 : {
2059 238 : GEN c = gel(list,i);
2060 238 : if (!valp(gel(c,2))) prod = RgX_mul(prod,c);
2061 : }
2062 98 : switch(degpol(prod))
2063 : {
2064 : GEN e0, e1, e2;
2065 0 : case 0:
2066 0 : dism = 0; break;
2067 7 : case 1:
2068 7 : e1 = gel(prod,3);
2069 7 : dism = 2*valp(e1); break;
2070 91 : case 2:
2071 91 : e0 = gel(prod,2);
2072 91 : e1 = gel(prod,3);
2073 91 : e2 = gel(prod,4);
2074 91 : dism = valp(gsub(gsqr(e1),gmul2n(gmul(e0,e2),2))); break;
2075 0 : default:
2076 0 : pari_err_BUG("genus2localred [padicfactors 2]");
2077 0 : dism = 0;
2078 : }
2079 98 : switch(t60/5+alpha-4)
2080 : {
2081 14 : case 0:
2082 14 : condp = Dmin-dism-1;
2083 14 : Ip->type = stack_sprintf("[IV-II{%ld}] page 175", dism);
2084 14 : Ip->neron = cyclic(3*dism+2); break;
2085 7 : case 1:
2086 7 : condp = Dmin-dism-10;
2087 7 : Ip->type = stack_sprintf("[II*-II*{%ld}] page 176",dism);
2088 7 : Ip->neron = groupH(dism+1); break;
2089 70 : case 2: case 3:
2090 70 : if (myval(RgX_coeff(polh,0),p) == 2)
2091 : {
2092 56 : if (Ip->tt>4) pari_err_BUG("genus2localred [tt 5]");
2093 56 : return tame(polh, t60, alpha, Dmin, I, Ip);
2094 : }
2095 14 : dism++;
2096 14 : indice = val[6]-(5*val[3]/2)-dism;
2097 14 : condp = Dmin-dism-indice-2;
2098 14 : Ip->type = stack_sprintf("[II{%ld-%ld}] page 182", dism,indice);
2099 14 : Ip->neron = both_odd(dism,indice)? dicyclic(2,2*dism): cyclic(4*dism);
2100 14 : break;
2101 7 : case 4:
2102 7 : condp = Dmin-dism-5;
2103 7 : Ip->type = stack_sprintf("[IV*-II{%ld}] page 175",dism+1);
2104 7 : Ip->neron = cyclic(3*dism+4); break;
2105 : }
2106 42 : break;
2107 721 : case 3:
2108 721 : if (!equaliu(p,3) || Ip->tt <= 4)
2109 490 : return tame(polh, t60, alpha, Dmin, I, Ip);
2110 231 : return labelm3(polh,t60,alpha,Dmin,I,Ip); /* p = 3 */
2111 0 : default: pari_err_BUG("genus2localred [switch lambda]");
2112 : }
2113 294 : if (condp < 2 || condp > get_maxc(p))
2114 0 : pari_err_BUG("genus2localred [conductor]");
2115 294 : return condp;
2116 : }
2117 :
2118 : static GEN
2119 1400 : hyperellintegralmodel(GEN PQ)
2120 : {
2121 : GEN D;
2122 1400 : PQ = Q_remove_denom(PQ, &D);
2123 1400 : if (!D) return PQ;
2124 14 : if (typ(PQ)==t_POL) return gmul(PQ,D);
2125 0 : if (typ(PQ) == t_VEC && lg(PQ) == 3)
2126 0 : return mkvec2(gmul(gel(PQ,1),D), gel(PQ,2));
2127 0 : pari_err_TYPE("hyperellintegralmodel",PQ);
2128 : return NULL; /* LCOV_EXCL_LINE */
2129 : }
2130 :
2131 : /* P,Q are ZX, study Y^2 + Q(X) Y = P(X) */
2132 : GEN
2133 1400 : genus2red(GEN PQ, GEN p)
2134 : {
2135 1400 : pari_sp av = avma;
2136 : struct igusa I;
2137 : GEN P, Q;
2138 : GEN j22, j42, j2j6, a0,a1,a2,a3,a4,a5,a6, V,polr,facto,factp, vecmini, cond;
2139 : long i, l, dd;
2140 1400 : PQ = hyperellminimalmodel(hyperellintegralmodel(PQ), NULL, p ? mkvec(p): p);
2141 1400 : P = gel(PQ,1);
2142 1400 : Q = gel(PQ,2);
2143 1400 : if (p && typ(p) != t_INT) pari_err_TYPE("genus2red", p);
2144 :
2145 1400 : polr = ZX_add(ZX_sqr(Q), gmul2n(P,2)); /* ZX */
2146 1400 : switch(degpol(polr))
2147 : {
2148 1400 : case 5: case 6: break;
2149 0 : default: pari_err_DOMAIN("genus2red","genus","!=", gen_2,mkvec2(P,Q));
2150 : }
2151 :
2152 1400 : RgX_to_03(polr, &a0,&a1,&a2,&a3);
2153 1400 : I.j10 = !signe(a0)? mulii(sqri(a1), ZX_disc(polr)): ZX_disc(polr);
2154 1400 : if (!signe(I.j10))
2155 0 : pari_err_DOMAIN("genus2red","genus","<",gen_2,mkvec2(P,Q));
2156 1400 : I.j10 = gmul2n(I.j10, -12); /* t_INT */
2157 :
2158 1400 : if (p == NULL)
2159 : {
2160 49 : facto = absZ_factor(I.j10);
2161 49 : factp = gel(facto,1);
2162 : }
2163 : else
2164 : {
2165 1351 : factp = mkcol(p);
2166 1351 : facto = mkmat2(factp, mkcol(gen_1));
2167 : }
2168 1400 : l = lg(factp);
2169 1400 : vecmini = cgetg(l, t_COL);
2170 2842 : for(i = 1; i<l; i++)
2171 : {
2172 1442 : GEN l = gel(factp,i), pm;
2173 1442 : if (i == 1 && absequaliu(l, 2)) { gel(vecmini,1) = gen_0; continue; }
2174 1421 : gel(vecmini,i) = pm = polymini(polr, l);
2175 1421 : polr = ZX_Q_mul(gel(pm,1), powiu(l, gel(pm,2)[3]));
2176 : }
2177 1400 : RgX_to_06(polr, &a0,&a1,&a2,&a3,&a4,&a5,&a6);
2178 1400 : I.j10 = !signe(a0)? mulii(sqri(a1), ZX_disc(polr)): ZX_disc(polr);
2179 1400 : I.j10 = gmul2n(I.j10,-12);
2180 :
2181 1400 : I.a0 = a0;
2182 1400 : I.A2 = apol2(a0,a1,a2);
2183 1400 : I.A3 = apol3(a0,a1,a2,a3);
2184 1400 : I.A5 = apol5(a0,a1,a2,a3,a4,a5);
2185 1400 : I.B2 = bpol2(a0,a1,a2,a3,a4);
2186 :
2187 1400 : I.j2 = igusaj2(a0,a1,a2,a3,a4,a5,a6);
2188 1400 : I.j4 = igusaj4(a0,a1,a2,a3,a4,a5,a6);
2189 1400 : I.i4 = gsub(gsqr(I.j2), gmulsg(24,I.j4));
2190 1400 : I.j6 = igusaj6(a0,a1,a2,a3,a4,a5,a6);
2191 1400 : j42 = gsqr(I.j4);
2192 1400 : j22 = gsqr(I.j2);
2193 1400 : j2j6 = gmul(I.j2,I.j6);
2194 1400 : I.j8 = gmul2n(gsub(j2j6,j42), -2);
2195 1400 : I.i12= gmul2n(gsub(gadd(gmul(j22,j42),gmulsg(36,gmul(j2j6,I.j4))),
2196 : gadd(gadd(gmulsg(32,gmul(j42,I.j4)),gmul(j2j6,j22)),gmulsg(108,gsqr(I.j6)))),-2);
2197 :
2198 2842 : for(i = 1; i < l; i++)
2199 1442 : gcoeff(facto,i,2) = stoi(Q_pval(I.j10, gel(factp,i)));
2200 1400 : dd = ZX_pval(polr,gen_2) & (~1); /* = 2 floor(val/2) */
2201 1400 : polr = gmul2n(polr, -dd);
2202 :
2203 1400 : V = cgetg(l, t_VEC);
2204 2842 : for (i = 1; i < l; i++)
2205 : {
2206 1442 : GEN q = gel(factp,i), red, N = NULL;
2207 : struct igusa_p Ip;
2208 1442 : long f = genus2localred(&I, &Ip, q, gel(vecmini,i));
2209 1442 : gcoeff(facto,i,2) = stoi(f);
2210 1442 : if (Ip.tame) Ip.type = stack_strcat("(tame) ", Ip.type);
2211 1442 : if (f >= 0)
2212 1421 : N = zv_snf(Ip.neron);
2213 1442 : if (DEBUGLEVEL)
2214 : {
2215 0 : if (!p) err_printf("p = %Ps\n", q);
2216 0 : err_printf("(potential) stable reduction: %Ps\n", Ip.stable);
2217 0 : if (f >= 0) {
2218 0 : err_printf("reduction at p: %s, %Ps", Ip.type, N);
2219 0 : err_printf(", f = %ld\n", f);
2220 : }
2221 : }
2222 1442 : red = f >= 0? mkvec2(strtoGENstr(Ip.type), N): cgetg(1, t_VEC);
2223 1442 : gel(V, i) = mkvec3(q, Ip.stable, red);
2224 : }
2225 1400 : if (p) V = gel(V,1);
2226 1400 : cond = factorback(facto);
2227 : /* remove denominator 2 coming from f = -1 in genuslocalred(, p = 2) */
2228 1400 : if (typ(cond) != t_INT) cond = gel(cond,1);
2229 1400 : return gerepilecopy(av, mkvec4(cond, facto, PQ, V));
2230 : }
|