DeRhamComputation/elementary_covers_of_supere...

83 KiB

class superelliptic:
    """Class of a superelliptic curve. Given a polynomial f(x) with coefficient field F, it constructs
       the curve y^m = f(x)"""
    def __init__(self, f, m):
        Rx = f.parent()
        x = Rx.gen()
        F = Rx.base()        
        Rx.<x> = PolynomialRing(F)
        Rxy.<x, y> = PolynomialRing(F, 2)
        Fxy = FractionField(Rxy)
        self.polynomial = Rx(f)
        self.exponent = m
        self.base_ring = F
        self.characteristic = F.characteristic()
        r = Rx(f).degree()
        delta = GCD(r, m)
        
    def __repr__(self):
        f = self.polynomial
        m = self.exponent
        F = self.base_ring
        return 'Superelliptic curve with the equation y^' + str(m) + ' = ' + str(f)+' over ' + str(F)

    #Auxilliary algorithm that returns the basis of holomorphic differentials
    #of the curve and (as a second argument) the list of pairs (i, j)
    #such that x^i dx/y^j is holomorphic.
    def basis_holomorphic_differentials_degree(self):
        f = self.polynomial
        m = self.exponent
        r = f.degree()
        delta = GCD(r, m)
        F = self.base_ring
        Rx.<x> = PolynomialRing(F)
        Rxy.<x, y> = PolynomialRing(F, 2)
        Fxy = FractionField(Rxy)
        #########basis of holomorphic differentials and de Rham

        basis_holo = []
        degrees0 = {}
        k = 0

        for j in range(1, m):
            for i in range(1, r):
                if (r*j - m*i >= delta):
                    basis_holo += [superelliptic_form(self, Fxy(x^(i-1)/y^j))]
                    degrees0[k] = (i-1, j)
                    k = k+1

        return(basis_holo, degrees0)
    
    #Returns the basis of holomorphic differentials using the previous algorithm.
    def holomorphic_differentials_basis(self):
        basis_holo, degrees0 = self.basis_holomorphic_differentials_degree()
        return basis_holo
    #Returns the list of pairs (i, j) such that x^i dx/y^j is holomorphic.
    def degrees_holomorphic_differentials(self):
        basis_holo, degrees0 = self.basis_holomorphic_differentials_degree()
        return degrees0

    def basis_de_rham_degrees(self):
        f = self.polynomial
        m = self.exponent
        r = f.degree()
        delta = GCD(r, m)
        F = self.base_ring
        Rx.<x> = PolynomialRing(F)
        Rxy.<x, y> = PolynomialRing(F, 2)
        Fxy = FractionField(Rxy)
        basis_holo = self.holomorphic_differentials_basis()
        basis = []
        #First g_X elements of basis are holomorphic differentials.
        for k in range(0, len(basis_holo)):
            basis += [superelliptic_cech(self, basis_holo[k], superelliptic_function(self, 0))]

        ## Next elements do not come from holomorphic differentials.
        t = len(basis)
        degrees0 = {}
        degrees1 = {}
        for j in range(1, m):
            for i in range(1, r):
                if (r*(m-j) - m*i >= delta): 
                    s = Rx(m-j)*Rx(x)*Rx(f.derivative()) - Rx(m)*Rx(i)*f
                    psi = Rx(cut(s, i))
                    basis += [superelliptic_cech(self, superelliptic_form(self, Fxy(psi/y^j)), superelliptic_function(self, Fxy(m*y^(m-j)/x^i)))]
                    degrees0[t] = (psi.degree(), j)
                    degrees1[t] = (-i, m-j)
                    t += 1
        return basis, degrees0, degrees1

    def de_rham_basis(self):
        basis, degrees0, degrees1 = self.basis_de_rham_degrees()
        return basis

    def degrees_de_rham0(self):
        basis, degrees0, degrees1 = self.basis_de_rham_degrees()
        return degrees0

    def degrees_de_rham1(self):
        basis, degrees0, degrees1 = self.basis_de_rham_degrees()
        return degrees1    
    
    def is_smooth(self):
        f = self.polynomial
        if f.discriminant() == 0:
            return 0
        return 1
    
    def genus(self):
        r = self.polynomial.degree()
        m = self.exponent
        delta = GCD(r, m)
        return 1/2*((r-1)*(m-1) - delta + 1)
    
    def verschiebung_matrix(self):
        basis = self.de_rham_basis()
        g = self.genus()
        p = self.characteristic
        F = self.base_ring
        M = matrix(F, 2*g, 2*g)
        for i in range(0, len(basis)):
            w = basis[i]
            v = w.verschiebung().coordinates()
            M[i, :] = v
        return M
    
    def frobenius_matrix(self):
        basis = self.de_rham_basis()
        g = self.genus()
        p = self.characteristic
        F = self.base_ring
        M = matrix(F, 2*g, 2*g)
        
        for i in range(0, len(basis)):
            w = basis[i]
            v = w.frobenius().coordinates()
            M[i, :] = v
        return M

    def cartier_matrix(self):
        basis = self.holomorphic_differentials_basis()
        g = self.genus()
        p = self.characteristic
        F = self.base_ring
        M = matrix(F, g, g)
        for i in range(0, len(basis)):
            w = basis[i]
            v = w.cartier().coordinates()
            M[i, :] = v
        return M     

#    def p_rank(self):
#        return self.cartier_matrix().rank()
    
    def a_number(self):
        g = C.genus()
        return g - self.cartier_matrix().rank()
    
    def final_type(self, test = 0):
        Fr = self.frobenius_matrix()
        V = self.verschiebung_matrix()
        p = self.characteristic
        return flag(Fr, V, p, test)

#Auxilliary. Given a superelliptic curve C : y^m = f(x) and a polynomial g(x, y)
#it replaces repeteadly all y^m's in g(x, y) by f(x). As a result
#you obtain \sum_{i = 0}^{m-1} y^i g_i(x).
def reduction(C, g):
    p = C.characteristic
    F = C.base_ring
    Rxy.<x, y> = PolynomialRing(F, 2)
    Fxy = FractionField(Rxy)
    f = C.polynomial
    r = f.degree()
    m = C.exponent
    g = Fxy(g)
    g1 = g.numerator()
    g2 = g.denominator()
    
    Rx.<x> = PolynomialRing(F)
    Fx = FractionField(Rx)
    FxRy.<y> = PolynomialRing(Fx)    
    (A, B, C) = xgcd(FxRy(g2), FxRy(y^m - f))
    g = FxRy(g1*B/A)
    
    while(g.degree(Rxy(y)) >= m):
        d = g.degree(Rxy(y))
        G = coff(g, d)
        i = floor(d/m)
        g = g - G*y^d + f^i * y^(d%m) *G
    
    return(FxRy(g))

#Auxilliary. Given a superelliptic curve C : y^m = f(x) and a polynomial g(x, y)
#it replaces repeteadly all y^m's in g(x, y) by f(x). As a result
#you obtain \sum_{i = 0}^{m-1} g_i(x)/y^i. This is needed for reduction of
#superelliptic forms.
def reduction_form(C, g):
    F = C.base_ring
    Rxy.<x, y> = PolynomialRing(F, 2)
    Fxy = FractionField(Rxy)
    f = C.polynomial
    r = f.degree()
    m = C.exponent
    g = reduction(C, g)

    g1 = Rxy(0)
    Rx.<x> = PolynomialRing(F)
    Fx = FractionField(Rx)
    FxRy.<y> = PolynomialRing(Fx)
    
    g = FxRy(g)
    for j in range(0, m):
        if j==0:
            G = coff(g, 0)
            g1 += FxRy(G)
        else:
            G = coff(g, j)
            g1 += Fxy(y^(j-m)*f*G)
    return(g1)

