Skip to content
This repository has been archived by the owner on Jan 30, 2023. It is now read-only.

Commit

Permalink
Merge remote-tracking branch 'origin/u/annasomoza/22173-shiodainv' in…
Browse files Browse the repository at this point in the history
…to local-merge-22173
  • Loading branch information
Anna Somoza committed Apr 11, 2019
2 parents 2863afa + a6a9331 commit e455b4d
Show file tree
Hide file tree
Showing 8 changed files with 177 additions and 5 deletions.
3 changes: 3 additions & 0 deletions src/doc/en/reference/references/index.rst
Original file line number Diff line number Diff line change
Expand Up @@ -3491,6 +3491,9 @@ REFERENCES:
Reduction*. Advances in Cryptology - EUROCRYPT '95. LNCS
Volume 921, 1995, pp 1-12.
.. [Sh1967] \T. Shioda. *On the graded ring of invariants of binary octavics* in
American J. of Math., 89(4):1022-1046, 1967.
.. [SHET2018] \O. Seker, P. Heggernes, T. Ekim, and Z. Caner Taskin.
*Generation of random chordal graphs using subtrees of a tree*,
:arxiv:`1810.13326v1`.
Expand Down
2 changes: 2 additions & 0 deletions src/sage/schemes/hyperelliptic_curves/all.py
Original file line number Diff line number Diff line change
Expand Up @@ -7,3 +7,5 @@
clebsch_invariants)
from .mestre import (Mestre_conic, HyperellipticCurve_from_invariants)
from . import monsky_washnitzer
from sage.misc.lazy_import import lazy_import
lazy_import('sage.schemes.hyperelliptic_curves.invariants', 'shioda_invariants')
6 changes: 5 additions & 1 deletion src/sage/schemes/hyperelliptic_curves/constructor.py
Original file line number Diff line number Diff line change
Expand Up @@ -5,6 +5,9 @@
- David Kohel (2006): initial version
- Sorina Ionica, Elisa Lorenzo Garcia, Anna Somoza (2017-01-11): Added
functionality for genus 3.
- Anna Somoza (2019-04): dynamic class creation
"""
Expand All @@ -24,6 +27,7 @@
from .hyperelliptic_rational_field import HyperellipticCurve_rational_field
from .hyperelliptic_padic_field import HyperellipticCurve_padic_field
from .hyperelliptic_g2 import HyperellipticCurve_g2
from .hyperelliptic_g3 import HyperellipticCurve_g3

from sage.rings.padics.all import is_pAdicField
from sage.rings.rational_field import is_RationalField
Expand Down Expand Up @@ -228,7 +232,7 @@ def HyperellipticCurve(f, h=0, names=None, PP=None, check_squarefree=True):
supercls = []
cls_name = []

genus_class = {2:HyperellipticCurve_g2}
genus_class = {2:HyperellipticCurve_g2, 3:HyperellipticCurve_g3}
if g in genus_class.keys():
supercls.append(genus_class[g])
cls_name.append("g"+str(g))
Expand Down
Original file line number Diff line number Diff line change
@@ -0,0 +1,20 @@
"""
Hyperelliptic curves of genus 3 over a finite field
"""
from __future__ import absolute_import

#*****************************************************************************
# Copyright (C) 2017 Sorina Ionica <sorina.ionica@gmail.com>
# 2017 Elisa Lorenzo Garcia <elisa.lorenzo@gmail.com>
# 2017 Anna Somoza <anna.somoza@upc.edu>
# Distributed under the terms of the GNU General Public License (GPL)
# http://www.gnu.org/licenses/
#*****************************************************************************

from . import hyperelliptic_finite_field, hyperelliptic_g3_generic

class HyperellipticCurve_g3_finite_field(
hyperelliptic_g3_generic.HyperellipticCurve_g3_generic,
hyperelliptic_finite_field.HyperellipticCurve_finite_field):
pass

54 changes: 54 additions & 0 deletions src/sage/schemes/hyperelliptic_curves/hyperelliptic_g3_generic.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,54 @@
r"""
Hyperelliptic curves of genus 3 over a general ring
AUTHORS:
- Sorina Ionica, Elisa Lorenzo Garcia, Anna Somoza (2017-01-11): Initial version
"""
from __future__ import absolute_import
#*****************************************************************************
# Copyright (C) 2017 Sorina Ionica <sorina.ionica@gmail.com>
# 2017 Elisa Lorenzo Garcia <elisa.lorenzo@gmail.com>
# 2017 Anna Somoza <anna.somoza@upc.edu>
# Distributed under the terms of the GNU General Public License (GPL)
# http://www.gnu.org/licenses/
#*****************************************************************************

from . import hyperelliptic_generic
from . import invariants


class HyperellipticCurve_g3_generic(hyperelliptic_generic.HyperellipticCurve_generic):
def shioda_invariants(self):
r"""
Return the Shioda invariants `(J2, J3, J4, J5, J6, J7, J8, J9, J10)` of Shioda, [Sh1967]_.
.. SEEALSO::
:meth:`sage.schemes.hyperelliptic_curves.invariants`
EXAMPLES::
sage: R.<x> = QQ[]
sage: f = x^7 - x^4 + 3
sage: C = HyperellipticCurve(f)
sage: C.shioda_invariants()
[8/35, -144/8575, 32/7203, -128/84035, 128/1058841, 10673116768/12353145,
149423644992/1008840175, 85384936192/1815912315, 398463049216/49433168575]
sage: C = HyperellipticCurve(x^8 + 3*x^2 + 2*x +1)
sage: C.shioda_invariants()
[32, 648/49, 512/3, 6912/49, -8192/9, -24064/49, -89293312/36015,
69632/21, 1441134592/108045]
sage: C = HyperellipticCurve(x*(x^6 + 1))
sage: C.shioda_invariants()
[-4, 0, 1/6, 0, 1/36, 0, 1/70, 0, 1/420]
sage: C = HyperellipticCurve(x**8+4*x**4+5, x)
sage: C.shioda_invariants()
[5728/35, 1901607/17150, 29738480/7203, 95245904/16807, -114656031488/1058841,
-1832866169368/12353145, -4626584012971696/3026520525, 1415110141852160/363182463,
713601762721839616/17795940687]
"""
f, h = self.hyperelliptic_polynomials()
return invariants.shioda_invariants(4*f + h**2)
Original file line number Diff line number Diff line change
@@ -0,0 +1,19 @@
"""
Hyperelliptic curves of genus 3 over a p-adic field
"""
from __future__ import absolute_import

#*****************************************************************************
# Copyright (C) 2017 Sorina Ionica <sorina.ionica@gmail.com>
# 2017 Elisa Lorenzo Garcia <elisa.lorenzo@gmail.com>
# 2017 Anna Somoza <anna.somoza@upc.edu>
# Distributed under the terms of the GNU General Public License (GPL)
# http://www.gnu.org/licenses/
#*****************************************************************************

from . import hyperelliptic_g3_generic, hyperelliptic_padic_field

class HyperellipticCurve_g3_padic_field(
hyperelliptic_g3_generic.HyperellipticCurve_g3_generic,
hyperelliptic_padic_field.HyperellipticCurve_padic_field):
pass
Original file line number Diff line number Diff line change
@@ -0,0 +1,19 @@
"""
Hyperelliptic curves of genus 3 over the rationals
"""
from __future__ import absolute_import

#*****************************************************************************
# Copyright (C) 2017 Sorina Ionica <sorina.ionica@gmail.com>
# 2017 Elisa Lorenzo Garcia <elisa.lorenzo@gmail.com>
# 2017 Anna Somoza <anna.somoza@upc.edu>
# Distributed under the terms of the GNU General Public License (GPL)
# http://www.gnu.org/licenses/
#*****************************************************************************

from . import hyperelliptic_g3_generic, hyperelliptic_rational_field

class HyperellipticCurve_g3_rational_field(
hyperelliptic_g3_generic.HyperellipticCurve_g3_generic,
hyperelliptic_rational_field.HyperellipticCurve_rational_field):
pass
59 changes: 55 additions & 4 deletions src/sage/schemes/hyperelliptic_curves/invariants.py
Original file line number Diff line number Diff line change
@@ -1,17 +1,18 @@
# -*- coding: utf-8 -*-
r"""
Compute invariants of quintics and sextics via 'Ueberschiebung'
Compute invariants of genus 2 and 3 hyperelliptic curves via 'Ueberschiebung'.
.. TODO::
* Implement invariants in small positive characteristic.
* Cardona-Quer and additional invariants for classifying automorphism groups.
AUTHOR:
- Nick Alexander
AUTHORS:
- Nick Alexander: Initial version
- Sorina Ionica, Elisa Lorenzo Garcia, Anna Somoza (2017-01-11): Added
functionality for genus 3.
"""
from sage.rings.all import ZZ
from sage.rings.all import PolynomialRing
Expand Down Expand Up @@ -369,6 +370,12 @@ def absolute_igusa_invariants_kohel(f):
`f` may be homogeneous in two variables or inhomogeneous in one.
REFERENCES:
.. [K] Kohel, David. ECHIDNA: Databases for Elliptic Curves
and Higher Dimensional Analogues.
Available at http://echidna.maths.usyd.edu.au/~kohel/dbs/
EXAMPLES::
sage: R.<x> = QQ[]
Expand Down Expand Up @@ -397,3 +404,47 @@ def absolute_igusa_invariants_kohel(f):
i2 = I2**3*I4/I10
i3 = I2**2*I6/I10
return (i1, i2, i3)

def shioda_invariants(f):
r"""
Given a octavic form `f`, return the 9 Shioda invariants, [Sh1967]_.
`f` may be homogeneous in two variables or inhomogeneous in one.
EXAMPLES::
sage: R.<x> = QQ[]
sage: f = x^7 - x^4 + 3
sage: from sage.schemes.hyperelliptic_curves.invariants import shioda_invariants
sage: shioda_invariants(f)
[1/70, -9/34300, 1/57624, -1/672280, 1/33882912, 333534899/6324810240,
2334744453/1033052339200, 333534907/1859494210560, 778248143/101239129241600]
sage: R.<x,y> = QQ[]
sage: shioda_invariants(x^8 + 3*x^2*y^6 + 2*x*y^7 + y^8)
[2, 81/392, 2/3, 27/196, -2/9, -47/1568, -174401/4609920, 17/1344,
703679/55319040]
"""
F = f
if f.parent().ngens() == 1:
f = PolynomialRing(f.parent().base_ring(), 1, f.parent().variable_name())(f)
x1, x2 = f.homogenize().parent().gens()
F = sum([ f[i]*x1**i*x2**(8-i) for i in range(9)])
H = Ueberschiebung(F,F,2)
g = Ueberschiebung(F,F,4)
k = Ueberschiebung(F,F,6)
h = Ueberschiebung(k,k,2)
m = Ueberschiebung(F,k,4)
n = Ueberschiebung(F,h,4)
p = Ueberschiebung(g,k,4)
q = Ueberschiebung(g,h,4)
J2 = Ueberschiebung(F,F,8)
J3 = Ueberschiebung(F,g,8)
J4 = Ueberschiebung(k,k,4)
J5 = Ueberschiebung(m,k,4)
J6 = Ueberschiebung(k,h,4)
J7 = Ueberschiebung(m,h,4)
J8 = Ueberschiebung(p,h,4)
J9 = Ueberschiebung(n,h,4)
J10 = Ueberschiebung(q,h,4)
JI = [J2,J3,J4,J5,J6,J7,J8,J9,J10]
return [j.constant_coefficient() for j in JI]

0 comments on commit e455b4d

Please sign in to comment.