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