#Class of rational functions on a superelliptic curve C. g = g(x, y) is a polynomial
#defining the function.
class superelliptic_function:
    def __init__(self, C, g):
        F = C.base_ring
        Rxy.<x, y> = PolynomialRing(F, 2)
        Fxy = FractionField(Rxy)
        f = C.polynomial
        r = f.degree()
        m = C.exponent
        
        self.curve = C
        g = reduction(C, g)
        self.function = g
        
    def __repr__(self):
        return str(self.function)
    
    def jth_component(self, j):
        g = self.function
        C = self.curve
        F = C.base_ring
        Rx.<x> = PolynomialRing(F)
        Fx.<x> = FractionField(Rx)
        FxRy.<y> = PolynomialRing(Fx)
        g = FxRy(g)
        return coff(g, j)
    
    def __add__(self, other):
        C = self.curve
        g1 = self.function
        g2 = other.function
        g = reduction(C, g1 + g2)
        return superelliptic_function(C, g)
    
    def __sub__(self, other):
        C = self.curve
        g1 = self.function
        g2 = other.function
        g = reduction(C, g1 - g2)
        return superelliptic_function(C, g)
    
    def __mul__(self, other):
        C = self.curve
        g1 = self.function
        g2 = other.function
        g = reduction(C, g1 * g2)
        return superelliptic_function(C, g)
    
    def __truediv__(self, other):
        C = self.curve
        g1 = self.function
        g2 = other.function
        g = reduction(C, g1 / g2)
        return superelliptic_function(C, g)

    def __pow__(self, exp):
        C = self.curve
        g = self.function
        return superelliptic_function(C, g^(exp))

    def diffn(self):
        C = self.curve
        f = C.polynomial
        m = C.exponent
        F = C.base_ring
        g = self.function
        Rxy.<x, y> = PolynomialRing(F, 2)
        Fxy = FractionField(Rxy)
        g = Fxy(g)
        A = g.derivative(x)
        B = g.derivative(y)*f.derivative(x)/(m*y^(m-1))
        return superelliptic_form(C, A+B)

    
    def expansion_at_infty(self, i = 0, prec=10):
        C = self.curve
        f = C.polynomial
        m = C.exponent
        F = C.base_ring
        Rx.<x> = PolynomialRing(F)
        f = Rx(f)
        Rt.<t> = LaurentSeriesRing(F, default_prec=prec)
        RptW.<W> = PolynomialRing(Rt)
        RptWQ = FractionField(RptW)
        Rxy.<x, y> = PolynomialRing(F)
        RxyQ = FractionField(Rxy)
        fct = self.function
        fct = RxyQ(fct)
        r = f.degree()
        delta, a, b = xgcd(m, r)
        b = -b
        M = m/delta
        R = r/delta
        while a<0:
            a += R
            b += M
        
        g = (x^r*f(x = 1/x))
        gW = RptWQ(g(x = t^M * W^b)) - W^(delta)
        ww = naive_hensel(gW, F, start = root_of_unity(F, delta)^i, prec = prec)
        xx = Rt(1/(t^M*ww^b))
        yy = 1/(t^R*ww^a)
        return Rt(fct(x = Rt(xx), y = Rt(yy)))
    
def naive_hensel(fct, F, start = 1, prec=10):
    Rt.<t> = LaurentSeriesRing(F, default_prec=prec)
    RtQ = FractionField(Rt)
    RptW.<W> = PolynomialRing(RtQ)
    RptWQ = FractionField(RptW)
    fct = RptWQ(fct)
    fct = RptW(numerator(fct))
    #return(fct)
    #while fct not in RptW:
    #    print(fct)
    #    fct *= W
    alpha = (fct.derivative())(W = start)
    w0 = Rt(start)
    i = 1
    while(i < prec):
        w0 = w0 - fct(W = w0)/alpha + O(t^(prec))
        i += 1
    return w0

class superelliptic_form:
    def __init__(self, C, g):
        F = C.base_ring
        Rxy.<x, y> = PolynomialRing(F, 2)
        Fxy = FractionField(Rxy)
        g = Fxy(reduction_form(C, g))
        self.form = g
        self.curve = C
        
    def __add__(self, other):
        C = self.curve
        g1 = self.form
        g2 = other.form
        g = reduction(C, g1 + g2)
        return superelliptic_form(C, g)
    
    def __sub__(self, other):
        C = self.curve
        g1 = self.form
        g2 = other.form
        g = reduction(C, g1 - g2)
        return superelliptic_form(C, g)
    
    def __repr__(self):
        g = self.form
        if len(str(g)) == 1:
            return str(g) + ' dx'
        return '('+str(g) + ') dx'

    def __rmul__(self, constant):
        C = self.curve
        omega = self.form
        return superelliptic_form(C, constant*omega)        
    
    def cartier(self):
        C = self.curve
        m = C.exponent
        p = C.characteristic
        f = C.polynomial
        F = C.base_ring
        Rx.<x> = PolynomialRing(F)
        Fx = FractionField(Rx)
        FxRy.<y> = PolynomialRing(Fx)
        Fxy = FractionField(FxRy)
        result = superelliptic_form(C, FxRy(0))
        mult_order = Integers(m)(p).multiplicative_order()
        M = Integer((p^(mult_order)-1)/m)
        
        for j in range(1, m):
            fct_j = self.jth_component(j)
            h = Rx(fct_j*f^(M*j))
            j1 = (p^(mult_order-1)*j)%m
            B = floor(p^(mult_order-1)*j/m)
            result += superelliptic_form(C, polynomial_part(p, h)/(f^B*y^(j1)))
        return result   
    
    def coordinates(self):
        C = self.curve
        F = C.base_ring
        m = C.exponent
        Rx.<x> = PolynomialRing(F)
        Fx = FractionField(Rx)
        FxRy.<y> = PolynomialRing(Fx)
        g = C.genus()
        degrees_holo = C.degrees_holomorphic_differentials()
        degrees_holo_inv = {b:a for a, b in degrees_holo.items()}
        basis = C.holomorphic_differentials_basis()
        
        for j in range(1, m):
            omega_j = Fx(self.jth_component(j))
            if omega_j != Fx(0):
                d = degree_of_rational_fctn(omega_j, F)
                index = degrees_holo_inv[(d, j)]
                a = coeff_of_rational_fctn(omega_j, F)
                a1 = coeff_of_rational_fctn(basis[index].jth_component(j), F)
                elt = self - (a/a1)*basis[index]
                return elt.coordinates() + a/a1*vector([F(i == index) for i in range(0, g)])
            
        return vector(g*[0])
    
    def jth_component(self, j):
        g = self.form
        C = self.curve
        F = C.base_ring
        Rx.<x> = PolynomialRing(F)
        Fx = FractionField(Rx)
        FxRy.<y> = PolynomialRing(Fx)
        Fxy = FractionField(FxRy)
        Ryinv.<y_inv> = PolynomialRing(Fx)
        g = Fxy(g)
        g = g(y = 1/y_inv)
        g = Ryinv(g)
        return coff(g, j)
    
    def is_regular_on_U0(self):
        C = self.curve
        F = C.base_ring
        m = C.exponent
        Rx.<x> = PolynomialRing(F)
        for j in range(1, m):
            if self.jth_component(j) not in Rx:
                return 0
        return 1
    
    def is_regular_on_Uinfty(self):
        C = self.curve
        F = C.base_ring
        m = C.exponent
        f = C.polynomial
        r = f.degree()
        delta = GCD(m, r)
        M = m/delta
        R = r/delta
        
        for j in range(1, m):
            A = self.jth_component(j)
            d = degree_of_rational_fctn(A, F)
            if(-d*M + j*R -(M+1)<0):
                return 0
        return 1

    def expansion_at_infty(self, i = 0, prec=10):
        g = self.form
        C = self.curve
        g = superelliptic_function(C, g)
        g = g.expansion_at_infty(i = i, prec=prec)
        x_series = superelliptic_function(C, x).expansion_at_infty(i = i, prec=prec)
        dx_series = x_series.derivative()
        return g*dx_series
            
