Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
8 changes: 4 additions & 4 deletions src/sage/schemes/elliptic_curves/ell_curve_isogeny.py
Original file line number Diff line number Diff line change
Expand Up @@ -79,7 +79,7 @@
# https://www.gnu.org/licenses/
# ****************************************************************************

from copy import copy
from copy import copy, deepcopy

from sage.structure.sequence import Sequence

Expand Down Expand Up @@ -1442,7 +1442,7 @@ def __neg__(self):
sage: negphi.rational_maps() # optional - sage.rings.number_field
((x^2 + (-a)*x - 2)/(x + (-a)), (-x^2*y + (2*a)*x*y - y)/(x^2 + (-2*a)*x - 1))
"""
output = copy(self)
output = deepcopy(self)
output._set_post_isomorphism(negation_morphism(output._codomain))
return output

Expand Down Expand Up @@ -3312,13 +3312,13 @@ def _composition_impl(left, right):
NotImplemented
"""
if isinstance(left, WeierstrassIsomorphism) and isinstance(right, EllipticCurveIsogeny):
result = copy(right)
result = deepcopy(right)
result._set_post_isomorphism(left)
return result

if isinstance(left, EllipticCurveIsogeny) and isinstance(right, WeierstrassIsomorphism):
assert isinstance(left, EllipticCurveIsogeny)
result = copy(left)
result = deepcopy(left)
result._set_pre_isomorphism(right)
return result

Expand Down
69 changes: 69 additions & 0 deletions src/sage/schemes/elliptic_curves/ell_field.py
Original file line number Diff line number Diff line change
Expand Up @@ -2034,3 +2034,72 @@ def compute_model(E, name):
return E.montgomery_model()

raise NotImplementedError(f'cannot compute {name} model')

def point_of_order(E, l):
r"""
Given an elliptic curve `E` over a finite field or a number field
and an integer `\ell \geq 1`, construct a point of order `\ell` on `E`,
possibly defined over an extension of the base field of `E`.

Currently only prime values of `\ell` are supported.

EXAMPLES::

sage: from sage.schemes.elliptic_curves.ell_field import point_of_order
sage: E = EllipticCurve(GF(101), [1,2,3,4,5])
sage: P = point_of_order(E, 5); P
(50*Y^5 + 48*Y^4 + 26*Y^3 + 37*Y^2 + 48*Y + 15 : 25*Y^5 + 31*Y^4 + 79*Y^3 + 39*Y^2 + 3*Y + 20 : 1)
sage: P.base_ring()
Finite Field in Y of size 101^6
sage: P.order()
5
sage: P.curve().a_invariants()
(1, 2, 3, 4, 5)

::

sage: from sage.schemes.elliptic_curves.ell_field import point_of_order
sage: E = EllipticCurve(QQ, [7,7])
sage: P = point_of_order(E, 3); P # random
(x : -Y : 1)
sage: P.base_ring()
Number Field in Y with defining polynomial Y^2 - x^3 - 7*x - 7 over its base field
sage: P.order()
3
sage: P.curve().a_invariants()
(0, 0, 0, 7, 7)
"""
# Construct the field extension defined by the given polynomial,
# in such a way that the result is recognized by Sage as a field.
def ffext(poly):
rng = poly.parent()
fld = rng.base_ring()
if fld in FiniteFields():
# Workaround: .extension() would return a PolynomialQuotientRing
# rather than another FiniteField.
return poly.splitting_field(rng.variable_name())
return fld.extension(poly, rng.variable_name())

l = ZZ(l)
if l == 1:
return E(0)

if not l.is_prime():
raise NotImplementedError('composite orders are currently unsupported')

xpoly = E.division_polynomial(l)
if xpoly.degree() < 1: # supersingular and l == p
raise ValueError('curve does not have any points of the specified order')

mu = xpoly.factor()[0][0]
FF = ffext(mu)
xx = mu.any_root(ring=FF, assume_squarefree=True)

Y = polygen(FF, 'Y')
ypoly = E.defining_polynomial()(xx, Y, 1)
if ypoly.is_irreducible():
FF = ffext(ypoly)
xx = FF(xx)

