{ "cells": [ { "cell_type": "code", "execution_count": 1, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ "class superelliptic:\n", " \"\"\"Class of a superelliptic curve. Given a polynomial f(x) with coefficient field F, it constructs\n", " the curve y^m = f(x)\"\"\"\n", " def __init__(self, f, m):\n", " Rx = f.parent()\n", " x = Rx.gen()\n", " F = Rx.base() \n", " Rx. = PolynomialRing(F)\n", " Rxy. = PolynomialRing(F, 2)\n", " Fxy = FractionField(Rxy)\n", " self.polynomial = Rx(f)\n", " self.exponent = m\n", " self.base_ring = F\n", " self.characteristic = F.characteristic()\n", " r = Rx(f).degree()\n", " delta = GCD(r, m)\n", " \n", " def __repr__(self):\n", " f = self.polynomial\n", " m = self.exponent\n", " F = self.base_ring\n", " return 'Superelliptic curve with the equation y^' + str(m) + ' = ' + str(f)+' over ' + str(F)\n", "\n", " #Auxilliary algorithm that returns the basis of holomorphic differentials\n", " #of the curve and (as a second argument) the list of pairs (i, j)\n", " #such that x^i dx/y^j is holomorphic.\n", " def basis_holomorphic_differentials_degree(self):\n", " f = self.polynomial\n", " m = self.exponent\n", " r = f.degree()\n", " delta = GCD(r, m)\n", " F = self.base_ring\n", " Rx. = PolynomialRing(F)\n", " Rxy. = PolynomialRing(F, 2)\n", " Fxy = FractionField(Rxy)\n", " #########basis of holomorphic differentials and de Rham\n", "\n", " basis_holo = []\n", " degrees0 = {}\n", " k = 0\n", "\n", " for j in range(1, m):\n", " for i in range(1, r):\n", " if (r*j - m*i >= delta):\n", " basis_holo += [superelliptic_form(self, Fxy(x^(i-1)/y^j))]\n", " degrees0[k] = (i-1, j)\n", " k = k+1\n", "\n", " return(basis_holo, degrees0)\n", " \n", " #Returns the basis of holomorphic differentials using the previous algorithm.\n", " def holomorphic_differentials_basis(self):\n", " basis_holo, degrees0 = self.basis_holomorphic_differentials_degree()\n", " return basis_holo\n", " #Returns the list of pairs (i, j) such that x^i dx/y^j is holomorphic.\n", " def degrees_holomorphic_differentials(self):\n", " basis_holo, degrees0 = self.basis_holomorphic_differentials_degree()\n", " return degrees0\n", "\n", " def basis_de_rham_degrees(self):\n", " f = self.polynomial\n", " m = self.exponent\n", " r = f.degree()\n", " delta = GCD(r, m)\n", " F = self.base_ring\n", " Rx. = PolynomialRing(F)\n", " Rxy. = PolynomialRing(F, 2)\n", " Fxy = FractionField(Rxy)\n", " basis_holo = self.holomorphic_differentials_basis()\n", " basis = []\n", " #First g_X elements of basis are holomorphic differentials.\n", " for k in range(0, len(basis_holo)):\n", " basis += [superelliptic_cech(self, basis_holo[k], superelliptic_function(self, 0))]\n", "\n", " ## Next elements do not come from holomorphic differentials.\n", " t = len(basis)\n", " degrees0 = {}\n", " degrees1 = {}\n", " for j in range(1, m):\n", " for i in range(1, r):\n", " if (r*(m-j) - m*i >= delta): \n", " s = Rx(m-j)*Rx(x)*Rx(f.derivative()) - Rx(m)*Rx(i)*f\n", " psi = Rx(cut(s, i))\n", " basis += [superelliptic_cech(self, superelliptic_form(self, Fxy(psi/y^j)), superelliptic_function(self, Fxy(m*y^(m-j)/x^i)))]\n", " degrees0[t] = (psi.degree(), j)\n", " degrees1[t] = (-i, m-j)\n", " t += 1\n", " return basis, degrees0, degrees1\n", "\n", " def de_rham_basis(self):\n", " basis, degrees0, degrees1 = self.basis_de_rham_degrees()\n", " return basis\n", "\n", " def degrees_de_rham0(self):\n", " basis, degrees0, degrees1 = self.basis_de_rham_degrees()\n", " return degrees0\n", "\n", " def degrees_de_rham1(self):\n", " basis, degrees0, degrees1 = self.basis_de_rham_degrees()\n", " return degrees1 \n", " \n", " def is_smooth(self):\n", " f = self.polynomial\n", " if f.discriminant() == 0:\n", " return 0\n", " return 1\n", " \n", " def genus(self):\n", " r = self.polynomial.degree()\n", " m = self.exponent\n", " delta = GCD(r, m)\n", " return 1/2*((r-1)*(m-1) - delta + 1)\n", " \n", " def verschiebung_matrix(self):\n", " basis = self.de_rham_basis()\n", " g = self.genus()\n", " p = self.characteristic\n", " F = self.base_ring\n", " M = matrix(F, 2*g, 2*g)\n", " for i in range(0, len(basis)):\n", " w = basis[i]\n", " v = w.verschiebung().coordinates()\n", " M[i, :] = v\n", " return M\n", " \n", " def frobenius_matrix(self):\n", " basis = self.de_rham_basis()\n", " g = self.genus()\n", " p = self.characteristic\n", " F = self.base_ring\n", " M = matrix(F, 2*g, 2*g)\n", " \n", " for i in range(0, len(basis)):\n", " w = basis[i]\n", " v = w.frobenius().coordinates()\n", " M[i, :] = v\n", " return M\n", "\n", " def cartier_matrix(self):\n", " basis = self.holomorphic_differentials_basis()\n", " g = self.genus()\n", " p = self.characteristic\n", " F = self.base_ring\n", " M = matrix(F, g, g)\n", " for i in range(0, len(basis)):\n", " w = basis[i]\n", " v = w.cartier().coordinates()\n", " M[i, :] = v\n", " return M \n", "\n", "# def p_rank(self):\n", "# return self.cartier_matrix().rank()\n", " \n", " def a_number(self):\n", " g = C.genus()\n", " return g - self.cartier_matrix().rank()\n", " \n", " def final_type(self, test = 0):\n", " Fr = self.frobenius_matrix()\n", " V = self.verschiebung_matrix()\n", " p = self.characteristic\n", " return flag(Fr, V, p, test)\n", "\n", "#Auxilliary. Given a superelliptic curve C : y^m = f(x) and a polynomial g(x, y)\n", "#it replaces repeteadly all y^m's in g(x, y) by f(x). As a result\n", "#you obtain \\sum_{i = 0}^{m-1} y^i g_i(x).\n", "def reduction(C, g):\n", " p = C.characteristic\n", " F = C.base_ring\n", " Rxy. = PolynomialRing(F, 2)\n", " Fxy = FractionField(Rxy)\n", " f = C.polynomial\n", " r = f.degree()\n", " m = C.exponent\n", " g = Fxy(g)\n", " g1 = g.numerator()\n", " g2 = g.denominator()\n", " \n", " Rx. = PolynomialRing(F)\n", " Fx = FractionField(Rx)\n", " FxRy. = PolynomialRing(Fx) \n", " (A, B, C) = xgcd(FxRy(g2), FxRy(y^m - f))\n", " g = FxRy(g1*B/A)\n", " \n", " while(g.degree(Rxy(y)) >= m):\n", " d = g.degree(Rxy(y))\n", " G = coff(g, d)\n", " i = floor(d/m)\n", " g = g - G*y^d + f^i * y^(d%m) *G\n", " \n", " return(FxRy(g))\n", "\n", "#Auxilliary. Given a superelliptic curve C : y^m = f(x) and a polynomial g(x, y)\n", "#it replaces repeteadly all y^m's in g(x, y) by f(x). As a result\n", "#you obtain \\sum_{i = 0}^{m-1} g_i(x)/y^i. This is needed for reduction of\n", "#superelliptic forms.\n", "def reduction_form(C, g):\n", " F = C.base_ring\n", " Rxy. = PolynomialRing(F, 2)\n", " Fxy = FractionField(Rxy)\n", " f = C.polynomial\n", " r = f.degree()\n", " m = C.exponent\n", " g = reduction(C, g)\n", "\n", " g1 = Rxy(0)\n", " Rx. = PolynomialRing(F)\n", " Fx = FractionField(Rx)\n", " FxRy. = PolynomialRing(Fx)\n", " \n", " g = FxRy(g)\n", " for j in range(0, m):\n", " if j==0:\n", " G = coff(g, 0)\n", " g1 += FxRy(G)\n", " else:\n", " G = coff(g, j)\n", " g1 += Fxy(y^(j-m)*f*G)\n", " return(g1)\n", "\n", "#Class of rational functions on a superelliptic curve C. g = g(x, y) is a polynomial\n", "#defining the function.\n", "class superelliptic_function:\n", " def __init__(self, C, g):\n", " F = C.base_ring\n", " Rxy. = PolynomialRing(F, 2)\n", " Fxy = FractionField(Rxy)\n", " f = C.polynomial\n", " r = f.degree()\n", " m = C.exponent\n", " \n", " self.curve = C\n", " g = reduction(C, g)\n", " self.function = g\n", " \n", " def __repr__(self):\n", " return str(self.function)\n", " \n", " def jth_component(self, j):\n", " g = self.function\n", " C = self.curve\n", " F = C.base_ring\n", " Rx. = PolynomialRing(F)\n", " Fx. = FractionField(Rx)\n", " FxRy. = PolynomialRing(Fx)\n", " g = FxRy(g)\n", " return coff(g, j)\n", " \n", " def __add__(self, other):\n", " C = self.curve\n", " g1 = self.function\n", " g2 = other.function\n", " g = reduction(C, g1 + g2)\n", " return superelliptic_function(C, g)\n", " \n", " def __sub__(self, other):\n", " C = self.curve\n", " g1 = self.function\n", " g2 = other.function\n", " g = reduction(C, g1 - g2)\n", " return superelliptic_function(C, g)\n", " \n", " def __mul__(self, other):\n", " C = self.curve\n", " g1 = self.function\n", " g2 = other.function\n", " g = reduction(C, g1 * g2)\n", " return superelliptic_function(C, g)\n", " \n", " def __truediv__(self, other):\n", " C = self.curve\n", " g1 = self.function\n", " g2 = other.function\n", " g = reduction(C, g1 / g2)\n", " return superelliptic_function(C, g)\n", "\n", " def __pow__(self, exp):\n", " C = self.curve\n", " g = self.function\n", " return superelliptic_function(C, g^(exp))\n", "\n", " def diffn(self):\n", " C = self.curve\n", " f = C.polynomial\n", " m = C.exponent\n", " F = C.base_ring\n", " g = self.function\n", " Rxy. = PolynomialRing(F, 2)\n", " Fxy = FractionField(Rxy)\n", " g = Fxy(g)\n", " A = g.derivative(x)\n", " B = g.derivative(y)*f.derivative(x)/(m*y^(m-1))\n", " return superelliptic_form(C, A+B)\n", "\n", " \n", " def expansion_at_infty(self, i = 0, prec=10):\n", " C = self.curve\n", " f = C.polynomial\n", " m = C.exponent\n", " F = C.base_ring\n", " Rx. = PolynomialRing(F)\n", " f = Rx(f)\n", " Rt. = LaurentSeriesRing(F, default_prec=prec)\n", " RptW. = PolynomialRing(Rt)\n", " RptWQ = FractionField(RptW)\n", " Rxy. = PolynomialRing(F)\n", " RxyQ = FractionField(Rxy)\n", " fct = self.function\n", " fct = RxyQ(fct)\n", " r = f.degree()\n", " delta, a, b = xgcd(m, r)\n", " b = -b\n", " M = m/delta\n", " R = r/delta\n", " while a<0:\n", " a += R\n", " b += M\n", " \n", " g = (x^r*f(x = 1/x))\n", " gW = RptWQ(g(x = t^M * W^b)) - W^(delta)\n", " ww = naive_hensel(gW, F, start = root_of_unity(F, delta)^i, prec = prec)\n", " xx = Rt(1/(t^M*ww^b))\n", " yy = 1/(t^R*ww^a)\n", " return Rt(fct(x = Rt(xx), y = Rt(yy)))\n", " \n", "def naive_hensel(fct, F, start = 1, prec=10):\n", " Rt. = LaurentSeriesRing(F, default_prec=prec)\n", " RtQ = FractionField(Rt)\n", " RptW. = PolynomialRing(RtQ)\n", " RptWQ = FractionField(RptW)\n", " fct = RptWQ(fct)\n", " fct = RptW(numerator(fct))\n", " #return(fct)\n", " #while fct not in RptW:\n", " # print(fct)\n", " # fct *= W\n", " alpha = (fct.derivative())(W = start)\n", " w0 = Rt(start)\n", " i = 1\n", " while(i < prec):\n", " w0 = w0 - fct(W = w0)/alpha + O(t^(prec))\n", " i += 1\n", " return w0\n", "\n", "class superelliptic_form:\n", " def __init__(self, C, g):\n", " F = C.base_ring\n", " Rxy. = PolynomialRing(F, 2)\n", " Fxy = FractionField(Rxy)\n", " g = Fxy(reduction_form(C, g))\n", " self.form = g\n", " self.curve = C\n", " \n", " def __add__(self, other):\n", " C = self.curve\n", " g1 = self.form\n", " g2 = other.form\n", " g = reduction(C, g1 + g2)\n", " return superelliptic_form(C, g)\n", " \n", " def __sub__(self, other):\n", " C = self.curve\n", " g1 = self.form\n", " g2 = other.form\n", " g = reduction(C, g1 - g2)\n", " return superelliptic_form(C, g)\n", " \n", " def __repr__(self):\n", " g = self.form\n", " if len(str(g)) == 1:\n", " return str(g) + ' dx'\n", " return '('+str(g) + ') dx'\n", "\n", " def __rmul__(self, constant):\n", " C = self.curve\n", " omega = self.form\n", " return superelliptic_form(C, constant*omega) \n", " \n", " def cartier(self):\n", " C = self.curve\n", " m = C.exponent\n", " p = C.characteristic\n", " f = C.polynomial\n", " F = C.base_ring\n", " Rx. = PolynomialRing(F)\n", " Fx = FractionField(Rx)\n", " FxRy. = PolynomialRing(Fx)\n", " Fxy = FractionField(FxRy)\n", " result = superelliptic_form(C, FxRy(0))\n", " mult_order = Integers(m)(p).multiplicative_order()\n", " M = Integer((p^(mult_order)-1)/m)\n", " \n", " for j in range(1, m):\n", " fct_j = self.jth_component(j)\n", " h = Rx(fct_j*f^(M*j))\n", " j1 = (p^(mult_order-1)*j)%m\n", " B = floor(p^(mult_order-1)*j/m)\n", " result += superelliptic_form(C, polynomial_part(p, h)/(f^B*y^(j1)))\n", " return result \n", " \n", " def coordinates(self):\n", " C = self.curve\n", " F = C.base_ring\n", " m = C.exponent\n", " Rx. = PolynomialRing(F)\n", " Fx = FractionField(Rx)\n", " FxRy. = PolynomialRing(Fx)\n", " g = C.genus()\n", " degrees_holo = C.degrees_holomorphic_differentials()\n", " degrees_holo_inv = {b:a for a, b in degrees_holo.items()}\n", " basis = C.holomorphic_differentials_basis()\n", " \n", " for j in range(1, m):\n", " omega_j = Fx(self.jth_component(j))\n", " if omega_j != Fx(0):\n", " d = degree_of_rational_fctn(omega_j, F)\n", " index = degrees_holo_inv[(d, j)]\n", " a = coeff_of_rational_fctn(omega_j, F)\n", " a1 = coeff_of_rational_fctn(basis[index].jth_component(j), F)\n", " elt = self - (a/a1)*basis[index]\n", " return elt.coordinates() + a/a1*vector([F(i == index) for i in range(0, g)])\n", " \n", " return vector(g*[0])\n", " \n", " def jth_component(self, j):\n", " g = self.form\n", " C = self.curve\n", " F = C.base_ring\n", " Rx. = PolynomialRing(F)\n", " Fx = FractionField(Rx)\n", " FxRy. = PolynomialRing(Fx)\n", " Fxy = FractionField(FxRy)\n", " Ryinv. = PolynomialRing(Fx)\n", " g = Fxy(g)\n", " g = g(y = 1/y_inv)\n", " g = Ryinv(g)\n", " return coff(g, j)\n", " \n", " def is_regular_on_U0(self):\n", " C = self.curve\n", " F = C.base_ring\n", " m = C.exponent\n", " Rx. = PolynomialRing(F)\n", " for j in range(1, m):\n", " if self.jth_component(j) not in Rx:\n", " return 0\n", " return 1\n", " \n", " def is_regular_on_Uinfty(self):\n", " C = self.curve\n", " F = C.base_ring\n", " m = C.exponent\n", " f = C.polynomial\n", " r = f.degree()\n", " delta = GCD(m, r)\n", " M = m/delta\n", " R = r/delta\n", " \n", " for j in range(1, m):\n", " A = self.jth_component(j)\n", " d = degree_of_rational_fctn(A, F)\n", " if(-d*M + j*R -(M+1)<0):\n", " return 0\n", " return 1\n", "\n", " def expansion_at_infty(self, i = 0, prec=10):\n", " g = self.form\n", " C = self.curve\n", " g = superelliptic_function(C, g)\n", " g = g.expansion_at_infty(i = i, prec=prec)\n", " x_series = superelliptic_function(C, x).expansion_at_infty(i = i, prec=prec)\n", " dx_series = x_series.derivative()\n", " return g*dx_series\n", " \n", "class superelliptic_cech:\n", " def __init__(self, C, omega, fct):\n", " self.omega0 = omega\n", " self.omega8 = omega - fct.diffn()\n", " self.f = fct\n", " self.curve = C\n", " \n", " def __add__(self, other):\n", " C = self.curve\n", " return superelliptic_cech(C, self.omega0 + other.omega0, self.f + other.f)\n", " \n", " def __sub__(self, other):\n", " C = self.curve\n", " return superelliptic_cech(C, self.omega0 - other.omega0, self.f - other.f)\n", "\n", " def __rmul__(self, constant):\n", " C = self.curve\n", " w1 = self.omega0.form\n", " f1 = self.f.function\n", " w2 = superelliptic_form(C, constant*w1)\n", " f2 = superelliptic_function(C, constant*f1)\n", " return superelliptic_cech(C, w2, f2) \n", " \n", " def __repr__(self):\n", " return \"(\" + str(self.omega0) + \", \" + str(self.f) + \", \" + str(self.omega8) + \")\" \n", " \n", " def verschiebung(self):\n", " C = self.curve\n", " omega = self.omega0\n", " F = C.base_ring\n", " Rx. = PolynomialRing(F)\n", " return superelliptic_cech(C, omega.cartier(), superelliptic_function(C, Rx(0)))\n", " \n", " def frobenius(self):\n", " C = self.curve\n", " fct = self.f.function\n", " p = C.characteristic\n", " Rx. = PolynomialRing(F)\n", " return superelliptic_cech(C, superelliptic_form(C, Rx(0)), superelliptic_function(C, fct^p))\n", "\n", " def coordinates(self):\n", " C = self.curve\n", " F = C.base_ring\n", " m = C.exponent\n", " Rx. = PolynomialRing(F)\n", " Fx = FractionField(Rx)\n", " FxRy. = PolynomialRing(Fx)\n", " g = C.genus()\n", " degrees_holo = C.degrees_holomorphic_differentials()\n", " degrees_holo_inv = {b:a for a, b in degrees_holo.items()}\n", " degrees0 = C.degrees_de_rham0()\n", " degrees0_inv = {b:a for a, b in degrees0.items()}\n", " degrees1 = C.degrees_de_rham1()\n", " degrees1_inv = {b:a for a, b in degrees1.items()}\n", " basis = C.de_rham_basis()\n", " \n", " omega = self.omega0\n", " fct = self.f\n", " \n", " if fct.function == Rx(0) and omega.form != Rx(0):\n", " for j in range(1, m):\n", " omega_j = Fx(omega.jth_component(j))\n", " if omega_j != Fx(0):\n", " d = degree_of_rational_fctn(omega_j, F)\n", " index = degrees_holo_inv[(d, j)]\n", " a = coeff_of_rational_fctn(omega_j, F)\n", " a1 = coeff_of_rational_fctn(basis[index].omega0.jth_component(j), F)\n", " elt = self - (a/a1)*basis[index]\n", " return elt.coordinates() + a/a1*vector([F(i == index) for i in range(0, 2*g)])\n", " \n", " for j in range(1, m):\n", " fct_j = Fx(fct.jth_component(j))\n", " if (fct_j != Rx(0)):\n", " d = degree_of_rational_fctn(fct_j, p)\n", " \n", " if (d, j) in degrees1.values():\n", " index = degrees1_inv[(d, j)]\n", " a = coeff_of_rational_fctn(fct_j, F)\n", " elt = self - (a/m)*basis[index]\n", " return elt.coordinates() + a/m*vector([F(i == index) for i in range(0, 2*g)])\n", " \n", " if d<0:\n", " a = coeff_of_rational_fctn(fct_j, F)\n", " h = superelliptic_function(C, FxRy(a*y^j*x^d))\n", " elt = superelliptic_cech(C, self.omega0, self.f - h)\n", " return elt.coordinates()\n", " \n", " if (fct_j != Rx(0)):\n", " G = superelliptic_function(C, y^j*x^d)\n", " a = coeff_of_rational_fctn(fct_j, F)\n", " elt =self - a*superelliptic_cech(C, diffn(G), G)\n", " return elt.coordinates()\n", "\n", " return vector(2*g*[0])\n", " \n", " def is_cocycle(self):\n", " w0 = self.omega0\n", " w8 = self.omega8\n", " fct = self.f\n", " if not w0.is_regular_on_U0() and not w8.is_regular_on_Uinfty():\n", " return('w0 & w8')\n", " if not w0.is_regular_on_U0():\n", " return('w0')\n", " if not w8.is_regular_on_Uinfty():\n", " return('w8')\n", " if w0.is_regular_on_U0() and w8.is_regular_on_Uinfty():\n", " return 1\n", " return 0\n", "\n", "#Auxilliary. If f = f1/f2 is a rational function, return deg f_1 - deg f_2.\n", "def degree_of_rational_fctn(f, F):\n", " Rx. = PolynomialRing(F)\n", " Fx = FractionField(Rx)\n", " f = Fx(f)\n", " f1 = f.numerator()\n", " f2 = f.denominator()\n", " d1 = f1.degree()\n", " d2 = f2.degree()\n", " return(d1 - d2)\n", "\n", "#Auxilliary. If f = f1/f2 is a rational function, return (leading coeff of f1)/(leading coeff of f2).\n", "def coeff_of_rational_fctn(f, F):\n", " Rx. = PolynomialRing(F)\n", " Fx = FractionField(Rx)\n", " f = Fx(f)\n", " if f == Rx(0):\n", " return 0\n", " f1 = f.numerator()\n", " f2 = f.denominator()\n", " d1 = f1.degree()\n", " d2 = f2.degree()\n", " a1 = f1.coefficients(sparse = false)[d1]\n", " a2 = f2.coefficients(sparse = false)[d2]\n", " return(a1/a2)\n", "\n", "#Auxilliary. Given polynomial f(x) and integer d, return\n", "#coefficient of x^d in f (and 0 is deg(f)= i+1} a_j x^{j - i -1}\n", "def cut(f, i):\n", " R = f.parent()\n", " coeff = f.coefficients(sparse = false)\n", " return sum(R(x^(j-i-1)) * coeff[j] for j in range(i+1, f.degree() + 1))\n", "\n", "def polynomial_part(p, h):\n", " F = GF(p)\n", " Rx. = PolynomialRing(F)\n", " h = Rx(h)\n", " result = Rx(0)\n", " for i in range(0, h.degree()+1):\n", " if (i%p) == p-1:\n", " power = Integer((i-(p-1))/p)\n", " result += Integer(h[i]) * x^(power) \n", " return result\n", "\n", "#Find delta-th root of unity in field F\n", "def root_of_unity(F, delta):\n", " Rx. = PolynomialRing(F)\n", " cyclotomic = x^(delta) - 1\n", " for root, a in cyclotomic.roots():\n", " powers = [root^d for d in delta.divisors() if d!= delta]\n", " if 1 not in powers:\n", " return root\n", " \n", "def preimage(U, V, M): #preimage of subspace U under M\n", " basis_preimage = M.right_kernel().basis()\n", " imageU = U.intersection(M.transpose().image())\n", " basis = imageU.basis()\n", " for v in basis:\n", " w = M.solve_right(v)\n", " basis_preimage = basis_preimage + [w]\n", " return V.subspace(basis_preimage)\n", "\n", "def image(U, V, M):\n", " basis = U.basis()\n", " basis_image = []\n", " for v in basis:\n", " basis_image += [M*v]\n", " return V.subspace(basis_image)\n", "\n", "def flag(F, V, p, test = 0):\n", " dim = F.dimensions()[0]\n", " space = VectorSpace(GF(p), dim)\n", " flag_subspaces = (dim+1)*[0]\n", " flag_used = (dim+1)*[0]\n", " final_type = (dim+1)*['?']\n", " \n", " flag_subspaces[dim] = space\n", " flag_used[dim] = 1\n", " \n", " \n", " while 1 in flag_used:\n", " index = flag_used.index(1)\n", " flag_used[index] = 0\n", " U = flag_subspaces[index]\n", " U_im = image(U, space, V)\n", " d_im = U_im.dimension()\n", " final_type[index] = d_im\n", " U_pre = preimage(U, space, F)\n", " d_pre = U_pre.dimension()\n", " \n", " if flag_subspaces[d_im] == 0:\n", " flag_subspaces[d_im] = U_im\n", " flag_used[d_im] = 1\n", " \n", " if flag_subspaces[d_pre] == 0:\n", " flag_subspaces[d_pre] = U_pre\n", " flag_used[d_pre] = 1\n", " \n", " if test == 1:\n", " print('(', final_type, ')')\n", " \n", " for i in range(0, dim+1):\n", " if final_type[i] == '?' and final_type[dim - i] != '?':\n", " i1 = dim - i\n", " final_type[i] = final_type[i1] - i1 + dim/2\n", " \n", " final_type[0] = 0\n", " \n", " for i in range(1, dim+1):\n", " if final_type[i] == '?':\n", " prev = final_type[i-1]\n", " if prev != '?' and prev in final_type[i+1:]:\n", " final_type[i] = prev\n", " \n", " for i in range(1, dim+1):\n", " if final_type[i] == '?':\n", " final_type[i] = min(final_type[i-1] + 1, dim/2)\n", " \n", " if is_final(final_type, dim/2):\n", " return final_type[1:dim/2 + 1]\n", " print('error:', final_type[1:dim/2 + 1])\n", " \n", "def is_final(final_type, dim):\n", " n = len(final_type)\n", " if final_type[0] != 0:\n", " return 0\n", " \n", " if final_type[n-1] != dim:\n", " return 0\n", " \n", " for i in range(1, n):\n", " if final_type[i] != final_type[i - 1] and final_type[i] != final_type[i - 1] + 1:\n", " return 0\n", " return 1" ] }, { "cell_type": "code", "execution_count": 2, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ "class as_cover:\n", " def __init__(self, C, list_of_fcts, prec = 10):\n", " self.quotient = C\n", " self.functions = list_of_fcts\n", " self.height = len(list_of_fcts)\n", " F = C.base_ring\n", " self.base_ring = F\n", " p = C.characteristic\n", " self.characteristic = p\n", " self.prec = prec\n", " f = C.polynomial\n", " m = C.exponent\n", " r = f.degree()\n", " delta = GCD(m, r)\n", " self.nb_of_pts_at_infty = delta\n", " Rxy. = PolynomialRing(F, 2)\n", " Rt. = LaurentSeriesRing(F, default_prec=prec)\n", "\n", " all_x_series = []\n", " all_y_series = []\n", " all_z_series = []\n", " all_dx_series = []\n", " all_jumps = []\n", "\n", " for i in range(delta):\n", " x_series = superelliptic_function(C, x).expansion_at_infty(i = i, prec=prec)\n", " y_series = superelliptic_function(C, y).expansion_at_infty(i = i, prec=prec)\n", " z_series = []\n", " jumps = []\n", " n = len(list_of_fcts)\n", " list_of_power_series = [g.expansion_at_infty(i = i, prec=prec) for g in list_of_fcts]\n", " for i in range(n):\n", " power_series = list_of_power_series[i]\n", " jump, correction, t_old, z = artin_schreier_transform(power_series, prec = prec)\n", " x_series = x_series(t = t_old)\n", " y_series = y_series(t = t_old)\n", " z_series = [zi(t = t_old) for zi in z_series]\n", " z_series += [z]\n", " jumps += [jump]\n", " list_of_power_series = [g(t = t_old) for g in list_of_power_series]\n", " \n", " all_jumps += [jumps]\n", " all_x_series += [x_series]\n", " all_y_series += [y_series]\n", " all_z_series += [z_series]\n", " all_dx_series += [x_series.derivative()]\n", " self.jumps = all_jumps\n", " self.x = all_x_series\n", " self.y = all_y_series\n", " self.z = all_z_series\n", " self.dx = all_dx_series\n", " \n", " def __repr__(self):\n", " n = self.height\n", " p = self.characteristic\n", " if n==1:\n", " return \"(Z/p)-cover of \" + str(self.quotient)+\" with the equation:\\n z^\" + str(p) + \" - z = \" + str(self.functions[0])\n", " \n", " result = \"(Z/p)^\"+str(self.height)+ \"-cover of \" + str(self.quotient)+\" with the equations:\\n\"\n", " for i in range(n):\n", " result += 'z' + str(i) + \"^\" + str(p) + \" - z\" + str(i) + \" = \" + str(self.functions[i]) + \"\\n\"\n", " return result\n", " \n", " def genus(self):\n", " jumps = self.jumps\n", " gY = self.quotient.genus()\n", " n = self.height\n", " delta = self.nb_of_pts_at_infty\n", " p = self.characteristic\n", " return p^n*gY + (p^n - 1)*(delta - 1) + sum(p^(n-j-1)*(jumps[i][j]-1)*(p-1)/2 for j in range(n) for i in range(delta))\n", " \n", " def exponent_of_different(self, i = 0):\n", " jumps = self.jumps\n", " n = self.height\n", " delta = self.nb_of_pts_at_infty\n", " p = self.characteristic\n", " return sum(p^(n-j-1)*(jumps[i][j]+1)*(p-1) for j in range(n))\n", "\n", " def exponent_of_different_prim(self, i = 0):\n", " jumps = self.jumps\n", " n = self.height\n", " delta = self.nb_of_pts_at_infty\n", " p = self.characteristic\n", " return sum(p^(n-j-1)*(jumps[i][j])*(p-1) for j in range(n))\n", " \n", " def holomorphic_differentials_basis(self, threshold = 8):\n", " from itertools import product\n", " x_series = self.x\n", " y_series = self.y\n", " z_series = self.z\n", " dx_series = self.dx\n", " delta = self.nb_of_pts_at_infty\n", " p = self.characteristic\n", " n = self.height\n", " prec = self.prec\n", " C = self.quotient\n", " F = self.base_ring\n", " m = C.exponent\n", " r = C.polynomial.degree()\n", " variable_names = 'x, y'\n", " for i in range(n):\n", " variable_names += ', z' + str(i)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " Rt. = LaurentSeriesRing(F, default_prec=prec)\n", " #Tworzymy zbiór S form z^i x^j y^k dx/y o waluacji >= waluacja z^(p-1)*dx/y\n", " S = []\n", " RQxyz = FractionField(Rxyz)\n", " pr = [list(GF(p)) for _ in range(n)]\n", " for i in range(0, threshold*r):\n", " for j in range(0, m):\n", " for k in product(*pr):\n", " eta = as_form(self, x^i * prod(z[i1]^(k[i1]) for i1 in range(n))/y^j)\n", " eta_exp = eta.expansion_at_infty()\n", " S += [(eta, eta_exp)]\n", "\n", " forms = holomorphic_combinations(S)\n", " \n", " for i in range(1, delta):\n", " forms = [(omega, omega.expansion_at_infty(i = i)) for omega in forms]\n", " forms = holomorphic_combinations(forms)\n", " \n", " if len(forms) < self.genus():\n", " print(\"I haven't found all forms.\")\n", " return holomorphic_differentials_basis(self, threshold = threshold + 1)\n", " if len(forms) > self.genus():\n", " print(\"Increase precision.\")\n", " return forms\n", "\n", " def at_most_poles(self, pole_order, threshold = 8):\n", " \"\"\" Find fcts with pole order in infty's at most pole_order. Threshold gives a bound on powers of x in the function.\n", " If you suspect that you haven't found all the functions, you may increase it.\"\"\"\n", " from itertools import product\n", " x_series = self.x\n", " y_series = self.y\n", " z_series = self.z\n", " delta = self.nb_of_pts_at_infty\n", " p = self.characteristic\n", " n = self.height\n", " prec = self.prec\n", " C = self.quotient\n", " F = self.base_ring\n", " m = C.exponent\n", " r = C.polynomial.degree()\n", " variable_names = 'x, y'\n", " for i in range(n):\n", " variable_names += ', z' + str(i)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " Rt. = LaurentSeriesRing(F, default_prec=prec)\n", " #Tworzymy zbiór S form z^i x^j y^k dx/y o waluacji >= waluacja z^(p-1)*dx/y\n", " S = []\n", " RQxyz = FractionField(Rxyz)\n", " pr = [list(GF(p)) for _ in range(n)]\n", " for i in range(0, threshold*r):\n", " for j in range(0, m):\n", " for k in product(*pr):\n", " eta = as_function(self, x^i * prod(z[i1]^(k[i1]) for i1 in range(n))*y^j)\n", " eta_exp = eta.expansion_at_infty()\n", " S += [(eta, eta_exp)]\n", "\n", " forms = holomorphic_combinations_fcts(S, pole_order)\n", "\n", " for i in range(1, delta):\n", " forms = [(omega, omega.expansion_at_infty(i = i)) for omega in forms]\n", " forms = holomorphic_combinations_fcts(forms, pole_order)\n", " \n", " return forms\n", " \n", " def magical_element(self, threshold = 8):\n", " list_of_elts = self.at_most_poles(self.exponent_of_different_prim(), threshold)\n", " result = []\n", " for a in list_of_elts:\n", " if a.trace().function != 0:\n", " result += [a]\n", " return result\n", "\n", " def pseudo_magical_element(self, threshold = 8):\n", " list_of_elts = self.at_most_poles(self.exponent_of_different(), threshold)\n", " result = []\n", " for a in list_of_elts:\n", " if a.trace().function != 0:\n", " result += [a]\n", " return result\n", "\n", " def at_most_poles_forms(self, pole_order, threshold = 8):\n", " \"\"\"Find forms with pole order in all the points at infty equat at most to pole_order. Threshold gives a bound on powers of x in the form.\n", " If you suspect that you haven't found all the functions, you may increase it.\"\"\"\n", " from itertools import product\n", " x_series = self.x\n", " y_series = self.y\n", " z_series = self.z\n", " delta = self.nb_of_pts_at_infty\n", " p = self.characteristic\n", " n = self.height\n", " prec = self.prec\n", " C = self.quotient\n", " F = self.base_ring\n", " m = C.exponent\n", " r = C.polynomial.degree()\n", " variable_names = 'x, y'\n", " for i in range(n):\n", " variable_names += ', z' + str(i)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " Rt. = LaurentSeriesRing(F, default_prec=prec)\n", " #Tworzymy zbiór S form z^i x^j y^k dx/y o waluacji >= waluacja z^(p-1)*dx/y\n", " S = []\n", " RQxyz = FractionField(Rxyz)\n", " pr = [list(GF(p)) for _ in range(n)]\n", " for i in range(0, threshold*r):\n", " for j in range(0, m):\n", " for k in product(*pr):\n", " eta = as_form(self, x^i * prod(z[i1]^(k[i1]) for i1 in range(n))/y^j)\n", " eta_exp = eta.expansion_at_infty()\n", " S += [(eta, eta_exp)]\n", "\n", " forms = holomorphic_combinations_forms(S, pole_order)\n", "\n", " for i in range(1, delta):\n", " forms = [(omega, omega.expansion_at_infty(i = i)) for omega in forms]\n", " forms = holomorphic_combinations_forms(forms, pole_order)\n", " \n", " return forms\n", "\n", "def holomorphic_combinations(S):\n", " \"\"\"Given a list S of pairs (form, corresponding Laurent series at some pt), find their combinations holomorphic at that pt.\"\"\"\n", " C_AS = S[0][0].curve\n", " p = C_AS.characteristic\n", " F = C_AS.base_ring\n", " prec = C_AS.prec\n", " Rt. = LaurentSeriesRing(F, default_prec=prec)\n", " RtQ = FractionField(Rt)\n", " minimal_valuation = min([g[1].valuation() for g in S])\n", " if minimal_valuation >= 0:\n", " return [s[0] for s in S]\n", " list_of_lists = [] #to będzie lista złożona z list współczynników część nieholomorficznych rozwinięcia form z S\n", " for eta, eta_exp in S:\n", " a = -minimal_valuation + eta_exp.valuation()\n", " list_coeffs = a*[0] + eta_exp.list() + (-minimal_valuation)*[0]\n", " list_coeffs = list_coeffs[:-minimal_valuation]\n", " list_of_lists += [list_coeffs]\n", " M = matrix(F, list_of_lists)\n", " V = M.kernel() #chcemy wyzerować części nieholomorficzne, biorąc kombinacje form z S\n", "\n", "\n", " # Sprawdzamy, jakim formom odpowiadają elementy V.\n", " forms = []\n", " for vec in V.basis():\n", " forma_holo = as_form(C_AS, 0)\n", " forma_holo_power_series = Rt(0)\n", " for vec_wspolrzedna, elt_S in zip(vec, S):\n", " eta = elt_S[0]\n", " #eta_exp = elt_S[1]\n", " forma_holo += vec_wspolrzedna*eta\n", " #forma_holo_power_series += vec_wspolrzedna*eta_exp\n", " forms += [forma_holo]\n", " return forms\n", "\n", "class as_function:\n", " def __init__(self, C, g):\n", " self.curve = C\n", " F = C.base_ring\n", " n = C.height\n", " variable_names = 'x, y'\n", " for i in range(n):\n", " variable_names += ', z' + str(i)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " self.function = RxyzQ(g)\n", " #self.function = as_reduction(AS, RxyzQ(g))\n", " \n", " def __repr__(self):\n", " return str(self.function)\n", "\n", " def __add__(self, other):\n", " C = self.curve\n", " g1 = self.function\n", " g2 = other.function\n", " return as_function(C, g1 + g2)\n", "\n", " def __sub__(self, other):\n", " C = self.curve\n", " g1 = self.function\n", " g2 = other.function\n", " return as_function(C, g1 - g2)\n", "\n", " def __rmul__(self, constant):\n", " C = self.curve\n", " g = self.function\n", " return as_function(C, constant*g)\n", " \n", " def __mul__(self, other):\n", " C = self.curve\n", " g1 = self.function\n", " g2 = other.function\n", " return as_function(C, g1*g2)\n", "\n", " def expansion_at_infty(self, i = 0):\n", " C = self.curve\n", " delta = C.nb_of_pts_at_infty\n", " F = C.base_ring\n", " x_series = C.x[i]\n", " y_series = C.y[i]\n", " z_series = C.z[i]\n", " n = C.height\n", " variable_names = 'x, y'\n", " for j in range(n):\n", " variable_names += ', z' + str(j)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " prec = C.prec\n", " Rt. = LaurentSeriesRing(F, default_prec=prec)\n", " g = self.function\n", " g = RxyzQ(g)\n", " sub_list = {x : x_series, y : y_series} | {z[j] : z_series[j] for j in range(n)}\n", " return g.substitute(sub_list)\n", "\n", " def group_action(self, ZN_tuple):\n", " C = self.curve\n", " n = C.height\n", " F = C.base_ring\n", " variable_names = 'x, y'\n", " for j in range(n):\n", " variable_names += ', z' + str(j)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " sub_list = {x : x, y : y} | {z[j] : z[j]+ZN_tuple[j] for j in range(n)}\n", " g = self.function\n", " return as_function(C, g.substitute(sub_list))\n", "\n", " def trace(self):\n", " C = self.curve\n", " C_super = C.quotient\n", " n = C.height\n", " F = C.base_ring\n", " variable_names = 'x, y'\n", " for j in range(n):\n", " variable_names += ', z' + str(j)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " g = self.function\n", " g = as_reduction(C, g)\n", " result = RxyzQ(0)\n", " g_num = Rxyz(numerator(g))\n", " g_den = Rxyz(denominator(g))\n", " z = prod(z[i] for i in range(n))^(p-1)\n", " for a in g_num.monomials():\n", " if (z.divides(a)):\n", " result += g_num.monomial_coefficient(a)*a/z\n", " result /= g_den\n", " result = as_reduction(C, result)\n", " Rxy. = PolynomialRing(F, 2)\n", " Qxy = FractionField(Rxy)\n", " return superelliptic_function(C_super, Qxy(result))\n", "\n", " def trace2(self):\n", " C = self.curve\n", " C_super = C.quotient\n", " n = C.height\n", " F = C.base_ring\n", " variable_names = 'x, y'\n", " for j in range(n):\n", " variable_names += ', z' + str(j)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " result = as_function(C, 0)\n", " for i in range(0, p):\n", " for j in range(0, p):\n", " result += self.group_action([i, j])\n", " result = result.function\n", " Rxy. = PolynomialRing(F, 2)\n", " Qxy = FractionField(Rxy)\n", " return superelliptic_function(C_super, Qxy(result)) \n", " \n", "class as_form:\n", " def __init__(self, C, g):\n", " self.curve = C\n", " n = C.height\n", " F = C.base_ring\n", " variable_names = 'x, y'\n", " for i in range(n):\n", " variable_names += ', z' + str(i)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " self.form = RxyzQ(g)\n", " \n", " def __repr__(self):\n", " return \"(\" + str(self.form)+\") * dx\"\n", " \n", " def expansion_at_infty(self, i = 0):\n", " C = self.curve\n", " delta = C.nb_of_pts_at_infty\n", " F = C.base_ring\n", " x_series = C.x[i]\n", " y_series = C.y[i]\n", " z_series = C.z[i]\n", " dx_series = C.dx[i]\n", " n = C.height\n", " variable_names = 'x, y'\n", " for j in range(n):\n", " variable_names += ', z' + str(j)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " prec = C.prec\n", " Rt. = LaurentSeriesRing(F, default_prec=prec)\n", " g = self.form\n", " sub_list = {x : x_series, y : y_series} | {z[j] : z_series[j] for j in range(n)}\n", " return g.substitute(sub_list)*dx_series\n", " \n", " def __add__(self, other):\n", " C = self.curve\n", " g1 = self.form\n", " g2 = other.form\n", " return as_form(C, g1 + g2)\n", " \n", " def __sub__(self, other):\n", " C = self.curve\n", " g1 = self.form\n", " g2 = other.form\n", " return as_form(C, g1 - g2)\n", " \n", " def __rmul__(self, constant):\n", " C = self.curve\n", " omega = self.form\n", " return as_form(C, constant*omega)\n", "\n", " def group_action(self, ZN_tuple):\n", " C = self.curve\n", " n = C.height\n", " F = C.base_ring\n", " variable_names = 'x, y'\n", " for j in range(n):\n", " variable_names += ', z' + str(j)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " sub_list = {x : x, y : y} | {z[j] : z[j]+ZN_tuple[j] for j in range(n)}\n", " g = self.form\n", " return as_form(C, g.substitute(sub_list))\n", "\n", " def coordinates(self, holo):\n", " \"\"\"Find coordinates of the given form self in terms of the basis forms in a list holo.\"\"\"\n", " C = self.curve\n", " n = C.height\n", " gC = C.genus()\n", " F = C.base_ring\n", " variable_names = 'x, y'\n", " for j in range(n):\n", " variable_names += ', z' + str(j)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " from sage.rings.polynomial.toy_variety import linear_representation\n", " return linear_representation(Rxyz(self.form), holo)\n", "\n", " def trace(self):\n", " C = self.curve\n", " C_super = C.quotient\n", " n = C.height\n", " F = C.base_ring\n", " variable_names = 'x, y'\n", " for j in range(n):\n", " variable_names += ', z' + str(j)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " g = self.form\n", " result = RxyzQ(0)\n", " g_num = Rxyz(numerator(g))\n", " g_den = Rxyz(denominator(g))\n", " z = prod(z[i] for i in range(n))^(p-1)\n", " for a in g_num.monomials():\n", " if (z.divides(a)):\n", " result += g_num.monomial_coefficient(a)*a/z\n", " result /= g_den\n", " Rxy. = PolynomialRing(F, 2)\n", " return superelliptic_form(C_super, Rxy(result))\n", "\n", " \n", " def trace2(self):\n", " C = self.curve\n", " C_super = C.quotient\n", " n = C.height\n", " F = C.base_ring\n", " variable_names = 'x, y'\n", " for j in range(n):\n", " variable_names += ', z' + str(j)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " result = as_form(C, 0)\n", " for i in range(0, p):\n", " for j in range(0, p):\n", " result += self.group_action([i, j])\n", " result = result.form\n", " Rxy. = PolynomialRing(F, 2)\n", " Qxy = FractionField(Rxy)\n", " return superelliptic_form(C_super, Qxy(result))\n", " \n", "# Given power_series, find its reverse (g with g \\circ power_series = id) with given precision\n", "\n", "def new_reverse(power_series, prec = 10):\n", " F = power_series.parent().base()\n", " Rt. = LaurentSeriesRing(F, default_prec=prec)\n", " RtQ = FractionField(Rt)\n", " power_series = RtQ(power_series)\n", " a = power_series.list()[0]\n", " g = 1/a*t\n", " n = 2\n", " while(n <= prec):\n", " aux = power_series(t = g) - t\n", " if aux.valuation() > n:\n", " b = 0\n", " else:\n", " b = aux.list()[0]\n", " g = g - b/a*t^n\n", " n += 1\n", " return g\n", "\n", "def artin_schreier_transform(power_series, prec = 10):\n", " \"\"\"Given a power_series, find correction such that power_series - (correction)^p +correction has valuation\n", " -jump non divisible by p. Also, express t (the variable) in terms of the uniformizer at infty on the curve\n", " z^p - z = power_series, where z = 1/t_new^(jump) and express z in terms of the new uniformizer.\"\"\"\n", " correction = 0\n", " F = power_series.parent().base()\n", " Rt. = LaurentSeriesRing(F, default_prec=prec)\n", " RtQ = FractionField(Rt)\n", " power_series = RtQ(power_series)\n", " if power_series.valuation() == +Infinity:\n", " return(0,0,t,0)\n", " while(power_series.valuation() % p == 0 and power_series.valuation() < 0):\n", " M = -power_series.valuation()/p\n", " coeff = power_series.list()[0] #wspolczynnik a_(-p) w f_AS\n", " correction += coeff.nth_root(p)*t^(-M)\n", " power_series = power_series - (coeff*t^(-p*M) - coeff.nth_root(p)*t^(-M))\n", " jump = max(-(power_series.valuation()), 0)\n", " try:\n", " T = ((power_series)^(-1)).nth_root(jump) #T is defined by power_series = 1/T^m\n", " except:\n", " print(\"no \", str(jump), \"-th root; divide by\", power_series.list()[0])\n", " return (jump, power_series.list()[0])\n", " T_rev = new_reverse(T, prec = prec)\n", " t_old = T_rev(t^p/(1 - t^((p-1)*jump)).nth_root(jump))\n", " z = 1/t^(jump) + Rt(correction)(t = t_old)\n", " return(jump, correction, t_old, z)\n", "\n", "\n", "def are_forms_linearly_dependent(set_of_forms):\n", " from sage.rings.polynomial.toy_variety import is_linearly_dependent\n", " C = set_of_forms[0].curve\n", " F = C.base_ring\n", " n = C.height\n", " variable_names = 'x, y'\n", " for i in range(n):\n", " variable_names += ', z' + str(i)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " denominators = prod(denominator(omega.form) for omega in set_of_forms)\n", " return is_linearly_dependent([Rxyz(denominators*omega.form) for omega in set_of_forms])\n", "\n", "def group_action_matrices(C_AS):\n", " F = C_AS.base_ring\n", " n = C_AS.height\n", " holo = C_AS.holomorphic_differentials_basis()\n", " holo_forms = [omega.form for omega in holo]\n", " denom = LCM([denominator(omega) for omega in holo_forms])\n", " variable_names = 'x, y'\n", " for j in range(n):\n", " variable_names += ', z' + str(j)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " holo_forms = [Rxyz(omega*denom) for omega in holo_forms]\n", " A = [[] for i in range(n)]\n", " for omega in holo:\n", " for i in range(n):\n", " ei = n*[0]\n", " ei[i] = 1\n", " omega1 = omega.group_action(ei)\n", " omega1 = denom * omega1\n", " v1 = omega1.coordinates(holo_forms)\n", " A[i] += [v1]\n", " for i in range(n):\n", " A[i] = matrix(F, A[i])\n", " A[i] = A[i].transpose()\n", " return A\n", "\n", "def group_action_matrices_log(C_AS):\n", " F = C_AS.base_ring\n", " n = C_AS.height\n", " holo = C_AS.at_most_poles_forms(1)\n", " holo_forms = [omega.form for omega in holo]\n", " denom = LCM([denominator(omega) for omega in holo_forms])\n", " variable_names = 'x, y'\n", " for j in range(n):\n", " variable_names += ', z' + str(j)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " holo_forms = [Rxyz(omega*denom) for omega in holo_forms]\n", " A = [[] for i in range(n)]\n", " for omega in holo:\n", " for i in range(n):\n", " ei = n*[0]\n", " ei[i] = 1\n", " omega1 = omega.group_action(ei)\n", " omega1 = denom * omega1\n", " v1 = omega1.coordinates(holo_forms)\n", " A[i] += [v1]\n", " for i in range(n):\n", " A[i] = matrix(F, A[i])\n", " A[i] = A[i].transpose()\n", " return A\n", "\n", "\n", "#given a set S of (form, corresponding Laurent series at some pt), find their combinations holomorphic at that pt\n", "def holomorphic_combinations_fcts(S, pole_order):\n", " C_AS = S[0][0].curve\n", " p = C_AS.characteristic\n", " F = C_AS.base_ring\n", " prec = C_AS.prec\n", " Rt. = LaurentSeriesRing(F, default_prec=prec)\n", " RtQ = FractionField(Rt)\n", " minimal_valuation = min([Rt(g[1]).valuation() for g in S])\n", " if minimal_valuation >= -pole_order:\n", " return [s[0] for s in S]\n", " list_of_lists = [] #to będzie lista złożona z list współczynników część nieholomorficznych rozwinięcia form z S\n", " for eta, eta_exp in S:\n", " a = -minimal_valuation + Rt(eta_exp).valuation()\n", " list_coeffs = a*[0] + Rt(eta_exp).list() + (-minimal_valuation)*[0]\n", " list_coeffs = list_coeffs[:-minimal_valuation - pole_order]\n", " list_of_lists += [list_coeffs]\n", " M = matrix(F, list_of_lists)\n", " V = M.kernel() #chcemy wyzerować części nieholomorficzne, biorąc kombinacje form z S\n", "\n", "\n", " # Sprawdzamy, jakim formom odpowiadają elementy V.\n", " forms = []\n", " for vec in V.basis():\n", " forma_holo = as_function(C_AS, 0)\n", " forma_holo_power_series = Rt(0)\n", " for vec_wspolrzedna, elt_S in zip(vec, S):\n", " eta = elt_S[0]\n", " #eta_exp = elt_S[1]\n", " forma_holo += vec_wspolrzedna*eta\n", " #forma_holo_power_series += vec_wspolrzedna*eta_exp\n", " forms += [forma_holo]\n", " return forms\n", "\n", "#given a set S of (form, corresponding Laurent series at some pt), find their combinations holomorphic at that pt\n", "def holomorphic_combinations_forms(S, pole_order):\n", " C_AS = S[0][0].curve\n", " p = C_AS.characteristic\n", " F = C_AS.base_ring\n", " prec = C_AS.prec\n", " Rt. = LaurentSeriesRing(F, default_prec=prec)\n", " RtQ = FractionField(Rt)\n", " minimal_valuation = min([Rt(g[1]).valuation() for g in S])\n", " if minimal_valuation >= -pole_order:\n", " return [s[0] for s in S]\n", " list_of_lists = [] #to będzie lista złożona z list współczynników część nieholomorficznych rozwinięcia form z S\n", " for eta, eta_exp in S:\n", " a = -minimal_valuation + Rt(eta_exp).valuation()\n", " list_coeffs = a*[0] + Rt(eta_exp).list() + (-minimal_valuation)*[0]\n", " list_coeffs = list_coeffs[:-minimal_valuation - pole_order]\n", " list_of_lists += [list_coeffs]\n", " M = matrix(F, list_of_lists)\n", " V = M.kernel() #chcemy wyzerować części nieholomorficzne, biorąc kombinacje form z S\n", "\n", "\n", " # Sprawdzamy, jakim formom odpowiadają elementy V.\n", " forms = []\n", " for vec in V.basis():\n", " forma_holo = as_form(C_AS, 0)\n", " forma_holo_power_series = Rt(0)\n", " for vec_wspolrzedna, elt_S in zip(vec, S):\n", " eta = elt_S[0]\n", " #eta_exp = elt_S[1]\n", " forma_holo += vec_wspolrzedna*eta\n", " #forma_holo_power_series += vec_wspolrzedna*eta_exp\n", " forms += [forma_holo]\n", " return forms\n", "\n", "#print only forms that are log at the branch pts, but not holomorphic\n", "def only_log_forms(C_AS):\n", " list1 = AS.at_most_poles_forms(0)\n", " list2 = AS.at_most_poles_forms(1)\n", " result = []\n", " for a in list2:\n", " if not(are_forms_linearly_dependent(list1 + result + [a])):\n", " result += [a]\n", " return result\n", "\n", "def magmathis(A, B, text = False):\n", " \"\"\"Find decomposition of Z/p^2-module given by matrices A, B into indecomposables using magma.\n", " If text = True, print the command for Magma. Else - return the output of Magma free.\"\"\"\n", " q = parent(A).base_ring().order()\n", " n = A.dimensions()[0]\n", " A = str(list(A))\n", " B = str(list(B))\n", " A = A.replace(\"(\", \"\")\n", " A = A.replace(\")\", \"\")\n", " B = B.replace(\"(\", \"\")\n", " B = B.replace(\")\", \"\")\n", " result = \"A := MatrixAlgebra;\"\n", " result += \"M := RModule(RSpace(GF(\"+str(q)+\"),\" + str(n) + \"), A);\"\n", " result += \"IndecomposableSummands(M);\"\n", " if text:\n", " return result\n", " print(magma_free(result))\n", " \n", " \n", "def as_reduction(AS, fct):\n", " n = AS.height\n", " F = AS.base_ring\n", " variable_names = 'x, y'\n", " for i in range(n):\n", " variable_names += ', z' + str(i)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " ff = AS.functions\n", " ff = [RxyzQ(F.function) for F in ff]\n", " fct = RxyzQ(fct)\n", " fct1 = numerator(fct)\n", " fct2 = denominator(fct)\n", " if fct2 != 1:\n", " return as_reduction(AS, fct1)/as_reduction(AS, fct2)\n", " \n", " result = RxyzQ(0)\n", " change = 0\n", " for a in fct1.monomials():\n", " degrees_zi = [a.degree(z[i]) for i in range(n)]\n", " d_div = [a.degree(z[i])//p for i in range(n)]\n", " if d_div != n*[0]:\n", " change = 1\n", " d_rem = [a.degree(z[i])%p for i in range(n)]\n", " monomial = fct1.coefficient(a)*x^(a.degree(x))*y^(a.degree(y))*prod(z[i]^(d_rem[i]) for i in range(n))*prod((z[i] + ff[i])^(d_div[i]) for i in range(n))\n", " result += RxyzQ(fct1.coefficient(a)*x^(a.degree(x))*y^(a.degree(y))*prod(z[i]^(d_rem[i]) for i in range(n))*prod((z[i] + ff[i])^(d_div[i]) for i in range(n)))\n", " \n", " \n", " if change == 0:\n", " return RxyzQ(result)\n", " else:\n", " print(fct, '\\n')\n", " return as_reduction(AS, RxyzQ(result))" ] }, { "cell_type": "markdown", "metadata": { "collapsed": false }, "source": [ "# Testy" ] }, { "cell_type": "code", "execution_count": 3, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "325\n" ] }, { "name": "stdout", "output_type": "stream", "text": [ "325\n" ] } ], "source": [ "p = 5\n", "m = 2\n", "Rx. = PolynomialRing(GF(p))\n", "f = x^3 + x^2 + 1\n", "C_super = superelliptic(f, m)\n", "Rxy. = PolynomialRing(GF(p), 2)\n", "fArS1 = superelliptic_function(C_super, y*x)\n", "fArS2 = superelliptic_function(C_super, y*x^2)\n", "fArS3 = superelliptic_function(C_super, y)\n", "AS1 = as_cover(C_super, [fArS1, fArS2, fArS3], prec=500)\n", "print(AS1.genus())\n", "AS2 = as_cover(C_super, [fArS2, fArS3, fArS1], prec=500)\n", "print(AS2.genus())" ] }, { "cell_type": "code", "execution_count": 0, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ "p = 5\n", "m = 2\n", "Rx. = PolynomialRing(GF(p))\n", "f = x^3 + x^2 + 1\n", "C_super = superelliptic(f, m)\n", "Rxy. = PolynomialRing(GF(p), 2)\n", "fArS1 = superelliptic_function(C_super, y*x)\n", "fArS2 = superelliptic_function(C_super, y*x^2)\n", "fArS3 = superelliptic_function(C_super, y)\n", "AS1 = as_cover(C_super, [fArS1, fArS2, fArS3], prec=1000)\n", "print(AS1.exponent_of_different())\n", "omega = as_form(AS1, 1/y)\n", "print(omega.expansion_at_infty())" ] }, { "cell_type": "markdown", "metadata": { "collapsed": false }, "source": [ "# Brudnopis" ] }, { "cell_type": "code", "execution_count": 8, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "[1,\n", " z1,\n", " z1^2,\n", " z0,\n", " z0*z1,\n", " -x*z0^2 + z0*z1^2,\n", " z0^2,\n", " z0^2*z1,\n", " z0^2*z1^2 + x*z0*z1 + x^2,\n", " y,\n", " -y*z0^2 + y*z1,\n", " y*z0,\n", " x,\n", " -x*z0^2 + x*z1,\n", " x*z0]" ] }, "execution_count": 8, "metadata": { }, "output_type": "execute_result" } ], "source": [ "AS.at_most_poles(14)" ] }, { "cell_type": "code", "execution_count": 30, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ ] }, { "cell_type": "code", "execution_count": 35, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "[((z0^2*z1^2 + x*z0*z1 + x^2)/y) * dx]" ] }, "execution_count": 35, "metadata": { }, "output_type": "execute_result" } ], "source": [ "only_log_forms(AS)" ] }, { "cell_type": "code", "execution_count": 162, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "I haven't found all forms.\n", "[1]\n" ] } ], "source": [ "dpp = AS.exponent_of_different_prim()\n", "print(at_most_poles(AS, dpp))" ] }, { "cell_type": "code", "execution_count": 66, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "((z0*z1^2 - z1^2 + (a + 1)*x)/y) * dx\n" ] } ], "source": [ "Rxy. = PolynomialRing(F, 4)\n", "omega = as_form(AS, z0^2*z1^2/y - (a+1)*z0*x/y)\n", "print(omega - omega.group_action([1, 0]))" ] }, { "cell_type": "code", "execution_count": 78, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "a*t^-10 + 2*t^-8 + a*t^-6 + (a + 1)*t^-2 + O(t^288)" ] }, "execution_count": 78, "metadata": { }, "output_type": "execute_result" } ], "source": [ "g = as_function(AS, z0^2*z1^2+z0*(a+1)*x)\n", "g.expansion_at_infty(i = 0)" ] }, { "cell_type": "code", "execution_count": 18, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ "p = 7\n", "m = 2\n", "F = GF(p)\n", "Rx. = PolynomialRing(F)\n", "f = x^3 + 1\n", "C_super = superelliptic(f, m)\n", "\n", "Rxy. = PolynomialRing(F, 2)\n", "f1 = superelliptic_function(C_super, x^2*y)\n", "f2 = superelliptic_function(C_super, x^3)\n", "AS = as_cover(C_super, [f1, f2], prec=1000)\n", "#print(AS.magical_element())\n", "#print(AS.pseudo_magical_element())" ] }, { "cell_type": "code", "execution_count": 19, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "True\n", "True\n", "True\n" ] }, { "name": "stdout", "output_type": "stream", "text": [ "\n" ] } ], "source": [ "#A, B = group_action_matrices2_log(AS1)\n", "A, B = group_action_matrices2(AS)\n", "n = A.dimensions()[0]\n", "print(A*B == B*A)\n", "print(A^p == identity_matrix(n))\n", "print(B^p == identity_matrix(n))\n", "magmathis(A, B, p, deg = 1)" ] }, { "cell_type": "code", "execution_count": 10, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "t^-5 + 2*t + O(t^5)\n", "t^-6 + 3 + O(t^4)\n" ] } ], "source": [ "print(f1.expansion_at_infty())\n", "print(f2.expansion_at_infty())" ] }, { "cell_type": "code", "execution_count": 48, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "[z0^2]" ] }, "execution_count": 48, "metadata": { }, "output_type": "execute_result" } ], "source": [ "AS.magical_element()" ] }, { "cell_type": "code", "execution_count": 51, "metadata": { "collapsed": false, "scrolled": true }, "outputs": [ { "data": { "text/plain": [ "68" ] }, "execution_count": 51, "metadata": { }, "output_type": "execute_result" } ], "source": [ "AS1.exponent_of_different_prim()" ] }, { "cell_type": "code", "execution_count": 39, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "[[5, 2], [5, 2]]" ] }, "execution_count": 39, "metadata": { }, "output_type": "execute_result" } ], "source": [ "AS1.jumps" ] }, { "cell_type": "code", "execution_count": 7, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ ] }, { "cell_type": "code", "execution_count": 34, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "[\n", "RModule of dimension 1 over GF(3^2),\n", "RModule of dimension 1 over GF(3^2),\n", "RModule of dimension 3 over GF(3^2),\n", "RModule of dimension 3 over GF(3^2),\n", "RModule of dimension 6 over GF(3^2),\n", "RModule of dimension 6 over GF(3^2),\n", "RModule of dimension 8 over GF(3^2)\n", "]\n" ] } ], "source": [ "magmathis(A, B, 3, deg = 2)" ] }, { "cell_type": "code", "execution_count": 21, "metadata": { "collapsed": false, "scrolled": true }, "outputs": [ { "name": "stderr", "output_type": "stream", "text": [ "WARNING: Some output was deleted.\n" ] } ], "source": [ "magmatxtx(A, B, p)" ] }, { "cell_type": "markdown", "metadata": { "collapsed": false }, "source": [ "A := MatrixAlgebra;\n", "M := RModule(RSpace(GF(3^1),14), A);\n", "lll := IndecomposableSummands(M);\n", "Morphism(lll[2], M);" ] }, { "cell_type": "code", "execution_count": 25, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "3" ] }, "execution_count": 25, "metadata": { }, "output_type": "execute_result" } ], "source": [ "parent(A).base_ring().order()" ] }, { "cell_type": "code", "execution_count": 4, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ "Rx. = PolynomialRing(GF(3))\n", "C = superelliptic(x^3 - x, 2)" ] }, { "cell_type": "code", "execution_count": 5, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "[((1/y) dx, 0, (1/y) dx), ((x/y) dx, 2/x*y, ((-1)/(x*y)) dx)]" ] }, "execution_count": 5, "metadata": { }, "output_type": "execute_result" } ], "source": [ "C.de_rham_basis()" ] }, { "cell_type": "code", "execution_count": 5, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ "Rxy. = PolynomialRing(GF(p), 2)\n", "omega = superelliptic_form(C, 1/y^3)\n", "ff = superelliptic_function(C, x)\n", "cycle = superelliptic_cech(C, omega, ff)" ] }, { "cell_type": "code", "execution_count": 7, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "(0, 0)" ] }, "execution_count": 7, "metadata": { }, "output_type": "execute_result" } ], "source": [ "cycle.coordinates()" ] }, { "cell_type": "code", "execution_count": 3, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ "def dual_elt(AS, zmag):\n", " p = AS.characteristic\n", " n = AS.height\n", " group_elts = [(i, j) for i in range(p) for j in range(p)]\n", " variable_names = 'x, y'\n", " for i in range(n):\n", " variable_names += ', z' + str(i)\n", " Rxyz = PolynomialRing(F, n+2, variable_names)\n", " x, y = Rxyz.gens()[:2]\n", " z = Rxyz.gens()[2:]\n", " RxyzQ = FractionField(Rxyz)\n", " M = matrix(RxyzQ, p^n, p^n)\n", " for i in range(p^n):\n", " for j in range(p^n):\n", " M[i, j] = (zmag.group_action(group_elts[i])*zmag.group_action(group_elts[j])).trace2()\n", " main_det = M.determinant()\n", " zvee = as_function(AS, 0)\n", " for i in range(p^n):\n", " Mprim = matrix(RxyzQ, M)\n", " Mprim[:, i] = vector([(j == 0) for j in range(p^2)])\n", " fi = Mprim.determinant()/main_det\n", " zvee += fi*zmag.group_action(group_elts[i])\n", " return zvee" ] }, { "cell_type": "code", "execution_count": 4, "metadata": { "collapsed": false, "scrolled": true }, "outputs": [ ], "source": [ "p = 5\n", "m = 1\n", "F = GF(p)\n", "Rx. = PolynomialRing(F)\n", "f = x\n", "C_super = superelliptic(f, m)\n", "\n", "Rxy. = PolynomialRing(F, 2)\n", "f1 = superelliptic_function(C_super, x^2)\n", "f2 = superelliptic_function(C_super, x^3)\n", "AS = as_cover(C_super, [f1, f2], prec=500)\n", "zmag = (AS.magical_element())[0]" ] }, { "cell_type": "code", "execution_count": 5, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ "n = 2\n", "variable_names = 'x, y'\n", "for i in range(n):\n", " variable_names += ', z' + str(i)\n", "Rxyz = PolynomialRing(F, n+2, variable_names)\n", "x, y = Rxyz.gens()[:2]\n", "z = Rxyz.gens()[2:]\n", "f3 = as_function(AS, z[0]^p - z[0])" ] }, { "cell_type": "code", "execution_count": 6, "metadata": { "collapsed": false, "scrolled": true }, "outputs": [ ], "source": [ "zdual = dual_elt(AS, zmag)" ] }, { "cell_type": "code", "execution_count": 8, "metadata": { "collapsed": false, "scrolled": true }, "outputs": [ { "data": { "text/plain": [ "1" ] }, "execution_count": 8, "metadata": { }, "output_type": "execute_result" } ], "source": [ "(zmag*(zdual.group_action([0, 0]))).trace2()" ] }, { "cell_type": "code", "execution_count": 5, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ "n = 2\n", "F = GF(5)\n", "variable_names = 'x, y'\n", "for i in range(n):\n", " variable_names += ', z' + str(i)\n", "Rxyz = PolynomialRing(F, n+2, variable_names)\n", "x, y = Rxyz.gens()[:2]\n", "z = Rxyz.gens()[2:]" ] }, { "cell_type": "code", "execution_count": 9, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ "def ith_component(omega, zvee, i):\n", " AS = omega.curve\n", " p = AS.characteristic\n", " group_elts = [(j1, j2) for j1 in range(p) for j2 in range(p)]\n", " z_vee_fct = zvee.group_action(group_elts[i]).function\n", " new_form = as_form(AS, z_vee_fct*omega.form)\n", " return new_form.trace2()" ] }, { "cell_type": "code", "execution_count": 12, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ "om = AS.holomorphic_differentials_basis()[4]" ] }, { "cell_type": "code", "execution_count": 13, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "(-z0^3*z1^2 - 2*x*z0^2*z1 + z1^4 + 2*x^2*z0) * dx" ] }, "execution_count": 13, "metadata": { }, "output_type": "execute_result" } ], "source": [ "om" ] }, { "cell_type": "code", "execution_count": 14, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ "def combination_components(omega, zmag, w):\n", " AS = omega.curve\n", " p = AS.characteristic\n", " group_elts = [(j1, j2) for j1 in range(p) for j2 in range(p)]\n", " zvee = dual_elt(AS, zmag)\n", " result = as_form(AS, 0)\n", " for i in range(p^2):\n", " omegai = ith_component(omega, zvee, i)\n", " aux_fct1 = w.group_action(group_elts[i]).function\n", " aux_fct2 = omegai.form\n", " result += as_form(AS, aux_fct1*aux_fct2)\n", " return result" ] }, { "cell_type": "code", "execution_count": 15, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "(-z0^3*z1^2 - 2*x*z0^2*z1 + z1^4 + 2*x^2*z0) * dx" ] }, "execution_count": 15, "metadata": { }, "output_type": "execute_result" } ], "source": [ "combination_components(om, zmag, zmag)" ] }, { "cell_type": "code", "execution_count": 16, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "(-z0^3*z1^2 - 2*x*z0^2*z1 + z1^4 + 2*x^2*z0) * dx" ] }, "execution_count": 16, "metadata": { }, "output_type": "execute_result" } ], "source": [ "om" ] }, { "cell_type": "code", "execution_count": 150, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "t^-68 + t^-48 + t^-40 + 2*t^-36 + 3*t^-28 + 3*t^-20 + 2*t^-16 + t^-12 + 3*t^-8 + 4*t^-4 + 4*t^4 + t^8 + 3*t^12 + 4*t^20 + t^24 + 3*t^28 + 3*t^32 + 4*t^36 + t^44 + 3*t^48 + t^52 + 2*t^56 + 4*t^60 + 3*t^64 + 2*t^72 + 3*t^76 + 2*t^88 + 2*t^92 + t^96 + O(t^100)" ] }, "execution_count": 150, "metadata": { }, "output_type": "execute_result" } ], "source": [ "zmag.expansion_at_infty()" ] }, { "cell_type": "code", "execution_count": 52, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "O(t^430)" ] }, "execution_count": 52, "metadata": { }, "output_type": "execute_result" } ], "source": [ "w1 = as_function(AS, z[0]^2*y^3/x - z[0]^2*y^2)\n", "w1.expansion_at_infty()" ] }, { "cell_type": "code", "execution_count": 48, "metadata": { "collapsed": false }, "outputs": [ { "data": { "text/plain": [ "4*t^-85 + t^-70 + t^-65 + t^-57 + 3*t^-50 + t^-45 + 3*t^-42 + t^-37 + t^-30 + 3*t^-25 + 3*t^-17 + t^-14 + t^-9 + t^-5 + 2*t^-2 + t^10 + t^11 + 3*t^15 + t^19 + 2*t^23 + 3*t^26 + 2*t^30 + t^38 + 2*t^43 + t^47 + 3*t^51 + 4*t^54 + t^55 + t^58 + t^66 + t^70 + t^75 + 4*t^79 + 4*t^83 + 2*t^86 + 2*t^90 + t^94 + 3*t^95 + 2*t^98 + 4*t^103 + 3*t^110 + 4*t^115 + t^122 + 4*t^123 + 3*t^126 + t^130 + 4*t^131 + 3*t^135 + 2*t^138 + 2*t^143 + 3*t^150 + 3*t^151 + 4*t^154 + 4*t^155 + 3*t^158 + 4*t^159 + t^166 + t^170 + t^175 + 3*t^178 + t^186 + 4*t^187 + 2*t^190 + 2*t^191 + 3*t^195 + t^198 + 3*t^206 + t^211 + 3*t^215 + 3*t^218 + 4*t^219 + 4*t^222 + t^223 + 4*t^230 + 3*t^234 + 4*t^235 + t^238 + t^247 + 4*t^250 + 4*t^251 + 3*t^254 + 2*t^258 + 3*t^262 + t^263 + 2*t^266 + t^270 + 2*t^275 + 2*t^279 + 2*t^283 + 4*t^286 + t^290 + 4*t^291 + 3*t^294 + 3*t^295 + t^298 + 4*t^303 + 3*t^310 + 2*t^311 + 2*t^315 + t^318 + 3*t^319 + 4*t^322 + 4*t^323 + 4*t^330 + 2*t^331 + 3*t^335 + t^338 + 4*t^343 + t^346 + 3*t^347 + 2*t^350 + t^351 + 3*t^354 + 3*t^355 + t^358 + t^363 + 4*t^370 + t^374 + 2*t^379 + 4*t^383 + 3*t^387 + 4*t^390 + 4*t^391 + 2*t^398 + t^402 + t^406 + t^410 + 3*t^411 + O(t^415)" ] }, "execution_count": 48, "metadata": { }, "output_type": "execute_result" } ], "source": [ ] }, { "cell_type": "code", "execution_count": 40, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ "om1 = combination_components(om, zmag, w1)" ] }, { "cell_type": "code", "execution_count": 43, "metadata": { "collapsed": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "(0) * dx\n", "O(t^474)\n" ] } ], "source": [ "print(om1)\n", "print(om1.expansion_at_infty())" ] }, { "cell_type": "code", "execution_count": 0, "metadata": { "collapsed": false }, "outputs": [ ], "source": [ ] } ], "metadata": { "kernelspec": { "display_name": "SageMath 9.7", "language": "sagemath", "metadata": { "cocalc": { "description": "Open-source mathematical software system", "priority": 10, "url": "https://www.sagemath.org/" } }, "name": "sage-9.7", "resource_dir": "/ext/jupyter/kernels/sage-9.7" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.10.5" } }, "nbformat": 4, "nbformat_minor": 4 }