class superelliptic_cech:
    def __init__(self, C, omega, fct):
        self.omega0 = omega
        self.omega8 = omega - fct.diffn()
        self.f = fct
        self.curve = C
        
    def __add__(self, other):
        C = self.curve
        return superelliptic_cech(C, self.omega0 + other.omega0, self.f + other.f)
    
    def __sub__(self, other):
        C = self.curve
        return superelliptic_cech(C, self.omega0 - other.omega0, self.f - other.f)

    def __rmul__(self, constant):
        C = self.curve
        w1 = self.omega0.form
        f1 = self.f.function
        w2 = superelliptic_form(C, constant*w1)
        f2 = superelliptic_function(C, constant*f1)
        return superelliptic_cech(C, w2, f2)    
    
    def __repr__(self):
        return "(" + str(self.omega0) + ", " + str(self.f) + ", " + str(self.omega8) + ")" 
    
    def verschiebung(self):
        C = self.curve
        omega = self.omega0
        F = C.base_ring
        Rx.<x> = PolynomialRing(F)
        return superelliptic_cech(C, omega.cartier(), superelliptic_function(C, Rx(0)))
    
    def frobenius(self):
        C = self.curve
        fct = self.f.function
        p = C.characteristic
        Rx.<x> = PolynomialRing(F)
        return superelliptic_cech(C, superelliptic_form(C, Rx(0)), superelliptic_function(C, fct^p))

    def coordinates(self):
        C = self.curve
        F = C.base_ring
        m = C.exponent
        Rx.<x> = PolynomialRing(F)
        Fx = FractionField(Rx)
        FxRy.<y> = PolynomialRing(Fx)
        g = C.genus()
        degrees_holo = C.degrees_holomorphic_differentials()
        degrees_holo_inv = {b:a for a, b in degrees_holo.items()}
        degrees0 = C.degrees_de_rham0()
        degrees0_inv = {b:a for a, b in degrees0.items()}
        degrees1 = C.degrees_de_rham1()
        degrees1_inv = {b:a for a, b in degrees1.items()}
        basis = C.de_rham_basis()
        
        omega = self.omega0
        fct = self.f
        
        if fct.function == Rx(0) and omega.form != Rx(0):
            for j in range(1, m):
                omega_j = Fx(omega.jth_component(j))
                if omega_j != Fx(0):
                    d = degree_of_rational_fctn(omega_j, F)
                    index = degrees_holo_inv[(d, j)]
                    a = coeff_of_rational_fctn(omega_j, F)
                    a1 = coeff_of_rational_fctn(basis[index].omega0.jth_component(j), F)
                    elt = self - (a/a1)*basis[index]
                    return elt.coordinates() + a/a1*vector([F(i == index) for i in range(0, 2*g)])
                    
        for j in range(1, m):
            fct_j = Fx(fct.jth_component(j))
            if (fct_j != Rx(0)):
                d = degree_of_rational_fctn(fct_j, p)
            
                if (d, j) in degrees1.values():
                    index = degrees1_inv[(d, j)]
                    a = coeff_of_rational_fctn(fct_j, F)
                    elt = self - (a/m)*basis[index]
                    return elt.coordinates() + a/m*vector([F(i == index) for i in range(0, 2*g)])
                
                if d<0:
                    a = coeff_of_rational_fctn(fct_j, F)
                    h = superelliptic_function(C, FxRy(a*y^j*x^d))
                    elt = superelliptic_cech(C, self.omega0, self.f - h)
                    return elt.coordinates()
            
                if (fct_j != Rx(0)):
                    G = superelliptic_function(C, y^j*x^d)
                    a = coeff_of_rational_fctn(fct_j, F)
                    elt =self - a*superelliptic_cech(C, diffn(G), G)
                    return elt.coordinates()

        return vector(2*g*[0])
    
    def is_cocycle(self):
        w0 = self.omega0
        w8 = self.omega8
        fct = self.f
        if not w0.is_regular_on_U0() and not w8.is_regular_on_Uinfty():
            return('w0 & w8')
        if not w0.is_regular_on_U0():
            return('w0')
        if not w8.is_regular_on_Uinfty():
            return('w8')
        if w0.is_regular_on_U0() and w8.is_regular_on_Uinfty():
            return 1
        return 0

#Auxilliary. If f = f1/f2 is a rational function, return deg f_1 - deg f_2.
def degree_of_rational_fctn(f, F):
    Rx.<x> = PolynomialRing(F)
    Fx = FractionField(Rx)
    f = Fx(f)
    f1 = f.numerator()
    f2 = f.denominator()
    d1 = f1.degree()
    d2 = f2.degree()
    return(d1 - d2)

#Auxilliary. If f = f1/f2 is a rational function, return (leading coeff of f1)/(leading coeff of f2).
def coeff_of_rational_fctn(f, F):
    Rx.<x> = PolynomialRing(F)
    Fx = FractionField(Rx)
    f = Fx(f)
    if f == Rx(0):
        return 0
    f1 = f.numerator()
    f2 = f.denominator()
    d1 = f1.degree()
    d2 = f2.degree()
    a1 = f1.coefficients(sparse = false)[d1]
    a2 = f2.coefficients(sparse = false)[d2]
    return(a1/a2)

#Auxilliary. Given polynomial f(x) and integer d, return
#coefficient of x^d in f (and 0 is deg(f)<d).
def coff(f, d):
    lista = f.coefficients(sparse = false)
    if len(lista) <= d:
        return 0
    return lista[d]

#Auxilliary. Given polynomial f(x) = \sum_i a_i x^i and integer i, return
#only \sum_{j >= i+1} a_j x^{j - i -1}
def cut(f, i):
    R = f.parent()
    coeff = f.coefficients(sparse = false)
    return sum(R(x^(j-i-1)) * coeff[j] for j in range(i+1, f.degree() + 1))

def polynomial_part(p, h):
    F = GF(p)
    Rx.<x> = PolynomialRing(F)
    h = Rx(h)
    result = Rx(0)
    for i in range(0, h.degree()+1):
        if (i%p) == p-1:
            power = Integer((i-(p-1))/p)
            result += Integer(h[i]) * x^(power)    
    return result

#Find delta-th root of unity in field F
def root_of_unity(F, delta):
    Rx.<x> = PolynomialRing(F)
    cyclotomic = x^(delta) - 1
    for root, a in cyclotomic.roots():
        powers = [root^d for d in delta.divisors() if d!= delta]
        if 1 not in powers:
            return root
        
def preimage(U, V, M): #preimage of subspace U under M
    basis_preimage = M.right_kernel().basis()
    imageU = U.intersection(M.transpose().image())
    basis = imageU.basis()
    for v in basis:
        w = M.solve_right(v)
        basis_preimage = basis_preimage + [w]
    return V.subspace(basis_preimage)

