Rational point sets on a Jacobian of a general hyperelliptic curve¶
AUTHORS:
Sabrina Kunzweiler, Gareth Ma, Giacomo Pope (2024): adapt to smooth model
- class sage.schemes.hyperelliptic_curves.jacobian_homset_generic.HyperellipticJacobianHomset(Y, X, **kwds)[source]¶
Bases:
SchemeHomset_pointsSet of rational points of the Jacobian.
- abelian_group()[source]¶
Return the group of rational points on this Jacobian as an
AdditiveAbelianGroupWrapperobject.EXAMPLES:
sage: x = polygen(GF(419^2)) sage: C = HyperellipticCurve(x^4 + x - 1) sage: C.is_split() True sage: C.jacobian().abelian_group() Additive abelian group isomorphic to Z/25193 + Z/7 embedded in Abelian group of points over Finite Field in z2 of size 419^2 on Jacobian of Hyperelliptic Curve over Finite Field in z2 of size 419^2 defined by y^2 = x^4 + x + 418
>>> from sage.all import * >>> x = polygen(GF(Integer(419)**Integer(2))) >>> C = HyperellipticCurve(x**Integer(4) + x - Integer(1)) >>> C.is_split() True >>> C.jacobian().abelian_group() Additive abelian group isomorphic to Z/25193 + Z/7 embedded in Abelian group of points over Finite Field in z2 of size 419^2 on Jacobian of Hyperelliptic Curve over Finite Field in z2 of size 419^2 defined by y^2 = x^4 + x + 418
sage: x = polygen(GF(419)) sage: C = HyperellipticCurve(x^5 + x) sage: C.is_ramified() True sage: C.jacobian().abelian_group() Additive abelian group isomorphic to Z/29346 + Z/6 embedded in Abelian group of points over Finite Field of size 419 on Jacobian of Hyperelliptic Curve over Finite Field of size 419 defined by y^2 = x^5 + x
[Python]>>> from sage.all import * >>> x = polygen(GF(Integer(419))) >>> C = HyperellipticCurve(x**Integer(5) + x) >>> C.is_ramified() True >>> C.jacobian().abelian_group() Additive abelian group isomorphic to Z/29346 + Z/6 embedded in Abelian group of points over Finite Field of size 419 on Jacobian of Hyperelliptic Curve over Finite Field of size 419 defined by y^2 = x^5 + x
sage: x = polygen(GF(419)) sage: C = HyperellipticCurve(-x^6 + x - 3) sage: C.is_inert() True sage: C.jacobian().abelian_group() Additive abelian group isomorphic to Z/174078 embedded in Abelian group of points over Finite Field of size 419 on Jacobian of Hyperelliptic Curve over Finite Field of size 419 defined by y^2 = 418*x^6 + x + 416
>>> from sage.all import * >>> x = polygen(GF(Integer(419))) >>> C = HyperellipticCurve(-x**Integer(6) + x - Integer(3)) >>> C.is_inert() True >>> C.jacobian().abelian_group() Additive abelian group isomorphic to Z/174078 embedded in Abelian group of points over Finite Field of size 419 on Jacobian of Hyperelliptic Curve over Finite Field of size 419 defined by y^2 = 418*x^6 + x + 416
sage: x = polygen(GF(419)) sage: C = HyperellipticCurve(11*x^6 + x) sage: C.is_inert() True sage: C.jacobian().abelian_group() Additive abelian group isomorphic to Z/420 + Z/420 embedded in Abelian group of points over Finite Field of size 419 on Jacobian of Hyperelliptic Curve over Finite Field of size 419 defined by y^2 = 11*x^6 + x
[Python]>>> from sage.all import * >>> x = polygen(GF(Integer(419))) >>> C = HyperellipticCurve(Integer(11)*x**Integer(6) + x) >>> C.is_inert() True >>> C.jacobian().abelian_group() Additive abelian group isomorphic to Z/420 + Z/420 embedded in Abelian group of points over Finite Field of size 419 on Jacobian of Hyperelliptic Curve over Finite Field of size 419 defined by y^2 = 11*x^6 + x
sage: x = polygen(GF(419)) sage: C = HyperellipticCurve(x^7 - x) sage: C.is_ramified() True sage: C.jacobian().abelian_group() Additive abelian group isomorphic to Z/420 + Z/420 + Z/210 + Z/2 embedded in Abelian group of points over Finite Field of size 419 on Jacobian of Hyperelliptic Curve over Finite Field of size 419 defined by y^2 = x^7 + 418*x
>>> from sage.all import * >>> x = polygen(GF(Integer(419))) >>> C = HyperellipticCurve(x**Integer(7) - x) >>> C.is_ramified() True >>> C.jacobian().abelian_group() Additive abelian group isomorphic to Z/420 + Z/420 + Z/210 + Z/2 embedded in Abelian group of points over Finite Field of size 419 on Jacobian of Hyperelliptic Curve over Finite Field of size 419 defined by y^2 = x^7 + 418*x
sage: x = polygen(GF(419)) sage: C = HyperellipticCurve(x^8 + x) sage: C.is_split() True sage: C.jacobian().abelian_group() Additive abelian group isomorphic to Z/420 + Z/420 + Z/420 embedded in Abelian group of points over Finite Field of size 419 on Jacobian of Hyperelliptic Curve over Finite Field of size 419 defined by y^2 = x^8 + x
[Python]>>> from sage.all import * >>> x = polygen(GF(Integer(419))) >>> C = HyperellipticCurve(x**Integer(8) + x) >>> C.is_split() True >>> C.jacobian().abelian_group() Additive abelian group isomorphic to Z/420 + Z/420 + Z/420 embedded in Abelian group of points over Finite Field of size 419 on Jacobian of Hyperelliptic Curve over Finite Field of size 419 defined by y^2 = x^8 + x
- cantor_composition(u1, v1, u2, v2)[source]¶
Return the Cantor composition of
(u1, v1)and(u2, v2).EXAMPLES:
sage: R.<x> = GF(13)[] sage: H = HyperellipticCurve(x^7 + x^5 + x + 1) sage: JF = Jacobian(H).point_homset() sage: (u1, v1) = (x^3 + 4*x^2, 10*x^2 + 7*x + 1) sage: (u2, v2) = (x^3 + 8*x^2 + 11*x + 2, x^2 + 9*x + 10) sage: JF.cantor_composition(u1, v1, u2, v2) (x^6 + 12*x^5 + 4*x^4 + 7*x^3 + 8*x^2, 5*x^5 + 2*x^4 + 12*x^2 + 7*x + 1)
>>> from sage.all import * >>> R = GF(Integer(13))['x']; (x,) = R._first_ngens(1) >>> H = HyperellipticCurve(x**Integer(7) + x**Integer(5) + x + Integer(1)) >>> JF = Jacobian(H).point_homset() >>> (u1, v1) = (x**Integer(3) + Integer(4)*x**Integer(2), Integer(10)*x**Integer(2) + Integer(7)*x + Integer(1)) >>> (u2, v2) = (x**Integer(3) + Integer(8)*x**Integer(2) + Integer(11)*x + Integer(2), x**Integer(2) + Integer(9)*x + Integer(10)) >>> JF.cantor_composition(u1, v1, u2, v2) (x^6 + 12*x^5 + 4*x^4 + 7*x^3 + 8*x^2, 5*x^5 + 2*x^4 + 12*x^2 + 7*x + 1)
- cantor_reduction(u0, v0)[source]¶
Apply one reduction step of Cantor’s algorithm to
(u0, v0).Note that, in general, several steps are necessary the representation of a reduced divisor.
EXAMPLES:
sage: R.<x> = GF(13)[] sage: H = HyperellipticCurve(x^7 + x^5 + x + 1) sage: g = H.genus() sage: JF = Jacobian(H).point_homset() sage: (u0, v0) = (x^6 + 12*x^5 + 4*x^4 + 7*x^3 + 8*x^2, 5*x^5 + 2*x^4 + 12*x^2 + 7*x + 1) sage: (u1, v1) = JF.cantor_reduction(u0, v0) sage: u1.degree() <= g False sage: (u2, v2) = JF.cantor_reduction(u1, v1) sage: u2.degree() <= g True
>>> from sage.all import * >>> R = GF(Integer(13))['x']; (x,) = R._first_ngens(1) >>> H = HyperellipticCurve(x**Integer(7) + x**Integer(5) + x + Integer(1)) >>> g = H.genus() >>> JF = Jacobian(H).point_homset() >>> (u0, v0) = (x**Integer(6) + Integer(12)*x**Integer(5) + Integer(4)*x**Integer(4) + Integer(7)*x**Integer(3) + Integer(8)*x**Integer(2), Integer(5)*x**Integer(5) + Integer(2)*x**Integer(4) + Integer(12)*x**Integer(2) + Integer(7)*x + Integer(1)) >>> (u1, v1) = JF.cantor_reduction(u0, v0) >>> u1.degree() <= g False >>> (u2, v2) = JF.cantor_reduction(u1, v1) >>> u2.degree() <= g True
Applying the reduction step to a reduced divisor might have unintended output, as is illustrated below.
sage: (u3, v3) = JF.cantor_reduction(u2, v2) sage: u3.degree() >= g True
- cardinality(extension_degree=1)[source]¶
Return \(|Jac(C) / \mathbb{F}_{q^n}|\).
EXAMPLES:
sage: R.<x> = GF(5)[] sage: H = HyperellipticCurve(x^6 + x + 1) sage: J = H.jacobian() sage: J(GF(5)).cardinality() 31 sage: J(GF(5^2)).cardinality() 961
>>> from sage.all import * >>> R = GF(Integer(5))['x']; (x,) = R._first_ngens(1) >>> H = HyperellipticCurve(x**Integer(6) + x + Integer(1)) >>> J = H.jacobian() >>> J(GF(Integer(5))).cardinality() 31 >>> J(GF(Integer(5)**Integer(2))).cardinality() 961
- count_points(n=1)[source]¶
Count the number of points of the Jacobian over all finite extensions of the base fields of degree less than or equal to n.
INPUT:
n– a positive integer
EXAMPLES:
sage: R.<x> = GF(5)[] sage: H = HyperellipticCurve(x^6 + x + 1) sage: J = H.jacobian() sage: J.count_points(10) == [J.change_ring(GF((5, k))).order() for k in range(1, 11)] True sage: J2 = J(GF((5, 2))) sage: J2.count_points(5) == J.count_points(10)[1::2] True
>>> from sage.all import * >>> R = GF(Integer(5))['x']; (x,) = R._first_ngens(1) >>> H = HyperellipticCurve(x**Integer(6) + x + Integer(1)) >>> J = H.jacobian() >>> J.count_points(Integer(10)) == [J.change_ring(GF((Integer(5), k))).order() for k in range(Integer(1), Integer(11))] True >>> J2 = J(GF((Integer(5), Integer(2)))) >>> J2.count_points(Integer(5)) == J.count_points(Integer(10))[Integer(1)::Integer(2)] True
- curve()[source]¶
On input the set of \(L\)-rational points of a Jacobian \(Jac(H)\) defined over \(K\), return the curve \(H\).
NOTE: The base field of \(H\) is not extended to \(L\).
EXAMPLES:
sage: R.<x> = QQ[] sage: K.<omega> = QQ.extension(x^2+x+1) sage: H = HyperellipticCurve(x^6-1) sage: JK = Jacobian(H)(K); JK Abelian group of points over Number Field in omega with defining polynomial x^2 + x + 1 on Jacobian of Hyperelliptic Curve over Rational Field defined by y^2 = x^6 - 1 sage: JK.curve() Hyperelliptic Curve over Rational Field defined by y^2 = x^6 - 1
>>> from sage.all import * >>> R = QQ['x']; (x,) = R._first_ngens(1) >>> K = QQ.extension(x**Integer(2)+x+Integer(1), names=('omega',)); (omega,) = K._first_ngens(1) >>> H = HyperellipticCurve(x**Integer(6)-Integer(1)) >>> JK = Jacobian(H)(K); JK Abelian group of points over Number Field in omega with defining polynomial x^2 + x + 1 on Jacobian of Hyperelliptic Curve over Rational Field defined by y^2 = x^6 - 1 >>> JK.curve() Hyperelliptic Curve over Rational Field defined by y^2 = x^6 - 1
- extended_curve()[source]¶
On input the set of \(L\)-rational points of a Jacobian \(Jac(H)\) defined over \(K\), return the curve \(H\) with base extended to \(L\).
EXAMPLES:
sage: R.<x> = QQ[] sage: K.<omega> = QQ.extension(x^2+x+1) sage: H = HyperellipticCurve(x^6-1) sage: JK = Jacobian(H)(K); JK Abelian group of points over Number Field in omega with defining polynomial x^2 + x + 1 on Jacobian of Hyperelliptic Curve over Rational Field defined by y^2 = x^6 - 1 sage: JK.extended_curve() Hyperelliptic Curve over Number Field in omega with defining polynomial x^2 + x + 1 defined by y^2 = x^6 - 1
>>> from sage.all import * >>> R = QQ['x']; (x,) = R._first_ngens(1) >>> K = QQ.extension(x**Integer(2)+x+Integer(1), names=('omega',)); (omega,) = K._first_ngens(1) >>> H = HyperellipticCurve(x**Integer(6)-Integer(1)) >>> JK = Jacobian(H)(K); JK Abelian group of points over Number Field in omega with defining polynomial x^2 + x + 1 on Jacobian of Hyperelliptic Curve over Rational Field defined by y^2 = x^6 - 1 >>> JK.extended_curve() Hyperelliptic Curve over Number Field in omega with defining polynomial x^2 + x + 1 defined by y^2 = x^6 - 1
- lift_u(u, all=False)[source]¶
Return one or all points with given \(u\)-coordinate.
This method is deterministic: it returns the same data each time when called with the same \(u\).
Currently only implemented for Jacobians over a finite field.
INPUT:
u– an element of the base ring of the Jacobianall– boolean (default:False); ifTrue, return a (possibly empty) list of all points with the given \(u\)-coordinate; ifFalse, return just one point, or raise aValueErrorif there are none.
OUTPUT:
A point on this Jacobian.
EXAMPLES:
sage: R.<x> = GF(1993)[] sage: H = HyperellipticCurve(x^5 + x + 1) sage: J = H.jacobian() sage: P = J.lift_u(x^2 + 42*x + 270); P (x^2 + 42*x + 270, 1837*x + 838)
>>> from sage.all import * >>> R = GF(Integer(1993))['x']; (x,) = R._first_ngens(1) >>> H = HyperellipticCurve(x**Integer(5) + x + Integer(1)) >>> J = H.jacobian() >>> P = J.lift_u(x**Integer(2) + Integer(42)*x + Integer(270)); P (x^2 + 42*x + 270, 1837*x + 838)
- order()[source]¶
Compute the order of the Jacobian.
EXAMPLES:
We compute the order of a superspecial hyperelliptic curve of genus 3:
sage: R.<x> = GF(7)[] sage: H = HyperellipticCurve(x^8 - 1) sage: J = H.jacobian() sage: J(GF(7)).order() == (7+1)^3 True sage: J(GF(7^2)).order() == (7+1)^6 True sage: R.<x> = QQ[] sage: H = HyperellipticCurve(x^8 - 1) sage: J = H.jacobian() sage: J.order() Traceback (most recent call last): ... NotImplementedError: order computation is only implemented for Jacobians over finite fields
>>> from sage.all import * >>> R = GF(Integer(7))['x']; (x,) = R._first_ngens(1) >>> H = HyperellipticCurve(x**Integer(8) - Integer(1)) >>> J = H.jacobian() >>> J(GF(Integer(7))).order() == (Integer(7)+Integer(1))**Integer(3) True >>> J(GF(Integer(7)**Integer(2))).order() == (Integer(7)+Integer(1))**Integer(6) True >>> R = QQ['x']; (x,) = R._first_ngens(1) >>> H = HyperellipticCurve(x**Integer(8) - Integer(1)) >>> J = H.jacobian() >>> J.order() Traceback (most recent call last): ... NotImplementedError: order computation is only implemented for Jacobians over finite fields
- point_to_mumford_coordinates(P)[source]¶
On input a point P, return the Mumford coordinates of (the affine part of) the divisor [P].
EXAMPLES:
sage: R.<x> = QQ[] sage: H = HyperellipticCurve(x^5 - 2*x^4 + 2*x^3 - x^2, 1) sage: P = H([2,3]); P (2 : 3 : 1) sage: JQ = H.jacobian()(QQ) sage: JQ.point_to_mumford_coordinates(P) (x - 2, 3)
>>> from sage.all import * >>> R = QQ['x']; (x,) = R._first_ngens(1) >>> H = HyperellipticCurve(x**Integer(5) - Integer(2)*x**Integer(4) + Integer(2)*x**Integer(3) - x**Integer(2), Integer(1)) >>> P = H([Integer(2),Integer(3)]); P (2 : 3 : 1) >>> JQ = H.jacobian()(QQ) >>> JQ.point_to_mumford_coordinates(P) (x - 2, 3)
- points()[source]¶
Return all points on this Jacobian \(J(K)\).
Warning
This code is not efficient at all.
EXAMPLES:
sage: R.<x> = GF(3)[] sage: H = HyperellipticCurve(x^7 + 2*x + 1) sage: J3 = H.jacobian()(GF(3)) sage: Pts = J3.points() sage: len(Pts) 94 sage: Pts[10] (x^2 + 2, x) sage: Pts[-1] (x^3 + 2*x^2 + 2*x + 1, 1)
>>> from sage.all import * >>> R = GF(Integer(3))['x']; (x,) = R._first_ngens(1) >>> H = HyperellipticCurve(x**Integer(7) + Integer(2)*x + Integer(1)) >>> J3 = H.jacobian()(GF(Integer(3))) >>> Pts = J3.points() >>> len(Pts) 94 >>> Pts[Integer(10)] (x^2 + 2, x) >>> Pts[-Integer(1)] (x^3 + 2*x^2 + 2*x + 1, 1)
- random_element(fast=False, *args, **kwargs)[source]¶
Returns a random element from the Jacobian. Distribution is NOT uniformly random.
INPUT:
fast– (boolean, defaultTrue) If set toTrue, a fast algorithm is used, but the output is NOT guaranteed to cover the entire Jacobian. See examples below. If set toFalse, a slower algorithm is used, but covers the entire Jacobian.
EXAMPLES:
sage: R.<x> = GF(5)[] sage: f = x^5 + 2*x^4 + 4*x^3 + x^2 + 4*x + 3 sage: H = HyperellipticCurve(f) sage: J = H.jacobian()
>>> from sage.all import * >>> R = GF(Integer(5))['x']; (x,) = R._first_ngens(1) >>> f = x**Integer(5) + Integer(2)*x**Integer(4) + Integer(4)*x**Integer(3) + x**Integer(2) + Integer(4)*x + Integer(3) >>> H = HyperellipticCurve(f) >>> J = H.jacobian()
This example demonstrates that the
fastalgorithm is not necessarily surjective on the rational points of the Jacobian:sage: JH = J.point_homset() sage: len(set(JH.random_element(fast=True) for _ in range(300))) 8 sage: len(set(JH.random_element(fast=False) for _ in range(300))) 16 sage: J.order() 16
[Python]>>> from sage.all import * >>> JH = J.point_homset() >>> len(set(JH.random_element(fast=True) for _ in range(Integer(300)))) 8 >>> len(set(JH.random_element(fast=False) for _ in range(Integer(300)))) 16 >>> J.order() 16
- zero(check=True)[source]¶
Return the zero element of this jacobian homset.
EXAMPLES:
sage: R.<x> = QQ[] sage: H = HyperellipticCurve(x^5 + 1) sage: JQ = H.jacobian()(QQ) sage: JQ.zero() (1, 0)
>>> from sage.all import * >>> R = QQ['x']; (x,) = R._first_ngens(1) >>> H = HyperellipticCurve(x**Integer(5) + Integer(1)) >>> JQ = H.jacobian()(QQ) >>> JQ.zero() (1, 0)