Points on elliptic curves

The base class EllipticCurvePoint currently provides little functionality of its own. Its derived class EllipticCurvePoint_field provides support for points on elliptic curves over general fields. The derived classes EllipticCurvePoint_number_field and EllipticCurvePoint_finite_field provide further support for points on curves over number fields (including the rational field \(\QQ\)) and over finite fields.

EXAMPLES:

An example over \(\QQ\):

sage: E = EllipticCurve('389a1')
sage: P = E(-1,1); P
(-1 : 1 : 1)
sage: Q = E(0,-1); Q
(0 : -1 : 1)
sage: P+Q
(4 : 8 : 1)
sage: P-Q
(1 : 0 : 1)
sage: 3*P-5*Q
(328/361 : -2800/6859 : 1)
>>> from sage.all import *
>>> E = EllipticCurve('389a1')
>>> P = E(-Integer(1),Integer(1)); P
(-1 : 1 : 1)
>>> Q = E(Integer(0),-Integer(1)); Q
(0 : -1 : 1)
>>> P+Q
(4 : 8 : 1)
>>> P-Q
(1 : 0 : 1)
>>> Integer(3)*P-Integer(5)*Q
(328/361 : -2800/6859 : 1)

An example over a number field:

sage: # needs sage.rings.number_field
sage: K.<i> = QuadraticField(-1)
sage: E = EllipticCurve(K, [1,0,0,0,-1])
sage: P = E(0,i); P
(0 : i : 1)
sage: P.order()
+Infinity
sage: 101*P - 100*P == P
True
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> K = QuadraticField(-Integer(1), names=('i',)); (i,) = K._first_ngens(1)
>>> E = EllipticCurve(K, [Integer(1),Integer(0),Integer(0),Integer(0),-Integer(1)])
>>> P = E(Integer(0),i); P
(0 : i : 1)
>>> P.order()
+Infinity
>>> Integer(101)*P - Integer(100)*P == P
True

An example over a finite field:

sage: # needs sage.rings.finite_rings
sage: K.<a> = GF((101,3))
sage: E = EllipticCurve(K, [1,0,0,0,-1])
sage: P = E(40*a^2 + 69*a + 84 , 58*a^2 + 73*a + 45)
sage: P.order()
1032210
sage: E.cardinality()
1032210
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> K = GF((Integer(101),Integer(3)), names=('a',)); (a,) = K._first_ngens(1)
>>> E = EllipticCurve(K, [Integer(1),Integer(0),Integer(0),Integer(0),-Integer(1)])
>>> P = E(Integer(40)*a**Integer(2) + Integer(69)*a + Integer(84) , Integer(58)*a**Integer(2) + Integer(73)*a + Integer(45))
>>> P.order()
1032210
>>> E.cardinality()
1032210

Arithmetic with a point over an extension of a finite field:

sage: # needs sage.rings.finite_rings
sage: k.<a> = GF((5,2))
sage: E = EllipticCurve(k,[1,0]); E
Elliptic Curve defined by y^2 = x^3 + x over Finite Field in a of size 5^2
sage: P = E([a,2*a+4])
sage: 5*P
(2*a + 3 : 2*a : 1)
sage: P*5
(2*a + 3 : 2*a : 1)
sage: P + P + P + P + P
(2*a + 3 : 2*a : 1)
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> k = GF((Integer(5),Integer(2)), names=('a',)); (a,) = k._first_ngens(1)
>>> E = EllipticCurve(k,[Integer(1),Integer(0)]); E
Elliptic Curve defined by y^2 = x^3 + x over Finite Field in a of size 5^2
>>> P = E([a,Integer(2)*a+Integer(4)])
>>> Integer(5)*P
(2*a + 3 : 2*a : 1)
>>> P*Integer(5)
(2*a + 3 : 2*a : 1)
>>> P + P + P + P + P
(2*a + 3 : 2*a : 1)

sage: F = Zmod(3)
sage: E = EllipticCurve(F, [1,0]);
sage: P = E([2,1])
sage: import sys
sage: n = sys.maxsize
sage: P*(n+1)-P*n == P
True
>>> from sage.all import *
>>> F = Zmod(Integer(3))
>>> E = EllipticCurve(F, [Integer(1),Integer(0)]);
>>> P = E([Integer(2),Integer(1)])
>>> import sys
>>> n = sys.maxsize
>>> P*(n+Integer(1))-P*n == P
True

Arithmetic over \(\ZZ/N\ZZ\) with composite \(N\) is supported. When an operation tries to invert a non-invertible element, a ZeroDivisionError is raised and a factorization of the modulus appears in the error message:

sage: N = 1715761513
sage: E = EllipticCurve(Integers(N), [3,-13])
sage: P = E(2,1)
sage: LCM([2..60])*P
Traceback (most recent call last):
...
ZeroDivisionError: Inverse of 26927 does not exist
(characteristic = 1715761513 = 26927*63719)
>>> from sage.all import *
>>> N = Integer(1715761513)
>>> E = EllipticCurve(Integers(N), [Integer(3),-Integer(13)])
>>> P = E(Integer(2),Integer(1))
>>> LCM((ellipsis_range(Integer(2),Ellipsis,Integer(60))))*P
Traceback (most recent call last):
...
ZeroDivisionError: Inverse of 26927 does not exist
(characteristic = 1715761513 = 26927*63719)

AUTHORS:

  • William Stein (2005) – Initial version

  • Robert Bradshaw et al….

  • John Cremona (Feb 2008) – Point counting and group structure for non-prime fields, Frobenius endomorphism and order, elliptic logs

  • John Cremona (Aug 2008) – Introduced EllipticCurvePoint_number_field class

  • Tobias Nagel, Michael Mardaus, John Cremona (Dec 2008) – \(p\)-adic elliptic logarithm over \(\QQ\)

  • David Hansen (Jan 2009) – Added weil_pairing function to EllipticCurvePoint_finite_field class

  • Mariah Lenox (March 2011) – Added tate_pairing and ate_pairing functions to EllipticCurvePoint_finite_field class

class sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint(X, v, check=True)[source]

Bases: AdditiveGroupElement, SchemeMorphism_point_projective_ring

A point on an elliptic curve.

curve()[source]

Return the curve that this point is on.

This is a synonym for scheme().

EXAMPLES:

sage: E = EllipticCurve('389a')
sage: P = E([-1, 1])
sage: P.curve()
Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field

sage: E = EllipticCurve(QQ, [1, 1])
sage: P = E(0, 1)
sage: P.scheme()
Elliptic Curve defined by y^2 = x^3 + x + 1 over Rational Field
sage: P.scheme() == P.curve()
True
sage: x = polygen(ZZ, 'x')
sage: K.<a> = NumberField(x^2 - 3,'a')                                      # needs sage.rings.number_field
sage: P = E.base_extend(K)(1, a)                                            # needs sage.rings.number_field
sage: P.scheme()                                                            # needs sage.rings.number_field
Elliptic Curve defined by y^2 = x^3 + x + 1 over
 Number Field in a with defining polynomial x^2 - 3
>>> from sage.all import *
>>> E = EllipticCurve('389a')
>>> P = E([-Integer(1), Integer(1)])
>>> P.curve()
Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field

>>> E = EllipticCurve(QQ, [Integer(1), Integer(1)])
>>> P = E(Integer(0), Integer(1))
>>> P.scheme()
Elliptic Curve defined by y^2 = x^3 + x + 1 over Rational Field
>>> P.scheme() == P.curve()
True
>>> x = polygen(ZZ, 'x')
>>> K = NumberField(x**Integer(2) - Integer(3),'a', names=('a',)); (a,) = K._first_ngens(1)# needs sage.rings.number_field
>>> P = E.base_extend(K)(Integer(1), a)                                            # needs sage.rings.number_field
>>> P.scheme()                                                            # needs sage.rings.number_field
Elliptic Curve defined by y^2 = x^3 + x + 1 over
 Number Field in a with defining polynomial x^2 - 3
class sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_field(curve, v, check=True)[source]

Bases: EllipticCurvePoint, SchemeMorphism_point_abelian_variety_field

A point on an elliptic curve over a field. The point has coordinates in the base field.

EXAMPLES:

sage: E = EllipticCurve('37a')
sage: E([0,0])
(0 : 0 : 1)
sage: E(0,0)               # brackets are optional
(0 : 0 : 1)
sage: E([GF(5)(0), 0])     # entries are coerced
(0 : 0 : 1)

sage: E(0.000, 0)
(0 : 0 : 1)

sage: E(1,0,0)
Traceback (most recent call last):
...
TypeError: Coordinates [1, 0, 0] do not define a point on
Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
>>> from sage.all import *
>>> E = EllipticCurve('37a')
>>> E([Integer(0),Integer(0)])
(0 : 0 : 1)
>>> E(Integer(0),Integer(0))               # brackets are optional
(0 : 0 : 1)
>>> E([GF(Integer(5))(Integer(0)), Integer(0)])     # entries are coerced
(0 : 0 : 1)

>>> E(RealNumber('0.000'), Integer(0))
(0 : 0 : 1)

>>> E(Integer(1),Integer(0),Integer(0))
Traceback (most recent call last):
...
TypeError: Coordinates [1, 0, 0] do not define a point on
Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field

sage: E = EllipticCurve([0,0,1,-1,0])
sage: S = E(QQ); S
Abelian group of points on
 Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field

sage: # needs sage.rings.number_field
sage: x = polygen(ZZ, 'x')
sage: K.<i> = NumberField(x^2 + 1)
sage: E = EllipticCurve(K, [0,1,0,-160,308])
sage: P = E(26, -120)
sage: Q = E(2+12*i, -36+48*i)
sage: P.order() == Q.order() == 4       # long time
True
sage: 2*P == 2*Q
False
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0),Integer(0),Integer(1),-Integer(1),Integer(0)])
>>> S = E(QQ); S
Abelian group of points on
 Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field

>>> # needs sage.rings.number_field
>>> x = polygen(ZZ, 'x')
>>> K = NumberField(x**Integer(2) + Integer(1), names=('i',)); (i,) = K._first_ngens(1)
>>> E = EllipticCurve(K, [Integer(0),Integer(1),Integer(0),-Integer(160),Integer(308)])
>>> P = E(Integer(26), -Integer(120))
>>> Q = E(Integer(2)+Integer(12)*i, -Integer(36)+Integer(48)*i)
>>> P.order() == Q.order() == Integer(4)       # long time
True
>>> Integer(2)*P == Integer(2)*Q
False

sage: K.<t> = FractionField(PolynomialRing(QQ,'t'))
sage: E = EllipticCurve([0,0,0,0,t^2])
sage: P = E(0,t)
sage: P, 2*P, 3*P
((0 : t : 1), (0 : -t : 1), (0 : 1 : 0))
>>> from sage.all import *
>>> K = FractionField(PolynomialRing(QQ,'t'), names=('t',)); (t,) = K._first_ngens(1)
>>> E = EllipticCurve([Integer(0),Integer(0),Integer(0),Integer(0),t**Integer(2)])
>>> P = E(Integer(0),t)
>>> P, Integer(2)*P, Integer(3)*P
((0 : t : 1), (0 : -t : 1), (0 : 1 : 0))
additive_order()[source]

Return the order of this point on the elliptic curve.

If the point is zero, returns 1, otherwise raise a NotImplementedError.

For curves over number fields and finite fields, see below.

Note

additive_order() is a synonym for order()

EXAMPLES:

sage: K.<t> = FractionField(PolynomialRing(QQ,'t'))
sage: E = EllipticCurve([0, 0, 0, -t^2, 0])
sage: P = E(t,0)
sage: P.order()
Traceback (most recent call last):
...
NotImplementedError: Computation of order of a point not implemented
over general fields.
sage: E(0).additive_order()
1
sage: E(0).order() == 1
True
>>> from sage.all import *
>>> K = FractionField(PolynomialRing(QQ,'t'), names=('t',)); (t,) = K._first_ngens(1)
>>> E = EllipticCurve([Integer(0), Integer(0), Integer(0), -t**Integer(2), Integer(0)])
>>> P = E(t,Integer(0))
>>> P.order()
Traceback (most recent call last):
...
NotImplementedError: Computation of order of a point not implemented
over general fields.
>>> E(Integer(0)).additive_order()
1
>>> E(Integer(0)).order() == Integer(1)
True
ate_pairing(Q, n, k, t, q=None)[source]

Return ate pairing of \(n\)-torsion points \(P\) (=self) and \(Q\).

Also known as the \(n\)-th modified ate pairing. \(P\) is \(GF(q)\)-rational, and \(Q\) must be an element of \(Ker(\pi-p)\), where \(\pi\) is the \(q\)-frobenius map (and hence \(Q\) is \(GF(q^k)\)-rational).

INPUT:

  • P (=self) – a point of order \(n\), in \(ker(\pi-1)\), where \(\pi\) is the \(q\)-Frobenius map (e.g., \(P\) is \(q\)-rational)

  • Q – a point of order \(n\) in \(ker(\pi-q)\)

  • n – the order of \(P\) and \(Q\)

  • k – the embedding degree

  • t – the trace of Frobenius of the curve over \(GF(q)\)

  • q – (default: None) the size of base field (the “big” field is \(GF(q^k)\)). \(q\) needs to be set only if its value cannot be deduced.

OUTPUT:

FiniteFieldElement in \(GF(q^k)\) – the ate pairing of \(P\) and \(Q\).

EXAMPLES:

An example with embedding degree 6:

sage: # needs sage.rings.finite_rings
sage: p = 7549; A = 0; B = 1; n = 157; k = 6; t = 14
sage: F = GF(p); E = EllipticCurve(F, [A, B])
sage: R.<x> = F[]; K.<a> = GF((p,k), modulus=x^k+2)
sage: EK = E.base_extend(K)
sage: P = EK(3050, 5371); Q = EK(6908*a^4, 3231*a^3)
sage: P.ate_pairing(Q, n, k, t)
6708*a^5 + 4230*a^4 + 4350*a^3 + 2064*a^2 + 4022*a + 6733
sage: s = Integer(randrange(1, n))
sage: (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)
True
sage: P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)^s
True
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> p = Integer(7549); A = Integer(0); B = Integer(1); n = Integer(157); k = Integer(6); t = Integer(14)
>>> F = GF(p); E = EllipticCurve(F, [A, B])
>>> R = F['x']; (x,) = R._first_ngens(1); K = GF((p,k), modulus=x**k+Integer(2), names=('a',)); (a,) = K._first_ngens(1)
>>> EK = E.base_extend(K)
>>> P = EK(Integer(3050), Integer(5371)); Q = EK(Integer(6908)*a**Integer(4), Integer(3231)*a**Integer(3))
>>> P.ate_pairing(Q, n, k, t)
6708*a^5 + 4230*a^4 + 4350*a^3 + 2064*a^2 + 4022*a + 6733
>>> s = Integer(randrange(Integer(1), n))
>>> (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)
True
>>> P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)**s
True

Another example with embedding degree 7 and positive trace:

sage: # needs sage.rings.finite_rings
sage: p = 2213; A = 1; B = 49; n = 1093; k = 7; t = 28
sage: F = GF(p); E = EllipticCurve(F, [A, B])
sage: R.<x> = F[]; K.<a> = GF((p,k), modulus=x^k+2)
sage: EK = E.base_extend(K)
sage: P = EK(1583, 1734)
sage: Qx = 1729*a^6+1767*a^5+245*a^4+980*a^3+1592*a^2+1883*a+722
sage: Qy = 1299*a^6+1877*a^5+1030*a^4+1513*a^3+1457*a^2+309*a+1636
sage: Q = EK(Qx, Qy)
sage: P.ate_pairing(Q, n, k, t)
1665*a^6 + 1538*a^5 + 1979*a^4 + 239*a^3 + 2134*a^2 + 2151*a + 654
sage: s = Integer(randrange(1, n))
sage: (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)
True
sage: P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)^s
True
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> p = Integer(2213); A = Integer(1); B = Integer(49); n = Integer(1093); k = Integer(7); t = Integer(28)
>>> F = GF(p); E = EllipticCurve(F, [A, B])
>>> R = F['x']; (x,) = R._first_ngens(1); K = GF((p,k), modulus=x**k+Integer(2), names=('a',)); (a,) = K._first_ngens(1)
>>> EK = E.base_extend(K)
>>> P = EK(Integer(1583), Integer(1734))
>>> Qx = Integer(1729)*a**Integer(6)+Integer(1767)*a**Integer(5)+Integer(245)*a**Integer(4)+Integer(980)*a**Integer(3)+Integer(1592)*a**Integer(2)+Integer(1883)*a+Integer(722)
>>> Qy = Integer(1299)*a**Integer(6)+Integer(1877)*a**Integer(5)+Integer(1030)*a**Integer(4)+Integer(1513)*a**Integer(3)+Integer(1457)*a**Integer(2)+Integer(309)*a+Integer(1636)
>>> Q = EK(Qx, Qy)
>>> P.ate_pairing(Q, n, k, t)
1665*a^6 + 1538*a^5 + 1979*a^4 + 239*a^3 + 2134*a^2 + 2151*a + 654
>>> s = Integer(randrange(Integer(1), n))
>>> (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)
True
>>> P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)**s
True

Another example with embedding degree 7 and negative trace:

sage: # needs sage.rings.finite_rings
sage: p = 2017; A = 1; B = 30; n = 29; k = 7; t = -70
sage: F = GF(p); E = EllipticCurve(F, [A, B])
sage: R.<x> = F[]; K.<a> = GF((p,k), modulus=x^k+2)
sage: EK = E.base_extend(K)
sage: P = EK(369, 716)
sage: Qx = 1226*a^6+1778*a^5+660*a^4+1791*a^3+1750*a^2+867*a+770
sage: Qy = 1764*a^6+198*a^5+1206*a^4+406*a^3+1200*a^2+273*a+1712
sage: Q = EK(Qx, Qy)
sage: P.ate_pairing(Q, n, k, t)
1794*a^6 + 1161*a^5 + 576*a^4 + 488*a^3 + 1950*a^2 + 1905*a + 1315
sage: s = Integer(randrange(1, n))
sage: (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)
True
sage: P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)^s
True
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> p = Integer(2017); A = Integer(1); B = Integer(30); n = Integer(29); k = Integer(7); t = -Integer(70)
>>> F = GF(p); E = EllipticCurve(F, [A, B])
>>> R = F['x']; (x,) = R._first_ngens(1); K = GF((p,k), modulus=x**k+Integer(2), names=('a',)); (a,) = K._first_ngens(1)
>>> EK = E.base_extend(K)
>>> P = EK(Integer(369), Integer(716))
>>> Qx = Integer(1226)*a**Integer(6)+Integer(1778)*a**Integer(5)+Integer(660)*a**Integer(4)+Integer(1791)*a**Integer(3)+Integer(1750)*a**Integer(2)+Integer(867)*a+Integer(770)
>>> Qy = Integer(1764)*a**Integer(6)+Integer(198)*a**Integer(5)+Integer(1206)*a**Integer(4)+Integer(406)*a**Integer(3)+Integer(1200)*a**Integer(2)+Integer(273)*a+Integer(1712)
>>> Q = EK(Qx, Qy)
>>> P.ate_pairing(Q, n, k, t)
1794*a^6 + 1161*a^5 + 576*a^4 + 488*a^3 + 1950*a^2 + 1905*a + 1315
>>> s = Integer(randrange(Integer(1), n))
>>> (s*P).ate_pairing(Q, n, k, t) == P.ate_pairing(s*Q, n, k, t)
True
>>> P.ate_pairing(s*Q, n, k, t) == P.ate_pairing(Q, n, k, t)**s
True

Using the same data, we show that the ate pairing is a power of the Tate pairing (see [HSV2006] end of section 3.1):

sage: # needs sage.rings.finite_rings
sage: c = (k*p^(k-1)).mod(n); T = t - 1
sage: N = gcd(T^k - 1, p^k - 1)
sage: s = Integer(N/n)
sage: L = Integer((T^k - 1)/N)
sage: M = (L*s*c.inverse_mod(n)).mod(n)
sage: P.ate_pairing(Q, n, k, t) == Q.tate_pairing(P, n, k)^M
True
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> c = (k*p**(k-Integer(1))).mod(n); T = t - Integer(1)
>>> N = gcd(T**k - Integer(1), p**k - Integer(1))
>>> s = Integer(N/n)
>>> L = Integer((T**k - Integer(1))/N)
>>> M = (L*s*c.inverse_mod(n)).mod(n)
>>> P.ate_pairing(Q, n, k, t) == Q.tate_pairing(P, n, k)**M
True

An example where we have to pass the base field size (and we again have agreement with the Tate pairing). Note that though \(Px\) is not \(F\)-rational, (it is the homomorphic image of an \(F\)-rational point) it is nonetheless in \(ker(\pi-1)\), and so is a legitimate input:

sage: # needs sage.rings.finite_rings
sage: q = 2^5; F.<a> = GF(q)
sage: n = 41; k = 4; t = -8
sage: E = EllipticCurve(F,[0,0,1,1,1])
sage: P = E(a^4 + 1, a^3)
sage: Fx.<b> = GF(q^k)
sage: Ex = EllipticCurve(Fx, [0,0,1,1,1])
sage: phi = Hom(F, Fx)(F.gen().minpoly().roots(Fx)[0][0])
sage: Px = Ex(phi(P.x()), phi(P.y()))
sage: 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)
sage: Qx = Ex(Qx[0]^q, Qx[1]^q) - Qx  # ensure Qx is in ker(pi - q)
sage: Px.ate_pairing(Qx, n, k, t)
Traceback (most recent call last):
...
ValueError: Unexpected field degree: set keyword argument q equal to
the size of the base field (big field is GF(q^4)).
sage: Px.ate_pairing(Qx, n, k, t, q)
b^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 + 1
sage: s = Integer(randrange(1, n))
sage: (s*Px).ate_pairing(Qx, n, k, t, q) == Px.ate_pairing(s*Qx, n, k, t, q)
True
sage: Px.ate_pairing(s*Qx, n, k, t, q) == Px.ate_pairing(Qx, n, k, t, q)^s
True
sage: c = (k*q^(k-1)).mod(n); T = t - 1
sage: N = gcd(T^k - 1, q^k - 1)
sage: s = Integer(N/n)
sage: L = Integer((T^k - 1)/N)
sage: M = (L*s*c.inverse_mod(n)).mod(n)
sage: Px.ate_pairing(Qx, n, k, t, q) == Qx.tate_pairing(Px, n, k, q)^M
True
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> q = Integer(2)**Integer(5); F = GF(q, names=('a',)); (a,) = F._first_ngens(1)
>>> n = Integer(41); k = Integer(4); t = -Integer(8)
>>> E = EllipticCurve(F,[Integer(0),Integer(0),Integer(1),Integer(1),Integer(1)])
>>> P = E(a**Integer(4) + Integer(1), a**Integer(3))
>>> Fx = GF(q**k, names=('b',)); (b,) = Fx._first_ngens(1)
>>> Ex = EllipticCurve(Fx, [Integer(0),Integer(0),Integer(1),Integer(1),Integer(1)])
>>> phi = Hom(F, Fx)(F.gen().minpoly().roots(Fx)[Integer(0)][Integer(0)])
>>> Px = Ex(phi(P.x()), phi(P.y()))
>>> Qx = Ex(b**Integer(19)+b**Integer(18)+b**Integer(16)+b**Integer(12)+b**Integer(10)+b**Integer(9)+b**Integer(8)+b**Integer(5)+b**Integer(3)+Integer(1),
...         b**Integer(18)+b**Integer(13)+b**Integer(10)+b**Integer(8)+b**Integer(5)+b**Integer(4)+b**Integer(3)+b)
>>> Qx = Ex(Qx[Integer(0)]**q, Qx[Integer(1)]**q) - Qx  # ensure Qx is in ker(pi - q)
>>> Px.ate_pairing(Qx, n, k, t)
Traceback (most recent call last):
...
ValueError: Unexpected field degree: set keyword argument q equal to
the size of the base field (big field is GF(q^4)).
>>> Px.ate_pairing(Qx, n, k, t, q)
b^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 + 1
>>> s = Integer(randrange(Integer(1), n))
>>> (s*Px).ate_pairing(Qx, n, k, t, q) == Px.ate_pairing(s*Qx, n, k, t, q)
True
>>> Px.ate_pairing(s*Qx, n, k, t, q) == Px.ate_pairing(Qx, n, k, t, q)**s
True
>>> c = (k*q**(k-Integer(1))).mod(n); T = t - Integer(1)
>>> N = gcd(T**k - Integer(1), q**k - Integer(1))
>>> s = Integer(N/n)
>>> L = Integer((T**k - Integer(1))/N)
>>> M = (L*s*c.inverse_mod(n)).mod(n)
>>> Px.ate_pairing(Qx, n, k, t, q) == Qx.tate_pairing(Px, n, k, q)**M
True

It is an error if \(Q\) is not in the kernel of \(\pi - p\), where \(\pi\) is the Frobenius automorphism:

sage: # needs sage.rings.finite_rings
sage: p = 29; A = 1; B = 0; n = 5; k = 2; t = 10
sage: F = GF(p); R.<x> = F[]
sage: E = EllipticCurve(F, [A, B]);
sage: K.<a> = GF((p,k), modulus=x^k+2); EK = E.base_extend(K)
sage: P = EK(13, 8); Q = EK(13, 21)
sage: P.ate_pairing(Q, n, k, t)
Traceback (most recent call last):
...
ValueError: Point (13 : 21 : 1) not in Ker(pi - q)
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> p = Integer(29); A = Integer(1); B = Integer(0); n = Integer(5); k = Integer(2); t = Integer(10)
>>> F = GF(p); R = F['x']; (x,) = R._first_ngens(1)
>>> E = EllipticCurve(F, [A, B]);
>>> K = GF((p,k), modulus=x**k+Integer(2), names=('a',)); (a,) = K._first_ngens(1); EK = E.base_extend(K)
>>> P = EK(Integer(13), Integer(8)); Q = EK(Integer(13), Integer(21))
>>> P.ate_pairing(Q, n, k, t)
Traceback (most recent call last):
...
ValueError: Point (13 : 21 : 1) not in Ker(pi - q)

It is also an error if \(P\) is not in the kernel os \(\pi - 1\):

sage: # needs sage.rings.finite_rings
sage: p = 29; A = 1; B = 0; n = 5; k = 2; t = 10
sage: F = GF(p); R.<x> = F[]
sage: E = EllipticCurve(F, [A, B]);
sage: K.<a> = GF((p,k), modulus=x^k+2); EK = E.base_extend(K)
sage: P = EK(14, 10*a); Q = EK(13, 21)
sage: P.ate_pairing(Q, n, k, t)
Traceback (most recent call last):
...
ValueError: This point (14 : 10*a : 1) is not in Ker(pi - 1)
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> p = Integer(29); A = Integer(1); B = Integer(0); n = Integer(5); k = Integer(2); t = Integer(10)
>>> F = GF(p); R = F['x']; (x,) = R._first_ngens(1)
>>> E = EllipticCurve(F, [A, B]);
>>> K = GF((p,k), modulus=x**k+Integer(2), names=('a',)); (a,) = K._first_ngens(1); EK = E.base_extend(K)
>>> P = EK(Integer(14), Integer(10)*a); Q = EK(Integer(13), Integer(21))
>>> P.ate_pairing(Q, n, k, t)
Traceback (most recent call last):
...
ValueError: This point (14 : 10*a : 1) is not in Ker(pi - 1)

Note

First defined in the paper of [HSV2006], the ate pairing can be computationally effective in those cases when the trace of the curve over the base field is significantly smaller than the expected value. This implementation is simply Miller’s algorithm followed by a naive exponentiation, and makes no claims towards efficiency.

AUTHORS:

  • Mariah Lenox (2011-03-08)

division_points(m, poly_only=False)[source]

Return a list of all points \(Q\) such that \(mQ=P\) where \(P\) = self.

Only points on the elliptic curve containing self and defined over the base field are included.

INPUT:

  • m – positive integer

  • poly_only – boolean (default: False); if True return polynomial whose roots give all possible \(x\)-coordinates of \(m\)-th roots of self

OUTPUT: a (possibly empty) list of solutions \(Q\) to \(mQ=P\), where \(P\) = self

EXAMPLES:

We find the five 5-torsion points on an elliptic curve:

sage: E = EllipticCurve('11a'); E
Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field
sage: P = E(0); P
(0 : 1 : 0)
sage: P.division_points(5)
[(0 : 1 : 0), (5 : -6 : 1), (5 : 5 : 1), (16 : -61 : 1), (16 : 60 : 1)]
>>> from sage.all import *
>>> E = EllipticCurve('11a'); E
Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field
>>> P = E(Integer(0)); P
(0 : 1 : 0)
>>> P.division_points(Integer(5))
[(0 : 1 : 0), (5 : -6 : 1), (5 : 5 : 1), (16 : -61 : 1), (16 : 60 : 1)]

Note above that 0 is included since [5]*0 = 0.

We create a curve of rank 1 with no torsion and do a consistency check:

sage: E = EllipticCurve('11a').quadratic_twist(-7)
sage: Q = E([44,-270])
sage: (4*Q).division_points(4)
[(44 : -270 : 1)]
>>> from sage.all import *
>>> E = EllipticCurve('11a').quadratic_twist(-Integer(7))
>>> Q = E([Integer(44),-Integer(270)])
>>> (Integer(4)*Q).division_points(Integer(4))
[(44 : -270 : 1)]

We create a curve over a non-prime finite field with group of order \(18\):

sage: # needs sage.rings.finite_rings
sage: k.<a> = GF((5,2))
sage: E = EllipticCurve(k, [1,2+a,3,4*a,2])
sage: P = E([3, 3*a+4])
sage: factor(E.order())
2 * 3^2
sage: P.order()
9
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> k = GF((Integer(5),Integer(2)), names=('a',)); (a,) = k._first_ngens(1)
>>> E = EllipticCurve(k, [Integer(1),Integer(2)+a,Integer(3),Integer(4)*a,Integer(2)])
>>> P = E([Integer(3), Integer(3)*a+Integer(4)])
>>> factor(E.order())
2 * 3^2
>>> P.order()
9

We find the \(1\)-division points as a consistency check – there is just one, of course:

sage: P.division_points(1)                                                  # needs sage.rings.finite_rings
[(3 : 3*a + 4 : 1)]
>>> from sage.all import *
>>> P.division_points(Integer(1))                                                  # needs sage.rings.finite_rings
[(3 : 3*a + 4 : 1)]

The point \(P\) has order coprime to 2 but divisible by 3, so:

sage: P.division_points(2)                                                  # needs sage.rings.finite_rings
[(2*a + 1 : 3*a + 4 : 1), (3*a + 1 : a : 1)]
>>> from sage.all import *
>>> P.division_points(Integer(2))                                                  # needs sage.rings.finite_rings
[(2*a + 1 : 3*a + 4 : 1), (3*a + 1 : a : 1)]

We check that each of the 2-division points works as claimed:

sage: [2*Q for Q in P.division_points(2)]                                   # needs sage.rings.finite_rings
[(3 : 3*a + 4 : 1), (3 : 3*a + 4 : 1)]
>>> from sage.all import *
>>> [Integer(2)*Q for Q in P.division_points(Integer(2))]                                   # needs sage.rings.finite_rings
[(3 : 3*a + 4 : 1), (3 : 3*a + 4 : 1)]

Some other checks:

sage: P.division_points(3)                                                  # needs sage.rings.finite_rings
[]
sage: P.division_points(4)                                                  # needs sage.rings.finite_rings
[(0 : 3*a + 2 : 1), (1 : 0 : 1)]
sage: P.division_points(5)                                                  # needs sage.rings.finite_rings
[(1 : 1 : 1)]
>>> from sage.all import *
>>> P.division_points(Integer(3))                                                  # needs sage.rings.finite_rings
[]
>>> P.division_points(Integer(4))                                                  # needs sage.rings.finite_rings
[(0 : 3*a + 2 : 1), (1 : 0 : 1)]
>>> P.division_points(Integer(5))                                                  # needs sage.rings.finite_rings
[(1 : 1 : 1)]

An example over a number field (see Issue #3383):

sage: # needs sage.rings.number_field
sage: E = EllipticCurve('19a1')
sage: x = polygen(ZZ, 'x')
sage: 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)
sage: EK = E.base_extend(K)
sage: E(0).division_points(3)
[(0 : 1 : 0), (5 : -10 : 1), (5 : 9 : 1)]
sage: EK(0).division_points(3)
[(0 : 1 : 0), (5 : 9 : 1), (5 : -10 : 1)]
sage: E(0).division_points(9)
[(0 : 1 : 0), (5 : -10 : 1), (5 : 9 : 1)]
sage: EK(0).division_points(9)
[(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)]
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> E = EllipticCurve('19a1')
>>> x = polygen(ZZ, 'x')
>>> K = NumberField(x**Integer(9) - Integer(3)*x**Integer(8) - Integer(4)*x**Integer(7) + Integer(16)*x**Integer(6) - Integer(3)*x**Integer(5)
...                     - Integer(21)*x**Integer(4) + Integer(5)*x**Integer(3) + Integer(7)*x**Integer(2) - Integer(7)*x + Integer(1), names=('t',)); (t,) = K._first_ngens(1)
>>> EK = E.base_extend(K)
>>> E(Integer(0)).division_points(Integer(3))
[(0 : 1 : 0), (5 : -10 : 1), (5 : 9 : 1)]
>>> EK(Integer(0)).division_points(Integer(3))
[(0 : 1 : 0), (5 : 9 : 1), (5 : -10 : 1)]
>>> E(Integer(0)).division_points(Integer(9))
[(0 : 1 : 0), (5 : -10 : 1), (5 : 9 : 1)]
>>> EK(Integer(0)).division_points(Integer(9))
[(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)]
has_finite_order()[source]

Return True if this point has finite additive order as an element of the group of points on this curve.

For fields other than number fields and finite fields, this is NotImplemented unless self.is_zero().

EXAMPLES:

sage: K.<t> = FractionField(PolynomialRing(QQ,'t'))
sage: E = EllipticCurve([0, 0, 0, -t^2, 0])
sage: P = E(0)
sage: P.has_finite_order()
True
sage: P = E(t,0)
sage: P.has_finite_order()
Traceback (most recent call last):
...
NotImplementedError: Computation of order of a point not implemented
over general fields.
sage: (2*P).is_zero()
True
>>> from sage.all import *
>>> K = FractionField(PolynomialRing(QQ,'t'), names=('t',)); (t,) = K._first_ngens(1)
>>> E = EllipticCurve([Integer(0), Integer(0), Integer(0), -t**Integer(2), Integer(0)])
>>> P = E(Integer(0))
>>> P.has_finite_order()
True
>>> P = E(t,Integer(0))
>>> P.has_finite_order()
Traceback (most recent call last):
...
NotImplementedError: Computation of order of a point not implemented
over general fields.
>>> (Integer(2)*P).is_zero()
True
has_infinite_order()[source]

Return True if this point has infinite additive order as an element of the group of points on this curve.

For fields other than number fields and finite fields, this is NotImplemented unless self.is_zero().

EXAMPLES:

sage: K.<t> = FractionField(PolynomialRing(QQ,'t'))
sage: E = EllipticCurve([0, 0, 0, -t^2, 0])
sage: P = E(0)
sage: P.has_infinite_order()
False
sage: P = E(t,0)
sage: P.has_infinite_order()
Traceback (most recent call last):
...
NotImplementedError: Computation of order of a point not implemented over general fields.
sage: (2*P).is_zero()
True
>>> from sage.all import *
>>> K = FractionField(PolynomialRing(QQ,'t'), names=('t',)); (t,) = K._first_ngens(1)
>>> E = EllipticCurve([Integer(0), Integer(0), Integer(0), -t**Integer(2), Integer(0)])
>>> P = E(Integer(0))
>>> P.has_infinite_order()
False
>>> P = E(t,Integer(0))
>>> P.has_infinite_order()
Traceback (most recent call last):
...
NotImplementedError: Computation of order of a point not implemented over general fields.
>>> (Integer(2)*P).is_zero()
True
has_order(n)[source]

Test if this point has order exactly \(n\).

INPUT:

ALGORITHM:

Compare a cached order if available, otherwise use sage.groups.generic.has_order().

EXAMPLES:

sage: E = EllipticCurve('26b1')
sage: P = E(1, 0)
sage: P.has_order(7)
True
sage: P._order
7
sage: P.has_order(7)
True
>>> from sage.all import *
>>> E = EllipticCurve('26b1')
>>> P = E(Integer(1), Integer(0))
>>> P.has_order(Integer(7))
True
>>> P._order
7
>>> P.has_order(Integer(7))
True

It also works with a Factorization object:

sage: E = EllipticCurve(GF(419), [1,0])
sage: P = E(-33, 8)
sage: P.has_order(factor(21))
True
sage: P._order
21
sage: P.has_order(factor(21))
True
>>> from sage.all import *
>>> E = EllipticCurve(GF(Integer(419)), [Integer(1),Integer(0)])
>>> P = E(-Integer(33), Integer(8))
>>> P.has_order(factor(Integer(21)))
True
>>> P._order
21
>>> P.has_order(factor(Integer(21)))
True

This method can be much faster than computing the order and comparing:

sage: # not tested -- timings are different each time
sage: p = 4 * prod(primes(3,377)) * 587 - 1
sage: E = EllipticCurve(GF(p), [1,0])
sage: %timeit P = E.random_point(); P.set_order(multiple=p+1)
72.4 ms ± 773 µs per loop (mean ± std. dev. of 7 runs, 1 loop each)
sage: %timeit P = E.random_point(); P.has_order(p+1)
32.8 ms ± 3.12 ms per loop (mean ± std. dev. of 7 runs, 10 loops each)
sage: fac = factor(p+1)
sage: %timeit P = E.random_point(); P.has_order(fac)
30.6 ms ± 3.48 ms per loop (mean ± std. dev. of 7 runs, 10 loops each)
>>> from sage.all import *
>>> # not tested -- timings are different each time
>>> p = Integer(4) * prod(primes(Integer(3),Integer(377))) * Integer(587) - Integer(1)
>>> E = EllipticCurve(GF(p), [Integer(1),Integer(0)])
>>> %timeit P = E.random_point(); P.set_order(multiple=p+Integer(1))
72.4 ms ± 773 µs per loop (mean ± std. dev. of 7 runs, 1 loop each)
>>> %timeit P = E.random_point(); P.has_order(p+Integer(1))
32.8 ms ± 3.12 ms per loop (mean ± std. dev. of 7 runs, 10 loops each)
>>> fac = factor(p+Integer(1))
>>> %timeit P = E.random_point(); P.has_order(fac)
30.6 ms ± 3.48 ms per loop (mean ± std. dev. of 7 runs, 10 loops each)

The order is cached once it has been confirmed once, and the cache is shared with order():

sage: # not tested -- timings are different each time
sage: P, = E.gens()
sage: delattr(P, '_order')
sage: %time P.has_order(p+1)
CPU times: user 83.6 ms, sys: 30 µs, total: 83.6 ms
Wall time: 83.8 ms
True
sage: %time P.has_order(p+1)
CPU times: user 31 µs, sys: 2 µs, total: 33 µs
Wall time: 37.9 µs
True
sage: %time P.order()
CPU times: user 11 µs, sys: 0 ns, total: 11 µs
Wall time: 16 µs
5326738796327623094747867617954605554069371494832722337612446642054009560026576537626892113026381253624626941643949444792662881241621373288942880288065660
sage: delattr(P, '_order')
sage: %time P.has_order(fac)
CPU times: user 68.6 ms, sys: 17 µs, total: 68.7 ms
Wall time: 68.7 ms
True
sage: %time P.has_order(fac)
CPU times: user 92 µs, sys: 0 ns, total: 92 µs
Wall time: 97.5 µs
True
sage: %time P.order()
CPU times: user 10 µs, sys: 1e+03 ns, total: 11 µs
Wall time: 14.5 µs
5326738796327623094747867617954605554069371494832722337612446642054009560026576537626892113026381253624626941643949444792662881241621373288942880288065660
>>> from sage.all import *
>>> # not tested -- timings are different each time
>>> P, = E.gens()
>>> delattr(P, '_order')
>>> %time P.has_order(p+Integer(1))
CPU times: user 83.6 ms, sys: 30 µs, total: 83.6 ms
Wall time: 83.8 ms
True
>>> %time P.has_order(p+Integer(1))
CPU times: user 31 µs, sys: 2 µs, total: 33 µs
Wall time: 37.9 µs
True
>>> %time P.order()
CPU times: user 11 µs, sys: 0 ns, total: 11 µs
Wall time: 16 µs
5326738796327623094747867617954605554069371494832722337612446642054009560026576537626892113026381253624626941643949444792662881241621373288942880288065660
>>> delattr(P, '_order')
>>> %time P.has_order(fac)
CPU times: user 68.6 ms, sys: 17 µs, total: 68.7 ms
Wall time: 68.7 ms
True
>>> %time P.has_order(fac)
CPU times: user 92 µs, sys: 0 ns, total: 92 µs
Wall time: 97.5 µs
True
>>> %time P.order()
CPU times: user 10 µs, sys: 1e+03 ns, total: 11 µs
Wall time: 14.5 µs
5326738796327623094747867617954605554069371494832722337612446642054009560026576537626892113026381253624626941643949444792662881241621373288942880288065660
is_divisible_by(m)[source]

Return True if there exists a point \(Q\) defined over the same field as self such that \(mQ\) == self.

INPUT:

  • m – positive integer

OUTPUT:

boolean; True if there is a solution, else False.

Warning

This function usually triggers the computation of the \(m\)-th division polynomial of the associated elliptic curve, which will be expensive if \(m\) is large, though it will be cached for subsequent calls with the same \(m\).

EXAMPLES:

sage: E = EllipticCurve('389a')
sage: Q = 5*E(0,0); Q
(-2739/1444 : -77033/54872 : 1)
sage: Q.is_divisible_by(4)
False
sage: Q.is_divisible_by(5)
True
>>> from sage.all import *
>>> E = EllipticCurve('389a')
>>> Q = Integer(5)*E(Integer(0),Integer(0)); Q
(-2739/1444 : -77033/54872 : 1)
>>> Q.is_divisible_by(Integer(4))
False
>>> Q.is_divisible_by(Integer(5))
True

A finite field example:

sage: E = EllipticCurve(GF(101), [23,34])
sage: E.cardinality().factor()
2 * 53
sage: Set([T.order() for T in E.points()])
{1, 106, 2, 53}
sage: len([T for T in E.points() if T.is_divisible_by(2)])
53
sage: len([T for T in E.points() if T.is_divisible_by(3)])
106
>>> from sage.all import *
>>> E = EllipticCurve(GF(Integer(101)), [Integer(23),Integer(34)])
>>> E.cardinality().factor()
2 * 53
>>> Set([T.order() for T in E.points()])
{1, 106, 2, 53}
>>> len([T for T in E.points() if T.is_divisible_by(Integer(2))])
53
>>> len([T for T in E.points() if T.is_divisible_by(Integer(3))])
106
is_finite_order()[source]

Return True if this point has finite additive order as an element of the group of points on this curve.

For fields other than number fields and finite fields, this is NotImplemented unless self.is_zero().

EXAMPLES:

sage: K.<t> = FractionField(PolynomialRing(QQ,'t'))
sage: E = EllipticCurve([0, 0, 0, -t^2, 0])
sage: P = E(0)
sage: P.has_finite_order()
True
sage: P = E(t,0)
sage: P.has_finite_order()
Traceback (most recent call last):
...
NotImplementedError: Computation of order of a point not implemented
over general fields.
sage: (2*P).is_zero()
True
>>> from sage.all import *
>>> K = FractionField(PolynomialRing(QQ,'t'), names=('t',)); (t,) = K._first_ngens(1)
>>> E = EllipticCurve([Integer(0), Integer(0), Integer(0), -t**Integer(2), Integer(0)])
>>> P = E(Integer(0))
>>> P.has_finite_order()
True
>>> P = E(t,Integer(0))
>>> P.has_finite_order()
Traceback (most recent call last):
...
NotImplementedError: Computation of order of a point not implemented
over general fields.
>>> (Integer(2)*P).is_zero()
True
order()[source]

Return the order of this point on the elliptic curve.

If the point is zero, returns 1, otherwise raise a NotImplementedError.

For curves over number fields and finite fields, see below.

Note

additive_order() is a synonym for order()

EXAMPLES:

sage: K.<t> = FractionField(PolynomialRing(QQ,'t'))
sage: E = EllipticCurve([0, 0, 0, -t^2, 0])
sage: P = E(t,0)
sage: P.order()
Traceback (most recent call last):
...
NotImplementedError: Computation of order of a point not implemented
over general fields.
sage: E(0).additive_order()
1
sage: E(0).order() == 1
True
>>> from sage.all import *
>>> K = FractionField(PolynomialRing(QQ,'t'), names=('t',)); (t,) = K._first_ngens(1)
>>> E = EllipticCurve([Integer(0), Integer(0), Integer(0), -t**Integer(2), Integer(0)])
>>> P = E(t,Integer(0))
>>> P.order()
Traceback (most recent call last):
...
NotImplementedError: Computation of order of a point not implemented
over general fields.
>>> E(Integer(0)).additive_order()
1
>>> E(Integer(0)).order() == Integer(1)
True
plot(**args)[source]

Plot this point on an elliptic curve.

INPUT:

  • **args – all arguments get passed directly onto the point plotting function

EXAMPLES:

sage: E = EllipticCurve('389a')
sage: P = E([-1,1])
sage: P.plot(pointsize=30, rgbcolor=(1,0,0))                                # needs sage.plot
Graphics object consisting of 1 graphics primitive
>>> from sage.all import *
>>> E = EllipticCurve('389a')
>>> P = E([-Integer(1),Integer(1)])
>>> P.plot(pointsize=Integer(30), rgbcolor=(Integer(1),Integer(0),Integer(0)))                                # needs sage.plot
Graphics object consisting of 1 graphics primitive
point_of_jacobian_of_curve()[source]

Return the point in the Jacobian of the curve.

The Jacobian is the one attached to the projective curve associated with this elliptic curve.

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: k.<a> = GF((5,2))
sage: E = EllipticCurve(k,[1,0]); E
Elliptic Curve defined by y^2 = x^3 + x over Finite Field in a of size 5^2
sage: E.order()
32
sage: P = E([a, 2*a + 4])
sage: P
(a : 2*a + 4 : 1)
sage: P.order()
8
sage: p = P.point_of_jacobian_of_curve()
sage: p
[Place (x + 4*a, y + 3*a + 1)]
sage: p.order()
8
sage: Q = 3*P
sage: q = Q.point_of_jacobian_of_curve()
sage: q == 3*p
True
sage: G = p.parent()
sage: G.order()
32
sage: G
Group of rational points of Jacobian over Finite Field in a of size 5^2 (Hess model)
sage: J = G.parent(); J
Jacobian of Projective Plane Curve over Finite Field in a of size 5^2
 defined by x^2*y + y^3 - x*z^2 (Hess model)
sage: J.curve() == E.affine_patch(2).projective_closure()
True
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> k = GF((Integer(5),Integer(2)), names=('a',)); (a,) = k._first_ngens(1)
>>> E = EllipticCurve(k,[Integer(1),Integer(0)]); E
Elliptic Curve defined by y^2 = x^3 + x over Finite Field in a of size 5^2
>>> E.order()
32
>>> P = E([a, Integer(2)*a + Integer(4)])
>>> P
(a : 2*a + 4 : 1)
>>> P.order()
8
>>> p = P.point_of_jacobian_of_curve()
>>> p
[Place (x + 4*a, y + 3*a + 1)]
>>> p.order()
8
>>> Q = Integer(3)*P
>>> q = Q.point_of_jacobian_of_curve()
>>> q == Integer(3)*p
True
>>> G = p.parent()
>>> G.order()
32
>>> G
Group of rational points of Jacobian over Finite Field in a of size 5^2 (Hess model)
>>> J = G.parent(); J
Jacobian of Projective Plane Curve over Finite Field in a of size 5^2
 defined by x^2*y + y^3 - x*z^2 (Hess model)
>>> J.curve() == E.affine_patch(Integer(2)).projective_closure()
True
set_order(value, multiple, check=None)[source]

Set the cached order of this point (i.e., the value of self._order) to the given value.

Alternatively, when multiple is given, this method will first run order_from_multiple() to determine the exact order from the given multiple of the point order, then cache the result.

Use this when you know a priori the order of this point, or a multiple of the order, to avoid a potentially expensive order calculation.

INPUT:

  • value – positive integer

  • multiple – positive integer; mutually exclusive with value

OUTPUT: none

EXAMPLES:

This example illustrates basic usage.

sage: E = EllipticCurve(GF(7), [0, 1])  # This curve has order 12
sage: G = E(5, 0)
sage: G.set_order(2)
sage: 2*G
(0 : 1 : 0)
sage: G = E(0, 6)
sage: G.set_order(multiple=12)
sage: G._order
3
>>> from sage.all import *
>>> E = EllipticCurve(GF(Integer(7)), [Integer(0), Integer(1)])  # This curve has order 12
>>> G = E(Integer(5), Integer(0))
>>> G.set_order(Integer(2))
>>> Integer(2)*G
(0 : 1 : 0)
>>> G = E(Integer(0), Integer(6))
>>> G.set_order(multiple=Integer(12))
>>> G._order
3

We now give a more interesting case, the NIST-P521 curve. Its order is too big to calculate with Sage, and takes a long time using other packages, so it is very useful here.

sage: # needs sage.rings.finite_rings
sage: p = 2^521 - 1
sage: prev_proof_state = proof.arithmetic()
sage: proof.arithmetic(False)  # turn off primality checking
sage: F = GF(p)
sage: A = p - 3
sage: B = 1093849038073734274511112390766805569936207598951683748994586394495953116150735016013708737573759623248592132296706313309438452531591012912142327488478985984
sage: q = 6864797660130609714981900799081393217269435300143305409394463459185543183397655394245057746333217197532963996371363321113864768612440380340372808892707005449
sage: E = EllipticCurve([F(A), F(B)])
sage: G = E.random_point()
sage: G.set_order(q)
sage: G.order() * G  # This takes practically no time.
(0 : 1 : 0)
sage: proof.arithmetic(prev_proof_state) # restore state
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> p = Integer(2)**Integer(521) - Integer(1)
>>> prev_proof_state = proof.arithmetic()
>>> proof.arithmetic(False)  # turn off primality checking
>>> F = GF(p)
>>> A = p - Integer(3)
>>> B = Integer(1093849038073734274511112390766805569936207598951683748994586394495953116150735016013708737573759623248592132296706313309438452531591012912142327488478985984)
>>> q = Integer(6864797660130609714981900799081393217269435300143305409394463459185543183397655394245057746333217197532963996371363321113864768612440380340372808892707005449)
>>> E = EllipticCurve([F(A), F(B)])
>>> G = E.random_point()
>>> G.set_order(q)
>>> G.order() * G  # This takes practically no time.
(0 : 1 : 0)
>>> proof.arithmetic(prev_proof_state) # restore state

Using .set_order() with a multiple= argument can be used to compute a point’s order significantly faster than calling order() if the point is already known to be \(m\)-torsion:

sage: F.<a> = GF((10007, 23))
sage: E = EllipticCurve(F, [9,9])
sage: n = E.order()
sage: m = 5 * 47 * 139 * 1427 * 2027 * 4831 * 275449 * 29523031
sage: assert m.divides(n)
sage: P = n/m * E.lift_x(6747+a)
sage: assert m * P == 0
sage: P.set_order(multiple=m)   # compute exact order
sage: factor(m // P.order())    # order is now cached
47 * 139
>>> from sage.all import *
>>> F = GF((Integer(10007), Integer(23)), names=('a',)); (a,) = F._first_ngens(1)
>>> E = EllipticCurve(F, [Integer(9),Integer(9)])
>>> n = E.order()
>>> m = Integer(5) * Integer(47) * Integer(139) * Integer(1427) * Integer(2027) * Integer(4831) * Integer(275449) * Integer(29523031)
>>> assert m.divides(n)
>>> P = n/m * E.lift_x(Integer(6747)+a)
>>> assert m * P == Integer(0)
>>> P.set_order(multiple=m)   # compute exact order
>>> factor(m // P.order())    # order is now cached
47 * 139

The algorithm used internally for this functionality is order_from_multiple(). Indeed, simply calling order() on P would take much longer since factoring n is fairly expensive:

sage: n == m * 6670822796985115651 * 441770032618665681677 * 9289973478285634606114927
True
>>> from sage.all import *
>>> n == m * Integer(6670822796985115651) * Integer(441770032618665681677) * Integer(9289973478285634606114927)
True

It is an error to pass a value equal to \(0\):

sage: # needs sage.rings.finite_rings
sage: E = EllipticCurve(GF(7), [0, 1])  # This curve has order 12
sage: G = E.random_point()
sage: G.set_order(0)
Traceback (most recent call last):
...
ValueError: Value 0 illegal for point order
sage: G.set_order(1000)
Traceback (most recent call last):
...
ValueError: Value 1000 illegal: outside max Hasse bound
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> E = EllipticCurve(GF(Integer(7)), [Integer(0), Integer(1)])  # This curve has order 12
>>> G = E.random_point()
>>> G.set_order(Integer(0))
Traceback (most recent call last):
...
ValueError: Value 0 illegal for point order
>>> G.set_order(Integer(1000))
Traceback (most recent call last):
...
ValueError: Value 1000 illegal: outside max Hasse bound

It is also very likely an error to pass a value which is not the actual order of this point. How unlikely is determined by the factorization of the actual order, and the actual group structure:

sage: E = EllipticCurve(GF(7), [0, 1])  # This curve has order 12
sage: G = E(5, 0)   # G has order 2
sage: G.set_order(11)
Traceback (most recent call last):
...
ValueError: Value 11 illegal: 11 * (5 : 0 : 1) is not the identity
>>> from sage.all import *
>>> E = EllipticCurve(GF(Integer(7)), [Integer(0), Integer(1)])  # This curve has order 12
>>> G = E(Integer(5), Integer(0))   # G has order 2
>>> G.set_order(Integer(11))
Traceback (most recent call last):
...
ValueError: Value 11 illegal: 11 * (5 : 0 : 1) is not the identity

However, set_order can be fooled. For instance, the order can be set to a multiple the actual order:

sage: E = EllipticCurve(GF(7), [0, 1])  # This curve has order 12
sage: G = E(5, 0)   # G has order 2
sage: G.set_order(8)
sage: G.order()
8
>>> from sage.all import *
>>> E = EllipticCurve(GF(Integer(7)), [Integer(0), Integer(1)])  # This curve has order 12
>>> G = E(Integer(5), Integer(0))   # G has order 2
>>> G.set_order(Integer(8))
>>> G.order()
8

AUTHORS:

  • Mariah Lenox (2011-02-16)

  • Lorenz Panny (2022): add multiple= option

tate_pairing(Q, n, k, q=None)[source]

Return Tate pairing of \(n\)-torsion point \(P = self\) and point \(Q\).

The value returned is \(f_{n,P}(Q)^e\) where \(f_{n,P}\) is a function with divisor \(n[P]-n[O].\). This is also known as the “modified Tate pairing”. It is a well-defined bilinear map.

INPUT:

  • P=self – elliptic curve point having order \(n\)

  • Q – elliptic curve point on same curve as \(P\) (can be any order)

  • n – positive integer; order of \(P\)

  • k – positive integer; embedding degree

  • q – positive integer; size of base field (the “big” field is \(GF(q^k)\). \(q\) needs to be set only if its value cannot be deduced.)

OUTPUT:

An \(n\)-th root of unity in the base field self.curve().base_field().

EXAMPLES:

A simple example, pairing a point with itself, and pairing a point with another rational point:

sage: p = 103; A = 1; B = 18; E = EllipticCurve(GF(p), [A, B])
sage: P = E(33, 91); n = P.order(); n
19
sage: k = GF(n)(p).multiplicative_order(); k
6
sage: P.tate_pairing(P, n, k)
1
sage: Q = E(87, 51)
sage: P.tate_pairing(Q, n, k)
1
sage: set_random_seed(35)
sage: P.tate_pairing(P, n, k)
1
>>> from sage.all import *
>>> p = Integer(103); A = Integer(1); B = Integer(18); E = EllipticCurve(GF(p), [A, B])
>>> P = E(Integer(33), Integer(91)); n = P.order(); n
19
>>> k = GF(n)(p).multiplicative_order(); k
6
>>> P.tate_pairing(P, n, k)
1
>>> Q = E(Integer(87), Integer(51))
>>> P.tate_pairing(Q, n, k)
1
>>> set_random_seed(Integer(35))
>>> P.tate_pairing(P, n, k)
1

We now let \(Q\) be a point on the same curve as above, but defined over the pairing extension field, and we also demonstrate the bilinearity of the pairing:

sage: # needs sage.rings.finite_rings
sage: K.<a> = GF((p,k))
sage: EK = E.base_extend(K); P = EK(P)
sage: Qx = 69*a^5 + 96*a^4 + 22*a^3 + 86*a^2 + 6*a + 35
sage: Qy = 34*a^5 + 24*a^4 + 16*a^3 + 41*a^2 + 4*a + 40
sage: Q = EK(Qx, Qy);
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> K = GF((p,k), names=('a',)); (a,) = K._first_ngens(1)
>>> EK = E.base_extend(K); P = EK(P)
>>> Qx = Integer(69)*a**Integer(5) + Integer(96)*a**Integer(4) + Integer(22)*a**Integer(3) + Integer(86)*a**Integer(2) + Integer(6)*a + Integer(35)
>>> Qy = Integer(34)*a**Integer(5) + Integer(24)*a**Integer(4) + Integer(16)*a**Integer(3) + Integer(41)*a**Integer(2) + Integer(4)*a + Integer(40)
>>> Q = EK(Qx, Qy);

Multiply by cofactor so Q has order n:

sage: # needs sage.rings.finite_rings
sage: h = 551269674; Q = h*Q
sage: P = EK(P); P.tate_pairing(Q, n, k)
24*a^5 + 34*a^4 + 3*a^3 + 69*a^2 + 86*a + 45
sage: s = Integer(randrange(1,n))
sage: ans1 = (s*P).tate_pairing(Q, n, k)
sage: ans2 = P.tate_pairing(s*Q, n, k)
sage: ans3 = P.tate_pairing(Q, n, k)^s
sage: ans1 == ans2 == ans3
True
sage: (ans1 != 1) and (ans1^n == 1)
True
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> h = Integer(551269674); Q = h*Q
>>> P = EK(P); P.tate_pairing(Q, n, k)
24*a^5 + 34*a^4 + 3*a^3 + 69*a^2 + 86*a + 45
>>> s = Integer(randrange(Integer(1),n))
>>> ans1 = (s*P).tate_pairing(Q, n, k)
>>> ans2 = P.tate_pairing(s*Q, n, k)
>>> ans3 = P.tate_pairing(Q, n, k)**s
>>> ans1 == ans2 == ans3
True
>>> (ans1 != Integer(1)) and (ans1**n == Integer(1))
True

Here is an example of using the Tate pairing to compute the Weil pairing (using the same data as above):

sage: # needs sage.rings.finite_rings
sage: e = Integer((p^k-1)/n); e
62844857712
sage: P.weil_pairing(Q, n)^e
94*a^5 + 99*a^4 + 29*a^3 + 45*a^2 + 57*a + 34
sage: P.tate_pairing(Q, n, k) == P._miller_(Q, n)^e
True
sage: Q.tate_pairing(P, n, k) == Q._miller_(P, n)^e
True
sage: P.tate_pairing(Q, n, k)/Q.tate_pairing(P, n, k)
94*a^5 + 99*a^4 + 29*a^3 + 45*a^2 + 57*a + 34
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> e = Integer((p**k-Integer(1))/n); e
62844857712
>>> P.weil_pairing(Q, n)**e
94*a^5 + 99*a^4 + 29*a^3 + 45*a^2 + 57*a + 34
>>> P.tate_pairing(Q, n, k) == P._miller_(Q, n)**e
True
>>> Q.tate_pairing(P, n, k) == Q._miller_(P, n)**e
True
>>> P.tate_pairing(Q, n, k)/Q.tate_pairing(P, n, k)
94*a^5 + 99*a^4 + 29*a^3 + 45*a^2 + 57*a + 34

An example where we have to pass the base field size (and we again have agreement with the Weil pairing):

sage: # needs sage.rings.finite_rings
sage: F.<a> = GF((2,5))
sage: E = EllipticCurve(F, [0,0,1,1,1])
sage: P = E(a^4 + 1, a^3)
sage: Fx.<b> = GF((2,4*5))
sage: Ex = EllipticCurve(Fx,[0,0,1,1,1])
sage: phi = Hom(F, Fx)(F.gen().minpoly().roots(Fx)[0][0])
sage: Px = Ex(phi(P.x()), phi(P.y()))
sage: 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)
sage: Px.tate_pairing(Qx, n=41, k=4)
Traceback (most recent call last):
...
ValueError: Unexpected field degree: set keyword argument q equal to
the size of the base field (big field is GF(q^4)).
sage: num = Px.tate_pairing(Qx, n=41, k=4, q=32); num
b^19 + b^14 + b^13 + b^12 + b^6 + b^4 + b^3
sage: den = Qx.tate_pairing(Px, n=41, k=4, q=32); den
b^19 + b^17 + b^16 + b^15 + b^14 + b^10 + b^6 + b^2 + 1
sage: e = Integer((32^4-1)/41); e
25575
sage: Px.weil_pairing(Qx, 41)^e == num/den
True
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> F = GF((Integer(2),Integer(5)), names=('a',)); (a,) = F._first_ngens(1)
>>> E = EllipticCurve(F, [Integer(0),Integer(0),Integer(1),Integer(1),Integer(1)])
>>> P = E(a**Integer(4) + Integer(1), a**Integer(3))
>>> Fx = GF((Integer(2),Integer(4)*Integer(5)), names=('b',)); (b,) = Fx._first_ngens(1)
>>> Ex = EllipticCurve(Fx,[Integer(0),Integer(0),Integer(1),Integer(1),Integer(1)])
>>> phi = Hom(F, Fx)(F.gen().minpoly().roots(Fx)[Integer(0)][Integer(0)])
>>> Px = Ex(phi(P.x()), phi(P.y()))
>>> Qx = Ex(b**Integer(19) + b**Integer(18) + b**Integer(16) + b**Integer(12) + b**Integer(10) + b**Integer(9) + b**Integer(8) + b**Integer(5) + b**Integer(3) + Integer(1),
...         b**Integer(18) + b**Integer(13) + b**Integer(10) + b**Integer(8) + b**Integer(5) + b**Integer(4) + b**Integer(3) + b)
>>> Px.tate_pairing(Qx, n=Integer(41), k=Integer(4))
Traceback (most recent call last):
...
ValueError: Unexpected field degree: set keyword argument q equal to
the size of the base field (big field is GF(q^4)).
>>> num = Px.tate_pairing(Qx, n=Integer(41), k=Integer(4), q=Integer(32)); num
b^19 + b^14 + b^13 + b^12 + b^6 + b^4 + b^3
>>> den = Qx.tate_pairing(Px, n=Integer(41), k=Integer(4), q=Integer(32)); den
b^19 + b^17 + b^16 + b^15 + b^14 + b^10 + b^6 + b^2 + 1
>>> e = Integer((Integer(32)**Integer(4)-Integer(1))/Integer(41)); e
25575
>>> Px.weil_pairing(Qx, Integer(41))**e == num/den
True

An example over a large base field:

sage: F = GF(65537^2, modulus=[3,46810,1], name='z2')
sage: F.inject_variables()
Defining z2
sage: E = EllipticCurve(F, [0,1])
sage: P = E(22, 28891)
sage: Q = E(-93, 40438*z2 + 31573)
sage: P.tate_pairing(Q, 7282, 2)
34585*z2 + 4063
>>> from sage.all import *
>>> F = GF(Integer(65537)**Integer(2), modulus=[Integer(3),Integer(46810),Integer(1)], name='z2')
>>> F.inject_variables()
Defining z2
>>> E = EllipticCurve(F, [Integer(0),Integer(1)])
>>> P = E(Integer(22), Integer(28891))
>>> Q = E(-Integer(93), Integer(40438)*z2 + Integer(31573))
>>> P.tate_pairing(Q, Integer(7282), Integer(2))
34585*z2 + 4063

ALGORITHM:

  • pari:elltatepairing computes the non-reduced tate pairing and the exponentiation is handled by Sage using user input for \(k\) (and optionally \(q\)).

AUTHORS:

  • Mariah Lenox (2011-03-07)

  • Giacomo Pope (2024): Use of PARI for the non-reduced Tate pairing

weil_pairing(Q, n, algorithm=None)[source]

Compute the Weil pairing of this point with another point \(Q\) on the same curve.

INPUT:

  • Q – another point on the same curve as self

  • n – integer \(n\) such that \(nP = nQ = (0:1:0)\), where \(P\) is self

  • algorithm – (default: None) choices are 'pari' and 'sage'. PARI is usually significantly faster, but it only works over finite fields. When None is given, a suitable algorithm is chosen automatically.

OUTPUT: an \(n\)-th root of unity in the base field of the curve

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: F.<a> = GF((2,5))
sage: E = EllipticCurve(F, [0,0,1,1,1])
sage: P = E(a^4 + 1, a^3)
sage: Fx.<b> = GF((2, 4*5))
sage: Ex = EllipticCurve(Fx, [0,0,1,1,1])
sage: phi = Hom(F, Fx)(F.gen().minpoly().roots(Fx)[0][0])
sage: Px = Ex(phi(P.x()), phi(P.y()))
sage: O = Ex(0)
sage: 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)
sage: Px.weil_pairing(Qx, 41) == b^19 + b^15 + b^9 + b^8 + b^6 + b^4 + b^3 + b^2 + 1
True
sage: Px.weil_pairing(17*Px, 41) == Fx(1)
True
sage: Px.weil_pairing(O, 41) == Fx(1)
True
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> F = GF((Integer(2),Integer(5)), names=('a',)); (a,) = F._first_ngens(1)
>>> E = EllipticCurve(F, [Integer(0),Integer(0),Integer(1),Integer(1),Integer(1)])
>>> P = E(a**Integer(4) + Integer(1), a**Integer(3))
>>> Fx = GF((Integer(2), Integer(4)*Integer(5)), names=('b',)); (b,) = Fx._first_ngens(1)
>>> Ex = EllipticCurve(Fx, [Integer(0),Integer(0),Integer(1),Integer(1),Integer(1)])
>>> phi = Hom(F, Fx)(F.gen().minpoly().roots(Fx)[Integer(0)][Integer(0)])
>>> Px = Ex(phi(P.x()), phi(P.y()))
>>> O = Ex(Integer(0))
>>> Qx = Ex(b**Integer(19) + b**Integer(18) + b**Integer(16) + b**Integer(12) + b**Integer(10) + b**Integer(9) + b**Integer(8) + b**Integer(5) + b**Integer(3) + Integer(1),
...         b**Integer(18) + b**Integer(13) + b**Integer(10) + b**Integer(8) + b**Integer(5) + b**Integer(4) + b**Integer(3) + b)
>>> Px.weil_pairing(Qx, Integer(41)) == b**Integer(19) + b**Integer(15) + b**Integer(9) + b**Integer(8) + b**Integer(6) + b**Integer(4) + b**Integer(3) + b**Integer(2) + Integer(1)
True
>>> Px.weil_pairing(Integer(17)*Px, Integer(41)) == Fx(Integer(1))
True
>>> Px.weil_pairing(O, Integer(41)) == Fx(Integer(1))
True

An error is raised if either point is not \(n\)-torsion:

sage: Px.weil_pairing(O, 40)                                                # needs sage.rings.finite_rings
Traceback (most recent call last):
...
ValueError: points must both be n-torsion
>>> from sage.all import *
>>> Px.weil_pairing(O, Integer(40))                                                # needs sage.rings.finite_rings
Traceback (most recent call last):
...
ValueError: points must both be n-torsion

