Square‑root Vélu algorithm for elliptic-curve isogenies

The square-root Vélu algorithm, also called the √élu algorithm, computes isogenies of elliptic curves in time O~() rather than naïvely O(), where is the degree.

The core idea is to reindex the points in the kernel subgroup in a baby-step-giant-step manner, then use fast resultant computations to evaluate “elliptic polynomials” (see FastEllipticPolynomial) in essentially square-root time.

Based on experiments with Sage version 9.7, the isogeny degree where EllipticCurveHom_velusqrt begins to outperform EllipticCurveIsogeny can be as low as 100, but is typically closer to 1000, depending on the exact situation.

REFERENCES: [BDLS2020]

EXAMPLES:

sage: from sage.schemes.elliptic_curves.hom_velusqrt import EllipticCurveHom_velusqrt
sage: E = EllipticCurve(GF(6666679), [5,5])
sage: K = E(9970, 1003793, 1)
sage: K.order()
10009
sage: phi = EllipticCurveHom_velusqrt(E, K)
sage: phi
Elliptic-curve isogeny (using square-root Vélu) of degree 10009:
  From: Elliptic Curve defined by y^2 = x^3 + 5*x + 5 over Finite Field of size 6666679
  To:   Elliptic Curve defined by y^2 = x^3 + 227975*x + 3596133 over Finite Field of size 6666679
sage: phi.codomain()
Elliptic Curve defined by y^2 = x^3 + 227975*x + 3596133 over Finite Field of size 6666679
>>> from sage.all import *
>>> from sage.schemes.elliptic_curves.hom_velusqrt import EllipticCurveHom_velusqrt
>>> E = EllipticCurve(GF(Integer(6666679)), [Integer(5),Integer(5)])
>>> K = E(Integer(9970), Integer(1003793), Integer(1))
>>> K.order()
10009
>>> phi = EllipticCurveHom_velusqrt(E, K)
>>> phi
Elliptic-curve isogeny (using square-root Vélu) of degree 10009:
  From: Elliptic Curve defined by y^2 = x^3 + 5*x + 5 over Finite Field of size 6666679
  To:   Elliptic Curve defined by y^2 = x^3 + 227975*x + 3596133 over Finite Field of size 6666679
>>> phi.codomain()
Elliptic Curve defined by y^2 = x^3 + 227975*x + 3596133 over Finite Field of size 6666679
from sage.schemes.elliptic_curves.hom_velusqrt import EllipticCurveHom_velusqrt
E = EllipticCurve(GF(6666679), [5,5])
K = E(9970, 1003793, 1)
K.order()
phi = EllipticCurveHom_velusqrt(E, K)
phi
phi.codomain()

Note that the isogeny is usually not identical to the one computed by EllipticCurveIsogeny:

sage: psi = EllipticCurveIsogeny(E, K)
sage: psi
Isogeny of degree 10009
    from Elliptic Curve defined by y^2 = x^3 + 5*x + 5 over Finite Field of size 6666679
    to Elliptic Curve defined by y^2 = x^3 + 5344836*x + 3950273 over Finite Field of size 6666679
>>> from sage.all import *
>>> psi = EllipticCurveIsogeny(E, K)
>>> psi
Isogeny of degree 10009
    from Elliptic Curve defined by y^2 = x^3 + 5*x + 5 over Finite Field of size 6666679
    to Elliptic Curve defined by y^2 = x^3 + 5344836*x + 3950273 over Finite Field of size 6666679
psi = EllipticCurveIsogeny(E, K)
psi

However, they are certainly separable isogenies with the same kernel and must therefore be equal up to post-isomorphism:

sage: isos = psi.codomain().isomorphisms(phi.codomain())
sage: sum(iso * psi == phi for iso in isos)
1
>>> from sage.all import *
>>> isos = psi.codomain().isomorphisms(phi.codomain())
>>> sum(iso * psi == phi for iso in isos)
1
isos = psi.codomain().isomorphisms(phi.codomain())
sum(iso * psi == phi for iso in isos)