def image(U, V, M):
    basis = U.basis()
    basis_image = []
    for v in basis:
        basis_image += [M*v]
    return V.subspace(basis_image)

def flag(F, V, p, test = 0):
    dim = F.dimensions()[0]
    space = VectorSpace(GF(p), dim)
    flag_subspaces = (dim+1)*[0]
    flag_used = (dim+1)*[0]
    final_type = (dim+1)*['?']
    
    flag_subspaces[dim] = space
    flag_used[dim] = 1
   
    
    while 1 in flag_used:
        index = flag_used.index(1)
        flag_used[index] = 0
        U = flag_subspaces[index]
        U_im = image(U, space, V)
        d_im = U_im.dimension()
        final_type[index] = d_im
        U_pre = preimage(U, space, F)
        d_pre = U_pre.dimension()
        
        if flag_subspaces[d_im] == 0:
            flag_subspaces[d_im] = U_im
            flag_used[d_im] = 1
        
        if flag_subspaces[d_pre] == 0:
            flag_subspaces[d_pre] = U_pre
            flag_used[d_pre] = 1
    
    if test == 1:
        print('(', final_type, ')')
    
    for i in range(0, dim+1):
        if final_type[i] == '?' and final_type[dim - i] != '?':
            i1 = dim - i
            final_type[i] = final_type[i1] - i1 + dim/2
    
    final_type[0] = 0
    
    for i in range(1, dim+1):
        if final_type[i] == '?':
            prev = final_type[i-1]
            if prev != '?' and prev in final_type[i+1:]:
                final_type[i] = prev
                
    for i in range(1, dim+1):
        if final_type[i] == '?':
            final_type[i] = min(final_type[i-1] + 1, dim/2)
    
    if is_final(final_type, dim/2):
        return final_type[1:dim/2 + 1]
    print('error:', final_type[1:dim/2 + 1])
    
def is_final(final_type, dim):
    n = len(final_type)
    if final_type[0] != 0:
        return 0
    
    if final_type[n-1] != dim:
        return 0
    
    for i in range(1, n):
        if final_type[i] != final_type[i - 1] and final_type[i] != final_type[i - 1] + 1:
            return 0
    return 1
class as_cover:
    def __init__(self, C, list_of_fcts, prec = 10):
        self.quotient = C
        self.functions = list_of_fcts
        self.height = len(list_of_fcts)
        F = C.base_ring
        self.base_ring = F
        p = C.characteristic
        self.characteristic = p
        self.prec = prec
        f = C.polynomial
        m = C.exponent
        r = f.degree()
        delta = GCD(m, r)
        self.nb_of_pts_at_infty = delta
        Rxy.<x, y> = PolynomialRing(F, 2)
        Rt.<t> = LaurentSeriesRing(F, default_prec=prec)

        all_x_series = []
        all_y_series = []
        all_z_series = []
        all_dx_series = []
        all_jumps = []

        for i in range(delta):
            x_series = superelliptic_function(C, x).expansion_at_infty(i = i, prec=prec)
            y_series = superelliptic_function(C, y).expansion_at_infty(i = i, prec=prec)
            z_series = []
            jumps = []
            n = len(list_of_fcts)
            list_of_power_series = [g.expansion_at_infty(i = i, prec=prec) for g in list_of_fcts]
            for i in range(n):
                power_series = list_of_power_series[i]
                jump, correction, t_old, z = artin_schreier_transform(power_series, prec = prec)
                x_series = x_series(t = t_old)
                y_series = y_series(t = t_old)
                z_series = [zi(t = t_old) for zi in z_series]
                z_series += [z]
                jumps += [jump]
                list_of_power_series = [g(t = t_old) for g in list_of_power_series]
            
            all_jumps += [jumps]
            all_x_series += [x_series]
            all_y_series += [y_series]
            all_z_series += [z_series]
            all_dx_series += [x_series.derivative()]
        self.jumps = all_jumps
        self.x = all_x_series
        self.y = all_y_series
        self.z = all_z_series
        self.dx = all_dx_series
        
    def __repr__(self):
        n = self.height
        p = self.characteristic
        if n==1:
            return "(Z/p)-cover of " + str(self.quotient)+" with the equation:\n z^" + str(p) + " - z = " + str(self.functions[0])
        
        result = "(Z/p)^"+str(self.height)+ "-cover of " + str(self.quotient)+" with the equations:\n"
        for i in range(n):
            result += 'z' + str(i) + "^" + str(p) + " - z" + str(i) + " = " + str(self.functions[i]) + "\n"
        return result
        
    def genus(self):
        jumps = self.jumps
        gY = self.quotient.genus()
        n = self.height
        delta = self.nb_of_pts_at_infty
        p = self.characteristic
        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))
    
    def exponent_of_different(self, i = 0):
        jumps = self.jumps
        n = self.height
        delta = self.nb_of_pts_at_infty
        p = self.characteristic
        return sum(p^(n-j-1)*(jumps[i][j]+1)*(p-1) for j in range(n))

    def exponent_of_different_prim(self, i = 0):
        jumps = self.jumps
        n = self.height
        delta = self.nb_of_pts_at_infty
        p = self.characteristic
        return sum(p^(n-j-1)*(jumps[i][j])*(p-1) for j in range(n))
    
    def holomorphic_differentials_basis(self, threshold = 8):
        from itertools import product
        x_series = self.x
        y_series = self.y
        z_series = self.z
        dx_series = self.dx
        delta = self.nb_of_pts_at_infty
        p = self.characteristic
        n = self.height
        prec = self.prec
        C = self.quotient
        F = self.base_ring
        m = C.exponent
        r = C.polynomial.degree()
        variable_names = 'x, y'
        for i in range(n):
            variable_names += ', z' + str(i)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        RxyzQ = FractionField(Rxyz)
        Rt.<t> = LaurentSeriesRing(F, default_prec=prec)
        #Tworzymy zbiór S form z^i x^j y^k dx/y o waluacji >= waluacja z^(p-1)*dx/y
        S = []
        RQxyz = FractionField(Rxyz)
        pr = [list(GF(p)) for _ in range(n)]
        for i in range(0, threshold*r):
            for j in range(0, m):
                for k in product(*pr):
                    eta = as_form(self, x^i * prod(z[i1]^(k[i1]) for i1 in range(n))/y^j)
                    eta_exp = eta.expansion_at_infty()
                    S += [(eta, eta_exp)]

        forms = holomorphic_combinations(S)
        
        for i in range(1, delta):
            forms = [(omega, omega.expansion_at_infty(i = i)) for omega in forms]
            forms = holomorphic_combinations(forms)
        
        if len(forms) < self.genus():
            print("I haven't found all forms.")
            return holomorphic_differentials_basis(self, threshold = threshold + 1)
        if len(forms) > self.genus():
            print("Increase precision.")
        return forms

    def at_most_poles(self, pole_order, threshold = 8):
        """ Find fcts with pole order in infty's at most pole_order. Threshold gives a bound on powers of x in the function.
           If you suspect that you haven't found all the functions, you may increase it."""
        from itertools import product
        x_series = self.x
        y_series = self.y
        z_series = self.z
        delta = self.nb_of_pts_at_infty
        p = self.characteristic
        n = self.height
        prec = self.prec
        C = self.quotient
        F = self.base_ring
        m = C.exponent
        r = C.polynomial.degree()
        variable_names = 'x, y'
        for i in range(n):
            variable_names += ', z' + str(i)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        RxyzQ = FractionField(Rxyz)
        Rt.<t> = LaurentSeriesRing(F, default_prec=prec)
        #Tworzymy zbiór S form z^i x^j y^k dx/y o waluacji >= waluacja z^(p-1)*dx/y
        S = []
        RQxyz = FractionField(Rxyz)
        pr = [list(GF(p)) for _ in range(n)]
        for i in range(0, threshold*r):
            for j in range(0, m):
                for k in product(*pr):
                    eta = as_function(self, x^i * prod(z[i1]^(k[i1]) for i1 in range(n))*y^j)
                    eta_exp = eta.expansion_at_infty()
                    S += [(eta, eta_exp)]

        forms = holomorphic_combinations_fcts(S, pole_order)

        for i in range(1, delta):
            forms = [(omega, omega.expansion_at_infty(i = i)) for omega in forms]
            forms = holomorphic_combinations_fcts(forms, pole_order)
                  
        return forms
    
    def magical_element(self, threshold = 8):
        list_of_elts = self.at_most_poles(self.exponent_of_different_prim(), threshold)
        result = []
        for a in list_of_elts:
            if a.trace().function != 0:
                result += [a]
        return result

    def pseudo_magical_element(self, threshold = 8):
        list_of_elts = self.at_most_poles(self.exponent_of_different(), threshold)
        result = []
        for a in list_of_elts:
            if a.trace().function != 0:
                result += [a]
        return result

    def at_most_poles_forms(self, pole_order, threshold = 8):
        """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.
           If you suspect that you haven't found all the functions, you may increase it."""
        from itertools import product
        x_series = self.x
        y_series = self.y
        z_series = self.z
        delta = self.nb_of_pts_at_infty
        p = self.characteristic
        n = self.height
        prec = self.prec
        C = self.quotient
        F = self.base_ring
        m = C.exponent
        r = C.polynomial.degree()
        variable_names = 'x, y'
        for i in range(n):
            variable_names += ', z' + str(i)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        RxyzQ = FractionField(Rxyz)
        Rt.<t> = LaurentSeriesRing(F, default_prec=prec)
        #Tworzymy zbiór S form z^i x^j y^k dx/y o waluacji >= waluacja z^(p-1)*dx/y
        S = []
        RQxyz = FractionField(Rxyz)
        pr = [list(GF(p)) for _ in range(n)]
        for i in range(0, threshold*r):
            for j in range(0, m):
                for k in product(*pr):
                    eta = as_form(self, x^i * prod(z[i1]^(k[i1]) for i1 in range(n))/y^j)
                    eta_exp = eta.expansion_at_infty()
                    S += [(eta, eta_exp)]

        forms = holomorphic_combinations_forms(S, pole_order)

        for i in range(1, delta):
            forms = [(omega, omega.expansion_at_infty(i = i)) for omega in forms]
            forms = holomorphic_combinations_forms(forms, pole_order)
                  
        return forms

