Path: blob/master/sage/schemes/elliptic_curves/ell_point.py
4128 views
# -*- coding: utf-8 -*-1r"""2Points on elliptic curves34The base class ``EllipticCurvePoint_field``, derived from5``AdditiveGroupElement``, provides support for points on elliptic6curves defined over general fields. The derived classes7``EllipticCurvePoint_number_field`` and8``EllipticCurvePoint_finite_field`` provide further support for point9on curves defined over number fields (including the rational field10`\QQ`) and over finite fields. Although there is no special11class for points over `\QQ`, there is currently greater12functionality implemented over `\QQ` than over other number13fields.1415The class ``EllipticCurvePoint``, which is based on16``SchemeMorphism_point_projective_ring``, currently has little extra17functionality.1819EXAMPLES:2021An example over `\QQ`::2223sage: E = EllipticCurve('389a1')24sage: P = E(-1,1); P25(-1 : 1 : 1)26sage: Q = E(0,-1); Q27(0 : -1 : 1)28sage: P+Q29(4 : 8 : 1)30sage: P-Q31(1 : 0 : 1)32sage: 3*P-5*Q33(328/361 : -2800/6859 : 1)3435An example over a number field::3637sage: K.<i> = QuadraticField(-1)38sage: E = EllipticCurve(K,[1,0,0,0,-1])39sage: P = E(0,i); P40(0 : i : 1)41sage: P.order()42+Infinity43sage: 101*P-100*P==P44True4546An example over a finite field::4748sage: K.<a> = GF(101^3)49sage: E = EllipticCurve(K,[1,0,0,0,-1])50sage: P = E(40*a^2 + 69*a + 84 , 58*a^2 + 73*a + 45)51sage: P.order()52103221053sage: E.cardinality()5410322105556Arithmetic with a point over an extension of a finite field::5758sage: k.<a> = GF(5^2)59sage: E = EllipticCurve(k,[1,0]); E60Elliptic Curve defined by y^2 = x^3 + x over Finite Field in a of size 5^261sage: P = E([a,2*a+4])62sage: 5*P63(2*a + 3 : 2*a : 1)64sage: P*565(2*a + 3 : 2*a : 1)66sage: P + P + P + P + P67(2*a + 3 : 2*a : 1)6869::7071sage: F = Zmod(3)72sage: E = EllipticCurve(F,[1,0]);73sage: P = E([2,1])74sage: import sys75sage: n = sys.maxint76sage: P*(n+1)-P*n == P77True7879Arithmetic over `\ZZ/N\ZZ` with composite `N` is supported. When an80operation tries to invert a non-invertible element, a81ZeroDivisionError is raised and a factorization of the modulus appears82in the error message::8384sage: N = 171576151385sage: E = EllipticCurve(Integers(N),[3,-13])86sage: P = E(2,1)87sage: LCM([2..60])*P88Traceback (most recent call last):89...90ZeroDivisionError: Inverse of 1520944668 does not exist (characteristic = 1715761513 = 26927*63719)919293AUTHORS:9495- William Stein (2005) -- Initial version9697- Robert Bradshaw et al....9899- John Cremona (Feb 2008) -- Point counting and group structure for100non-prime fields, Frobenius endomorphism and order, elliptic logs101102- John Cremona (Aug 2008) -- Introduced ``EllipticCurvePoint_number_field`` class103104- Tobias Nagel, Michael Mardaus, John Cremona (Dec 2008) -- `p`-adic elliptic logarithm over `\QQ`105106- David Hansen (Jan 2009) -- Added ``weil_pairing`` function to ``EllipticCurvePoint_finite_field`` class107108- Mariah Lenox (March 2011) -- Added ``tate_pairing`` and ``ate_pairing``109functions to ``EllipticCurvePoint_finite_field`` class110111"""112113#*****************************************************************************114# Copyright (C) 2005 William Stein <[email protected]>115#116# Distributed under the terms of the GNU General Public License (GPL)117# as published by the Free Software Foundation; either version 2 of118# the License, or (at your option) any later version.119# http://www.gnu.org/licenses/120#*****************************************************************************121122123import math124125from sage.structure.element import AdditiveGroupElement126import sage.plot.all as plot127128from sage.rings.padics.factory import Qp129from sage.rings.padics.precision_error import PrecisionError130131import sage.rings.all as rings132from sage.groups.all import AbelianGroup133import sage.groups.generic as generic134from sage.libs.pari.all import pari, PariError135from sage.structure.sequence import Sequence136137from sage.schemes.plane_curves.projective_curve import Hasse_bounds138from sage.schemes.generic.morphism import (SchemeMorphism_point_projective_ring,139SchemeMorphism_point_abelian_variety_field,140is_SchemeMorphism)141142from constructor import EllipticCurve143144oo = rings.infinity # infinity145146class EllipticCurvePoint(SchemeMorphism_point_projective_ring):147"""148A point on an elliptic curve.149"""150def __cmp__(self, other):151"""152Standard comparison function for points on elliptic curves, to153allow sorting and equality testing.154155EXAMPLES:156sage: E=EllipticCurve(QQ,[1,1])157sage: P=E(0,1)158sage: P.order()159+Infinity160sage: Q=P+P161sage: P==Q162False163sage: Q+Q == 4*P164True165"""166assert isinstance(other, (int, long, rings.Integer)) and other == 0167if self.is_zero():168return 0169else:170return -1171172173class EllipticCurvePoint_field(SchemeMorphism_point_abelian_variety_field):174"""175A point on an elliptic curve over a field. The point has coordinates176in the base field.177178EXAMPLES::179180sage: E = EllipticCurve('37a')181sage: E([0,0])182(0 : 0 : 1)183sage: E(0,0) # brackets are optional184(0 : 0 : 1)185sage: E([GF(5)(0), 0]) # entries are coerced186(0 : 0 : 1)187188sage: E(0.000, 0)189(0 : 0 : 1)190191sage: E(1,0,0)192Traceback (most recent call last):193...194TypeError: Coordinates [1, 0, 0] do not define a point on195Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field196197::198199sage: E = EllipticCurve([0,0,1,-1,0])200sage: S = E(QQ); S201Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field202203sage: K.<i>=NumberField(x^2+1)204sage: E=EllipticCurve(K,[0,1,0,-160,308])205sage: P=E(26,-120)206sage: Q=E(2+12*i,-36+48*i)207sage: P.order() == Q.order() == 4 # long time (3s)208True209sage: 2*P==2*Q210False211212::213214sage: K.<t>=FractionField(PolynomialRing(QQ,'t'))215sage: E=EllipticCurve([0,0,0,0,t^2])216sage: P=E(0,t)217sage: P,2*P,3*P218((0 : t : 1), (0 : -t : 1), (0 : 1 : 0))219220221TESTS::222223sage: loads(S.dumps()) == S224True225sage: E = EllipticCurve('37a')226sage: P = E(0,0); P227(0 : 0 : 1)228sage: loads(P.dumps()) == P229True230sage: T = 100*P231sage: loads(T.dumps()) == T232True233234Test pickling an elliptic curve that has known points on it::235236sage: e = EllipticCurve([0, 0, 1, -1, 0]); g = e.gens(); loads(dumps(e)) == e237True238"""239def __init__(self, curve, v, check=True):240"""241Constructor for a point on an elliptic curve.242243INPUT:244245- curve -- an elliptic curve246- v -- data determining a point (another point, the integer2470, or a tuple of coordinates)248249EXAMPLE::250251sage: E = EllipticCurve('43a')252sage: P = E([2, -4, 2]); P253(1 : -2 : 1)254sage: P == E([1,-2])255True256sage: P = E(0); P257(0 : 1 : 0)258sage: P=E(2, -4, 2); P259(1 : -2 : 1)260"""261point_homset = curve.point_homset()262if is_SchemeMorphism(v) or isinstance(v,EllipticCurvePoint_field):263v = list(v)264elif v == 0:265# some of the code assumes that E(0) has integral entries266# irregardless of the base ring...267#R = self.base_ring()268#v = (R.zero(),R.one(),R.zero())269v = (0,1,0)270if check:271# mostly from SchemeMorphism_point_projective_field272d = point_homset.codomain().ambient_space().ngens()273if not isinstance(v,(list,tuple)):274raise TypeError, \275"Argument v (= %s) must be a scheme point, list, or tuple."%str(v)276if len(v) != d and len(v) != d-1:277raise TypeError, "v (=%s) must have %s components"%(v, d)278v = Sequence(v, point_homset.value_ring())279if len(v) == d-1: # very common special case280v.append(v.universe()(1))281282n = len(v)283all_zero = True284for i in range(n):285c = v[n-1-i]286if c:287all_zero = False288if c == 1:289break290for j in range(n-i):291v[j] /= c292break293if all_zero:294raise ValueError, "%s does not define a valid point since all entries are 0"%repr(v)295296x, y, z = v297if z == 0:298test = x299else:300a1, a2, a3, a4, a6 = curve.ainvs()301test = y**2 + (a1*x+a3)*y - (((x+a2)*x+a4)*x+a6)302if not test == 0:303raise TypeError, "Coordinates %s do not define a point on %s"%(list(v),curve)304305SchemeMorphism_point_abelian_variety_field.__init__(self, point_homset, v, check=False)306#AdditiveGroupElement.__init__(self, point_homset)307308def _repr_(self):309"""310Return a string representation of this point.311312EXAMPLE::313314sage: E = EllipticCurve('39a')315sage: P = E([-2, 1, 1])316sage: P._repr_()317'(-2 : 1 : 1)'318"""319return self.codomain().ambient_space()._repr_generic_point(self._coords)320321def _latex_(self):322"""323Return a LaTeX representation of this point.324325EXAMPLE::326327sage: E = EllipticCurve('40a')328sage: P = E([3, 0])329sage: P._latex_()330'\\left(3 : 0 : 1\\right)'331"""332return self.codomain().ambient_space()._latex_generic_point(self._coords)333334def __getitem__(self, n):335"""336Return the n'th coordinate of this point.337338EXAMPLE::339340sage: E = EllipticCurve('42a')341sage: P = E([-17, -51, 17])342sage: [P[i] for i in [2,1,0]]343[1, -3, -1]344"""345return self._coords[n]346347def __iter__(self):348"""349Return the coordinates of this point as a list.350351EXAMPLE::352353sage: E = EllipticCurve('37a')354sage: list(E([0,0]))355[0, 0, 1]356"""357return iter(self._coords)358359def __tuple__(self):360"""361Return the coordinates of this point as a tuple.362363EXAMPLE::364365sage: E = EllipticCurve('44a')366sage: P = E([1, -2, 1])367sage: P.__tuple__()368(1, -2, 1)369"""370return tuple(self._coords) # Warning: _coords is a list!371372def __cmp__(self, other):373"""374Comparison function for points to allow sorting and equality testing.375376EXAMPLES::377378sage: E = EllipticCurve('45a')379sage: P = E([2, -1, 1])380sage: P == E(0)381False382sage: P+P == E(0)383True384"""385if not isinstance(other, EllipticCurvePoint_field):386try:387other = self.codomain().ambient_space()(other)388except TypeError:389return -1390return cmp(self._coords, other._coords)391392def _pari_(self):393r"""394Converts this point to PARI format.395396EXAMPLES::397398sage: E = EllipticCurve([0,0,0,3,0])399sage: O = E(0)400sage: P = E.point([1,2])401sage: O._pari_()402[0]403sage: P._pari_()404[1, 2]405406The following implicitly calls O._pari_() and P._pari_()::407408sage: pari(E).elladd(O,P)409[1, 2]410411TESTS::412413Try the same over a finite field::414415sage: E = EllipticCurve(GF(11), [0,0,0,3,0])416sage: O = E(0)417sage: P = E.point([1,2])418sage: O._pari_()419[0]420sage: P._pari_()421[Mod(1, 11), Mod(2, 11)]422sage: pari(E).elladd(O,P)423[Mod(1, 11), Mod(2, 11)]424"""425if self[2]:426return pari([self[0]/self[2], self[1]/self[2]])427else:428return pari([0])429430def scheme(self):431"""432Return the scheme of this point, i.e., the curve it is on.433This is synonymous with curve() which is perhaps more434intuitive.435436EXAMPLES::437438sage: E=EllipticCurve(QQ,[1,1])439sage: P=E(0,1)440sage: P.scheme()441Elliptic Curve defined by y^2 = x^3 + x + 1 over Rational Field442sage: P.scheme() == P.curve()443True444sage: K.<a>=NumberField(x^2-3,'a')445sage: P=E.base_extend(K)(1,a)446sage: P.scheme()447Elliptic Curve defined by y^2 = x^3 + x + 1 over Number Field in a with defining polynomial x^2 - 3448"""449#The following text is just not true: it applies to the class450#EllipticCurvePoint, which appears to be never used, but does451#not apply to EllipticCurvePoint_field which is simply derived452#from AdditiveGroupElement.453#454#"Technically, points on curves in Sage are scheme maps from455# the domain Spec(F) where F is the base field of the curve to456# the codomain which is the curve. See also domain() and457# codomain()."458459return self.codomain()460461def domain(self):462"""463Return the domain of this point, which is `Spec(F)` where `F` is464the field of definition.465466EXAMPLES::467468sage: E=EllipticCurve(QQ,[1,1])469sage: P=E(0,1)470sage: P.domain()471Spectrum of Rational Field472sage: K.<a>=NumberField(x^2-3,'a')473sage: P=E.base_extend(K)(1,a)474sage: P.domain()475Spectrum of Number Field in a with defining polynomial x^2 - 3476"""477return self.parent().domain()478479def codomain(self):480"""481Return the codomain of this point, which is the curve it is482on. Synonymous with curve() which is perhaps more intuitive.483484EXAMPLES::485486sage: E=EllipticCurve(QQ,[1,1])487sage: P=E(0,1)488sage: P.domain()489Spectrum of Rational Field490sage: K.<a>=NumberField(x^2-3,'a')491sage: P=E.base_extend(K)(1,a)492sage: P.codomain()493Elliptic Curve defined by y^2 = x^3 + x + 1 over Number Field in a with defining polynomial x^2 - 3494sage: P.codomain() == P.curve()495True496"""497return self.parent().codomain()498499def order(self):500r"""501Return the order of this point on the elliptic curve.502503If the point is zero, returns 1, otherwise raise a504NotImplementedError.505506For curves over number fields and finite fields, see below.507508.. note::509510:meth:`additive_order` is a synonym for :meth:`order`511512EXAMPLE::513514sage: K.<t>=FractionField(PolynomialRing(QQ,'t'))515sage: E=EllipticCurve([0,0,0,-t^2,0])516sage: P=E(t,0)517sage: P.order()518Traceback (most recent call last):519...520NotImplementedError: Computation of order of a point not implemented over general fields.521sage: E(0).additive_order()5221523sage: E(0).order() == 1524True525526"""527if self.is_zero():528return rings.Integer(1)529raise NotImplementedError, "Computation of order of a point not implemented over general fields."530531additive_order = order532533def curve(self):534"""535Return the curve that this point is on.536537EXAMPLES::538539sage: E = EllipticCurve('389a')540sage: P = E([-1,1])541sage: P.curve()542Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field543"""544return self.scheme()545546def __nonzero__(self):547"""548Return True if this is not the zero point on the curve.549550EXAMPLES::551552sage: E = EllipticCurve('37a')553sage: P = E(0); P554(0 : 1 : 0)555sage: P.is_zero()556True557sage: P = E.gens()[0]558sage: P.is_zero()559False560"""561return bool(self[2])562563def has_finite_order(self):564"""565Return True if this point has finite additive order as an element566of the group of points on this curve.567568For fields other than number fields and finite fields, this is569NotImplemented unless self.is_zero().570571EXAMPLES::572573sage: K.<t>=FractionField(PolynomialRing(QQ,'t'))574sage: E=EllipticCurve([0,0,0,-t^2,0])575sage: P = E(0)576sage: P.has_finite_order()577True578sage: P=E(t,0)579sage: P.has_finite_order()580Traceback (most recent call last):581...582NotImplementedError: Computation of order of a point not implemented over general fields.583sage: (2*P).is_zero()584True585"""586if self.is_zero(): return True587return self.order() != oo588589is_finite_order = has_finite_order # for backward compatibility590591def has_infinite_order(self):592"""593Return True if this point has infinite additive order as an element594of the group of points on this curve.595596For fields other than number fields and finite fields, this is597NotImplemented unless self.is_zero().598599EXAMPLES::600601sage: K.<t>=FractionField(PolynomialRing(QQ,'t'))602sage: E=EllipticCurve([0,0,0,-t^2,0])603sage: P = E(0)604sage: P.has_infinite_order()605False606sage: P=E(t,0)607sage: P.has_infinite_order()608Traceback (most recent call last):609...610NotImplementedError: Computation of order of a point not implemented over general fields.611sage: (2*P).is_zero()612True613"""614if self.is_zero(): return False615return self.order() == oo616617def plot(self, **args):618"""619Plot this point on an elliptic curve.620621INPUT:622623- ``**args`` -- all arguments get passed directly onto the point624plotting function.625626EXAMPLES::627628sage: E = EllipticCurve('389a')629sage: P = E([-1,1])630sage: P.plot(pointsize=30, rgbcolor=(1,0,0))631"""632if self.is_zero():633return plot.text("$\\infty$", (-3,3), **args)634635else:636return plot.point((self[0], self[1]), **args)637638def _add_(self, right):639"""640Add self to right.641642EXAMPLES::643644sage: E = EllipticCurve('389a')645sage: P = E([-1,1]); Q = E([0,0])646sage: P + Q647(1 : 0 : 1)648sage: P._add_(Q) == P + Q649True650651Example to show that bug \#4820 is fixed::652653sage: [type(c) for c in 2*EllipticCurve('37a1').gen(0)]654[<type 'sage.rings.rational.Rational'>,655<type 'sage.rings.rational.Rational'>,656<type 'sage.rings.rational.Rational'>]657"""658# Use Prop 7.1.7 of Cohen "A Course in Computational Algebraic Number Theory"659if self.is_zero():660return right661if right.is_zero():662return self663E = self.curve()664a1, a2, a3, a4, a6 = E.ainvs()665x1, y1 = self[0], self[1]666x2, y2 = right[0], right[1]667if x1 == x2 and y1 == -y2 - a1*x2 - a3:668return E(0) # point at infinity669670if x1 == x2 and y1 == y2:671try:672m = (3*x1*x1 + 2*a2*x1 + a4 - a1*y1) / (2*y1 + a1*x1 + a3)673except ZeroDivisionError:674R = E.base_ring()675if R.is_finite():676N = R.characteristic()677from sage.rings.all import ZZ678N1 = N.gcd(ZZ(2*y1 + a1*x1 + a3))679N2 = N//N1680raise ZeroDivisionError, "Inverse of %s does not exist (characteristic = %s = %s*%s)"%(2*y1 + a1*x1 + a3, N,N1,N2)681else:682raise ZeroDivisionError, "Inverse of %s does not exist"%(2*y1 + a1*x1 + a3)683else:684try:685m = (y1-y2)/(x1-x2)686except ZeroDivisionError:687R = E.base_ring()688if R.is_finite():689N = R.characteristic()690from sage.rings.all import ZZ691N1 = N.gcd(ZZ(x1-x2))692N2 = N//N1693raise ZeroDivisionError, "Inverse of %s does not exist (characteristic = %s = %s*%s)"%(x1-x2, N,N1,N2)694else:695raise ZeroDivisionError, "Inverse of %s does not exist"%(x1-x2)696697x3 = -x1 - x2 - a2 + m*(m+a1)698y3 = -y1 - a3 - a1*x3 + m*(x1-x3)699# See \#4820 for why we need to coerce 1 into the base ring here:700return E.point([x3, y3, E.base_ring()(1)], check=False)701702def _sub_(self, right):703"""704Subtract right from self.705706EXAMPLES::707708sage: E = EllipticCurve('389a')709sage: P = E([-1,1]); Q = E([0,0])710sage: P - Q711(4 : 8 : 1)712sage: P - Q == P._sub_(Q)713True714sage: (P - Q) + Q715(-1 : 1 : 1)716sage: P717(-1 : 1 : 1)718"""719return self + (-right)720721def __neg__(self):722"""723Return the additive inverse of this point.724725EXAMPLES::726727sage: E = EllipticCurve('389a')728sage: P = E([-1,1])729sage: Q = -P; Q730(-1 : -2 : 1)731sage: Q + P732(0 : 1 : 0)733734Example to show that bug \#4820 is fixed::735736sage: [type(c) for c in -EllipticCurve('37a1').gen(0)]737[<type 'sage.rings.rational.Rational'>,738<type 'sage.rings.rational.Rational'>,739<type 'sage.rings.rational.Rational'>]740"""741if self.is_zero():742return self743E, x, y = self.curve(), self[0], self[1]744# See \#4820 for why we need to coerce 1 into the base ring here:745return E.point([x, -y - E.a1()*x - E.a3(), E.base_ring()(1)], check=False)746747def xy(self):748"""749Return the `x` and `y` coordinates of this point, as a 2-tuple.750If this is the point at infinity a ZeroDivisionError is raised.751752EXAMPLES::753754sage: E = EllipticCurve('389a')755sage: P = E([-1,1])756sage: P.xy()757(-1, 1)758sage: Q = E(0); Q759(0 : 1 : 0)760sage: Q.xy()761Traceback (most recent call last):762...763ZeroDivisionError: Rational division by zero764"""765if self[2] == 1:766return self[0], self[1]767else:768return self[0]/self[2], self[1]/self[2]769770def is_divisible_by(self, m):771"""772Return True if there exists a point `Q` defined over the same773field as self such that `mQ` == self.774775INPUT:776777- ``m`` -- a positive integer.778779OUTPUT:780781(bool) -- True if there is a solution, else False.782783.. warning::784785This function usually triggers the computation of the786`m`-th division polynomial of the associated elliptic787curve, which will be expensive if `m` is large, though it788will be cached for subsequent calls with the same `m`.789790EXAMPLES::791792sage: E = EllipticCurve('389a')793sage: Q = 5*E(0,0); Q794(-2739/1444 : -77033/54872 : 1)795sage: Q.is_divisible_by(4)796False797sage: Q.is_divisible_by(5)798True799800A finite field example::801802sage: E = EllipticCurve(GF(101),[23,34])803sage: E.cardinality().factor()8042 * 53805sage: Set([T.order() for T in E.points()])806{1, 106, 2, 53}807sage: len([T for T in E.points() if T.is_divisible_by(2)])80853809sage: len([T for T in E.points() if T.is_divisible_by(3)])810106811812TESTS:813814This shows that the bug reported at #10076 is fixed::815816sage: K = QuadraticField(8,'a')817sage: E = EllipticCurve([K(0),0,0,-1,0])818sage: P = E([-1,0])819sage: P.is_divisible_by(2)820False821sage: P.division_points(2)822[]823824Note that it is not sufficient to test that825``self.division_points(m,poly_only=True)`` has roots::826827sage: P.division_points(2, poly_only=True).roots()828[(1/2*a - 1, 1), (-1/2*a - 1, 1)]829830sage: tor = E.torsion_points(); len(tor)8318832sage: [T.order() for T in tor]833[4, 2, 4, 2, 4, 2, 4, 1]834sage: all([T.is_divisible_by(3) for T in tor])835True836sage: Set([T for T in tor if T.is_divisible_by(2)])837{(0 : 1 : 0), (1 : 0 : 1)}838sage: Set([2*T for T in tor])839{(0 : 1 : 0), (1 : 0 : 1)}840841842"""843# Coerce the input m to an integer844m = rings.Integer(m)845846# Check for trivial cases of m = 1, -1 and 0.847if m == 1 or m == -1:848return True849if m == 0:850return self == 0 # then m*self=self for all m!851m = m.abs()852853# Now the following line would of course be correct, but we854# work harder to be more efficient:855# return len(self.division_points(m)) > 0856857P = self858859# If P has finite order n and gcd(m,n)=1 then the result is860# True. However, over general fields computing P.order() is861# not implemented.862863try:864n = P.order()865if not n == oo:866if m.gcd(n)==1:867return True868except NotImplementedError:869pass870871P_is_2_torsion = (P==-P)872g = P.division_points(m, poly_only=True)873874if not P_is_2_torsion:875# In this case deg(g)=m^2, and each root in K lifts to two876# points Q,-Q both in E(K), of which exactly one is a877# solution. So we just check the existence of roots:878return len(g.roots())>0879880# Now 2*P==0881882if m%2==1:883return True # P itself is a solution when m is odd884885# Now m is even and 2*P=0. Roots of g in K may or may not886# lift to solutions in E(K), so we fall back to the default.887# Note that division polynomials are cached so this is not888# inefficient:889890return len(self.division_points(m)) > 0891892def division_points(self, m, poly_only=False):893r"""894Return a list of all points `Q` such that `mQ=P` where `P` = self.895896Only points on the elliptic curve containing self and defined897over the base field are included.898899INPUT:900901- ``m`` -- a positive integer902- ``poly_only`` -- bool (default: False); if True return polynomial whose roots give all possible `x`-coordinates of `m`-th roots of self.903904OUTPUT:905906(list) -- a (possibly empty) list of solutions `Q` to `mQ=P`, where `P` = self.907908EXAMPLES:909910We find the five 5-torsion points on an elliptic curve::911912sage: E = EllipticCurve('11a'); E913Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field914sage: P = E(0); P915(0 : 1 : 0)916sage: P.division_points(5)917[(0 : 1 : 0), (5 : -6 : 1), (5 : 5 : 1), (16 : -61 : 1), (16 : 60 : 1)]918919Note above that 0 is included since [5]*0 = 0.920921We create a curve of rank 1 with no torsion and do a consistency check::922923sage: E = EllipticCurve('11a').quadratic_twist(-7)924sage: Q = E([44,-270])925sage: (4*Q).division_points(4)926[(44 : -270 : 1)]927928We create a curve over a non-prime finite field with group of order `18`::929930sage: k.<a> = GF(25)931sage: E = EllipticCurve(k, [1,2+a,3,4*a,2])932sage: P = E([3,3*a+4])933sage: factor(E.order())9342 * 3^2935sage: P.order()9369937938We find the `1`-division points as a consistency check -- there939is just one, of course::940941sage: P.division_points(1)942[(3 : 3*a + 4 : 1)]943944The point `P` has order coprime to 2 but divisible by 3, so::945946sage: P.division_points(2)947[(2*a + 1 : 3*a + 4 : 1), (3*a + 1 : a : 1)]948949We check that each of the 2-division points works as claimed::950951sage: [2*Q for Q in P.division_points(2)]952[(3 : 3*a + 4 : 1), (3 : 3*a + 4 : 1)]953954Some other checks::955956sage: P.division_points(3)957[]958sage: P.division_points(4)959[(0 : 3*a + 2 : 1), (1 : 0 : 1)]960sage: P.division_points(5)961[(1 : 1 : 1)]962963An example over a number field (see trac #3383)::964965sage: E = EllipticCurve('19a1')966sage: K.<t> = NumberField(x^9-3*x^8-4*x^7+16*x^6-3*x^5-21*x^4+5*x^3+7*x^2-7*x+1)967sage: EK = E.base_extend(K)968sage: E(0).division_points(3)969[(0 : 1 : 0), (5 : -10 : 1), (5 : 9 : 1)]970sage: EK(0).division_points(3)971[(0 : 1 : 0), (5 : 9 : 1), (5 : -10 : 1)]972sage: E(0).division_points(9)973[(0 : 1 : 0), (5 : -10 : 1), (5 : 9 : 1)]974sage: EK(0).division_points(9)975[(0 : 1 : 0), (5 : 9 : 1), (5 : -10 : 1), (-150/121*t^8 + 414/121*t^7 + 1481/242*t^6 - 2382/121*t^5 - 103/242*t^4 + 629/22*t^3 - 367/242*t^2 - 1307/121*t + 625/121 : 35/484*t^8 - 133/242*t^7 + 445/242*t^6 - 799/242*t^5 + 373/484*t^4 + 113/22*t^3 - 2355/484*t^2 - 753/242*t + 1165/484 : 1), (-150/121*t^8 + 414/121*t^7 + 1481/242*t^6 - 2382/121*t^5 - 103/242*t^4 + 629/22*t^3 - 367/242*t^2 - 1307/121*t + 625/121 : -35/484*t^8 + 133/242*t^7 - 445/242*t^6 + 799/242*t^5 - 373/484*t^4 - 113/22*t^3 + 2355/484*t^2 + 753/242*t - 1649/484 : 1), (-1383/484*t^8 + 970/121*t^7 + 3159/242*t^6 - 5211/121*t^5 + 37/484*t^4 + 654/11*t^3 - 909/484*t^2 - 4831/242*t + 6791/484 : 927/121*t^8 - 5209/242*t^7 - 8187/242*t^6 + 27975/242*t^5 - 1147/242*t^4 - 1729/11*t^3 + 1566/121*t^2 + 12873/242*t - 10871/242 : 1), (-1383/484*t^8 + 970/121*t^7 + 3159/242*t^6 - 5211/121*t^5 + 37/484*t^4 + 654/11*t^3 - 909/484*t^2 - 4831/242*t + 6791/484 : -927/121*t^8 + 5209/242*t^7 + 8187/242*t^6 - 27975/242*t^5 + 1147/242*t^4 + 1729/11*t^3 - 1566/121*t^2 - 12873/242*t + 10629/242 : 1), (-4793/484*t^8 + 6791/242*t^7 + 10727/242*t^6 - 18301/121*t^5 + 2347/484*t^4 + 2293/11*t^3 - 7311/484*t^2 - 17239/242*t + 26767/484 : 30847/484*t^8 - 21789/121*t^7 - 34605/121*t^6 + 117164/121*t^5 - 10633/484*t^4 - 29437/22*t^3 + 39725/484*t^2 + 55428/121*t - 176909/484 : 1), (-4793/484*t^8 + 6791/242*t^7 + 10727/242*t^6 - 18301/121*t^5 + 2347/484*t^4 + 2293/11*t^3 - 7311/484*t^2 - 17239/242*t + 26767/484 : -30847/484*t^8 + 21789/121*t^7 + 34605/121*t^6 - 117164/121*t^5 + 10633/484*t^4 + 29437/22*t^3 - 39725/484*t^2 - 55428/121*t + 176425/484 : 1)]976977"""978# Coerce the input m to an integer979m = rings.Integer(m)980# Check for trivial cases of m = 1, -1 and 0.981if m == 1 or m == -1:982return [self]983if m == 0:984if self == 0: # then every point Q is a solution, but...985return [self]986else:987return []988989# ans will contain the list of division points.990ans = []991992# We compute a polynomial g whose roots give all possible x993# coordinates of the m-division points. The number of994# solutions (over the algebraic closure) is m^2, assuming that995# the characteristic does not divide m.996997E = self.curve()998P = self999nP = -P1000P_is_2_torsion = (P==nP)10011002# If self is the 0, then self is a solution, and the correct1003# poly is the m'th division polynomial1004if P == 0:1005ans.append(P)1006g = E.division_polynomial(m)1007else:1008# The poly g here is 0 at x(Q) iff x(m*Q) = x(P).1009g = E._multiple_x_numerator(m) - P[0]*E._multiple_x_denominator(m)10101011# When 2*P=0, then -Q is a solution iff Q is. For even m,1012# no 2-torsion point is a solution, so that g is the1013# square of a polynomial g1 of degree m^2/2, and each root1014# of g1 leads to a pair of solutions Q, -Q to m*Q=P. For1015# odd m, P itself is the only 2-torsion solution, so g has1016# the form (x-x(P))*g1(x)^2 where g1 has degree (m^2-1)/21017# and each root of g1 leads to a pair Q, -Q.10181019if P_is_2_torsion:1020if m%2==0:1021# This computes g.sqrt() which is not implemented1022g = g.gcd(g.derivative())*g.leading_coefficient().sqrt()10231024# When 2*P!=0, then for each solution Q to m*Q=P, -Q is1025# not a solution (and points of order 2 are not1026# solutions). Hence the roots of g are distinct and each1027# gives rise to precisely one solution Q.10281029else:1030g0 = g.variables()[0] - P[0]1031g = g // g01032g = g.gcd(g.derivative())*g.leading_coefficient().sqrt()1033g = g0*g10341035if poly_only:1036return g10371038for x in g.roots(multiplicities=False):1039if E.is_x_coord(x):1040# Make a point on the curve with this x coordinate.1041Q = E.lift_x(x)1042nQ = -Q1043mQ = m*Q1044# if P==-P then Q works iff -Q works, so we include1045# both unless they are equal:1046if P_is_2_torsion:1047if mQ == P:1048ans.append(Q)1049if nQ != Q:1050ans.append(nQ)1051else:1052# P is not 2-torsion so at most one of Q, -Q works1053# and we must try both:1054if mQ == P:1055ans.append(Q)1056elif mQ == nP:1057ans.append(nQ)10581059# Finally, sort and return1060ans.sort()1061return ans10621063def _divide_out(self,p):1064r"""1065Return `(Q,k)` where `p^kQ` == self and `Q` cannot be divided by `p`.10661067..WARNING:10681069It is up to the caller to make sure that this does not loop1070endlessly. It is used in1071``EllipticCurve_generic._p_primary_torsion_basis()``, when1072self will always have (finite) order which is a power of `p`,1073so that the order of `Q` increases by a factor of `p` at each1074stage.10751076Since it will clearly be in danger of looping when1077self.is_zero(), this case is caught, but otherwise caveat1078user.10791080EXAMPLES::10811082sage: E = EllipticCurve('37a1')1083sage: P = E([0, 0])1084sage: R = 12*P1085sage: R._divide_out(2)1086((-1 : -1 : 1), 2)1087sage: R._divide_out(3)1088((2 : -3 : 1), 1)1089sage: R._divide_out(5)1090((1357/841 : 28888/24389 : 1), 0)1091sage: R._divide_out(12)1092Traceback (most recent call last):1093...1094ValueError: p (=12) should be prime.1095"""1096p = rings.Integer(p)1097if not p.is_prime():1098raise ValueError, "p (=%s) should be prime."%p10991100if self.is_zero():1101raise ValueError, "self must not be 0."11021103k=0; Q=self1104pts = Q.division_points(p)1105while len(pts) > 0:1106Q = pts[0]1107k += 11108pts = Q.division_points(p)1109return (Q,k)11101111def set_order(self, value):1112r"""1113Set the value of self._order to value.11141115Use this when you know a priori the order of this point to avoid a1116potentially expensive order calculation.11171118INPUT:11191120- ``value`` - positive Integer11211122OUTPUT:11231124None11251126EXAMPLES:11271128This example illustrates basic usage.11291130::11311132sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61133sage: G = E(5, 0)1134sage: G.set_order(2)1135sage: 2*G1136(0 : 1 : 0)11371138We now give a more interesting case, the NIST-P521 curve. Its1139order is too big to calculate with SAGE, and takes a long time1140using other packages, so it is very useful here.11411142::11431144sage: p = 2^521 - 11145sage: prev_proof_state = proof.arithmetic()1146sage: proof.arithmetic(False) # turn off primality checking1147sage: F = GF(p)1148sage: A = p - 31149sage: B = 10938490380737342745111123907668055699362075989516837489945863944959531161507350160137087375737596232485921322967063133094384525315910129121423274884789859841150sage: q = 68647976601306097149819007990813932172694353001433054093944634591855431833976553942450577463332171975329639963713633211138647686124403803403728088927070054491151sage: E = EllipticCurve([F(A), F(B)])1152sage: G = E.random_point()1153sage: G.set_order(q)1154sage: G.order() * G # This takes practically no time.1155(0 : 1 : 0)1156sage: proof.arithmetic(prev_proof_state) # restore state11571158It is an error to pass a `value` equal to `0`::11591160sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61161sage: G = E.random_point()1162sage: G.set_order(0)1163Traceback (most recent call last):1164...1165ValueError: Value 0 illegal for point order1166sage: G.set_order(1000)1167Traceback (most recent call last):1168...1169ValueError: Value 1000 illegal: outside max Hasse bound11701171It is also very likely an error to pass a value which is not the actual1172order of this point. How unlikely is determined by the factorization of1173the actual order, and the actual group structure::11741175sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61176sage: G = E(5, 0) # G has order 21177sage: G.set_order(11)1178Traceback (most recent call last):1179...1180ValueError: Value 11 illegal: 11 * (5 : 0 : 1) is not the identity11811182However, set_order can be fooled, though it's not likely in "real cases1183of interest". For instance, the order can be set to a multiple the1184actual order::11851186sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61187sage: G = E(5, 0) # G has order 21188sage: G.set_order(8)1189sage: G.order()1190811911192NOTES:11931194The implementation is based of the fact that orders of elliptic curve1195points are cached in the (pseudo-private) _order slot.11961197AUTHORS:11981199- Mariah Lenox (2011-02-16)1200"""1201E = self.curve()1202q = E.base_ring().order()1203O = E(0)1204if value == 0:1205raise ValueError('Value 0 illegal for point order')1206if value == 1:1207if self == O:1208self._order = 11209return1210else:1211raise ValueError('Value 1 illegal order for non-identity')1212low, hi = Hasse_bounds(q)1213if value > hi:1214raise ValueError('Value %s illegal: outside max Hasse bound'%value)1215if value * self != O:1216raise ValueError('Value %s illegal: %s * %s is not the identity'%(value, value, self))1217if (value - 1) * self == O:1218raise ValueError('Value %s illegal: %s * %s is the identity'%(value, value-1, self))1219self._order = value12201221############################## end ################################12221223def _line_(self,R,Q):1224r"""1225Computes the value at `Q` of a straight line through points self and `R`.12261227INPUT:12281229- ``R, Q`` -- points on self.curve() with ``Q`` nonzero.12301231OUTPUT:12321233An element of the base field self.curve().base_field().1234A ValueError is raised if ``Q`` is zero.12351236EXAMPLES::12371238sage: F.<a>=GF(2^5)1239sage: E=EllipticCurve(F,[0,0,1,1,1])1240sage: P = E(a^4 + 1, a^3)1241sage: Q = E(a^4, a^4 + a^3)1242sage: O = E(0)1243sage: P._line_(P,-2*P) == 01244True1245sage: P._line_(Q,-(P+Q)) == 01246True1247sage: O._line_(O,Q) == F(1)1248True1249sage: P._line_(O,Q) == a^4 - a^4 + 11250True1251sage: P._line_(13*P,Q) == a^41252True1253sage: P._line_(P,Q) == a^4 + a^3 + a^2 + 11254True12551256See trac #7116::12571258sage: P._line_ (Q,O)1259Traceback (most recent call last):1260...1261ValueError: Q must be nonzero.12621263..NOTES:12641265This function is used in _miller_ algorithm.12661267AUTHOR:12681269- David Hansen (2009-01-25)1270"""1271if Q.is_zero():1272raise ValueError, "Q must be nonzero."12731274if self.is_zero() or R.is_zero():1275if self == R:1276return self.curve().base_field().one_element()1277if self.is_zero():1278return Q[0] - R[0]1279if R.is_zero():1280return Q[0] - self[0]1281elif self != R:1282if self[0] == R[0]:1283return Q[0] - self[0]1284else:1285l = (R[1] - self[1])/(R[0] - self[0])1286return Q[1] - self[1] - l * (Q[0] - self[0])1287else:1288a1, a2, a3, a4, a6 = self.curve().a_invariants()1289numerator = (3*self[0]**2 + 2*a2*self[0] + a4 - a1*self[1])1290denominator = (2*self[1] + a1*self[0] + a3)1291if denominator == 0:1292return Q[0] - self[0]1293else:1294l = numerator/denominator1295return Q[1] - self[1] - l * (Q[0] - self[0])12961297def _miller_(self,Q,n):1298r"""1299Return the value at `Q` of the rational function `f_{n,P}`, where the divisor of `f_{n,P}` is `n[P]-[nP]-(n-1)[O]`.13001301INPUT:13021303- ``Q`` -- a nonzero point on self.curve().13041305- ``n`` -- an nonzero integer. If `n<0` then return `Q`1306evaluated at `1/(v_{nP}*f_{n,P})` (used in the ate pairing).13071308OUTPUT:13091310An element in the base field self.curve().base_field()1311A ValueError is raised if `Q` is zero.13121313EXAMPLES::13141315sage: F.<a>=GF(2^5)1316sage: E=EllipticCurve(F,[0,0,1,1,1])1317sage: P = E(a^4 + 1, a^3)1318sage: Fx.<b>=GF(2^(4*5))1319sage: Ex=EllipticCurve(Fx,[0,0,1,1,1])1320sage: phi=Hom(F,Fx)(F.gen().minpoly().roots(Fx)[0][0])1321sage: Px=Ex(phi(P.xy()[0]),phi(P.xy()[1]))1322sage: Qx = Ex(b^19 + b^18 + b^16 + b^12 + b^10 + b^9 + b^8 + b^5 + b^3 + 1, b^18 + b^13 + b^10 + b^8 + b^5 + b^4 + b^3 + b)1323sage: Px._miller_(Qx,41) == b^17 + b^13 + b^12 + b^9 + b^8 + b^6 + b^4 + 11324True1325sage: Qx._miller_(Px,41) == b^13 + b^10 + b^8 + b^7 + b^6 + b^51326True1327sage: P._miller_(E(0),41)1328Traceback (most recent call last):1329...1330ValueError: Q must be nonzero.13311332An example of even order::13331334sage: F.<a> = GF(19^4)1335sage: E = EllipticCurve(F,[-1,0])1336sage: P = E(15*a^3 + 17*a^2 + 14*a + 13,16*a^3 + 7*a^2 + a + 18)1337sage: Q = E(10*a^3 + 16*a^2 + 4*a + 2, 6*a^3 + 4*a^2 + 3*a + 2)1338sage: x=P.weil_pairing(Q,360)1339sage: x^360 == F(1)1340True13411342You can use the _miller_ function on linearly dependent points, but with the risk of a dividing with zero::13431344sage: Px._miller_(2*Px,41)1345Traceback (most recent call last):1346...1347ZeroDivisionError: division by zero in finite field.13481349A small example of embedding degree 6::13501351sage: q = 401; F = GF(q); a = 146; b = 400; k = 61352sage: E = EllipticCurve([F(a), F(b)])1353sage: R.<x> = F[]; K.<a> = GF(q^k, modulus=x^6 + 4*x^4 + 115*x^3 + 81*x^2 + 51*x + 3)1354sage: EK = E.base_extend(K)1355sage: P = E([F(338), F(227)])1356sage: Q_x = 333*a^5 + 391*a^4 + 160*a^3 + 335*a^2 + 71*a + 931357sage: Q_y = 343*a^5 + 273*a^4 + 26*a^3 + 342*a^2 + 340*a + 2101358sage: Q = EK([Q_x, Q_y])1359sage: P._miller_(Q, 127)1360371*a^5 + 39*a^4 + 355*a^3 + 233*a^2 + 20*a + 27513611362A series of small examples and small torsions. We start with1363`n=1`, which is trivial: the function is the constant13641::13651366sage: E = EllipticCurve([GF(7)(0), 2])1367sage: P = E(5, 1); Q = E(0, 3); I = E(0)1368sage: I._miller_(P, 1)136911370sage: I._miller_(Q, 1)1371113721373A two-torsion example. In this case `f_{n,P}(Q) = x_Q - x_P`::13741375sage: E = EllipticCurve([GF(7)(-1), 0])1376sage: P = E(0,0); Q = E(1, 0);1377sage: P._miller_(P, 2)137801379sage: Q._miller_(Q, 2)138001381sage: P._miller_(Q, 2)138211383sage: Q._miller_(P, 2)1384613851386A three-torsion example::13871388sage: E = EllipticCurve([GF(7)(0), 2])1389sage: P = E(5, 1); Q = E(0, 3);1390sage: P._miller_(Q, 3)1391413921393A 4-torsion example::13941395sage: E = EllipticCurve([GF(7)(-1), 0])1396sage: P = E(5, 1); Q = E(4, 2)1397sage: P._miller_(Q, 4)1398313991400A 5-torsion example::14011402sage: E = EllipticCurve([GF(7)(-1), 4])1403sage: P = E(4, 1); Q = E(6, 5)1404sage: P._miller_(Q, 5)1405114061407A 6-torsion example::14081409sage: E = EllipticCurve([GF(7)(3), 1])1410sage: P = E(5, 1); Q = E(3, 3)1411sage: P._miller_(Q, 6)1412514131414An example which is part of an ate pairing calculation. The trace of1415the curve is negative, so it should exercise the `n<0` case in the1416code::14171418sage: p = 2017; A = 1; B = 30; r = 29; t = -70; k = 7;1419sage: F = GF(p); R.<x> = F[]1420sage: E = EllipticCurve([F(A), F(B)]); P = E(369, 716)1421sage: K.<a> = GF(p^k, modulus=x^k+2); EK = E.base_extend(K)1422sage: Qx = 1226*a^6 + 1778*a^5 + 660*a^4 + 1791*a^3 + 1750*a^2 + 867*a + 7701423sage: Qy = 1764*a^6 + 198*a^5 + 1206*a^4 + 406*a^3 + 1200*a^2 + 273*a + 17121424sage: Q = EK(Qx, Qy)1425sage: Q._miller_(P, t-1)14261311*a^6 + 1362*a^5 + 1177*a^4 + 807*a^3 + 1331*a^2 + 1530*a + 193114271428ALGORITHM:14291430Double-and-add.14311432REFERENCES:14331434- [Mil04] Victor S. Miller, "The Weil pairing, and its efficient calculation", J. Cryptol., 17(4):235-261, 200414351436AUTHORS:14371438- David Hansen (2009-01-25)1439- Mariah Lenox (2011-03-07) -- Added more doctests and support for1440negative n.1441"""1442if Q.is_zero():1443raise ValueError, "Q must be nonzero."1444if n.is_zero():1445raise ValueError, "n must be nonzero."1446n_is_negative = False1447if n < 0:1448n = n.abs()1449n_is_negative = True14501451one = self.curve().base_field().one_element()1452t = one1453V = self1454S = 2*V1455nbin = n.bits()1456i = n.nbits() - 21457while i > -1:1458S = 2*V1459ell = V._line_(V,Q)1460vee = S._line_(-S,Q)1461t = (t**2)*(ell/vee)1462V = S1463if nbin[i] == 1:1464S = V+self1465ell = V._line_(self,Q)1466vee = S._line_(-S,Q)1467t=t*(ell/vee)1468V = S1469i=i-11470if n_is_negative:1471vee = V._line_(-V, Q)1472t = 1/(t*vee)1473return t14741475def weil_pairing(self, Q, n):1476r"""1477Compute the Weil pairing of self and `Q` using Miller's algorithm.14781479INPUT:14801481- ``Q`` -- a point on self.curve().14821483- ``n`` -- an integer `n` such that `nP = nQ = (0:1:0)` where `P` = self.14841485OUTPUT:14861487An `n`'th root of unity in the base field self.curve().base_field()14881489EXAMPLES::14901491sage: F.<a>=GF(2^5)1492sage: E=EllipticCurve(F,[0,0,1,1,1])1493sage: P = E(a^4 + 1, a^3)1494sage: Fx.<b>=GF(2^(4*5))1495sage: Ex=EllipticCurve(Fx,[0,0,1,1,1])1496sage: phi=Hom(F,Fx)(F.gen().minpoly().roots(Fx)[0][0])1497sage: Px=Ex(phi(P.xy()[0]),phi(P.xy()[1]))1498sage: O = Ex(0)1499sage: Qx = Ex(b^19 + b^18 + b^16 + b^12 + b^10 + b^9 + b^8 + b^5 + b^3 + 1, b^18 + b^13 + b^10 + b^8 + b^5 + b^4 + b^3 + b)1500sage: Px.weil_pairing(Qx,41) == b^19 + b^15 + b^9 + b^8 + b^6 + b^4 + b^3 + b^2 + 11501True1502sage: Px.weil_pairing(17*Px,41) == Fx(1)1503True1504sage: Px.weil_pairing(O,41) == Fx(1)1505True15061507An error is raised if either point is not n-torsion::15081509sage: Px.weil_pairing(O,40)1510Traceback (most recent call last):1511...1512ValueError: points must both be n-torsion15131514A larger example (see trac \#4964)::15151516sage: P,Q = EllipticCurve(GF(19^4,'a'),[-1,0]).gens()1517sage: P.order(), Q.order()1518(360, 360)1519sage: z = P.weil_pairing(Q,360)1520sage: z.multiplicative_order()152136015221523An example over a number field::15241525sage: P,Q = EllipticCurve('11a1').change_ring(CyclotomicField(5)).torsion_subgroup().gens() # long time (10s)1526sage: P,Q = (P.element(), Q.element()) # long time1527sage: (P.order(),Q.order()) # long time1528(5, 5)1529sage: P.weil_pairing(Q,5) # long time1530zeta5^21531sage: Q.weil_pairing(P,5) # long time1532zeta5^315331534ALGORITHM:15351536Implemented using Proposition 8 in [Mil04]. The value 1 is1537returned for linearly dependent input points. This condition1538is caught via a DivisionByZeroError, since the use of a1539discrete logarithm test for linear dependence, is much to slow1540for large `n`.15411542REFERENCES:15431544[Mil04] Victor S. Miller, "The Weil pairing, and its efficient1545calculation", J. Cryptol., 17(4):235-261, 200415461547AUTHOR:15481549- David Hansen (2009-01-25)1550"""1551P = self1552E = P.curve()15531554if not Q.curve() is E:1555raise ValueError, "points must both be on the same curve"15561557# Test if P, Q are both in E[n]1558if not ((n*P).is_zero() and (n*Q).is_zero()):1559raise ValueError, "points must both be n-torsion"15601561one = E.base_field().one_element()15621563# Case where P = Q1564if P == Q:1565return one15661567# Case where P = O or Q = O1568if P.is_zero() or Q.is_zero():1569return one15701571# The non-trivial case P != Q15721573# Note (2010-12-29): The following code block should not be1574# used. It attempts to reduce the pairing calculation to order1575# d = gcd(|P|,|Q|) instead of order n, but this approach is1576# counterproductive, since calculating |P| and |Q| is much1577# slower than calculating the pairing [BGN05].1578#1579# [BGN05] D. Boneh, E. Goh, and K. Nissim, "Evaluating 2-DNF1580# Formulas on Ciphertexts", TCC 2005, LNCS 3378, pp. 325-341.1581#1582# Reduction to order d = gcd(|P|,|Q|); value is a d'th root of unity1583# try:1584# nP = P.order()1585# except AttributeError:1586# nP = generic.order_from_multiple(P,n,operation='+')1587# try:1588# nQ = Q.order()1589# except AttributeError:1590# nQ = generic.order_from_multiple(Q,n,operation='+')1591# d = arith.gcd(nP,nQ)1592# if d==1:1593# return one1594#1595# P = (nP//d)*P # order d1596# Q = (nQ//d)*Q # order d1597# n = d15981599try:1600return ((-1)**n.test_bit(0))*(P._miller_(Q,n)/Q._miller_(P,n))1601except ZeroDivisionError, detail:1602return one16031604def tate_pairing(self, Q, n, k, q=None):1605r"""1606Return Tate pairing of `n`-torsion point `P = self` and point `Q`.16071608The value returned is `f_{n,P}(Q)^e` where `f_{n,P}` is a function with1609divisor `n[P]-n[O].`. This is also known as the "modified Tate1610pairing". It is a well-defined bilinear map.16111612INPUT:16131614- ``P=self`` -- Elliptic curve point having order n16151616- ``Q`` -- Elliptic curve point on same curve as P (can be any order)16171618- ``n`` -- positive integer: order of P16191620- ``k`` -- positive integer: embedding degree16211622- ``q`` -- positive integer: size of base field (the "big"1623field is `GF(q^k)`. `q` needs to be set only if its value1624cannot be deduced.)16251626OUTPUT:16271628An `n`'th root of unity in the base field self.curve().base_field()16291630EXAMPLES:16311632A simple example, pairing a point with itself, and pairing a point with1633another rational point::16341635sage: p = 103; A = 1; B = 18; E = EllipticCurve(GF(p), [A, B])1636sage: P = E(33, 91); n = P.order(); n1637191638sage: k = GF(n)(p).multiplicative_order(); k163961640sage: P.tate_pairing(P, n, k)164111642sage: Q = E(87, 51)1643sage: P.tate_pairing(Q, n, k)164411645sage: set_random_seed(35)1646sage: P.tate_pairing(P,n,k)1647116481649We now let Q be a point on the same curve as above, but defined over1650the pairing extension field, and we also demonstrate the bilinearity of1651the pairing::16521653sage: K.<a> = GF(p^k)1654sage: EK = E.base_extend(K); P = EK(P)1655sage: Qx = 69*a^5 + 96*a^4 + 22*a^3 + 86*a^2 + 6*a + 351656sage: Qy = 34*a^5 + 24*a^4 + 16*a^3 + 41*a^2 + 4*a + 401657sage: Q = EK(Qx, Qy);1658sage: # multiply by cofactor so Q has order n:1659sage: h = 551269674; Q = h*Q1660sage: P = EK(P); P.tate_pairing(Q, n, k)166124*a^5 + 34*a^4 + 3*a^3 + 69*a^2 + 86*a + 451662sage: s = Integer(randrange(1,n))1663sage: ans1 = (s*P).tate_pairing(Q, n, k)1664sage: ans2 = P.tate_pairing(s*Q, n, k)1665sage: ans3 = P.tate_pairing(Q, n, k)^s1666sage: ans1 == ans2 == ans31667True1668sage: (ans1 != 1) and (ans1^n == 1)1669True16701671Here is an example of using the Tate pairing to compute the Weil1672pairing (using the same data as above)::16731674sage: e = Integer((p^k-1)/n); e1675628448577121676sage: P.weil_pairing(Q, n)^e167794*a^5 + 99*a^4 + 29*a^3 + 45*a^2 + 57*a + 341678sage: P.tate_pairing(Q, n, k) == P._miller_(Q, n)^e1679True1680sage: Q.tate_pairing(P, n, k) == Q._miller_(P, n)^e1681True1682sage: P.tate_pairing(Q, n, k)/Q.tate_pairing(P, n, k)168394*a^5 + 99*a^4 + 29*a^3 + 45*a^2 + 57*a + 3416841685An example where we have to pass the base field size (and we again have1686agreement with the Weil pairing)::16871688sage: F.<a>=GF(2^5)1689sage: E=EllipticCurve(F,[0,0,1,1,1])1690sage: P = E(a^4 + 1, a^3)1691sage: Fx.<b>=GF(2^(4*5))1692sage: Ex=EllipticCurve(Fx,[0,0,1,1,1])1693sage: phi=Hom(F,Fx)(F.gen().minpoly().roots(Fx)[0][0])1694sage: Px=Ex(phi(P.xy()[0]),phi(P.xy()[1]))1695sage: Qx = Ex(b^19+b^18+b^16+b^12+b^10+b^9+b^8+b^5+b^3+1, b^18+b^13+b^10+b^8+b^5+b^4+b^3+b)1696sage: Px.tate_pairing(Qx, n=41, k=4)1697Traceback (most recent call last):1698...1699ValueError: Unexpected field degree: set keyword argument q equal to the size of the base field (big field is GF(q^4)).1700sage: num = Px.tate_pairing(Qx, n=41, k=4, q=32); num1701b^19 + b^14 + b^13 + b^12 + b^6 + b^4 + b^31702sage: den = Qx.tate_pairing(Px, n=41, k=4, q=32); den1703b^19 + b^17 + b^16 + b^15 + b^14 + b^10 + b^6 + b^2 + 11704sage: e = Integer((32^4-1)/41); e1705255751706sage: Px.weil_pairing(Qx, 41)^e == num/den1707True17081709NOTES:17101711This function uses Miller's algorithm, followed by a naive1712exponentiation. It does not do anything fancy. In the case1713that there is an issue with `Q` being on one of the lines1714generated in the `r*P` calculation, `Q` is offset by a random1715point `R` and P.tate_pairing(Q+R,n,k)/P.tate_pairing(R,n,k)1716is returned.17171718AUTHORS:17191720- Mariah Lenox (2011-03-07)1721"""1722P = self1723E = P.curve()17241725if not Q.curve() is E:1726raise ValueError, "Points must both be on the same curve"17271728K = E.base_ring()1729d = K.degree()1730if q == None:1731if d == 1:1732q = K.order()1733elif d == k:1734q = K.base_ring().order()1735else:1736raise ValueError("Unexpected field degree: set keyword argument q equal to the size of the base field (big field is GF(q^%s))."%k)17371738if n*P != E(0):1739raise ValueError('This point is not of order n=%s' %n)17401741# In small cases, or in the case of pairing an element with1742# itself, Q could be on one of the lines in the Miller1743# algorithm. If this happens we try again, with an offset of a1744# random point.1745try:1746ret = self._miller_(Q, n)1747e = rings.Integer((q**k - 1)/n)1748ret = ret**e1749except (ZeroDivisionError, ValueError):1750R = E.random_point()1751ret = self.tate_pairing(Q + R, n, k)/self.tate_pairing(R, n, k)1752return ret17531754def ate_pairing(self, Q, n, k, t, q=None):1755r"""1756Return ate pairing of `n`-torsion points `P=self` and `Q`.17571758Also known as the `n`-th modified ate pairing. `P` is `GF(q)`-rational,1759and `Q` must be an element of `Ker(\pi-p)`, where `\pi` is the1760`q`-frobenius map (and hence `Q` is `GF(q^k)`-rational).17611762INPUT:17631764- ``P=self`` -- a point of order `n`, in `ker(\pi-1)`, where1765`\pi` is the `q`-Frobenius map (e.g., `P` is `q-rational`).17661767- ``Q`` -- a point of order `n` in `ker(\pi-q)`17681769- ``n`` -- the order of `P` and `Q`.17701771- ``k`` -- the embedding degree.17721773- ``t`` -- the trace of Frobenius of the curve over `GF(q)`.17741775- ``q`` -- (default:None) the size of base field (the "big"1776field is `GF(q^k)`). `q` needs to be set only if its value1777cannot be deduced.17781779OUTPUT:17801781FiniteFieldElement in `GF(q^k)` -- the ate pairing of `P` and `Q`.17821783EXAMPLES:17841785An example with embedding degree 6::17861787sage: p = 7549; A = 0; B = 1; n = 157; k = 6; t = 141788sage: F = GF(p); E = EllipticCurve(F, [A, B])1789sage: R.<x> = F[]; K.<a> = GF(p^k, modulus=x^k+2)1790sage: EK = E.base_extend(K)1791sage: P = EK(3050, 5371); Q = EK(6908*a^4, 3231*a^3)1792sage: P.ate_pairing(Q, n, k, t)17936708*a^5 + 4230*a^4 + 4350*a^3 + 2064*a^2 + 4022*a + 67331794sage: s = Integer(randrange(1, n))1795sage: (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)1796True1797sage: P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)^s1798True17991800Another example with embedding degree 7 and positive trace::18011802sage: p = 2213; A = 1; B = 49; n = 1093; k = 7; t = 281803sage: F = GF(p); E = EllipticCurve(F, [A, B])1804sage: R.<x> = F[]; K.<a> = GF(p^k, modulus=x^k+2)1805sage: EK = E.base_extend(K)1806sage: P = EK(1583, 1734)1807sage: Qx = 1729*a^6+1767*a^5+245*a^4+980*a^3+1592*a^2+1883*a+7221808sage: Qy = 1299*a^6+1877*a^5+1030*a^4+1513*a^3+1457*a^2+309*a+16361809sage: Q = EK(Qx, Qy)1810sage: P.ate_pairing(Q, n, k, t)18111665*a^6 + 1538*a^5 + 1979*a^4 + 239*a^3 + 2134*a^2 + 2151*a + 6541812sage: s = Integer(randrange(1, n))1813sage: (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)1814True1815sage: P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)^s1816True18171818Another example with embedding degree 7 and negative trace::18191820sage: p = 2017; A = 1; B = 30; n = 29; k = 7; t = -701821sage: F = GF(p); E = EllipticCurve(F, [A, B])1822sage: R.<x> = F[]; K.<a> = GF(p^k, modulus=x^k+2)1823sage: EK = E.base_extend(K)1824sage: P = EK(369, 716)1825sage: Qx = 1226*a^6+1778*a^5+660*a^4+1791*a^3+1750*a^2+867*a+7701826sage: Qy = 1764*a^6+198*a^5+1206*a^4+406*a^3+1200*a^2+273*a+17121827sage: Q = EK(Qx, Qy)1828sage: P.ate_pairing(Q, n, k, t)18291794*a^6 + 1161*a^5 + 576*a^4 + 488*a^3 + 1950*a^2 + 1905*a + 13151830sage: s = Integer(randrange(1, n))1831sage: (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)1832True1833sage: P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)^s1834True18351836Using the same data, we show that the ate pairing is a power of the1837Tate pairing (see [HSV]_ end of section 3.1)::18381839sage: c = (k*p^(k-1)).mod(n); T = t - 11840sage: N = gcd(T^k - 1, p^k - 1)1841sage: s = Integer(N/n)1842sage: L = Integer((T^k - 1)/N)1843sage: M = (L*s*c.inverse_mod(n)).mod(n)1844sage: P.ate_pairing(Q, n, k, t) == Q.tate_pairing(P, n, k)^M1845True18461847An example where we have to pass the base field size (and we again have1848agreement with the Tate pairing). Note that though `Px` is not1849`F`-rational, (it is the homomorphic image of an `F`-rational point) it1850is nonetheless in `ker(\pi-1)`, and so is a legitimate input::18511852sage: q = 2^5; F.<a>=GF(q)1853sage: n = 41; k = 4; t = -81854sage: E=EllipticCurve(F,[0,0,1,1,1])1855sage: P = E(a^4 + 1, a^3)1856sage: Fx.<b>=GF(q^k)1857sage: Ex=EllipticCurve(Fx,[0,0,1,1,1])1858sage: phi=Hom(F,Fx)(F.gen().minpoly().roots(Fx)[0][0])1859sage: Px=Ex(phi(P.xy()[0]),phi(P.xy()[1]))1860sage: Qx = Ex(b^19+b^18+b^16+b^12+b^10+b^9+b^8+b^5+b^3+1, b^18+b^13+b^10+b^8+b^5+b^4+b^3+b)1861sage: Qx = Ex(Qx[0]^q, Qx[1]^q) - Qx # ensure Qx is in ker(pi - q)1862sage: Px.ate_pairing(Qx, n, k, t)1863Traceback (most recent call last):1864...1865ValueError: Unexpected field degree: set keyword argument q equal to the size of the base field (big field is GF(q^4)).1866sage: Px.ate_pairing(Qx, n, k, t, q)1867b^19 + b^18 + b^17 + b^16 + b^15 + b^14 + b^13 + b^12 + b^11 + b^9 + b^8 + b^5 + b^4 + b^2 + b + 11868sage: s = Integer(randrange(1, n))1869sage: (s*Px).ate_pairing(Qx, n, k, t, q) == Px.ate_pairing(s*Qx, n, k, t, q)1870True1871sage: Px.ate_pairing(s*Qx, n, k, t, q) == Px.ate_pairing(Qx, n, k, t, q)^s1872True1873sage: c = (k*q^(k-1)).mod(n); T = t - 11874sage: N = gcd(T^k - 1, q^k - 1)1875sage: s = Integer(N/n)1876sage: L = Integer((T^k - 1)/N)1877sage: M = (L*s*c.inverse_mod(n)).mod(n)1878sage: Px.ate_pairing(Qx, n, k, t, q) == Qx.tate_pairing(Px, n, k, q)^M1879True18801881It is an error if `Q` is not in the kernel of `\pi - p`, where `\pi` is1882the Frobenius automorphism::18831884sage: p = 29; A = 1; B = 0; n = 5; k = 2; t = 101885sage: F = GF(p); R.<x> = F[]1886sage: E = EllipticCurve(F, [A, B]);1887sage: K.<a> = GF(p^k, modulus=x^k+2); EK = E.base_extend(K)1888sage: P = EK(13, 8); Q = EK(13, 21)1889sage: P.ate_pairing(Q, n, k, t)1890Traceback (most recent call last):1891...1892ValueError: Point (13 : 21 : 1) not in Ker(pi - q)18931894It is also an error if `P` is not in the kernel os `\pi - 1`::18951896sage: p = 29; A = 1; B = 0; n = 5; k = 2; t = 101897sage: F = GF(p); R.<x> = F[]1898sage: E = EllipticCurve(F, [A, B]);1899sage: K.<a> = GF(p^k, modulus=x^k+2); EK = E.base_extend(K)1900sage: P = EK(14, 10*a); Q = EK(13, 21)1901sage: P.ate_pairing(Q, n, k, t)1902Traceback (most recent call last):1903...1904ValueError: This point (14 : 10*a : 1) is not in Ker(pi - 1)19051906NOTES:19071908First defined in the paper of [HSV]_, the ate pairing can be1909computationally effective in those cases when the trace of the curve1910over the base field is significantly smaller than the expected1911value. This implementation is simply Miller's algorithm followed by a1912naive exponentiation, and makes no claims towards efficiency.19131914REFERENCES:19151916.. [HSV] Hess, Smart, Vercauteren, "The Eta Pairing Revisited",1917IEEE Trans. Information Theory, 52(10): 4595-4602, 2006.19181919AUTHORS:19201921- Mariah Lenox (2011-03-08)1922"""1923P = self1924# check for same curve1925E = P.curve()1926O = E(0)1927if not Q.curve() is E:1928raise ValueError, "Points must both be on the same curve"19291930# set q to be the order of the base field1931K = E.base_ring()1932F = K.base_ring()1933if q == None:1934d = K.degree()1935if d == k:1936q = K.base_ring().order()1937else:1938raise ValueError("Unexpected field degree: set keyword argument q equal to the size of the base field (big field is GF(q^%s))." %k)19391940# check order of P1941if n*P != O:1942raise ValueError('This point %s is not of order n=%s' % (P,n))19431944# check for P in kernel pi - 1:1945piP = E(P[0]**q, P[1]**q)1946if piP - P != O:1947raise ValueError('This point %s is not in Ker(pi - 1)' % P)19481949# check for Q in kernel pi - q:1950piQ = E(Q[0]**q, Q[1]**q)1951if piQ - q*Q != O:1952raise ValueError('Point %s not in Ker(pi - q)' % Q)19531954T = t-11955ret = Q._miller_(P, T)1956e = rings.Integer((q**k - 1)/n)1957ret = ret**e1958return ret19591960class EllipticCurvePoint_number_field(EllipticCurvePoint_field):1961"""1962A point on an elliptic curve over a number field.19631964Most of the functionality is derived from the parent class1965``EllipticCurvePoint_field``. In addition we have support for the1966order of a point, and heights (currently only implemented over1967`\QQ`).19681969EXAMPLES::19701971sage: E = EllipticCurve('37a')1972sage: E([0,0])1973(0 : 0 : 1)1974sage: E(0,0) # brackets are optional1975(0 : 0 : 1)1976sage: E([GF(5)(0), 0]) # entries are coerced1977(0 : 0 : 1)19781979sage: E(0.000, 0)1980(0 : 0 : 1)19811982sage: E(1,0,0)1983Traceback (most recent call last):1984...1985TypeError: Coordinates [1, 0, 0] do not define a point on1986Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field19871988::19891990sage: E = EllipticCurve([0,0,1,-1,0])1991sage: S = E(QQ); S1992Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field19931994TESTS::19951996sage: loads(S.dumps()) == S1997True1998sage: P = E(0,0); P1999(0 : 0 : 1)2000sage: loads(P.dumps()) == P2001True2002sage: T = 100*P2003sage: loads(T.dumps()) == T2004True20052006Test pickling an elliptic curve that has known points on it::20072008sage: e = EllipticCurve([0, 0, 1, -1, 0]); g = e.gens(); loads(dumps(e)) == e2009True2010"""20112012def order(self):2013r"""2014Return the order of this point on the elliptic curve.20152016If the point has infinite order, returns +Infinity. For2017curves defined over `\QQ`, we call PARI; over other2018number fields we implement the function here.20192020.. note::20212022:meth:`additive_order` is a synonym for :meth:`order`20232024EXAMPLES::20252026sage: E = EllipticCurve([0,0,1,-1,0])2027sage: P = E([0,0]); P2028(0 : 0 : 1)2029sage: P.order()2030+Infinity20312032::20332034sage: E = EllipticCurve([0,1])2035sage: P = E([-1,0])2036sage: P.order()203722038sage: P.additive_order()2039220402041"""2042try:2043return self._order2044except AttributeError:2045pass20462047if self.is_zero():2048self._order = rings.Integer(1)2049return self._order20502051E = self.curve()20522053# Special code for curves over Q, calling PARI2054try:2055n = int(E.pari_curve().ellorder(self))2056if n == 0: n = oo2057self._order = n2058return n2059except PariError:2060pass20612062# Get the torsion order if known, else a bound on (multiple2063# of) the order. We do not compute the torsion if it is not2064# already known, since computing the bound is faster (and is2065# also cached).20662067try:2068N = E._torsion_order2069except AttributeError:2070N = E._torsion_bound()20712072# Now self is a torsion point iff it is killed by N:2073if not (N*self).is_zero():2074self._order = oo2075return self._order20762077# Finally we find the exact order using the generic code:2078self._order = generic.order_from_multiple(self,N,operation='+')2079return self._order20802081additive_order = order20822083def has_finite_order(self):2084"""2085Return True iff this point has finite order on the elliptic curve.20862087EXAMPLES::20882089sage: E = EllipticCurve([0,0,1,-1,0])2090sage: P = E([0,0]); P2091(0 : 0 : 1)2092sage: P.has_finite_order()2093False20942095::20962097sage: E = EllipticCurve([0,1])2098sage: P = E([-1,0])2099sage: P.has_finite_order()2100True2101"""2102if self.is_zero(): return True2103return self.order() != oo21042105def has_infinite_order(self):2106r"""2107Return True iff this point has infinite order on the elliptic curve.21082109EXAMPLES::21102111sage: E = EllipticCurve([0,0,1,-1,0])2112sage: P = E([0,0]); P2113(0 : 0 : 1)2114sage: P.has_infinite_order()2115True21162117::21182119sage: E = EllipticCurve([0,1])2120sage: P = E([-1,0])2121sage: P.has_infinite_order()2122False2123"""2124if self.is_zero(): return False2125return self.order() == oo21262127def is_on_identity_component(self, embedding=None):2128r"""2129Returns True iff this point is on the identity component of2130its curve with respect to a given (real or complex) embedding.21312132INPUT:21332134- ``self`` -- a point on a curve over any ordered field (e.g. `\QQ`)21352136- ``embedding`` -- an embedding from the base_field of the2137point's curve into `\RR` or `\CC`; if None (the2138default) it uses the first embedding of the base_field into2139`\RR` if any, else the first embedding into `\CC`.21402141OUTPUT:21422143(bool) -- True iff the point is on the identity component of2144the curve. (If the point is zero then the result is True.)21452146EXAMPLES:21472148For `K=\QQ` there is no need to specify an embedding::21492150sage: E=EllipticCurve('5077a1')2151sage: [E.lift_x(x).is_on_identity_component() for x in range(-3,5)]2152[False, False, False, False, False, True, True, True]21532154An example over a field with two real embeddings::21552156sage: L.<a> = QuadraticField(2)2157sage: E=EllipticCurve(L,[0,1,0,a,a])2158sage: P=E(-1,0)2159sage: [P.is_on_identity_component(e) for e in L.embeddings(RR)]2160[False, True]21612162We can check this as follows::21632164sage: [e(E.discriminant())>0 for e in L.embeddings(RR)]2165[True, False]2166sage: e = L.embeddings(RR)[0]2167sage: E1 = EllipticCurve(RR,[e(ai) for ai in E.ainvs()])2168sage: e1,e2,e3 = E1.two_division_polynomial().roots(RR,multiplicities=False)2169sage: e1 < e2 < e3 and e(P[0]) < e32170True21712172"""2173if self.is_zero(): # trivial case2174return True21752176e = embedding2177# It is also trivially true if we have a complex embedding2178if not e is None:2179if not rings.is_RealField(e.codomain()):2180return True21812182# find a suitable embedding if none was supplied:2183E = self.curve()2184K = E.base_field()2185if e is None:2186try:2187e = K.embeddings(rings.RealField())[0]2188except IndexError:2189e = K.embeddings(rings.ComplexField())[0]21902191# If there is only one component, the result is True:2192if not rings.is_RealField(e.codomain()): # complex embedding2193return True2194if e(E.discriminant()) < 0: # only one component2195return True21962197# Now we have a real embedding and two components and have to work:2198gx = E.two_division_polynomial()2199gxd = gx.derivative()2200gxdd = gxd.derivative()2201return ( e(gxd(self[0])) > 0 and e(gxdd(self[0])) > 0)22022203def has_good_reduction(self, P=None):2204r"""2205Returns True iff this point has good reduction modulo a prime.22062207INPUT:22082209- ``P`` -- a prime of the base_field of the point's curve, or2210None (default)22112212OUTPUT:22132214(bool) If a prime `P` of the base field is specified, returns2215True iff the point has good reduction at `P`; otherwise,2216return true if the point has god reduction at all primes in2217the support of the discriminant of this model.22182219EXAMPLES::22202221sage: E = EllipticCurve('990e1')2222sage: P = E.gen(0); P2223(15 : 51 : 1)2224sage: [E.has_good_reduction(p) for p in [2,3,5,7]]2225[False, False, False, True]2226sage: [P.has_good_reduction(p) for p in [2,3,5,7]]2227[True, False, True, True]2228sage: [E.tamagawa_exponent(p) for p in [2,3,5,7]]2229[2, 2, 1, 1]2230sage: [(2*P).has_good_reduction(p) for p in [2,3,5,7]]2231[True, True, True, True]2232sage: P.has_good_reduction()2233False2234sage: (2*P).has_good_reduction()2235True2236sage: (3*P).has_good_reduction()2237False22382239::22402241sage: K.<i> = NumberField(x^2+1)2242sage: E = EllipticCurve(K,[0,1,0,-160,308])2243sage: P = E(26,-120)2244sage: E.discriminant().support()2245[Fractional ideal (i + 1),2246Fractional ideal (-i - 2),2247Fractional ideal (i - 2),2248Fractional ideal (3)]2249sage: [E.tamagawa_exponent(p) for p in E.discriminant().support()]2250[1, 4, 4, 4]2251sage: P.has_good_reduction()2252False2253sage: (2*P).has_good_reduction()2254False2255sage: (4*P).has_good_reduction()2256True22572258TESTS:22592260An example showing that #8498 is fixed::22612262sage: E = EllipticCurve('11a1')2263sage: K.<t> = NumberField(x^2+47)2264sage: EK = E.base_extend(K)2265sage: T = EK(5,5)2266sage: P = EK(-2, -1/2*t - 1/2)2267sage: p = K.ideal(11)2268sage: T.has_good_reduction(p)2269False2270sage: P.has_good_reduction(p)2271True2272"""2273if self.is_zero(): # trivial case2274return True22752276E = self.curve()2277if P is None:2278return all([self.has_good_reduction(P) for P in E.discriminant().support()])2279K = E.base_field()2280from sage.schemes.elliptic_curves.ell_local_data import check_prime2281P = check_prime(K,P)22822283# If the curve has good reduction at P, the result is True:2284t = E.local_data(P).bad_reduction_type()2285if t is None:2286return True22872288# Make sure the curve is integral and locally minimal at P:2289Emin = E.local_minimal_model(P)2290urst = E.isomorphism_to(Emin)2291Q = urst(self)22922293# Scale the homogeneous coordinates of the point to be primitive:2294xyz = list(Q)2295e = min([c.valuation(P) for c in xyz])2296if e !=0:2297if K is rings.QQ:2298pi = P2299else:2300pi = K.uniformizer(P)2301pie = pi**e2302xyz = [c/pie for c in xyz]23032304# Evaluate the partial derivatives at the point to see if they2305# are zero mod P23062307# See #8498: sometimes evaluating F's derivatives at xyz2308# returns a constant polynomial instead of a constant23092310F = Emin.defining_polynomial()2311for v in F.variables():2312c = (F.derivative(v))(xyz)2313try:2314val = c.valuation(P)2315except AttributeError:2316val = c.constant_coefficient().valuation(P)2317if val == 0:2318return True2319return False23202321def reduction(self,p):2322"""2323This finds the reduction of a point `P` on the elliptic curve modulo the prime `p`.23242325INPUT:23262327- ``self`` -- A point on an elliptic curve.23282329- ``p`` -- a prime number23302331OUTPUT:23322333The point reduced to be a point on the elliptic curve modulo `p`.23342335EXAMPLES::23362337sage: E = EllipticCurve([1,2,3,4,0])2338sage: P = E(0,0)2339sage: P.reduction(5)2340(0 : 0 : 1)2341sage: Q = E(98,931)2342sage: Q.reduction(5)2343(3 : 1 : 1)2344sage: Q.reduction(5).curve() == E.reduction(5)2345True23462347::23482349sage: F.<a> = NumberField(x^2+5)2350sage: E = EllipticCurve(F,[1,2,3,4,0])2351sage: Q = E(98,931)2352sage: Q.reduction(a)2353(3 : 1 : 1)2354sage: Q.reduction(11)2355(10 : 7 : 1)23562357::23582359sage: F.<a> = NumberField(x^3+x^2+1)2360sage: E = EllipticCurve(F,[a,2])2361sage: P = E(a,1)2362sage: P.reduction(F.ideal(5))2363(abar : 1 : 1)2364sage: P.reduction(F.ideal(a^2-4*a-2))2365(abar : 1 : 1)23662367"""2368P = self2369E = P.curve()2370return E.reduction(p)(P)23712372def height(self, precision=None):2373"""2374The Neron-Tate canonical height of the point.23752376INPUT:23772378- ``self`` -- a point on a curve over a number field23792380- ``precision`` -- (int or None (default)): the precision in2381bits of the result (default real precision if None)23822383OUTPUT:23842385The rational number 0, or a nonzero real field element.23862387The returned height is normalized to be independant of the base field.2388Fixing this, there are two normalizations used in the literature,2389one of which is double the other. We use the larger of the two,2390which is the one appropriate for the BSD conjecture. This is consistant2391with [Cre] and double that of [Sil].23922393REFERENCES:23942395- [Cre] John Cremona, Algorithms for modular elliptic curves,2396Cambridge University Press, 1997.23972398- [Sil] Silverman, Joseph H. The arithmetic of elliptic curves.2399Second edition. Graduate Texts in Mathematics, 106. Springer, 2009.24002401EXAMPLES::24022403sage: E = EllipticCurve('11a'); E2404Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field2405sage: P = E([5,5]); P2406(5 : 5 : 1)2407sage: P.height()240802409sage: Q = 5*P2410sage: Q.height()2411024122413::24142415sage: E = EllipticCurve('37a'); E2416Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field2417sage: P = E([0,0])2418sage: P.height()24190.05111140823996882420sage: P.order()2421+Infinity2422sage: E.regulator()24230.0511114082399688...24242425sage: def naive_height(P):2426... return log(RR(max(abs(P[0].numerator()), abs(P[0].denominator()))))2427sage: for n in [1..10]:2428... print naive_height(2^n*P)/4^n24290.00000000000000024300.043321698784996624310.050294934763565624320.051100633561864524330.051100783479961224340.051101366615246624350.051103419990774324360.051110649290647124370.051111408154108224380.051111408154118024392440::24412442sage: E = EllipticCurve('4602a1'); E2443Elliptic Curve defined by y^2 + x*y = x^3 + x^2 - 37746035*x - 89296920339 over Rational Field2444sage: x = 779859224589749492468582291959451034715902445sage: y = 195752602300153137022613790221516759619651571089202635945452232446sage: d = 22540207618847822432447sage: E([ x / d^2, y / d^3 ]).height()244886.740656138127524492450::24512452sage: E = EllipticCurve([17, -60, -120, 0, 0]); E2453Elliptic Curve defined by y^2 + 17*x*y - 120*y = x^3 - 60*x^2 over Rational Field2454sage: E([30, -90]).height()2455024562457sage: E = EllipticCurve('389a1'); E2458Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field2459sage: [P,Q] = [E(-1,1),E(0,-1)]2460sage: P.height(precision=100)24610.686667083305586585723552102952462sage: (3*Q).height(precision=100)/Q.height(precision=100)24639.00000000000000000000000000002464sage: _.parent()2465Real Field with 100 bits of precision24662467Canonical heights over number fields are implemented as well::24682469sage: R.<x> = QQ[]2470sage: K.<a> = NumberField(x^3-2)2471sage: E = EllipticCurve([a, 4]); E2472Elliptic Curve defined by y^2 = x^3 + a*x + 4 over Number Field in a with defining polynomial x^3 - 22473sage: P = E((0,2))2474sage: P.height()24750.8104630965859252476sage: P.height(precision=100)24770.810463096585925368639918105772478sage: P.height(precision=200)24790.810463096585925368639918105768651588961302864171558323780862480sage: (2*P).height() / P.height()24814.000000000000002482sage: (100*P).height() / P.height()248310000.000000000024842485Some consistency checks::24862487sage: E = EllipticCurve('5077a1')2488sage: P = E([-2,3,1])2489sage: P.height()24901.3685725053539324912492sage: EK = E.change_ring(QuadraticField(-3,'a'))2493sage: PK = EK([-2,3,1])2494sage: PK.height()24951.3685725053539324962497sage: K.<i> = NumberField(x^2+1)2498sage: E = EllipticCurve(K, [0,0,4,6*i,0])2499sage: Q = E.lift_x(-9/4); Q2500(-9/4 : 27/8*i - 4 : 1)2501sage: Q.height()25022.695185600179092503sage: (15*Q).height() / Q.height()2504225.00000000000025052506sage: E = EllipticCurve('37a')2507sage: P = E([0,-1])2508sage: P.height()25090.05111140823996882510sage: K.<a> = QuadraticField(-7)2511sage: ED = E.quadratic_twist(-7)2512sage: Q = E.isomorphism_to(ED.change_ring(K))(P); Q2513(0 : -7/2*a - 1/2 : 1)2514sage: Q.height()25150.05111140823996882516sage: Q.height(precision=100)25170.05111140823996884023588609975525182519An example to show that the bug at \#5252 is fixed::25202521sage: E = EllipticCurve([1, -1, 1, -2063758701246626370773726978, 32838647793306133075103747085833809114881])2522sage: P = E([-30987785091199, 258909576181697016447])2523sage: P.height()252425.86031706754622525sage: P.height(precision=100)252625.8603170675461907438688407412527sage: P.height(precision=250)252825.8603170675461907438688407407351103230988729038444162155771710417835725132529sage: P.height(precision=500)253025.860317067546190743868840740735110323098872903844416215577171041783572512955113057088981328179215727850763990997211285601919023612536291419545232172025312532An example to show that the bug at \#8319 is fixed (correct height when the curve is not minimal)::25332534sage: E = EllipticCurve([-5580472329446114952805505804593498080000,-157339733785368110382973689903536054787700497223306368000000])2535sage: xP = 204885147732879546487576840131729064308289385547094673627174585676211859152978311600/236255019070579481322622171889836812048569076577531784154303612536sage: P = E.lift_x(xP)2537sage: P.height()2538157.4325985167542539sage: Q = 2*P2540sage: Q.height() # long time (4s)2541629.7303940670162542sage: Q.height()-4*P.height() # long time25430.00000000000000025442545"""2546if self.has_finite_order():2547return rings.QQ(0)25482549if precision is None:2550precision = rings.RealField().precision()25512552try:2553height = self.__height2554if height.prec() >= precision:2555return rings.RealField(precision)(height)2556except AttributeError:2557pass25582559if self.curve().base_ring() is rings.QQ:2560E = self.curve()2561Emin = E.minimal_model()2562iso = E.isomorphism_to(Emin)2563P = iso(self)2564h = Emin.pari_curve(prec=precision).ellheight(P, precision=precision)2565height = rings.RealField(precision)(h)2566else:2567from sage.misc.stopgap import stopgap2568stopgap("Computation of heights on elliptic curves over fields other than Q can return very imprecise results.", 12509)2569height = (self.nonarchimedian_local_height(prec=precision)2570+ self.archimedian_local_height(prec=precision))25712572self.__height = height2573return height257425752576def archimedian_local_height(self, v=None, prec=None):2577"""2578Computes the local height of self at the archimedian place `v`.25792580If `v` is None, returns the (weighted) sum of all archimedian2581contributions to the height.25822583The normalization is taken to be independant of the base field,2584but twice that in the paper. Note also that this local height depends2585on the model of the curve.25862587INPUT:25882589- ``v`` -- a real or complex embedding, or None25902591- ``prec`` -- the precision of the computation. If None, the precision2592is deduced from v.25932594ALGORITHM:25952596See section 4 of Silverman, J. Computing Heights on Elliptic Curves.2597Mathematics of Computation, Vol. 51, No. 183 (Jul., 1988), pp. 339-358259825992600EXAMPLES:26012602Examples 1, 2, and 3 from the above paper::26032604sage: K.<a> = QuadraticField(-2)2605sage: E = EllipticCurve(K, [0,-1,1,0,0]); E2606Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 over Number Field in a with defining polynomial x^2 + 22607sage: P = E.lift_x(2+a); P2608(a + 2 : 2*a + 1 : 1)2609sage: P.archimedian_local_height(K.places(prec=170)[0]) / 226100.4575477328752327673621121074142365434657602981469526112612sage: K.<i> = NumberField(x^2+1)2613sage: E = EllipticCurve(K, [0,0,4,6*i,0]); E2614Elliptic Curve defined by y^2 + 4*y = x^3 + 6*i*x over Number Field in i with defining polynomial x^2 + 12615sage: P = E((0,0))2616sage: P.archimedian_local_height(K.places()[0]) / 226170.51018499516237326182619sage: Q = E.lift_x(-9/4); Q2620(-9/4 : 27/8*i - 4 : 1)2621sage: Q.archimedian_local_height(K.places()[0]) / 226220.6544456195296002623"""2624if self.has_finite_order():2625return QQ(0)26262627if v is None:2628K = self.curve().base_ring()2629def local_degree(v):2630"""2631Computes [ K_v : Q_v ]2632"""2633return 2 - int(v.im_gens()[0] in rings.RR)2634return sum(local_degree(v) * self.archimedian_local_height(v, prec) for v in K.places(prec=prec)) / K.degree()26352636if prec is None:2637prec = v.codomain().prec()2638E = self.curve()2639b2, b4, b6, b8 = [v(b) for b in E.b_invariants()]2640H = max(4, abs(b2), 2*abs(b4), 2*abs(b6), abs(b8))2641nterms = int(math.ceil(0.5*prec + 0.5 + 3*math.log(7+4*math.log(H)/3)/4 + math.log(max(1, ~abs(v(E.discriminant()))))/3))2642b2p = b2 - 122643b4p = b4 - b2 + 62644b6p = b6 - 2*b4 + b2 - 42645b8p = b8 - 3*b6 + 3*b4 - b2 + 326462647t = rings.polygen(v.codomain())2648fw = 4*t + b2 * t**2 + 2*b4 * t**3 + b6 * t**42649fz = 1 - b4 * t**2 - 2*b6 * t**3 - b8 * t**42650fwp = 4*t + b2p * t**2 + 2*b4p * t**3 + b6p * t**42651fzp = 1 - b4p * t**2 - 2*b6p * t**3 - b8p * t**426522653x = v(self[0])2654if abs(x) >= .5:2655t = 1/x2656beta = True2657else:2658t = 1/(x+1)2659beta = False2660lam = -abs(t).log()2661mu = 02662four_to_n = rings.QQ(1)26632664for n in range(nterms):2665if beta:2666w = fw(t)2667z = fz(t)2668if abs(w) <= 2 * abs(z):2669mu += four_to_n * abs(z).log()2670t = w/z2671else:2672mu += four_to_n * abs(z+w).log()2673t = w/(z+w)2674beta = not beta2675else:2676w = fwp(t)2677z = fzp(t)2678if abs(w) <= 2 * abs(z):2679mu += four_to_n * abs(z).log()2680t = w/z2681else:2682mu += four_to_n * abs(z-w).log()2683t = w/(z-w)2684beta = not beta2685four_to_n >>= 22686return lam + mu/426872688def nonarchimedian_local_height(self, v=None, prec=None, weighted=False, is_minimal=None):2689"""2690Computes the local height of self at the non-archimedian place `v`.26912692If `v` is None, returns the (weighted) sum of all non-archimedian2693contributions to the height.26942695The normalization is taken to be independant of the base field,2696but twice that in the paper. Note also that this local height depends2697on the model of the curve.26982699INPUT:27002701- ``v`` -- a non-archimedean place of the base field of the curve,2702or None, in which case the total nonarchimedian contribution2703is returned27042705- ``prec`` -- working precision, or None in which case the height2706is returned symbolically27072708ALGORITHM:27092710See section 5 of Silverman, J. Computing Heights on Elliptic Curves.2711Mathematics of Computation, Vol. 51, No. 183 (Jul., 1988), pp. 339-35827122713EXAMPLES:27142715Examples 2 and 3 from the above paper::27162717sage: K.<i> = NumberField(x^2+1)2718sage: E = EllipticCurve(K, [0,0,4,6*i,0]); E2719Elliptic Curve defined by y^2 + 4*y = x^3 + 6*i*x over Number Field in i with defining polynomial x^2 + 12720sage: P = E((0,0))2721sage: P.nonarchimedian_local_height(K.ideal(i+1))2722-1/2*log(2)2723sage: P.nonarchimedian_local_height(K.ideal(3))272402725sage: P.nonarchimedian_local_height(K.ideal(1-2*i))2726027272728sage: Q = E.lift_x(-9/4); Q2729(-9/4 : 27/8*i - 4 : 1)2730sage: Q.nonarchimedian_local_height(K.ideal(1+i))27312*log(2)2732sage: Q.nonarchimedian_local_height(K.ideal(3))273302734sage: Q.nonarchimedian_local_height(K.ideal(1-2*i))273502736sage: Q.nonarchimedian_local_height()27371/2*log(16)27382739TESTS::27402741sage: Q.nonarchimedian_local_height(prec=100)27421.38629436111989061883446424292743sage: (3*Q).nonarchimedian_local_height()27441/2*log(75923153929839865104)27452746sage: F.<a> = NumberField(x^4 + 2*x^3 + 19*x^2 + 18*x + 288)2747sage: F.ring_of_integers().gens()2748[1, 5/6*a^3 + 1/6*a, 1/6*a^3 + 1/6*a^2, a^3]2749sage: F.class_number()2750122751sage: E = EllipticCurve('37a').change_ring(F)2752sage: P = E((-a^2/6 - a/6 - 1, a)); P2753(-1/6*a^2 - 1/6*a - 1 : a : 1)2754sage: P[0].is_integral()2755True2756sage: P.nonarchimedian_local_height()275702758"""2759if self.has_finite_order():2760return rings.QQ(0)27612762if prec:2763log = lambda x: rings.RealField(prec)(x).log()2764else:2765from sage.functions.log import log27662767if v is None:2768D = self.curve().discriminant()2769K = self.curve().base_ring()2770factorD = K.factor(D)2771if self[0] == 0:2772c = K.ideal(1)2773else:2774c = K.ideal(self[0]).denominator()2775# The last sum is for bad primes that divide c where the model is not minimal.2776return (log(c.norm())2777+ sum(self.nonarchimedian_local_height(v, prec, True, e < 12) for v,e in factorD if not v.divides(c))2778+ sum(self.nonarchimedian_local_height(v, prec, True) - c.valuation(v) * log(v.norm()) for v,e in factorD if e >= 12 and v.divides(c))2779) / K.degree()27802781if is_minimal:2782E = self.curve()2783P = self2784offset = 02785else:2786E = self.curve().local_minimal_model(v)2787P = self.curve().isomorphism_to(E)(self)2788# Silverman's normalization is not invariant under change of model,2789# but it all cancels out in the global height.2790offset = (self.curve().discriminant()/E.discriminant()).valuation(v)27912792a1, a2, a3, a4, a6 = E.a_invariants()2793b2, b4, b6, b8 = E.b_invariants()2794c4 = E.c4()2795x,y = P.xy()2796D = E.discriminant()2797N = D.valuation(v)2798A = (3*x**2 + 2*a2*x + a4 - a1*y).valuation(v)2799B = (2*y+a1*x+a3).valuation(v)2800C = (3*x**4 + b2*x**3 + 3*b4*x**2 + 3*b6*x + b8).valuation(v)2801if A <= 0 or B <= 0:2802r = max(0, -x.valuation(v))2803elif c4.valuation(v) == 0:2804n = min(B, N/2)2805r = -n*(N-n)/N2806elif C > 3*B:2807r = -2*B/32808else:2809r = -C/42810r -= offset/62811if not r:2812return rings.QQ(0)2813else:2814if E.base_ring() is rings.QQ:2815Nv = rings.ZZ(v)2816else:2817Nv = v.norm()2818if not weighted:2819r /= v.ramification_index() * v.residue_class_degree()2820return r * log(Nv)28212822def elliptic_logarithm(self, embedding=None, precision=100, algorithm='pari'):2823r"""2824Returns the elliptic logarithm of this elliptic curve point.28252826An embedding of the base field into `\RR` or `\CC` (with2827arbitrary precision) may be given; otherwise the first real2828embedding is used (with the specified precision) if any, else2829the first complex embedding.28302831INPUT:28322833- ``embedding``: an embedding of the base field into `\RR` or `\CC`28342835- ``precision``: a positive integer (default 100) setting the2836number of bits of precision for the computation28372838- ``algorithm``: either 'pari' (default for real embeddings)2839to use PARI's ``ellpointtoz{}``, or 'sage' for a native2840implementation. Ignored for complex embeddings.28412842ALGORITHM:28432844See [Co2] Cohen H., A Course in Computational Algebraic Number2845Theory GTM 138, Springer 1996 for the case of real embeddings,2846and Cremona, J.E. and Thongjunthug , T. 2010 for the complex2847case.28482849AUTHORS:28502851- Michael Mardaus (2008-07),2852- Tobias Nagel (2008-07) -- original version from [Co2].2853- John Cremona (2008-07) -- revision following eclib code.2854- John Cremona (2010-03) -- implementation for complex embeddings.28552856EXAMPLES::28572858sage: E = EllipticCurve('389a')2859sage: E.discriminant() > 02860True2861sage: P = E([-1,1])2862sage: P.is_on_identity_component ()2863False2864sage: P.elliptic_logarithm (precision=96)28650.4793482501902193161295330101 + 0.985868850775824102211203849...*I2866sage: Q=E([3,5])2867sage: Q.is_on_identity_component()2868True2869sage: Q.elliptic_logarithm (precision=96)28701.93112827154255944248858522028712872An example with negative discriminant, and a torsion point::28732874sage: E = EllipticCurve('11a1')2875sage: E.discriminant() < 02876True2877sage: P = E([16,-61])2878sage: P.elliptic_logarithm(precision=70)28790.253841860855910684342880sage: E.period_lattice().real_period(prec=70) / P.elliptic_logarithm(precision=70)28815.000000000000000000028822883A larger example. The default algorithm uses PARI and makes2884sure the result has the requested precision::28852886sage: E = EllipticCurve([1, 0, 1, -85357462, 303528987048]) #18074g12887sage: P = E([4458713781401/835903744, -64466909836503771/24167649046528, 1])2888sage: P.elliptic_logarithm() # 100 bits28890.2765620401410706146407620309728902891The native algorithm 'sage' used to have trouble with2892precision in this example, but no longer::28932894sage: P.elliptic_logarithm(algorithm='sage') # 100 bits28950.2765620401410706146407620309728962897This shows that the bug reported at \#4901 has been fixed::28982899sage: E = EllipticCurve("4390c2")2900sage: P = E(683762969925/44944,-565388972095220019/9528128)2901sage: P.elliptic_logarithm()29020.000256387258865202253531989325292903sage: P.elliptic_logarithm(precision=64)29040.0002563872588652022542905sage: P.elliptic_logarithm(precision=65)29060.00025638725886520225352907sage: P.elliptic_logarithm(precision=128)29080.000256387258865202253531989325286664274122909sage: P.elliptic_logarithm(precision=129)29100.000256387258865202253531989325286664274122911sage: P.elliptic_logarithm(precision=256)29120.00025638725886520225353198932528666427411683880083463700150051421280096109363732913sage: P.elliptic_logarithm(precision=257)29140.0002563872588652022535319893252866642741168388008346370015005142128009610936373029152916Examples over number fields::29172918sage: K.<a> = NumberField(x^3-2)2919sage: embs = K.embeddings(CC)2920sage: E = EllipticCurve([0,1,0,a,a])2921sage: Ls = [E.period_lattice(e) for e in embs]2922sage: [L.real_flag for L in Ls]2923[0, 0, -1]2924sage: P = E(-1,0) # order 22925sage: [L.elliptic_logarithm(P) for L in Ls]2926[-1.73964256006716 - 1.07861534489191*I, -0.363756518406398 - 1.50699412135253*I, 1.90726488608927]29272928sage: E = EllipticCurve([-a^2 - a - 1, a^2 + a])2929sage: Ls = [E.period_lattice(e) for e in embs]2930sage: pts = [E(2*a^2 - a - 1 , -2*a^2 - 2*a + 6 ), E(-2/3*a^2 - 1/3 , -4/3*a - 2/3 ), E(5/4*a^2 - 1/2*a , -a^2 - 1/4*a + 9/4 ), E(2*a^2 + 3*a + 4 , -7*a^2 - 10*a - 12 )]2931sage: [[L.elliptic_logarithm(P) for P in pts] for L in Ls]2932[[0.250819591818930 - 0.411963479992219*I, -0.290994550611374 - 1.37239400324105*I, -0.693473752205595 - 2.45028458830342*I, -0.151659609775291 - 1.48985406505459*I], [1.33444787667954 - 1.50889756650544*I, 0.792633734249234 - 0.548467043256610*I, 0.390154532655013 + 0.529423541805758*I, 0.931968675085317 - 0.431006981443071*I], [1.14758249500109 + 0.853389664016075*I, 2.59823462472518 + 0.853389664016075*I, 1.75372176444709, 0.303069634723001]]29332934::29352936sage: K.<i> = QuadraticField(-1)2937sage: E = EllipticCurve([0,0,0,9*i-10,21-i])2938sage: emb = K.embeddings(CC)[1]2939sage: L = E.period_lattice(emb)2940sage: P = E(2-i,4+2*i)2941sage: L.elliptic_logarithm(P,prec=100)29420.70448375537782208460499649302 - 0.79246725643650979858266018068*I29432944"""2945from sage.rings.number_field.number_field import refine_embedding2946from sage.rings.all import RealField, ComplexField, is_RealField, QQ29472948# Check the trivial case:29492950C = ComplexField(precision)2951if self.is_zero():2952return C.zero()29532954# find a suitable embedding if none was supplied:29552956E = self.curve()2957K = E.base_field()2958rational = (K is QQ)2959emb = embedding29602961if emb is None:2962emb = K.embeddings(RealField(precision))2963if len(emb)>0:2964emb = emb[0]2965else:2966emb = K.embeddings(ComplexField(precision))[0]2967else:2968# Get the precision of the supplied embedding2969prec = emb.codomain().precision()2970# if the precision parameter is greater, refine the embedding:2971if precision > prec:2972emb = refine_embedding(emb,precision)29732974L = E.period_lattice(emb)29752976if algorithm == 'sage' or not is_RealField(emb.codomain):2977return L.elliptic_logarithm(self,precision)29782979if algorithm <> 'pari':2980raise ValueError, "algorithm must be either 'pari' or 'sage'"29812982# From now on emb() is a real embedding of K into2983# RealField(precision). We interface with the PARI library.29842985x, y = self.xy()2986if rational: # work with exact coordinates2987E_work = E2988pt_pari = pari([x,y])2989else: # use the embedding to get real coordinates2990ai = [emb(a) for a in E.a_invariants()]2991E_work = EllipticCurve(ai) # defined over RR2992pt_pari = pari([emb(x), emb(y)])2993working_prec = precision2994E_pari = E_work.pari_curve(prec=working_prec)2995log_pari = E_pari.ellpointtoz(pt_pari, precision=working_prec)29962997while prec_words_to_bits(log_pari.precision()) < precision:2998# result is not precise enough, re-compute with double2999# precision. if the base field is not QQ, this3000# requires modifying the precision of the embedding,3001# the curve, and the point3002working_prec = 2*working_prec3003if not rational:3004emb = refine_embedding(emb, working_prec)3005ai = [emb(a) for a in E.a_invariants()]3006E_work = EllipticCurve(ai) # defined over RR3007pt_pari = pari([emb(x), emb(y)])3008E_pari = E_work.pari_curve(prec=working_prec)3009log_pari = E_pari.ellpointtoz(pt_pari, precision=working_prec)30103011# normalization step3012r, i = C(log_pari)3013wR, wI = L.basis(prec=precision)3014k = (r/wR).floor()3015if k:3016r -= k*wR3017if self.is_on_identity_component(emb):3018return C(r)3019# Now there are two components and P is on the non-identity one3020return C(r)+C(wI/2)30213022def padic_elliptic_logarithm(self, p, absprec=20):3023r"""3024Computes the `p`-adic elliptic logarithm of this point.30253026INPUT:30273028``p`` - integer: a prime ``absprec`` - integer (default: 20):3029the initial `p`-adic absolute precision of the computation30303031OUTPUT:30323033The `p`-adic elliptic logarithm of self, with precision ``absprec``.30343035AUTHORS:30363037- Tobias Nagel3038- Michael Mardaus3039- John Cremona30403041ALGORITHM:30423043For points in the formal group (i.e. not integral at `p`) we3044take the ``log()`` function from the formal groups module and3045evaluate it at `-x/y`. Otherwise we first multiply the point3046to get into the formal group, and divide the result3047afterwards.30483049TODO:30503051See comments at trac \#4805. Currently the absolute precision3052of the result may be less than the given value of absprec, and3053error-handling is imperfect.30543055EXAMPLES::30563057sage: E = EllipticCurve([0,1,1,-2,0])3058sage: E(0).padic_elliptic_logarithm(3)305903060sage: P = E(0,0)3061sage: P.padic_elliptic_logarithm(3)30622 + 2*3 + 3^3 + 2*3^7 + 3^8 + 3^9 + 3^11 + 3^15 + 2*3^17 + 3^18 + O(3^19)3063sage: P.padic_elliptic_logarithm(3).lift()30646602575223065sage: P = E(-11/9,28/27)3066sage: [(2*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(20)] # long time (3s)3067[2 + O(2^19), 2 + O(3^20), 2 + O(5^19), 2 + O(7^19), 2 + O(11^19), 2 + O(13^19), 2 + O(17^19), 2 + O(19^19)]3068sage: [(3*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(12)] # long time (2s)3069[1 + 2 + O(2^19), 3 + 3^20 + O(3^21), 3 + O(5^19), 3 + O(7^19), 3 + O(11^19)]3070sage: [(5*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(12)] # long time (2s)3071[1 + 2^2 + O(2^19), 2 + 3 + O(3^20), 5 + O(5^19), 5 + O(7^19), 5 + O(11^19)]30723073An example which arose during reviewing #4741::30743075sage: E = EllipticCurve('794a1')3076sage: P = E(-1,2)3077sage: P.padic_elliptic_logarithm(2) # default precision=2030782^4 + 2^5 + 2^6 + 2^8 + 2^9 + 2^13 + 2^14 + 2^15 + O(2^16)3079sage: P.padic_elliptic_logarithm(2, absprec=30)30802^4 + 2^5 + 2^6 + 2^8 + 2^9 + 2^13 + 2^14 + 2^15 + 2^22 + 2^23 + 2^24 + O(2^26)3081sage: P.padic_elliptic_logarithm(2, absprec=40)30822^4 + 2^5 + 2^6 + 2^8 + 2^9 + 2^13 + 2^14 + 2^15 + 2^22 + 2^23 + 2^24 + 2^28 + 2^29 + 2^31 + 2^34 + O(2^35)3083"""3084if not p.is_prime():3085raise ValueError,'p must be prime'3086debug = False # True3087if debug:3088print "P=",self,"; p=",p," with precision ",precision3089E = self.curve()3090Q_p = Qp(p, absprec)3091if self.has_finite_order():3092return Q_p(0)3093while True:3094try:3095Ep = E.change_ring(Q_p)3096P = Ep(self)3097x,y = P.xy()3098break3099except (PrecisionError, ArithmeticError, ZeroDivisionError):3100absprec *=23101Q_p = Qp(p, absprec)3102if debug:3103print "x,y=",(x,y)3104f = 1 # f will be such that f*P is in the formal group E^1(Q_p)3105if x.valuation() >=0: # P is not in E^13106if not self.has_good_reduction(p): # P is not in E^03107n = E.tamagawa_exponent(p) # n*P has good reduction at p3108if debug:3109print "Tamagawa exponent = =",n3110f = n3111P = n*P # lies in E^03112if debug:3113print "P=",P3114try:3115x,y = P.xy()3116except ZeroDivisionError:3117raise ValueError, "Insufficient precision in p-adic_elliptic_logarithm()"3118if debug:3119print "x,y=",(x,y)3120if x.valuation() >=0: # P is still not in E^13121t = E.local_data(p).bad_reduction_type()3122if t is None:3123m = E.reduction(p).abelian_group().exponent()3124else:3125m = p - t3126if debug:3127print "mod p exponent = =",m3128# now m*(n*P) reduces to the identity mod p, so is in E^1(Q_p)3129f *= m3130P = m*P # lies in E^13131try:3132x,y = P.xy()3133except ZeroDivisionError:3134raise ValueError, "Insufficient precision in p-adic_elliptic_logarithm()"3135if debug:3136print "f=",f3137print "x,y=",(x,y)3138vx = x.valuation()3139vy = y.valuation()3140v = vx-vy3141if not ( v>0 and vx==-2*v and vy == -3*v):3142raise ValueError, "Insufficient precision in p-adic_elliptic_logarithm()"3143try:3144t = -x/y3145except (ZeroDivisionError, PrecisionError):3146raise ValueError, "Insufficient precision in p-adic_elliptic_logarithm()"3147if debug:3148print "t=",t,", with valuation ",v3149phi = Ep.formal().log(prec=1+absprec//v)3150return phi(t)/f315131523153class EllipticCurvePoint_finite_field(EllipticCurvePoint_field):3154r"""3155Class for elliptic curve points over finite fields.3156"""3157def _magma_init_(self, magma):3158"""3159Return a string representation of self that ``MAGMA`` can3160use for input.31613162EXAMPLE::31633164sage: E = EllipticCurve(GF(17), [1,-1])3165sage: P = E([13, 4])3166sage: P._magma_init_(magma) # optional - magma3167'EllipticCurve([_sage_ref...|GF(17)!0,GF(17)!0,GF(17)!0,GF(17)!1,GF(17)!16])![13,4]'3168"""3169E = self.curve()._magma_init_(magma)3170x,y = self.xy()3171return "%s![%s,%s]"%(E,x,y)31723173def discrete_log(self, Q, ord=None):3174r"""3175Returns discrete log of `Q` with respect to `P` =self.31763177INPUT:31783179- ``Q`` (point) -- another point on the same curve as self.31803181- ``ord`` (integer or None (default)) -- the order of self.31823183OUTPUT:31843185(integer) -- The discrete log of `Q` with respect to `P`, which is an3186integer `m` with `0\le m<o(P)` such that `mP=Q`, if one3187exists. A ValueError is raised if there is no solution.31883189.. note::31903191The order of self is computed if not supplied.31923193AUTHOR:31943195- John Cremona. Adapted to use generic functions 2008-04-05.31963197EXAMPLE::31983199sage: F = GF(3^6,'a')3200sage: a = F.gen()3201sage: E = EllipticCurve([0,1,1,a,a])3202sage: E.cardinality()32037623204sage: A = E.abelian_group()3205sage: P = A.gen(0).element()3206sage: Q=400*P3207sage: P.discrete_log(Q)32084003209"""3210if ord==None: ord=self.order()3211try:3212return generic.discrete_log(Q,self,ord,operation='+')3213except:3214raise ValueError, "ECDLog problem has no solution"321532163217def order(self):3218r"""3219Return the order of this point on the elliptic curve.32203221ALGORITHM:32223223Use generic functions from ``sage.groups.generic``. If the3224group order is known, use ``order_from_multiple()``, otherwise3225use ``order_from_bounds()`` with the Hasse bounds for the base3226field. In the latter case, we might find that we have a3227generator for the group, in which case it is cached.32283229We do not cause the group order to be calculated when not3230known, since this function is used in determining the group3231order via computation of several random points and their3232orders.32333234.. note::32353236:meth:`additive_order` is a synonym for :meth:`order`32373238AUTHOR:32393240- John Cremona, 2008-02-10, adapted 2008-04-05 to use generic functions.32413242EXAMPLES::32433244sage: k.<a> = GF(5^5)3245sage: E = EllipticCurve(k,[2,4]); E3246Elliptic Curve defined by y^2 = x^3 + 2*x + 4 over Finite Field in a of size 5^53247sage: P = E(3*a^4 + 3*a , 2*a + 1 )3248sage: P.order()324932273250sage: Q = E(0,2)3251sage: Q.order()325273253sage: Q.additive_order()3254732553256In the next example, the cardinality of E will be computed (using SEA) and cached::32573258sage: p=next_prime(2^150)3259sage: E=EllipticCurve(GF(p),[1,1])3260sage: P=E(831623307675610677632782670796608848711856078, 42295786042873366706573292533588638217232964)3261sage: P.order()326214272476927059598810582625452724743006282814483263sage: P.order()==E.cardinality()3264True326532663267"""3268try:3269return self._order3270except AttributeError:3271pass3272if self.is_zero():3273return rings.Integer(1)3274E = self.curve()3275K = E.base_ring()3276from sage.schemes.plane_curves.projective_curve import Hasse_bounds3277bounds = Hasse_bounds(K.order())32783279try:3280M = E._order3281try:3282plist = E._prime_factors_of_order3283except:3284plist = M.prime_divisors()3285E._prime_factors_of_order = plist3286N = generic.order_from_multiple(self,M,plist,operation='+')3287except:3288if K.is_prime_field():3289M = E.cardinality() # computed and cached3290plist = M.prime_divisors()3291E._prime_factors_of_order = plist3292N = generic.order_from_multiple(self,M,plist,operation='+')3293else:3294N = generic.order_from_bounds(self,bounds,operation='+')32953296if 2*N>bounds[1]: # then we have a generator, so cache this3297try:3298dummy = E._order3299except:3300E._order = N3301try:3302dummy = E.__abelian_group3303except:3304E.__abelian_group = AbelianGroup([N]), (self,)33053306self._order = N3307return self._order33083309additive_order = order33103311331233133314