Path: blob/master/src/sage/schemes/elliptic_curves/ell_point.py
8821 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.1112The class ``EllipticCurvePoint``, which is based on13``SchemeMorphism_point_projective_ring``, currently has little extra14functionality.1516EXAMPLES:1718An example over `\QQ`::1920sage: E = EllipticCurve('389a1')21sage: P = E(-1,1); P22(-1 : 1 : 1)23sage: Q = E(0,-1); Q24(0 : -1 : 1)25sage: P+Q26(4 : 8 : 1)27sage: P-Q28(1 : 0 : 1)29sage: 3*P-5*Q30(328/361 : -2800/6859 : 1)3132An example over a number field::3334sage: K.<i> = QuadraticField(-1)35sage: E = EllipticCurve(K,[1,0,0,0,-1])36sage: P = E(0,i); P37(0 : i : 1)38sage: P.order()39+Infinity40sage: 101*P-100*P==P41True4243An example over a finite field::4445sage: K.<a> = GF(101^3)46sage: E = EllipticCurve(K,[1,0,0,0,-1])47sage: P = E(40*a^2 + 69*a + 84 , 58*a^2 + 73*a + 45)48sage: P.order()49103221050sage: E.cardinality()5110322105253Arithmetic with a point over an extension of a finite field::5455sage: k.<a> = GF(5^2)56sage: E = EllipticCurve(k,[1,0]); E57Elliptic Curve defined by y^2 = x^3 + x over Finite Field in a of size 5^258sage: P = E([a,2*a+4])59sage: 5*P60(2*a + 3 : 2*a : 1)61sage: P*562(2*a + 3 : 2*a : 1)63sage: P + P + P + P + P64(2*a + 3 : 2*a : 1)6566::6768sage: F = Zmod(3)69sage: E = EllipticCurve(F,[1,0]);70sage: P = E([2,1])71sage: import sys72sage: n = sys.maxint73sage: P*(n+1)-P*n == P74True7576Arithmetic over `\ZZ/N\ZZ` with composite `N` is supported. When an77operation tries to invert a non-invertible element, a78ZeroDivisionError is raised and a factorization of the modulus appears79in the error message::8081sage: N = 171576151382sage: E = EllipticCurve(Integers(N),[3,-13])83sage: P = E(2,1)84sage: LCM([2..60])*P85Traceback (most recent call last):86...87ZeroDivisionError: Inverse of 1520944668 does not exist (characteristic = 1715761513 = 26927*63719)888990AUTHORS:9192- William Stein (2005) -- Initial version9394- Robert Bradshaw et al....9596- John Cremona (Feb 2008) -- Point counting and group structure for97non-prime fields, Frobenius endomorphism and order, elliptic logs9899- John Cremona (Aug 2008) -- Introduced ``EllipticCurvePoint_number_field`` class100101- Tobias Nagel, Michael Mardaus, John Cremona (Dec 2008) -- `p`-adic elliptic logarithm over `\QQ`102103- David Hansen (Jan 2009) -- Added ``weil_pairing`` function to ``EllipticCurvePoint_finite_field`` class104105- Mariah Lenox (March 2011) -- Added ``tate_pairing`` and ``ate_pairing``106functions to ``EllipticCurvePoint_finite_field`` class107108"""109110#*****************************************************************************111# Copyright (C) 2005 William Stein <[email protected]>112#113# Distributed under the terms of the GNU General Public License (GPL)114# as published by the Free Software Foundation; either version 2 of115# the License, or (at your option) any later version.116# http://www.gnu.org/licenses/117#*****************************************************************************118119120import math121122import sage.plot.all as plot123124from sage.rings.padics.factory import Qp125from sage.rings.padics.precision_error import PrecisionError126127import sage.rings.all as rings128from sage.rings.real_mpfr import is_RealField129from sage.rings.all import ZZ130from sage.groups.all import AbelianGroup131import sage.groups.generic as generic132from sage.libs.pari.pari_instance import pari, prec_words_to_bits133from sage.structure.sequence import Sequence134135from sage.schemes.plane_curves.projective_curve import Hasse_bounds136from sage.schemes.projective.projective_point import (SchemeMorphism_point_projective_ring,137SchemeMorphism_point_abelian_variety_field)138from sage.schemes.generic.morphism import is_SchemeMorphism139140from constructor import EllipticCurve141from sage.misc.superseded import deprecated_function_alias142143oo = rings.infinity # infinity144145146class 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.154155.. NOTE::156157``__eq__`` and ``__ne__`` are implemented in158SchemeMorphism_point_projective_ring159160EXAMPLES:161sage: E=EllipticCurve(QQ,[1,1])162sage: P=E(0,1)163sage: P.order()164+Infinity165sage: Q=P+P166sage: P==Q167False168sage: Q+Q == 4*P169True170"""171assert isinstance(other, (int, long, rings.Integer)) and other == 0172if self.is_zero():173return 0174else:175return -1176177178class EllipticCurvePoint_field(SchemeMorphism_point_abelian_variety_field):179"""180A point on an elliptic curve over a field. The point has coordinates181in the base field.182183EXAMPLES::184185sage: E = EllipticCurve('37a')186sage: E([0,0])187(0 : 0 : 1)188sage: E(0,0) # brackets are optional189(0 : 0 : 1)190sage: E([GF(5)(0), 0]) # entries are coerced191(0 : 0 : 1)192193sage: E(0.000, 0)194(0 : 0 : 1)195196sage: E(1,0,0)197Traceback (most recent call last):198...199TypeError: Coordinates [1, 0, 0] do not define a point on200Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field201202::203204sage: E = EllipticCurve([0,0,1,-1,0])205sage: S = E(QQ); S206Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field207208sage: K.<i>=NumberField(x^2+1)209sage: E=EllipticCurve(K,[0,1,0,-160,308])210sage: P=E(26,-120)211sage: Q=E(2+12*i,-36+48*i)212sage: P.order() == Q.order() == 4 # long time (3s)213True214sage: 2*P==2*Q215False216217::218219sage: K.<t>=FractionField(PolynomialRing(QQ,'t'))220sage: E=EllipticCurve([0,0,0,0,t^2])221sage: P=E(0,t)222sage: P,2*P,3*P223((0 : t : 1), (0 : -t : 1), (0 : 1 : 0))224225226TESTS::227228sage: loads(S.dumps()) == S229True230sage: E = EllipticCurve('37a')231sage: P = E(0,0); P232(0 : 0 : 1)233sage: loads(P.dumps()) == P234True235sage: T = 100*P236sage: loads(T.dumps()) == T237True238239Test pickling an elliptic curve that has known points on it::240241sage: e = EllipticCurve([0, 0, 1, -1, 0]); g = e.gens(); loads(dumps(e)) == e242True243"""244def __init__(self, curve, v, check=True):245"""246Constructor for a point on an elliptic curve.247248INPUT:249250- curve -- an elliptic curve251- v -- data determining a point (another point, the integer2520, or a tuple of coordinates)253254EXAMPLE::255256sage: E = EllipticCurve('43a')257sage: P = E([2, -4, 2]); P258(1 : -2 : 1)259sage: P == E([1,-2])260True261sage: P = E(0); P262(0 : 1 : 0)263sage: P=E(2, -4, 2); P264(1 : -2 : 1)265"""266point_homset = curve.point_homset()267if is_SchemeMorphism(v) or isinstance(v, EllipticCurvePoint_field):268v = list(v)269elif v == 0:270# some of the code assumes that E(0) has integral entries271# irregardless of the base ring...272#R = self.base_ring()273#v = (R.zero(),R.one(),R.zero())274v = (0, 1, 0)275if check:276# mostly from SchemeMorphism_point_projective_field277d = point_homset.codomain().ambient_space().ngens()278if not isinstance(v, (list, tuple)):279raise TypeError("Argument v (= %s) must be a scheme point, list, or tuple." % str(v))280if len(v) != d and len(v) != d-1:281raise TypeError("v (=%s) must have %s components" % (v, d))282v = Sequence(v, point_homset.value_ring())283if len(v) == d-1: # very common special case284v.append(v.universe()(1))285286n = len(v)287all_zero = True288for i in range(n):289c = v[n-1-i]290if c:291all_zero = False292if c == 1:293break294for j in range(n-i):295v[j] /= c296break297if all_zero:298raise ValueError("%s does not define a valid point "299"since all entries are 0" % repr(v))300301x, y, z = v302if z == 0:303test = x304else:305a1, a2, a3, a4, a6 = curve.ainvs()306test = y**2 + (a1*x+a3)*y - (((x+a2)*x+a4)*x+a6)307if not test == 0:308raise TypeError("Coordinates %s do not define a point on %s" % (list(v), curve))309310SchemeMorphism_point_abelian_variety_field.__init__(self, point_homset, v, check=False)311#AdditiveGroupElement.__init__(self, point_homset)312313def _repr_(self):314"""315Return a string representation of this point.316317EXAMPLE::318319sage: E = EllipticCurve('39a')320sage: P = E([-2, 1, 1])321sage: P._repr_()322'(-2 : 1 : 1)'323"""324return self.codomain().ambient_space()._repr_generic_point(self._coords)325326def _latex_(self):327"""328Return a LaTeX representation of this point.329330EXAMPLE::331332sage: E = EllipticCurve('40a')333sage: P = E([3, 0])334sage: P._latex_()335'\\left(3 : 0 : 1\\right)'336"""337return self.codomain().ambient_space()._latex_generic_point(self._coords)338339def __getitem__(self, n):340"""341Return the n'th coordinate of this point.342343EXAMPLE::344345sage: E = EllipticCurve('42a')346sage: P = E([-17, -51, 17])347sage: [P[i] for i in [2,1,0]]348[1, -3, -1]349"""350return self._coords[n]351352def __iter__(self):353"""354Return the coordinates of this point as a list.355356EXAMPLE::357358sage: E = EllipticCurve('37a')359sage: list(E([0,0]))360[0, 0, 1]361"""362return iter(self._coords)363364def __tuple__(self):365"""366Return the coordinates of this point as a tuple.367368EXAMPLE::369370sage: E = EllipticCurve('44a')371sage: P = E([1, -2, 1])372sage: P.__tuple__()373(1, -2, 1)374"""375return tuple(self._coords) # Warning: _coords is a list!376377def __cmp__(self, other):378"""379Comparison function for points to allow sorting and equality testing.380381.. NOTE::382383``__eq__`` and ``__ne__`` are implemented in384SchemeMorphism_point_projective_field385386EXAMPLES::387388sage: E = EllipticCurve('45a')389sage: P = E([2, -1, 1])390sage: P == E(0)391False392sage: P+P == E(0)393True394"""395if not isinstance(other, EllipticCurvePoint_field):396try:397other = self.codomain().ambient_space()(other)398except TypeError:399return -1400return cmp(self._coords, other._coords)401402def _pari_(self):403r"""404Converts this point to PARI format.405406EXAMPLES::407408sage: E = EllipticCurve([0,0,0,3,0])409sage: O = E(0)410sage: P = E.point([1,2])411sage: O._pari_()412[0]413sage: P._pari_()414[1, 2]415416The following implicitly calls O._pari_() and P._pari_()::417418sage: pari(E).elladd(O,P)419[1, 2]420421TESTS::422423Try the same over a finite field::424425sage: E = EllipticCurve(GF(11), [0,0,0,3,0])426sage: O = E(0)427sage: P = E.point([1,2])428sage: O._pari_()429[0]430sage: P._pari_()431[Mod(1, 11), Mod(2, 11)]432433We no longer need to explicitly call ``pari(O)`` and ``pari(P)``434after :trac:`11868`::435436sage: pari(E).elladd(O, P)437[Mod(1, 11), Mod(2, 11)]438"""439if self[2]:440return pari([self[0]/self[2], self[1]/self[2]])441else:442return pari([0])443444def scheme(self):445"""446Return the scheme of this point, i.e., the curve it is on.447This is synonymous with :meth:`curve` which is perhaps more448intuitive.449450EXAMPLES::451452sage: E=EllipticCurve(QQ,[1,1])453sage: P=E(0,1)454sage: P.scheme()455Elliptic Curve defined by y^2 = x^3 + x + 1 over Rational Field456sage: P.scheme() == P.curve()457True458sage: K.<a>=NumberField(x^2-3,'a')459sage: P=E.base_extend(K)(1,a)460sage: P.scheme()461Elliptic Curve defined by y^2 = x^3 + x + 1 over Number Field in a with defining polynomial x^2 - 3462"""463#The following text is just not true: it applies to the class464#EllipticCurvePoint, which appears to be never used, but does465#not apply to EllipticCurvePoint_field which is simply derived466#from AdditiveGroupElement.467#468#"Technically, points on curves in Sage are scheme maps from469# the domain Spec(F) where F is the base field of the curve to470# the codomain which is the curve. See also domain() and471# codomain()."472473return self.codomain()474475def domain(self):476"""477Return the domain of this point, which is `Spec(F)` where `F` is478the field of definition.479480EXAMPLES::481482sage: E=EllipticCurve(QQ,[1,1])483sage: P=E(0,1)484sage: P.domain()485Spectrum of Rational Field486sage: K.<a>=NumberField(x^2-3,'a')487sage: P=E.base_extend(K)(1,a)488sage: P.domain()489Spectrum of Number Field in a with defining polynomial x^2 - 3490"""491return self.parent().domain()492493def codomain(self):494"""495Return the codomain of this point, which is the curve it is496on. Synonymous with :meth:`curve` which is perhaps more intuitive.497498EXAMPLES::499500sage: E=EllipticCurve(QQ,[1,1])501sage: P=E(0,1)502sage: P.domain()503Spectrum of Rational Field504sage: K.<a>=NumberField(x^2-3,'a')505sage: P=E.base_extend(K)(1,a)506sage: P.codomain()507Elliptic Curve defined by y^2 = x^3 + x + 1 over Number Field in a with defining polynomial x^2 - 3508sage: P.codomain() == P.curve()509True510"""511return self.parent().codomain()512513def order(self):514r"""515Return the order of this point on the elliptic curve.516517If the point is zero, returns 1, otherwise raise a518NotImplementedError.519520For curves over number fields and finite fields, see below.521522.. NOTE::523524:meth:`additive_order` is a synonym for :meth:`order`525526EXAMPLE::527528sage: K.<t>=FractionField(PolynomialRing(QQ,'t'))529sage: E=EllipticCurve([0,0,0,-t^2,0])530sage: P=E(t,0)531sage: P.order()532Traceback (most recent call last):533...534NotImplementedError: Computation of order of a point not implemented over general fields.535sage: E(0).additive_order()5361537sage: E(0).order() == 1538True539540"""541if self.is_zero():542return rings.Integer(1)543raise NotImplementedError("Computation of order of a point "544"not implemented over general fields.")545546additive_order = order547548def curve(self):549"""550Return the curve that this point is on.551552EXAMPLES::553554sage: E = EllipticCurve('389a')555sage: P = E([-1,1])556sage: P.curve()557Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field558"""559return self.scheme()560561def __nonzero__(self):562"""563Return True if this is not the zero point on the curve.564565EXAMPLES::566567sage: E = EllipticCurve('37a')568sage: P = E(0); P569(0 : 1 : 0)570sage: P.is_zero()571True572sage: P = E.gens()[0]573sage: P.is_zero()574False575"""576return bool(self[2])577578def has_finite_order(self):579"""580Return True if this point has finite additive order as an element581of the group of points on this curve.582583For fields other than number fields and finite fields, this is584NotImplemented unless self.is_zero().585586EXAMPLES::587588sage: K.<t>=FractionField(PolynomialRing(QQ,'t'))589sage: E=EllipticCurve([0,0,0,-t^2,0])590sage: P = E(0)591sage: P.has_finite_order()592True593sage: P=E(t,0)594sage: P.has_finite_order()595Traceback (most recent call last):596...597NotImplementedError: Computation of order of a point not implemented over general fields.598sage: (2*P).is_zero()599True600"""601if self.is_zero():602return True603return self.order() != oo604605is_finite_order = has_finite_order # for backward compatibility606607def has_infinite_order(self):608"""609Return True if this point has infinite additive order as an element610of the group of points on this curve.611612For fields other than number fields and finite fields, this is613NotImplemented unless self.is_zero().614615EXAMPLES::616617sage: K.<t>=FractionField(PolynomialRing(QQ,'t'))618sage: E=EllipticCurve([0,0,0,-t^2,0])619sage: P = E(0)620sage: P.has_infinite_order()621False622sage: P=E(t,0)623sage: P.has_infinite_order()624Traceback (most recent call last):625...626NotImplementedError: Computation of order of a point not implemented over general fields.627sage: (2*P).is_zero()628True629"""630if self.is_zero():631return False632return self.order() == oo633634def plot(self, **args):635"""636Plot this point on an elliptic curve.637638INPUT:639640- ``**args`` -- all arguments get passed directly onto the point641plotting function.642643EXAMPLES::644645sage: E = EllipticCurve('389a')646sage: P = E([-1,1])647sage: P.plot(pointsize=30, rgbcolor=(1,0,0))648"""649if self.is_zero():650return plot.text("$\\infty$", (-3, 3), **args)651652else:653return plot.point((self[0], self[1]), **args)654655def _add_(self, right):656"""657Add self to right.658659EXAMPLES::660661sage: E = EllipticCurve('389a')662sage: P = E([-1,1]); Q = E([0,0])663sage: P + Q664(1 : 0 : 1)665sage: P._add_(Q) == P + Q666True667668Example to show that bug :trac:`4820` is fixed::669670sage: [type(c) for c in 2*EllipticCurve('37a1').gen(0)]671[<type 'sage.rings.rational.Rational'>,672<type 'sage.rings.rational.Rational'>,673<type 'sage.rings.rational.Rational'>]674"""675# Use Prop 7.1.7 of Cohen "A Course in Computational Algebraic676# Number Theory"677if self.is_zero():678return right679if right.is_zero():680return self681E = self.curve()682a1, a2, a3, a4, a6 = E.ainvs()683x1, y1 = self[0], self[1]684x2, y2 = right[0], right[1]685if x1 == x2 and y1 == -y2 - a1*x2 - a3:686return E(0) # point at infinity687688if x1 == x2 and y1 == y2:689try:690m = (3*x1*x1 + 2*a2*x1 + a4 - a1*y1) / (2*y1 + a1*x1 + a3)691except ZeroDivisionError:692R = E.base_ring()693if R.is_finite():694N = R.characteristic()695N1 = N.gcd(ZZ(2*y1 + a1*x1 + a3))696N2 = N//N1697raise ZeroDivisionError("Inverse of %s does not exist (characteristic = %s = %s*%s)" % (2*y1 + a1*x1 + a3, N, N1, N2))698else:699raise ZeroDivisionError("Inverse of %s does not exist" % (2*y1 + a1*x1 + a3))700else:701try:702m = (y1-y2)/(x1-x2)703except ZeroDivisionError:704R = E.base_ring()705if R.is_finite():706N = R.characteristic()707from sage.rings.all import ZZ708N1 = N.gcd(ZZ(x1-x2))709N2 = N//N1710raise ZeroDivisionError("Inverse of %s does not exist (characteristic = %s = %s*%s)" % (x1-x2, N, N1, N2))711else:712raise ZeroDivisionError("Inverse of %s does not exist" % (x1-x2))713714x3 = -x1 - x2 - a2 + m*(m+a1)715y3 = -y1 - a3 - a1*x3 + m*(x1-x3)716# See trac #4820 for why we need to coerce 1 into the base ring here:717return E.point([x3, y3, E.base_ring()(1)], check=False)718719def _sub_(self, right):720"""721Subtract right from self.722723EXAMPLES::724725sage: E = EllipticCurve('389a')726sage: P = E([-1,1]); Q = E([0,0])727sage: P - Q728(4 : 8 : 1)729sage: P - Q == P._sub_(Q)730True731sage: (P - Q) + Q732(-1 : 1 : 1)733sage: P734(-1 : 1 : 1)735"""736return self + (-right)737738def __neg__(self):739"""740Return the additive inverse of this point.741742EXAMPLES::743744sage: E = EllipticCurve('389a')745sage: P = E([-1,1])746sage: Q = -P; Q747(-1 : -2 : 1)748sage: Q + P749(0 : 1 : 0)750751Example to show that bug :trac:`4820` is fixed::752753sage: [type(c) for c in -EllipticCurve('37a1').gen(0)]754[<type 'sage.rings.rational.Rational'>,755<type 'sage.rings.rational.Rational'>,756<type 'sage.rings.rational.Rational'>]757"""758if self.is_zero():759return self760E, x, y = self.curve(), self[0], self[1]761# See trac #4820 for why we need to coerce 1 into the base ring here:762return E.point([x, -y - E.a1()*x - E.a3(), E.base_ring()(1)], check=False)763764def xy(self):765"""766Return the `x` and `y` coordinates of this point, as a 2-tuple.767If this is the point at infinity a ZeroDivisionError is raised.768769EXAMPLES::770771sage: E = EllipticCurve('389a')772sage: P = E([-1,1])773sage: P.xy()774(-1, 1)775sage: Q = E(0); Q776(0 : 1 : 0)777sage: Q.xy()778Traceback (most recent call last):779...780ZeroDivisionError: Rational division by zero781"""782if self[2] == 1:783return self[0], self[1]784else:785return self[0]/self[2], self[1]/self[2]786787def is_divisible_by(self, m):788"""789Return True if there exists a point `Q` defined over the same790field as self such that `mQ` == self.791792INPUT:793794- ``m`` -- a positive integer.795796OUTPUT:797798(bool) -- True if there is a solution, else False.799800.. WARNING::801802This function usually triggers the computation of the803`m`-th division polynomial of the associated elliptic804curve, which will be expensive if `m` is large, though it805will be cached for subsequent calls with the same `m`.806807EXAMPLES::808809sage: E = EllipticCurve('389a')810sage: Q = 5*E(0,0); Q811(-2739/1444 : -77033/54872 : 1)812sage: Q.is_divisible_by(4)813False814sage: Q.is_divisible_by(5)815True816817A finite field example::818819sage: E = EllipticCurve(GF(101),[23,34])820sage: E.cardinality().factor()8212 * 53822sage: Set([T.order() for T in E.points()])823{1, 106, 2, 53}824sage: len([T for T in E.points() if T.is_divisible_by(2)])82553826sage: len([T for T in E.points() if T.is_divisible_by(3)])827106828829TESTS:830831This shows that the bug reported at :trac:`10076` is fixed::832833sage: K = QuadraticField(8,'a')834sage: E = EllipticCurve([K(0),0,0,-1,0])835sage: P = E([-1,0])836sage: P.is_divisible_by(2)837False838sage: P.division_points(2)839[]840841Note that it is not sufficient to test that842``self.division_points(m,poly_only=True)`` has roots::843844sage: P.division_points(2, poly_only=True).roots()845[(1/2*a - 1, 1), (-1/2*a - 1, 1)]846847sage: tor = E.torsion_points(); len(tor)8488849sage: [T.order() for T in tor]850[2, 4, 4, 2, 1, 2, 4, 4]851sage: all([T.is_divisible_by(3) for T in tor])852True853sage: Set([T for T in tor if T.is_divisible_by(2)])854{(0 : 1 : 0), (1 : 0 : 1)}855sage: Set([2*T for T in tor])856{(0 : 1 : 0), (1 : 0 : 1)}857858859"""860# Coerce the input m to an integer861m = rings.Integer(m)862863# Check for trivial cases of m = 1, -1 and 0.864if m == 1 or m == -1:865return True866if m == 0:867return self == 0 # then m*self=self for all m!868m = m.abs()869870# Now the following line would of course be correct, but we871# work harder to be more efficient:872# return len(self.division_points(m)) > 0873874P = self875876# If P has finite order n and gcd(m,n)=1 then the result is877# True. However, over general fields computing P.order() is878# not implemented.879880try:881n = P.order()882if not n == oo:883if m.gcd(n) == 1:884return True885except NotImplementedError:886pass887888P_is_2_torsion = (P == -P)889g = P.division_points(m, poly_only=True)890891if not P_is_2_torsion:892# In this case deg(g)=m^2, and each root in K lifts to two893# points Q,-Q both in E(K), of which exactly one is a894# solution. So we just check the existence of roots:895return len(g.roots()) > 0896897# Now 2*P==0898899if m % 2 == 1:900return True # P itself is a solution when m is odd901902# Now m is even and 2*P=0. Roots of g in K may or may not903# lift to solutions in E(K), so we fall back to the default.904# Note that division polynomials are cached so this is not905# inefficient:906907return len(self.division_points(m)) > 0908909def division_points(self, m, poly_only=False):910r"""911Return a list of all points `Q` such that `mQ=P` where `P` = self.912913Only points on the elliptic curve containing self and defined914over the base field are included.915916INPUT:917918- ``m`` -- a positive integer919920- ``poly_only`` -- bool (default: False); if True return921polynomial whose roots give all possible `x`-coordinates of922`m`-th roots of self.923924OUTPUT:925926(list) -- a (possibly empty) list of solutions `Q` to `mQ=P`,927where `P` = self.928929EXAMPLES:930931We find the five 5-torsion points on an elliptic curve::932933sage: E = EllipticCurve('11a'); E934Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field935sage: P = E(0); P936(0 : 1 : 0)937sage: P.division_points(5)938[(0 : 1 : 0), (5 : -6 : 1), (5 : 5 : 1), (16 : -61 : 1), (16 : 60 : 1)]939940Note above that 0 is included since [5]*0 = 0.941942We create a curve of rank 1 with no torsion and do a consistency check::943944sage: E = EllipticCurve('11a').quadratic_twist(-7)945sage: Q = E([44,-270])946sage: (4*Q).division_points(4)947[(44 : -270 : 1)]948949We create a curve over a non-prime finite field with group of950order `18`::951952sage: k.<a> = GF(25)953sage: E = EllipticCurve(k, [1,2+a,3,4*a,2])954sage: P = E([3,3*a+4])955sage: factor(E.order())9562 * 3^2957sage: P.order()9589959960We find the `1`-division points as a consistency check -- there961is just one, of course::962963sage: P.division_points(1)964[(3 : 3*a + 4 : 1)]965966The point `P` has order coprime to 2 but divisible by 3, so::967968sage: P.division_points(2)969[(2*a + 1 : 3*a + 4 : 1), (3*a + 1 : a : 1)]970971We check that each of the 2-division points works as claimed::972973sage: [2*Q for Q in P.division_points(2)]974[(3 : 3*a + 4 : 1), (3 : 3*a + 4 : 1)]975976Some other checks::977978sage: P.division_points(3)979[]980sage: P.division_points(4)981[(0 : 3*a + 2 : 1), (1 : 0 : 1)]982sage: P.division_points(5)983[(1 : 1 : 1)]984985An example over a number field (see :trac:`3383`)::986987sage: E = EllipticCurve('19a1')988sage: 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)989sage: EK = E.base_extend(K)990sage: E(0).division_points(3)991[(0 : 1 : 0), (5 : -10 : 1), (5 : 9 : 1)]992sage: EK(0).division_points(3)993[(0 : 1 : 0), (5 : 9 : 1), (5 : -10 : 1)]994sage: E(0).division_points(9)995[(0 : 1 : 0), (5 : -10 : 1), (5 : 9 : 1)]996sage: EK(0).division_points(9)997[(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)]998999"""1000# Coerce the input m to an integer1001m = rings.Integer(m)1002# Check for trivial cases of m = 1, -1 and 0.1003if m == 1 or m == -1:1004return [self]1005if m == 0:1006if self == 0: # then every point Q is a solution, but...1007return [self]1008else:1009return []10101011# ans will contain the list of division points.1012ans = []10131014# We compute a polynomial g whose roots give all possible x1015# coordinates of the m-division points. The number of1016# solutions (over the algebraic closure) is m^2, assuming that1017# the characteristic does not divide m.10181019E = self.curve()1020P = self1021nP = -P1022P_is_2_torsion = (P == nP)10231024# If self is the 0, then self is a solution, and the correct1025# poly is the m'th division polynomial1026if P == 0:1027ans.append(P)1028g = E.division_polynomial(m)1029else:1030# The poly g here is 0 at x(Q) iff x(m*Q) = x(P).1031g = E._multiple_x_numerator(m) - P[0]*E._multiple_x_denominator(m)10321033# When 2*P=0, then -Q is a solution iff Q is. For even m,1034# no 2-torsion point is a solution, so that g is the1035# square of a polynomial g1 of degree m^2/2, and each root1036# of g1 leads to a pair of solutions Q, -Q to m*Q=P. For1037# odd m, P itself is the only 2-torsion solution, so g has1038# the form (x-x(P))*g1(x)^2 where g1 has degree (m^2-1)/21039# and each root of g1 leads to a pair Q, -Q.10401041if P_is_2_torsion:1042if m % 2 == 0:1043# This computes g.sqrt() which is not implemented1044g = g.gcd(g.derivative())*g.leading_coefficient().sqrt()10451046# When 2*P!=0, then for each solution Q to m*Q=P, -Q is1047# not a solution (and points of order 2 are not1048# solutions). Hence the roots of g are distinct and each1049# gives rise to precisely one solution Q.10501051else:1052g0 = g.variables()[0] - P[0]1053g = g // g01054g = g.gcd(g.derivative())*g.leading_coefficient().sqrt()1055g = g0*g10561057if poly_only:1058return g10591060for x in g.roots(multiplicities=False):1061if E.is_x_coord(x):1062# Make a point on the curve with this x coordinate.1063Q = E.lift_x(x)1064nQ = -Q1065mQ = m*Q1066# if P==-P then Q works iff -Q works, so we include1067# both unless they are equal:1068if P_is_2_torsion:1069if mQ == P:1070ans.append(Q)1071if nQ != Q:1072ans.append(nQ)1073else:1074# P is not 2-torsion so at most one of Q, -Q works1075# and we must try both:1076if mQ == P:1077ans.append(Q)1078elif mQ == nP:1079ans.append(nQ)10801081# Finally, sort and return1082ans.sort()1083return ans10841085def _divide_out(self, p):1086r"""1087Return `(Q,k)` where `p^kQ` == self and `Q` cannot be divided by `p`.10881089..WARNING:10901091It is up to the caller to make sure that this does not loop1092endlessly. It is used in1093``EllipticCurve_generic._p_primary_torsion_basis()``, when1094self will always have (finite) order which is a power of `p`,1095so that the order of `Q` increases by a factor of `p` at each1096stage.10971098Since it will clearly be in danger of looping when1099self.is_zero(), this case is caught, but otherwise caveat1100user.11011102EXAMPLES::11031104sage: E = EllipticCurve('37a1')1105sage: P = E([0, 0])1106sage: R = 12*P1107sage: R._divide_out(2)1108((-1 : -1 : 1), 2)1109sage: R._divide_out(3)1110((2 : -3 : 1), 1)1111sage: R._divide_out(5)1112((1357/841 : 28888/24389 : 1), 0)1113sage: R._divide_out(12)1114Traceback (most recent call last):1115...1116ValueError: p (=12) should be prime.1117"""1118p = rings.Integer(p)1119if not p.is_prime():1120raise ValueError("p (=%s) should be prime." % p)11211122if self.is_zero():1123raise ValueError("self must not be 0.")11241125k = 01126Q = self1127pts = Q.division_points(p)1128while len(pts) > 0:1129Q = pts[0]1130k += 11131pts = Q.division_points(p)1132return (Q, k)11331134def set_order(self, value):1135r"""1136Set the value of self._order to value.11371138Use this when you know a priori the order of this point to avoid a1139potentially expensive order calculation.11401141INPUT:11421143- ``value`` - positive Integer11441145OUTPUT:11461147``None``11481149EXAMPLES:11501151This example illustrates basic usage.11521153::11541155sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61156sage: G = E(5, 0)1157sage: G.set_order(2)1158sage: 2*G1159(0 : 1 : 0)11601161We now give a more interesting case, the NIST-P521 curve. Its1162order is too big to calculate with SAGE, and takes a long time1163using other packages, so it is very useful here.11641165::11661167sage: p = 2^521 - 11168sage: prev_proof_state = proof.arithmetic()1169sage: proof.arithmetic(False) # turn off primality checking1170sage: F = GF(p)1171sage: A = p - 31172sage: B = 10938490380737342745111123907668055699362075989516837489945863944959531161507350160137087375737596232485921322967063133094384525315910129121423274884789859841173sage: q = 68647976601306097149819007990813932172694353001433054093944634591855431833976553942450577463332171975329639963713633211138647686124403803403728088927070054491174sage: E = EllipticCurve([F(A), F(B)])1175sage: G = E.random_point()1176sage: G.set_order(q)1177sage: G.order() * G # This takes practically no time.1178(0 : 1 : 0)1179sage: proof.arithmetic(prev_proof_state) # restore state11801181It is an error to pass a `value` equal to `0`::11821183sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61184sage: G = E.random_point()1185sage: G.set_order(0)1186Traceback (most recent call last):1187...1188ValueError: Value 0 illegal for point order1189sage: G.set_order(1000)1190Traceback (most recent call last):1191...1192ValueError: Value 1000 illegal: outside max Hasse bound11931194It is also very likely an error to pass a value which is not the actual1195order of this point. How unlikely is determined by the factorization of1196the actual order, and the actual group structure::11971198sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61199sage: G = E(5, 0) # G has order 21200sage: G.set_order(11)1201Traceback (most recent call last):1202...1203ValueError: Value 11 illegal: 11 * (5 : 0 : 1) is not the identity12041205However, set_order can be fooled, though it's not likely in "real cases1206of interest". For instance, the order can be set to a multiple the1207actual order::12081209sage: E = EllipticCurve(GF(7), [0, 1]) # This curve has order 61210sage: G = E(5, 0) # G has order 21211sage: G.set_order(8)1212sage: G.order()1213812141215NOTES:12161217The implementation is based of the fact that orders of elliptic curve1218points are cached in the (pseudo-private) _order slot.12191220AUTHORS:12211222- Mariah Lenox (2011-02-16)1223"""1224E = self.curve()1225q = E.base_ring().order()1226O = E(0)1227if value == 0:1228raise ValueError('Value 0 illegal for point order')1229if value == 1:1230if self == O:1231self._order = 11232return1233else:1234raise ValueError('Value 1 illegal order for non-identity')1235low, hi = Hasse_bounds(q)1236if value > hi:1237raise ValueError('Value %s illegal: outside max Hasse bound' % value)1238if value * self != O:1239raise ValueError('Value %s illegal: %s * %s is not the identity' % (value, value, self))1240if (value - 1) * self == O:1241raise ValueError('Value %s illegal: %s * %s is the identity' % (value, value-1, self))1242self._order = value12431244############################## end ################################12451246def _line_(self, R, Q):1247r"""1248Computes the value at `Q` of a straight line through points1249self and `R`.12501251INPUT:12521253- ``R, Q`` -- points on self.curve() with ``Q`` nonzero.12541255OUTPUT:12561257An element of the base field self.curve().base_field().1258A ValueError is raised if ``Q`` is zero.12591260EXAMPLES::12611262sage: F.<a>=GF(2^5)1263sage: E=EllipticCurve(F,[0,0,1,1,1])1264sage: P = E(a^4 + 1, a^3)1265sage: Q = E(a^4, a^4 + a^3)1266sage: O = E(0)1267sage: P._line_(P,-2*P) == 01268True1269sage: P._line_(Q,-(P+Q)) == 01270True1271sage: O._line_(O,Q) == F(1)1272True1273sage: P._line_(O,Q) == a^4 - a^4 + 11274True1275sage: P._line_(13*P,Q) == a^41276True1277sage: P._line_(P,Q) == a^4 + a^3 + a^2 + 11278True12791280See :trac:`7116`::12811282sage: P._line_ (Q,O)1283Traceback (most recent call last):1284...1285ValueError: Q must be nonzero.12861287..NOTES:12881289This function is used in _miller_ algorithm.12901291AUTHOR:12921293- David Hansen (2009-01-25)1294"""1295if Q.is_zero():1296raise ValueError("Q must be nonzero.")12971298if self.is_zero() or R.is_zero():1299if self == R:1300return self.curve().base_field().one_element()1301if self.is_zero():1302return Q[0] - R[0]1303if R.is_zero():1304return Q[0] - self[0]1305elif self != R:1306if self[0] == R[0]:1307return Q[0] - self[0]1308else:1309l = (R[1] - self[1])/(R[0] - self[0])1310return Q[1] - self[1] - l * (Q[0] - self[0])1311else:1312a1, a2, a3, a4, a6 = self.curve().a_invariants()1313numerator = (3*self[0]**2 + 2*a2*self[0] + a4 - a1*self[1])1314denominator = (2*self[1] + a1*self[0] + a3)1315if denominator == 0:1316return Q[0] - self[0]1317else:1318l = numerator/denominator1319return Q[1] - self[1] - l * (Q[0] - self[0])13201321def _miller_(self, Q, n):1322r"""1323Return the value at `Q` of the rational function `f_{n,P}`, where the1324divisor of `f_{n,P}` is `n[P]-[nP]-(n-1)[O]`.13251326INPUT:13271328- ``Q`` -- a nonzero point on self.curve().13291330- ``n`` -- an nonzero integer. If `n<0` then return `Q`1331evaluated at `1/(v_{nP}*f_{n,P})` (used in the ate pairing).13321333OUTPUT:13341335An element in the base field self.curve().base_field()1336A ValueError is raised if `Q` is zero.13371338EXAMPLES::13391340sage: F.<a>=GF(2^5)1341sage: E=EllipticCurve(F,[0,0,1,1,1])1342sage: P = E(a^4 + 1, a^3)1343sage: Fx.<b>=GF(2^(4*5))1344sage: Ex=EllipticCurve(Fx,[0,0,1,1,1])1345sage: phi=Hom(F,Fx)(F.gen().minpoly().roots(Fx)[0][0])1346sage: Px=Ex(phi(P.xy()[0]),phi(P.xy()[1]))1347sage: 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)1348sage: Px._miller_(Qx,41) == b^17 + b^13 + b^12 + b^9 + b^8 + b^6 + b^4 + 11349True1350sage: Qx._miller_(Px,41) == b^13 + b^10 + b^8 + b^7 + b^6 + b^51351True1352sage: P._miller_(E(0),41)1353Traceback (most recent call last):1354...1355ValueError: Q must be nonzero.13561357An example of even order::13581359sage: F.<a> = GF(19^4)1360sage: E = EllipticCurve(F,[-1,0])1361sage: P = E(15*a^3 + 17*a^2 + 14*a + 13,16*a^3 + 7*a^2 + a + 18)1362sage: Q = E(10*a^3 + 16*a^2 + 4*a + 2, 6*a^3 + 4*a^2 + 3*a + 2)1363sage: x=P.weil_pairing(Q,360)1364sage: x^360 == F(1)1365True13661367You can use the _miller_ function on linearly dependent1368points, but with the risk of a dividing with zero::13691370sage: Px._miller_(2*Px,41)1371Traceback (most recent call last):1372...1373ZeroDivisionError: division by zero in finite field.13741375A small example of embedding degree 6::13761377sage: q = 401; F = GF(q); a = 146; b = 400; k = 61378sage: E = EllipticCurve([F(a), F(b)])1379sage: R.<x> = F[]; K.<a> = GF(q^k, modulus=x^6 + 4*x^4 + 115*x^3 + 81*x^2 + 51*x + 3)1380sage: EK = E.base_extend(K)1381sage: P = E([F(338), F(227)])1382sage: Q_x = 333*a^5 + 391*a^4 + 160*a^3 + 335*a^2 + 71*a + 931383sage: Q_y = 343*a^5 + 273*a^4 + 26*a^3 + 342*a^2 + 340*a + 2101384sage: Q = EK([Q_x, Q_y])1385sage: P._miller_(Q, 127)1386371*a^5 + 39*a^4 + 355*a^3 + 233*a^2 + 20*a + 27513871388A series of small examples and small torsions. We start with1389`n=1`, which is trivial: the function is the constant13901::13911392sage: E = EllipticCurve([GF(7)(0), 2])1393sage: P = E(5, 1); Q = E(0, 3); I = E(0)1394sage: I._miller_(P, 1)139511396sage: I._miller_(Q, 1)1397113981399A two-torsion example. In this case `f_{n,P}(Q) = x_Q - x_P`::14001401sage: E = EllipticCurve([GF(7)(-1), 0])1402sage: P = E(0,0); Q = E(1, 0);1403sage: P._miller_(P, 2)140401405sage: Q._miller_(Q, 2)140601407sage: P._miller_(Q, 2)140811409sage: Q._miller_(P, 2)1410614111412A three-torsion example::14131414sage: E = EllipticCurve([GF(7)(0), 2])1415sage: P = E(5, 1); Q = E(0, 3);1416sage: P._miller_(Q, 3)1417414181419A 4-torsion example::14201421sage: E = EllipticCurve([GF(7)(-1), 0])1422sage: P = E(5, 1); Q = E(4, 2)1423sage: P._miller_(Q, 4)1424314251426A 5-torsion example::14271428sage: E = EllipticCurve([GF(7)(-1), 4])1429sage: P = E(4, 1); Q = E(6, 5)1430sage: P._miller_(Q, 5)1431114321433A 6-torsion example::14341435sage: E = EllipticCurve([GF(7)(3), 1])1436sage: P = E(5, 1); Q = E(3, 3)1437sage: P._miller_(Q, 6)1438514391440An example which is part of an ate pairing calculation. The trace of1441the curve is negative, so it should exercise the `n<0` case in the1442code::14431444sage: p = 2017; A = 1; B = 30; r = 29; t = -70; k = 7;1445sage: F = GF(p); R.<x> = F[]1446sage: E = EllipticCurve([F(A), F(B)]); P = E(369, 716)1447sage: K.<a> = GF(p^k, modulus=x^k+2); EK = E.base_extend(K)1448sage: Qx = 1226*a^6 + 1778*a^5 + 660*a^4 + 1791*a^3 + 1750*a^2 + 867*a + 7701449sage: Qy = 1764*a^6 + 198*a^5 + 1206*a^4 + 406*a^3 + 1200*a^2 + 273*a + 17121450sage: Q = EK(Qx, Qy)1451sage: Q._miller_(P, t-1)14521311*a^6 + 1362*a^5 + 1177*a^4 + 807*a^3 + 1331*a^2 + 1530*a + 193114531454ALGORITHM:14551456Double-and-add. See also [Mil04]_.14571458AUTHORS:14591460- David Hansen (2009-01-25)1461- Mariah Lenox (2011-03-07) -- Added more doctests and support for1462negative n.1463"""1464if Q.is_zero():1465raise ValueError("Q must be nonzero.")1466if n.is_zero():1467raise ValueError("n must be nonzero.")1468n_is_negative = False1469if n < 0:1470n = n.abs()1471n_is_negative = True14721473one = self.curve().base_field().one_element()1474t = one1475V = self1476S = 2*V1477nbin = n.bits()1478i = n.nbits() - 21479while i > -1:1480S = 2*V1481ell = V._line_(V, Q)1482vee = S._line_(-S, Q)1483t = (t**2)*(ell/vee)1484V = S1485if nbin[i] == 1:1486S = V+self1487ell = V._line_(self, Q)1488vee = S._line_(-S, Q)1489t = t*(ell/vee)1490V = S1491i = i-11492if n_is_negative:1493vee = V._line_(-V, Q)1494t = 1/(t*vee)1495return t14961497def weil_pairing(self, Q, n):1498r"""1499Compute the Weil pairing of self and `Q` using Miller's algorithm.15001501INPUT:15021503- ``Q`` -- a point on self.curve().15041505- ``n`` -- an integer `n` such that `nP = nQ = (0:1:0)` where1506`P` = self.15071508OUTPUT:15091510An `n`'th root of unity in the base field self.curve().base_field()15111512EXAMPLES::15131514sage: F.<a>=GF(2^5)1515sage: E=EllipticCurve(F,[0,0,1,1,1])1516sage: P = E(a^4 + 1, a^3)1517sage: Fx.<b>=GF(2^(4*5))1518sage: Ex=EllipticCurve(Fx,[0,0,1,1,1])1519sage: phi=Hom(F,Fx)(F.gen().minpoly().roots(Fx)[0][0])1520sage: Px=Ex(phi(P.xy()[0]),phi(P.xy()[1]))1521sage: O = Ex(0)1522sage: 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)1523sage: Px.weil_pairing(Qx,41) == b^19 + b^15 + b^9 + b^8 + b^6 + b^4 + b^3 + b^2 + 11524True1525sage: Px.weil_pairing(17*Px,41) == Fx(1)1526True1527sage: Px.weil_pairing(O,41) == Fx(1)1528True15291530An error is raised if either point is not n-torsion::15311532sage: Px.weil_pairing(O,40)1533Traceback (most recent call last):1534...1535ValueError: points must both be n-torsion15361537A larger example (see :trac:`4964`)::15381539sage: P,Q = EllipticCurve(GF(19^4,'a'),[-1,0]).gens()1540sage: P.order(), Q.order()1541(360, 360)1542sage: z = P.weil_pairing(Q,360)1543sage: z.multiplicative_order()154436015451546An example over a number field::15471548sage: P,Q = EllipticCurve('11a1').change_ring(CyclotomicField(5)).torsion_subgroup().gens() # long time (10s)1549sage: P,Q = (P.element(), Q.element()) # long time1550sage: (P.order(),Q.order()) # long time1551(5, 5)1552sage: P.weil_pairing(Q,5) # long time1553zeta5^21554sage: Q.weil_pairing(P,5) # long time1555zeta5^315561557ALGORITHM:15581559Implemented using Proposition 8 in [Mil04]_. The value 1 is1560returned for linearly dependent input points. This condition1561is caught via a DivisionByZeroError, since the use of a1562discrete logarithm test for linear dependence, is much too slow1563for large `n`.15641565REFERENCES:15661567.. [Mil04] Victor S. Miller, "The Weil pairing, and its1568efficient calculation", J. Cryptol., 17(4):235-261, 200415691570AUTHOR:15711572- David Hansen (2009-01-25)1573"""1574P = self1575E = P.curve()15761577if not Q.curve() is E:1578raise ValueError("points must both be on the same curve")15791580# Test if P, Q are both in E[n]1581if not ((n*P).is_zero() and (n*Q).is_zero()):1582raise ValueError("points must both be n-torsion")15831584one = E.base_field().one_element()15851586# Case where P = Q1587if P == Q:1588return one15891590# Case where P = O or Q = O1591if P.is_zero() or Q.is_zero():1592return one15931594# The non-trivial case P != Q15951596# Note (2010-12-29): The following code block should not be1597# used. It attempts to reduce the pairing calculation to order1598# d = gcd(|P|,|Q|) instead of order n, but this approach is1599# counterproductive, since calculating |P| and |Q| is much1600# slower than calculating the pairing [BGN05].1601#1602# [BGN05] D. Boneh, E. Goh, and K. Nissim, "Evaluating 2-DNF1603# Formulas on Ciphertexts", TCC 2005, LNCS 3378, pp. 325-341.1604#1605# Reduction to order d = gcd(|P|,|Q|); value is a d'th root of unity1606# try:1607# nP = P.order()1608# except AttributeError:1609# nP = generic.order_from_multiple(P,n,operation='+')1610# try:1611# nQ = Q.order()1612# except AttributeError:1613# nQ = generic.order_from_multiple(Q,n,operation='+')1614# d = arith.gcd(nP,nQ)1615# if d==1:1616# return one1617#1618# P = (nP//d)*P # order d1619# Q = (nQ//d)*Q # order d1620# n = d16211622try:1623return ((-1)**n.test_bit(0))*(P._miller_(Q, n)/Q._miller_(P, n))1624except ZeroDivisionError:1625return one16261627def tate_pairing(self, Q, n, k, q=None):1628r"""1629Return Tate pairing of `n`-torsion point `P = self` and point `Q`.16301631The value returned is `f_{n,P}(Q)^e` where `f_{n,P}` is a function with1632divisor `n[P]-n[O].`. This is also known as the "modified Tate1633pairing". It is a well-defined bilinear map.16341635INPUT:16361637- ``P=self`` -- Elliptic curve point having order n16381639- ``Q`` -- Elliptic curve point on same curve as P (can be any order)16401641- ``n`` -- positive integer: order of P16421643- ``k`` -- positive integer: embedding degree16441645- ``q`` -- positive integer: size of base field (the "big"1646field is `GF(q^k)`. `q` needs to be set only if its value1647cannot be deduced.)16481649OUTPUT:16501651An `n`'th root of unity in the base field self.curve().base_field()16521653EXAMPLES:16541655A simple example, pairing a point with itself, and pairing a point with1656another rational point::16571658sage: p = 103; A = 1; B = 18; E = EllipticCurve(GF(p), [A, B])1659sage: P = E(33, 91); n = P.order(); n1660191661sage: k = GF(n)(p).multiplicative_order(); k166261663sage: P.tate_pairing(P, n, k)166411665sage: Q = E(87, 51)1666sage: P.tate_pairing(Q, n, k)166711668sage: set_random_seed(35)1669sage: P.tate_pairing(P,n,k)1670116711672We now let Q be a point on the same curve as above, but defined over1673the pairing extension field, and we also demonstrate the bilinearity of1674the pairing::16751676sage: K.<a> = GF(p^k)1677sage: EK = E.base_extend(K); P = EK(P)1678sage: Qx = 69*a^5 + 96*a^4 + 22*a^3 + 86*a^2 + 6*a + 351679sage: Qy = 34*a^5 + 24*a^4 + 16*a^3 + 41*a^2 + 4*a + 401680sage: Q = EK(Qx, Qy);1681sage: # multiply by cofactor so Q has order n:1682sage: h = 551269674; Q = h*Q1683sage: P = EK(P); P.tate_pairing(Q, n, k)168424*a^5 + 34*a^4 + 3*a^3 + 69*a^2 + 86*a + 451685sage: s = Integer(randrange(1,n))1686sage: ans1 = (s*P).tate_pairing(Q, n, k)1687sage: ans2 = P.tate_pairing(s*Q, n, k)1688sage: ans3 = P.tate_pairing(Q, n, k)^s1689sage: ans1 == ans2 == ans31690True1691sage: (ans1 != 1) and (ans1^n == 1)1692True16931694Here is an example of using the Tate pairing to compute the Weil1695pairing (using the same data as above)::16961697sage: e = Integer((p^k-1)/n); e1698628448577121699sage: P.weil_pairing(Q, n)^e170094*a^5 + 99*a^4 + 29*a^3 + 45*a^2 + 57*a + 341701sage: P.tate_pairing(Q, n, k) == P._miller_(Q, n)^e1702True1703sage: Q.tate_pairing(P, n, k) == Q._miller_(P, n)^e1704True1705sage: P.tate_pairing(Q, n, k)/Q.tate_pairing(P, n, k)170694*a^5 + 99*a^4 + 29*a^3 + 45*a^2 + 57*a + 3417071708An example where we have to pass the base field size (and we again have1709agreement with the Weil pairing)::17101711sage: F.<a>=GF(2^5)1712sage: E=EllipticCurve(F,[0,0,1,1,1])1713sage: P = E(a^4 + 1, a^3)1714sage: Fx.<b>=GF(2^(4*5))1715sage: Ex=EllipticCurve(Fx,[0,0,1,1,1])1716sage: phi=Hom(F,Fx)(F.gen().minpoly().roots(Fx)[0][0])1717sage: Px=Ex(phi(P.xy()[0]),phi(P.xy()[1]))1718sage: 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)1719sage: Px.tate_pairing(Qx, n=41, k=4)1720Traceback (most recent call last):1721...1722ValueError: Unexpected field degree: set keyword argument q equal to the size of the base field (big field is GF(q^4)).1723sage: num = Px.tate_pairing(Qx, n=41, k=4, q=32); num1724b^19 + b^14 + b^13 + b^12 + b^6 + b^4 + b^31725sage: den = Qx.tate_pairing(Px, n=41, k=4, q=32); den1726b^19 + b^17 + b^16 + b^15 + b^14 + b^10 + b^6 + b^2 + 11727sage: e = Integer((32^4-1)/41); e1728255751729sage: Px.weil_pairing(Qx, 41)^e == num/den1730True17311732NOTES:17331734This function uses Miller's algorithm, followed by a naive1735exponentiation. It does not do anything fancy. In the case1736that there is an issue with `Q` being on one of the lines1737generated in the `r*P` calculation, `Q` is offset by a random1738point `R` and P.tate_pairing(Q+R,n,k)/P.tate_pairing(R,n,k)1739is returned.17401741AUTHORS:17421743- Mariah Lenox (2011-03-07)1744"""1745P = self1746E = P.curve()17471748if not Q.curve() is E:1749raise ValueError("Points must both be on the same curve")17501751K = E.base_ring()1752d = K.degree()1753if q is None:1754if d == 1:1755q = K.order()1756elif d == k:1757q = K.base_ring().order()1758else:1759raise ValueError("Unexpected field degree: set keyword argument q equal to the size of the base field (big field is GF(q^%s))." % k)17601761if n*P != E(0):1762raise ValueError('This point is not of order n=%s' % n)17631764# In small cases, or in the case of pairing an element with1765# itself, Q could be on one of the lines in the Miller1766# algorithm. If this happens we try again, with an offset of a1767# random point.1768try:1769ret = self._miller_(Q, n)1770e = rings.Integer((q**k - 1)/n)1771ret = ret**e1772except (ZeroDivisionError, ValueError):1773R = E.random_point()1774ret = self.tate_pairing(Q + R, n, k)/self.tate_pairing(R, n, k)1775return ret17761777def ate_pairing(self, Q, n, k, t, q=None):1778r"""1779Return ate pairing of `n`-torsion points `P=self` and `Q`.17801781Also known as the `n`-th modified ate pairing. `P` is `GF(q)`-rational,1782and `Q` must be an element of `Ker(\pi-p)`, where `\pi` is the1783`q`-frobenius map (and hence `Q` is `GF(q^k)`-rational).17841785INPUT:17861787- ``P=self`` -- a point of order `n`, in `ker(\pi-1)`, where1788`\pi` is the `q`-Frobenius map (e.g., `P` is `q-rational`).17891790- ``Q`` -- a point of order `n` in `ker(\pi-q)`17911792- ``n`` -- the order of `P` and `Q`.17931794- ``k`` -- the embedding degree.17951796- ``t`` -- the trace of Frobenius of the curve over `GF(q)`.17971798- ``q`` -- (default:None) the size of base field (the "big"1799field is `GF(q^k)`). `q` needs to be set only if its value1800cannot be deduced.18011802OUTPUT:18031804FiniteFieldElement in `GF(q^k)` -- the ate pairing of `P` and `Q`.18051806EXAMPLES:18071808An example with embedding degree 6::18091810sage: p = 7549; A = 0; B = 1; n = 157; k = 6; t = 141811sage: F = GF(p); E = EllipticCurve(F, [A, B])1812sage: R.<x> = F[]; K.<a> = GF(p^k, modulus=x^k+2)1813sage: EK = E.base_extend(K)1814sage: P = EK(3050, 5371); Q = EK(6908*a^4, 3231*a^3)1815sage: P.ate_pairing(Q, n, k, t)18166708*a^5 + 4230*a^4 + 4350*a^3 + 2064*a^2 + 4022*a + 67331817sage: s = Integer(randrange(1, n))1818sage: (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)1819True1820sage: P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)^s1821True18221823Another example with embedding degree 7 and positive trace::18241825sage: p = 2213; A = 1; B = 49; n = 1093; k = 7; t = 281826sage: F = GF(p); E = EllipticCurve(F, [A, B])1827sage: R.<x> = F[]; K.<a> = GF(p^k, modulus=x^k+2)1828sage: EK = E.base_extend(K)1829sage: P = EK(1583, 1734)1830sage: Qx = 1729*a^6+1767*a^5+245*a^4+980*a^3+1592*a^2+1883*a+7221831sage: Qy = 1299*a^6+1877*a^5+1030*a^4+1513*a^3+1457*a^2+309*a+16361832sage: Q = EK(Qx, Qy)1833sage: P.ate_pairing(Q, n, k, t)18341665*a^6 + 1538*a^5 + 1979*a^4 + 239*a^3 + 2134*a^2 + 2151*a + 6541835sage: s = Integer(randrange(1, n))1836sage: (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)1837True1838sage: P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)^s1839True18401841Another example with embedding degree 7 and negative trace::18421843sage: p = 2017; A = 1; B = 30; n = 29; k = 7; t = -701844sage: F = GF(p); E = EllipticCurve(F, [A, B])1845sage: R.<x> = F[]; K.<a> = GF(p^k, modulus=x^k+2)1846sage: EK = E.base_extend(K)1847sage: P = EK(369, 716)1848sage: Qx = 1226*a^6+1778*a^5+660*a^4+1791*a^3+1750*a^2+867*a+7701849sage: Qy = 1764*a^6+198*a^5+1206*a^4+406*a^3+1200*a^2+273*a+17121850sage: Q = EK(Qx, Qy)1851sage: P.ate_pairing(Q, n, k, t)18521794*a^6 + 1161*a^5 + 576*a^4 + 488*a^3 + 1950*a^2 + 1905*a + 13151853sage: s = Integer(randrange(1, n))1854sage: (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)1855True1856sage: P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)^s1857True18581859Using the same data, we show that the ate pairing is a power of the1860Tate pairing (see [HSV]_ end of section 3.1)::18611862sage: c = (k*p^(k-1)).mod(n); T = t - 11863sage: N = gcd(T^k - 1, p^k - 1)1864sage: s = Integer(N/n)1865sage: L = Integer((T^k - 1)/N)1866sage: M = (L*s*c.inverse_mod(n)).mod(n)1867sage: P.ate_pairing(Q, n, k, t) == Q.tate_pairing(P, n, k)^M1868True18691870An example where we have to pass the base field size (and we again have1871agreement with the Tate pairing). Note that though `Px` is not1872`F`-rational, (it is the homomorphic image of an `F`-rational point) it1873is nonetheless in `ker(\pi-1)`, and so is a legitimate input::18741875sage: q = 2^5; F.<a>=GF(q)1876sage: n = 41; k = 4; t = -81877sage: E=EllipticCurve(F,[0,0,1,1,1])1878sage: P = E(a^4 + 1, a^3)1879sage: Fx.<b>=GF(q^k)1880sage: Ex=EllipticCurve(Fx,[0,0,1,1,1])1881sage: phi=Hom(F,Fx)(F.gen().minpoly().roots(Fx)[0][0])1882sage: Px=Ex(phi(P.xy()[0]),phi(P.xy()[1]))1883sage: 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)1884sage: Qx = Ex(Qx[0]^q, Qx[1]^q) - Qx # ensure Qx is in ker(pi - q)1885sage: Px.ate_pairing(Qx, n, k, t)1886Traceback (most recent call last):1887...1888ValueError: Unexpected field degree: set keyword argument q equal to the size of the base field (big field is GF(q^4)).1889sage: Px.ate_pairing(Qx, n, k, t, q)1890b^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 + 11891sage: s = Integer(randrange(1, n))1892sage: (s*Px).ate_pairing(Qx, n, k, t, q) == Px.ate_pairing(s*Qx, n, k, t, q)1893True1894sage: Px.ate_pairing(s*Qx, n, k, t, q) == Px.ate_pairing(Qx, n, k, t, q)^s1895True1896sage: c = (k*q^(k-1)).mod(n); T = t - 11897sage: N = gcd(T^k - 1, q^k - 1)1898sage: s = Integer(N/n)1899sage: L = Integer((T^k - 1)/N)1900sage: M = (L*s*c.inverse_mod(n)).mod(n)1901sage: Px.ate_pairing(Qx, n, k, t, q) == Qx.tate_pairing(Px, n, k, q)^M1902True19031904It is an error if `Q` is not in the kernel of `\pi - p`, where `\pi` is1905the Frobenius automorphism::19061907sage: p = 29; A = 1; B = 0; n = 5; k = 2; t = 101908sage: F = GF(p); R.<x> = F[]1909sage: E = EllipticCurve(F, [A, B]);1910sage: K.<a> = GF(p^k, modulus=x^k+2); EK = E.base_extend(K)1911sage: P = EK(13, 8); Q = EK(13, 21)1912sage: P.ate_pairing(Q, n, k, t)1913Traceback (most recent call last):1914...1915ValueError: Point (13 : 21 : 1) not in Ker(pi - q)19161917It is also an error if `P` is not in the kernel os `\pi - 1`::19181919sage: p = 29; A = 1; B = 0; n = 5; k = 2; t = 101920sage: F = GF(p); R.<x> = F[]1921sage: E = EllipticCurve(F, [A, B]);1922sage: K.<a> = GF(p^k, modulus=x^k+2); EK = E.base_extend(K)1923sage: P = EK(14, 10*a); Q = EK(13, 21)1924sage: P.ate_pairing(Q, n, k, t)1925Traceback (most recent call last):1926...1927ValueError: This point (14 : 10*a : 1) is not in Ker(pi - 1)19281929NOTES:19301931First defined in the paper of [HSV]_, the ate pairing can be1932computationally effective in those cases when the trace of the curve1933over the base field is significantly smaller than the expected1934value. This implementation is simply Miller's algorithm followed by a1935naive exponentiation, and makes no claims towards efficiency.19361937REFERENCES:19381939.. [HSV] Hess, Smart, Vercauteren, "The Eta Pairing Revisited",1940IEEE Trans. Information Theory, 52(10): 4595-4602, 2006.19411942AUTHORS:19431944- Mariah Lenox (2011-03-08)1945"""1946P = self1947# check for same curve1948E = P.curve()1949O = E(0)1950if not Q.curve() is E:1951raise ValueError("Points must both be on the same curve")19521953# set q to be the order of the base field1954K = E.base_ring()1955if q is None:1956d = K.degree()1957if d == k:1958q = K.base_ring().order()1959else:1960raise ValueError("Unexpected field degree: set keyword argument q equal to the size of the base field (big field is GF(q^%s))." % k)19611962# check order of P1963if n*P != O:1964raise ValueError('This point %s is not of order n=%s' % (P, n))19651966# check for P in kernel pi - 1:1967piP = E(P[0]**q, P[1]**q)1968if piP - P != O:1969raise ValueError('This point %s is not in Ker(pi - 1)' % P)19701971# check for Q in kernel pi - q:1972piQ = E(Q[0]**q, Q[1]**q)1973if piQ - q*Q != O:1974raise ValueError('Point %s not in Ker(pi - q)' % Q)19751976T = t-11977ret = Q._miller_(P, T)1978e = rings.Integer((q**k - 1)/n)1979ret = ret**e1980return ret198119821983class EllipticCurvePoint_number_field(EllipticCurvePoint_field):1984"""1985A point on an elliptic curve over a number field.19861987Most of the functionality is derived from the parent class1988``EllipticCurvePoint_field``. In addition we have support for1989orders, heights, reduction modulo primes, and elliptic logarithms.19901991EXAMPLES::19921993sage: E = EllipticCurve('37a')1994sage: E([0,0])1995(0 : 0 : 1)1996sage: E(0,0) # brackets are optional1997(0 : 0 : 1)1998sage: E([GF(5)(0), 0]) # entries are coerced1999(0 : 0 : 1)20002001sage: E(0.000, 0)2002(0 : 0 : 1)20032004sage: E(1,0,0)2005Traceback (most recent call last):2006...2007TypeError: Coordinates [1, 0, 0] do not define a point on2008Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field20092010::20112012sage: E = EllipticCurve([0,0,1,-1,0])2013sage: S = E(QQ); S2014Abelian group of points on Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field20152016TESTS::20172018sage: loads(S.dumps()) == S2019True2020sage: P = E(0,0); P2021(0 : 0 : 1)2022sage: loads(P.dumps()) == P2023True2024sage: T = 100*P2025sage: loads(T.dumps()) == T2026True20272028Test pickling an elliptic curve that has known points on it::20292030sage: e = EllipticCurve([0, 0, 1, -1, 0]); g = e.gens(); loads(dumps(e)) == e2031True2032"""20332034def order(self):2035r"""2036Return the order of this point on the elliptic curve.20372038If the point has infinite order, returns +Infinity. For2039curves defined over `\QQ`, we call PARI; over other2040number fields we implement the function here.20412042.. NOTE::20432044:meth:`additive_order` is a synonym for :meth:`order`20452046EXAMPLES::20472048sage: E = EllipticCurve([0,0,1,-1,0])2049sage: P = E([0,0]); P2050(0 : 0 : 1)2051sage: P.order()2052+Infinity20532054::20552056sage: E = EllipticCurve([0,1])2057sage: P = E([-1,0])2058sage: P.order()205922060sage: P.additive_order()2061220622063"""2064try:2065return self._order2066except AttributeError:2067pass20682069if self.is_zero():2070self._order = rings.Integer(1)2071return self._order20722073E = self.curve()20742075# Special code for curves over Q, calling PARI2076from sage.libs.pari.all import PariError2077try:2078n = int(E.pari_curve().ellorder(self))2079if n == 0:2080n = oo2081self._order = n2082return n2083except PariError:2084pass20852086# Get the torsion order if known, else a bound on (multiple2087# of) the order. We do not compute the torsion if it is not2088# already known, since computing the bound is faster (and is2089# also cached).20902091try:2092N = E._torsion_order2093except AttributeError:2094N = E._torsion_bound()20952096# Now self is a torsion point iff it is killed by N:2097if not (N*self).is_zero():2098self._order = oo2099return self._order21002101# Finally we find the exact order using the generic code:2102self._order = generic.order_from_multiple(self, N, operation='+')2103return self._order21042105additive_order = order21062107def has_finite_order(self):2108"""2109Return True iff this point has finite order on the elliptic curve.21102111EXAMPLES::21122113sage: E = EllipticCurve([0,0,1,-1,0])2114sage: P = E([0,0]); P2115(0 : 0 : 1)2116sage: P.has_finite_order()2117False21182119::21202121sage: E = EllipticCurve([0,1])2122sage: P = E([-1,0])2123sage: P.has_finite_order()2124True2125"""2126if self.is_zero():2127return True2128return self.order() != oo21292130def has_infinite_order(self):2131r"""2132Return True iff this point has infinite order on the elliptic curve.21332134EXAMPLES::21352136sage: E = EllipticCurve([0,0,1,-1,0])2137sage: P = E([0,0]); P2138(0 : 0 : 1)2139sage: P.has_infinite_order()2140True21412142::21432144sage: E = EllipticCurve([0,1])2145sage: P = E([-1,0])2146sage: P.has_infinite_order()2147False2148"""2149if self.is_zero():2150return False2151return self.order() == oo21522153def is_on_identity_component(self, embedding=None):2154r"""2155Returns True iff this point is on the identity component of2156its curve with respect to a given (real or complex) embedding.21572158INPUT:21592160- ``self`` -- a point on a curve over any ordered field (e.g. `\QQ`)21612162- ``embedding`` -- an embedding from the base_field of the2163point's curve into `\RR` or `\CC`; if ``None`` (the2164default) it uses the first embedding of the base_field into2165`\RR` if any, else the first embedding into `\CC`.21662167OUTPUT:21682169(bool) -- ``True`` iff the point is on the identity component of2170the curve. (If the point is zero then the result is True.)21712172EXAMPLES:21732174For `K=\QQ` there is no need to specify an embedding::21752176sage: E=EllipticCurve('5077a1')2177sage: [E.lift_x(x).is_on_identity_component() for x in range(-3,5)]2178[False, False, False, False, False, True, True, True]21792180An example over a field with two real embeddings::21812182sage: L.<a> = QuadraticField(2)2183sage: E=EllipticCurve(L,[0,1,0,a,a])2184sage: P=E(-1,0)2185sage: [P.is_on_identity_component(e) for e in L.embeddings(RR)]2186[False, True]21872188We can check this as follows::21892190sage: [e(E.discriminant())>0 for e in L.embeddings(RR)]2191[True, False]2192sage: e = L.embeddings(RR)[0]2193sage: E1 = EllipticCurve(RR,[e(ai) for ai in E.ainvs()])2194sage: e1,e2,e3 = E1.two_division_polynomial().roots(RR,multiplicities=False)2195sage: e1 < e2 < e3 and e(P[0]) < e32196True21972198"""2199if self.is_zero(): # trivial case2200return True22012202e = embedding2203# It is also trivially true if we have a complex embedding2204if not e is None:2205if not is_RealField(e.codomain()):2206return True22072208# find a suitable embedding if none was supplied:2209E = self.curve()2210K = E.base_field()2211if e is None:2212try:2213e = K.embeddings(rings.RealField())[0]2214except IndexError:2215e = K.embeddings(rings.ComplexField())[0]22162217# If there is only one component, the result is True:2218if not is_RealField(e.codomain()): # complex embedding2219return True2220if e(E.discriminant()) < 0: # only one component2221return True22222223# Now we have a real embedding and two components and have to work:2224gx = E.two_division_polynomial()2225gxd = gx.derivative()2226gxdd = gxd.derivative()2227return (e(gxd(self[0])) > 0 and e(gxdd(self[0])) > 0)22282229def has_good_reduction(self, P=None):2230r"""2231Returns True iff this point has good reduction modulo a prime.22322233INPUT:22342235- ``P`` -- a prime of the base_field of the point's curve, or2236``None`` (default)22372238OUTPUT:22392240(bool) If a prime `P` of the base field is specified, returns2241True iff the point has good reduction at `P`; otherwise,2242return true if the point has god reduction at all primes in2243the support of the discriminant of this model.22442245EXAMPLES::22462247sage: E = EllipticCurve('990e1')2248sage: P = E.gen(0); P2249(15 : 51 : 1)2250sage: [E.has_good_reduction(p) for p in [2,3,5,7]]2251[False, False, False, True]2252sage: [P.has_good_reduction(p) for p in [2,3,5,7]]2253[True, False, True, True]2254sage: [E.tamagawa_exponent(p) for p in [2,3,5,7]]2255[2, 2, 1, 1]2256sage: [(2*P).has_good_reduction(p) for p in [2,3,5,7]]2257[True, True, True, True]2258sage: P.has_good_reduction()2259False2260sage: (2*P).has_good_reduction()2261True2262sage: (3*P).has_good_reduction()2263False22642265::22662267sage: K.<i> = NumberField(x^2+1)2268sage: E = EllipticCurve(K,[0,1,0,-160,308])2269sage: P = E(26,-120)2270sage: E.discriminant().support()2271[Fractional ideal (i + 1),2272Fractional ideal (-i - 2),2273Fractional ideal (i - 2),2274Fractional ideal (3)]2275sage: [E.tamagawa_exponent(p) for p in E.discriminant().support()]2276[1, 4, 4, 4]2277sage: P.has_good_reduction()2278False2279sage: (2*P).has_good_reduction()2280False2281sage: (4*P).has_good_reduction()2282True22832284TESTS:22852286An example showing that :trac:`8498` is fixed::22872288sage: E = EllipticCurve('11a1')2289sage: K.<t> = NumberField(x^2+47)2290sage: EK = E.base_extend(K)2291sage: T = EK(5,5)2292sage: P = EK(-2, -1/2*t - 1/2)2293sage: p = K.ideal(11)2294sage: T.has_good_reduction(p)2295False2296sage: P.has_good_reduction(p)2297True2298"""2299if self.is_zero(): # trivial case2300return True23012302E = self.curve()2303if P is None:2304return all([self.has_good_reduction(Pr) for Pr in E.discriminant().support()])2305K = E.base_field()2306from sage.schemes.elliptic_curves.ell_local_data import check_prime2307P = check_prime(K, P)23082309# If the curve has good reduction at P, the result is True:2310t = E.local_data(P).bad_reduction_type()2311if t is None:2312return True23132314# Make sure the curve is integral and locally minimal at P:2315Emin = E.local_minimal_model(P)2316urst = E.isomorphism_to(Emin)2317Q = urst(self)23182319# Scale the homogeneous coordinates of the point to be primitive:2320xyz = list(Q)2321e = min([c.valuation(P) for c in xyz])2322if e != 0:2323if K is rings.QQ:2324pi = P2325else:2326pi = K.uniformizer(P)2327pie = pi**e2328xyz = [c/pie for c in xyz]23292330# Evaluate the partial derivatives at the point to see if they2331# are zero mod P23322333# See #8498: sometimes evaluating F's derivatives at xyz2334# returns a constant polynomial instead of a constant23352336F = Emin.defining_polynomial()2337for v in F.variables():2338c = (F.derivative(v))(xyz)2339try:2340val = c.valuation(P)2341except AttributeError:2342val = c.constant_coefficient().valuation(P)2343if val == 0:2344return True2345return False23462347def reduction(self, p):2348"""2349This finds the reduction of a point `P` on the elliptic curve2350modulo the prime `p`.23512352INPUT:23532354- ``self`` -- A point on an elliptic curve.23552356- ``p`` -- a prime number23572358OUTPUT:23592360The point reduced to be a point on the elliptic curve modulo `p`.23612362EXAMPLES::23632364sage: E = EllipticCurve([1,2,3,4,0])2365sage: P = E(0,0)2366sage: P.reduction(5)2367(0 : 0 : 1)2368sage: Q = E(98,931)2369sage: Q.reduction(5)2370(3 : 1 : 1)2371sage: Q.reduction(5).curve() == E.reduction(5)2372True23732374::23752376sage: F.<a> = NumberField(x^2+5)2377sage: E = EllipticCurve(F,[1,2,3,4,0])2378sage: Q = E(98,931)2379sage: Q.reduction(a)2380(3 : 1 : 1)2381sage: Q.reduction(11)2382(10 : 7 : 1)23832384::23852386sage: F.<a> = NumberField(x^3+x^2+1)2387sage: E = EllipticCurve(F,[a,2])2388sage: P = E(a,1)2389sage: P.reduction(F.ideal(5))2390(abar : 1 : 1)2391sage: P.reduction(F.ideal(a^2-4*a-2))2392(abar : 1 : 1)23932394"""2395P = self2396E = P.curve()2397return E.reduction(p)(P)23982399def height(self, precision=None, normalised=True, algorithm='pari'):2400"""2401Return the Néron-Tate canonical height of the point.24022403INPUT:24042405- ``self`` -- a point on an elliptic curve over a number field2406`K`.24072408- ``precision`` -- positive integer, or None (default). The2409precision in bits of the result. If None, the default real2410precision is used.24112412- ``normalised`` -- boolean. If True (default), the height is2413normalised to be invariant under extension of `K`. If False,2414return this normalised height multiplied by the degree of2415`K`.24162417- ``algorithm`` -- string: either 'pari' (default) or 'sage'.2418If 'pari' and the base field is `\QQ`, use the PARI library2419function; otherwise use the Sage implementation.24202421OUTPUT:24222423The rational number 0, or a non-negative real number.24242425There are two normalisations used in the literature, one of2426which is double the other. We use the larger of the two, which2427is the one appropriate for the BSD conjecture. This is2428consistent with [Cre]_ and double that of [SilBook]_.24292430See :wikipedia:`Néron-Tate height`24312432REFERENCES:24332434.. [Cre] John Cremona, Algorithms for Modular Elliptic Curves.2435Cambridge University Press, 1997.24362437.. [Sil1988] Joseph H. Silverman, Computing heights on2438elliptic curves. Mathematics of Computation, Vol. 51,2439No. 183 (Jul., 1988), pp. 339-358.24402441.. [SilBook] Joseph H. Silverman, The Arithmetic of Elliptic2442Curves. Second edition. Graduate Texts in Mathematics, 106.2443Springer, 2009.24442445EXAMPLES::24462447sage: E = EllipticCurve('11a'); E2448Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field2449sage: P = E([5,5]); P2450(5 : 5 : 1)2451sage: P.height()245202453sage: Q = 5*P2454sage: Q.height()2455024562457::24582459sage: E = EllipticCurve('37a'); E2460Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field2461sage: P = E([0,0])2462sage: P.height()24630.05111140823996882464sage: P.order()2465+Infinity2466sage: E.regulator()24670.0511114082399688...24682469sage: def naive_height(P):2470... return log(RR(max(abs(P[0].numerator()), abs(P[0].denominator()))))2471sage: for n in [1..10]:2472... print naive_height(2^n*P)/4^n24730.00000000000000024740.043321698784996624750.050294934763565624760.051100633561864524770.051100783479961224780.051101366615246624790.051103419990774324800.051110649290647124810.051111408154108224820.051111408154118024832484::24852486sage: E = EllipticCurve('4602a1'); E2487Elliptic Curve defined by y^2 + x*y = x^3 + x^2 - 37746035*x - 89296920339 over Rational Field2488sage: x = 779859224589749492468582291959451034715902489sage: y = 195752602300153137022613790221516759619651571089202635945452232490sage: d = 22540207618847822432491sage: E([ x / d^2, y / d^3 ]).height()249286.740656138127524932494::24952496sage: E = EllipticCurve([17, -60, -120, 0, 0]); E2497Elliptic Curve defined by y^2 + 17*x*y - 120*y = x^3 - 60*x^2 over Rational Field2498sage: E([30, -90]).height()2499025002501sage: E = EllipticCurve('389a1'); E2502Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field2503sage: [P,Q] = [E(-1,1),E(0,-1)]2504sage: P.height(precision=100)25050.686667083305586585723552102952506sage: (3*Q).height(precision=100)/Q.height(precision=100)25079.00000000000000000000000000002508sage: _.parent()2509Real Field with 100 bits of precision25102511Canonical heights over number fields are implemented as well::25122513sage: R.<x> = QQ[]2514sage: K.<a> = NumberField(x^3-2)2515sage: E = EllipticCurve([a, 4]); E2516Elliptic Curve defined by y^2 = x^3 + a*x + 4 over Number Field in a with defining polynomial x^3 - 22517sage: P = E((0,2))2518sage: P.height()25190.8104630965859252520sage: P.height(precision=100)25210.810463096585925368639918105772522sage: P.height(precision=200)25230.810463096585925368639918105768651588961302864171558323780862524sage: (2*P).height() / P.height()25254.000000000000002526sage: (100*P).height() / P.height()252710000.000000000025282529Setting normalised=False multiplies the height by the degree of `K`::25302531sage: E = EllipticCurve('37a')2532sage: P = E([0,0])2533sage: P.height()25340.05111140823996882535sage: P.height(normalised=False)25360.05111140823996882537sage: K.<z> = CyclotomicField(5)2538sage: EK = E.change_ring(K)2539sage: PK = EK([0,0])2540sage: PK.height()25410.05111140823996882542sage: PK.height(normalised=False)25430.20444563295987525442545Some consistency checks::25462547sage: E = EllipticCurve('5077a1')2548sage: P = E([-2,3,1])2549sage: P.height()25501.3685725053539325512552sage: EK = E.change_ring(QuadraticField(-3,'a'))2553sage: PK = EK([-2,3,1])2554sage: PK.height()25551.3685725053539325562557sage: K.<i> = NumberField(x^2+1)2558sage: E = EllipticCurve(K, [0,0,4,6*i,0])2559sage: Q = E.lift_x(-9/4); Q2560(-9/4 : -27/8*i : 1)2561sage: Q.height()25622.695185600179092563sage: (15*Q).height() / Q.height()2564225.00000000000025652566sage: E = EllipticCurve('37a')2567sage: P = E([0,-1])2568sage: P.height()25690.05111140823996882570sage: K.<a> = QuadraticField(-7)2571sage: ED = E.quadratic_twist(-7)2572sage: Q = E.isomorphism_to(ED.change_ring(K))(P); Q2573(0 : -7/2*a - 1/2 : 1)2574sage: Q.height()25750.05111140823996882576sage: Q.height(precision=100)25770.05111140823996884023588609975725782579An example to show that the bug at :trac:`5252` is fixed::25802581sage: E = EllipticCurve([1, -1, 1, -2063758701246626370773726978, 32838647793306133075103747085833809114881])2582sage: P = E([-30987785091199, 258909576181697016447])2583sage: P.height()258425.86031706754622585sage: P.height(precision=100)258625.8603170675461907438688407412587sage: P.height(precision=250)258825.8603170675461907438688407407351103230988729038444162155771710417835725132589sage: P.height(precision=500)259025.860317067546190743868840740735110323098872903844416215577171041783572512955113057088981328179215727850763990997211285601919023612536291419545232172025912592sage: P.height(precision=100) == P.non_archimedean_local_height(prec=100)+P.archimedean_local_height(prec=100)2593True25942595An example to show that the bug at :trac:`8319` is fixed (correct height when the curve is not minimal)::25962597sage: E = EllipticCurve([-5580472329446114952805505804593498080000,-157339733785368110382973689903536054787700497223306368000000])2598sage: xP = 204885147732879546487576840131729064308289385547094673627174585676211859152978311600/236255019070579481322622171889836812048569076577531784154303612599sage: P = E.lift_x(xP)2600sage: P.height()2601157.4325985167542602sage: Q = 2*P2603sage: Q.height() # long time (4s)2604629.7303940670162605sage: Q.height()-4*P.height() # long time26060.00000000000000026072608An example to show that the bug at :trac:`12509` is fixed (precision issues)::26092610sage: x = polygen(QQ)2611sage: K.<a> = NumberField(x^2-x-1)2612sage: v = [0, a + 1, 1, 28665*a - 46382, 2797026*a - 4525688]2613sage: E = EllipticCurve(v)2614sage: P = E([72*a - 509/5, -682/25*a - 434/25])2615sage: P.height()26161.388777116887272617sage: (2*P).height()/P.height()26184.000000000000002619sage: (2*P).height(precision=100)/P.height(precision=100)26204.00000000000000000000000000002621sage: (2*P).height(precision=1000)/P.height(precision=1000)26224.0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000026232624This shows that the bug reported at :trac:`13951` has been fixed::26252626sage: E = EllipticCurve([0,17])2627sage: P1 = E(2,5)2628sage: P1.height()26291.062481376525282630sage: F = E.change_ring(QuadraticField(-3,'a'))2631sage: P2 = F([2,5])2632sage: P2.height()26331.062481376525282634"""2635if self.has_finite_order():2636return rings.QQ(0)26372638E = self.curve()2639K = E.base_ring()26402641if precision is None:2642precision = rings.RealField().precision()26432644known_prec = -12645try:2646height = self.__height2647known_prec = height.prec()2648if known_prec > precision:2649height = rings.RealField(precision)(height)2650except AttributeError:2651pass26522653if known_prec < precision:2654if algorithm == 'pari' and K is rings.QQ:2655Emin = E.minimal_model()2656iso = E.isomorphism_to(Emin)2657P = iso(self)2658h = Emin.pari_curve(prec=precision).ellheight(P, precision=precision)2659height = rings.RealField(precision)(h)2660else:2661height = (self.non_archimedean_local_height(prec=precision)2662+ self.archimedean_local_height(prec=precision))26632664# The cached height is the one that is independent of the base field.2665self.__height = height2666if not normalised:2667height *= K.degree()2668return height26692670def archimedean_local_height(self, v=None, prec=None, weighted=False):2671"""2672Compute the local height of self at the archimedean place `v`.26732674INPUT:26752676- ``self`` -- a point on an elliptic curve over a number field2677`K`.26782679- ``v`` -- a real or complex embedding of K, or None (default).2680If `v` is a real or complex embedding, return the local2681height of self at `v`. If `v` is None, return the total2682archimedean contribution to the global height.26832684- ``prec`` -- integer, or None (default). The precision of the2685computation. If None, the precision is deduced from v.26862687- ``weighted`` -- boolean. If False (default), the height is2688normalised to be invariant under extension of `K`. If True,2689return this normalised height multiplied by the local degree2690if `v` is a single place, or by the degree of `K` if `v` is2691None.26922693OUTPUT:26942695A real number. The normalisation is twice that in Silverman's2696paper [Sil1988]_. Note that this local height depends on the2697model of the curve.26982699ALGORITHM:27002701See [Sil1988]_, Section 4.27022703EXAMPLES:27042705Examples 1, 2, and 3 from [Sil1988]_::27062707sage: K.<a> = QuadraticField(-2)2708sage: E = EllipticCurve(K, [0,-1,1,0,0]); E2709Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 over Number Field in a with defining polynomial x^2 + 22710sage: P = E.lift_x(2+a); P2711(a + 2 : 2*a + 1 : 1)2712sage: P.archimedean_local_height(K.places(prec=170)[0]) / 227130.4575477328752327673621121074142365434657602981469527142715sage: K.<i> = NumberField(x^2+1)2716sage: E = EllipticCurve(K, [0,0,4,6*i,0]); E2717Elliptic Curve defined by y^2 + 4*y = x^3 + 6*i*x over Number Field in i with defining polynomial x^2 + 12718sage: P = E((0,0))2719sage: P.archimedean_local_height(K.places()[0]) / 227200.51018499516237327212722sage: Q = E.lift_x(-9/4); Q2723(-9/4 : -27/8*i : 1)2724sage: Q.archimedean_local_height(K.places()[0]) / 227250.65444561952960027262727An example over the rational numbers::27282729sage: E = EllipticCurve([0, 0, 0, -36, 0])2730sage: P = E([-3, 9])2731sage: P.archimedean_local_height()27321.9872381635077327332734Local heights of torsion points can be non-zero (unlike the2735global height)::27362737sage: K.<i> = QuadraticField(-1)2738sage: E = EllipticCurve([0, 0, 0, K(1), 0])2739sage: P = E(i, 0)2740sage: P.archimedean_local_height()27410.34657359027997327422743TESTS:27442745See :trac:`12509`::27462747sage: x = polygen(QQ)2748sage: K.<a> = NumberField(x^2-x-1)2749sage: v = [0, a + 1, 1, 28665*a - 46382, 2797026*a - 4525688]2750sage: E = EllipticCurve(v)2751sage: P = E([72*a - 509/5, -682/25*a - 434/25])2752sage: P.archimedean_local_height()2753-0.220660795546827849218336274693027542755"""2756E = self.curve()2757K = E.base_ring()27582759if v is None:2760if K is rings.QQ:2761v = K.embeddings(rings.RR)[0]2762h = self.archimedean_local_height(v, prec)2763else:2764r1, r2 = K.signature()2765pl = K.places()2766h = (sum(self.archimedean_local_height(pl[i], prec, weighted=False)2767for i in range(r1))2768+ 2 * sum(self.archimedean_local_height(pl[i], prec, weighted=False)2769for i in range(r1, r1 + r2)))2770if not weighted:2771h /= K.degree()2772return h27732774from sage.rings.number_field.number_field import refine_embedding2775prec_v = v.codomain().prec()2776if prec is None:2777prec = prec_v2778if K is rings.QQ:2779vv = K.embeddings(rings.RealField(max(2*prec, prec_v)))[0]2780else:2781vv = refine_embedding(v, 2*prec) # vv.prec() = max(2*prec, prec_v)2782b2, b4, b6, b8 = [vv(b) for b in E.b_invariants()]2783H = max(4, abs(b2), 2*abs(b4), 2*abs(b6), abs(b8))2784# The following comes from Silverman Theorem 4.2. Silverman2785# uses decimal precision d, so his term (5/3)d =2786# (5/3)*(log(2)/log(10))*prec = 0.5017*prec, which we round2787# up. The rest of the expression was wrongly transcribed in2788# Sage versions <5.6 (see #12509).2789nterms = int(math.ceil(0.51*prec + 0.5 + 3*math.log(7+(4*math.log(H)+math.log(max(1, ~abs(v(E.discriminant())))))/3)/4))2790b2p = b2 - 122791b4p = b4 - b2 + 62792b6p = b6 - 2*b4 + b2 - 42793b8p = b8 - 3*b6 + 3*b4 - b2 + 327942795fz = lambda T: 1 - T**2 * (b4 + T*(2*b6 + T*b8))2796fzp = lambda T: 1 - T**2 * (b4p + T*(2*b6p + T*b8p))2797fw = lambda T: T*(4 + T*(b2 + T*(2*b4 + T*b6)))2798fwp = lambda T: T*(4 + T*(b2p + T*(2*b4p + T*b6p)))27992800x = vv(self[0])2801if abs(x) >= .5:2802t = 1/x2803beta = True2804else:2805t = 1/(x+1)2806beta = False2807lam = -t.abs().log()2808mu = 02809four_to_n = rings.QQ(1)28102811for n in range(nterms):2812if beta:2813w = fw(t)2814z = fz(t)2815if abs(w) <= 2 * abs(z):2816mu += four_to_n * z.abs().log()2817t = w/z2818else:2819mu += four_to_n * (z+w).abs().log()2820t = w/(z+w)2821beta = not beta2822else:2823w = fwp(t)2824z = fzp(t)2825if abs(w) <= 2 * abs(z):2826mu += four_to_n * z.abs().log()2827t = w/z2828else:2829mu += four_to_n * (z-w).abs().log()2830t = w/(z-w)2831beta = not beta2832four_to_n >>= 22833h = rings.RealField(prec)(lam + mu/4)2834if weighted and not v.im_gens()[0] in rings.RR:2835h *= 22836return h28372838archimedian_local_height = deprecated_function_alias(13951, archimedean_local_height)28392840def non_archimedean_local_height(self, v=None, prec=None,2841weighted=False, is_minimal=None):2842"""2843Compute the local height of self at the non-archimedean place `v`.28442845INPUT:28462847- ``self`` -- a point on an elliptic curve over a number field2848`K`.28492850- ``v`` -- a non-archimedean place of `K`, or None (default).2851If `v` is a non-archimedean place, return the local height2852of self at `v`. If `v` is None, return the total2853non-archimedean contribution to the global height.28542855- ``prec`` -- integer, or None (default). The precision of the2856computation. If None, the height is returned symbolically.28572858- ``weighted`` -- boolean. If False (default), the height is2859normalised to be invariant under extension of `K`. If True,2860return this normalised height multiplied by the local degree2861if `v` is a single place, or by the degree of `K` if `v` is2862None.28632864OUTPUT:28652866A real number. The normalisation is twice that in Silverman's2867paper [Sil1988]_. Note that this local height depends on the2868model of the curve.28692870ALGORITHM:28712872See [Sil1988]_, Section 5.28732874EXAMPLES:28752876Examples 2 and 3 from [Sil1988]_::28772878sage: K.<i> = NumberField(x^2+1)2879sage: E = EllipticCurve(K, [0,0,4,6*i,0]); E2880Elliptic Curve defined by y^2 + 4*y = x^3 + 6*i*x over Number Field in i with defining polynomial x^2 + 12881sage: P = E((0,0))2882sage: P.non_archimedean_local_height(K.ideal(i+1))2883-1/2*log(2)2884sage: P.non_archimedean_local_height(K.ideal(3))288502886sage: P.non_archimedean_local_height(K.ideal(1-2*i))2887028882889sage: Q = E.lift_x(-9/4); Q2890(-9/4 : -27/8*i : 1)2891sage: Q.non_archimedean_local_height(K.ideal(1+i))28922*log(2)2893sage: Q.non_archimedean_local_height(K.ideal(3))289402895sage: Q.non_archimedean_local_height(K.ideal(1-2*i))289602897sage: Q.non_archimedean_local_height()28981/2*log(16)28992900An example over the rational numbers::29012902sage: E = EllipticCurve([0, 0, 0, -36, 0])2903sage: P = E([-3, 9])2904sage: P.non_archimedean_local_height()2905-log(3)29062907Local heights of torsion points can be non-zero (unlike the2908global height)::29092910sage: K.<i> = QuadraticField(-1)2911sage: E = EllipticCurve([0, 0, 0, K(1), 0])2912sage: P = E(i, 0)2913sage: P.non_archimedean_local_height()2914-1/2*log(2)29152916TESTS::29172918sage: Q.non_archimedean_local_height(prec=100)29191.38629436111989061883446424292920sage: (3*Q).non_archimedean_local_height()29211/2*log(75923153929839865104)29222923sage: F.<a> = NumberField(x^4 + 2*x^3 + 19*x^2 + 18*x + 288)2924sage: F.ring_of_integers().gens()2925[1, 5/6*a^3 + 1/6*a, 1/6*a^3 + 1/6*a^2, a^3]2926sage: F.class_number()2927122928sage: E = EllipticCurve('37a').change_ring(F)2929sage: P = E((-a^2/6 - a/6 - 1, a)); P2930(-1/6*a^2 - 1/6*a - 1 : a : 1)2931sage: P[0].is_integral()2932True2933sage: P.non_archimedean_local_height()2934029352936This shows that the bug reported at :trac:`13951` has been fixed::29372938sage: E = EllipticCurve([0,17])2939sage: P = E(2,5)2940sage: P.non_archimedean_local_height(2)2941-2/3*log(2)2942"""2943if prec:2944log = lambda x: rings.RealField(prec)(x).log()2945else:2946from sage.functions.log import log29472948if v is None:2949D = self.curve().discriminant()2950K = self.curve().base_ring()2951if K is rings.QQ:2952factorD = D.factor()2953if self[0] == 0:2954c = 12955else:2956c = self[0].denominator()2957# The last sum is for bad primes that divide c where2958# the model is not minimal.2959h = (log(c)2960+ sum(self.non_archimedean_local_height(p, prec, weighted=True, is_minimal=(e < 12))2961for p,e in factorD if not p.divides(c))2962+ sum(self.non_archimedean_local_height(p, prec, weighted=True)2963- c.valuation(p) * log(p)2964for p,e in factorD if e >= 12 and p.divides(c)))2965else:2966factorD = K.factor(D)2967if self[0] == 0:2968c = K.ideal(1)2969else:2970c = K.ideal(self[0]).denominator()2971# The last sum is for bad primes that divide c where2972# the model is not minimal.2973h = (log(c.norm())2974+ sum(self.non_archimedean_local_height(v, prec, weighted=True, is_minimal=(e < 12))2975for v,e in factorD if not v.divides(c))2976+ sum(self.non_archimedean_local_height(v, prec, weighted=True)2977- c.valuation(v) * log(v.norm())2978for v,e in factorD if e >= 12 and v.divides(c)))2979if not weighted:2980h /= K.degree()2981return h29822983if is_minimal:2984E = self.curve()2985P = self2986offset = 02987else:2988E = self.curve().local_minimal_model(v)2989P = self.curve().isomorphism_to(E)(self)2990# Silverman's normalization is not invariant under change of model,2991# but it all cancels out in the global height.2992offset = (self.curve().discriminant()/E.discriminant()).valuation(v)29932994a1, a2, a3, a4, a6 = E.a_invariants()2995b2, b4, b6, b8 = E.b_invariants()2996c4 = E.c4()2997x, y = P.xy()2998D = E.discriminant()2999N = D.valuation(v)3000A = (3*x**2 + 2*a2*x + a4 - a1*y).valuation(v)3001B = (2*y+a1*x+a3).valuation(v)3002C = (3*x**4 + b2*x**3 + 3*b4*x**2 + 3*b6*x + b8).valuation(v)3003if A <= 0 or B <= 0:3004r = max(0, -x.valuation(v))3005elif c4.valuation(v) == 0:3006n = min(B, N/2)3007r = -n*(N-n)/N3008elif C >= 3*B:3009r = -2*B/33010else:3011r = -C/43012r -= offset/63013if not r:3014return rings.QQ(0)3015else:3016if E.base_ring() is rings.QQ:3017Nv = rings.ZZ(v)3018else:3019Nv = v.norm()3020if not weighted:3021r /= v.ramification_index() * v.residue_class_degree()3022return r * log(Nv)30233024nonarchimedian_local_height = deprecated_function_alias(13951, non_archimedean_local_height)30253026def elliptic_logarithm(self, embedding=None, precision=100,3027algorithm='pari'):3028r"""3029Returns the elliptic logarithm of this elliptic curve point.30303031An embedding of the base field into `\RR` or `\CC` (with3032arbitrary precision) may be given; otherwise the first real3033embedding is used (with the specified precision) if any, else3034the first complex embedding.30353036INPUT:30373038- ``embedding``: an embedding of the base field into `\RR` or `\CC`30393040- ``precision``: a positive integer (default 100) setting the3041number of bits of precision for the computation30423043- ``algorithm``: either 'pari' (default for real embeddings)3044to use PARI's ``ellpointtoz{}``, or 'sage' for a native3045implementation. Ignored for complex embeddings.30463047ALGORITHM:30483049See [Co2] Cohen H., A Course in Computational Algebraic Number3050Theory GTM 138, Springer 1996 for the case of real embeddings,3051and Cremona, J.E. and Thongjunthug , T. 2010 for the complex3052case.30533054AUTHORS:30553056- Michael Mardaus (2008-07),3057- Tobias Nagel (2008-07) -- original version from [Co2].3058- John Cremona (2008-07) -- revision following eclib code.3059- John Cremona (2010-03) -- implementation for complex embeddings.30603061EXAMPLES::30623063sage: E = EllipticCurve('389a')3064sage: E.discriminant() > 03065True3066sage: P = E([-1,1])3067sage: P.is_on_identity_component ()3068False3069sage: P.elliptic_logarithm (precision=96)30700.4793482501902193161295330101 + 0.985868850775824102211203849...*I3071sage: Q=E([3,5])3072sage: Q.is_on_identity_component()3073True3074sage: Q.elliptic_logarithm (precision=96)30751.93112827154255944248858522030763077An example with negative discriminant, and a torsion point::30783079sage: E = EllipticCurve('11a1')3080sage: E.discriminant() < 03081True3082sage: P = E([16,-61])3083sage: P.elliptic_logarithm(precision=70)30840.253841860855910684343085sage: E.period_lattice().real_period(prec=70) / P.elliptic_logarithm(precision=70)30865.000000000000000000030873088A larger example. The default algorithm uses PARI and makes3089sure the result has the requested precision::30903091sage: E = EllipticCurve([1, 0, 1, -85357462, 303528987048]) #18074g13092sage: P = E([4458713781401/835903744, -64466909836503771/24167649046528, 1])3093sage: P.elliptic_logarithm() # 100 bits30940.2765620401410706146407620309730953096The native algorithm 'sage' used to have trouble with3097precision in this example, but no longer::30983099sage: P.elliptic_logarithm(algorithm='sage') # 100 bits31000.2765620401410706146407620309731013102This shows that the bug reported at :trac:`4901` has been fixed::31033104sage: E = EllipticCurve("4390c2")3105sage: P = E(683762969925/44944,-565388972095220019/9528128)3106sage: P.elliptic_logarithm()31070.000256387258865202253531989325293108sage: P.elliptic_logarithm(precision=64)31090.0002563872588652022543110sage: P.elliptic_logarithm(precision=65)31110.00025638725886520225353112sage: P.elliptic_logarithm(precision=128)31130.000256387258865202253531989325286664274123114sage: P.elliptic_logarithm(precision=129)31150.000256387258865202253531989325286664274123116sage: P.elliptic_logarithm(precision=256)31170.00025638725886520225353198932528666427411683880083463700150051421280096109363733118sage: P.elliptic_logarithm(precision=257)31190.0002563872588652022535319893252866642741168388008346370015005142128009610936373031203121Examples over number fields::31223123sage: K.<a> = NumberField(x^3-2)3124sage: embs = K.embeddings(CC)3125sage: E = EllipticCurve([0,1,0,a,a])3126sage: Ls = [E.period_lattice(e) for e in embs]3127sage: [L.real_flag for L in Ls]3128[0, 0, -1]3129sage: P = E(-1,0) # order 23130sage: [L.elliptic_logarithm(P) for L in Ls]3131[-1.73964256006716 - 1.07861534489191*I, -0.363756518406398 - 1.50699412135253*I, 1.90726488608927]31323133sage: E = EllipticCurve([-a^2 - a - 1, a^2 + a])3134sage: Ls = [E.period_lattice(e) for e in embs]3135sage: 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 )]3136sage: [[L.elliptic_logarithm(P) for P in pts] for L in Ls]3137[[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]]31383139::31403141sage: K.<i> = QuadraticField(-1)3142sage: E = EllipticCurve([0,0,0,9*i-10,21-i])3143sage: emb = K.embeddings(CC)[1]3144sage: L = E.period_lattice(emb)3145sage: P = E(2-i,4+2*i)3146sage: L.elliptic_logarithm(P,prec=100)31470.70448375537782208460499649302 - 0.79246725643650979858266018068*I31483149"""3150from sage.rings.number_field.number_field import refine_embedding3151from sage.rings.all import RealField, ComplexField, QQ31523153# Check the trivial case:31543155C = ComplexField(precision)3156if self.is_zero():3157return C.zero()31583159# find a suitable embedding if none was supplied:31603161E = self.curve()3162K = E.base_field()3163rational = (K is QQ)3164emb = embedding31653166if emb is None:3167emb = K.embeddings(RealField(precision))3168if len(emb) > 0:3169emb = emb[0]3170else:3171emb = K.embeddings(ComplexField(precision))[0]3172else:3173# Get the precision of the supplied embedding3174prec = emb.codomain().precision()3175# if the precision parameter is greater, refine the embedding:3176if precision > prec:3177emb = refine_embedding(emb, precision)31783179L = E.period_lattice(emb)31803181if algorithm == 'sage' or not is_RealField(emb.codomain):3182return L.elliptic_logarithm(self, precision)31833184if algorithm != 'pari':3185raise ValueError("algorithm must be either 'pari' or 'sage'")31863187# From now on emb() is a real embedding of K into3188# RealField(precision). We interface with the PARI library.31893190x, y = self.xy()3191if rational: # work with exact coordinates3192E_work = E3193pt_pari = pari([x, y])3194else: # use the embedding to get real coordinates3195ai = [emb(a) for a in E.a_invariants()]3196E_work = EllipticCurve(ai) # defined over RR3197pt_pari = pari([emb(x), emb(y)])3198working_prec = precision3199E_pari = E_work.pari_curve(prec=working_prec)3200log_pari = E_pari.ellpointtoz(pt_pari, precision=working_prec)32013202while prec_words_to_bits(log_pari.precision()) < precision:3203# result is not precise enough, re-compute with double3204# precision. if the base field is not QQ, this3205# requires modifying the precision of the embedding,3206# the curve, and the point3207working_prec = 2*working_prec3208if not rational:3209emb = refine_embedding(emb, working_prec)3210ai = [emb(a) for a in E.a_invariants()]3211E_work = EllipticCurve(ai) # defined over RR3212pt_pari = pari([emb(x), emb(y)])3213E_pari = E_work.pari_curve(prec=working_prec)3214log_pari = E_pari.ellpointtoz(pt_pari, precision=working_prec)32153216# normalization step3217r, i = C(log_pari)3218wR, wI = L.basis(prec=precision)3219k = (r/wR).floor()3220if k:3221r -= k*wR3222if self.is_on_identity_component(emb):3223return C(r)3224# Now there are two components and P is on the non-identity one3225return C(r)+C(wI/2)32263227def padic_elliptic_logarithm(self, p, absprec=20):3228r"""3229Computes the `p`-adic elliptic logarithm of this point.32303231INPUT:32323233``p`` - integer: a prime ``absprec`` - integer (default: 20):3234the initial `p`-adic absolute precision of the computation32353236OUTPUT:32373238The `p`-adic elliptic logarithm of self, with precision ``absprec``.32393240AUTHORS:32413242- Tobias Nagel3243- Michael Mardaus3244- John Cremona32453246ALGORITHM:32473248For points in the formal group (i.e. not integral at `p`) we3249take the ``log()`` function from the formal groups module and3250evaluate it at `-x/y`. Otherwise we first multiply the point3251to get into the formal group, and divide the result3252afterwards.32533254.. TODO::32553256See comments at :trac:`4805`. Currently the absolute3257precision of the result may be less than the given value3258of absprec, and error-handling is imperfect.32593260EXAMPLES::32613262sage: E = EllipticCurve([0,1,1,-2,0])3263sage: E(0).padic_elliptic_logarithm(3)326403265sage: P = E(0,0)3266sage: P.padic_elliptic_logarithm(3)32672 + 2*3 + 3^3 + 2*3^7 + 3^8 + 3^9 + 3^11 + 3^15 + 2*3^17 + 3^18 + O(3^19)3268sage: P.padic_elliptic_logarithm(3).lift()32696602575223270sage: P = E(-11/9,28/27)3271sage: [(2*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(20)] # long time (3s)3272[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)]3273sage: [(3*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(12)] # long time (2s)3274[1 + 2 + O(2^19), 3 + 3^20 + O(3^21), 3 + O(5^19), 3 + O(7^19), 3 + O(11^19)]3275sage: [(5*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(12)] # long time (2s)3276[1 + 2^2 + O(2^19), 2 + 3 + O(3^20), 5 + O(5^19), 5 + O(7^19), 5 + O(11^19)]32773278An example which arose during reviewing :trac:`4741`::32793280sage: E = EllipticCurve('794a1')3281sage: P = E(-1,2)3282sage: P.padic_elliptic_logarithm(2) # default precision=2032832^4 + 2^5 + 2^6 + 2^8 + 2^9 + 2^13 + 2^14 + 2^15 + O(2^16)3284sage: P.padic_elliptic_logarithm(2, absprec=30)32852^4 + 2^5 + 2^6 + 2^8 + 2^9 + 2^13 + 2^14 + 2^15 + 2^22 + 2^23 + 2^24 + O(2^26)3286sage: P.padic_elliptic_logarithm(2, absprec=40)32872^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)3288"""3289if not p.is_prime():3290raise ValueError('p must be prime')3291debug = False # True3292if debug:3293print "P=", self, "; p=", p, " with precision ", absprec3294E = self.curve()3295Q_p = Qp(p, absprec)3296if self.has_finite_order():3297return Q_p(0)3298while True:3299try:3300Ep = E.change_ring(Q_p)3301P = Ep(self)3302x, y = P.xy()3303break3304except (PrecisionError, ArithmeticError, ZeroDivisionError):3305absprec *= 23306Q_p = Qp(p, absprec)3307if debug:3308print "x,y=", (x, y)3309f = 1 # f will be such that f*P is in the formal group E^1(Q_p)3310if x.valuation() >= 0: # P is not in E^13311if not self.has_good_reduction(p): # P is not in E^03312n = E.tamagawa_exponent(p) # n*P has good reduction at p3313if debug:3314print "Tamagawa exponent = =", n3315f = n3316P = n*P # lies in E^03317if debug:3318print "P=", P3319try:3320x, y = P.xy()3321except ZeroDivisionError:3322raise ValueError("Insufficient precision in "3323"p-adic_elliptic_logarithm()")3324if debug:3325print "x,y=", (x, y)3326if x.valuation() >= 0: # P is still not in E^13327t = E.local_data(p).bad_reduction_type()3328if t is None:3329m = E.reduction(p).abelian_group().exponent()3330else:3331m = p - t3332if debug:3333print "mod p exponent = =", m3334# now m*(n*P) reduces to the identity mod p, so is3335# in E^1(Q_p)3336f *= m3337P = m*P # lies in E^13338try:3339x, y = P.xy()3340except ZeroDivisionError:3341raise ValueError("Insufficient precision in "3342"p-adic_elliptic_logarithm()")3343if debug:3344print "f=", f3345print "x,y=", (x, y)3346vx = x.valuation()3347vy = y.valuation()3348v = vx-vy3349if not (v > 0 and vx == -2*v and vy == -3*v):3350raise ValueError("Insufficient precision in "3351"p-adic_elliptic_logarithm()")3352try:3353t = -x/y3354except (ZeroDivisionError, PrecisionError):3355raise ValueError("Insufficient precision in "3356"p-adic_elliptic_logarithm()")3357if debug:3358print "t=", t, ", with valuation ", v3359phi = Ep.formal().log(prec=1+absprec//v)3360return phi(t)/f336133623363class EllipticCurvePoint_finite_field(EllipticCurvePoint_field):3364r"""3365Class for elliptic curve points over finite fields.3366"""3367def _magma_init_(self, magma):3368"""3369Return a string representation of self that ``MAGMA`` can3370use for input.33713372EXAMPLE::33733374sage: E = EllipticCurve(GF(17), [1,-1])3375sage: P = E([13, 4])3376sage: P._magma_init_(magma) # optional - magma3377'EllipticCurve([_sage_ref...|GF(17)!0,GF(17)!0,GF(17)!0,GF(17)!1,GF(17)!16])![13,4]'3378"""3379E = self.curve()._magma_init_(magma)3380x, y = self.xy()3381return "%s![%s,%s]" % (E, x, y)33823383def discrete_log(self, Q, ord=None):3384r"""3385Returns discrete log of `Q` with respect to `P` =self.33863387INPUT:33883389- ``Q`` (point) -- another point on the same curve as self.33903391- ``ord`` (integer or ``None`` (default)) -- the order of self.33923393OUTPUT:33943395(integer) -- The discrete log of `Q` with respect to `P`, which is an3396integer `m` with `0\le m<o(P)` such that `mP=Q`, if one3397exists. A ValueError is raised if there is no solution.33983399.. NOTE::34003401The order of self is computed if not supplied.34023403AUTHOR:34043405- John Cremona. Adapted to use generic functions 2008-04-05.34063407EXAMPLE::34083409sage: F = GF(3^6,'a')3410sage: a = F.gen()3411sage: E = EllipticCurve([0,1,1,a,a])3412sage: E.cardinality()34137623414sage: A = E.abelian_group()3415sage: P = A.gen(0).element()3416sage: Q = 400*P3417sage: P.discrete_log(Q)34184003419"""3420if ord is None:3421ord = self.order()3422try:3423return generic.discrete_log(Q, self, ord, operation='+')3424except Exception:3425raise ValueError("ECDLog problem has no solution")34263427def order(self):3428r"""3429Return the order of this point on the elliptic curve.34303431ALGORITHM:34323433Use generic functions from :mod:`sage.groups.generic`. If the3434group order is known, use ``order_from_multiple()``, otherwise3435use ``order_from_bounds()`` with the Hasse bounds for the base3436field. In the latter case, we might find that we have a3437generator for the group, in which case it is cached.34383439We do not cause the group order to be calculated when not3440known, since this function is used in determining the group3441order via computation of several random points and their3442orders.34433444.. NOTE::34453446:meth:`additive_order` is a synonym for :meth:`order`34473448AUTHOR:34493450- John Cremona, 2008-02-10, adapted 2008-04-05 to use generic functions.34513452EXAMPLES::34533454sage: k.<a> = GF(5^5)3455sage: E = EllipticCurve(k,[2,4]); E3456Elliptic Curve defined by y^2 = x^3 + 2*x + 4 over Finite Field in a of size 5^53457sage: P = E(3*a^4 + 3*a , 2*a + 1 )3458sage: P.order()345932273460sage: Q = E(0,2)3461sage: Q.order()346273463sage: Q.additive_order()3464734653466In the next example, the cardinality of E will be computed3467(using SEA) and cached::34683469sage: p=next_prime(2^150)3470sage: E=EllipticCurve(GF(p),[1,1])3471sage: P=E(831623307675610677632782670796608848711856078, 42295786042873366706573292533588638217232964)3472sage: P.order()347314272476927059598810582625452724743006282814483474sage: P.order()==E.cardinality()3475True3476"""3477try:3478return self._order3479except AttributeError:3480pass3481if self.is_zero():3482return rings.Integer(1)3483E = self.curve()3484K = E.base_ring()3485from sage.schemes.plane_curves.projective_curve import Hasse_bounds3486bounds = Hasse_bounds(K.order())34873488try:3489M = E._order3490try:3491plist = E._prime_factors_of_order3492except AttributeError:3493plist = M.prime_divisors()3494E._prime_factors_of_order = plist3495N = generic.order_from_multiple(self, M, plist, operation='+')3496except Exception:3497if K.is_prime_field():3498M = E.cardinality() # computed and cached3499plist = M.prime_divisors()3500E._prime_factors_of_order = plist3501N = generic.order_from_multiple(self, M, plist, operation='+')3502else:3503N = generic.order_from_bounds(self, bounds, operation='+')35043505if 2*N > bounds[1]: # then we have a generator, so cache this3506if not hasattr(E, '_order'):3507E._order = N3508if not hasattr(E, '__abelian_group'):3509E.__abelian_group = AbelianGroup([N]), (self, )35103511self._order = N3512return self._order35133514additive_order = order351535163517