def holomorphic_combinations(S):
    """Given a list S of pairs (form, corresponding Laurent series at some pt), find their combinations holomorphic at that pt."""
    C_AS = S[0][0].curve
    p = C_AS.characteristic
    F = C_AS.base_ring
    prec = C_AS.prec
    Rt.<t> = LaurentSeriesRing(F, default_prec=prec)
    RtQ = FractionField(Rt)
    minimal_valuation = min([g[1].valuation() for g in S])
    if minimal_valuation >= 0:
        return [s[0] for s in S]
    list_of_lists = [] #to będzie lista złożona z list współczynników część nieholomorficznych rozwinięcia form z S
    for eta, eta_exp in S:
        a = -minimal_valuation + eta_exp.valuation()
        list_coeffs = a*[0] + eta_exp.list() + (-minimal_valuation)*[0]
        list_coeffs = list_coeffs[:-minimal_valuation]
        list_of_lists += [list_coeffs]
    M = matrix(F, list_of_lists)
    V = M.kernel() #chcemy wyzerować części nieholomorficzne, biorąc kombinacje form z S


    # Sprawdzamy, jakim formom odpowiadają elementy V.
    forms = []
    for vec in V.basis():
        forma_holo = as_form(C_AS, 0)
        forma_holo_power_series = Rt(0)
        for vec_wspolrzedna, elt_S in zip(vec, S):
            eta = elt_S[0]
            #eta_exp = elt_S[1]
            forma_holo += vec_wspolrzedna*eta
            #forma_holo_power_series += vec_wspolrzedna*eta_exp
        forms += [forma_holo]
    return forms

class as_function:
    def __init__(self, C, g):
        self.curve = C
        F = C.base_ring
        n = C.height
        variable_names = 'x, y'
        for i in range(n):
            variable_names += ', z' + str(i)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        RxyzQ = FractionField(Rxyz)
        self.function = RxyzQ(g)
        #self.function = as_reduction(AS, RxyzQ(g))
        
    def __repr__(self):
        return str(self.function)

    def __add__(self, other):
        C = self.curve
        g1 = self.function
        g2 = other.function
        return as_function(C, g1 + g2)

    def __sub__(self, other):
        C = self.curve
        g1 = self.function
        g2 = other.function
        return as_function(C, g1 - g2)

    def __rmul__(self, constant):
        C = self.curve
        g = self.function
        return as_function(C, constant*g)
    
    def __mul__(self, other):
        C = self.curve
        g1 = self.function
        g2 = other.function
        return as_function(C, g1*g2)

    def expansion_at_infty(self, i = 0):
        C = self.curve
        delta = C.nb_of_pts_at_infty
        F = C.base_ring
        x_series = C.x[i]
        y_series = C.y[i]
        z_series = C.z[i]
        n = C.height
        variable_names = 'x, y'
        for j in range(n):
            variable_names += ', z' + str(j)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        RxyzQ = FractionField(Rxyz)
        prec = C.prec
        Rt.<t> = LaurentSeriesRing(F, default_prec=prec)
        g = self.function
        g = RxyzQ(g)
        sub_list = {x : x_series, y : y_series} | {z[j] : z_series[j] for j in range(n)}
        return g.substitute(sub_list)

    def group_action(self, ZN_tuple):
        C = self.curve
        n = C.height
        F = C.base_ring
        variable_names = 'x, y'
        for j in range(n):
            variable_names += ', z' + str(j)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        RxyzQ = FractionField(Rxyz)
        sub_list = {x : x, y : y} | {z[j] : z[j]+ZN_tuple[j] for j in range(n)}
        g = self.function
        return as_function(C, g.substitute(sub_list))

    def trace(self):
        C = self.curve
        C_super = C.quotient
        n = C.height
        F = C.base_ring
        variable_names = 'x, y'
        for j in range(n):
            variable_names += ', z' + str(j)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        RxyzQ = FractionField(Rxyz)
        g = self.function
        g = as_reduction(C, g)
        result = RxyzQ(0)
        g_num = Rxyz(numerator(g))
        g_den = Rxyz(denominator(g))
        z = prod(z[i] for i in range(n))^(p-1)
        for a in g_num.monomials():
            if (z.divides(a)):
                result += g_num.monomial_coefficient(a)*a/z
        result /= g_den
        result = as_reduction(C, result)
        Rxy.<x, y> = PolynomialRing(F, 2)
        Qxy = FractionField(Rxy)
        return superelliptic_function(C_super, Qxy(result))

    def trace2(self):
        C = self.curve
        C_super = C.quotient
        n = C.height
        F = C.base_ring
        variable_names = 'x, y'
        for j in range(n):
            variable_names += ', z' + str(j)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        RxyzQ = FractionField(Rxyz)
        result = as_function(C, 0)
        for i in range(0, p):
            for j in range(0, p):
                result += self.group_action([i, j])
        result = result.function
        Rxy.<x, y> = PolynomialRing(F, 2)
        Qxy = FractionField(Rxy)
        return superelliptic_function(C_super, Qxy(result))    
    
