Path: blob/master/sage/schemes/elliptic_curves/descent_two_isogeny.pyx
4158 views
"""1Descent on elliptic curves over QQ with a 2-isogeny.2"""34#*****************************************************************************5# Copyright (C) 2009 Robert L. Miller <[email protected]>6#7# Distributed under the terms of the GNU General Public License (GPL)8# http://www.gnu.org/licenses/9#*****************************************************************************1011from sage.rings.all import ZZ12from sage.rings.polynomial.polynomial_ring import polygen13cdef object x_ZZ = polygen(ZZ)14from sage.rings.polynomial.real_roots import real_roots15from sage.rings.arith import prime_divisors16from sage.misc.all import walltime, cputime17from sage.libs.pari.gen import pari18from sage.all import ntl1920from sage.rings.integer cimport Integer2122include "../../ext/cdefs.pxi"23include "../../ext/interrupt.pxi"24include "../../libs/flint/fmpz_poly.pxi"2526from sage.libs.flint.zmod_poly cimport *, zmod_poly_t27from sage.libs.flint.long_extras cimport *, factor_t28from sage.libs.ratpoints cimport ratpoints_mpz_exists_only2930cdef int N_RES_CLASSES_BSD = 103132cdef unsigned long ui0 = <unsigned long>033cdef unsigned long ui1 = <unsigned long>134cdef unsigned long ui2 = <unsigned long>235cdef unsigned long ui3 = <unsigned long>336cdef unsigned long ui4 = <unsigned long>437cdef unsigned long ui8 = <unsigned long>83839cdef unsigned long valuation(mpz_t a, mpz_t p):40"""41Return the number of times p divides a.42"""43cdef mpz_t aa44cdef unsigned long v45mpz_init(aa)46v = mpz_remove(aa,a,p)47mpz_clear(aa)48return v4950def test_valuation(a, p):51"""52Doctest function for cdef long valuation(mpz_t, mpz_t).5354EXAMPLE::5556sage: from sage.schemes.elliptic_curves.descent_two_isogeny import test_valuation as tv57sage: for i in [1..20]:58... print '%10s'%factor(i), tv(i,2), tv(i,3), tv(i,5)591 0 0 0602 1 0 0613 0 1 0622^2 2 0 0635 0 0 1642 * 3 1 1 0657 0 0 0662^3 3 0 0673^2 0 2 0682 * 5 1 0 16911 0 0 0702^2 * 3 2 1 07113 0 0 0722 * 7 1 0 0733 * 5 0 1 1742^4 4 0 07517 0 0 0762 * 3^2 1 2 07719 0 0 0782^2 * 5 2 0 17980"""81cdef Integer A = Integer(a)82cdef Integer P = Integer(p)83return valuation(A.value, P.value)8485cdef int padic_square(mpz_t a, mpz_t p):86"""87Test if a is a p-adic square.88"""89cdef unsigned long v90cdef mpz_t aa91cdef int result9293if mpz_sgn(a) == 0: return 19495v = valuation(a,p)96if v&(ui1): return 09798mpz_init_set(aa,a)99while v:100v -= 1101mpz_divexact(aa, aa, p)102if mpz_cmp_ui(p, ui2)==0:103result = ( mpz_fdiv_ui(aa,ui8)==ui1 )104else:105result = ( mpz_legendre(aa,p)==1 )106mpz_clear(aa)107return result108109def test_padic_square(a, p):110"""111Doctest function for cdef int padic_square(mpz_t, unsigned long).112113EXAMPLE::114115sage: from sage.schemes.elliptic_curves.descent_two_isogeny import test_padic_square as ps116sage: for i in [1..300]:117... for p in prime_range(100):118... if not Qp(p)(i).is_square()==bool(ps(i,p)):119... print i, p120121"""122cdef Integer A = Integer(a)123cdef Integer P = Integer(p)124return padic_square(A.value, P.value)125126cdef int lemma6(mpz_t a, mpz_t b, mpz_t c, mpz_t d, mpz_t e,127mpz_t x, mpz_t p, unsigned long nu):128"""129Implements Lemma 6 of BSD's "Notes on elliptic curves, I" for odd p.130131Returns -1 for insoluble, 0 for undecided, +1 for soluble.132"""133cdef mpz_t g_of_x, g_prime_of_x134cdef unsigned long lambd, mu135cdef int result = -1136137mpz_init(g_of_x)138mpz_mul(g_of_x, a, x)139mpz_add(g_of_x, g_of_x, b)140mpz_mul(g_of_x, g_of_x, x)141mpz_add(g_of_x, g_of_x, c)142mpz_mul(g_of_x, g_of_x, x)143mpz_add(g_of_x, g_of_x, d)144mpz_mul(g_of_x, g_of_x, x)145mpz_add(g_of_x, g_of_x, e)146147if padic_square(g_of_x, p):148mpz_clear(g_of_x)149return +1 # soluble150151mpz_init_set(g_prime_of_x, x)152mpz_mul(g_prime_of_x, a, x)153mpz_mul_ui(g_prime_of_x, g_prime_of_x, ui4)154mpz_addmul_ui(g_prime_of_x, b, ui3)155mpz_mul(g_prime_of_x, g_prime_of_x, x)156mpz_addmul_ui(g_prime_of_x, c, ui2)157mpz_mul(g_prime_of_x, g_prime_of_x, x)158mpz_add(g_prime_of_x, g_prime_of_x, d)159160lambd = valuation(g_of_x, p)161if mpz_sgn(g_prime_of_x)==0:162if lambd >= 2*nu: result = 0 # undecided163else:164mu = valuation(g_prime_of_x, p)165if lambd > 2*mu: result = +1 # soluble166elif lambd >= 2*nu and mu >= nu: result = 0 # undecided167168mpz_clear(g_prime_of_x)169mpz_clear(g_of_x)170return result171172cdef int lemma7(mpz_t a, mpz_t b, mpz_t c, mpz_t d, mpz_t e,173mpz_t x, mpz_t p, unsigned long nu):174"""175Implements Lemma 7 of BSD's "Notes on elliptic curves, I" for p=2.176177Returns -1 for insoluble, 0 for undecided, +1 for soluble.178"""179cdef mpz_t g_of_x, g_prime_of_x, g_of_x_odd_part180cdef unsigned long lambd, mu, g_of_x_odd_part_mod_4181cdef int result = -1182183mpz_init(g_of_x)184mpz_mul(g_of_x, a, x)185mpz_add(g_of_x, g_of_x, b)186mpz_mul(g_of_x, g_of_x, x)187mpz_add(g_of_x, g_of_x, c)188mpz_mul(g_of_x, g_of_x, x)189mpz_add(g_of_x, g_of_x, d)190mpz_mul(g_of_x, g_of_x, x)191mpz_add(g_of_x, g_of_x, e)192193if padic_square(g_of_x, p):194mpz_clear(g_of_x)195return +1 # soluble196197mpz_init_set(g_prime_of_x, x)198mpz_mul(g_prime_of_x, a, x)199mpz_mul_ui(g_prime_of_x, g_prime_of_x, ui4)200mpz_addmul_ui(g_prime_of_x, b, ui3)201mpz_mul(g_prime_of_x, g_prime_of_x, x)202mpz_addmul_ui(g_prime_of_x, c, ui2)203mpz_mul(g_prime_of_x, g_prime_of_x, x)204mpz_add(g_prime_of_x, g_prime_of_x, d)205206lambd = valuation(g_of_x, p)207mpz_init_set(g_of_x_odd_part, g_of_x)208while mpz_even_p(g_of_x_odd_part):209mpz_divexact_ui(g_of_x_odd_part, g_of_x_odd_part, ui2)210g_of_x_odd_part_mod_4 = mpz_fdiv_ui(g_of_x_odd_part, ui4)211if mpz_sgn(g_prime_of_x)==0:212if lambd >= 2*nu: result = 0 # undecided213elif lambd == 2*nu-2 and g_of_x_odd_part_mod_4==1:214result = 0 # undecided215else:216mu = valuation(g_prime_of_x, p)217if lambd > 2*mu: result = +1 # soluble218elif nu > mu:219if lambd >= mu+nu: result = +1 # soluble220elif lambd+1 == mu+nu and lambd&ui1==0:221result = +1 # soluble222elif lambd+2 == mu+nu and lambd&ui1==0 and g_of_x_odd_part_mod_4==1:223result = +1 # soluble224else: # nu <= mu225if lambd >= 2*nu: result = 0 # undecided226elif lambd+2 == 2*nu and g_of_x_odd_part_mod_4==1:227result = 0 # undecided228229mpz_clear(g_prime_of_x)230mpz_clear(g_of_x)231return result232233cdef int Zp_soluble_BSD(mpz_t a, mpz_t b, mpz_t c, mpz_t d, mpz_t e,234mpz_t x_k, mpz_t p, unsigned long k):235"""236Uses the approach of BSD's "Notes on elliptic curves, I" to test for237solubility of y^2 == ax^4 + bx^3 + cx^2 + dx + e over Zp, with238x=x_k (mod p^k).239"""240# returns solubility of y^2 = ax^4 + bx^3 + cx^2 + dx + e241# in Zp with x=x_k (mod p^k)242cdef int code243cdef unsigned long t244cdef mpz_t s245246if mpz_cmp_ui(p, ui2)==0:247code = lemma7(a,b,c,d,e,x_k,p,k)248else:249code = lemma6(a,b,c,d,e,x_k,p,k)250if code == 1:251return 1252if code == -1:253return 0254255# now code == 0256t = 0257mpz_init(s)258while code == 0 and mpz_cmp_ui(p, t) > 0 and t < N_RES_CLASSES_BSD:259mpz_pow_ui(s, p, k)260mpz_mul_ui(s, s, t)261mpz_add(s, s, x_k)262code = Zp_soluble_BSD(a,b,c,d,e,s,p,k+1)263t += 1264mpz_clear(s)265return code266267cdef bint Zp_soluble_siksek(mpz_t a, mpz_t b, mpz_t c, mpz_t d, mpz_t e,268mpz_t pp, unsigned long pp_ui,269zmod_poly_factor_t f_factzn, zmod_poly_t f,270fmpz_poly_t f1, fmpz_poly_t linear,271double pp_ui_inv):272"""273Uses the approach of Algorithm 5.3.1 of Siksek's thesis to test for274solubility of y^2 == ax^4 + bx^3 + cx^2 + dx + e over Zp.275"""276cdef unsigned long v_min, v277cdef unsigned long roots[4]278cdef int i, j, has_roots, has_single_roots279cdef bint result280281cdef mpz_t aa, bb, cc, dd, ee282cdef mpz_t aaa, bbb, ccc, ddd, eee283cdef unsigned long qq284cdef unsigned long rr, ss285cdef mpz_t tt286287# Step 0: divide out all common p from the quartic288v_min = valuation(a, pp)289if mpz_cmp_ui(b, ui0) != 0:290v = valuation(b, pp)291if v < v_min: v_min = v292if mpz_cmp_ui(c, ui0) != 0:293v = valuation(c, pp)294if v < v_min: v_min = v295if mpz_cmp_ui(d, ui0) != 0:296v = valuation(d, pp)297if v < v_min: v_min = v298if mpz_cmp_ui(e, ui0) != 0:299v = valuation(e, pp)300if v < v_min: v_min = v301for 0 <= v < v_min:302mpz_divexact(a, a, pp)303mpz_divexact(b, b, pp)304mpz_divexact(c, c, pp)305mpz_divexact(d, d, pp)306mpz_divexact(e, e, pp)307308if not v_min%2:309# Step I in Alg. 5.3.1 of Siksek's thesis310zmod_poly_set_coeff_ui(f, 0, mpz_fdiv_ui(e, pp_ui))311zmod_poly_set_coeff_ui(f, 1, mpz_fdiv_ui(d, pp_ui))312zmod_poly_set_coeff_ui(f, 2, mpz_fdiv_ui(c, pp_ui))313zmod_poly_set_coeff_ui(f, 3, mpz_fdiv_ui(b, pp_ui))314zmod_poly_set_coeff_ui(f, 4, mpz_fdiv_ui(a, pp_ui))315316result = 0317(<zmod_poly_factor_struct *>f_factzn)[0].num_factors = 0 # reset data struct318qq = zmod_poly_factor(f_factzn, f)319for i from 0 <= i < f_factzn.num_factors:320if f_factzn.exponents[i]&1:321result = 1322break323if result == 0 and z_legendre_precomp(qq, pp_ui, pp_ui_inv) == 1:324result = 1325if result:326return 1327328zmod_poly_zero(f)329zmod_poly_set_coeff_ui(f, 0, ui1)330for i from 0 <= i < f_factzn.num_factors:331for j from 0 <= j < (f_factzn.exponents[i]>>1):332zmod_poly_mul(f, f, f_factzn.factors[i])333334(<zmod_poly_factor_struct *>f_factzn)[0].num_factors = 0 # reset data struct335zmod_poly_factor(f_factzn, f)336has_roots = 0337j = 0338for i from 0 <= i < f_factzn.num_factors:339if zmod_poly_degree(f_factzn.factors[i]) == 1 and 0 != zmod_poly_get_coeff_ui(f_factzn.factors[i], 1):340has_roots = 1341roots[j] = pp_ui - zmod_poly_get_coeff_ui(f_factzn.factors[i], 0)342j += 1343if not has_roots:344return 0345346i = zmod_poly_degree(f)347mpz_init(aaa)348mpz_init(bbb)349mpz_init(ccc)350mpz_init(ddd)351mpz_init(eee)352353if i == 0: # g == 1354mpz_set(aaa, a)355mpz_set(bbb, b)356mpz_set(ccc, c)357mpz_set(ddd, d)358mpz_sub_ui(eee, e, qq)359elif i == 1: # g == x + rr360mpz_set(aaa, a)361mpz_set(bbb, b)362mpz_sub_ui(ccc, c, qq)363rr = zmod_poly_get_coeff_ui(f, 0)364ss = rr*qq365mpz_set(ddd,d)366mpz_sub_ui(ddd, ddd, ss*ui2)367mpz_set(eee,e)368mpz_sub_ui(eee, eee, ss*rr)369elif i == 2: # g == x^2 + rr*x + ss370mpz_sub_ui(aaa, a, qq)371rr = zmod_poly_get_coeff_ui(f, 1)372mpz_init(tt)373mpz_set_ui(tt, rr*qq)374mpz_set(bbb,b)375mpz_submul_ui(bbb, tt, ui2)376mpz_set(ccc,c)377mpz_submul_ui(ccc, tt, rr)378ss = zmod_poly_get_coeff_ui(f, 0)379mpz_set_ui(tt, ss*qq)380mpz_set(eee,e)381mpz_submul_ui(eee, tt, ss)382mpz_mul_ui(tt, tt, ui2)383mpz_sub(ccc, ccc, tt)384mpz_set(ddd,d)385mpz_submul_ui(ddd, tt, rr)386mpz_clear(tt)387mpz_divexact(aaa, aaa, pp)388mpz_divexact(bbb, bbb, pp)389mpz_divexact(ccc, ccc, pp)390mpz_divexact(ddd, ddd, pp)391mpz_divexact(eee, eee, pp)392# now aaa,bbb,ccc,ddd,eee represents h(x)393394result = 0395mpz_init(tt)396for i from 0 <= i < j:397mpz_mul_ui(tt, aaa, roots[i])398mpz_add(tt, tt, bbb)399mpz_mul_ui(tt, tt, roots[i])400mpz_add(tt, tt, ccc)401mpz_mul_ui(tt, tt, roots[i])402mpz_add(tt, tt, ddd)403mpz_mul_ui(tt, tt, roots[i])404mpz_add(tt, tt, eee)405# tt == h(r) mod p406mpz_mod(tt, tt, pp)407if mpz_sgn(tt) == 0:408fmpz_poly_zero(f1)409fmpz_poly_zero(linear)410fmpz_poly_set_coeff_mpz(f1, 0, e)411fmpz_poly_set_coeff_mpz(f1, 1, d)412fmpz_poly_set_coeff_mpz(f1, 2, c)413fmpz_poly_set_coeff_mpz(f1, 3, b)414fmpz_poly_set_coeff_mpz(f1, 4, a)415fmpz_poly_set_coeff_ui(linear, 0, roots[i])416fmpz_poly_set_coeff_mpz(linear, 1, pp)417fmpz_poly_compose(f1, f1, linear)418fmpz_poly_scalar_div_ui(f1, f1, pp_ui)419fmpz_poly_scalar_div_ui(f1, f1, pp_ui)420mpz_init(aa)421mpz_init(bb)422mpz_init(cc)423mpz_init(dd)424mpz_init(ee)425fmpz_poly_get_coeff_mpz(aa, f1, 4)426fmpz_poly_get_coeff_mpz(bb, f1, 3)427fmpz_poly_get_coeff_mpz(cc, f1, 2)428fmpz_poly_get_coeff_mpz(dd, f1, 1)429fmpz_poly_get_coeff_mpz(ee, f1, 0)430result = Zp_soluble_siksek(aa, bb, cc, dd, ee, pp, pp_ui, f_factzn, f, f1, linear, pp_ui_inv)431mpz_clear(aa)432mpz_clear(bb)433mpz_clear(cc)434mpz_clear(dd)435mpz_clear(ee)436if result == 1:437break438mpz_clear(aaa)439mpz_clear(bbb)440mpz_clear(ccc)441mpz_clear(ddd)442mpz_clear(eee)443mpz_clear(tt)444return result445else:446# Step II in Alg. 5.3.1 of Siksek's thesis447zmod_poly_set_coeff_ui(f, 0, mpz_fdiv_ui(e, pp_ui))448zmod_poly_set_coeff_ui(f, 1, mpz_fdiv_ui(d, pp_ui))449zmod_poly_set_coeff_ui(f, 2, mpz_fdiv_ui(c, pp_ui))450zmod_poly_set_coeff_ui(f, 3, mpz_fdiv_ui(b, pp_ui))451zmod_poly_set_coeff_ui(f, 4, mpz_fdiv_ui(a, pp_ui))452(<zmod_poly_factor_struct *>f_factzn)[0].num_factors = 0 # reset data struct453zmod_poly_factor(f_factzn, f)454has_roots = 0455has_single_roots = 0456j = 0457for i from 0 <= i < f_factzn.num_factors:458if zmod_poly_degree(f_factzn.factors[i]) == 1 and 0 != zmod_poly_get_coeff_ui(f_factzn.factors[i], 1):459has_roots = 1460if f_factzn.exponents[i] == 1:461has_single_roots = 1462break463roots[j] = pp_ui - zmod_poly_get_coeff_ui(f_factzn.factors[i], 0)464j += 1465466if not has_roots: return 0467if has_single_roots: return 1468469result = 0470if j > 0:471mpz_init(aa)472mpz_init(bb)473mpz_init(cc)474mpz_init(dd)475mpz_init(ee)476for i from 0 <= i < j:477fmpz_poly_zero(f1)478fmpz_poly_zero(linear)479fmpz_poly_set_coeff_mpz(f1, 0, e)480fmpz_poly_set_coeff_mpz(f1, 1, d)481fmpz_poly_set_coeff_mpz(f1, 2, c)482fmpz_poly_set_coeff_mpz(f1, 3, b)483fmpz_poly_set_coeff_mpz(f1, 4, a)484fmpz_poly_set_coeff_ui(linear, 0, roots[i])485fmpz_poly_set_coeff_mpz(linear, 1, pp)486fmpz_poly_compose(f1, f1, linear)487fmpz_poly_scalar_div_ui(f1, f1, pp_ui)488fmpz_poly_get_coeff_mpz(aa, f1, 4)489fmpz_poly_get_coeff_mpz(bb, f1, 3)490fmpz_poly_get_coeff_mpz(cc, f1, 2)491fmpz_poly_get_coeff_mpz(dd, f1, 1)492fmpz_poly_get_coeff_mpz(ee, f1, 0)493result = Zp_soluble_siksek(aa, bb, cc, dd, ee, pp, pp_ui, f_factzn, f, f1, linear, pp_ui_inv)494if result == 1:495break496if j > 0:497mpz_clear(aa)498mpz_clear(bb)499mpz_clear(cc)500mpz_clear(dd)501mpz_clear(ee)502return result503504cdef bint Zp_soluble_siksek_large_p(mpz_t a, mpz_t b, mpz_t c, mpz_t d, mpz_t e, mpz_t pp,505fmpz_poly_t f1, fmpz_poly_t linear):506"""507Uses the approach of Algorithm 5.3.1 of Siksek's thesis to test for508solubility of y^2 == ax^4 + bx^3 + cx^2 + dx + e over Zp.509"""510cdef unsigned long v_min, v511cdef mpz_t roots[4]512cdef int i, j, has_roots, has_single_roots513cdef bint result514515cdef mpz_t aa, bb, cc, dd, ee516cdef mpz_t aaa, bbb, ccc, ddd, eee517cdef mpz_t qq, rr, ss, tt518cdef Integer A,B,C,D,E,P519520# Step 0: divide out all common p from the quartic521v_min = valuation(a, pp)522if mpz_cmp_ui(b, ui0) != 0:523v = valuation(b, pp)524if v < v_min: v_min = v525if mpz_cmp_ui(c, ui0) != 0:526v = valuation(c, pp)527if v < v_min: v_min = v528if mpz_cmp_ui(d, ui0) != 0:529v = valuation(d, pp)530if v < v_min: v_min = v531if mpz_cmp_ui(e, ui0) != 0:532v = valuation(e, pp)533if v < v_min: v_min = v534for 0 <= v < v_min:535mpz_divexact(a, a, pp)536mpz_divexact(b, b, pp)537mpz_divexact(c, c, pp)538mpz_divexact(d, d, pp)539mpz_divexact(e, e, pp)540541if not v_min%2:542# Step I in Alg. 5.3.1 of Siksek's thesis543A = Integer(0); B = Integer(0); C = Integer(0); D = Integer(0); E = Integer(0); P = Integer(0)544mpz_set(A.value, a); mpz_set(B.value, b); mpz_set(C.value, c); mpz_set(D.value, d); mpz_set(E.value, e); mpz_set(P.value, pp)545f = ntl.ZZ_pX([E,D,C,B,A], P)546f /= ntl.ZZ_pX([A], P) # now f is monic, and we are done with A,B,C,D,E547mpz_set(qq, A.value) # qq is the leading coefficient of the polynomial548f_factzn = f.factor()549result = 0550for factor, exponent in f_factzn:551if exponent&1:552result = 1553break554if result == 0 and mpz_legendre(qq, pp) == 1:555result = 1556if result:557return 1558559f = ntl.ZZ_pX([1], P)560for factor, exponent in f_factzn:561for j from 0 <= j < (exponent/2):562f *= factor563564f /= f.leading_coefficient()565f_factzn = f.factor()566567has_roots = 0568j = 0569for factor, exponent in f_factzn:570if factor.degree() == 1:571has_roots = 1572A = P - Integer(factor[0])573mpz_set(roots[j], A.value)574j += 1575if not has_roots:576return 0577578i = f.degree()579mpz_init(aaa)580mpz_init(bbb)581mpz_init(ccc)582mpz_init(ddd)583mpz_init(eee)584585if i == 0: # g == 1586mpz_set(aaa, a)587mpz_set(bbb, b)588mpz_set(ccc, c)589mpz_set(ddd, d)590mpz_sub(eee, e, qq)591elif i == 1: # g == x + rr592mpz_set(aaa, a)593mpz_set(bbb, b)594mpz_sub(ccc, c, qq)595A = Integer(f[0])596mpz_set(rr, A.value)597mpz_mul(ss, rr, qq)598mpz_set(ddd,d)599mpz_sub(ddd, ddd, ss)600mpz_sub(ddd, ddd, ss)601mpz_set(eee,e)602mpz_mul(ss, ss, rr)603mpz_sub(eee, eee, ss)604mpz_divexact(ss, ss, rr)605elif i == 2: # g == x^2 + rr*x + ss606mpz_sub(aaa, a, qq)607A = Integer(f[1])608mpz_set(rr, A.value)609mpz_init(tt)610mpz_mul(tt, rr, qq)611mpz_set(bbb,b)612mpz_submul_ui(bbb, tt, ui2)613mpz_set(ccc,c)614mpz_submul(ccc, tt, rr)615A = Integer(f[0])616mpz_set(ss, A.value)617mpz_mul(tt, ss, qq)618mpz_set(eee,e)619mpz_submul(eee, tt, ss)620mpz_mul_ui(tt, tt, ui2)621mpz_sub(ccc, ccc, tt)622mpz_set(ddd,d)623mpz_submul(ddd, tt, rr)624mpz_clear(tt)625mpz_divexact(aaa, aaa, pp)626mpz_divexact(bbb, bbb, pp)627mpz_divexact(ccc, ccc, pp)628mpz_divexact(ddd, ddd, pp)629mpz_divexact(eee, eee, pp)630# now aaa,bbb,ccc,ddd,eee represents h(x)631632result = 0633mpz_init(tt)634for i from 0 <= i < j:635mpz_mul(tt, aaa, roots[i])636mpz_add(tt, tt, bbb)637mpz_mul(tt, tt, roots[i])638mpz_add(tt, tt, ccc)639mpz_mul(tt, tt, roots[i])640mpz_add(tt, tt, ddd)641mpz_mul(tt, tt, roots[i])642mpz_add(tt, tt, eee)643# tt == h(r) mod p644mpz_mod(tt, tt, pp)645if mpz_sgn(tt) == 0:646fmpz_poly_zero(f1)647fmpz_poly_zero(linear)648fmpz_poly_set_coeff_mpz(f1, 0, e)649fmpz_poly_set_coeff_mpz(f1, 1, d)650fmpz_poly_set_coeff_mpz(f1, 2, c)651fmpz_poly_set_coeff_mpz(f1, 3, b)652fmpz_poly_set_coeff_mpz(f1, 4, a)653fmpz_poly_set_coeff_mpz(linear, 0, roots[i])654fmpz_poly_set_coeff_mpz(linear, 1, pp)655fmpz_poly_compose(f1, f1, linear)656fmpz_poly_scalar_div_mpz(f1, f1, pp)657fmpz_poly_scalar_div_mpz(f1, f1, pp)658659mpz_init(aa)660mpz_init(bb)661mpz_init(cc)662mpz_init(dd)663mpz_init(ee)664fmpz_poly_get_coeff_mpz(aa, f1, 4)665fmpz_poly_get_coeff_mpz(bb, f1, 3)666fmpz_poly_get_coeff_mpz(cc, f1, 2)667fmpz_poly_get_coeff_mpz(dd, f1, 1)668fmpz_poly_get_coeff_mpz(ee, f1, 0)669result = Zp_soluble_siksek_large_p(aa, bb, cc, dd, ee, pp, f1, linear)670mpz_clear(aa)671mpz_clear(bb)672mpz_clear(cc)673mpz_clear(dd)674mpz_clear(ee)675if result == 1:676break677mpz_clear(aaa)678mpz_clear(bbb)679mpz_clear(ccc)680mpz_clear(ddd)681mpz_clear(eee)682mpz_clear(tt)683return result684else:685# Step II in Alg. 5.3.1 of Siksek's thesis686A = Integer(0); B = Integer(0); C = Integer(0); D = Integer(0); E = Integer(0); P = Integer(0)687mpz_set(A.value, a); mpz_set(B.value, b); mpz_set(C.value, c); mpz_set(D.value, d); mpz_set(E.value, e); mpz_set(P.value, pp)688f = ntl.ZZ_pX([E,D,C,B,A], P)689f /= ntl.ZZ_pX([A], P) # now f is monic690f_factzn = f.factor()691692has_roots = 0693has_single_roots = 0694j = 0695for factor, exponent in f_factzn:696if factor.degree() == 1:697has_roots = 1698if exponent == 1:699has_single_roots = 1700break701A = P - Integer(factor[0])702mpz_set(roots[j], A.value)703j += 1704705if not has_roots: return 0706if has_single_roots: return 1707708result = 0709if j > 0:710mpz_init(aa)711mpz_init(bb)712mpz_init(cc)713mpz_init(dd)714mpz_init(ee)715for i from 0 <= i < j:716fmpz_poly_zero(f1)717fmpz_poly_zero(linear)718fmpz_poly_set_coeff_mpz(f1, 0, e)719fmpz_poly_set_coeff_mpz(f1, 1, d)720fmpz_poly_set_coeff_mpz(f1, 2, c)721fmpz_poly_set_coeff_mpz(f1, 3, b)722fmpz_poly_set_coeff_mpz(f1, 4, a)723fmpz_poly_set_coeff_mpz(linear, 0, roots[i])724fmpz_poly_set_coeff_mpz(linear, 1, pp)725fmpz_poly_compose(f1, f1, linear)726fmpz_poly_scalar_div_mpz(f1, f1, pp)727fmpz_poly_get_coeff_mpz(aa, f1, 4)728fmpz_poly_get_coeff_mpz(bb, f1, 3)729fmpz_poly_get_coeff_mpz(cc, f1, 2)730fmpz_poly_get_coeff_mpz(dd, f1, 1)731fmpz_poly_get_coeff_mpz(ee, f1, 0)732result = Zp_soluble_siksek_large_p(aa, bb, cc, dd, ee, pp, f1, linear)733if result == 1:734break735if j > 0:736mpz_clear(aa)737mpz_clear(bb)738mpz_clear(cc)739mpz_clear(dd)740mpz_clear(ee)741return result742743cdef bint Qp_soluble_siksek(mpz_t A, mpz_t B, mpz_t C, mpz_t D, mpz_t E,744mpz_t p, unsigned long P,745zmod_poly_factor_t f_factzn, fmpz_poly_t f1,746fmpz_poly_t linear):747"""748Uses Samir Siksek's thesis results to determine whether the quartic is749locally soluble at p.750"""751cdef int result = 0752cdef mpz_t a,b,c,d,e753cdef zmod_poly_t f754cdef double pp_ui_inv = z_precompute_inverse(P)755zmod_poly_init(f, P)756757mpz_init_set(a,A)758mpz_init_set(b,B)759mpz_init_set(c,C)760mpz_init_set(d,D)761mpz_init_set(e,E)762763if Zp_soluble_siksek(a,b,c,d,e,p,P,f_factzn, f, f1, linear, pp_ui_inv):764result = 1765else:766mpz_set(a,A)767mpz_set(b,B)768mpz_set(c,C)769mpz_set(d,D)770mpz_set(e,E)771if Zp_soluble_siksek(e,d,c,b,a,p,P,f_factzn, f, f1, linear, pp_ui_inv):772result = 1773774mpz_clear(a)775mpz_clear(b)776mpz_clear(c)777mpz_clear(d)778mpz_clear(e)779zmod_poly_clear(f)780return result781782cdef bint Qp_soluble_siksek_large_p(mpz_t A, mpz_t B, mpz_t C, mpz_t D, mpz_t E,783mpz_t p, fmpz_poly_t f1, fmpz_poly_t linear):784"""785Uses Samir Siksek's thesis results to determine whether the quartic is786locally soluble at p, when p is bigger than wordsize, and we can't use787FLINT.788"""789cdef int result = 0790cdef mpz_t a,b,c,d,e791792mpz_init_set(a,A)793mpz_init_set(b,B)794mpz_init_set(c,C)795mpz_init_set(d,D)796mpz_init_set(e,E)797798if Zp_soluble_siksek_large_p(a,b,c,d,e,p,f1,linear):799result = 1800else:801mpz_set(a,A)802mpz_set(b,B)803mpz_set(c,C)804mpz_set(d,D)805mpz_set(e,E)806if Zp_soluble_siksek_large_p(e,d,c,b,a,p,f1,linear):807result = 1808809mpz_clear(a)810mpz_clear(b)811mpz_clear(c)812mpz_clear(d)813mpz_clear(e)814return result815816cdef bint Qp_soluble_BSD(mpz_t a, mpz_t b, mpz_t c, mpz_t d, mpz_t e, mpz_t p):817"""818Uses the original test of Birch and Swinnerton-Dyer to test for local819solubility of the quartic at p.820"""821cdef mpz_t zero822cdef int result = 0823mpz_init_set_ui(zero, ui0)824if Zp_soluble_BSD(a,b,c,d,e,zero,p,0):825result = 1826elif Zp_soluble_BSD(e,d,c,b,a,zero,p,1):827result = 1828mpz_clear(zero)829return result830831cdef bint Qp_soluble(mpz_t a, mpz_t b, mpz_t c, mpz_t d, mpz_t e, mpz_t p):832"""833Try the BSD approach for a few residue classes and if no solution is found,834switch to Siksek to try to prove insolubility.835"""836cdef int bsd_sol, sik_sol837cdef unsigned long pp838cdef fmpz_poly_t f1, linear839cdef zmod_poly_factor_t f_factzn840bsd_sol = Qp_soluble_BSD(a, b, c, d, e, p)841if mpz_cmp_ui(p,N_RES_CLASSES_BSD)>0 and not bsd_sol:842fmpz_poly_init(f1)843fmpz_poly_init(linear)844if mpz_fits_ulong_p(p):845zmod_poly_factor_init(f_factzn)846pp = mpz_get_ui(p)847sik_sol = Qp_soluble_siksek(a,b,c,d,e,p,pp,f_factzn,f1,linear)848zmod_poly_factor_clear(f_factzn)849else:850sik_sol = Qp_soluble_siksek_large_p(a,b,c,d,e,p,f1,linear)851fmpz_poly_clear(f1)852fmpz_poly_clear(linear)853else:854sik_sol = bsd_sol855return sik_sol856857def test_qpls(a,b,c,d,e,p):858"""859Testing function for Qp_soluble.860861EXAMPLE:862sage: from sage.schemes.elliptic_curves.descent_two_isogeny import test_qpls as tq863sage: tq(1,2,3,4,5,7)8641865866"""867cdef Integer A,B,C,D,E,P868cdef int i, result869cdef mpz_t aa,bb,cc,dd,ee,pp870A=Integer(a); B=Integer(b); C=Integer(c); D=Integer(d); E=Integer(e); P=Integer(p)871mpz_init_set(aa, A.value)872mpz_init_set(bb, B.value)873mpz_init_set(cc, C.value)874mpz_init_set(dd, D.value)875mpz_init_set(ee, E.value)876mpz_init_set(pp, P.value)877result = Qp_soluble(aa, bb, cc, dd, ee, pp)878mpz_clear(aa)879mpz_clear(bb)880mpz_clear(cc)881mpz_clear(dd)882mpz_clear(ee)883mpz_clear(pp)884return result885886cdef int everywhere_locally_soluble(mpz_t a, mpz_t b, mpz_t c, mpz_t d, mpz_t e) except -1:887"""888Returns whether the quartic has local solutions at all primes p.889"""890cdef Integer A,B,C,D,E,Delta,p891cdef mpz_t mpz_2892A=Integer(0); B=Integer(0); C=Integer(0); D=Integer(0); E=Integer(0)893mpz_set(A.value, a); mpz_set(B.value, b); mpz_set(C.value, c); mpz_set(D.value, d); mpz_set(E.value, e);894f = (((A*x_ZZ + B)*x_ZZ + C)*x_ZZ + D)*x_ZZ + E895896# RR soluble:897if mpz_sgn(a)!=1:898if len(real_roots(f)) == 0:899return 0900901# Q2 soluble:902mpz_init_set_ui(mpz_2, ui2)903if not Qp_soluble(a,b,c,d,e,mpz_2):904mpz_clear(mpz_2)905return 0906mpz_clear(mpz_2)907908# Odd finite primes909Delta = f.discriminant()910for p in prime_divisors(Delta):911if p == 2: continue912if not Qp_soluble(a,b,c,d,e,p.value): return 0913914return 1915916def test_els(a,b,c,d,e):917"""918Doctest function for cdef int everywhere_locally_soluble(mpz_t, mpz_t, mpz_t, mpz_t, mpz_t).919920EXAMPLE::921922sage: from sage.schemes.elliptic_curves.descent_two_isogeny import test_els923sage: from sage.libs.ratpoints import ratpoints924sage: for _ in range(1000):925... a,b,c,d,e = randint(1,1000), randint(1,1000), randint(1,1000), randint(1,1000), randint(1,1000)926... if len(ratpoints([e,d,c,b,a], 1000)) > 0:927... try:928... if not test_els(a,b,c,d,e):929... print "This never happened", a,b,c,d,e930... except ValueError:931... continue932933"""934cdef Integer A,B,C,D,E,Delta935A=Integer(a); B=Integer(b); C=Integer(c); D=Integer(d); E=Integer(e)936return everywhere_locally_soluble(A.value, B.value, C.value, D.value, E.value)937938cdef int count(mpz_t c_mpz, mpz_t d_mpz, mpz_t *p_list, unsigned long p_list_len,939int global_limit_small, int global_limit_large,940int verbosity, bint selmer_only, mpz_t n1, mpz_t n2) except -1:941"""942Count the number of els/gls quartic 2-covers of E.943"""944cdef unsigned long n_primes, i945cdef bint found_global_points, els, check_negs, verbose = (verbosity > 4)946cdef Integer a_Int, c_Int, e_Int947cdef mpz_t c_sq_mpz, d_prime_mpz948cdef mpz_t n_divisors, j949cdef mpz_t *coeffs_ratp950951952mpz_init(c_sq_mpz)953mpz_mul(c_sq_mpz, c_mpz, c_mpz)954mpz_init_set(d_prime_mpz, c_sq_mpz)955mpz_submul_ui(d_prime_mpz, d_mpz, ui4)956check_negs = 0957if mpz_sgn(d_prime_mpz) > 0:958if mpz_sgn(c_mpz) >= 0 or mpz_cmp(c_sq_mpz, d_prime_mpz) <= 0:959check_negs = 1960mpz_clear(c_sq_mpz)961mpz_clear(d_prime_mpz)962963964# Set up coefficient array, and static variables965cdef mpz_t *coeffs = <mpz_t *> sage_malloc(5 * sizeof(mpz_t))966for i from 0 <= i <= 4:967mpz_init(coeffs[i])968mpz_set_ui(coeffs[1], ui0) #969mpz_set(coeffs[2], c_mpz) # These never change970mpz_set_ui(coeffs[3], ui0) #971972if not selmer_only:973# allocate space for ratpoints974coeffs_ratp = <mpz_t *> sage_malloc(5 * sizeof(mpz_t))975for i from 0 <= i <= 4:976mpz_init(coeffs_ratp[i])977978# Get prime divisors, and put them in an mpz_t array979# (this block, by setting check_negs, takes care of980# local solubility over RR)981cdef mpz_t *p_div_d_mpz = <mpz_t *> sage_malloc((p_list_len+1) * sizeof(mpz_t))982n_primes = 0983for i from 0 <= i < p_list_len:984if mpz_divisible_p(d_mpz, p_list[i]):985mpz_init(p_div_d_mpz[n_primes])986mpz_set(p_div_d_mpz[n_primes], p_list[i])987n_primes += 1988if check_negs:989mpz_init(p_div_d_mpz[n_primes])990mpz_set_si(p_div_d_mpz[n_primes], -1)991n_primes += 1992mpz_init_set_ui(n_divisors, ui1)993mpz_mul_2exp(n_divisors, n_divisors, n_primes)994# if verbosity > 3:995# print '\nDivisors of d which may lead to RR-soluble quartics:', p_div_d996997mpz_init_set_ui(j, ui0)998if not selmer_only:999mpz_set_ui(n1, ui0)1000mpz_set_ui(n2, ui0)1001while mpz_cmp(j, n_divisors) < 0:1002mpz_set_ui(coeffs[4], ui1)1003for i from 0 <= i < n_primes:1004if mpz_tstbit(j, i):1005mpz_mul(coeffs[4], coeffs[4], p_div_d_mpz[i])1006if verbosity > 3:1007a_Int = Integer(0); mpz_set(a_Int.value, coeffs[4])1008print '\nSquarefree divisor:', a_Int1009mpz_divexact(coeffs[0], d_mpz, coeffs[4])1010found_global_points = 01011if not selmer_only:1012if verbose:1013print "\nCalling ratpoints for small point search"1014for i from 0 <= i <= 4:1015mpz_set(coeffs_ratp[i], coeffs[i])1016sig_on()1017found_global_points = ratpoints_mpz_exists_only(coeffs_ratp, global_limit_small, 4, verbose)1018sig_off()1019if found_global_points:1020if verbosity > 2:1021a_Int = Integer(0); mpz_set(a_Int.value, coeffs[4])1022c_Int = Integer(0); mpz_set(c_Int.value, coeffs[2])1023e_Int = Integer(0); mpz_set(e_Int.value, coeffs[0])1024print 'Found small global point, quartic (%d,%d,%d,%d,%d)'%(a_Int,0,c_Int,0,e_Int)1025mpz_add_ui(n1, n1, ui1)1026mpz_add_ui(n2, n2, ui1)1027if verbose:1028print "\nDone calling ratpoints for small point search"1029if not found_global_points:1030# Test whether the quartic is everywhere locally soluble:1031els = 11032for i from 0 <= i < p_list_len:1033if not Qp_soluble(coeffs[4], coeffs[3], coeffs[2], coeffs[1], coeffs[0], p_list[i]):1034els = 01035break1036if els:1037if verbosity > 2:1038a_Int = Integer(0); mpz_set(a_Int.value, coeffs[4])1039c_Int = Integer(0); mpz_set(c_Int.value, coeffs[2])1040e_Int = Integer(0); mpz_set(e_Int.value, coeffs[0])1041print 'ELS without small global points, quartic (%d,%d,%d,%d,%d)'%(a_Int,0,c_Int,0,e_Int)1042mpz_add_ui(n2, n2, ui1)1043if not selmer_only:1044if verbose:1045print "\nCalling ratpoints for large point search"1046for i from 0 <= i <= 4:1047mpz_set(coeffs_ratp[i], coeffs[i])1048sig_on()1049found_global_points = ratpoints_mpz_exists_only(coeffs_ratp, global_limit_large, 4, verbose)1050sig_off()1051if found_global_points:1052if verbosity > 2:1053print ' -- Found large global point.'1054mpz_add_ui(n1, n1, ui1)1055if verbose:1056print "\nDone calling ratpoints for large point search"1057mpz_add_ui(j, j, ui1)1058if not selmer_only:1059for i from 0 <= i <= 4:1060mpz_clear(coeffs_ratp[i])1061sage_free(coeffs_ratp)1062mpz_clear(j)1063for i from 0 <= i < n_primes:1064mpz_clear(p_div_d_mpz[i])1065sage_free(p_div_d_mpz)1066mpz_clear(n_divisors)1067for i from 0 <= i <= 4:1068mpz_clear(coeffs[i])1069sage_free(coeffs)1070return 010711072def two_descent_by_two_isogeny(E,1073int global_limit_small = 10,1074int global_limit_large = 10000,1075int verbosity = 0,1076bint selmer_only = 0, bint proof = 1):1077"""1078Given an elliptic curve E with a two-isogeny phi : E --> E' and dual isogeny1079phi', runs a two-isogeny descent on E, returning n1, n2, n1' and n2'. Here1080n1 is the number of quartic covers found with a rational point, and n2 is1081the number which are ELS.10821083EXAMPLES::10841085sage: from sage.schemes.elliptic_curves.descent_two_isogeny import two_descent_by_two_isogeny1086sage: E = EllipticCurve('14a')1087sage: n1, n2, n1_prime, n2_prime = two_descent_by_two_isogeny(E)1088sage: log(n1,2) + log(n1_prime,2) - 2 # the rank108901090sage: E = EllipticCurve('65a')1091sage: n1, n2, n1_prime, n2_prime = two_descent_by_two_isogeny(E)1092sage: log(n1,2) + log(n1_prime,2) - 2 # the rank109311094sage: x,y = var('x,y')1095sage: E = EllipticCurve(y^2 == x^3 + x^2 - 25*x + 39)1096sage: n1, n2, n1_prime, n2_prime = two_descent_by_two_isogeny(E)1097sage: log(n1,2) + log(n1_prime,2) - 2 # the rank109821099sage: E = EllipticCurve(y^2 + x*y + y == x^3 - 131*x + 558)1100sage: n1, n2, n1_prime, n2_prime = two_descent_by_two_isogeny(E)1101sage: log(n1,2) + log(n1_prime,2) - 2 # the rank1102311031104Using the verbosity option::11051106sage: E = EllipticCurve('14a')1107sage: two_descent_by_two_isogeny(E, verbosity=1)11082-isogeny1109Results:11102 <= #E(Q)/phi'(E'(Q)) <= 211112 <= #E'(Q)/phi(E(Q)) <= 21112#Sel^(phi')(E'/Q) = 21113#Sel^(phi)(E/Q) = 211141 <= #Sha(E'/Q)[phi'] <= 111151 <= #Sha(E/Q)[phi] <= 111161 <= #Sha(E/Q)[2], #Sha(E'/Q)[2] <= 111170 <= rank of E(Q) = rank of E'(Q) <= 01118(2, 2, 2, 2)11191120Handling curves whose discriminants involve larger than wordsize primes::11211122sage: E = EllipticCurve('14a')1123sage: E = E.quadratic_twist(next_prime(10^20))1124sage: E1125Elliptic Curve defined by y^2 = x^3 + x^2 + 716666666666666667225666666666666666775672*x - 391925925925925926384240370370370370549019837037037037060249356 over Rational Field1126sage: E.discriminant().factor()1127-1 * 2^18 * 7^3 * 100000000000000000039^61128sage: log(100000000000000000039.0, 2.0)112966.438...1130sage: n1, n2, n1_prime, n2_prime = two_descent_by_two_isogeny(E)1131sage: log(n1,2) + log(n1_prime,2) - 2 # the rank1132011331134TESTS:11351136Here we contrive an example to demonstrate that a keyboard interrupt1137is caught. Here we let `E` be the smallest optimal curve with two-torsion1138and nontrivial `Sha[2]`. This ensures that the two-descent will be looking1139for rational points which do not exist, and by setting global_limit_large1140to a very high bound, it will still be working when we simulate a ``CTRL-C``::11411142sage: from sage.schemes.elliptic_curves.descent_two_isogeny import two_descent_by_two_isogeny1143sage: import sage.tests.interrupt1144sage: E = EllipticCurve('960d'); E1145Elliptic Curve defined by y^2 = x^3 - x^2 - 900*x - 10098 over Rational Field1146sage: E.sha().an()114741148sage: try:1149... sage.tests.interrupt.interrupt_after_delay(1000)1150... two_descent_by_two_isogeny(E, global_limit_large=10^8)1151... except KeyboardInterrupt:1152... print "Caught!"1153Caught!1154"""1155cdef Integer a1, a2, a3, a4, a6, s2, s4, s61156cdef Integer c, d, x01157cdef list x_list1158assert E.torsion_order()%2==0, 'Need rational two-torsion for isogeny descent.'1159if verbosity > 0:1160print '\n2-isogeny'1161if verbosity > 1:1162print '\nchanging coordinates'1163a1 = Integer(E.a1())1164a2 = Integer(E.a2())1165a3 = Integer(E.a3())1166a4 = Integer(E.a4())1167a6 = Integer(E.a6())1168if a1==0 and a3==0:1169s2=a2; s4=a4; s6=a61170else:1171s2=a1*a1+4*a2; s4=8*(a1*a3+2*a4); s6=16*(a3*a3+4*a6)1172f = ((x_ZZ + s2)*x_ZZ + s4)*x_ZZ + s61173x_list = f.roots() # over ZZ -- use FLINT directly?1174x0 = x_list[0][0]1175c = 3*x0+s2; d = (c+s2)*x0+s41176return two_descent_by_two_isogeny_work(c, d,1177global_limit_small, global_limit_large, verbosity, selmer_only, proof)11781179def two_descent_by_two_isogeny_work(Integer c, Integer d,1180int global_limit_small = 10, int global_limit_large = 10000,1181int verbosity = 0, bint selmer_only = 0, bint proof = 1):1182"""1183Do all the work in doing a two-isogeny descent.11841185EXAMPLES::11861187sage: from sage.schemes.elliptic_curves.descent_two_isogeny import two_descent_by_two_isogeny_work1188sage: n1, n2, n1_prime, n2_prime = two_descent_by_two_isogeny_work(13,128)1189sage: log(n1,2) + log(n1_prime,2) - 2 # the rank119001191sage: n1, n2, n1_prime, n2_prime = two_descent_by_two_isogeny_work(1,-16)1192sage: log(n1,2) + log(n1_prime,2) - 2 # the rank119311194sage: n1, n2, n1_prime, n2_prime = two_descent_by_two_isogeny_work(10,8)1195sage: log(n1,2) + log(n1_prime,2) - 2 # the rank119621197sage: n1, n2, n1_prime, n2_prime = two_descent_by_two_isogeny_work(85,320)1198sage: log(n1,2) + log(n1_prime,2) - 2 # the rank1199312001201"""1202cdef mpz_t c_mpz, d_mpz, c_prime_mpz, d_prime_mpz1203cdef mpz_t *p_list_mpz1204cdef unsigned long i, j, p, p_list_len1205cdef Integer P, n1, n2, n1_prime, n2_prime, c_prime, d_prime1206cdef object PO1207cdef bint found, too_big, d_neg, d_prime_neg1208cdef factor_t fact1209cdef list primes1210mpz_init_set(c_mpz, c.value) #1211mpz_init_set(d_mpz, d.value) #1212mpz_init(c_prime_mpz) #1213mpz_init(d_prime_mpz) #1214mpz_mul_si(c_prime_mpz, c_mpz, -2)1215mpz_mul(d_prime_mpz, c_mpz, c_mpz)1216mpz_submul_ui(d_prime_mpz, d_mpz, ui4)12171218d_neg = 01219d_prime_neg = 01220if mpz_sgn(d_mpz) < 0:1221d_neg = 11222mpz_neg(d_mpz, d_mpz)1223if mpz_sgn(d_prime_mpz) < 0:1224d_prime_neg = 11225mpz_neg(d_prime_mpz, d_prime_mpz)1226if mpz_fits_ulong_p(d_mpz) and mpz_fits_ulong_p(d_prime_mpz):1227# Factor very quickly using FLINT.1228p_list_mpz = <mpz_t *> sage_malloc(20 * sizeof(mpz_t))1229mpz_init_set_ui(p_list_mpz[0], ui2)1230p_list_len = 11231z_factor(&fact, mpz_get_ui(d_mpz), proof)1232for i from 0 <= i < fact.num:1233p = fact.p[i]1234if p != ui2:1235mpz_init_set_ui(p_list_mpz[p_list_len], p)1236p_list_len += 11237z_factor(&fact, mpz_get_ui(d_prime_mpz), proof)1238for i from 0 <= i < fact.num:1239p = fact.p[i]1240found = 01241for j from 0 <= j < p_list_len:1242if mpz_cmp_ui(p_list_mpz[j], p)==0:1243found = 11244break1245if not found:1246mpz_init_set_ui(p_list_mpz[p_list_len], p)1247p_list_len += 11248else:1249# Factor more slowly using Pari via Python.1250d = Integer(0)1251mpz_set(d.value, d_mpz)1252primes = list(pari(d).factor()[0])1253d_prime = Integer(0)1254mpz_set(d_prime.value, d_prime_mpz)1255for PO in pari(d_prime).factor()[0]:1256if PO not in primes:1257primes.append(PO)1258P = Integer(2)1259if P not in primes: primes.append(P)1260p_list_len = len(primes)1261p_list_mpz = <mpz_t *> sage_malloc(p_list_len * sizeof(mpz_t))1262for i from 0 <= i < p_list_len:1263P = Integer(primes[i])1264mpz_init_set(p_list_mpz[i], P.value)1265if d_neg:1266mpz_neg(d_mpz, d_mpz)1267if d_prime_neg:1268mpz_neg(d_prime_mpz, d_prime_mpz)12691270if verbosity > 1:1271c_prime = -2*c1272d_prime = c*c-4*d1273print '\nnew curve is y^2 == x( x^2 + (%d)x + (%d) )'%(int(c),int(d))1274print 'new isogenous curve is' + \1275' y^2 == x( x^2 + (%d)x + (%d) )'%(int(c_prime),int(d_prime))12761277n1 = Integer(0); n2 = Integer(0)1278n1_prime = Integer(0); n2_prime = Integer(0)1279count(c.value, d.value, p_list_mpz, p_list_len,1280global_limit_small, global_limit_large, verbosity, selmer_only,1281n1.value, n2.value)1282count(c_prime_mpz, d_prime_mpz, p_list_mpz, p_list_len,1283global_limit_small, global_limit_large, verbosity, selmer_only,1284n1_prime.value, n2_prime.value)12851286for i from 0 <= i < p_list_len:1287mpz_clear(p_list_mpz[i])1288sage_free(p_list_mpz)12891290if verbosity > 0:1291print "\nResults:"1292print n1, "<= #E(Q)/phi'(E'(Q)) <=", n21293print n1_prime, "<= #E'(Q)/phi(E(Q)) <=", n2_prime1294print "#Sel^(phi')(E'/Q) =", n21295print "#Sel^(phi)(E/Q) =", n2_prime1296print "1 <= #Sha(E'/Q)[phi'] <=", n2/n11297print "1 <= #Sha(E/Q)[phi] <=", n2_prime/n1_prime1298print "1 <= #Sha(E/Q)[2], #Sha(E'/Q)[2] <=", (n2_prime/n1_prime)*(n2/n1)1299a = Integer(n1*n1_prime).log(Integer(2))1300e = Integer(n2*n2_prime).log(Integer(2))1301print a - 2, "<= rank of E(Q) = rank of E'(Q) <=", e - 213021303return n1, n2, n1_prime, n2_prime13041305130613071308