-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathdifferent.py
416 lines (317 loc) · 13.9 KB
/
different.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
r"""
Computing the different of a finite extension of discretely valued fields
=========================================================================
Let `(K, v)` be a discretely valued field, `L/K` a finite separable field
extension and `w` an extension of `v` to `L` (which means that `w|_K=v`).
We assume that `v` and `w` are nontrivial and take values in
`\mathbb{Q}\cup\{\infty\}`; we do not assume that either of `v` and `w` are
normalized.
Let `\mathcal{o}_v` denote the valuation ring of `v` and `\mathcal{O}`
the integral closure of `\mathcal{o}_v` in `L`. Then `\mathcal{O}/\mathcal{o}_v`
is a finite ring extension, `\mathcal{O}` is a semilocal ring whose maximale
ideals correspond to the extensions of `v` to `L`, and the valuation ring
`\mathcal{O}_w` of `w` is the localization of `\mathcal{O}` at the maximal ideal
corresponding to `w`.
The different `\mathfrak{D}_{\mathcal{O}/\mathcal{o}_v}` is an ideal of
`\mathcal{O}`, defined in [Neukirch]_, Definition III.2.1. The *different* of
the extension `w/v` is defined as the ideal
MATH::
\mathfrak{D}_{w/v} := \mathfrak{D}_{\mathcal{O}/\mathcal{o}_v}\mathcal{O}_w
of `\mathcal{O}_w`. Frequently, we will mean by the different of `w/v` the
*valuation* of the different,
MATH::
\mathfrak{d}_{w/v} := w(\mathfrak{D}_{w/v}).
Note that this value depends on how the valuation `v` is normalized!
In this module we implement a function ``different_of_finite_extension(v, w)``
which computes the valuation of the different `\mathfrak{d}_{w/v}`, assuming
that `L`, the domain of `w`, is given as a simple extension of `K`, the domain
of `v`. The computation uses mostly valuation theory, in particular the methods of
MacLane [MacLane]_, implemented into Sage by Julian Rüth [Rueth]_.
REFERENCES:
.. [MacLane] \S. MacLane, *A construction for absolute values in polynomial rings.*
Trans. Amer. Math. Soc, 40(3):363–395, 1936.
.. [Neukirch] \J. Neukirch, *Algebraische Zahlentheorie*
.. [Rueth] \J. Rüth, *Models of curves and valuations.* PhD thesis, Universität Ulm, 2015.
EXAMPLES:
sage: from regular_models.different import different_of_finite_extension
sage: R.<x> = QQ[]
sage: L.<alpha> = QQ.extension(x^4+2*x^2+4)
sage: v = QQ.valuation(2)
sage: w = v.extension(L)
sage: different_of_finite_extension(v, w)
3/2
sage: K.<x> = FunctionField(QQ)
sage: v = K.valuation(GaussValuation(K._ring, QQ.valuation(2)))
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^4 + 2*x^2*y^2 + x^2)
sage: w = v.extension(L)
sage: different_of_finite_extension(v, w)
2
AUTHORS:
- Stefan Wewers (2020): initial version
"""
# *****************************************************************************
# Copyright (C) 2020 Stefan Wewers <[email protected]>
#
# This program is free software: you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation, either version 2 of the License, or
# (at your option) any later version.
# https://www.gnu.org/licenses/
# *****************************************************************************
def different_of_finite_extension(v, w):
r""" Return the valuation of the different of this extension of valued fields.
INPUT:
- v -- a discrete valuation on a field `K`
- w -- an extension of `v` to a finite separable field extension `L/K`
OUTPUT: the valuation of the different, `w(\mathcal{D}_{w/v})`.
.. NOTE::
The valuation `w` may be given simply as a map from L to the usual
codomain of a discrete valuation. For instance, `w` may be realized
as a composition of a ring homomorphism from `L` to some field `M` and
a discrete valuation on `M`.
If `w` is internally realized as a pseudovaluation on the polynomial ring
`K[x]`, then the valualtion `v` is actually not needed as a parameter.
But if it is not, `v` must be given as an actual discrete valuation; in
particular, it needs to have a method ``mac_lane_approximants``.
"""
if hasattr(w, "_from_base"):
# if w is induced by an automorphism of its domain, replace w by its base valuation
w = w._base_valuation
K = v.domain()
L = w.domain()
assert K == L.base_field(), "the domain of v must be the base field of the domain of w"
# check that L/K is finite and separable!
if hasattr(L, 'defining_polynomial'):
f = L.defining_polynomial()
else:
f = L.polynomial()
assert f.derivative() != 0, "L/K must be separable"
wt = _good_approximation(w, v)
phi = wt.phi()
L = K.extension(phi, 'beta')
w = v.extension(L)
# We replace L by the extension generated by phi.
# Now there is a unique extension of v to L/K, and its completion
# is isomorphic to the completion of the original extension L at w.
# Hence we can use the new extension L/K to compute the different of w/v.
k_w = w.residue_field()
e = w.value_group().index(v.value_group())
if e == L.degree() or k_w.characteristic() == 0 or k_w.is_finite():
# If e = deg(g) then k_w/k_v is trivial.
# If k_w has characteristik 0 or is finite, then k_w/k_v is separable.
# We can then find a generator t of OO_w over OO_v, using the trick
# from the proof of [Neukirch], Lemma II.10.4.
pi = w.uniformizer()
tb = k_w.gen()
gb = _minimal_polynomial(tb, v.residue_field())
g = gb.map_coefficients(v.lift, K)
t = w.lift(tb)
if w(g(t)) == w(pi):
g = t.minpoly()
return w(g.derivative()(t))
else:
t = t + pi
g = t.minpoly()
return w(g.derivative()(t))
else:
b = _integral_bases(v, w)
n = len(b)
from sage.matrix.constructor import matrix
A = matrix(K, n)
for i in range(n):
for j in range(n):
A[i, j] = (b[i]*b[j]).trace()
return v(A.det())/n
def _integral_bases(v, w):
r""" Return an integral basis of OO_w over OO_v.
INPUT:
- ``v`` -- a discrete valuation on a field `K`
- ``w`` -- an extension of `v` to a finite separable field extension `L` of `K`
It is assumed that `w` is the only extension of `v` to `L`.
OUTPUT: a list of elements of `L` which is an integral basis of the valuation
ring `\mathcal{O}_w`, as a free `\mathcal{O}_v`-module.
"""
K = v.domain()
L = w.domain()
assert K is L.base_field()
n = L.degree()
V, from_V, to_V = L.free_module()
# our assumptions imply that OO_w is generated over OO_v by pi and t:
pi = w.uniformizer()
kw = w.residue_field()
if hasattr(kw, "primitive_element"):
tb = kw.primitive_element()
else:
tb = kw.gen()
t = w.lift(tb)
generators = []
for i in range(n):
for j in range(n):
generators.append(to_V(pi**i * t**j))
from regular_models.RR_spaces.lattices import DVR_Lattice
M = DVR_Lattice(v, generators)
basis = M.basis()
return [from_V(m) for m in basis]
def _good_approximation(w, v):
r""" Return a good approximation of this extension of a discrete valuation.
INPUT:
- ``w`` -- a discrete valuation on a field `L`, which is a finite extension
of a field `K`
- ``v`` -- the restriction of `w` to `K`
OUTPUT: a discrete valuation on the polynomial ring over `K`, which
is a 'good' approximation of `w`.
Let `\alpha` denote the generator of the extension `L/K`, and `f\in K[x]`
the minimal polynomial of `\alpha`. Also, let `v` denote the restriction
of `w` to `K`. The finitely many extension of `v` to `L` correspond to the
irreducible factors of `f` over the completion of `K` with respect to `v`.
Write `g` for the factor corresponding to the given extension `w`.
An *approximation* of `w` is a discrete valuation `\tilde{w}` on `K[x]` with
the following properties:
- the restriction of `\tilde{w}` to `K` is equal to `v`
- `\tilde{w}(x)\geq 0`
- `\tilde{w}(h) \leq w(h)`, for all `h\in K[x]`
An approximation `\tilde{w}` is called *selective* if it is not an
approximation for any other extension of `v` to `L` than `w`.
Let `\tilde{w}` be a selective approximation of `w`. It is of the form
MATH::
\tilde{w} = [w_0, w_1(\phi_1)=\lambda_1,\ldots,w_n(\phi_n)=\lambda_n].
Write `\phi:=\phi_n`. The approximation `\tilde{w}` is called *good* if
- `\phi` has the same degree as the factor `g` of `f`, and
- `\phi` has a root `\beta` (in the algebraic closure of the completion of
`K`) which is closer to `\alpha` then any root of `f` different from `\alpha`
It then follows from Krasner's Lemma that `\phi` generates the same extension
of the completion of `K` than `g`, i.e.
MATH::
\hat{K}[\alpha]\cong\hat{K}[\beta].
"""
from sage.geometry.newton_polygon import NewtonPolygon
L = w.domain()
alpha = L.gen()
f = alpha.minpoly()
x = f.parent().gen()
F = f(alpha + x).shift(-1)
np = NewtonPolygon([(i, w(F[i])) for i in range(F.degree()+1)])
mu = -np.slopes()[0]
wt = _some_approximation(w, v)
f = wt.domain()(f)
assert hasattr(wt, "phi"), "wt = {}, L = {}".format(wt, L)
while w(wt.phi()(alpha)) <= mu*wt.phi().degree():
wt = wt.mac_lane_step(f)[0]
"""
try:
wt = wt.mac_lane_step(f)[0]
except:
print("wt = ", wt)
print("domain = ", wt.domain())
raise ValueError()
"""
return wt
def _some_approximation(w, v):
r""" Return an approximation of this extension of a discrete valuation.
INPUT:
- ``w`` -- a discrete valuation on a field `L`, which is a finite extension
of a field `K`
- ``v`` -- the restriction of `w` to `K`
OUTPUT: a discrete valuation on the polynomial ring over `K`, which
is a selective approximation of `w` (see the docstring for
``good_approximation`` for the definition of these terms).
"""
# very likely, w is implemented as a limit valuation, in which case we
# already have access to a selective approximation
if hasattr(w, '_base_valuation'):
wb = w._base_valuation
if hasattr(wb, '_approximation'):
return wb._approximation
else:
pass
if hasattr(wb, '_base_valuation'):
wb = wb._base_valuation
if hasattr(wb, "_approximation"):
return wb._approximation
# if the above does not work, e.g. if w is defined as a composition of a
# valuation and a ring homomorphism, we compute selective approximations for all
# extensions of v=w|_K to L and chose the one which is an approximation of w
L = w.domain()
alpha = L.gen()
f = alpha.minpoly()
V = v.mac_lane_approximants(f, require_incomparability=True, require_maximal_degree=True)
for wt in V:
if _is_approximation(wt, w):
return wt
def _is_approximation(wt, w):
r""" Check whether `\tilde{w}` is an approximation of `w`.
INPUT:
``wt`` -- a discrete valuation on a polynomial ring over a field `K`
``w`` -- a discrete valuation on a finite field extension `L` of `K`
OUTPUT: ``True`` if `\tilde{w}` is an approximation of `w`, ``False`` otherwise.
We *assume* that the restrictions of `w` and `\tilde{w}` to `K` are equal;
this is not checked!
We also assume that `L=K[\alpha]` is a simple extension, with generator
`\alpha`, which is `w`-integral, i.e. `w(\alpha)\geq 0`. Then we may view
`w` as a pseudovaluation on `K[x]` via the evaluation map
MATH::
K[x] \to L, \quad h \mapsto h(\alpha).
By definition, `\wt` is an approximation of `w` if the following conditions
hold:
- the restrictions of `w` and `\tilde{w}` to `K` are equal,
- `\tilde{w}(x) \geq 0`,
- `\tilde{w}(g) \leq w(h)`, for all `g\in K[x]`.
As already remarked, we assume the first condition. The second condition is
easy to check. To check the third condition, we write `\tilde{w}` as an
*inductive valuation* (see [MacLane]_, [Rueth]_),
MATH::
\tilde{w} = [w_0, w_1(\phi_1)=\lambda_1,\ldots,w_n(\phi_n)=\lambda_n].
Then by [Rueth]_, Theorem 5.65, the third condition holds if and only if
MATH::
`w(\phi_n) \geq \lambda_n`.
"""
R = wt.domain()
x = R.gen()
if wt(x) < 0:
return False
K = R.base_ring()
L = w.domain()
alpha = L.gen()
assert L.base_field() == K, "the domains of wt and w must have the same base fields"
return w(wt.phi()(alpha)) >= wt.mu()
def _minimal_polynomial(alpha, K):
r""" Return the minimal polynomial of `\alpha` over the field `k`.
INPUT:
- ``alpha`` -- an elemenent of a field `L`
- ``K`` -- a subfield of `L` over which `alpha` is algebraic
OUTPUT: the minimal polynomial of `\alpha` over `K`.
"""
L = alpha.parent()
assert K.is_subring(L), "K must be a subfield of L"
if L.base_ring() == K:
return alpha.minpoly()
else:
A = _matrix_of_extension_element(alpha, K)
f = A.minimal_polynomial()
assert f(alpha) == 0, "Error!"
return f
def _matrix_of_extension_element(alpha, K):
r""" Return the matrix corresponding to this element of the field extension.
"""
L = alpha.parent()
M = L.base_ring()
if M == K:
return alpha.matrix()
assert K.is_subring(M), "K must be a subring of M"
A = alpha.matrix()
n = A.nrows()
D = {}
for i in range(n):
for j in range(n):
D[(i, j)] = _matrix_of_extension_element(A[i, j], K)
m = D[(0, 0)].nrows()
N = n*m
from sage.matrix.constructor import matrix
B = matrix(K, N, N)
for i in range(n):
for j in range(n):
for k in range(m):
for ell in range(m):
B[i*m+k, j*m+ell] = D[(i, j)][k, ell]
return B