class as_form:
    def __init__(self, C, g):
        self.curve = C
        n = C.height
        F = C.base_ring
        variable_names = 'x, y'
        for i in range(n):
            variable_names += ', z' + str(i)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        RxyzQ = FractionField(Rxyz)
        self.form = RxyzQ(g)
        
    def __repr__(self):
        return "(" + str(self.form)+") * dx"
    
    def expansion_at_infty(self, i = 0):
        C = self.curve
        delta = C.nb_of_pts_at_infty
        F = C.base_ring
        x_series = C.x[i]
        y_series = C.y[i]
        z_series = C.z[i]
        dx_series = C.dx[i]
        n = C.height
        variable_names = 'x, y'
        for j in range(n):
            variable_names += ', z' + str(j)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        RxyzQ = FractionField(Rxyz)
        prec = C.prec
        Rt.<t> = LaurentSeriesRing(F, default_prec=prec)
        g = self.form
        sub_list = {x : x_series, y : y_series} | {z[j] : z_series[j] for j in range(n)}
        return g.substitute(sub_list)*dx_series
    
    def __add__(self, other):
        C = self.curve
        g1 = self.form
        g2 = other.form
        return as_form(C, g1 + g2)
    
    def __sub__(self, other):
        C = self.curve
        g1 = self.form
        g2 = other.form
        return as_form(C, g1 - g2)
    
    def __rmul__(self, constant):
        C = self.curve
        omega = self.form
        return as_form(C, constant*omega)

    def group_action(self, ZN_tuple):
        C = self.curve
        n = C.height
        F = C.base_ring
        variable_names = 'x, y'
        for j in range(n):
            variable_names += ', z' + str(j)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        RxyzQ = FractionField(Rxyz)
        sub_list = {x : x, y : y} | {z[j] : z[j]+ZN_tuple[j] for j in range(n)}
        g = self.form
        return as_form(C, g.substitute(sub_list))

    def coordinates(self, holo):
        """Find coordinates of the given form self in terms of the basis forms in a list holo."""
        C = self.curve
        n = C.height
        gC = C.genus()
        F = C.base_ring
        variable_names = 'x, y'
        for j in range(n):
            variable_names += ', z' + str(j)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        from sage.rings.polynomial.toy_variety import linear_representation
        return linear_representation(Rxyz(self.form), holo)

    def trace(self):
        C = self.curve
        C_super = C.quotient
        n = C.height
        F = C.base_ring
        variable_names = 'x, y'
        for j in range(n):
            variable_names += ', z' + str(j)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        RxyzQ = FractionField(Rxyz)
        g = self.form
        result = RxyzQ(0)
        g_num = Rxyz(numerator(g))
        g_den = Rxyz(denominator(g))
        z = prod(z[i] for i in range(n))^(p-1)
        for a in g_num.monomials():
            if (z.divides(a)):
                result += g_num.monomial_coefficient(a)*a/z
        result /= g_den
        Rxy.<x, y> = PolynomialRing(F, 2)
        return superelliptic_form(C_super, Rxy(result))

    
    def trace2(self):
        C = self.curve
        C_super = C.quotient
        n = C.height
        F = C.base_ring
        variable_names = 'x, y'
        for j in range(n):
            variable_names += ', z' + str(j)
        Rxyz = PolynomialRing(F, n+2, variable_names)
        x, y = Rxyz.gens()[:2]
        z = Rxyz.gens()[2:]
        RxyzQ = FractionField(Rxyz)
        result = as_form(C, 0)
        for i in range(0, p):
            for j in range(0, p):
                result += self.group_action([i, j])
        result = result.form
        Rxy.<x, y> = PolynomialRing(F, 2)
        Qxy = FractionField(Rxy)
        return superelliptic_form(C_super, Qxy(result))
    
# Given power_series, find its reverse (g with g \circ power_series = id) with given precision

def new_reverse(power_series, prec = 10):
    F = power_series.parent().base()
    Rt.<t> = LaurentSeriesRing(F, default_prec=prec)
    RtQ = FractionField(Rt)
    power_series = RtQ(power_series)
    a = power_series.list()[0]
    g = 1/a*t
    n = 2
    while(n <= prec):
        aux = power_series(t = g) - t
        if aux.valuation() > n:
            b = 0
        else:
            b = aux.list()[0]
        g = g - b/a*t^n
        n += 1
    return g

def artin_schreier_transform(power_series, prec = 10):
    """Given a power_series, find correction such that power_series - (correction)^p +correction has valuation
       -jump non divisible by p. Also, express t (the variable) in terms of the uniformizer at infty on the curve
       z^p - z = power_series, where z = 1/t_new^(jump) and express z in terms of the new uniformizer."""
    correction = 0
    F = power_series.parent().base()
    Rt.<t> = LaurentSeriesRing(F, default_prec=prec)
    RtQ = FractionField(Rt)
    power_series = RtQ(power_series)
    if power_series.valuation() == +Infinity:
        return(0,0,t,0)
    while(power_series.valuation() % p == 0 and power_series.valuation() < 0):
        M = -power_series.valuation()/p
        coeff = power_series.list()[0]   #wspolczynnik a_(-p) w f_AS
        correction += coeff.nth_root(p)*t^(-M)
        power_series = power_series - (coeff*t^(-p*M) - coeff.nth_root(p)*t^(-M))
    jump = max(-(power_series.valuation()), 0)
    try:
        T = ((power_series)^(-1)).nth_root(jump) #T is defined by power_series = 1/T^m
    except:
        print("no ", str(jump), "-th root; divide by", power_series.list()[0])
        return (jump, power_series.list()[0])
    T_rev = new_reverse(T, prec = prec)
    t_old = T_rev(t^p/(1 - t^((p-1)*jump)).nth_root(jump))
    z = 1/t^(jump) + Rt(correction)(t = t_old)
    return(jump, correction, t_old, z)


def are_forms_linearly_dependent(set_of_forms):
    from sage.rings.polynomial.toy_variety import is_linearly_dependent
    C = set_of_forms[0].curve
    F = C.base_ring
    n = C.height
    variable_names = 'x, y'
    for i in range(n):
        variable_names += ', z' + str(i)
    Rxyz = PolynomialRing(F, n+2, variable_names)
    denominators = prod(denominator(omega.form) for omega in set_of_forms)
    return is_linearly_dependent([Rxyz(denominators*omega.form) for omega in set_of_forms])

def group_action_matrices(C_AS):
    F = C_AS.base_ring
    n = C_AS.height
    holo = C_AS.holomorphic_differentials_basis()
    holo_forms = [omega.form for omega in holo]
    denom = LCM([denominator(omega) for omega in holo_forms])
    variable_names = 'x, y'
    for j in range(n):
        variable_names += ', z' + str(j)
    Rxyz = PolynomialRing(F, n+2, variable_names)
    x, y = Rxyz.gens()[:2]
    z = Rxyz.gens()[2:]
    holo_forms = [Rxyz(omega*denom) for omega in holo_forms]
    A = [[] for i in range(n)]
    for omega in holo:
        for i in range(n):
            ei = n*[0]
            ei[i] = 1
            omega1 = omega.group_action(ei)
            omega1 = denom * omega1
            v1 = omega1.coordinates(holo_forms)
            A[i] += [v1]
    for i in range(n):
        A[i] = matrix(F, A[i])
        A[i] = A[i].transpose()
    return A

