Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
Download

Testing latest pari + WASM + node.js... and it works?! Wow.

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