Implemented splitting for quasi-indecomposable bundles

This commit is contained in:
Montessinos Mickael Gerard Bernard 2024-03-12 23:18:32 +02:00
parent 2809d9c8c5
commit 7c59a58b5e
5 changed files with 188 additions and 54 deletions

View File

@ -171,21 +171,22 @@ def wedderburn_malcev_basis(A):
imat = mat**-1
#We may now begin running the algorithm in earnest.
basis = [A(c) for c in mat.columns()]
table = [imat * c.left_matrix() * mat for c in basis[:r]]
table = [imat * c.left_matrix().transpose() * mat for c in basis[:r]]
for first, stop in zip(rad_power_steps[:-1], rad_power_steps[1:]):
b_pairs = list(itertools.product(enumerate(basis[:r]), repeat=2))
d_tuples = [[0]*i + [d] + [0]*(r-1-i)
for i in range(r) for d in basis[first:stop]]
#(first, stop) marks the positions of the vectors generating V_i.
eq = matrix([sum([list(mat[:,first:stop].solve_right(
eq = matrix([sum([list(mat.solve_right(
(bi*ds[j]
+ ds[i]*bj
- sum([table[i][j,s]*d for (s, d) in enumerate(ds)])).vector()))
- sum([table[i][s,j]*d
for (s, d) in enumerate(ds)])).vector())[first:stop])
for (i, bi),(j, bj) in b_pairs], [])
for ds in d_tuples])
target = vector(sum([list(mat[:,first:stop].solve_right(
(sum([table[i][j,s]*b for s, b in enumerate(basis[:r])])
- bi*bj).vector()))
target = vector(sum([list(mat.solve_right(
(sum([table[i][s,j]*b for s, b in enumerate(basis[:r])])
- bi*bj).vector())[first:stop])
for (i, bi), (j, bj) in b_pairs], []))
sol = eq.solve_left(target)
sol_chunks = [sol[i: i + stop - first]
@ -357,7 +358,7 @@ def zero_div_from_min_poly(a, f):
return((facto[0][0]**facto[0][1])(a),
prod([(fa[0]**fa[1])(a) for fa in facto[1:]]))
if facto[0][1] > 1:
return((facto[0][0])(a), (facto[0]**(facto[0][1]-1))(a))
return((facto[0][0])(a), (facto[0][0]**(facto[0][1]-1))(a))
return None

View File

@ -94,7 +94,10 @@ def atiyah_bundle(field, rank, degree, base=None):
sage: F.<x> = FunctionField(GF(11))
sage: R.<y> = F[]
sage: K.<y> = F.extension(y^2 - x^3 - x)
sage: base = VectorBundle(K, K.places_finite()[0].divisor())
sage: base = VectorBundle(
....: K,
....: K.places_finite()[0].divisor()
....: - K.places_infinite()[0].divisor())
sage: E = atiyah_bundle(K, 5, 3, base)
sage: E.rank()
5
@ -133,6 +136,7 @@ def atiyah_bundle(field, rank, degree, base=None):
line_bundle = VectorBundle(field, divisor)
for _ in range(starting_rank - 1):
result = result.extension_by_global_sections()
result = result.tensor_product(base)
if degree > 0:
result = result.tensor_product(line_bundle)
for op, reps in reversed(plan):

View File

@ -248,6 +248,9 @@ def infinite_order_xgcd(ideals):
sage: all([a[i] in primes[i] for i in range(2)])
True
"""
s = len(ideals)
non_zero_indices = [i for i, ideal in enumerate(ideals) if ideal != 0]
ideals = [ideal for i, ideal in enumerate(ideals) if ideal != 0]
order_basis = ideals[0].ring().basis()
if order_basis[0] != 1:
raise ValueError('The first element of the basis of the order should'
@ -265,7 +268,11 @@ def infinite_order_xgcd(ideals):
if not (H/den).is_one():
raise ValueError("The ideals should be coprime.")
v = U[:,0].list()
return [sum([ideals_bases[i][j]*v[n*i+j] for j in range(n)]) for i in range(k)]
coefs = [sum([ideals_bases[i][j]*v[n*i+j] for j in range(n)]) for i in range(k)]
res = [0]*s
for i, c in zip(non_zero_indices, coefs):
res[i] = c
return res
def infinite_approximation(places,valuations,residues):
@ -878,9 +885,9 @@ def pseudo_hermite_form(ideals, mat, include_zero_cols=True, transformation=Fals
j = k-1
for i in range(n-1, -1, -1):
#Check zero
if all([h[i,m] == 0 for m in range(j)]):
if all([h[i,m] == 0 for m in range(j+1)]):
continue
m = [h[i,m] == 0 for m in range(j, -1, -1)].index(False)
m = [h[i,m] == 0 for m in range(j+1)].index(False)
h[:,m], h[:,j] = h[:,j], h[:,m]
U[:,m], U[:,j] = U[:,j], U[:,m]
h_ideals[m], h_ideals[j] = h_ideals[j], h_ideals[m]
@ -950,18 +957,19 @@ def hermite_form_infinite_polymod(mat, include_zero_cols=True, transformation=Fa
U = identity_matrix(K,k)
j = k-1
for i in range(n-1, -1, -1):
if all([h[i,m] == 0 for m in range(j)]):
if all([h[i,m] == 0 for m in range(j+1)]):
continue
min_vals = [min([h[i,m].valuation(place) for m in range(j+1)])
for place in places]
gcd = prod([pi**m for pi, m in zip(pis, min_vals)])
#put gcd on the diagonal
ideals = [O.ideal(h[i,m]/gcd) for m in range(j+1)]
ideals = [O.ideal(h[i,m]/gcd) if not h[i,m].is_zero() else 0
for m in range(j+1)]
coefs = infinite_order_xgcd(ideals)
m = [c == 0 for c in coefs].index(False)
U[:,j], U[:, m] = gcd*sum([c/h[i,m]*U[:,m]
ell = [c == 0 for c in coefs].index(False)
U[:,j], U[:, ell] = gcd*sum([(c/h[i,ell])*U[:,m]
for m,c in enumerate(coefs)]), U[:, j]
h[:,j], h[:, m] = gcd*sum([c/h[i,m]*h[:,m]
h[:,j], h[:, ell] = gcd*sum([(c/h[i,ell])*h[:,m]
for m,c in enumerate(coefs)]), h[:, j]
assert(all([u in O for u in U.list()]))
#eliminate coefficients left of diagonal

View File

@ -304,8 +304,12 @@ class HomBundle(VectorBundle):
sage: all([all([a in O_infinity for a in (x**-2 * mat * g_infinite).list()]) for mat in h0])
True
"""
if self._h0 is not None:
return self._h0
h0 = super().h0()
return [self._vector_to_matrix(v) for v in h0]
h0 = [self._vector_to_matrix(v) for v in h0]
self._h0 = h0
return h0
def coordinates_in_h0(self, mat):
r"""
@ -339,7 +343,8 @@ class HomBundle(VectorBundle):
`H^0(\mathrm{self})`.
That is, a vector bundle `V` together with an injective morphism of `V`
into ``self.codomain``
into ``self.codomain`` such that the image of `V` in ``self.codomain``
is also the image of ``f``.
EXAMPLES ::
@ -390,11 +395,11 @@ class EndBundle(HomBundle):
"""
def __init__(self, bundle):
HomBundle.__init__(self, bundle, bundle)
self._A = None
def __repr__(self):
return "Endomorphism bundle of %s" % (self._domain)
@cached_method
def global_algebra(self):
r"""
Return ``self.h0()`` as a k-algebra in which computations may be done.
@ -413,6 +418,8 @@ class EndBundle(HomBundle):
Finite-dimensional algebra of degree 4 over Finite Field of size 7
sage: h0 = end.h0()
"""
if self._A is not None:
return self._A
k = self._function_field.constant_base_field()
category = Algebras(k).FiniteDimensional().WithBasis().Associative()
basis = self.h0()
@ -425,6 +432,7 @@ class EndBundle(HomBundle):
category = category)
to_a = lambda mat : algebra(self.coordinates_in_h0(mat))
from_a = lambda a : self.h0_from_vector(a.vector())
self._A = (algebra, to_a, from_a)
return algebra, to_a, from_a
@cached_method

View File

@ -78,8 +78,9 @@ from sage.arith.functions import lcm
from sage.arith.misc import integer_ceil
from sage.functions.log import logb, log
from sage.matrix.constructor import matrix
from sage.matrix.special import block_matrix, elementary_matrix\
, identity_matrix
from sage.matrix.special import (block_matrix, elementary_matrix,
identity_matrix, diagonal_matrix,
zero_matrix, block_diagonal_matrix)
from sage.matrix.matrix_space import MatrixSpace
from sage.rings.function_field.ideal import FunctionFieldIdeal
from sage.rings.function_field.function_field_rational\
@ -138,12 +139,13 @@ class VectorBundle(SageObject):
Vector bundle of rank 1 over Function field in y defined by y^2 + 2*x^3 + 2*x
"""
def __init__(self,function_field, ideals,g_finite=None,g_infinite=None):
def __init__(self,function_field, ideals,g_finite=None,g_infinite=None, check=True):
if g_finite is None or g_infinite is None:
self._line_bundle_from_divisor(function_field, ideals)
self._line_bundle_from_divisor(function_field, ideals, check=check)
else:
self._vector_bundle_from_data(function_field, ideals,
g_finite,g_infinite)
g_finite,g_infinite, check)
self._h0 = None
self._h0_matrix = None
self._h0_Kp = None
self._h0_vs = None
@ -167,13 +169,14 @@ class VectorBundle(SageObject):
self._function_field,
)
def _line_bundle_from_divisor(self,function_field,divisor):
def _line_bundle_from_divisor(self,function_field,divisor, check=True):
r"""
Build a line bundle from a divisor
"""
if not function_field == divisor.parent().function_field():
raise ValueError('The divisor should be defined over the '
+ 'function field.')
if check:
if not function_field == divisor.parent().function_field():
raise ValueError('The divisor should be defined over the '
+ 'function field.')
self._function_field = function_field
couples = divisor.list()
finite_part = [c for c in couples if not c[0].is_infinite_place()]
@ -191,7 +194,7 @@ class VectorBundle(SageObject):
self._g_infinite = matrix(function_field,[[pi]])
def _vector_bundle_from_data(self,function_field,
ideals,g_finite,g_infinite):
ideals,g_finite,g_infinite, check=True):
r"""
Construct a vector bundle from data.
"""
@ -210,22 +213,23 @@ class VectorBundle(SageObject):
g_finite.change_ring(function_field)
g_infinite.change_ring(function_field)
r = len(ideals)
if (g_finite.nrows() != r
or g_finite.ncols() != r
or g_infinite.nrows() != r
or g_infinite.ncols() != r):
raise ValueError('The length of the ideal list must equal the \
size of the basis matrices')
if not g_finite.is_invertible() or not g_infinite.is_invertible():
raise ValueError('The basis matrices must be invertible')
if not all([isinstance(I,FunctionFieldIdeal)
for I in ideals]):
raise TypeError('The second argument must be a list of \
FunctionFieldIdeals.')
if not all([I.base_ring() == function_field.maximal_order()
for I in ideals]):
raise ValueError('All ideals must have the maximal order of\
function_field as base ring.')
if check:
if (g_finite.nrows() != r
or g_finite.ncols() != r
or g_infinite.nrows() != r
or g_infinite.ncols() != r):
raise ValueError('The length of the ideal list must equal'
+ ' the size of the basis matrices')
if not g_finite.is_invertible() or not g_infinite.is_invertible():
raise ValueError('The basis matrices must be invertible')
if not all([isinstance(I,FunctionFieldIdeal)
for I in ideals]):
raise TypeError('The second argument must be a list of \
FunctionFieldIdeals.')
if not all([I.base_ring() == function_field.maximal_order()
for I in ideals]):
raise ValueError('All ideals must have the maximal order of\
function_field as base ring.')
self._function_field = function_field
self._ideals = ideals
self._g_finite = g_finite
@ -379,10 +383,11 @@ class VectorBundle(SageObject):
I = prod(self._ideals)
determinant_finite = self._g_finite.determinant()
determinant_infinite = self._g_infinite.determinant()
return VectorBundle(self._function_field
,I
,determinant_finite
,determinant_infinite)
return VectorBundle(self._function_field,
I,
determinant_finite,
determinant_infinite,
check=False)
def degree(self):
r"""
@ -554,7 +559,8 @@ class VectorBundle(SageObject):
ideals = self._ideals + other._ideals
g_finite = block_matrix([[self._g_finite,0],[0,other._g_finite]])
g_infinite = block_matrix([[self._g_infinite,0],[0,other._g_infinite]])
return VectorBundle(self._function_field, ideals,g_finite,g_infinite)
return VectorBundle(self._function_field, ideals,
g_finite, g_infinite, check=False)
def _direct_sum_rec(self,acc,n):
r"""
@ -612,7 +618,8 @@ class VectorBundle(SageObject):
ideals = [I * J for I in self._ideals for J in other._ideals]
g_finite = self._g_finite.tensor_product(other._g_finite)
g_infinite = self._g_infinite.tensor_product(other._g_infinite)
return VectorBundle(self._function_field, ideals,g_finite,g_infinite)
return VectorBundle(self._function_field, ideals,
g_finite, g_infinite, check=False)
def _tensor_power_aux(self,acc,n):
r"""
@ -668,7 +675,8 @@ class VectorBundle(SageObject):
"""
O = K.maximal_order()
ideals = [O.ideal(I.gens()) for I in self._ideals]
return VectorBundle(K, ideals,self._g_finite,self._g_infinite)
return VectorBundle(K, ideals,self._g_finite,
self._g_infinite, check=False)
def restriction(self):
r"""
@ -721,7 +729,8 @@ class VectorBundle(SageObject):
for c in gen_infinite])
g_infinite = matrix([sum([a.list() for a in collumn],[])
for collumn in g_infinite]).transpose()
return VectorBundle(F, ideals,g_finite,g_infinite)
return VectorBundle(F, ideals,g_finite,
g_infinite, check=False)
def _h0_rational(self):
r"""
@ -855,6 +864,8 @@ class VectorBundle(SageObject):
sage: len(L.h0())
1
"""
if self._h0 is not None:
return self._h0
if isinstance(self._function_field,RationalFunctionField):
#Compute restriction to normalize the coefficient ideals.
return self.restriction()._h0_rational()
@ -867,6 +878,7 @@ class VectorBundle(SageObject):
h0.append(vector([sum([y**j * v[i*deg + j]
for j in range(deg)])
for i in range(self.rank())]))
self._h0 = h0
return h0
@cached_method
@ -1108,6 +1120,106 @@ class VectorBundle(SageObject):
if mat.is_unit():
return P
def _isomorphism_indecomposable(self, other):
r"""
Return an isomorphism from self to other if it exists.
Assumes that self is indecomposable.
TESTS::
sage: from vector_bundle import atiyah_bundle, canonical_bundle
sage: F.<x> = FunctionField(GF(7))
sage: R.<y> = F[]
sage: K.<y> = F.extension(y^2 - x^3 - x)
sage: V = atiyah_bundle(K, 2, 0)
sage: W = atiyah_bundle(K, 2, 0,canonical_bundle(K))
sage: isom = V._isomorphism_indecomposable(W)
sage: isom is not None
True
sage: hom_1 = V.hom(W)
sage: hom_2 = W.hom(V)
sage: isom == hom_1.h0_from_vector(hom_1.coordinates_in_h0(isom))
True
sage: isom^-1 == hom_2.h0_from_vector(hom_2.coordinates_in_h0(isom^-1))
True
"""
r = self.rank()
if other.rank() != r:
return None
V = self.direct_sum(other)
End = V.end()
A, to_A, from_A = End.global_algebra()
(S, to_S, from_S), factors = End._global_algebra_split()
if len(factors) > 1:
return None
split = factors[0]
if split.M.nrows() != 2:
return None
K = self._function_field
id_self = split.to_M(to_S(to_A(
diagonal_matrix(K, [1]*r + [0]*r))))
im_self = id_self.transpose().echelon_form().rows()
im_self = [r for r in im_self if not r.is_zero()]
id_other = split.to_M(to_S(to_A(
diagonal_matrix(K, [0]*r + [1]*r))))
im_other = id_other.transpose().echelon_form().rows()
im_other = [r for r in im_other if not r.is_zero()]
P = matrix(im_self + im_other).transpose()
F = split._K
s = split._n // 2
id_s = identity_matrix(F, s)
zero_mat = zero_matrix(F, s)
isom = P * block_matrix([[zero_mat, zero_mat], [id_s, zero_mat]]) * P**-1
isom = from_A(from_S(split.from_M(isom)))
return isom[r:,:r]
def _indecomposable_power_split(self):
r"""
Check if self is a direct sum of copies of an indecomposable bundle.
If so, return this indecomposable `L` and an isomorphism from
``L.direct_sum_repeat(s)`` to ``self``.
TESTS ::
sage: #long time (25 seconds)
sage: from vector_bundle import atiyah_bundle
sage: F.<x> = FunctionField(GF(7))
sage: R.<y> = F[]
sage: K.<y> = F.extension(y^2 - x^3 - x)
sage: V = atiyah_bundle(K, 2, 0)
sage: W = V.direct_sum_repeat(2)
sage: T = matrix(K, 4, 4,
....: [x+1, 4, 3*x+4, 6*x+6,
....: 4*x+5, 6*x+5, 2*x+1, 4*x+3,
....: 3*x+1, 5*x, 3*x+1, x+6,
....: 5*x+4, 6*x, 5*x+2, 6*x+6])
sage: W = W.apply_isomorphism(T)
sage: ind, s, isom = W._indecomposable_power_split()
sage: s
2
sage: V._isomorphism_indecomposable(ind) is not None
True
"""
End = self.end()
A, to_A, from_A = End.global_algebra()
(S, to_S, from_S), factors = End._global_algebra_split()
if len(factors) > 1:
return None, None
split = factors[0]
K = split._K
s = split._n
idems = [diagonal_matrix(K, [0]*i + [1] + [0]*(s-i-1))
for i in range(s)]
idems = [from_A(from_S(split.from_M(idem))) for idem in idems]
images = [End.image(idem) for idem in idems]
factor = images[0][0]
isoms = [factor._isomorphism_indecomposable(im[0]) for im in images[1:]]
phis = ([images[0][1]]
+ [im[1]*isom for im, isom in zip(images[1:], isoms)])
return (factor,
len(phis),
block_matrix(self._function_field, [phis]))
def isomorphism_to(self, other):
r"""
Return an isomorphism from self to other if it exists and None otherwise
@ -1148,4 +1260,5 @@ class VectorBundle(SageObject):
return VectorBundle(self._function_field,
self._ideals,
isom * self._g_finite,
isom * self._g_infinite)
isom * self._g_infinite,
check=False)