def group_action_matrices_log(C_AS):
    F = C_AS.base_ring
    n = C_AS.height
    holo = C_AS.at_most_poles_forms(1)
    holo_forms = [omega.form for omega in holo]
    denom = LCM([denominator(omega) for omega in holo_forms])
    variable_names = 'x, y'
    for j in range(n):
        variable_names += ', z' + str(j)
    Rxyz = PolynomialRing(F, n+2, variable_names)
    x, y = Rxyz.gens()[:2]
    z = Rxyz.gens()[2:]
    holo_forms = [Rxyz(omega*denom) for omega in holo_forms]
    A = [[] for i in range(n)]
    for omega in holo:
        for i in range(n):
            ei = n*[0]
            ei[i] = 1
            omega1 = omega.group_action(ei)
            omega1 = denom * omega1
            v1 = omega1.coordinates(holo_forms)
            A[i] += [v1]
    for i in range(n):
        A[i] = matrix(F, A[i])
        A[i] = A[i].transpose()
    return A


#given a set S of (form, corresponding Laurent series at some pt), find their combinations holomorphic at that pt
def holomorphic_combinations_fcts(S, pole_order):
    C_AS = S[0][0].curve
    p = C_AS.characteristic
    F = C_AS.base_ring
    prec = C_AS.prec
    Rt.<t> = LaurentSeriesRing(F, default_prec=prec)
    RtQ = FractionField(Rt)
    minimal_valuation = min([Rt(g[1]).valuation() for g in S])
    if minimal_valuation >= -pole_order:
        return [s[0] for s in S]
    list_of_lists = [] #to będzie lista złożona z list współczynników część nieholomorficznych rozwinięcia form z S
    for eta, eta_exp in S:
        a = -minimal_valuation + Rt(eta_exp).valuation()
        list_coeffs = a*[0] + Rt(eta_exp).list() + (-minimal_valuation)*[0]
        list_coeffs = list_coeffs[:-minimal_valuation - pole_order]
        list_of_lists += [list_coeffs]
    M = matrix(F, list_of_lists)
    V = M.kernel() #chcemy wyzerować części nieholomorficzne, biorąc kombinacje form z S


    # Sprawdzamy, jakim formom odpowiadają elementy V.
    forms = []
    for vec in V.basis():
        forma_holo = as_function(C_AS, 0)
        forma_holo_power_series = Rt(0)
        for vec_wspolrzedna, elt_S in zip(vec, S):
            eta = elt_S[0]
            #eta_exp = elt_S[1]
            forma_holo += vec_wspolrzedna*eta
            #forma_holo_power_series += vec_wspolrzedna*eta_exp
        forms += [forma_holo]
    return forms

#given a set S of (form, corresponding Laurent series at some pt), find their combinations holomorphic at that pt
def holomorphic_combinations_forms(S, pole_order):
    C_AS = S[0][0].curve
    p = C_AS.characteristic
    F = C_AS.base_ring
    prec = C_AS.prec
    Rt.<t> = LaurentSeriesRing(F, default_prec=prec)
    RtQ = FractionField(Rt)
    minimal_valuation = min([Rt(g[1]).valuation() for g in S])
    if minimal_valuation >= -pole_order:
        return [s[0] for s in S]
    list_of_lists = [] #to będzie lista złożona z list współczynników część nieholomorficznych rozwinięcia form z S
    for eta, eta_exp in S:
        a = -minimal_valuation + Rt(eta_exp).valuation()
        list_coeffs = a*[0] + Rt(eta_exp).list() + (-minimal_valuation)*[0]
        list_coeffs = list_coeffs[:-minimal_valuation - pole_order]
        list_of_lists += [list_coeffs]
    M = matrix(F, list_of_lists)
    V = M.kernel() #chcemy wyzerować części nieholomorficzne, biorąc kombinacje form z S


    # Sprawdzamy, jakim formom odpowiadają elementy V.
    forms = []
    for vec in V.basis():
        forma_holo = as_form(C_AS, 0)
        forma_holo_power_series = Rt(0)
        for vec_wspolrzedna, elt_S in zip(vec, S):
            eta = elt_S[0]
            #eta_exp = elt_S[1]
            forma_holo += vec_wspolrzedna*eta
            #forma_holo_power_series += vec_wspolrzedna*eta_exp
        forms += [forma_holo]
    return forms

#print only forms that are log at the branch pts, but not holomorphic
def only_log_forms(C_AS):
    list1 = AS.at_most_poles_forms(0)
    list2 = AS.at_most_poles_forms(1)
    result = []
    for a in list2:
        if not(are_forms_linearly_dependent(list1 + result + [a])):
            result += [a]
    return result

def magmathis(A, B, text = False):
    """Find decomposition of Z/p^2-module given by matrices A, B into indecomposables using magma.
    If text = True, print the command for Magma. Else - return the output of Magma free."""
    q = parent(A).base_ring().order()
    n = A.dimensions()[0]
    A = str(list(A))
    B = str(list(B))
    A = A.replace("(", "")
    A = A.replace(")", "")
    B = B.replace("(", "")
    B = B.replace(")", "")
    result = "A := MatrixAlgebra<GF("+str(q) + "),"+ str(n) + "|"
    result += A + "," + B
    result += ">;"
    result += "M := RModule(RSpace(GF("+str(q)+")," + str(n) + "), A);"
    result += "IndecomposableSummands(M);"
    if text:
        return result
    print(magma_free(result))
    
    
def as_reduction(AS, fct):
    n = AS.height
    F = AS.base_ring
    variable_names = 'x, y'
    for i in range(n):
        variable_names += ', z' + str(i)
    Rxyz = PolynomialRing(F, n+2, variable_names)
    x, y = Rxyz.gens()[:2]
    z = Rxyz.gens()[2:]
    RxyzQ = FractionField(Rxyz)
    ff = AS.functions
    ff = [RxyzQ(F.function) for F in ff]
    fct = RxyzQ(fct)
    fct1 = numerator(fct)
    fct2 = denominator(fct)
    if fct2 != 1:
        return as_reduction(AS, fct1)/as_reduction(AS, fct2)
    
    result = RxyzQ(0)
    change = 0
    for a in fct1.monomials():
        degrees_zi = [a.degree(z[i]) for i in range(n)]
        d_div = [a.degree(z[i])//p for i in range(n)]
        if d_div != n*[0]:
            change = 1
        d_rem = [a.degree(z[i])%p for i in range(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))
        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)))
        
        
    if change == 0:
        return RxyzQ(result)
    else:
        print(fct, '\n')
        return as_reduction(AS, RxyzQ(result))

Testy

p = 5
m = 2
Rx.<x> = PolynomialRing(GF(p))
f = x^3 + x^2 + 1
C_super = superelliptic(f, m)
Rxy.<x, y> = PolynomialRing(GF(p), 2)
fArS1 = superelliptic_function(C_super, y*x)
fArS2 = superelliptic_function(C_super, y*x^2)
fArS3 = superelliptic_function(C_super, y)
AS1 = as_cover(C_super, [fArS1, fArS2, fArS3], prec=500)
print(AS1.genus())
AS2 = as_cover(C_super, [fArS2, fArS3, fArS1], prec=500)
print(AS2.genus())
325
325
p = 5
m = 2
Rx.<x> = PolynomialRing(GF(p))
f = x^3 + x^2 + 1
C_super = superelliptic(f, m)
Rxy.<x, y> = PolynomialRing(GF(p), 2)
fArS1 = superelliptic_function(C_super, y*x)
fArS2 = superelliptic_function(C_super, y*x^2)
fArS3 = superelliptic_function(C_super, y)
AS1 = as_cover(C_super, [fArS1, fArS2, fArS3], prec=1000)
print(AS1.exponent_of_different())
omega = as_form(AS1, 1/y)
print(omega.expansion_at_infty())