Just like EllipticCurveIsogeny, the constructor supports a model keyword argument:

sage: E = EllipticCurve(GF(6666679), [1,1])
sage: K = E(9091, 517864)
sage: phi = EllipticCurveHom_velusqrt(E, K, model='montgomery')
sage: phi
Elliptic-curve isogeny (using square-root Vélu) of degree 2999:
  From: Elliptic Curve defined by y^2 = x^3 + x + 1 over Finite Field of size 6666679
  To:   Elliptic Curve defined by y^2 = x^3 + 1559358*x^2 + x over Finite Field of size 6666679
>>> from sage.all import *
>>> E = EllipticCurve(GF(Integer(6666679)), [Integer(1),Integer(1)])
>>> K = E(Integer(9091), Integer(517864))
>>> phi = EllipticCurveHom_velusqrt(E, K, model='montgomery')
>>> phi
Elliptic-curve isogeny (using square-root Vélu) of degree 2999:
  From: Elliptic Curve defined by y^2 = x^3 + x + 1 over Finite Field of size 6666679
  To:   Elliptic Curve defined by y^2 = x^3 + 1559358*x^2 + x over Finite Field of size 6666679
E = EllipticCurve(GF(6666679), [1,1])
K = E(9091, 517864)
phi = EllipticCurveHom_velusqrt(E, K, model='montgomery')
phi

Internally, EllipticCurveHom_velusqrt works on short Weierstraß curves, but it performs the conversion automatically:

sage: E = EllipticCurve(GF(101), [1,2,3,4,5])
sage: K = E(1, 2)
sage: K.order()
37
sage: EllipticCurveHom_velusqrt(E, K)
Elliptic-curve isogeny (using square-root Vélu) of degree 37:
  From: Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Finite Field of size 101
  To:   Elliptic Curve defined by y^2 = x^3 + 66*x + 86 over Finite Field of size 101
>>> from sage.all import *
>>> E = EllipticCurve(GF(Integer(101)), [Integer(1),Integer(2),Integer(3),Integer(4),Integer(5)])
>>> K = E(Integer(1), Integer(2))
>>> K.order()
37
>>> EllipticCurveHom_velusqrt(E, K)
Elliptic-curve isogeny (using square-root Vélu) of degree 37:
  From: Elliptic Curve defined by y^2 + x*y + 3*y = x^3 + 2*x^2 + 4*x + 5 over Finite Field of size 101
  To:   Elliptic Curve defined by y^2 = x^3 + 66*x + 86 over Finite Field of size 101
E = EllipticCurve(GF(101), [1,2,3,4,5])
K = E(1, 2)
K.order()
EllipticCurveHom_velusqrt(E, K)

However, this does imply not all elliptic curves are supported. Curves without a short Weierstraß model exist in characteristics 2 and 3:

sage: F.<t> = GF(3^3)
sage: E = EllipticCurve(F, [1,1,1,1,1])
sage: P = E(t^2+2, 1)
sage: P.order()
19
sage: EllipticCurveHom_velusqrt(E, P)
Traceback (most recent call last):
...
NotImplementedError: only implemented for curves having a short Weierstrass model
>>> from sage.all import *
>>> F = GF(Integer(3)**Integer(3), names=('t',)); (t,) = F._first_ngens(1)
>>> E = EllipticCurve(F, [Integer(1),Integer(1),Integer(1),Integer(1),Integer(1)])
>>> P = E(t**Integer(2)+Integer(2), Integer(1))
>>> P.order()
19
>>> EllipticCurveHom_velusqrt(E, P)
Traceback (most recent call last):
...
NotImplementedError: only implemented for curves having a short Weierstrass model
F.<t> = GF(3^3)
E = EllipticCurve(F, [1,1,1,1,1])
P = E(t^2+2, 1)
P.order()
EllipticCurveHom_velusqrt(E, P)

Furthermore, the implementation is restricted to finite fields, since this appears to be the most relevant application for the square-root Vélu algorithm:

sage: E = EllipticCurve('26b1')
sage: P = E(1,0)
sage: P.order()
7
sage: EllipticCurveHom_velusqrt(E, P)
Traceback (most recent call last):
...
NotImplementedError: only implemented for elliptic curves over finite fields
>>> from sage.all import *
>>> E = EllipticCurve('26b1')
>>> P = E(Integer(1),Integer(0))
>>> P.order()
7
>>> EllipticCurveHom_velusqrt(E, P)
Traceback (most recent call last):
...
NotImplementedError: only implemented for elliptic curves over finite fields
E = EllipticCurve('26b1')
P = E(1,0)
P.order()
EllipticCurveHom_velusqrt(E, P)

Note

Some of the methods inherited from EllipticCurveHom compute data whose size is linear in the degree; this includes kernel polynomial and rational maps. In consequence, those methods cannot possibly run in the otherwise advertised square-root complexity, as merely storing the result already takes linear time.

AUTHORS:

  • Lorenz Panny (2022)

class sage.schemes.elliptic_curves.hom_velusqrt.EllipticCurveHom_velusqrt(E, P, *, codomain=None, model=None, Q=None)[source]

Bases: EllipticCurveHom

This class implements separable odd-degree isogenies of elliptic curves over finite fields using the square-root Vélu algorithm.

The complexity is O~() base-field operations, where is the degree.

REFERENCES: [BDLS2020]

INPUT:

  • E – an elliptic curve over a finite field

  • P – a point on E of odd order 9

  • codomain – codomain elliptic curve (optional)

  • model – string (optional); input to compute_model()

  • Q – a point on E outside P, or None

EXAMPLES:

sage: from sage.schemes.elliptic_curves.hom_velusqrt import EllipticCurveHom_velusqrt
sage: F.<t> = GF(10009^3)
sage: E = EllipticCurve(F, [t,t])
sage: K = E(2154*t^2 + 5711*t + 2899, 7340*t^2 + 4653*t + 6935)
sage: phi = EllipticCurveHom_velusqrt(E, K); phi
Elliptic-curve isogeny (using square-root Vélu) of degree 601:
  From: Elliptic Curve defined by y^2 = x^3 + t*x + t over Finite Field in t of size 10009^3
  To:   Elliptic Curve defined by y^2 = x^3 + (263*t^2+3173*t+4759)*x + (3898*t^2+6111*t+9443) over Finite Field in t of size 10009^3
sage: phi(K)
(0 : 1 : 0)
sage: P = E(2, 3163*t^2 + 7293*t + 5999)
sage: phi(P)
(6085*t^2 + 855*t + 8720 : 8078*t^2 + 9889*t + 6030 : 1)
sage: Q = E(6, 5575*t^2 + 6607*t + 9991)
sage: phi(Q)
(626*t^2 + 9749*t + 1291 : 5931*t^2 + 8549*t + 3111 : 1)
sage: phi(P + Q)
(983*t^2 + 4894*t + 4072 : 5047*t^2 + 9325*t + 336 : 1)
sage: phi(P) + phi(Q)
(983*t^2 + 4894*t + 4072 : 5047*t^2 + 9325*t + 336 : 1)
>>> from sage.all import *
>>> from sage.schemes.elliptic_curves.hom_velusqrt import EllipticCurveHom_velusqrt
>>> F = GF(Integer(10009)**Integer(3), names=('t',)); (t,) = F._first_ngens(1)
>>> E = EllipticCurve(F, [t,t])
>>> K = E(Integer(2154)*t**Integer(2) + Integer(5711)*t + Integer(2899), Integer(7340)*t**Integer(2) + Integer(4653)*t + Integer(6935))
>>> phi = EllipticCurveHom_velusqrt(E, K); phi
Elliptic-curve isogeny (using square-root Vélu) of degree 601:
  From: Elliptic Curve defined by y^2 = x^3 + t*x + t over Finite Field in t of size 10009^3
  To:   Elliptic Curve defined by y^2 = x^3 + (263*t^2+3173*t+4759)*x + (3898*t^2+6111*t+9443) over Finite Field in t of size 10009^3
