# SAGEMATH module: superelliptic curves and their p-group covers ## Usage The main file is *init.sage*. In order to use it, navigate to the folder with the module (*cd yourpath/DeRhamComputation*) and type: ```load('init.sage')``` The main two "packages" are intended for: - superelliptic curves, - $(\mathbb Z/p)^n$-covers of superelliptic curves. See below and the file *examples.sage* for examples. Also, folders with name *tests* contain tests, which might be useful for learning. ## Superelliptic curves In order to define a superelliptic curve $C : y^4 = x^6 + 1$ over the finite field with 25 elements, use the following commands: ``` F. = GF(25, 'a') Rx. = PolynomialRing(F) f = x^6 + 1 C = superelliptic(f, 4) ``` The class $C$ has an optional argument *prec*, which gives the precision of precomputed expansions at infinity of the functions of the curve $C$. Note that curve of the form $y^m = f(x)$ has $\delta := GCD(\deg f, m)$ points at infinity and that $f(x)$ must be separable in order for $C$ to be smooth. There are three auxilliary classes: superelliptic_function (for functions defined on superelliptic curves), superelliptic_form (for forms defined on superelliptic curves) and superelliptic_cech (for cech cocycles for the de Rham cohomology on superelliptic curves). For example, we can define the function $x + 2y + 1$ on our curve $C$ like this: ``` Rxy. = PolynomialRing(F, 2) fct = superelliptic_function(C, x + 2*y + 1) ``` or simpler: ``` fct = C.x + 2*C.y + C.one ``` Similarly, in order to define the form $\omega = y \cdot dx$ we may use: ``` omega = superelliptic_form(C, y) ``` or simpler: ``` omega = C.y * C.dx ``` The cech cocycles are given as triples: $$ (\omega_0, f, \omega_{\infty}), $$ where $\omega_0$ is a form regular on $U_0$ (i.e. on the affine curve $y^m = f(x)$), $\omega_{\infty}$ is a form regular on $U_{\infty}$, the affine curve containing the points at infinity (explicitly given by $w^{\delta} = g(v^M \cdot w^b)$, $g(x) = x^{\deg f} \cdot f(1/x)$, $\delta := GCD(m, \deg f)$, $br - am = \delta$, $M := m/\delta$) and $f$ is a function regular on $U_0 \cap U_{\infty}$ such that $\omega_0 - \omega_{\infty} = df$. See e.g. [Section 2 in article of Kock and Tait](https://arxiv.org/pdf/1709.03422.pdf). In order to access the arguments omega_0, f, omega_{\infty} of a cocyle *eta* we use the arguments *eta.omega0*, *eta.f*, *eta.omega8* respectively. Thus, let us check that the cocycle condition omega_0 - omega_{\infty} = df is satisfied for an exemplary cocycle: ``` eta = C.de_rham_basis()[-1] # we pick one of the forms in the de Rham basis of C print(eta.omega0 - eta.omega8 == eta.f.diffn()) ``` The module allows to compute the basis of of holomorphic differential forms: ``` print(C.holomorphic_differentials_basis()) ``` One may also compute the coordinates of a given holomorphic differential form. On default, the coordinates are computed with respect to *C.holomorphic_differentials_basis()*. One may also give a basis as an optional argument. Note that this speeds up computation, since the basis is not calculated several times. ``` omega = (2*C.y^2 - C.y + C.one)/C.y^3 * C.dx print(omega.coordinates()) basis = C.holomorphic_differentials_basis() print(omega.coordinates(basis = basis)) ``` The method *expansion_at_infty()* allows to compute the Laurent expansion of a given function at a place at infinity. The parameter *place* is optional. It is a number from 0 to $\delta - 1$, giving a place at infinity in which the expansion should be computed. ``` print(omega.expansion_at_infty(place=0)) print(omega.expansion_at_infty(place=1)) ``` One can check valuation of form/function at given place at infinity, using *valuation()* method. ## $p$-group covers of superelliptic curves This module allows to define $p$\-group covers of superelliptic curves in characteristic $p$ that are **ramified over the points of infinity**. For now the following $p$\-groups are possible: $\mathbb Z/p^n$, $(\mathbb Z/p)^n$, $E(p^3)$ \(the Heisenberg group mod $p$ of order $p^3$\), $Q_8$ \(the Heisenberg group $8$\). We define now a $(\mathbb Z/3)^2$ cover of curve $C : y^2 = x^3 + x$, given by the equations $z_0^3 - z_0 = x^2 y$, $z_1^3 - z_1 = x^3$. ``` F = GF(3) Rx. = PolynomialRing(F) f = x^3 + x C = superelliptic(f, 2) f1 = C.x^2*C.y f2 = C.x^3 AS = elementary_cover([f1, f2], prec=1000) ``` Similarly, one defines $\mathbb Z/p^n$, $E(p^3)$ and $Q_8$-covers: ``` sage: F = GF(3) sage: Rx. = PolynomialRing(F) sage: P1 = superelliptic(x, 1) sage: AS = witt_cover([P1.x, P1.x^2], prec = 600) sage: AS; (Z/p^n)-cover of Superelliptic curve with the equation y^1 = x over Finite Field of size 3 with the equations: z0^p - z0 = (x) z1^p - z1 = -z0^7 + z0^5 + (x^2) sage: AS = heisenberg_cover([P1.x^2, P1.x, P1.x], prec = 600) sage: AS; (E(p^3))-cover of Superelliptic curve with the equation y^1 = x over Finite Field of size 3 with the equations: z0^p - z0 = (x^2) z1^p - z1 = (x) z2^p - z2 = z0*(x) - z1*(x) + (x) ``` Note that defining a cover may take quite a long time, since several parameters are computed. Again *prec* parameter is optional and is required to compute some parameters of the cover. Note that the functions f1, f2 **must be polynomials in x and y** so that AS has ramification points at infinity. Similarly, the are classes _as\_function, as\_form, as\_cech_ and one can write _AS.x, AS.dx_, etc. There are also methods _holomorphic\_differentials\_basis\(\)_, _de\_rham\_basis\(\)_, _coordinates\(\)_, _expansion\_at\_infty\(\)_, *valuation()* etc. Note that some functions \(e.g. _holomorphic\_differential\_basis_\) have optional _threshold_ parameter. Increase it in case of problems. In order to compute the group action of $(\mathbb Z/p)^n$ on a given function/form/cocycle, use *group_action()*, e.g. ``` G = AS.group() #p-group acting on curve print(G.elts()) omega = AS.holomorphic_differentials_basis()[1] print(omega.group_action([1, 0])) #group action by element [1, 0] print(omega.group_action([0, 1])) #group action by element [0, 1] ``` In order to compute the matrices of the action, use *group_action_matrices_holo* and *group_action_matrices_dR*: ``` p = 3 F = GF(3) Rx. = PolynomialRing(F) f = x^3 + x C = superelliptic(f, 2) f1 = C.x^2*C.y f2 = C.x^3 AS = elementary_cover([f1, f2], prec=1000) A, B = AS.group_action_matrices_holo() n = A.dimensions()[0] #Let us check that they commute and are of order p: print(A*B == B*A) print(A^p == identity_matrix(n)) print(B^p == identity_matrix(n)) ``` One can decompose it into indecomposable $(\mathbb Z/p)^2$-modules, using *magma_module_decomposition*: ``` print(magma_module_decomposition(A, B)) ``` Note that this won't work for large genus of AS, as it uses free Magma with limited input. You may however use it with argument *text=True*, to obtain Magma command on output. One can also look for magical elements: ``` print(AS.magical_element()) ``` ## Polydifferential forms on abelian covers For any $(\mathbb Z/p)^n$\-cover as above, one can define a polydifferential form (i.e. a section of $\Omega^{\otimes n}$) as follows: ``` F = GF(3) Rx. = PolynomialRing(F) f = x^3 + x C = superelliptic(f, 2) AS = as_cover(C, [C.x^2*C.y, C.x^3], prec=1000) omega = as_polyform(AS.x^5, 3) # the first argument is a function on AS, the second is the multiplicity of polyform print(omega) print(omega.expansion_at_infty()) ``` Using the command *holo_polydifferentials_basis* one may compute the basis of $H^0(\Omega^{\otimes n})$: ``` p = 5 F = GF(p) Rx. = PolynomialRing(F) C = superelliptic(x, 1) AS = as_cover(C, [C.x^8], prec = 200) print(AS.holo_polydifferentials_basis(2, threshold = 15)) #we increase the threshold if needed, see below ``` The class *as_symmetric_product_forms* may be used to define an element of $\textrm{Sym}^n \, H^0(\Omega_X)$. The command *as_symmetric_power_basis* returns a basis of $\textrm{Sym}^n \, H^0(\Omega_X)$. ``` p = 5 F = GF(p) Rx. = PolynomialRing(F) C = superelliptic(x, 1) AS = as_cover(C, [C.x^8], prec = 200) omega = as_symmetric_product_forms([[2, AS.dx, AS.z[0]*AS.dx], [-1, AS.x*AS.dx, AS.z[0]*AS.dx]]) #note that the argument is a list of lists (first coefficient, then forms that occur in the given simple tensor) print(omega) print(as_symmetric_power_basis(AS, 2)) ``` The method *canonical_ideal* computes the elements of $\textrm{Sym}^n \, H^0(\Omega_X)$ that are in the kernel of the multiplication map with codomain in $H^0(\Omega_X^{\otimes n})$ (i.e. the n-th homogeneous part of the canonical ideal). The method *canonical_ideal_polynomials* computes the corresponding polynomials and *group_action_canonical_ideal* the matrices of the group action on the n-th homogeneous part of the canonical ideal. ``` p = 5 F = GF(p) Rx. = PolynomialRing(F) C = superelliptic(x, 1) AS = as_cover(C, [C.x^8], prec = 200) print(AS.canonical_ideal(2, threshold = 15)) print(AS.canonical_ideal_polynomials(2, threshold = 15)) print(AS.group_action_canonical_ideal(2, threshold = 15)) ``` ## Quaternion covers Some of the above methods are also implemented for quaternion covers in characteristic $2$. Those are defined by the equations $z_0^2 + z_0 = f_0, z_1^2 + z_1 = f_1, z_2^2 + z_2 = f_2 + z_0f_0 + z_1 (f_0 + f_1)$. The arguments of *quaternion_cover* are: the covered superelliptic curve $C$ and the functions $f_0$, $f_1$, $f_2$ on $C$. ``` p = 2 F. = GF(p^2) Rx. = PolynomialRing(F) C = superelliptic(x, 1) Q = quaternion_cover(C, [C.x^3, a*C.x^3, 0*C.x], prec = 300) print(Q) print(Q.genus()) print(Q.holomorphic_differentials_basis()) print(Q.canonical_ideal_polynomials(2)) Qi, Qj = QuaternionGroup().gens() omega = Q.z[2]*Q.dx print(omega.group_action(Qi), omega.group_action(Qj)) ``` ## Common errors: 1. *Increase precision.* - Increase the *prec* argument of the curve. 1. *I haven't found all forms, only x of y* - Increase threshold when computing a basis. 1. *no 12 -th root; divide by 2* - when defining AS cover, one needs to compute roots of some numbers. This error means that a number is not in the field. You can either enlarge the base field, or divide one of the functions by given number and study the modified curve. 1. *unsupported operand parent(s) for %: 'The Infinity Ring' and 'The Infinity Ring'* - One of the power series turned out to be zero. Probably the AS cover that you've given is not connected (for example it is of the form $z_0^p - z_0 = f^p - f$).