Brudnopis

AS.at_most_poles(14)
[1,
 z1,
 z1^2,
 z0,
 z0*z1,
 -x*z0^2 + z0*z1^2,
 z0^2,
 z0^2*z1,
 z0^2*z1^2 + x*z0*z1 + x^2,
 y,
 -y*z0^2 + y*z1,
 y*z0,
 x,
 -x*z0^2 + x*z1,
 x*z0]
only_log_forms(AS)
[((z0^2*z1^2 + x*z0*z1 + x^2)/y) * dx]
dpp = AS.exponent_of_different_prim()
print(at_most_poles(AS, dpp))
I haven't found all forms.
[1]
Rxy.<x, y, z0, z1> = PolynomialRing(F, 4)
omega = as_form(AS, z0^2*z1^2/y - (a+1)*z0*x/y)
print(omega - omega.group_action([1, 0]))
((z0*z1^2 - z1^2 + (a + 1)*x)/y) * dx
g = as_function(AS, z0^2*z1^2+z0*(a+1)*x)
g.expansion_at_infty(i = 0)
a*t^-10 + 2*t^-8 + a*t^-6 + (a + 1)*t^-2 + O(t^288)
p = 7
m = 2
F = GF(p)
Rx.<x> = PolynomialRing(F)
f = x^3 + 1
C_super = superelliptic(f, m)

Rxy.<x, y> = PolynomialRing(F, 2)
f1 = superelliptic_function(C_super, x^2*y)
f2 = superelliptic_function(C_super, x^3)
AS = as_cover(C_super, [f1, f2], prec=1000)
#print(AS.magical_element())
#print(AS.pseudo_magical_element())
#A, B = group_action_matrices2_log(AS1)
A, B = group_action_matrices2(AS)
n = A.dimensions()[0]
print(A*B == B*A)
print(A^p == identity_matrix(n))
print(B^p == identity_matrix(n))
magmathis(A, B, p, deg = 1)
True
True
True

print(f1.expansion_at_infty())
print(f2.expansion_at_infty())
t^-5 + 2*t + O(t^5)
t^-6 + 3 + O(t^4)
AS.magical_element()
[z0^2]
AS1.exponent_of_different_prim()
68
AS1.jumps
[[5, 2], [5, 2]]
magmathis(A, B, 3, deg = 2)
[
RModule of dimension 1 over GF(3^2),
RModule of dimension 1 over GF(3^2),
RModule of dimension 3 over GF(3^2),
RModule of dimension 3 over GF(3^2),
RModule of dimension 6 over GF(3^2),
RModule of dimension 6 over GF(3^2),
RModule of dimension 8 over GF(3^2)
]
magmatxtx(A, B, p)
WARNING: Some output was deleted.

A := MatrixAlgebra<GF(3^1),14|[1, 2, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 1, 2, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 1, 1],[1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1]>; M := RModule(RSpace(GF(3^1),14), A); lll := IndecomposableSummands(M); Morphism(lll[2], M);

parent(A).base_ring().order()
3
Rx.<x> = PolynomialRing(GF(3))
C = superelliptic(x^3 - x, 2)
C.de_rham_basis()
[((1/y) dx, 0, (1/y) dx), ((x/y) dx, 2/x*y, ((-1)/(x*y)) dx)]
Rxy.<x, y> = PolynomialRing(GF(p), 2)
omega = superelliptic_form(C, 1/y^3)
ff = superelliptic_function(C, x)
cycle = superelliptic_cech(C, omega, ff)
cycle.coordinates()
(0, 0)
def dual_elt(AS, zmag):
    p = AS.characteristic
    n = AS.height
    group_elts = [(i, j) for i in range(p) for j in range(p)]
    variable_names = 'x, y'
    for i in range(n):
        variable_names += ', z' + str(i)
    Rxyz = PolynomialRing(F, n+2, variable_names)
    x, y = Rxyz.gens()[:2]
    z = Rxyz.gens()[2:]
    RxyzQ = FractionField(Rxyz)
    M = matrix(RxyzQ, p^n, p^n)
    for i in range(p^n):
        for j in range(p^n):
            M[i, j] = (zmag.group_action(group_elts[i])*zmag.group_action(group_elts[j])).trace2()
    main_det = M.determinant()
    zvee = as_function(AS, 0)
    for i in range(p^n):
        Mprim = matrix(RxyzQ, M)
        Mprim[:, i] = vector([(j == 0) for j in range(p^2)])
        fi = Mprim.determinant()/main_det
        zvee += fi*zmag.group_action(group_elts[i])
    return zvee
p = 5
m = 1
F = GF(p)
Rx.<x> = PolynomialRing(F)
f = x
C_super = superelliptic(f, m)

Rxy.<x, y> = PolynomialRing(F, 2)
f1 = superelliptic_function(C_super, x^2)
f2 = superelliptic_function(C_super, x^3)
AS = as_cover(C_super, [f1, f2], prec=500)
zmag = (AS.magical_element())[0]
n = 2
variable_names = 'x, y'
for i in range(n):
    variable_names += ', z' + str(i)
Rxyz = PolynomialRing(F, n+2, variable_names)
x, y = Rxyz.gens()[:2]
z = Rxyz.gens()[2:]
f3 = as_function(AS, z[0]^p - z[0])
zdual = dual_elt(AS, zmag)
(zmag*(zdual.group_action([0, 0]))).trace2()
1
n = 2
F = GF(5)
variable_names = 'x, y'
for i in range(n):
    variable_names += ', z' + str(i)
Rxyz = PolynomialRing(F, n+2, variable_names)
x, y = Rxyz.gens()[:2]
z = Rxyz.gens()[2:]
def ith_component(omega, zvee, i):
    AS = omega.curve
    p = AS.characteristic
    group_elts = [(j1, j2) for j1 in range(p) for j2 in range(p)]
    z_vee_fct = zvee.group_action(group_elts[i]).function
    new_form = as_form(AS, z_vee_fct*omega.form)
    return new_form.trace2()
om = AS.holomorphic_differentials_basis()[4]
om
(-z0^3*z1^2 - 2*x*z0^2*z1 + z1^4 + 2*x^2*z0) * dx
def combination_components(omega, zmag, w):
    AS = omega.curve
    p = AS.characteristic
    group_elts = [(j1, j2) for j1 in range(p) for j2 in range(p)]
    zvee = dual_elt(AS, zmag)
    result = as_form(AS, 0)
    for i in range(p^2):
        omegai = ith_component(omega, zvee, i)
        aux_fct1 = w.group_action(group_elts[i]).function
        aux_fct2 = omegai.form
        result += as_form(AS, aux_fct1*aux_fct2)
    return result
combination_components(om, zmag, zmag)
(-z0^3*z1^2 - 2*x*z0^2*z1 + z1^4 + 2*x^2*z0) * dx
om
(-z0^3*z1^2 - 2*x*z0^2*z1 + z1^4 + 2*x^2*z0) * dx
zmag.expansion_at_infty()
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)
w1 = as_function(AS, z[0]^2*y^3/x - z[0]^2*y^2)
w1.expansion_at_infty()
O(t^430)
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)
om1 = combination_components(om, zmag, w1)
print(om1)
print(om1.expansion_at_infty())
(0) * dx
O(t^474)