rozniczkowanie form drw zrobione (?)

This commit is contained in:
jgarnek 2023-02-24 12:51:49 +00:00
parent 856d74212b
commit 22872e266b
6 changed files with 1466 additions and 69 deletions

File diff suppressed because one or more lines are too long

View File

@ -1,9 +1,9 @@
p = 2
m = 1
p = 3
m = 2
F = GF(p)
Rx.<x> = PolynomialRing(F)
f = x
f = x^3 - x
C = superelliptic(f, m)
xx = C.x
AS = as_cover(C, [xx^5, xx^5 + xx^3])
print(AS.magical_element())
a = superelliptic_drw_form(C.one, 0*C.dx, 0*C.x)
b = a+a+a+a+a+a+a+a+a
print(b)

View File

@ -0,0 +1,30 @@
def patch(C):
if C.exponent != 2:
raise ValueError("Not implemented yet!")
Fxy, Rxy, x, y = C.fct_field
F = C.base_ring
Rx.<x> = PolynomialRing(F)
f = C.polynomial
g = C.genus()
f_star = Rx(x^(2*g+2)*f(1/x))
return superelliptic(f_star, 2)
def second_patch(argument):
C = argument.curve
C1 = patch(C)
Fxy, Rxy, x, y = C.fct_field
g = C.genus()
if isinstance(argument, superelliptic_function):
fct = Fxy(argument.function)
fct1 = Fxy(fct.subs({x : 1/x, y : y/x^(g+1)}))
return superelliptic_function(C1, fct1)
if isinstance(argument, superelliptic_form):
fct = Fxy(argument.form)
fct1 = Fxy(fct.subs({x : 1/x, y : y/x^(g+1)}))
fct1 *= -x^(-2)
return superelliptic_form(C1, fct1)
def lift_form_to_drw(omega):
A, B = regular_form(omega)
A, B = A.change_ring(QQ), B.change_ring(QQ)
print("%s dx + %s dy"%(A, B))

View File