A larger example (see Issue #4964):

sage: # needs sage.rings.finite_rings
sage: P, Q = EllipticCurve(GF((19,4),'a'), [-1,0]).gens()
sage: P.order(), Q.order()
(360, 360)
sage: z = P.weil_pairing(Q, 360)
sage: z.multiplicative_order()
360
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> P, Q = EllipticCurve(GF((Integer(19),Integer(4)),'a'), [-Integer(1),Integer(0)]).gens()
>>> P.order(), Q.order()
(360, 360)
>>> z = P.weil_pairing(Q, Integer(360))
>>> z.multiplicative_order()
360

Another larger example:

sage: F = GF(65537^2, modulus=[3,-1,1], name='a')
sage: F.inject_variables()
Defining a
sage: E = EllipticCurve(F, [0,1])
sage: P = E(22, 28891)
sage: Q = E(-93, 2728*a + 64173)
sage: P.weil_pairing(Q, 7282, algorithm='sage')
53278*a + 36700
>>> from sage.all import *
>>> F = GF(Integer(65537)**Integer(2), modulus=[Integer(3),-Integer(1),Integer(1)], name='a')
>>> F.inject_variables()
Defining a
>>> E = EllipticCurve(F, [Integer(0),Integer(1)])
>>> P = E(Integer(22), Integer(28891))
>>> Q = E(-Integer(93), Integer(2728)*a + Integer(64173))
>>> P.weil_pairing(Q, Integer(7282), algorithm='sage')
53278*a + 36700

An example over a number field:

sage: # needs sage.rings.number_field
sage: E = EllipticCurve('11a1').change_ring(CyclotomicField(5))
sage: P, Q = E.torsion_subgroup().gens()
sage: P, Q = (P.element(), Q.element())
sage: (P.order(), Q.order())
(5, 5)
sage: P.weil_pairing(Q, 5)
zeta5^2
sage: Q.weil_pairing(P, 5)
zeta5^3
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> E = EllipticCurve('11a1').change_ring(CyclotomicField(Integer(5)))
>>> P, Q = E.torsion_subgroup().gens()
>>> P, Q = (P.element(), Q.element())
>>> (P.order(), Q.order())
(5, 5)
>>> P.weil_pairing(Q, Integer(5))
zeta5^2
>>> Q.weil_pairing(P, Integer(5))
zeta5^3

ALGORITHM:

  • For algorithm='pari': pari:ellweilpairing.

  • For algorithm='sage': Implemented using Proposition 8 in [Mil2004]. The value 1 is returned for linearly dependent input points. This condition is caught via a ZeroDivisionError, since the use of a discrete logarithm test for linear dependence is much too slow for large \(n\).

AUTHORS:

  • David Hansen (2009-01-25)

  • Lorenz Panny (2022): algorithm='pari'

x()[source]

Return the \(x\) coordinate of this point, as an element of the base field. If this is the point at infinity, a ZeroDivisionError is raised.

EXAMPLES:

sage: E = EllipticCurve('389a')
sage: P = E([-1,1])
sage: P.x()
-1
sage: Q = E(0); Q
(0 : 1 : 0)
sage: Q.x()
Traceback (most recent call last):
...
ZeroDivisionError: rational division by zero
>>> from sage.all import *
>>> E = EllipticCurve('389a')
>>> P = E([-Integer(1),Integer(1)])
>>> P.x()
-1
>>> Q = E(Integer(0)); Q
(0 : 1 : 0)
>>> Q.x()
Traceback (most recent call last):
...
ZeroDivisionError: rational division by zero
xy()[source]

Return the \(x\) and \(y\) coordinates of this point, as a 2-tuple. If this is the point at infinity, a ZeroDivisionError is raised.

EXAMPLES:

sage: E = EllipticCurve('389a')
sage: P = E([-1,1])
sage: P.xy()
(-1, 1)
sage: Q = E(0); Q
(0 : 1 : 0)
sage: Q.xy()
Traceback (most recent call last):
...
ZeroDivisionError: rational division by zero
>>> from sage.all import *
>>> E = EllipticCurve('389a')
>>> P = E([-Integer(1),Integer(1)])
>>> P.xy()
(-1, 1)
>>> Q = E(Integer(0)); Q
(0 : 1 : 0)
>>> Q.xy()
Traceback (most recent call last):
...
ZeroDivisionError: rational division by zero
y()[source]

Return the \(y\) coordinate of this point, as an element of the base field. If this is the point at infinity, a ZeroDivisionError is raised.

EXAMPLES:

sage: E = EllipticCurve('389a')
sage: P = E([-1,1])
sage: P.y()
1
sage: Q = E(0); Q
(0 : 1 : 0)
sage: Q.y()
Traceback (most recent call last):
...
ZeroDivisionError: rational division by zero
>>> from sage.all import *
>>> E = EllipticCurve('389a')
>>> P = E([-Integer(1),Integer(1)])
>>> P.y()
1
>>> Q = E(Integer(0)); Q
(0 : 1 : 0)
>>> Q.y()
Traceback (most recent call last):
...
ZeroDivisionError: rational division by zero
class sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_finite_field(curve, v, check=True)[source]

Bases: EllipticCurvePoint_field

Class for elliptic curve points over finite fields.

additive_order()[source]

Return the order of this point on the elliptic curve.

ALGORITHM: Use PARI function pari:ellorder.

Note

additive_order() is a synonym for order()

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: k.<a> = GF((5,5))
sage: E = EllipticCurve(k,[2,4]); E
Elliptic Curve defined by y^2 = x^3 + 2*x + 4 over Finite Field in a of size 5^5
sage: P = E(3*a^4 + 3*a, 2*a + 1)
sage: P.order()
3227
sage: Q = E(0,2)
sage: Q.order()
7
sage: Q.additive_order()
7
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> k = GF((Integer(5),Integer(5)), names=('a',)); (a,) = k._first_ngens(1)
>>> E = EllipticCurve(k,[Integer(2),Integer(4)]); E
Elliptic Curve defined by y^2 = x^3 + 2*x + 4 over Finite Field in a of size 5^5
>>> P = E(Integer(3)*a**Integer(4) + Integer(3)*a, Integer(2)*a + Integer(1))
>>> P.order()
3227
>>> Q = E(Integer(0),Integer(2))
>>> Q.order()
7
>>> Q.additive_order()
7

sage: # needs sage.rings.finite_rings
sage: p = next_prime(2^150)
sage: E = EllipticCurve(GF(p), [1,1])
sage: P = E(831623307675610677632782670796608848711856078,
....:       42295786042873366706573292533588638217232964)
sage: P.order()
1427247692705959881058262545272474300628281448
sage: P.order() == E.cardinality()
True
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> p = next_prime(Integer(2)**Integer(150))
>>> E = EllipticCurve(GF(p), [Integer(1),Integer(1)])
>>> P = E(Integer(831623307675610677632782670796608848711856078),
...       Integer(42295786042873366706573292533588638217232964))
>>> P.order()
1427247692705959881058262545272474300628281448
>>> P.order() == E.cardinality()
True

The next example has \(j(E)=0\):

sage: # needs sage.rings.finite_rings
sage: p = 33554501
sage: F.<u> = GF((p,2))
sage: E = EllipticCurve(F, [0,1])
sage: E.j_invariant()
0
sage: P = E.random_point()
sage: P.order()  # random
16777251
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> p = Integer(33554501)
>>> F = GF((p,Integer(2)), names=('u',)); (u,) = F._first_ngens(1)
>>> E = EllipticCurve(F, [Integer(0),Integer(1)])
>>> E.j_invariant()
0
>>> P = E.random_point()
>>> P.order()  # random
16777251

Similarly when \(j(E)=1728\):

sage: # needs sage.rings.finite_rings
sage: p = 33554473
sage: F.<u> = GF((p,2))
sage: E = EllipticCurve(F, [1,0])
sage: E.j_invariant()
1728
sage: P = E.random_point()
sage: P.order()  # random
46912611635760
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> p = Integer(33554473)
>>> F = GF((p,Integer(2)), names=('u',)); (u,) = F._first_ngens(1)
>>> E = EllipticCurve(F, [Integer(1),Integer(0)])
>>> E.j_invariant()
1728
>>> P = E.random_point()
>>> P.order()  # random
46912611635760
discrete_log(Q)[source]

Legacy version of log() with its arguments swapped.

Note that this method uses the opposite argument ordering of all other logarithm methods in Sage; see Issue #37150.

EXAMPLES:

sage: E = EllipticCurve(j=GF(101)(5))
sage: P, = E.gens()
sage: (2*P).log(P)
2
sage: (2*P).discrete_log(P)
doctest:warning ...
DeprecationWarning: The syntax P.discrete_log(Q) ... Please update your code. ...
45
sage: P.discrete_log(2*P)
2
>>> from sage.all import *
>>> E = EllipticCurve(j=GF(Integer(101))(Integer(5)))
>>> P, = E.gens()
>>> (Integer(2)*P).log(P)
2
>>> (Integer(2)*P).discrete_log(P)
doctest:warning ...
DeprecationWarning: The syntax P.discrete_log(Q) ... Please update your code. ...
45
>>> P.discrete_log(Integer(2)*P)
2
has_finite_order()[source]

Return True if this point has finite additive order as an element of the group of points on this curve.

Since the base field is finite, the answer will always be True.

EXAMPLES:

sage: E = EllipticCurve(GF(7), [1,3])
sage: P = E.points()[3]
sage: P.has_finite_order()
True
>>> from sage.all import *
>>> E = EllipticCurve(GF(Integer(7)), [Integer(1),Integer(3)])
>>> P = E.points()[Integer(3)]
>>> P.has_finite_order()
True
log(base)[source]

Return the discrete logarithm of this point to the given base. In other words, return an integer \(x\) such that \(xP = Q\) where \(P\) is base and \(Q\) is this point.

If base is a list or tuple of two points, then this function solves a two-dimensional discrete logarithm: Given \((P_1,P_2)\) in base, it returns a tuple of integers \((x,y)\) such that \([x]P_1 + [y]P_2 = Q\), where \(Q\) is this point.

A ValueError is raised if there is no solution.

ALGORITHM:

To compute the actual logarithm, pari:elllog is called.

However, elllog() does not guarantee termination if \(Q\) is not a multiple of \(P\), so we first need to check subgroup membership. This is done as follows:

  • Let \(n\) denote the order of \(P\). First check that \(nQ\) equals the point at infinity (and hence the order of \(Q\) divides \(n\)).

  • If the curve order \(\#E\) has been cached, check whether \(\gcd(n^2, \#E) = n\). If this holds, the curve has cyclic \(n\)-torsion, hence all points whose order divides \(n\) must be multiples of \(P\) and we are done.

  • Otherwise (if this test is inconclusive), check that the Weil pairing of \(P\) and \(Q\) is trivial.

For anomalous curves with \(\#E = p\), the padic_elliptic_logarithm() function is called.

For two-dimensional logarithms, we first compute the Weil pairings of \(Q\) with \(P_1\) and \(P_2\) and their logarithms relative to the pairing of \(P_1\) and \(P_2\); this allows reading off \(x\) and \(y\) modulo the part of the order where \(P_1\) and \(P_2\) are independent. Modulo the remaining part of the order the logarithm ends up being effectively one-dimensional, so we can reduce the problem to the basic one-dimensional case and finally recombine the results.

INPUT:

  • base – another point or sequence of two points on the same curve as self

OUTPUT:

(integer) – The discrete logarithm of \(Q\) with respect to \(P\), which is an integer \(x\) with \(0\le x<\mathrm{ord}(P)\) such that \(xP=Q\), if one exists. In the case of two points \(P_1,P_2\), two integers \(x,y\) with \(0\le x<\mathrm{ord}(P_1)\) and \(0\le y<\mathrm{ord}(P_2)\) such that \([x]P_1 + [y]P_2 = Q\).

AUTHORS:

  • John Cremona. Adapted to use generic functions 2008-04-05.

  • Lorenz Panny (2022): switch to PARI.

  • Lorenz Panny (2024): the two-dimensional case.

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: F = GF((3,6),'a')
sage: a = F.gen()
sage: E = EllipticCurve([0,1,1,a,a])
sage: E.cardinality()
762
sage: P = E.gens()[0]
sage: Q = 400*P
sage: Q.log(P)
400
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> F = GF((Integer(3),Integer(6)),'a')
>>> a = F.gen()
>>> E = EllipticCurve([Integer(0),Integer(1),Integer(1),a,a])
>>> E.cardinality()
762
>>> P = E.gens()[Integer(0)]
>>> Q = Integer(400)*P
>>> Q.log(P)
400

sage: # needs sage.rings.finite_rings
sage: F = GF((5, 60), 'a')
sage: E = EllipticCurve(F, [1, 1])
sage: E.abelian_group()
Additive abelian group isomorphic to Z/194301464603136995341424045476456938000 + Z/4464 embedded in Abelian group of points on Elliptic Curve defined by y^2 = x^3 + x + 1 over Finite Field in a of size 5^60
sage: P, Q = E.gens()  # cached generators from .abelian_group()
sage: T = 1234567890987654321*P + 1337*Q
sage: T.log([P, Q])
(1234567890987654321, 1337)
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> F = GF((Integer(5), Integer(60)), 'a')
>>> E = EllipticCurve(F, [Integer(1), Integer(1)])
>>> E.abelian_group()
Additive abelian group isomorphic to Z/194301464603136995341424045476456938000 + Z/4464 embedded in Abelian group of points on Elliptic Curve defined by y^2 = x^3 + x + 1 over Finite Field in a of size 5^60
>>> P, Q = E.gens()  # cached generators from .abelian_group()
>>> T = Integer(1234567890987654321)*P + Integer(1337)*Q
>>> T.log([P, Q])
(1234567890987654321, 1337)
order()[source]

Return the order of this point on the elliptic curve.

ALGORITHM: Use PARI function pari:ellorder.

Note

additive_order() is a synonym for order()

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: k.<a> = GF((5,5))
sage: E = EllipticCurve(k,[2,4]); E
Elliptic Curve defined by y^2 = x^3 + 2*x + 4 over Finite Field in a of size 5^5
sage: P = E(3*a^4 + 3*a, 2*a + 1)
sage: P.order()
3227
sage: Q = E(0,2)
sage: Q.order()
7
sage: Q.additive_order()
7
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> k = GF((Integer(5),Integer(5)), names=('a',)); (a,) = k._first_ngens(1)
>>> E = EllipticCurve(k,[Integer(2),Integer(4)]); E
Elliptic Curve defined by y^2 = x^3 + 2*x + 4 over Finite Field in a of size 5^5
>>> P = E(Integer(3)*a**Integer(4) + Integer(3)*a, Integer(2)*a + Integer(1))
>>> P.order()
3227
>>> Q = E(Integer(0),Integer(2))
>>> Q.order()
7
>>> Q.additive_order()
7

sage: # needs sage.rings.finite_rings
sage: p = next_prime(2^150)
sage: E = EllipticCurve(GF(p), [1,1])
sage: P = E(831623307675610677632782670796608848711856078,
....:       42295786042873366706573292533588638217232964)
sage: P.order()
1427247692705959881058262545272474300628281448
sage: P.order() == E.cardinality()
True
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> p = next_prime(Integer(2)**Integer(150))
>>> E = EllipticCurve(GF(p), [Integer(1),Integer(1)])
>>> P = E(Integer(831623307675610677632782670796608848711856078),
...       Integer(42295786042873366706573292533588638217232964))
>>> P.order()
1427247692705959881058262545272474300628281448
>>> P.order() == E.cardinality()
True

The next example has \(j(E)=0\):

sage: # needs sage.rings.finite_rings
sage: p = 33554501
sage: F.<u> = GF((p,2))
sage: E = EllipticCurve(F, [0,1])
sage: E.j_invariant()
0
sage: P = E.random_point()
sage: P.order()  # random
16777251
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> p = Integer(33554501)
>>> F = GF((p,Integer(2)), names=('u',)); (u,) = F._first_ngens(1)
>>> E = EllipticCurve(F, [Integer(0),Integer(1)])
>>> E.j_invariant()
0
>>> P = E.random_point()
>>> P.order()  # random
16777251

Similarly when \(j(E)=1728\):

sage: # needs sage.rings.finite_rings
sage: p = 33554473
sage: F.<u> = GF((p,2))
sage: E = EllipticCurve(F, [1,0])
sage: E.j_invariant()
1728
sage: P = E.random_point()
sage: P.order()  # random
46912611635760
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> p = Integer(33554473)
>>> F = GF((p,Integer(2)), names=('u',)); (u,) = F._first_ngens(1)
>>> E = EllipticCurve(F, [Integer(1),Integer(0)])
>>> E.j_invariant()
1728
>>> P = E.random_point()
>>> P.order()  # random
46912611635760
padic_elliptic_logarithm(Q, p)[source]

Return the discrete logarithm of \(Q\) to base \(P\) = self, that is, an integer \(x\) such that \(xP = Q\) only for anomalous curves.

ALGORITHM:

Discrete logarithm computed as in [Sma1999] with a loop to avoid the canonical lift.

INPUT:

  • Q – point; another point on the same curve as self

  • p – integer; a prime equal to the order of the curve

OUTPUT:

(integer) – The discrete logarithm of \(Q\) with respect to \(P\), which is an integer \(x\) with \(0\le x<\mathrm{ord}(P)\) such that \(xP=Q\).

AUTHORS:

  • Sylvain Pelissier (2022) based on Samuel Neves code.

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: p = 235322474717419
sage: b = 8856682
sage: E = EllipticCurve(GF(p), [0, b])
sage: P = E(200673830421813, 57025307876612)
sage: Q = E(40345734829479, 211738132651297)
sage: x = P.padic_elliptic_logarithm(Q, p)                                  # needs sage.rings.padics
sage: x * P == Q                                                            # needs sage.rings.padics
True
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> p = Integer(235322474717419)
>>> b = Integer(8856682)
>>> E = EllipticCurve(GF(p), [Integer(0), b])
>>> P = E(Integer(200673830421813), Integer(57025307876612))
>>> Q = E(Integer(40345734829479), Integer(211738132651297))
>>> x = P.padic_elliptic_logarithm(Q, p)                                  # needs sage.rings.padics
>>> x * P == Q                                                            # needs sage.rings.padics
True
class sage.schemes.elliptic_curves.ell_point.EllipticCurvePoint_number_field(curve, v, check=True)[source]

Bases: EllipticCurvePoint_field

A point on an elliptic curve over a number field.

Most of the functionality is derived from the parent class EllipticCurvePoint_field. In addition we have support for orders, heights, reduction modulo primes, and elliptic logarithms.

EXAMPLES:

sage: E = EllipticCurve('37a')
sage: E([0,0])
(0 : 0 : 1)
sage: E(0,0)               # brackets are optional
(0 : 0 : 1)
sage: E([GF(5)(0), 0])     # entries are coerced
(0 : 0 : 1)

sage: E(0.000, 0)
(0 : 0 : 1)

sage: E(1,0,0)
Traceback (most recent call last):
...
TypeError: Coordinates [1, 0, 0] do not define a point on
Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
>>> from sage.all import *
>>> E = EllipticCurve('37a')
>>> E([Integer(0),Integer(0)])
(0 : 0 : 1)
>>> E(Integer(0),Integer(0))               # brackets are optional
(0 : 0 : 1)
>>> E([GF(Integer(5))(Integer(0)), Integer(0)])     # entries are coerced
(0 : 0 : 1)

>>> E(RealNumber('0.000'), Integer(0))
(0 : 0 : 1)

>>> E(Integer(1),Integer(0),Integer(0))
Traceback (most recent call last):
...
TypeError: Coordinates [1, 0, 0] do not define a point on
Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field

sage: E = EllipticCurve([0,0,1,-1,0])
sage: S = E(QQ); S
Abelian group of points on
Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0),Integer(0),Integer(1),-Integer(1),Integer(0)])
>>> S = E(QQ); S
Abelian group of points on
Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
additive_order()[source]

Return the order of this point on the elliptic curve.