>>> phi(K)
(0 : 1 : 0)
>>> P = E(Integer(2), Integer(3163)*t**Integer(2) + Integer(7293)*t + Integer(5999))
>>> phi(P)
(6085*t^2 + 855*t + 8720 : 8078*t^2 + 9889*t + 6030 : 1)
>>> Q = E(Integer(6), Integer(5575)*t**Integer(2) + Integer(6607)*t + Integer(9991))
>>> phi(Q)
(626*t^2 + 9749*t + 1291 : 5931*t^2 + 8549*t + 3111 : 1)
>>> phi(P + Q)
(983*t^2 + 4894*t + 4072 : 5047*t^2 + 9325*t + 336 : 1)
>>> phi(P) + phi(Q)
(983*t^2 + 4894*t + 4072 : 5047*t^2 + 9325*t + 336 : 1)
from sage.schemes.elliptic_curves.hom_velusqrt import EllipticCurveHom_velusqrt
F.<t> = GF(10009^3)
E = EllipticCurve(F, [t,t])
K = E(2154*t^2 + 5711*t + 2899, 7340*t^2 + 4653*t + 6935)
phi = EllipticCurveHom_velusqrt(E, K); phi
phi(K)
P = E(2, 3163*t^2 + 7293*t + 5999)
phi(P)
Q = E(6, 5575*t^2 + 6607*t + 9991)
phi(Q)
phi(P + Q)
phi(P) + phi(Q)
dual()[source]

Return the dual of this square-root Vélu isogeny as an EllipticCurveHom.

Note

The dual is computed by EllipticCurveIsogeny, hence it does not benefit from the square-root Vélu speedup.

EXAMPLES:

sage: E = EllipticCurve(GF(101^2), [1, 1, 1, 1, 1])
sage: K = E.cardinality() // 11 * E.gens()[0]
sage: phi = E.isogeny(K, algorithm='velusqrt'); phi
Elliptic-curve isogeny (using square-root Vélu) of degree 11:
  From: Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 + x + 1 over Finite Field in z2 of size 101^2
  To:   Elliptic Curve defined by y^2 = x^3 + 39*x + 40 over Finite Field in z2 of size 101^2
sage: phi.dual()
Isogeny of degree 11 from Elliptic Curve defined by y^2 = x^3 + 39*x + 40 over Finite Field in z2 of size 101^2 to Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 + x + 1 over Finite Field in z2 of size 101^2
sage: phi.dual() * phi == phi.domain().scalar_multiplication(11)
True
sage: phi * phi.dual() == phi.codomain().scalar_multiplication(11)
True
>>> from sage.all import *
>>> E = EllipticCurve(GF(Integer(101)**Integer(2)), [Integer(1), Integer(1), Integer(1), Integer(1), Integer(1)])
>>> K = E.cardinality() // Integer(11) * E.gens()[Integer(0)]
>>> phi = E.isogeny(K, algorithm='velusqrt'); phi
Elliptic-curve isogeny (using square-root Vélu) of degree 11:
  From: Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 + x + 1 over Finite Field in z2 of size 101^2
  To:   Elliptic Curve defined by y^2 = x^3 + 39*x + 40 over Finite Field in z2 of size 101^2
>>> phi.dual()
Isogeny of degree 11 from Elliptic Curve defined by y^2 = x^3 + 39*x + 40 over Finite Field in z2 of size 101^2 to Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 + x + 1 over Finite Field in z2 of size 101^2
>>> phi.dual() * phi == phi.domain().scalar_multiplication(Integer(11))
True
>>> phi * phi.dual() == phi.codomain().scalar_multiplication(Integer(11))
True
E = EllipticCurve(GF(101^2), [1, 1, 1, 1, 1])
K = E.cardinality() // 11 * E.gens()[0]
phi = E.isogeny(K, algorithm='velusqrt'); phi
phi.dual()
phi.dual() * phi == phi.domain().scalar_multiplication(11)
phi * phi.dual() == phi.codomain().scalar_multiplication(11)
inseparable_degree()[source]