@ -1,39 +1,105 @@
class superelliptic_witt:
def __init__(self, C, t, f0, f1):
''' Define Witt function on C of the form [t] + f0([x], [y]) + V(f1). '''
self.curve = C
def __init__(self, t, f):
''' Define Witt function on C of the form [t] + V(f). '''
self.curve = t.curve
C = t.curve
p = C.characteristic
self.t = t #superelliptic_function
self.f0 = reduce_rational_fct(f0, p) #polynomial/rational function over Z/p^2
self.f1 = f1 #superelliptic_function
self.f = f #superelliptic_function
def __repr__(self):
f0 = self.f0
f1 = self.f1
f = self.f
t = self.t
return "[" + str(t) + "] + " + str(f0).replace("x", "[x]").replace("y", "[y]") + " + V(" + str(f1) + ")"
if f.function == 0:
return "[" + str(t) + "]"
if t.function == 0:
return "V(" + str(f) + ")"
return "[" + str(t) + "] + V(" + str(f) + ")"
def __neg__(self):
f = self.f
t = self.t
return superelliptic_witt(-t, -f)
def __add__(self, other):
C = self.curve
return superelliptic_witt(C, self.t + other.t, self.f0 + other.f0, self.f1 + other.f1)
def teichmuller_representation(self):
'''Represents as [f] + V(g), i.e. f0 = 0.'''
C = self.curve
Fxy, Rxy, x, y = self.fct_field
F = C.base_ring
function = Rxy(self.f0)
if self.f0 == 0:
return self
M = Rxy(function.monomials()[0])
a = F(function.monomial_coefficient(M))
fct1 = fct - superelliptic_function(C, a*M)
function1 = fct1.function
return teichmuller(fct1) + superelliptic_witt(C, (a.lift())^p*M.change_ring(QQ), superelliptic_function(C, function1^2*a*M + function1*(a*M)^2))
second_coor = 0*C.x
X = self.t
Y = other.t
for i in range(1, p):
second_coor -= binomial_prim(p, i)*X^i*Y^(p-i)
return superelliptic_witt(self.t + other.t, self.f + other.f + second_coor)
def antiteichmuller_representation(self):
'''Represents as f([x], [y]) + V(g), i.e. teichmuller part is zero.'''
return 0
def __sub__(self, other):
return self + (-other)
def __rmul__(self, other):
p = self.curve.characteristic
if other in ZZ:
if other == 0:
return superelliptic_witt(0*C.x, 0*C.x)
if other > 0:
return self + (other-1)*self
if other < 0:
return (-other)*(-self)
if other in QQ:
other_integer = Integers(p^2)(other)^(-1)
return other_integer*self
def __mul__(self, other):
if isinstance(other, superelliptic_witt):
t1 = self.t
f1 = self.f
p = self.curve.characteristic
t2 = other.t
f2 = other.f
return superelliptic_witt(t1*t2, t1^p*f2 + t2^p*f1)
if isinstance(other, superelliptic_drw_form):
h1 = other.h1
h2 = other.h2
omega = other.omega
p = other.curve.characteristic
t = self.t
f = self.f
aux_form = t^p*omega - h2*t^(p-1)*t.diffn() + f*h1^p*(C.x)^(p-1)*C.dx
return superelliptic_drw_form(t*h1, aux_form, t^p*h2)
def __eq__(self, other):
return self.t == other.t and self.f == other.f
def diffn(self):
C = self.curve
t = self.t
f = self.f
fC = C.polynomial
F = C.base_ring
Rxy.<x, y> = PolynomialRing(F, 2)
if t.function == 0:
return superelliptic_drw_form(0*C.x, 0*C.dx, f)
t_polynomial = t.function
num = t_polynomial.numerator()
den = t_polynomial.denominator()
num_t_fct = superelliptic_function(C, num)
den_t_fct = superelliptic_function(C, den)
inv_den_t_fct = superelliptic_function(C, 1/den)
if den != 1:
# d([N/D] + V(f)) = [1/D]*d([N]) - [N]*[D^(-2)]*d([D]) + dV(f)
return ((den_t_fct)^(-1)).teichmuller()*num_t_fct.teichmuller().diffn() - ((den_t_fct)^(-2)).teichmuller()*num_t_fct.teichmuller()*den_t_fct.teichmuller().diffn() + superelliptic_drw_form(0*C.x, 0*C.dx, f)
t_polynomial = Rxy(t_polynomial)
fC = superelliptic_function(C, fC)
fC = fC.teichmuller()
dy_w = 1/2* ((C.y)^(-1)).teichmuller()*auxilliary_derivative(fC)
M = t_polynomial.monomials()[0]
a = t_polynomial.monomial_coefficient(M)
#[P] = [aM] + Q, where Q = ([P] - [aM])
aM_fct = superelliptic_function(C, a*M)
Q = self - aM_fct.teichmuller()
exp_x = M.exponents()[0][0]
exp_y = M.exponents()[0][1]
return Q.diffn() + exp_x*superelliptic_drw_form(aM_fct/C.x, 0*C.dx, 0*C.x) + exp_y*(aM_fct/C.y).teichmuller()*dy_w
def binomial_prim(p, i):
return binomial(p, i)/p
def reduce_rational_fct(fct, p):
Rxy.<x, y> = PolynomialRing(QQ)
@ -53,58 +119,127 @@ def reduce_rational_fct(fct, p):
def teichmuller(fct):
C = fct.curve
Fxy, Rxy, x, y = C.fct_field
F = C.base_ring
function = Rxy(fct.function)
if function == 0:
return superelliptic_witt(C, 0, 0*C.x)
M = Rxy(function.monomials()[0])
a = F(function.monomial_coefficient(M))
fct1 = fct - superelliptic_function(C, a*M)
function1 = fct1.function
return teichmuller(fct1) + superelliptic_witt(C, (a.lift())^p*M.change_ring(QQ), superelliptic_function(C, function1^2*a*M + function1*(a*M)^2))
return superelliptic_witt(fct, 0*C.x)
superelliptic_function.teichmuller = teichmuller
#dy = [f(x)]'/2*y dx
#[f1 + M] = [f1] + [M] + V(cos)
#d[f1 + M] = d[f1] + d[M] + dV(f1*M)
#M = b x^a
#d[M] = a*[b x^(a-1)]
def auxilliary_derivative(P):
'''Return "derivative" of P, where P depends only on x. '''
P0 = P.t.function
P1 = P.f.function
C = P.curve
F = C.base_ring
Rx.<x> = PolynomialRing(F)
P0 = Rx(P0)
P1 = Rx(P1)
if P0 == 0:
return superelliptic_drw_form(0*C.x, 0*C.dx, P.f)
M = P0.monomials()[0]
a = P0.monomial_coefficient(M)
#[P] = [aM] + Q, where Q = ([P] - [aM])
aM_fct = superelliptic_function(C, a*M)
Q = P - aM_fct.teichmuller()
exp = M.exponents()[0]
return auxilliary_derivative(Q) + exp*superelliptic_drw_form(aM_fct/C.x, 0*C.dx, 0*C.x)
class superelliptic_drw_form:
def __init__(self, C, omega_x, omega_y, omega, h):
'''Form [omega_x] d[x] + [omega_y] d[y] + V(omega) + dV([h])'''
self.curve = C
self.omega_x = omega_x
self.omega_y = omega_y
def __init__(self, h1, omega, h2):
'''Form [h1] d[x] + V(omega) + dV([h])'''
self.curve = h1.curve
self.h1 = h1
self.omega = omega
self.h = h
self.h2 = h2
def __eq__(self, other):
eq1 = (self.omega1 == self.omega1)
eq1 = (self.h1 == self.h1)
try:
H = (self.h - other.h).pthroot()
H = (self.h2 - other.h2).pth_root()
except:
return False
eq2 = (self.omega2 - other.omega2).cartier() - H.diffn()
eq2 = (self.omega - other.omega).cartier() - H.diffn()
if eq1 and eq2:
return True
return False
def __repr__(self):
C = self.curve
omega_x = self.omega_x
omega_y = self.omega_y
h = self.h
return str(omega_x) + "] dx + [" + str(omega_x.form) + "] d[x] " + "+ V(" + str(omega2) + ") + dV([" + str(h) +"])"
h1 = self.h1
omega = self.omega
h2 = self.h2
result = ""
if h1.function != 0:
result += "[" + str(h1) + "] d[x]"
if h1.function !=0 and omega.form != 0:
result += " + "
if omega.form != 0:
result += "V(" + str(omega) + ")"
if h2.function !=0 and omega.form != 0:
result += " + "
if h2.function != 0:
result += "dV([" + str(h2) +"])"
if result == "":
result += "0"
return result
def __rmul__(self, other):
h1 = self.h1
h2 = self.h2
omega = self.omega
p = self.curve.characteristic
if other in ZZ:
if other == 0:
return superelliptic_drw_form(0*C.x, 0*C.dx, 0*C.x)
if other > 0:
return self + (other-1)*self
if other < 0:
return (-other)*(-self)
if other in QQ:
other_integer = Integers(p^2)(other)^(-1)
return other_integer*self
t = other.t
f = other.f
aux_form = t^p*omega - h2*t^(p-1)*t.diffn() + f*h1^p*(C.x)^(p-1)*C.dx
return superelliptic_drw_form(t*h1, aux_form, t^p*h2)
def __neg__(self):
C = self.curve
h1 = self.h1
h2 = self.h2
omega = self.omega
return superelliptic_drw_form(-h1, -omega, -h2)
def __add__(self, other):
C = self.curve
h1 = self.h1
h2 = self.h2
omega = self.omega
H1 = other.h1
H2 = other.h2
OMEGA = other.omega
aux = (teichmuller(h1) + teichmuller(H1))*superelliptic_drw_form(C.one, 0*C.dx, 0*C.x)
h1_aux = aux.h1
h2_aux = aux.h2
omega_aux = aux.omega
return superelliptic_drw_form(h1_aux, omega + OMEGA + omega_aux, h2 + H2 + h2_aux)
def frobenius(self):
C = self.curve
h1 = self.h1
h2 = self.h2
p = C.characteristic
return h1^p*C.x^(p-1)*C.dx + h2.diffn()
def mult_by_p(omega):
C = omega.curve
fct = omega.form
Fxy, Rxy, x, y = C.fct_field
omega2 = superelliptic_form(C, fct^p * x^(p-1))
result = superelliptic_drw_form(C, 0*C.dx, omega2, 0*C.x)
omega = superelliptic_form(C, fct^p * x^(p-1))
result = superelliptic_drw_form(C, 0*C.dx, omega, 0*C.x)
return result
def basis_W2Omega(C):
basis = C.holomorphic_differentials_basis()
result = []
for omega in basis:
result += [mult_by_p(omega)]
image_of_cartier = []
return result

View File

@ -21,6 +21,11 @@ class superelliptic_form:
g = reduction(C, g1 - g2)
return superelliptic_form(C, g)
def __neg__(self):
C = self.curve
g = self.form
return superelliptic_form(C, -g)
def __repr__(self):
g = self.form
if len(str(g)) == 1:
@ -32,6 +37,9 @@ class superelliptic_form:
omega = self.form
return superelliptic_form(C, constant*omega)
def __eq__(self, other):
return self.form == other.form
def cartier(self):
'''Computes Cartier operator on the form. Idea: y^m = f(x) -> y^(p^r - 1) = f(x)^M, where r = ord_p(m),
M = (p^r - 1)/m. Thus h(x)/y^j dx = h(x) f(x)^(M*j)/y^(p^r * j) dx. Thus C(h(x)/y^j dx) = 1/y^(p^(r-1)*j) C(h(x) f(x)^(M*j) dx).'''

View File

@ -38,6 +38,11 @@ class superelliptic_function:
g = reduction(C, g1 + g2)
return superelliptic_function(C, g)
def __neg__(self):
C = self.curve
g = self.function
return superelliptic_function(C, -g)
def __sub__(self, other):
C = self.curve
g1 = self.function