Path: blob/master/src/sage/schemes/elliptic_curves/constructor.py
8820 views
"""1Elliptic curve constructor23AUTHORS:45- William Stein (2005): Initial version67- John Cremona (2008-01): EllipticCurve(j) fixed for all cases8"""910#*****************************************************************************11# Copyright (C) 2005 William Stein <[email protected]>12#13# Distributed under the terms of the GNU General Public License (GPL)14#15# This code is distributed in the hope that it will be useful,16# but WITHOUT ANY WARRANTY; without even the implied warranty of17# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU18# General Public License for more details.19#20# The full text of the GPL is available at:21#22# http://www.gnu.org/licenses/23#*****************************************************************************242526import sage.rings.all as rings2728from sage.rings.finite_rings.integer_mod_ring import is_IntegerModRing29from sage.rings.rational_field import is_RationalField30from sage.rings.polynomial.multi_polynomial_ring import is_MPolynomialRing31from sage.rings.finite_rings.constructor import is_FiniteField32from sage.rings.number_field.number_field import is_NumberField33from sage.rings.polynomial.multi_polynomial_element import is_MPolynomial34from sage.rings.ring import is_Ring35from sage.rings.ring_element import is_RingElement3637from sage.categories.fields import Fields38_Fields = Fields()3940from sage.structure.sequence import Sequence41from sage.structure.element import parent42from sage.symbolic.ring import SR43from sage.symbolic.expression import is_SymbolicEquation444546def EllipticCurve(x=None, y=None, j=None, minimal_twist=True):47r"""48Construct an elliptic curve.4950In Sage, an elliptic curve is always specified by its a-invariants5152.. math::5354y^2 + a_1 xy + a_3 y = x^3 + a_2 x^2 + a_4 x + a_6.5556INPUT:5758There are several ways to construct an elliptic curve:5960- ``EllipticCurve([a1,a2,a3,a4,a6])``: Elliptic curve with given61a-invariants. The invariants are coerced into the parent of the62first element. If all are integers, they are coerced into the63rational numbers.6465- ``EllipticCurve([a4,a6])``: Same as above, but `a_1=a_2=a_3=0`.6667- ``EllipticCurve(label)``: Returns the elliptic curve over Q from68the Cremona database with the given label. The label is a69string, such as ``"11a"`` or ``"37b2"``. The letters in the70label *must* be lower case (Cremona's new labeling).7172- ``EllipticCurve(R, [a1,a2,a3,a4,a6])``: Create the elliptic73curve over ``R`` with given a-invariants. Here ``R`` can be an74arbitrary ring. Note that addition need not be defined.7576- ``EllipticCurve(j=j0)`` or ``EllipticCurve_from_j(j0)``: Return77an elliptic curve with j-invariant ``j0``.7879- ``EllipticCurve(polynomial)``: Read off the a-invariants from80the polynomial coefficients, see81:func:`EllipticCurve_from_Weierstrass_polynomial`.8283In each case above where the input is a list of length 2 or 5, one84can instead give a 2 or 5-tuple instead.8586EXAMPLES:8788We illustrate creating elliptic curves::8990sage: EllipticCurve([0,0,1,-1,0])91Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field9293We create a curve from a Cremona label::9495sage: EllipticCurve('37b2')96Elliptic Curve defined by y^2 + y = x^3 + x^2 - 1873*x - 31833 over Rational Field97sage: EllipticCurve('5077a')98Elliptic Curve defined by y^2 + y = x^3 - 7*x + 6 over Rational Field99sage: EllipticCurve('389a')100Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field101102Old Cremona labels are allowed::103104sage: EllipticCurve('2400FF')105Elliptic Curve defined by y^2 = x^3 + x^2 + 2*x + 8 over Rational Field106107Unicode labels are allowed::108109sage: EllipticCurve(u'389a')110Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field111112We create curves over a finite field as follows::113114sage: EllipticCurve([GF(5)(0),0,1,-1,0])115Elliptic Curve defined by y^2 + y = x^3 + 4*x over Finite Field of size 5116sage: EllipticCurve(GF(5), [0, 0,1,-1,0])117Elliptic Curve defined by y^2 + y = x^3 + 4*x over Finite Field of size 5118119Elliptic curves over `\ZZ/N\ZZ` with `N` prime are of type120"elliptic curve over a finite field"::121122sage: F = Zmod(101)123sage: EllipticCurve(F, [2, 3])124Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Ring of integers modulo 101125sage: E = EllipticCurve([F(2), F(3)])126sage: type(E)127<class 'sage.schemes.elliptic_curves.ell_finite_field.EllipticCurve_finite_field_with_category'>128sage: E.category()129Category of schemes over Ring of integers modulo 101130131In contrast, elliptic curves over `\ZZ/N\ZZ` with `N` composite132are of type "generic elliptic curve"::133134sage: F = Zmod(95)135sage: EllipticCurve(F, [2, 3])136Elliptic Curve defined by y^2 = x^3 + 2*x + 3 over Ring of integers modulo 95137sage: E = EllipticCurve([F(2), F(3)])138sage: type(E)139<class 'sage.schemes.elliptic_curves.ell_generic.EllipticCurve_generic_with_category'>140sage: E.category()141Category of schemes over Ring of integers modulo 95142143The following is a curve over the complex numbers::144145sage: E = EllipticCurve(CC, [0,0,1,-1,0])146sage: E147Elliptic Curve defined by y^2 + 1.00000000000000*y = x^3 + (-1.00000000000000)*x over Complex Field with 53 bits of precision148sage: E.j_invariant()1492988.97297297297150151We can also create elliptic curves by giving the Weierstrass equation::152153sage: x, y = var('x,y')154sage: EllipticCurve(y^2 + y == x^3 + x - 9)155Elliptic Curve defined by y^2 + y = x^3 + x - 9 over Rational Field156157sage: R.<x,y> = GF(5)[]158sage: EllipticCurve(x^3 + x^2 + 2 - y^2 - y*x)159Elliptic Curve defined by y^2 + x*y = x^3 + x^2 + 2 over Finite Field of size 5160161We can explicitly specify the `j`-invariant::162163sage: E = EllipticCurve(j=1728); E; E.j_invariant(); E.label()164Elliptic Curve defined by y^2 = x^3 - x over Rational Field1651728166'32a2'167168sage: E = EllipticCurve(j=GF(5)(2)); E; E.j_invariant()169Elliptic Curve defined by y^2 = x^3 + x + 1 over Finite Field of size 51702171172See :trac:`6657` ::173174sage: EllipticCurve(GF(144169),j=1728)175Elliptic Curve defined by y^2 = x^3 + x over Finite Field of size 144169176177By default, when a rational value of `j` is given, the constructed178curve is a minimal twist (minimal conductor for curves with that179`j`-invariant). This can be changed by setting the optional180parameter ``minimal_twist``, which is True by default, to False::181182183sage: EllipticCurve(j=100)184Elliptic Curve defined by y^2 = x^3 + x^2 + 3392*x + 307888 over Rational Field185sage: E =EllipticCurve(j=100); E186Elliptic Curve defined by y^2 = x^3 + x^2 + 3392*x + 307888 over Rational Field187sage: E.conductor()18833129800189sage: E.j_invariant()190100191sage: E =EllipticCurve(j=100, minimal_twist=False); E192Elliptic Curve defined by y^2 = x^3 + 488400*x - 530076800 over Rational Field193sage: E.conductor()194298168200195sage: E.j_invariant()196100197198Without this option, constructing the curve could take a long time199since both `j` and `j-1728` have to be factored to compute the200minimal twist (see :trac:`13100`)::201202sage: E = EllipticCurve_from_j(2^256+1,minimal_twist=False)203sage: E.j_invariant() == 2^256+1204True205206TESTS::207208sage: R = ZZ['u', 'v']209sage: EllipticCurve(R, [1,1])210Elliptic Curve defined by y^2 = x^3 + x + 1 over Multivariate Polynomial Ring in u, v211over Integer Ring212213We create a curve and a point over QQbar (see #6879)::214215sage: E = EllipticCurve(QQbar,[0,1])216sage: E(0)217(0 : 1 : 0)218sage: E.base_field()219Algebraic Field220221sage: E = EllipticCurve(RR,[1,2]); E; E.base_field()222Elliptic Curve defined by y^2 = x^3 + 1.00000000000000*x + 2.00000000000000 over Real Field with 53 bits of precision223Real Field with 53 bits of precision224sage: EllipticCurve(CC,[3,4]); E; E.base_field()225Elliptic Curve defined by y^2 = x^3 + 3.00000000000000*x + 4.00000000000000 over Complex Field with 53 bits of precision226Elliptic Curve defined by y^2 = x^3 + 1.00000000000000*x + 2.00000000000000 over Real Field with 53 bits of precision227Real Field with 53 bits of precision228sage: E = EllipticCurve(QQbar,[5,6]); E; E.base_field()229Elliptic Curve defined by y^2 = x^3 + 5*x + 6 over Algebraic Field230Algebraic Field231232See :trac:`6657` ::233234sage: EllipticCurve(3,j=1728)235Traceback (most recent call last):236...237ValueError: First parameter (if present) must be a ring when j is specified238239sage: EllipticCurve(GF(5),j=3/5)240Traceback (most recent call last):241...242ValueError: First parameter must be a ring containing 3/5243244If the universe of the coefficients is a general field, the object245constructed has type EllipticCurve_field. Otherwise it is246EllipticCurve_generic. See :trac:`9816` ::247248sage: E = EllipticCurve([QQbar(1),3]); E249Elliptic Curve defined by y^2 = x^3 + x + 3 over Algebraic Field250sage: type(E)251<class 'sage.schemes.elliptic_curves.ell_field.EllipticCurve_field_with_category'>252253sage: E = EllipticCurve([RR(1),3]); E254Elliptic Curve defined by y^2 = x^3 + 1.00000000000000*x + 3.00000000000000 over Real Field with 53 bits of precision255sage: type(E)256<class 'sage.schemes.elliptic_curves.ell_field.EllipticCurve_field_with_category'>257258sage: E = EllipticCurve([i,i]); E259Elliptic Curve defined by y^2 = x^3 + I*x + I over Symbolic Ring260sage: type(E)261<class 'sage.schemes.elliptic_curves.ell_field.EllipticCurve_field_with_category'>262sage: E.category()263Category of schemes over Symbolic Ring264sage: SR in Fields()265True266267sage: F = FractionField(PolynomialRing(QQ,'t'))268sage: t = F.gen()269sage: E = EllipticCurve([t,0]); E270Elliptic Curve defined by y^2 = x^3 + t*x over Fraction Field of Univariate Polynomial Ring in t over Rational Field271sage: type(E)272<class 'sage.schemes.elliptic_curves.ell_field.EllipticCurve_field_with_category'>273sage: E.category()274Category of schemes over Fraction Field of Univariate Polynomial Ring in t over Rational Field275276See :trac:`12517`::277278sage: E = EllipticCurve([1..5])279sage: EllipticCurve(E.a_invariants())280Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Rational Field281282See :trac:`11773`::283284sage: E = EllipticCurve()285Traceback (most recent call last):286...287TypeError: invalid input to EllipticCurve constructor288289"""290import ell_generic, ell_field, ell_finite_field, ell_number_field, ell_rational_field, ell_padic_field # here to avoid circular includes291292if j is not None:293if not x is None:294if is_Ring(x):295try:296j = x(j)297except (ZeroDivisionError, ValueError, TypeError):298raise ValueError, "First parameter must be a ring containing %s"%j299else:300raise ValueError, "First parameter (if present) must be a ring when j is specified"301return EllipticCurve_from_j(j, minimal_twist)302303if x is None:304raise TypeError, "invalid input to EllipticCurve constructor"305306if is_SymbolicEquation(x):307x = x.lhs() - x.rhs()308309if parent(x) is SR:310x = x._polynomial_(rings.QQ['x', 'y'])311312if is_MPolynomial(x):313if y is None:314return EllipticCurve_from_Weierstrass_polynomial(x)315else:316return EllipticCurve_from_cubic(x, y, morphism=False)317318if is_Ring(x):319if is_RationalField(x):320return ell_rational_field.EllipticCurve_rational_field(x, y)321elif is_FiniteField(x) or (is_IntegerModRing(x) and x.characteristic().is_prime()):322return ell_finite_field.EllipticCurve_finite_field(x, y)323elif rings.is_pAdicField(x):324return ell_padic_field.EllipticCurve_padic_field(x, y)325elif is_NumberField(x):326return ell_number_field.EllipticCurve_number_field(x, y)327elif x in _Fields:328return ell_field.EllipticCurve_field(x, y)329return ell_generic.EllipticCurve_generic(x, y)330331if isinstance(x, unicode):332x = str(x)333334if isinstance(x, basestring):335return ell_rational_field.EllipticCurve_rational_field(x)336337if is_RingElement(x) and y is None:338raise TypeError, "invalid input to EllipticCurve constructor"339340if not isinstance(x, (list, tuple)):341raise TypeError, "invalid input to EllipticCurve constructor"342343x = Sequence(x)344if not (len(x) in [2,5]):345raise ValueError, "sequence of coefficients must have length 2 or 5"346R = x.universe()347348if isinstance(x[0], (rings.Rational, rings.Integer, int, long)):349return ell_rational_field.EllipticCurve_rational_field(x, y)350351elif is_NumberField(R):352return ell_number_field.EllipticCurve_number_field(x, y)353354elif rings.is_pAdicField(R):355return ell_padic_field.EllipticCurve_padic_field(x, y)356357elif is_FiniteField(R) or (is_IntegerModRing(R) and R.characteristic().is_prime()):358return ell_finite_field.EllipticCurve_finite_field(x, y)359360elif R in _Fields:361return ell_field.EllipticCurve_field(x, y)362363return ell_generic.EllipticCurve_generic(x, y)364365366367def EllipticCurve_from_Weierstrass_polynomial(f):368"""369Return the elliptic curve defined by a cubic in (long) Weierstrass370form.371372INPUT:373374- ``f`` -- a inhomogeneous cubic polynomial in long Weierstrass375form.376377OUTPUT:378379The elliptic curve defined by it.380381EXAMPLES::382383sage: R.<x,y> = QQ[]384sage: f = y^2 + 1*x*y + 3*y - (x^3 + 2*x^2 + 4*x + 6)385sage: EllipticCurve(f)386Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 6 over Rational Field387sage: EllipticCurve(f).a_invariants()388(1, 2, 3, 4, 6)389390The polynomial ring may have extra variables as long as they391do not occur in the polynomial itself::392393sage: R.<x,y,z,w> = QQ[]394sage: EllipticCurve(-y^2 + x^3 + 1)395Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field396sage: EllipticCurve(-x^2 + y^3 + 1)397Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field398sage: EllipticCurve(-w^2 + z^3 + 1)399Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field400401TESTS::402403sage: from sage.schemes.elliptic_curves.constructor import EllipticCurve_from_Weierstrass_polynomial404sage: EllipticCurve_from_Weierstrass_polynomial(-w^2 + z^3 + 1)405Elliptic Curve defined by y^2 = x^3 + 1 over Rational Field406"""407R = f.parent()408cubic_variables = [ x for x in R.gens() if f.degree(x) == 3 ]409quadratic_variables = [ y for y in R.gens() if f.degree(y) == 2 ]410try:411x = cubic_variables[0]412y = quadratic_variables[0]413except IndexError:414raise ValueError('polynomial is not in long Weierstrass form')415416a1 = a2 = a3 = a4 = a6 = 0417x3 = y2 = None418for coeff, mon in f:419if mon == x**3:420x3 = coeff421elif mon == x**2:422a2 = coeff423elif mon == x:424a4 = coeff425elif mon == 1:426a6 = coeff427elif mon == y**2:428y2 = -coeff429elif mon == x*y:430a1 = -coeff431elif mon == y:432a3 = -coeff433else:434raise ValueError('polynomial is not in long Weierstrass form')435436if x3 != y2:437raise ValueError('the coefficient of x^3 and -y^2 must be the same')438elif x3 != 1:439a1, a2, a3, a4, a6 = a1/x3, a2/x3, a3/x3, a4/x3, a6/x3440return EllipticCurve([a1, a2, a3, a4, a6])441442443def EllipticCurve_from_c4c6(c4, c6):444"""445Return an elliptic curve with given `c_4` and446`c_6` invariants.447448EXAMPLES::449450sage: E = EllipticCurve_from_c4c6(17, -2005)451sage: E452Elliptic Curve defined by y^2 = x^3 - 17/48*x + 2005/864 over Rational Field453sage: E.c_invariants()454(17, -2005)455"""456try:457K = c4.parent()458except AttributeError:459K = rings.RationalField()460if K not in _Fields:461K = K.fraction_field()462return EllipticCurve([-K(c4)/K(48), -K(c6)/K(864)])463464465def EllipticCurve_from_j(j, minimal_twist=True):466"""467Return an elliptic curve with given `j`-invariant.468469INPUT:470471- ``j`` -- an element of some field.472473- ``minimal_twist`` (boolean, default True) -- If True and ``j`` is in `\QQ`, the curve returned is a474minimal twist, i.e. has minimal conductor. If `j` is not in `\QQ` this parameter is ignored.475476OUTPUT:477478An elliptic curve with `j`-invariant `j`.479480EXAMPLES::481482sage: E = EllipticCurve_from_j(0); E; E.j_invariant(); E.label()483Elliptic Curve defined by y^2 + y = x^3 over Rational Field4840485'27a3'486487sage: E = EllipticCurve_from_j(1728); E; E.j_invariant(); E.label()488Elliptic Curve defined by y^2 = x^3 - x over Rational Field4891728490'32a2'491492sage: E = EllipticCurve_from_j(1); E; E.j_invariant()493Elliptic Curve defined by y^2 + x*y = x^3 + 36*x + 3455 over Rational Field4941495496The ``minimal_twist`` parameter (ignored except over `\QQ` and497True by default) controls whether or not a minimal twist is498computed::499500sage: EllipticCurve_from_j(100)501Elliptic Curve defined by y^2 = x^3 + x^2 + 3392*x + 307888 over Rational Field502sage: _.conductor()50333129800504sage: EllipticCurve_from_j(100, minimal_twist=False)505Elliptic Curve defined by y^2 = x^3 + 488400*x - 530076800 over Rational Field506sage: _.conductor()507298168200508509Since computing the minimal twist requires factoring both `j` and510`j-1728` the following example would take a long time without511setting ``minimal_twist`` to False::512513sage: E = EllipticCurve_from_j(2^256+1,minimal_twist=False)514sage: E.j_invariant() == 2^256+1515True516517"""518try:519K = j.parent()520except AttributeError:521K = rings.RationalField()522if K not in _Fields:523K = K.fraction_field()524525char=K.characteristic()526if char==2:527if j == 0:528return EllipticCurve(K, [ 0, 0, 1, 0, 0 ])529else:530return EllipticCurve(K, [ 1, 0, 0, 0, 1/j ])531if char == 3:532if j==0:533return EllipticCurve(K, [ 0, 0, 0, 1, 0 ])534else:535return EllipticCurve(K, [ 0, j, 0, 0, -j**2 ])536537if K is rings.RationalField():538# we construct the minimal twist, i.e. the curve with minimal539# conductor with this j_invariant:540if j == 0:541return EllipticCurve(K, [ 0, 0, 1, 0, 0 ]) # 27a3542if j == 1728:543return EllipticCurve(K, [ 0, 0, 0, -1, 0 ]) # 32a2544545if not minimal_twist:546k=j-1728547return EllipticCurve(K, [0,0,0,-3*j*k, -2*j*k**2])548549n = j.numerator()550m = n-1728*j.denominator()551a4 = -3*n*m552a6 = -2*n*m**2553554# Now E=[0,0,0,a4,a6] has j-invariant j=n/d555from sage.sets.set import Set556for p in Set(n.prime_divisors()+m.prime_divisors()):557e = min(a4.valuation(p)//2,a6.valuation(p)//3)558if e>0:559p = p**e560a4 /= p**2561a6 /= p**3562563# Now E=[0,0,0,a4,a6] is minimal at all p != 2,3564tw = [-1,2,-2,3,-3,6,-6]565E1 = EllipticCurve([0,0,0,a4,a6])566Elist = [E1] + [E1.quadratic_twist(t) for t in tw]567crv_cmp = lambda E,F: cmp(E.conductor(),F.conductor())568Elist.sort(cmp=crv_cmp)569return Elist[0]570571# defaults for all other fields:572if j == 0:573return EllipticCurve(K, [ 0, 0, 0, 0, 1 ])574if j == 1728:575return EllipticCurve(K, [ 0, 0, 0, 1, 0 ])576k=j-1728577return EllipticCurve(K, [0,0,0,-3*j*k, -2*j*k**2])578579580def EllipticCurve_from_cubic(F, P, morphism=True):581r"""582Construct an elliptic curve from a ternary cubic with a rational point.583584If you just want the Weierstrass form and are not interested in585the morphism then it is easier to use586:func:`~sage.schemes.elliptic_curves.jacobian.Jacobian`587instead. This will construct the same elliptic curve but you don't588have to supply the point ``P``.589590INPUT:591592- ``F`` -- a homogeneous cubic in three variables with rational593coefficients, as a polynomial ring element, defining a smooth594plane cubic curve.595596- ``P`` -- a 3-tuple `(x,y,z)` defining a projective point on the597curve `F=0`. Need not be a flex, but see caveat on output.598599- ``morphism`` -- boolean (default: ``True``). Whether to return600the morphism or just the elliptic curve.601602OUTPUT:603604An elliptic curve in long Weierstrass form isomorphic to the curve605`F=0`.606607If ``morphism=True`` is passed, then a birational equivalence608between F and the Weierstrass curve is returned. If the point609happens to be a flex, then this is an isomorphism.610611EXAMPLES:612613First we find that the Fermat cubic is isomorphic to the curve614with Cremona label 27a1::615616sage: R.<x,y,z> = QQ[]617sage: cubic = x^3+y^3+z^3618sage: P = [1,-1,0]619sage: E = EllipticCurve_from_cubic(cubic, P, morphism=False); E620Elliptic Curve defined by y^2 + 2*x*y + 1/3*y = x^3 - x^2 - 1/3*x - 1/27 over Rational Field621sage: E.cremona_label()622'27a1'623sage: EllipticCurve_from_cubic(cubic, [0,1,-1], morphism=False).cremona_label()624'27a1'625sage: EllipticCurve_from_cubic(cubic, [1,0,-1], morphism=False).cremona_label()626'27a1'627628Next we find the minimal model and conductor of the Jacobian of the629Selmer curve::630631sage: R.<a,b,c> = QQ[]632sage: cubic = a^3+b^3+60*c^3633sage: P = [1,-1,0]634sage: E = EllipticCurve_from_cubic(cubic, P, morphism=False); E635Elliptic Curve defined by y^2 + 2*x*y + 20*y = x^3 - x^2 - 20*x - 400/3 over Rational Field636sage: E.minimal_model()637Elliptic Curve defined by y^2 = x^3 - 24300 over Rational Field638sage: E.conductor()63924300640641We can also get the birational equivalence to and from the642Weierstrass form. We start with an example where ``P`` is a flex643and the equivalence is an isomorphism::644645sage: f = EllipticCurve_from_cubic(cubic, P, morphism=True)646sage: f647Scheme morphism:648From: Closed subscheme of Projective Space of dimension 2 over Rational Field defined by:649a^3 + b^3 + 60*c^3650To: Elliptic Curve defined by y^2 + 2*x*y + 20*y = x^3 - x^2 - 20*x - 400/3651over Rational Field652Defn: Defined on coordinates by sending (a : b : c) to653(-c : -b + c : 1/20*a + 1/20*b)654655sage: finv = f.inverse(); finv656Scheme morphism:657From: Elliptic Curve defined by y^2 + 2*x*y + 20*y = x^3 - x^2 - 20*x - 400/3658over Rational Field659To: Closed subscheme of Projective Space of dimension 2 over Rational Field defined by:660a^3 + b^3 + 60*c^3661Defn: Defined on coordinates by sending (x : y : z) to662(x + y + 20*z : -x - y : -x)663664We verify that `f` maps the chosen point `P=(1,-1,0)` on the cubic665to the origin of the elliptic curve::666667sage: f([1,-1,0])668(0 : 1 : 0)669sage: finv([0,1,0])670(-1 : 1 : 0)671672To verify the output, we plug in the polynomials to check that673this indeed transforms the cubic into Weierstrass form::674675sage: cubic(finv.defining_polynomials()) * finv.post_rescaling()676-x^3 + x^2*z + 2*x*y*z + y^2*z + 20*x*z^2 + 20*y*z^2 + 400/3*z^3677678sage: E.defining_polynomial()(f.defining_polynomials()) * f.post_rescaling()679a^3 + b^3 + 60*c^3680681If the point is not a flex then the cubic can not be transformed682to a Weierstrass equation by a linear transformation. The general683birational transformation is quadratic::684685sage: cubic = a^3+7*b^3+64*c^3686sage: P = [2,2,-1]687sage: f = EllipticCurve_from_cubic(cubic, P, morphism=True)688sage: E = f.codomain(); E689Elliptic Curve defined by y^2 - 722*x*y - 21870000*y = x^3690+ 23579*x^2 over Rational Field691sage: E.minimal_model()692Elliptic Curve defined by y^2 + y = x^3 - 331 over Rational Field693694sage: f695Scheme morphism:696From: Closed subscheme of Projective Space of dimension 2 over Rational Field defined by:697a^3 + 7*b^3 + 64*c^3698To: Elliptic Curve defined by y^2 - 722*x*y - 21870000*y =699x^3 + 23579*x^2 over Rational Field700Defn: Defined on coordinates by sending (a : b : c) to701(-5/112896*a^2 - 17/40320*a*b - 1/1280*b^2 - 29/35280*a*c702- 13/5040*b*c - 4/2205*c^2 :703-4055/112896*a^2 - 4787/40320*a*b - 91/1280*b^2 - 7769/35280*a*c704- 1993/5040*b*c - 724/2205*c^2 :7051/4572288000*a^2 + 1/326592000*a*b + 1/93312000*b^2 + 1/142884000*a*c706+ 1/20412000*b*c + 1/17860500*c^2)707708sage: finv = f.inverse(); finv709Scheme morphism:710From: Elliptic Curve defined by y^2 - 722*x*y - 21870000*y =711x^3 + 23579*x^2 over Rational Field712To: Closed subscheme of Projective Space of dimension 2 over Rational Field defined by:713a^3 + 7*b^3 + 64*c^3714Defn: Defined on coordinates by sending (x : y : z) to715(2*x^2 + 227700*x*z - 900*y*z :7162*x^2 - 32940*x*z + 540*y*z :717-x^2 - 56520*x*z - 180*y*z)718719sage: cubic(finv.defining_polynomials()) * finv.post_rescaling()720-x^3 - 23579*x^2*z - 722*x*y*z + y^2*z - 21870000*y*z^2721722sage: E.defining_polynomial()(f.defining_polynomials()) * f.post_rescaling()723a^3 + 7*b^3 + 64*c^3724725TESTS::726727sage: R.<x,y,z> = QQ[]728sage: cubic = x^2*y + 4*x*y^2 + x^2*z + 8*x*y*z + 4*y^2*z + 9*x*z^2 + 9*y*z^2729sage: EllipticCurve_from_cubic(cubic, [1,-1,1], morphism=False)730Elliptic Curve defined by y^2 - 882*x*y - 2560000*y = x^3 - 127281*x^2 over Rational Field731"""732import sage.matrix.all as matrix733734# check the input735R = F.parent()736if not is_MPolynomialRing(R):737raise TypeError('equation must be a polynomial')738if R.ngens() != 3:739raise TypeError('equation must be a polynomial in three variables')740if not F.is_homogeneous():741raise TypeError('equation must be a homogeneous polynomial')742K = F.parent().base_ring()743try:744P = [K(c) for c in P]745except TypeError:746raise TypeError('cannot convert %s into %s'%(P,K))747if F(P) != 0:748raise ValueError('%s is not a point on %s'%(P,F))749if len(P) != 3:750raise TypeError('%s is not a projective point'%P)751x, y, z = R.gens()752753# First case: if P = P2 then P is a flex754P2 = chord_and_tangent(F, P)755if are_projectively_equivalent(P, P2, base_ring=K):756# find the tangent to F in P757dx = K(F.derivative(x)(P))758dy = K(F.derivative(y)(P))759dz = K(F.derivative(z)(P))760# find a second point Q on the tangent line but not on the cubic761for tangent in [[dy, -dx, K.zero()], [dz, K.zero(), -dx], [K.zero(), -dz, dx]]:762tangent = projective_point(tangent)763Q = [tangent[0]+P[0], tangent[1]+P[1], tangent[2]+P[2]]764F_Q = F(Q)765if F_Q != 0: # At most one further point may accidentally be on the cubic766break767assert F_Q != 0768# pick linearly independent third point769for third_point in [(1,0,0), (0,1,0), (0,0,1)]:770M = matrix.matrix(K, [Q, P, third_point]).transpose()771if M.is_invertible():772break773F2 = R(M.act_on_polynomial(F))774# scale and dehomogenise775a = K(F2.coefficient(x**3))776F3 = F2/a777b = K(F3.coefficient(y*y*z))778S = rings.PolynomialRing(K, 'x,y,z')779# elliptic curve coordinates780X, Y, Z = S.gen(0), S.gen(1), S(-1/b)*S.gen(2)781F4 = F3(X, Y, Z)782E = EllipticCurve(F4.subs(z=1))783if not morphism:784return E785inv_defining_poly = [ M[i,0]*X + M[i,1]*Y + M[i,2]*Z for i in range(3) ]786inv_post = -1/a787M = M.inverse()788trans_x, trans_y, trans_z = [ M[i,0]*x + M[i,1]*y + M[i,2]*z for i in range(3) ]789fwd_defining_poly = [trans_x, trans_y, -b*trans_z]790fwd_post = -a791792# Second case: P is not a flex, then P, P2, P3 are different793else:794P3 = chord_and_tangent(F, P2)795# send P, P2, P3 to (1:0:0), (0:1:0), (0:0:1) respectively796M = matrix.matrix(K, [P, P2, P3]).transpose()797F2 = M.act_on_polynomial(F)798# substitute x = U^2, y = V*W, z = U*W, and rename (x,y,z)=(U,V,W)799F3 = F2.substitute({x:x**2, y:y*z, z:x*z}) // (x**2*z)800# scale and dehomogenise801a = K(F3.coefficient(x**3))802F4 = F3/a803b = K(F4.coefficient(y*y*z))804# change to a polynomial in only two variables805S = rings.PolynomialRing(K, 'x,y,z')806# elliptic curve coordinates807X, Y, Z = S.gen(0), S.gen(1), S(-1/b)*S.gen(2)808F5 = F4(X, Y, Z)809E = EllipticCurve(F5.subs(z=1))810if not morphism:811return E812inv_defining_poly = [ M[i,0]*X*X + M[i,1]*Y*Z + M[i,2]*X*Z for i in range(3) ]813inv_post = -1/a/(X**2)/Z814M = M.inverse()815trans_x, trans_y, trans_z = [816(M[i,0]*x + M[i,1]*y + M[i,2]*z) for i in range(3) ]817fwd_defining_poly = [ trans_x*trans_z, trans_x*trans_y, -b*trans_z*trans_z ]818fwd_post = -a/(trans_x*trans_z*trans_z)819820# Construct the morphism821from sage.schemes.projective.projective_space import ProjectiveSpace822P2 = ProjectiveSpace(2, K, names=map(str, R.gens()))823cubic = P2.subscheme(F)824from sage.schemes.elliptic_curves.weierstrass_transform import \825WeierstrassTransformationWithInverse826return WeierstrassTransformationWithInverse(827cubic, E, fwd_defining_poly, fwd_post, inv_defining_poly, inv_post)828829830def chord_and_tangent(F, P):831"""832Use the chord and tangent method to get another point on a cubic.833834INPUT:835836- ``F`` -- a homogeneous cubic in three variables with rational837coefficients, as a polynomial ring element, defining a smooth838plane cubic curve.839840- ``P`` -- a 3-tuple `(x,y,z)` defining a projective point on the841curve `F=0`.842843OUTPUT:844845Another point satisfying the equation ``F``.846847EXAMPLES::848849sage: R.<x,y,z> = QQ[]850sage: from sage.schemes.elliptic_curves.constructor import chord_and_tangent851sage: F = x^3+y^3+60*z^3852sage: chord_and_tangent(F, [1,-1,0])853[1, -1, 0]854855sage: F = x^3+7*y^3+64*z^3856sage: p0 = [2,2,-1]857sage: p1 = chord_and_tangent(F, p0); p1858[-5, 3, -1]859sage: p2 = chord_and_tangent(F, p1); p2860[1265, -183, -314]861862TESTS::863864sage: F(p2)8650866sage: map(type, p2)867[<type 'sage.rings.rational.Rational'>,868<type 'sage.rings.rational.Rational'>,869<type 'sage.rings.rational.Rational'>]870"""871# check the input872R = F.parent()873if not is_MPolynomialRing(R):874raise TypeError('equation must be a polynomial')875if R.ngens() != 3:876raise TypeError('%s is not a polynomial in three variables'%F)877if not F.is_homogeneous():878raise TypeError('%s is not a homogeneous polynomial'%F)879x, y, z = R.gens()880if len(P) != 3:881raise TypeError('%s is not a projective point'%P)882K = R.base_ring()883try:884P = [K(c) for c in P]885except TypeError:886raise TypeError('cannot coerce %s into %s'%(P,K))887if F(P) != 0:888raise ValueError('%s is not a point on %s'%(P,F))889890# find the tangent to F in P891dx = K(F.derivative(x)(P))892dy = K(F.derivative(y)(P))893dz = K(F.derivative(z)(P))894# if dF/dy(P) = 0, change variables so that dF/dy != 0895if dy == 0:896if dx != 0:897g = F.substitute({x:y, y:x})898Q = [P[1], P[0], P[2]]899R = chord_and_tangent(g, Q)900return [R[1], R[0], R[2]]901elif dz != 0:902g = F.substitute({y:z, z:y})903Q = [P[0], P[2], P[1]]904R = chord_and_tangent(g, Q)905return [R[0], R[2], R[1]]906else:907raise ValueError('%s is singular at %s'%(F, P))908909# t will be our choice of parmeter of the tangent plane910# dx*(x-P[0]) + dy*(y-P[1]) + dz*(z-P[2])911# through the point P912t = rings.PolynomialRing(K, 't').gen(0)913Ft = F(dy*t+P[0], -dx*t+P[1], P[2])914if Ft == 0: # (dy, -dx, 0) is projectively equivalent to P915# then (0, -dz, dy) is not projectively equivalent to P916g = F.substitute({x:z, z:x})917Q = [P[2], P[1], P[0]]918R = chord_and_tangent(g, Q)919return [R[2], R[1], R[0]]920# Ft has a double zero at t=0 by construction, which we now remove921Ft = Ft // t**2922923# first case: the third point is at t=infinity924if Ft.is_constant():925return projective_point([dy, -dx, 0])926# second case: the third point is at finite t927else:928assert Ft.degree() == 1929t0 = Ft.roots()[0][0]930return projective_point([dy*t0+P[0], -dx*t0+P[1], P[2]])931932933def projective_point(p):934"""935Return equivalent point with denominators removed936937INPUT:938939- ``P``, ``Q`` -- list/tuple of projective coordinates.940941OUTPUT:942943List of projective coordinates.944945EXAMPLES::946947sage: from sage.schemes.elliptic_curves.constructor import projective_point948sage: projective_point([4/5, 6/5, 8/5])949[2, 3, 4]950sage: F = GF(11)951sage: projective_point([F(4), F(8), F(2)])952[4, 8, 2]953"""954try:955p_gcd = rings.integer.GCD_list([x.numerator() for x in p])956p_lcm = rings.integer.LCM_list([x.denominator() for x in p])957except AttributeError:958return p959scale = p_lcm / p_gcd960return [scale * x for x in p]961962963def are_projectively_equivalent(P, Q, base_ring):964"""965Test whether ``P`` and ``Q`` are projectively equivalent.966967INPUT:968969- ``P``, ``Q`` -- list/tuple of projective coordinates.970971- ``base_ring`` -- the base ring.972973OUTPUT:974975Boolean.976977EXAMPLES::978979sage: from sage.schemes.elliptic_curves.constructor import are_projectively_equivalent980sage: are_projectively_equivalent([0,1,2,3], [0,1,2,2], base_ring=QQ)981False982sage: are_projectively_equivalent([0,1,2,3], [0,2,4,6], base_ring=QQ)983True984"""985from sage.matrix.constructor import matrix986return matrix(base_ring, [P, Q]).rank() < 2987988989def EllipticCurve_from_plane_curve(C, P):990"""991Deprecated way to construct an elliptic curve.992993Use :meth:`~sage.schemes.elliptic_curves.jacobian.Jacobian` instead.994995EXAMPLES::996997sage: R.<x,y,z> = QQ[]998sage: C = Curve(x^3+y^3+z^3)999sage: P = C(1,-1,0)1000sage: E = EllipticCurve_from_plane_curve(C,P); E # long time (3s on sage.math, 2013)1001doctest:...: DeprecationWarning: use Jacobian(C) instead1002See http://trac.sagemath.org/3416 for details.1003Elliptic Curve defined by y^2 = x^3 - 27/4 over Rational Field1004"""1005from sage.misc.superseded import deprecation1006deprecation(3416, 'use Jacobian(C) instead')1007# Note: this function never used the rational point1008from sage.schemes.elliptic_curves.jacobian import Jacobian1009return Jacobian(C)101010111012def EllipticCurves_with_good_reduction_outside_S(S=[], proof=None, verbose=False):1013r"""1014Returns a sorted list of all elliptic curves defined over `Q`1015with good reduction outside the set `S` of primes.10161017INPUT:10181019- ``S`` - list of primes (default: empty list).10201021- ``proof`` - True/False (default True): the MW basis for1022auxiliary curves will be computed with this proof flag.10231024- ``verbose`` - True/False (default False): if True, some details1025of the computation will be output.10261027.. note::10281029Proof flag: The algorithm used requires determining all1030S-integral points on several auxiliary curves, which in turn1031requires the computation of their generators. This is not1032always possible (even in theory) using current knowledge.10331034The value of this flag is passed to the function which1035computes generators of various auxiliary elliptic curves, in1036order to find their S-integral points. Set to False if the1037default (True) causes warning messages, but note that you can1038then not rely on the set of curves returned being1039complete.10401041EXAMPLES::10421043sage: EllipticCurves_with_good_reduction_outside_S([])1044[]1045sage: elist = EllipticCurves_with_good_reduction_outside_S([2])1046sage: elist1047[Elliptic Curve defined by y^2 = x^3 + 4*x over Rational Field,1048Elliptic Curve defined by y^2 = x^3 - x over Rational Field,1049...1050Elliptic Curve defined by y^2 = x^3 - x^2 - 13*x + 21 over Rational Field]1051sage: len(elist)1052241053sage: ', '.join([e.label() for e in elist])1054'32a1, 32a2, 32a3, 32a4, 64a1, 64a2, 64a3, 64a4, 128a1, 128a2, 128b1, 128b2, 128c1, 128c2, 128d1, 128d2, 256a1, 256a2, 256b1, 256b2, 256c1, 256c2, 256d1, 256d2'10551056Without ``Proof=False``, this example gives two warnings::10571058sage: elist = EllipticCurves_with_good_reduction_outside_S([11],proof=False) # long time (14s on sage.math, 2011)1059sage: len(elist) # long time1060121061sage: ', '.join([e.label() for e in elist]) # long time1062'11a1, 11a2, 11a3, 121a1, 121a2, 121b1, 121b2, 121c1, 121c2, 121d1, 121d2, 121d3'10631064sage: elist = EllipticCurves_with_good_reduction_outside_S([2,3]) # long time (26s on sage.math, 2011)1065sage: len(elist) # long time10667521067sage: max([e.conductor() for e in elist]) # long time1068622081069sage: [N.factor() for N in Set([e.conductor() for e in elist])] # long time1070[2^7,10712^8,10722^3 * 3^4,10732^2 * 3^3,10742^8 * 3^4,10752^4 * 3^4,10762^3 * 3,10772^7 * 3,10782^3 * 3^5,10793^3,10802^8 * 3,10812^5 * 3^4,10822^4 * 3,10832 * 3^4,10842^2 * 3^2,10852^6 * 3^4,10862^6,10872^7 * 3^2,10882^4 * 3^5,10892^4 * 3^3,10902 * 3^3,10912^6 * 3^3,10922^6 * 3,10932^5,10942^2 * 3^4,10952^3 * 3^2,10962^5 * 3,10972^7 * 3^4,10982^2 * 3^5,10992^8 * 3^2,11002^5 * 3^2,11012^7 * 3^5,11022^8 * 3^5,11032^3 * 3^3,11042^8 * 3^3,11052^5 * 3^5,11062^4 * 3^2,11072 * 3^5,11082^5 * 3^3,11092^6 * 3^5,11102^7 * 3^3,11113^5,11122^6 * 3^2]1113"""1114from ell_egros import (egros_from_jlist, egros_get_j)1115return egros_from_jlist(egros_get_j(S, proof=proof, verbose=verbose), S)111611171118