Return the inseparable degree of this square-root Vélu isogeny.

Since EllipticCurveHom_velusqrt only implements separable isogenies, this method always returns one.

kernel_polynomial()[source]

Return the kernel polynomial of this square-root Vélu isogeny.

Note

The data returned by this method has size linear in the degree.

EXAMPLES:

sage: E = EllipticCurve(GF(65537^2,'a'), [5,5])
sage: K = E.cardinality()//31 * E.gens()[0]
sage: phi = E.isogeny(K, algorithm='velusqrt')
sage: h = phi.kernel_polynomial(); h
x^15 + 21562*x^14 + 8571*x^13 + 20029*x^12 + 1775*x^11 + 60402*x^10 + 17481*x^9 + 46543*x^8 + 46519*x^7 + 18590*x^6 + 36554*x^5 + 36499*x^4 + 48857*x^3 + 3066*x^2 + 23264*x + 53937
sage: h == E.isogeny(K).kernel_polynomial()
True
sage: h(K.x())
0
>>> from sage.all import *
>>> E = EllipticCurve(GF(Integer(65537)**Integer(2),'a'), [Integer(5),Integer(5)])
>>> K = E.cardinality()//Integer(31) * E.gens()[Integer(0)]
>>> phi = E.isogeny(K, algorithm='velusqrt')
>>> h = phi.kernel_polynomial(); h
x^15 + 21562*x^14 + 8571*x^13 + 20029*x^12 + 1775*x^11 + 60402*x^10 + 17481*x^9 + 46543*x^8 + 46519*x^7 + 18590*x^6 + 36554*x^5 + 36499*x^4 + 48857*x^3 + 3066*x^2 + 23264*x + 53937
>>> h == E.isogeny(K).kernel_polynomial()
True
>>> h(K.x())
0
E = EllipticCurve(GF(65537^2,'a'), [5,5])
K = E.cardinality()//31 * E.gens()[0]
phi = E.isogeny(K, algorithm='velusqrt')
h = phi.kernel_polynomial(); h
h == E.isogeny(K).kernel_polynomial()
h(K.x())
rational_maps()[source]

Return the pair of explicit rational maps of this square-root Vélu isogeny as fractions of bivariate polynomials in x and y.

Note

The data returned by this method has size linear in the degree.

EXAMPLES:

sage: E = EllipticCurve(GF(101^2), [1, 1, 1, 1, 1])
sage: K = (E.cardinality() // 11) * E.gens()[0]
sage: phi = E.isogeny(K, algorithm='velusqrt'); phi
Elliptic-curve isogeny (using square-root Vélu) of degree 11:
  From: Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 + x + 1 over Finite Field in z2 of size 101^2
  To:   Elliptic Curve defined by y^2 = x^3 + 39*x + 40 over Finite Field in z2 of size 101^2
sage: phi.rational_maps()
((-17*x^11 - 34*x^10 - 36*x^9 + ... - 29*x^2 - 25*x - 25)/(x^10 + 10*x^9 + 19*x^8 - ... + x^2 + 47*x + 24),
 (-3*x^16 - 6*x^15*y - 48*x^15 + ... - 49*x - 9*y + 46)/(x^15 + 15*x^14 - 35*x^13 - ... + 3*x^2 - 45*x + 47))
>>> from sage.all import *
>>> E = EllipticCurve(GF(Integer(101)**Integer(2)), [Integer(1), Integer(1), Integer(1), Integer(1), Integer(1)])
>>> K = (E.cardinality() // Integer(11)) * E.gens()[Integer(0)]
>>> phi = E.isogeny(K, algorithm='velusqrt'); phi
Elliptic-curve isogeny (using square-root Vélu) of degree 11:
  From: Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 + x + 1 over Finite Field in z2 of size 101^2
  To:   Elliptic Curve defined by y^2 = x^3 + 39*x + 40 over Finite Field in z2 of size 101^2
>>> phi.rational_maps()
((-17*x^11 - 34*x^10 - 36*x^9 + ... - 29*x^2 - 25*x - 25)/(x^10 + 10*x^9 + 19*x^8 - ... + x^2 + 47*x + 24),
 (-3*x^16 - 6*x^15*y - 48*x^15 + ... - 49*x - 9*y + 46)/(x^15 + 15*x^14 - 35*x^13 - ... + 3*x^2 - 45*x + 47))
E = EllipticCurve(GF(101^2), [1, 1, 1, 1, 1])
K = (E.cardinality() // 11) * E.gens()[0]
phi = E.isogeny(K, algorithm='velusqrt'); phi
phi.rational_maps()
scaling_factor()[source]

Return the Weierstrass scaling factor associated to this square-root Vélu isogeny.

The scaling factor is the constant u (in the base field) such that φω2=uω1, where φ:E1E2 is this isogeny and ωi are the standard Weierstrass differentials on Ei defined by dx/(2y+a1x+a3).

EXAMPLES:

sage: E = EllipticCurve(GF(101^2), [1, 1, 1, 1, 1])
sage: K = (E.cardinality() // 11) * E.gens()[0]
sage: phi = E.isogeny(K, algorithm='velusqrt', model='montgomery'); phi
Elliptic-curve isogeny (using square-root Vélu) of degree 11:
  From: Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 + x + 1 over Finite Field in z2 of size 101^2
  To:   Elliptic Curve defined by y^2 = x^3 + 61*x^2 + x over Finite Field in z2 of size 101^2
sage: phi.scaling_factor()
55
sage: phi.scaling_factor() == phi.formal()[1]
True
>>> from sage.all import *
>>> E = EllipticCurve(GF(Integer(101)**Integer(2)), [Integer(1), Integer(1), Integer(1), Integer(1), Integer(1)])
>>> K = (E.cardinality() // Integer(11)) * E.gens()[Integer(0)]
>>> phi = E.isogeny(K, algorithm='velusqrt', model='montgomery'); phi
Elliptic-curve isogeny (using square-root Vélu) of degree 11:
  From: Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 + x + 1 over Finite Field in z2 of size 101^2
  To:   Elliptic Curve defined by y^2 = x^3 + 61*x^2 + x over Finite Field in z2 of size 101^2
>>> phi.scaling_factor()
55
>>> phi.scaling_factor() == phi.formal()[Integer(1)]
True
E = EllipticCurve(GF(101^2), [1, 1, 1, 1, 1])
K = (E.cardinality() // 11) * E.gens()[0]
phi = E.isogeny(K, algorithm='velusqrt', model='montgomery'); phi
phi.scaling_factor()
phi.scaling_factor() == phi.formal()[1]
x_rational_map()[source]

Return the x-coordinate rational map of this square-root Vélu isogeny as a univariate rational function in x.

Note

The data returned by this method has size linear in the degree.

EXAMPLES:

sage: E = EllipticCurve(GF(101^2), [1, 1, 1, 1, 1])
sage: K = (E.cardinality() // 11) * E.gens()[0]
sage: phi = E.isogeny(K, algorithm='velusqrt'); phi
Elliptic-curve isogeny (using square-root Vélu) of degree 11:
  From: Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 + x + 1 over Finite Field in z2 of size 101^2
  To:   Elliptic Curve defined by y^2 = x^3 + 39*x + 40 over Finite Field in z2 of size 101^2
sage: phi.x_rational_map()
(84*x^11 + 67*x^10 + 65*x^9 + ... + 72*x^2 + 76*x + 76)/(x^10 + 10*x^9 + 19*x^8 + ... + x^2 + 47*x + 24)
sage: phi.x_rational_map() == phi.rational_maps()[0]
True
>>> from sage.all import *
>>> E = EllipticCurve(GF(Integer(101)**Integer(2)), [Integer(1), Integer(1), Integer(1), Integer(1), Integer(1)])
>>> K = (E.cardinality() // Integer(11)) * E.gens()[Integer(0)]
>>> phi = E.isogeny(K, algorithm='velusqrt'); phi
Elliptic-curve isogeny (using square-root Vélu) of degree 11:
  From: Elliptic Curve defined by y^2 + x*y + y = x^3 + x^2 + x + 1 over Finite Field in z2 of size 101^2
  To:   Elliptic Curve defined by y^2 = x^3 + 39*x + 40 over Finite Field in z2 of size 101^2
>>> phi.x_rational_map()
(84*x^11 + 67*x^10 + 65*x^9 + ... + 72*x^2 + 76*x + 76)/(x^10 + 10*x^9 + 19*x^8 + ... + x^2 + 47*x + 24)
>>> phi.x_rational_map() == phi.rational_maps()[Integer(0)]
True
E = EllipticCurve(GF(101^2), [1, 1, 1, 1, 1])
K = (E.cardinality() // 11) * E.gens()[0]
phi = E.isogeny(K, algorithm='velusqrt'); phi
phi.x_rational_map()
phi.x_rational_map() == phi.rational_maps()[0]
class sage.schemes.elliptic_curves.hom_velusqrt.FastEllipticPolynomial(E, n, P, Q=None)[source]

Bases: object

A class to represent and evaluate an elliptic polynomial, and optionally its derivative, in essentially square-root time.

The elliptic polynomials computed by this class are of the form

hS(Z)=iS(Zx(Q+[i]P))

where P is a point of odd order n5 and Q is either None, in which case it is assumed to be , or an arbitrary point which is not a multiple of P.

The index set S is chosen as follows:

  • If Q is given, then S={0,1,2,3,...,n1}.

  • If Q is omitted, then S={1,3,5,...,n2}. Note that in this case, h{1,2,3,...,n1} can be computed as hS2 since n is odd.

INPUT:

  • E – an elliptic curve in short Weierstraß form

  • n – an odd integer 5

  • P – a point on E

  • Q – a point on E, or None

ALGORITHM: [BDLS2020], Algorithm 2

Note

Currently only implemented for short Weierstraß curves.

EXAMPLES:

sage: from sage.schemes.elliptic_curves.hom_velusqrt import FastEllipticPolynomial
sage: E = EllipticCurve(GF(71), [5,5])
sage: P = E(4, 35)
sage: hP = FastEllipticPolynomial(E, P.order(), P); hP
Fast elliptic polynomial prod(Z - x(i*P) for i in range(1,n,2)) with n = 19, P = (4 : 35 : 1)
sage: hP(7)
19
sage: prod(7 - (i*P).x() for i in range(1,P.order(),2))
19
>>> from sage.all import *
>>> from sage.schemes.elliptic_curves.hom_velusqrt import FastEllipticPolynomial
>>> E = EllipticCurve(GF(Integer(71)), [Integer(5),Integer(5)])
>>> P = E(Integer(4), Integer(35))
>>> hP = FastEllipticPolynomial(E, P.order(), P); hP
Fast elliptic polynomial prod(Z - x(i*P) for i in range(1,n,2)) with n = 19, P = (4 : 35 : 1)
>>> hP(Integer(7))
19
>>> prod(Integer(7) - (i*P).x() for i in range(Integer(1),P.order(),Integer(2)))
19
from sage.schemes.elliptic_curves.hom_velusqrt import FastEllipticPolynomial
E = EllipticCurve(GF(71), [5,5])
P = E(4, 35)
hP = FastEllipticPolynomial(E, P.order(), P); hP
hP(7)
prod(7 - (i*P).x() for i in range(1,P.order(),2))

Passing Q changes the index set:

sage: Q = E(0, 17)
sage: hPQ = FastEllipticPolynomial(E, P.order(), P, Q)
sage: hPQ(7)
58
sage: prod(7 - (Q+i*P).x() for i in range(P.order()))
58
>>> from sage.all import *
>>> Q = E(Integer(0), Integer(17))
>>> hPQ = FastEllipticPolynomial(E, P.order(), P, Q)
>>> hPQ(Integer(7))
58
>>> prod(Integer(7) - (Q+i*P).x() for i in range(P.order()))
58
Q = E(0, 17)
hPQ = FastEllipticPolynomial(E, P.order(), P, Q)
hPQ(7)
prod(7 - (Q+i*P).x() for i in range(P.order()))

The call syntax has an optional keyword argument derivative, which will make the function return the pair (hS(α),hS(α)) instead of just hS(α):

sage: hP(7, derivative=True)
(19, 15)
sage: R.<Z> = E.base_field()[]
sage: HP = prod(Z - (i*P).x() for i in range(1,P.order(),2))
sage: HP
Z^9 + 16*Z^8 + 57*Z^7 + 6*Z^6 + 45*Z^5 + 31*Z^4 + 46*Z^3 + 10*Z^2 + 28*Z + 41
sage: HP(7)
19
sage: HP.derivative()(7)
15
>>> from sage.all import *
>>> hP(Integer(7), derivative=True)
(19, 15)
>>> R = E.base_field()['Z']; (Z,) = R._first_ngens(1)
>>> HP = prod(Z - (i*P).x() for i in range(Integer(1),P.order(),Integer(2)))
>>> HP
Z^9 + 16*Z^8 + 57*Z^7 + 6*Z^6 + 45*Z^5 + 31*Z^4 + 46*Z^3 + 10*Z^2 + 28*Z + 41
>>> HP(Integer(7))
19
>>> HP.derivative()(Integer(7))
15
hP(7, derivative=True)
R.<Z> = E.base_field()[]
HP = prod(Z - (i*P).x() for i in range(1,P.order(),2))
HP
HP(7)
HP.derivative()(7)

sage: hPQ(7, derivative=True)
(58, 62)
sage: R.<Z> = E.base_field()[]
sage: HPQ = prod(Z - (Q+i*P).x() for i in range(P.order()))
sage: HPQ
Z^19 + 53*Z^18 + 67*Z^17 + 39*Z^16 + 56*Z^15 + 32*Z^14 + 44*Z^13 + 6*Z^12 + 27*Z^11 + 29*Z^10 + 38*Z^9 + 48*Z^8 + 38*Z^7 + 43*Z^6 + 21*Z^5 + 25*Z^4 + 33*Z^3 + 49*Z^2 + 60*Z
sage: HPQ(7)
58
sage: HPQ.derivative()(7)
62
>>> from sage.all import *
>>> hPQ(Integer(7), derivative=True)
(58, 62)
>>> R = E.base_field()['Z']; (Z,) = R._first_ngens(1)
>>> HPQ = prod(Z - (Q+i*P).x() for i in range(P.order()))
>>> HPQ
Z^19 + 53*Z^18 + 67*Z^17 + 39*Z^16 + 56*Z^15 + 32*Z^14 + 44*Z^13 + 6*Z^12 + 27*Z^11 + 29*Z^10 + 38*Z^9 + 48*Z^8 + 38*Z^7 + 43*Z^6 + 21*Z^5 + 25*Z^4 + 33*Z^3 + 49*Z^2 + 60*Z
>>> HPQ(Integer(7))
58
>>> HPQ.derivative()(Integer(7))
62
hPQ(7, derivative=True)
R.<Z> = E.base_field()[]
HPQ = prod(Z - (Q+i*P).x() for i in range(P.order()))
HPQ
HPQ(7)
HPQ.derivative()(7)

The input can be an element of any algebra over the base ring:

sage: R.<T> = GF(71)[]
sage: S.<t> = R.quotient(T^2)
sage: hP(7 + t)
15*t + 19
>>> from sage.all import *
>>> R = GF(Integer(71))['T']; (T,) = R._first_ngens(1)
>>> S = R.quotient(T**Integer(2), names=('t',)); (t,) = S._first_ngens(1)
>>> hP(Integer(7) + t)
15*t + 19
R.<T> = GF(71)[]
S.<t> = R.quotient(T^2)
hP(7 + t)