EE = E.change_ring(FF)
return EE.lift_x(xx)
202 changes: 199 additions & 3 deletions src/sage/schemes/elliptic_curves/hom.py
Original file line number Diff line number Diff line change
Expand Up @@ -23,6 +23,8 @@
the common :class:`EllipticCurveHom` interface.

- Lorenz Panny (2022): :meth:`~EllipticCurveHom.matrix_on_subgroup`

- Lorenz Panny (2023): :meth:`~EllipticCurveHom.trace`, :meth:`~EllipticCurveHom.characteristic_polynomial`
"""
from sage.misc.cachefunc import cached_method
from sage.structure.richcmp import richcmp_not_equal, richcmp, op_EQ, op_NE
Expand All @@ -31,6 +33,7 @@

from sage.arith.misc import integer_floor

from sage.rings.integer_ring import ZZ
from sage.rings.finite_rings import finite_field_base
from sage.rings.number_field import number_field_base

Expand Down Expand Up @@ -283,6 +286,95 @@ def degree(self):
except AttributeError:
raise NotImplementedError('children must implement')

@cached_method
def trace(self):
r"""
Return the trace of this elliptic-curve morphism, which must
be an endomorphism.

ALGORITHM: :func:`compute_trace_generic`

EXAMPLES::

sage: E = EllipticCurve(QQ, [42, 42])
sage: m5 = E.scalar_multiplication(5)
sage: m5.trace()
10

::

sage: E = EllipticCurve(GF(71^2), [45, 45])
sage: P = E.lift_x(27)
sage: P.order()
71
sage: tau = E.isogeny(P, codomain=E)
sage: tau.trace()
-1

TESTS:

Make sure the cached value of the trace is not accidentally
copied on composition with automorphisms::

sage: aut = E.automorphisms()[1] # [-1]
sage: (aut * tau).trace()
1

It also works for more complicated :class:`EllipticCurveHom`
children::

sage: tau = E.isogeny(P, codomain=E, algorithm='velusqrt')
sage: tau.trace()
-1

Check that negation commutes with taking the trace::

sage: (-tau).trace()
1
"""
F = self.domain().base_field()
if F.characteristic().is_zero():
d = self.degree()
s = self.scaling_factor()
return ZZ(s + d/s)
return compute_trace_generic(self)

def characteristic_polynomial(self):
r"""
Return the characteristic polynomial of this elliptic-curve
morphism, which must be an endomorphism.

.. SEEALSO::

- :meth:`degree`
- :meth:`trace`

EXAMPLES::

sage: E = EllipticCurve(QQ, [42, 42])
sage: m5 = E.scalar_multiplication(5)
sage: m5.characteristic_polynomial()
x^2 - 10*x + 25

::

sage: E = EllipticCurve(GF(71), [42, 42])
sage: pi = E.frobenius_endomorphism()
sage: pi.characteristic_polynomial()
x^2 - 8*x + 71
sage: E.frobenius().charpoly()
x^2 - 8*x + 71

TESTS::

sage: m5.characteristic_polynomial().parent()
Univariate Polynomial Ring in x over Integer Ring
sage: pi.characteristic_polynomial().parent()
Univariate Polynomial Ring in x over Integer Ring
"""
R = ZZ['x']
return R([self.degree(), -self.trace(), 1])

def kernel_polynomial(self):
r"""
Return the kernel polynomial of this elliptic-curve morphism.
Expand Down Expand Up @@ -692,7 +784,7 @@ def __hash__(self):
sage: EllipticCurveIsogeny(E,X^3-13*X^2-58*X+503,check=False)
Isogeny of degree 7 from Elliptic Curve defined by y^2 + x*y = x^3 - x^2 - 107*x + 552 over Rational Field to Elliptic Curve defined by y^2 + x*y = x^3 - x^2 - 5252*x - 178837 over Rational Field
"""
return hash((self.domain(), self.codomain(), self.kernel_polynomial()))
return hash((self.domain(), self.codomain(), self.kernel_polynomial(), self.scaling_factor()))

def as_morphism(self):
r"""
Expand Down Expand Up @@ -970,8 +1062,6 @@ def find_post_isomorphism(phi, psi):
raise ValueError('codomains not isomorphic')