If the point has infinite order, returns +Infinity. For curves defined over \(\QQ\), we call PARI; over other number fields we implement the function here.

Note

additive_order() is a synonym for order()

EXAMPLES:

sage: E = EllipticCurve([0,0,1,-1,0])
sage: P = E([0,0]); P
(0 : 0 : 1)
sage: P.order()
+Infinity
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0),Integer(0),Integer(1),-Integer(1),Integer(0)])
>>> P = E([Integer(0),Integer(0)]); P
(0 : 0 : 1)
>>> P.order()
+Infinity

sage: E = EllipticCurve([0,1])
sage: P = E([-1,0])
sage: P.order()
2
sage: P.additive_order()
2
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0),Integer(1)])
>>> P = E([-Integer(1),Integer(0)])
>>> P.order()
2
>>> P.additive_order()
2
archimedean_local_height(v=None, prec=None, weighted=False)[source]

Compute the local height of self at the archimedean place \(v\).

INPUT:

  • self – a point on an elliptic curve over a number field \(K\)

  • v – a real or complex embedding of K, or None (default). If \(v\) is a real or complex embedding, return the local height of self at \(v\). If \(v\) is None, return the total archimedean contribution to the global height.

  • prec – integer or None (default). The precision of the computation. If None, the precision is deduced from \(v\)

  • weighted – boolean. If False (default), the height is normalised to be invariant under extension of \(K\). If True, return this normalised height multiplied by the local degree if \(v\) is a single place, or by the degree of \(K\) if \(v\) is None.

OUTPUT:

A real number. The normalisation is twice that in Silverman’s paper [Sil1988]. Note that this local height depends on the model of the curve.

ALGORITHM:

See [Sil1988], Section 4.

EXAMPLES:

Examples 1, 2, and 3 from [Sil1988]:

sage: # needs sage.rings.number_field
sage: K.<a> = QuadraticField(-2)
sage: E = EllipticCurve(K, [0,-1,1,0,0]); E
Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 over Number Field
 in a with defining polynomial x^2 + 2 with a = 1.414213562373095?*I
sage: P = E.lift_x(2 + a); P
(a + 2 : -2*a - 2 : 1)
sage: P.archimedean_local_height(K.places(prec=170)[0]) / 2
0.45754773287523276736211210741423654346576029814695

sage: # needs sage.rings.number_field
sage: x = polygen(ZZ, 'x')
sage: K.<i> = NumberField(x^2 + 1)
sage: E = EllipticCurve(K, [0,0,4,6*i,0]); E
Elliptic Curve defined by y^2 + 4*y = x^3 + 6*i*x
 over Number Field in i with defining polynomial x^2 + 1
sage: P = E((0,0))
sage: P.archimedean_local_height(K.places()[0]) / 2
0.510184995162373

sage: Q = E.lift_x(-9/4); Q                                                 # needs sage.rings.number_field
(-9/4 : 27/8*i - 4 : 1)
sage: Q.archimedean_local_height(K.places()[0]) / 2                         # needs sage.rings.number_field
0.654445619529600
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> K = QuadraticField(-Integer(2), names=('a',)); (a,) = K._first_ngens(1)
>>> E = EllipticCurve(K, [Integer(0),-Integer(1),Integer(1),Integer(0),Integer(0)]); E
Elliptic Curve defined by y^2 + y = x^3 + (-1)*x^2 over Number Field
 in a with defining polynomial x^2 + 2 with a = 1.414213562373095?*I
>>> P = E.lift_x(Integer(2) + a); P
(a + 2 : -2*a - 2 : 1)
>>> P.archimedean_local_height(K.places(prec=Integer(170))[Integer(0)]) / Integer(2)
0.45754773287523276736211210741423654346576029814695

>>> # needs sage.rings.number_field
>>> x = polygen(ZZ, 'x')
>>> K = NumberField(x**Integer(2) + Integer(1), names=('i',)); (i,) = K._first_ngens(1)
>>> E = EllipticCurve(K, [Integer(0),Integer(0),Integer(4),Integer(6)*i,Integer(0)]); E
Elliptic Curve defined by y^2 + 4*y = x^3 + 6*i*x
 over Number Field in i with defining polynomial x^2 + 1
>>> P = E((Integer(0),Integer(0)))
>>> P.archimedean_local_height(K.places()[Integer(0)]) / Integer(2)
0.510184995162373

>>> Q = E.lift_x(-Integer(9)/Integer(4)); Q                                                 # needs sage.rings.number_field
(-9/4 : 27/8*i - 4 : 1)
>>> Q.archimedean_local_height(K.places()[Integer(0)]) / Integer(2)                         # needs sage.rings.number_field
0.654445619529600

An example over the rational numbers:

sage: E = EllipticCurve([0, 0, 0, -36, 0])
sage: P = E([-3, 9])
sage: P.archimedean_local_height()
1.98723816350773
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0), Integer(0), Integer(0), -Integer(36), Integer(0)])
>>> P = E([-Integer(3), Integer(9)])
>>> P.archimedean_local_height()
1.98723816350773

Local heights of torsion points can be nonzero (unlike the global height):

sage: # needs sage.rings.number_field
sage: K.<i> = QuadraticField(-1)
sage: E = EllipticCurve([0, 0, 0, K(1), 0])
sage: P = E(i, 0)
sage: P.archimedean_local_height()
0.346573590279973
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> K = QuadraticField(-Integer(1), names=('i',)); (i,) = K._first_ngens(1)
>>> E = EllipticCurve([Integer(0), Integer(0), Integer(0), K(Integer(1)), Integer(0)])
>>> P = E(i, Integer(0))
>>> P.archimedean_local_height()
0.346573590279973
elliptic_logarithm(embedding=None, precision=100, algorithm='pari')[source]

Return the elliptic logarithm of this elliptic curve point.

An embedding of the base field into \(\RR\) or \(\CC\) (with arbitrary precision) may be given; otherwise the first real embedding is used (with the specified precision) if any, else the first complex embedding.

INPUT:

  • embedding – an embedding of the base field into \(\RR\) or \(\CC\)

  • precision – a positive integer (default: 100) setting the number of bits of precision for the computation

  • algorithm – either 'pari' (default: for real embeddings) to use PARI’s pari:ellpointtoz, or 'sage' for a native implementation. Ignored for complex embeddings.

ALGORITHM:

See [Coh1993] for the case of real embeddings, and Cremona, J.E. and Thongjunthug, T. 2010 for the complex case.

AUTHORS:

  • Michael Mardaus (2008-07),

  • Tobias Nagel (2008-07) – original version from [Coh1993].

  • John Cremona (2008-07) – revision following eclib code.

  • John Cremona (2010-03) – implementation for complex embeddings.

EXAMPLES:

sage: E = EllipticCurve('389a')
sage: E.discriminant() > 0
True
sage: P = E([-1,1])
sage: P.is_on_identity_component ()
False
sage: P.elliptic_logarithm (precision=96)
0.4793482501902193161295330101 + 0.985868850775824102211203849...*I
sage: Q = E([3,5])
sage: Q.is_on_identity_component()
True
sage: Q.elliptic_logarithm (precision=96)
1.931128271542559442488585220
>>> from sage.all import *
>>> E = EllipticCurve('389a')
>>> E.discriminant() > Integer(0)
True
>>> P = E([-Integer(1),Integer(1)])
>>> P.is_on_identity_component ()
False
>>> P.elliptic_logarithm (precision=Integer(96))
0.4793482501902193161295330101 + 0.985868850775824102211203849...*I
>>> Q = E([Integer(3),Integer(5)])
>>> Q.is_on_identity_component()
True
>>> Q.elliptic_logarithm (precision=Integer(96))
1.931128271542559442488585220

An example with negative discriminant, and a torsion point:

sage: E = EllipticCurve('11a1')
sage: E.discriminant() < 0
True
sage: P = E([16,-61])
sage: P.elliptic_logarithm(precision=70)
0.25384186085591068434
sage: E.period_lattice().real_period(prec=70) / P.elliptic_logarithm(precision=70)
5.0000000000000000000
>>> from sage.all import *
>>> E = EllipticCurve('11a1')
>>> E.discriminant() < Integer(0)
True
>>> P = E([Integer(16),-Integer(61)])
>>> P.elliptic_logarithm(precision=Integer(70))
0.25384186085591068434
>>> E.period_lattice().real_period(prec=Integer(70)) / P.elliptic_logarithm(precision=Integer(70))
5.0000000000000000000

A larger example. The default algorithm uses PARI and makes sure the result has the requested precision:

sage: E = EllipticCurve([1, 0, 1, -85357462, 303528987048]) #18074g1
sage: P = E([4458713781401/835903744, -64466909836503771/24167649046528, 1])
sage: P.elliptic_logarithm()  # 100 bits
0.27656204014107061464076203097
>>> from sage.all import *
>>> E = EllipticCurve([Integer(1), Integer(0), Integer(1), -Integer(85357462), Integer(303528987048)]) #18074g1
>>> P = E([Integer(4458713781401)/Integer(835903744), -Integer(64466909836503771)/Integer(24167649046528), Integer(1)])
>>> P.elliptic_logarithm()  # 100 bits
0.27656204014107061464076203097

The native algorithm 'sage' used to have trouble with precision in this example, but no longer:

sage: P.elliptic_logarithm(algorithm='sage')  # 100 bits
0.27656204014107061464076203097
>>> from sage.all import *
>>> P.elliptic_logarithm(algorithm='sage')  # 100 bits
0.27656204014107061464076203097

This shows that the bug reported at Issue #4901 has been fixed:

sage: E = EllipticCurve("4390c2")
sage: P = E(683762969925/44944,-565388972095220019/9528128)
sage: P.elliptic_logarithm()
0.00025638725886520225353198932529
sage: P.elliptic_logarithm(precision=64)
0.000256387258865202254
sage: P.elliptic_logarithm(precision=65)
0.0002563872588652022535
sage: P.elliptic_logarithm(precision=128)
0.00025638725886520225353198932528666427412
sage: P.elliptic_logarithm(precision=129)
0.00025638725886520225353198932528666427412
sage: P.elliptic_logarithm(precision=256)
0.0002563872588652022535319893252866642741168388008346370015005142128009610936373
sage: P.elliptic_logarithm(precision=257)
0.00025638725886520225353198932528666427411683880083463700150051421280096109363730
>>> from sage.all import *
>>> E = EllipticCurve("4390c2")
>>> P = E(Integer(683762969925)/Integer(44944),-Integer(565388972095220019)/Integer(9528128))
>>> P.elliptic_logarithm()
0.00025638725886520225353198932529
>>> P.elliptic_logarithm(precision=Integer(64))
0.000256387258865202254
>>> P.elliptic_logarithm(precision=Integer(65))
0.0002563872588652022535
>>> P.elliptic_logarithm(precision=Integer(128))
0.00025638725886520225353198932528666427412
>>> P.elliptic_logarithm(precision=Integer(129))
0.00025638725886520225353198932528666427412
>>> P.elliptic_logarithm(precision=Integer(256))
0.0002563872588652022535319893252866642741168388008346370015005142128009610936373
>>> P.elliptic_logarithm(precision=Integer(257))
0.00025638725886520225353198932528666427411683880083463700150051421280096109363730

Examples over number fields:

sage: # needs sage.rings.number_field
sage: x = polygen(ZZ, 'x')
sage: K.<a> = NumberField(x^3 - 2)
sage: embs = K.embeddings(CC)
sage: E = EllipticCurve([0,1,0,a,a])
sage: Ls = [E.period_lattice(e) for e in embs]
sage: [L.real_flag for L in Ls]
[0, 0, -1]
sage: P = E(-1,0)  # order 2
sage: [L.elliptic_logarithm(P) for L in Ls]
[-1.73964256006716 - 1.07861534489191*I,
 -0.363756518406398 - 1.50699412135253*I, 1.90726488608927]

sage: # needs sage.rings.number_field
sage: E = EllipticCurve([-a^2 - a - 1, a^2 + a])
sage: Ls = [E.period_lattice(e) for e in embs]
sage: 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 )]
sage: [[L.elliptic_logarithm(P) for P in pts] for L in Ls]
[[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]]
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> x = polygen(ZZ, 'x')
>>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1)
>>> embs = K.embeddings(CC)
>>> E = EllipticCurve([Integer(0),Integer(1),Integer(0),a,a])
>>> Ls = [E.period_lattice(e) for e in embs]
>>> [L.real_flag for L in Ls]
[0, 0, -1]
>>> P = E(-Integer(1),Integer(0))  # order 2
>>> [L.elliptic_logarithm(P) for L in Ls]
[-1.73964256006716 - 1.07861534489191*I,
 -0.363756518406398 - 1.50699412135253*I, 1.90726488608927]

>>> # needs sage.rings.number_field
>>> E = EllipticCurve([-a**Integer(2) - a - Integer(1), a**Integer(2) + a])
>>> Ls = [E.period_lattice(e) for e in embs]
>>> pts = [E(Integer(2)*a**Integer(2) - a - Integer(1) , -Integer(2)*a**Integer(2) - Integer(2)*a + Integer(6) ),
...        E(-Integer(2)/Integer(3)*a**Integer(2) - Integer(1)/Integer(3) , -Integer(4)/Integer(3)*a - Integer(2)/Integer(3) ),
...        E(Integer(5)/Integer(4)*a**Integer(2) - Integer(1)/Integer(2)*a , -a**Integer(2) - Integer(1)/Integer(4)*a + Integer(9)/Integer(4) ),
...        E(Integer(2)*a**Integer(2) + Integer(3)*a + Integer(4) , -Integer(7)*a**Integer(2) - Integer(10)*a - Integer(12) )]
>>> [[L.elliptic_logarithm(P) for P in pts] for L in Ls]
[[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]]

sage: # needs sage.rings.number_field
sage: K.<i> = QuadraticField(-1)
sage: E = EllipticCurve([0,0,0,9*i-10,21-i])
sage: emb = K.embeddings(CC)[1]
sage: L = E.period_lattice(emb)
sage: P = E(2-i, 4+2*i)
sage: L.elliptic_logarithm(P, prec=100)
0.70448375537782208460499649302 - 0.79246725643650979858266018068*I
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> K = QuadraticField(-Integer(1), names=('i',)); (i,) = K._first_ngens(1)
>>> E = EllipticCurve([Integer(0),Integer(0),Integer(0),Integer(9)*i-Integer(10),Integer(21)-i])
>>> emb = K.embeddings(CC)[Integer(1)]
>>> L = E.period_lattice(emb)
>>> P = E(Integer(2)-i, Integer(4)+Integer(2)*i)
>>> L.elliptic_logarithm(P, prec=Integer(100))
0.70448375537782208460499649302 - 0.79246725643650979858266018068*I
has_finite_order()[source]

Return True iff this point has finite order on the elliptic curve.

EXAMPLES:

sage: E = EllipticCurve([0,0,1,-1,0])
sage: P = E([0,0]); P
(0 : 0 : 1)
sage: P.has_finite_order()
False
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0),Integer(0),Integer(1),-Integer(1),Integer(0)])
>>> P = E([Integer(0),Integer(0)]); P
(0 : 0 : 1)
>>> P.has_finite_order()
False

sage: E = EllipticCurve([0,1])
sage: P = E([-1,0])
sage: P.has_finite_order()
True
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0),Integer(1)])
>>> P = E([-Integer(1),Integer(0)])
>>> P.has_finite_order()
True
has_good_reduction(P=None)[source]

Return True iff this point has good reduction modulo a prime.

INPUT:

  • P – a prime of the base_field of the point’s curve, or None (default)

OUTPUT:

boolean; if a prime \(P\) of the base field is specified, returns True iff the point has good reduction at \(P\); otherwise, return True if the point has god reduction at all primes in the support of the discriminant of this model.

EXAMPLES:

sage: E = EllipticCurve('990e1')
sage: P = E.gen(0); P
(15 : 51 : 1)
sage: [E.has_good_reduction(p) for p in [2,3,5,7]]
[False, False, False, True]
sage: [P.has_good_reduction(p) for p in [2,3,5,7]]
[True, False, True, True]
sage: [E.tamagawa_exponent(p) for p in [2,3,5,7]]
[2, 2, 1, 1]
sage: [(2*P).has_good_reduction(p) for p in [2,3,5,7]]
[True, True, True, True]
sage: P.has_good_reduction()
False
sage: (2*P).has_good_reduction()
True
sage: (3*P).has_good_reduction()
False
>>> from sage.all import *
>>> E = EllipticCurve('990e1')
>>> P = E.gen(Integer(0)); P
(15 : 51 : 1)
>>> [E.has_good_reduction(p) for p in [Integer(2),Integer(3),Integer(5),Integer(7)]]
[False, False, False, True]
>>> [P.has_good_reduction(p) for p in [Integer(2),Integer(3),Integer(5),Integer(7)]]
[True, False, True, True]
>>> [E.tamagawa_exponent(p) for p in [Integer(2),Integer(3),Integer(5),Integer(7)]]
[2, 2, 1, 1]
>>> [(Integer(2)*P).has_good_reduction(p) for p in [Integer(2),Integer(3),Integer(5),Integer(7)]]
[True, True, True, True]
>>> P.has_good_reduction()
False
>>> (Integer(2)*P).has_good_reduction()
True
>>> (Integer(3)*P).has_good_reduction()
False

sage: # needs sage.rings.number_field
sage: x = polygen(ZZ, 'x')
sage: K.<i> = NumberField(x^2 + 1)
sage: E = EllipticCurve(K, [0,1,0,-160,308])
sage: P = E(26, -120)
sage: E.discriminant().support()
[Fractional ideal (i + 1),
 Fractional ideal (-i - 2),
 Fractional ideal (2*i + 1),
 Fractional ideal (3)]
