Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/schemes/plane_curves/affine_curve.py
4158 views
1
"""
2
Affine plane curves over a general ring
3
4
AUTHORS:
5
6
- William Stein (2005-11-13)
7
8
- David Joyner (2005-11-13)
9
10
- David Kohel (2006-01)
11
"""
12
13
#*****************************************************************************
14
# Copyright (C) 2005 William Stein <[email protected]>
15
#
16
# Distributed under the terms of the GNU General Public License (GPL)
17
#
18
# The full text of the GPL is available at:
19
#
20
# http://www.gnu.org/licenses/
21
#*****************************************************************************
22
23
from sage.interfaces.all import singular
24
25
from sage.misc.all import add
26
27
from sage.rings.all import degree_lowest_rational_function
28
29
from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing
30
31
from sage.schemes.generic.affine_space import is_AffineSpace
32
from sage.schemes.generic.algebraic_scheme import AlgebraicScheme_subscheme_affine
33
34
35
from curve import Curve_generic
36
37
class AffineSpaceCurve_generic(Curve_generic, AlgebraicScheme_subscheme_affine):
38
def _repr_type(self):
39
return "Affine Space"
40
41
def __init__(self, A, X):
42
if not is_AffineSpace(A):
43
raise TypeError, "A (=%s) must be an affine space"%A
44
Curve_generic.__init__(self, A, X)
45
d = self.dimension()
46
if d != 1:
47
raise ValueError, "defining equations (=%s) define a scheme of dimension %s != 1"%(X,d)
48
49
class AffineCurve_generic(Curve_generic):
50
def __init__(self, A, f):
51
P = f.parent()
52
if not (is_AffineSpace(A) and A.dimension != 2):
53
raise TypeError, "Argument A (= %s) must be an affine plane."%A
54
Curve_generic.__init__(self, A, [f])
55
56
def _repr_type(self):
57
return "Affine"
58
59
def divisor_of_function(self, r):
60
"""
61
Return the divisor of a function on a curve.
62
63
INPUT: r is a rational function on X
64
65
OUTPUT:
66
67
68
- ``list`` - The divisor of r represented as a list of
69
coefficients and points. (TODO: This will change to a more
70
structural output in the future.)
71
72
73
EXAMPLES::
74
75
sage: F = GF(5)
76
sage: P2 = AffineSpace(2, F, names = 'xy')
77
sage: R = P2.coordinate_ring()
78
sage: x, y = R.gens()
79
sage: f = y^2 - x^9 - x
80
sage: C = Curve(f)
81
sage: K = FractionField(R)
82
sage: r = 1/x
83
sage: C.divisor_of_function(r) # todo: not implemented (broken)
84
[[-1, (0, 0, 1)]]
85
sage: r = 1/x^3
86
sage: C.divisor_of_function(r) # todo: not implemented (broken)
87
[[-3, (0, 0, 1)]]
88
"""
89
F = self.base_ring()
90
f = self.defining_polynomial()
91
pts = self.places_on_curve()
92
numpts = len(pts)
93
R = f.parent()
94
x,y = R.gens()
95
R0 = PolynomialRing(F,3,names = [str(x),str(y),"t"])
96
vars0 = R0.gens()
97
t = vars0[2]
98
divf = []
99
for pt0 in pts:
100
if pt0[2] != F(0):
101
lcs = self.local_coordinates(pt0,5)
102
yt = lcs[1]
103
xt = lcs[0]
104
ldg = degree_lowest_rational_function(r(xt,yt),t)
105
if ldg[0] != 0:
106
divf.append([ldg[0],pt0])
107
return divf
108
109
def local_coordinates(self, pt, n):
110
r"""
111
Return local coordinates to precision n at the given point.
112
113
Behaviour is flaky - some choices of `n` are worst that
114
others.
115
116
117
INPUT:
118
119
120
- ``pt`` - an F-rational point on X which is not a
121
point of ramification for the projection (x,y) - x.
122
123
- ``n`` - the number of terms desired
124
125
126
OUTPUT: x = x0 + t y = y0 + power series in t
127
128
EXAMPLES::
129
130
sage: F = GF(5)
131
sage: pt = (2,3)
132
sage: R = PolynomialRing(F,2, names = ['x','y'])
133
sage: x,y = R.gens()
134
sage: f = y^2-x^9-x
135
sage: C = Curve(f)
136
sage: C.local_coordinates(pt, 9)
137
[t + 2, -2*t^12 - 2*t^11 + 2*t^9 + t^8 - 2*t^7 - 2*t^6 - 2*t^4 + t^3 - 2*t^2 - 2]
138
"""
139
f = self.defining_polynomial()
140
R = f.parent()
141
F = self.base_ring()
142
p = F.characteristic()
143
x0 = F(pt[0])
144
y0 = F(pt[1])
145
astr = ["a"+str(i) for i in range(1,2*n)]
146
x,y = R.gens()
147
R0 = PolynomialRing(F,2*n+2,names = [str(x),str(y),"t"]+astr)
148
vars0 = R0.gens()
149
t = vars0[2]
150
yt = y0*t**0+add([vars0[i]*t**(i-2) for i in range(3,2*n+2)])
151
xt = x0+t
152
ft = f(xt,yt)
153
S = singular
154
S.eval('ring s = '+str(p)+','+str(R0.gens())+',lp;')
155
S.eval('poly f = '+str(ft) + ';')
156
c = S('coeffs(%s, t)'%ft)
157
N = int(c.size())
158
b = ["%s[%s,1],"%(c.name(), i) for i in range(2,N/2-4)]
159
b = ''.join(b)
160
b = b[:len(b)-1] # to cut off the trailing comma
161
cmd = 'ideal I = '+b
162
S.eval(cmd)
163
S.eval('short=0') # print using *'s and ^'s.
164
c = S.eval('slimgb(I)')
165
d = c.split("=")
166
d = d[1:]
167
d[len(d)-1] += "\n"
168
e = [x[:x.index("\n")] for x in d]
169
vals = []
170
for x in e:
171
for y in vars0:
172
if str(y) in x:
173
if len(x.replace(str(y),"")) != 0:
174
i = x.find("-")
175
if i>0:
176
vals.append([eval(x[1:i]),x[:i],F(eval(x[i+1:]))])
177
i = x.find("+")
178
if i>0:
179
vals.append([eval(x[1:i]),x[:i],-F(eval(x[i+1:]))])
180
else:
181
vals.append([eval(str(y)[1:]),str(y),F(0)])
182
vals.sort()
183
k = len(vals)
184
v = [x0+t,y0+add([vals[i][2]*t**(i+1) for i in range(k)])]
185
return v
186
187
def plot(self, *args, **kwds):
188
"""
189
Plot the real points on this affine plane curve.
190
191
INPUT:
192
193
194
- ``self`` - an affine plane curve
195
196
- ``*args`` - optional tuples (variable, minimum, maximum) for
197
plotting dimensions
198
199
- ``**kwds`` - optional keyword arguments passed on to
200
``implicit_plot``
201
202
203
EXAMPLES:
204
205
A cuspidal curve::
206
207
sage: R.<x, y> = QQ[]
208
sage: C = Curve(x^3 - y^2)
209
sage: C.plot()
210
211
A 5-nodal curve of degree 11. This example also illustrates
212
some of the optional arguments::
213
214
sage: R.<x, y> = ZZ[]
215
sage: C = Curve(32*x^2 - 2097152*y^11 + 1441792*y^9 - 360448*y^7 + 39424*y^5 - 1760*y^3 + 22*y - 1)
216
sage: C.plot((x, -1, 1), (y, -1, 1), plot_points=400)
217
218
A line over `\mathbf{RR}`::
219
220
sage: R.<x, y> = RR[]
221
sage: C = Curve(R(y - sqrt(2)*x))
222
sage: C.plot()
223
"""
224
I = self.defining_ideal()
225
return I.plot(*args, **kwds)
226
227
class AffineCurve_finite_field(AffineCurve_generic):
228
def rational_points(self, algorithm="enum"):
229
r"""
230
Return sorted list of all rational points on this curve.
231
232
Use *very* naive point enumeration to find all rational points on
233
this curve over a finite field.
234
235
EXAMPLE::
236
237
sage: A, (x,y) = AffineSpace(2,GF(9,'a')).objgens()
238
sage: C = Curve(x^2 + y^2 - 1)
239
sage: C
240
Affine Curve over Finite Field in a of size 3^2 defined by x0^2 + x1^2 - 1
241
sage: C.rational_points()
242
[(0, 1), (0, 2), (1, 0), (2, 0), (a + 1, a + 1), (a + 1, 2*a + 2), (2*a + 2, a + 1), (2*a + 2, 2*a + 2)]
243
"""
244
f = self.defining_polynomial()
245
R = f.parent()
246
K = R.base_ring()
247
points = []
248
for x in K:
249
for y in K:
250
if f(x,y) == 0:
251
points.append(self((x,y)))
252
points.sort()
253
return points
254
255
256
class AffineCurve_prime_finite_field(AffineCurve_finite_field):
257
# CHECK WHAT ASSUMPTIONS ARE MADE REGARDING AFFINE VS. PROJECTIVE MODELS!!!
258
# THIS IS VERY DIRTY STILL -- NO DATASTRUCTURES FOR DIVISORS.
259
260
def riemann_roch_basis(self,D):
261
"""
262
Interfaces with Singular's BrillNoether command.
263
264
INPUT:
265
266
267
- ``self`` - a plane curve defined by a polynomial eqn f(x,y)
268
= 0 over a prime finite field F = GF(p) in 2 variables x,y
269
representing a curve X: f(x,y) = 0 having n F-rational
270
points (see the Sage function places_on_curve)
271
272
- ``D`` - an n-tuple of integers
273
`(d1, ..., dn)` representing the divisor
274
`Div = d1*P1+...+dn*Pn`, where
275
`X(F) = \{P1,...,Pn\}`.
276
*The ordering is that dictated by places_on_curve.*
277
278
279
OUTPUT: basis of L(Div)
280
281
EXAMPLE::
282
283
sage: R = PolynomialRing(GF(5),2,names = ["x","y"])
284
sage: x, y = R.gens()
285
sage: f = y^2 - x^9 - x
286
sage: C = Curve(f)
287
sage: D = [6,0,0,0,0,0]
288
sage: C.riemann_roch_basis(D)
289
[1, (y^2*z^4 - x*z^5)/x^6, (y^2*z^5 - x*z^6)/x^7, (y^2*z^6 - x*z^7)/x^8]
290
"""
291
f = self.defining_polynomial()
292
R = f.parent()
293
F = self.base_ring()
294
p = F.characteristic()
295
Dstr = str(tuple(D))
296
G = singular(','.join([str(x) for x in D]), type='intvec')
297
298
singular.LIB('brnoeth.lib')
299
300
S = singular.ring(p, R.gens(), 'lp')
301
fsing = singular(str(f))
302
303
X = fsing.Adj_div()
304
P = singular.NSplaces(1, X)
305
T = P[1][2]
306
T.set_ring()
307
308
LG = G.BrillNoether(P)
309
310
dim = len(LG)
311
basis = [(LG[i][1], LG[i][2]) for i in range(1,dim+1)]
312
x, y, z = PolynomialRing(F, 3, names = ["x","y","z"]).gens()
313
V = []
314
for g in basis:
315
T.set_ring() # necessary...
316
V.append(eval(g[0].sage_polystring())/eval(g[1].sage_polystring()))
317
return V
318
319
320
def rational_points(self, algorithm="enum"):
321
r"""
322
Return sorted list of all rational points on this curve.
323
324
INPUT:
325
326
327
- ``algorithm`` - string:
328
329
+ ``'enum'`` - straightforward enumeration
330
331
+ ``'bn'`` - via Singular's Brill-Noether package.
332
333
+ ``'all'`` - use all implemented algorithms and
334
verify that they give the same answer, then return it
335
336
337
.. note::
338
339
The Brill-Noether package does not always work. When it
340
fails a RuntimeError exception is raised.
341
342
EXAMPLE::
343
344
sage: x, y = (GF(5)['x,y']).gens()
345
sage: f = y^2 - x^9 - x
346
sage: C = Curve(f); C
347
Affine Curve over Finite Field of size 5 defined by -x^9 + y^2 - x
348
sage: C.rational_points(algorithm='bn')
349
[(0, 0), (2, 2), (2, 3), (3, 1), (3, 4)]
350
sage: C = Curve(x - y + 1)
351
sage: C.rational_points()
352
[(0, 1), (1, 2), (2, 3), (3, 4), (4, 0)]
353
354
We compare Brill-Noether and enumeration::
355
356
sage: x, y = (GF(17)['x,y']).gens()
357
sage: C = Curve(x^2 + y^5 + x*y - 19)
358
sage: v = C.rational_points(algorithm='bn')
359
sage: w = C.rational_points(algorithm='enum')
360
sage: len(v)
361
20
362
sage: v == w
363
True
364
"""
365
if algorithm == "enum":
366
367
return AffineCurve_finite_field.rational_points(self, algorithm="enum")
368
369
elif algorithm == "bn":
370
f = self.defining_polynomial()._singular_()
371
singular = f.parent()
372
singular.lib('brnoeth')
373
try:
374
X1 = f.Adj_div()
375
except (TypeError, RuntimeError), s:
376
raise RuntimeError, str(s) + "\n\n ** Unable to use the Brill-Noether Singular package to compute all points (see above)."
377
378
X2 = singular.NSplaces(1, X1)
379
R = X2[5][1][1]
380
singular.set_ring(R)
381
382
# We use sage_flattened_str_list since iterating through
383
# the entire list through the sage/singular interface directly
384
# would involve hundreds of calls to singular, and timing issues
385
# with the expect interface could crop up. Also, this is vastly
386
# faster (and more robust).
387
v = singular('POINTS').sage_flattened_str_list()
388
pnts = [self(int(v[3*i]), int(v[3*i+1])) for i in range(len(v)/3) if int(v[3*i+2])!=0]
389
# remove multiple points
390
pnts = list(set(pnts))
391
pnts.sort()
392
return pnts
393
394
elif algorithm == "all":
395
396
S_enum = self.rational_points(algorithm = "enum")
397
S_bn = self.rational_points(algorithm = "bn")
398
if S_enum != S_bn:
399
raise RuntimeError, "Bug in rational_points -- different algorithms give different answers for curve %s!"%self
400
return S_enum
401
402
else:
403
raise ValueError, "No algorithm '%s' known"%algorithm
404
405