F = E.base_ring()
from sage.rings.finite_rings import finite_field_base
from sage.rings.number_field import number_field_base

if isinstance(F, finite_field_base.FiniteField):
while len(isos) > 1:
Expand Down Expand Up @@ -1006,3 +1096,109 @@ def find_post_isomorphism(phi, psi):

# found no suitable isomorphism -- either doesn't exist or a bug
raise ValueError('isogenies not equal up to post-isomorphism')


def compute_trace_generic(phi):
r"""
Compute the trace of the given elliptic-curve endomorphism.

ALGORITHM: Simple variant of Schoof's algorithm.
For enough small primes `\ell`, we find an order-`\ell` point `P`
on `E` and use a discrete-logarithm calculation to find the unique
scalar `t_\ell \in \{0,...,\ell-1\}` such that
`\varphi^2(P)+[\deg(\varphi)]P = [t_\ell]\varphi(P)`.
Then `t_\ell` equals the trace of `\varphi` modulo `\ell`, which
can therefore be recovered using the Chinese remainder theorem.

EXAMPLES:

It works over finite fields::

sage: from sage.schemes.elliptic_curves.hom import compute_trace_generic
sage: E = EllipticCurve(GF(31337), [1,1])
sage: compute_trace_generic(E.frobenius_endomorphism())
314

It works over `\QQ`::

sage: from sage.schemes.elliptic_curves.hom import compute_trace_generic
sage: E = EllipticCurve(QQ, [1,2,3,4,5])
sage: dbl = E.scalar_multiplication(2)
sage: compute_trace_generic(dbl)
4

It works over number fields (for a CM curve)::

sage: from sage.schemes.elliptic_curves.hom import compute_trace_generic
sage: x = polygen(QQ)
sage: K.<t> = NumberField(5*x^2 - 2*x + 1)
sage: E = EllipticCurve(K, [1,0])
sage: phi = E.isogeny([t,0,1], codomain=E) # phi = 2 + i
sage: compute_trace_generic(phi)
4

TESTS:

Check on random elliptic curves over finite fields that
the result for Frobenius matches
:meth:`~sage.schemes.elliptic_curves.ell_finite_field.EllipticCurve_finite_field.trace_of_frobenius`::

sage: from sage.schemes.elliptic_curves.hom import compute_trace_generic
sage: p = random_prime(10^3)
sage: e = randrange(1, ceil(log(10^5,p)))
sage: F.<t> = GF((p, e))
sage: E = choice(EllipticCurve(j=F.random_element()).twists())
sage: pi = E.frobenius_endomorphism()
sage: compute_trace_generic(pi) == E.trace_of_frobenius()
True

Check that the nonexistence of `p`-torsion for supersingular curves
does not cause trouble::

sage: from sage.schemes.elliptic_curves.hom import compute_trace_generic
sage: E = EllipticCurve(GF(5), [0,1])
sage: E.division_polynomial(5)
4
sage: m7 = E.scalar_multiplication(7)
sage: compute_trace_generic(-m7)
-14
"""
from sage.rings.finite_rings.integer_mod import Mod
from sage.groups.generic import discrete_log
from sage.sets.primes import Primes
from sage.schemes.elliptic_curves.ell_field import point_of_order

E = phi.domain()
if phi.codomain() != E:
raise ValueError('trace only makes sense for endomorphisms')

d = phi.degree()

M = 4 * d.isqrt() + 1 # |trace| <= 2 sqrt(deg)
tr = Mod(0,1)

F = E.base_field()
p = F.characteristic()
if p:
s = phi.scaling_factor()
if s:
tr = Mod(ZZ(s + d/s), p)

for l in Primes():
if tr.modulus() >= M:
break
try:
P = point_of_order(E, l)
except ValueError:
continue # supersingular and l == p

Q = phi._eval(P)
if not Q: # we learn nothing when P lies in the kernel
continue
R = phi._eval(Q)
t = discrete_log(R + d*P, Q, ord=l, operation='+')
# assert not R - t*Q + d*P

tr = tr.crt(Mod(t, l))

return tr.lift_centered()