sage: [E.tamagawa_exponent(p) for p in E.discriminant().support()]
[1, 4, 4, 4]
sage: P.has_good_reduction()
False
sage: (2*P).has_good_reduction()
False
sage: (4*P).has_good_reduction()
True
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> x = polygen(ZZ, 'x')
>>> K = NumberField(x**Integer(2) + Integer(1), names=('i',)); (i,) = K._first_ngens(1)
>>> E = EllipticCurve(K, [Integer(0),Integer(1),Integer(0),-Integer(160),Integer(308)])
>>> P = E(Integer(26), -Integer(120))
>>> E.discriminant().support()
[Fractional ideal (i + 1),
 Fractional ideal (-i - 2),
 Fractional ideal (2*i + 1),
 Fractional ideal (3)]
>>> [E.tamagawa_exponent(p) for p in E.discriminant().support()]
[1, 4, 4, 4]
>>> P.has_good_reduction()
False
>>> (Integer(2)*P).has_good_reduction()
False
>>> (Integer(4)*P).has_good_reduction()
True
has_infinite_order()[source]

Return True iff this point has infinite order on the elliptic curve.

EXAMPLES:

sage: E = EllipticCurve([0,0,1,-1,0])
sage: P = E([0,0]); P
(0 : 0 : 1)
sage: P.has_infinite_order()
True
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0),Integer(0),Integer(1),-Integer(1),Integer(0)])
>>> P = E([Integer(0),Integer(0)]); P
(0 : 0 : 1)
>>> P.has_infinite_order()
True

sage: E = EllipticCurve([0,1])
sage: P = E([-1,0])
sage: P.has_infinite_order()
False
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0),Integer(1)])
>>> P = E([-Integer(1),Integer(0)])
>>> P.has_infinite_order()
False
height(precision=None, normalised=True, algorithm='pari')[source]

Return the Néron-Tate canonical height of the point.

INPUT:

  • self – a point on an elliptic curve over a number field \(K\)

  • precision – positive integer, or None (default). The precision in bits of the result. If None, the default real precision is used.

  • normalised – boolean. If True (default), the height is normalised to be invariant under extension of \(K\). If False, return this normalised height multiplied by the degree of \(K\).

  • algorithm – string; either 'pari' (default) or 'sage'. If 'pari' and the base field is \(\QQ\), use the PARI library function; otherwise use the Sage implementation.

OUTPUT:

The rational number 0, or a nonnegative real number.

There are two normalisations used in the literature, one of which is double the other. We use the larger of the two, which is the one appropriate for the BSD conjecture. This is consistent with [Cre1997] and double that of [Sil2009].

See Wikipedia article Néron-Tate height.

Note

The correct height to use for the regulator in the BSD formula is the non-normalised height.

EXAMPLES:

sage: E = EllipticCurve('11a'); E
Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field
sage: P = E([5,5]); P
(5 : 5 : 1)
sage: P.height()
0
sage: Q = 5*P
sage: Q.height()
0
>>> from sage.all import *
>>> E = EllipticCurve('11a'); E
Elliptic Curve defined by y^2 + y = x^3 - x^2 - 10*x - 20 over Rational Field
>>> P = E([Integer(5),Integer(5)]); P
(5 : 5 : 1)
>>> P.height()
0
>>> Q = Integer(5)*P
>>> Q.height()
0

sage: E = EllipticCurve('37a'); E
Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
sage: P = E([0,0])
sage: P.height()
0.0511114082399688
sage: P.order()
+Infinity
sage: E.regulator()
0.0511114082399688...

sage: def naive_height(P):
....:     return log(RR(max(abs(P[0].numerator()), abs(P[0].denominator()))))
sage: for n in [1..10]:
....:     print(naive_height(2^n*P)/4^n)
0.000000000000000
0.0433216987849966
0.0502949347635656
0.0511006335618645
0.0511007834799612
0.0511013666152466
0.0511034199907743
0.0511106492906471
0.0511114081541082
0.0511114081541180
>>> from sage.all import *
>>> E = EllipticCurve('37a'); E
Elliptic Curve defined by y^2 + y = x^3 - x over Rational Field
>>> P = E([Integer(0),Integer(0)])
>>> P.height()
0.0511114082399688
>>> P.order()
+Infinity
>>> E.regulator()
0.0511114082399688...

>>> def naive_height(P):
...     return log(RR(max(abs(P[Integer(0)].numerator()), abs(P[Integer(0)].denominator()))))
>>> for n in (ellipsis_range(Integer(1),Ellipsis,Integer(10))):
...     print(naive_height(Integer(2)**n*P)/Integer(4)**n)
0.000000000000000
0.0433216987849966
0.0502949347635656
0.0511006335618645
0.0511007834799612
0.0511013666152466
0.0511034199907743
0.0511106492906471
0.0511114081541082
0.0511114081541180

sage: E = EllipticCurve('4602a1'); E
Elliptic Curve defined by y^2 + x*y  = x^3 + x^2 - 37746035*x - 89296920339
over Rational Field
sage: x = 77985922458974949246858229195945103471590
sage: y = 19575260230015313702261379022151675961965157108920263594545223
sage: d = 2254020761884782243
sage: E([ x / d^2,  y / d^3 ]).height()
86.7406561381275
>>> from sage.all import *
>>> E = EllipticCurve('4602a1'); E
Elliptic Curve defined by y^2 + x*y  = x^3 + x^2 - 37746035*x - 89296920339
over Rational Field
>>> x = Integer(77985922458974949246858229195945103471590)
>>> y = Integer(19575260230015313702261379022151675961965157108920263594545223)
>>> d = Integer(2254020761884782243)
>>> E([ x / d**Integer(2),  y / d**Integer(3) ]).height()
86.7406561381275

sage: E = EllipticCurve([17, -60, -120, 0, 0]); E
Elliptic Curve defined by y^2 + 17*x*y - 120*y = x^3 - 60*x^2 over Rational Field
sage: E([30, -90]).height()
0

sage: E = EllipticCurve('389a1'); E
Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field
sage: P, Q = E(-1,1), E(0,-1)
sage: P.height(precision=100)
0.68666708330558658572355210295
sage: (3*Q).height(precision=100)/Q.height(precision=100)
9.0000000000000000000000000000
sage: _.parent()
Real Field with 100 bits of precision
>>> from sage.all import *
>>> E = EllipticCurve([Integer(17), -Integer(60), -Integer(120), Integer(0), Integer(0)]); E
Elliptic Curve defined by y^2 + 17*x*y - 120*y = x^3 - 60*x^2 over Rational Field
>>> E([Integer(30), -Integer(90)]).height()
0

>>> E = EllipticCurve('389a1'); E
Elliptic Curve defined by y^2 + y = x^3 + x^2 - 2*x over Rational Field
>>> P, Q = E(-Integer(1),Integer(1)), E(Integer(0),-Integer(1))
>>> P.height(precision=Integer(100))
0.68666708330558658572355210295
>>> (Integer(3)*Q).height(precision=Integer(100))/Q.height(precision=Integer(100))
9.0000000000000000000000000000
>>> _.parent()
Real Field with 100 bits of precision

Canonical heights over number fields are implemented as well:

sage: R.<x> = QQ[]
sage: K.<a> = NumberField(x^3 - 2)                                          # needs sage.rings.number_field
sage: E = EllipticCurve([a, 4]); E                                          # needs sage.rings.number_field
Elliptic Curve defined by y^2 = x^3 + a*x + 4
 over Number Field in a with defining polynomial x^3 - 2
sage: P = E((0,2))                                                          # needs sage.rings.number_field
sage: P.height()                                                            # needs sage.rings.number_field
0.810463096585925
sage: P.height(precision=100)                                               # needs sage.rings.number_field
0.81046309658592536863991810577
sage: P.height(precision=200)                                               # needs sage.rings.number_field
0.81046309658592536863991810576865158896130286417155832378086
sage: (2*P).height() / P.height()                                           # needs sage.rings.number_field
4.00000000000000
sage: (100*P).height() / P.height()                                         # needs sage.rings.number_field
10000.0000000000
>>> from sage.all import *
>>> R = QQ['x']; (x,) = R._first_ngens(1)
>>> K = NumberField(x**Integer(3) - Integer(2), names=('a',)); (a,) = K._first_ngens(1)# needs sage.rings.number_field
>>> E = EllipticCurve([a, Integer(4)]); E                                          # needs sage.rings.number_field
Elliptic Curve defined by y^2 = x^3 + a*x + 4
 over Number Field in a with defining polynomial x^3 - 2
>>> P = E((Integer(0),Integer(2)))                                                          # needs sage.rings.number_field
>>> P.height()                                                            # needs sage.rings.number_field
0.810463096585925
>>> P.height(precision=Integer(100))                                               # needs sage.rings.number_field
0.81046309658592536863991810577
>>> P.height(precision=Integer(200))                                               # needs sage.rings.number_field
0.81046309658592536863991810576865158896130286417155832378086
>>> (Integer(2)*P).height() / P.height()                                           # needs sage.rings.number_field
4.00000000000000
>>> (Integer(100)*P).height() / P.height()                                         # needs sage.rings.number_field
10000.0000000000

Setting normalised=False multiplies the height by the degree of \(K\):

sage: E = EllipticCurve('37a')
sage: P = E([0,0])
sage: P.height()
0.0511114082399688
sage: P.height(normalised=False)
0.0511114082399688
sage: K.<z> = CyclotomicField(5)                                            # needs sage.rings.number_field
sage: EK = E.change_ring(K)                                                 # needs sage.rings.number_field
sage: PK = EK([0,0])                                                        # needs sage.rings.number_field
sage: PK.height()                                                           # needs sage.rings.number_field
0.0511114082399688
sage: PK.height(normalised=False)                                           # needs sage.rings.number_field
0.204445632959875
>>> from sage.all import *
>>> E = EllipticCurve('37a')
>>> P = E([Integer(0),Integer(0)])
>>> P.height()
0.0511114082399688
>>> P.height(normalised=False)
0.0511114082399688
>>> K = CyclotomicField(Integer(5), names=('z',)); (z,) = K._first_ngens(1)# needs sage.rings.number_field
>>> EK = E.change_ring(K)                                                 # needs sage.rings.number_field
>>> PK = EK([Integer(0),Integer(0)])                                                        # needs sage.rings.number_field
>>> PK.height()                                                           # needs sage.rings.number_field
0.0511114082399688
>>> PK.height(normalised=False)                                           # needs sage.rings.number_field
0.204445632959875

Some consistency checks:

sage: E = EllipticCurve('5077a1')
sage: P = E([-2,3,1])
sage: P.height()
1.36857250535393

sage: EK = E.change_ring(QuadraticField(-3,'a'))                            # needs sage.rings.number_field
sage: PK = EK([-2,3,1])                                                     # needs sage.rings.number_field
sage: PK.height()                                                           # needs sage.rings.number_field
1.36857250535393

sage: # needs sage.rings.number_field
sage: K.<i> = NumberField(x^2 + 1)
sage: E = EllipticCurve(K, [0,0,4,6*i,0])
sage: Q = E.lift_x(-9/4); Q
(-9/4 : 27/8*i - 4 : 1)
sage: Q.height()
2.69518560017909
sage: (15*Q).height() / Q.height()
225.000000000000

sage: E = EllipticCurve('37a')
sage: P = E([0,-1])
sage: P.height()
0.0511114082399688
sage: K.<a> = QuadraticField(-7)                                            # needs sage.rings.number_field
sage: ED = E.quadratic_twist(-7)                                            # needs sage.rings.number_field
sage: Q = E.isomorphism_to(ED.change_ring(K))(P); Q                         # needs sage.rings.number_field
(0 : -7/2*a - 1/2 : 1)
sage: Q.height()                                                            # needs sage.rings.number_field
0.0511114082399688
sage: Q.height(precision=100)                                               # needs sage.rings.number_field
0.051111408239968840235886099757
>>> from sage.all import *
>>> E = EllipticCurve('5077a1')
>>> P = E([-Integer(2),Integer(3),Integer(1)])
>>> P.height()
1.36857250535393

>>> EK = E.change_ring(QuadraticField(-Integer(3),'a'))                            # needs sage.rings.number_field
>>> PK = EK([-Integer(2),Integer(3),Integer(1)])                                                     # needs sage.rings.number_field
>>> PK.height()                                                           # needs sage.rings.number_field
1.36857250535393

>>> # needs sage.rings.number_field
>>> K = NumberField(x**Integer(2) + Integer(1), names=('i',)); (i,) = K._first_ngens(1)
>>> E = EllipticCurve(K, [Integer(0),Integer(0),Integer(4),Integer(6)*i,Integer(0)])
>>> Q = E.lift_x(-Integer(9)/Integer(4)); Q
(-9/4 : 27/8*i - 4 : 1)
>>> Q.height()
2.69518560017909
>>> (Integer(15)*Q).height() / Q.height()
225.000000000000

>>> E = EllipticCurve('37a')
>>> P = E([Integer(0),-Integer(1)])
>>> P.height()
0.0511114082399688
>>> K = QuadraticField(-Integer(7), names=('a',)); (a,) = K._first_ngens(1)# needs sage.rings.number_field
>>> ED = E.quadratic_twist(-Integer(7))                                            # needs sage.rings.number_field
>>> Q = E.isomorphism_to(ED.change_ring(K))(P); Q                         # needs sage.rings.number_field
(0 : -7/2*a - 1/2 : 1)
>>> Q.height()                                                            # needs sage.rings.number_field
0.0511114082399688
>>> Q.height(precision=Integer(100))                                               # needs sage.rings.number_field
0.051111408239968840235886099757

An example to show that the bug at Issue #5252 is fixed:

sage: E = EllipticCurve([1, -1, 1, -2063758701246626370773726978, 32838647793306133075103747085833809114881])
sage: P = E([-30987785091199, 258909576181697016447])
sage: P.height()
25.8603170675462
sage: P.height(precision=100)
25.860317067546190743868840741
sage: P.height(precision=250)
25.860317067546190743868840740735110323098872903844416215577171041783572513
sage: P.height(precision=500)
25.8603170675461907438688407407351103230988729038444162155771710417835725129551130570889813281792157278507639909972112856019190236125362914195452321720

sage: P.height(precision=100) == P.non_archimedean_local_height(prec=100)+P.archimedean_local_height(prec=100)
True
>>> from sage.all import *
>>> E = EllipticCurve([Integer(1), -Integer(1), Integer(1), -Integer(2063758701246626370773726978), Integer(32838647793306133075103747085833809114881)])
>>> P = E([-Integer(30987785091199), Integer(258909576181697016447)])
>>> P.height()
25.8603170675462
>>> P.height(precision=Integer(100))
25.860317067546190743868840741
>>> P.height(precision=Integer(250))
25.860317067546190743868840740735110323098872903844416215577171041783572513
>>> P.height(precision=Integer(500))
25.8603170675461907438688407407351103230988729038444162155771710417835725129551130570889813281792157278507639909972112856019190236125362914195452321720

>>> P.height(precision=Integer(100)) == P.non_archimedean_local_height(prec=Integer(100))+P.archimedean_local_height(prec=Integer(100))
True

An example to show that the bug at Issue #8319 is fixed (correct height when the curve is not minimal):

sage: E = EllipticCurve([-5580472329446114952805505804593498080000,-157339733785368110382973689903536054787700497223306368000000])
sage: xP = 204885147732879546487576840131729064308289385547094673627174585676211859152978311600/23625501907057948132262217188983681204856907657753178415430361
sage: P = E.lift_x(xP)
sage: P.height()
157.432598516754
sage: Q = 2*P
sage: Q.height() # long time (4s)
629.730394067016
sage: Q.height()-4*P.height() # long time
0.000000000000000
>>> from sage.all import *
>>> E = EllipticCurve([-Integer(5580472329446114952805505804593498080000),-Integer(157339733785368110382973689903536054787700497223306368000000)])
>>> xP = Integer(204885147732879546487576840131729064308289385547094673627174585676211859152978311600)/Integer(23625501907057948132262217188983681204856907657753178415430361)
>>> P = E.lift_x(xP)
>>> P.height()
157.432598516754
>>> Q = Integer(2)*P
>>> Q.height() # long time (4s)
629.730394067016
>>> Q.height()-Integer(4)*P.height() # long time
0.000000000000000

An example to show that the bug at Issue #12509 is fixed (precision issues):

sage: # needs sage.rings.number_field
sage: x = polygen(QQ)
sage: K.<a> = NumberField(x^2 - x - 1)
sage: v = [0, a + 1, 1, 28665*a - 46382, 2797026*a - 4525688]
sage: E = EllipticCurve(v)
sage: P = E([72*a - 509/5,  -682/25*a - 434/25])
sage: P.height()
1.38877711688727
sage: (2*P).height()/P.height()
4.00000000000000
sage: (2*P).height(precision=100)/P.height(precision=100)
4.0000000000000000000000000000
sage: (2*P).height(precision=1000)/P.height(precision=1000)
4.00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> x = polygen(QQ)
>>> K = NumberField(x**Integer(2) - x - Integer(1), names=('a',)); (a,) = K._first_ngens(1)
>>> v = [Integer(0), a + Integer(1), Integer(1), Integer(28665)*a - Integer(46382), Integer(2797026)*a - Integer(4525688)]
>>> E = EllipticCurve(v)
>>> P = E([Integer(72)*a - Integer(509)/Integer(5),  -Integer(682)/Integer(25)*a - Integer(434)/Integer(25)])
>>> P.height()
1.38877711688727
>>> (Integer(2)*P).height()/P.height()
4.00000000000000
>>> (Integer(2)*P).height(precision=Integer(100))/P.height(precision=Integer(100))
4.0000000000000000000000000000
>>> (Integer(2)*P).height(precision=Integer(1000))/P.height(precision=Integer(1000))
4.00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000

This shows that the bug reported at Issue #13951 has been fixed:

sage: E = EllipticCurve([0,17])
sage: P1 = E(2,5)
sage: P1.height()
1.06248137652528
sage: F = E.change_ring(QuadraticField(-3, 'a'))                            # needs sage.rings.number_field
sage: P2 = F([2,5])                                                         # needs sage.rings.number_field
sage: P2.height()                                                           # needs sage.rings.number_field
1.06248137652528
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0),Integer(17)])
>>> P1 = E(Integer(2),Integer(5))
>>> P1.height()
1.06248137652528
>>> F = E.change_ring(QuadraticField(-Integer(3), 'a'))                            # needs sage.rings.number_field
>>> P2 = F([Integer(2),Integer(5)])                                                         # needs sage.rings.number_field
>>> P2.height()                                                           # needs sage.rings.number_field
1.06248137652528

This shows that the bug reported at Issue #36834 (incorrect value when the model is not integral) has been fixed:

sage: # needs sage.rings.number_field
sage: K.<a> = NumberField(x^2 - 84131656042917)
sage: E = EllipticCurve(K, [0, 0, 0, -5482707841/48, -244634179112639/864])
sage: P = E(349189/12, 1/2*a)
sage: P.height()
10.4560438181991
sage: [(n*P).height()/P.height() for n in [2,3,4,5]]
[4.00000000000000, 9.00000000000000, 16.0000000000000, 25.0000000000000]
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> K = NumberField(x**Integer(2) - Integer(84131656042917), names=('a',)); (a,) = K._first_ngens(1)
>>> E = EllipticCurve(K, [Integer(0), Integer(0), Integer(0), -Integer(5482707841)/Integer(48), -Integer(244634179112639)/Integer(864)])
>>> P = E(Integer(349189)/Integer(12), Integer(1)/Integer(2)*a)
>>> P.height()
10.4560438181991
>>> [(n*P).height()/P.height() for n in [Integer(2),Integer(3),Integer(4),Integer(5)]]
[4.00000000000000, 9.00000000000000, 16.0000000000000, 25.0000000000000]
is_on_identity_component(embedding=None)[source]

Return True iff this point is on the identity component of its curve with respect to a given (real or complex) embedding.

INPUT:

  • self – a point on a curve over any ordered field (e.g. \(\QQ\))

  • embedding – an embedding from the base_field of the point’s curve into \(\RR\) or \(\CC\); if None (the default) it uses the first embedding of the base_field into \(\RR\) if any, else the first embedding into \(\CC\).

OUTPUT:

boolean; True iff the point is on the identity component of the curve. (If the point is zero then the result is True.)

EXAMPLES:

For \(K=\QQ\) there is no need to specify an embedding:

sage: E = EllipticCurve('5077a1')
sage: [E.lift_x(x).is_on_identity_component() for x in srange(-3,5)]
[False, False, False, False, False, True, True, True]
>>> from sage.all import *
>>> E = EllipticCurve('5077a1')
>>> [E.lift_x(x).is_on_identity_component() for x in srange(-Integer(3),Integer(5))]
[False, False, False, False, False, True, True, True]

An example over a field with two real embeddings:

sage: # needs sage.rings.number_field
sage: L.<a> = QuadraticField(2)
sage: E = EllipticCurve(L, [0,1,0,a,a])
sage: P = E(-1,0)
sage: [P.is_on_identity_component(e) for e in L.embeddings(RR)]
[False, True]
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> L = QuadraticField(Integer(2), names=('a',)); (a,) = L._first_ngens(1)
>>> E = EllipticCurve(L, [Integer(0),Integer(1),Integer(0),a,a])
>>> P = E(-Integer(1),Integer(0))
>>> [P.is_on_identity_component(e) for e in L.embeddings(RR)]
[False, True]

We can check this as follows:

sage: # needs sage.rings.number_field
sage: [e(E.discriminant()) > 0 for e in L.embeddings(RR)]
[True, False]
sage: e = L.embeddings(RR)[0]
sage: E1 = EllipticCurve(RR, [e(ai) for ai in E.ainvs()])
sage: e1, e2, e3 = E1.two_division_polynomial().roots(RR,
....:                                                 multiplicities=False)
sage: e1 < e2 < e3 and e(P[0]) < e3
True
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> [e(E.discriminant()) > Integer(0) for e in L.embeddings(RR)]
[True, False]
>>> e = L.embeddings(RR)[Integer(0)]
>>> E1 = EllipticCurve(RR, [e(ai) for ai in E.ainvs()])
>>> e1, e2, e3 = E1.two_division_polynomial().roots(RR,
...                                                 multiplicities=False)
>>> e1 < e2 < e3 and e(P[Integer(0)]) < e3
True
non_archimedean_local_height(v=None, prec=None, weighted=False, is_minimal=None)[source]

Compute the local height of self at non-archimedean places.

INPUT:

  • self – a point on an elliptic curve over a number field \(K\)

  • v – a non-archimedean place of \(K\), or None (default). If \(v\) is a non-archimedean place, return the local height of self at \(v\). If \(v\) is None, return the total non-archimedean contribution to the global height.

  • prec – integer; or None (default). The precision of the computation. If None, the height is returned symbolically

  • weighted – boolean. If False (default), the height is normalised to be invariant under extension of \(K\). If True, return this normalised height multiplied by the local degree if \(v\) is a single place, or by the degree of \(K\) if \(v\) is None.

  • is_minimal – boolean, or None (default). Ignored when v is None (default) or True. Otherwise, when the place v is specified: if True, the model is assumed to be both locally integral and a local minimal model; if None (default) or False, a local minimal model is computed and the computation is done on that model.

OUTPUT:

A real number. The normalisation is twice that in Silverman’s paper [Sil1988]. Note that this local height depends on the model of the curve.

ALGORITHM:

See [Sil1988], Section 5.

EXAMPLES:

Examples 2 and 3 from [Sil1988]:

sage: # needs sage.rings.number_field
sage: x = polygen(ZZ, 'x')
sage: K.<i> = NumberField(x^2 + 1)
sage: E = EllipticCurve(K, [0,0,4,6*i,0]); E
Elliptic Curve defined by y^2 + 4*y = x^3 + 6*i*x
 over Number Field in i with defining polynomial x^2 + 1
sage: P = E((0,0))
sage: P.non_archimedean_local_height(K.ideal(i+1))
-1/2*log(2)
sage: P.non_archimedean_local_height(K.ideal(3))
0
sage: P.non_archimedean_local_height(K.ideal(1-2*i))
0

sage: # needs sage.rings.number_field
sage: Q = E.lift_x(-9/4); Q
(-9/4 : 27/8*i - 4 : 1)
sage: Q.non_archimedean_local_height(K.ideal(1+i))
2*log(2)
sage: Q.non_archimedean_local_height(K.ideal(3))
0
sage: Q.non_archimedean_local_height(K.ideal(1-2*i))
0
sage: Q.non_archimedean_local_height()
2*log(2)
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> x = polygen(ZZ, 'x')
>>> K = NumberField(x**Integer(2) + Integer(1), names=('i',)); (i,) = K._first_ngens(1)
>>> E = EllipticCurve(K, [Integer(0),Integer(0),Integer(4),Integer(6)*i,Integer(0)]); E
Elliptic Curve defined by y^2 + 4*y = x^3 + 6*i*x
 over Number Field in i with defining polynomial x^2 + 1
>>> P = E((Integer(0),Integer(0)))
>>> P.non_archimedean_local_height(K.ideal(i+Integer(1)))
-1/2*log(2)
>>> P.non_archimedean_local_height(K.ideal(Integer(3)))
0
>>> P.non_archimedean_local_height(K.ideal(Integer(1)-Integer(2)*i))
0

>>> # needs sage.rings.number_field
>>> Q = E.lift_x(-Integer(9)/Integer(4)); Q
(-9/4 : 27/8*i - 4 : 1)
>>> Q.non_archimedean_local_height(K.ideal(Integer(1)+i))
2*log(2)
>>> Q.non_archimedean_local_height(K.ideal(Integer(3)))
0
>>> Q.non_archimedean_local_height(K.ideal(Integer(1)-Integer(2)*i))
0
>>> Q.non_archimedean_local_height()
2*log(2)

An example over the rational numbers:

sage: E = EllipticCurve([0, 0, 0, -36, 0])
sage: P = E([-3, 9])
sage: P.non_archimedean_local_height()
-log(3)
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0), Integer(0), Integer(0), -Integer(36), Integer(0)])
>>> P = E([-Integer(3), Integer(9)])
>>> P.non_archimedean_local_height()
-log(3)

Local heights of torsion points can be nonzero (unlike the global height):

sage: # needs sage.rings.number_field
sage: K.<i> = QuadraticField(-1)
sage: E = EllipticCurve([0, 0, 0, K(1), 0])
sage: P = E(i, 0)
sage: P.non_archimedean_local_height()
-1/2*log(2)
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> K = QuadraticField(-Integer(1), names=('i',)); (i,) = K._first_ngens(1)
>>> E = EllipticCurve([Integer(0), Integer(0), Integer(0), K(Integer(1)), Integer(0)])
>>> P = E(i, Integer(0))
>>> P.non_archimedean_local_height()
-1/2*log(2)
order()[source]

Return the order of this point on the elliptic curve.

If the point has infinite order, returns +Infinity. For curves defined over \(\QQ\), we call PARI; over other number fields we implement the function here.

Note

additive_order() is a synonym for order()

EXAMPLES:

sage: E = EllipticCurve([0,0,1,-1,0])
sage: P = E([0,0]); P
(0 : 0 : 1)
sage: P.order()
+Infinity
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0),Integer(0),Integer(1),-Integer(1),Integer(0)])
>>> P = E([Integer(0),Integer(0)]); P
(0 : 0 : 1)
>>> P.order()
+Infinity

sage: E = EllipticCurve([0,1])
sage: P = E([-1,0])
sage: P.order()
2
sage: P.additive_order()
2
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0),Integer(1)])
>>> P = E([-Integer(1),Integer(0)])
>>> P.order()
2
>>> P.additive_order()
2
padic_elliptic_logarithm(p, absprec=20)[source]

Compute the \(p\)-adic elliptic logarithm of this point.

INPUT:

  • p – integer; a prime

  • absprec – integer (default: 20); the initial \(p\)-adic absolute precision of the computation

OUTPUT:

The \(p\)-adic elliptic logarithm of self, with precision absprec.

AUTHORS:

  • Tobias Nagel

  • Michael Mardaus

  • John Cremona

ALGORITHM:

For points in the formal group (i.e. not integral at \(p\)) we take the log() function from the formal groups module and evaluate it at \(-x/y\). Otherwise we first multiply the point to get into the formal group, and divide the result afterwards.

Todo

See comments at Issue #4805. Currently the absolute precision of the result may be less than the given value of absprec, and error-handling is imperfect.

EXAMPLES:

sage: E = EllipticCurve([0,1,1,-2,0])
sage: E(0).padic_elliptic_logarithm(3)                                      # needs sage.rings.padics
0
sage: P = E(0, 0)                                                           # needs sage.rings.padics
sage: P.padic_elliptic_logarithm(3)                                         # needs sage.rings.padics
2 + 2*3 + 3^3 + 2*3^7 + 3^8 + 3^9 + 3^11 + 3^15 + 2*3^17 + 3^18 + O(3^19)
sage: P.padic_elliptic_logarithm(3).lift()                                  # needs sage.rings.padics
660257522
sage: P = E(-11/9, 28/27)                                                   # needs sage.rings.padics
sage: [(2*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(20)]  # long time, needs sage.rings.padics
[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)]
sage: [(3*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(12)]  # long time, needs sage.rings.padics
[1 + 2 + O(2^19), 3 + 3^20 + O(3^21), 3 + O(5^19), 3 + O(7^19), 3 + O(11^19)]
sage: [(5*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(12)]  # long time, needs sage.rings.padics
[1 + 2^2 + O(2^19), 2 + 3 + O(3^20), 5 + O(5^19), 5 + O(7^19), 5 + O(11^19)]
>>> from sage.all import *
>>> E = EllipticCurve([Integer(0),Integer(1),Integer(1),-Integer(2),Integer(0)])
>>> E(Integer(0)).padic_elliptic_logarithm(Integer(3))                                      # needs sage.rings.padics
0
>>> P = E(Integer(0), Integer(0))                                                           # needs sage.rings.padics
>>> P.padic_elliptic_logarithm(Integer(3))                                         # needs sage.rings.padics
2 + 2*3 + 3^3 + 2*3^7 + 3^8 + 3^9 + 3^11 + 3^15 + 2*3^17 + 3^18 + O(3^19)
>>> P.padic_elliptic_logarithm(Integer(3)).lift()                                  # needs sage.rings.padics
660257522
>>> P = E(-Integer(11)/Integer(9), Integer(28)/Integer(27))                                                   # needs sage.rings.padics
>>> [(Integer(2)*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(Integer(20))]  # long time, needs sage.rings.padics
[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)]
>>> [(Integer(3)*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(Integer(12))]  # long time, needs sage.rings.padics
[1 + 2 + O(2^19), 3 + 3^20 + O(3^21), 3 + O(5^19), 3 + O(7^19), 3 + O(11^19)]
>>> [(Integer(5)*P).padic_elliptic_logarithm(p)/P.padic_elliptic_logarithm(p) for p in prime_range(Integer(12))]  # long time, needs sage.rings.padics
[1 + 2^2 + O(2^19), 2 + 3 + O(3^20), 5 + O(5^19), 5 + O(7^19), 5 + O(11^19)]

An example which arose during reviewing Issue #4741:

sage: E = EllipticCurve('794a1')
sage: P = E(-1,2)
sage: P.padic_elliptic_logarithm(2)  # default precision=20                 # needs sage.rings.padics
2^4 + 2^5 + 2^6 + 2^8 + 2^9 + 2^13 + 2^14 + 2^15 + O(2^16)
sage: P.padic_elliptic_logarithm(2, absprec=30)                             # needs sage.rings.padics
2^4 + 2^5 + 2^6 + 2^8 + 2^9 + 2^13 + 2^14 + 2^15 + 2^22 + 2^23 + 2^24 + O(2^26)
sage: P.padic_elliptic_logarithm(2, absprec=40)                             # needs sage.rings.padics
2^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)
>>> from sage.all import *
>>> E = EllipticCurve('794a1')
>>> P = E(-Integer(1),Integer(2))
>>> P.padic_elliptic_logarithm(Integer(2))  # default precision=20                 # needs sage.rings.padics
2^4 + 2^5 + 2^6 + 2^8 + 2^9 + 2^13 + 2^14 + 2^15 + O(2^16)
>>> P.padic_elliptic_logarithm(Integer(2), absprec=Integer(30))                             # needs sage.rings.padics
2^4 + 2^5 + 2^6 + 2^8 + 2^9 + 2^13 + 2^14 + 2^15 + 2^22 + 2^23 + 2^24 + O(2^26)
>>> P.padic_elliptic_logarithm(Integer(2), absprec=Integer(40))                             # needs sage.rings.padics
2^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)
reduction(p)[source]

This finds the reduction of a point \(P\) on the elliptic curve modulo the prime \(p\).

INPUT:

  • self – a point on an elliptic curve

  • p – a prime number

OUTPUT: the point reduced to be a point on the elliptic curve modulo \(p\)

EXAMPLES:

sage: E = EllipticCurve([1,2,3,4,0])
sage: P = E(0,0)
sage: P.reduction(5)
(0 : 0 : 1)
sage: Q = E(98,931)
sage: Q.reduction(5)
(3 : 1 : 1)
sage: Q.reduction(5).curve() == E.reduction(5)
True
>>> from sage.all import *
>>> E = EllipticCurve([Integer(1),Integer(2),Integer(3),Integer(4),Integer(0)])
>>> P = E(Integer(0),Integer(0))
>>> P.reduction(Integer(5))
(0 : 0 : 1)
>>> Q = E(Integer(98),Integer(931))
>>> Q.reduction(Integer(5))
(3 : 1 : 1)
>>> Q.reduction(Integer(5)).curve() == E.reduction(Integer(5))
True

sage: # needs sage.rings.number_field
sage: x = polygen(ZZ, 'x')
sage: F.<a> = NumberField(x^2 + 5)
sage: E = EllipticCurve(F, [1,2,3,4,0])
sage: Q = E(98, 931)
sage: Q.reduction(a)
(3 : 1 : 1)
sage: Q.reduction(11)
(10 : 7 : 1)
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> x = polygen(ZZ, 'x')
>>> F = NumberField(x**Integer(2) + Integer(5), names=('a',)); (a,) = F._first_ngens(1)
>>> E = EllipticCurve(F, [Integer(1),Integer(2),Integer(3),Integer(4),Integer(0)])
>>> Q = E(Integer(98), Integer(931))
>>> Q.reduction(a)
(3 : 1 : 1)
>>> Q.reduction(Integer(11))
(10 : 7 : 1)

sage: # needs sage.rings.number_field
sage: F.<a> = NumberField(x^3 + x^2 + 1)
sage: E = EllipticCurve(F, [a,2])
sage: P = E(a, 1)
sage: P.reduction(F.ideal(5))
(abar : 1 : 1)
sage: P.reduction(F.ideal(a^2 - 4*a - 2))
(abar : 1 : 1)
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> F = NumberField(x**Integer(3) + x**Integer(2) + Integer(1), names=('a',)); (a,) = F._first_ngens(1)
>>> E = EllipticCurve(F, [a,Integer(2)])
>>> P = E(a, Integer(1))
>>> P.reduction(F.ideal(Integer(5)))
(abar : 1 : 1)
>>> P.reduction(F.ideal(a**Integer(2) - Integer(4)*a - Integer(2)))
(abar : 1 : 1)