Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagesmc
Path: blob/master/src/sage/schemes/affine/affine_rational_point.py
8820 views
1
r"""
2
Enumeration of rational points on affine schemes
3
4
Naive algorithms for enumerating rational points over `\QQ` or finite fields
5
over for general schemes.
6
7
.. WARNING::
8
9
Incorrect results and infinite loops may occur if using a wrong function.
10
(For instance using an affine function for a projective scheme or a finite
11
field function for a scheme defined over an infinite field.)
12
13
EXAMPLES:
14
15
Affine, over `\QQ`::
16
17
sage: from sage.schemes.affine.affine_rational_point import enum_affine_rational_field
18
sage: A.<x,y,z> = AffineSpace(3,QQ)
19
sage: S = A.subscheme([2*x-3*y])
20
sage: enum_affine_rational_field(S,2)
21
[(0, 0, -2), (0, 0, -1), (0, 0, -1/2), (0, 0, 0),
22
(0, 0, 1/2), (0, 0, 1), (0, 0, 2)]
23
24
Affine over a finite field::
25
26
sage: from sage.schemes.affine.affine_rational_point import enum_affine_finite_field
27
sage: A.<w,x,y,z> = AffineSpace(4,GF(2))
28
sage: enum_affine_finite_field(A(GF(2)))
29
[(0, 0, 0, 0), (0, 0, 0, 1), (0, 0, 1, 0), (0, 0, 1, 1), (0, 1, 0, 0),
30
(0, 1, 0, 1), (0, 1, 1, 0), (0, 1, 1, 1), (1, 0, 0, 0), (1, 0, 0, 1),
31
(1, 0, 1, 0), (1, 0, 1, 1), (1, 1, 0, 0), (1, 1, 0, 1), (1, 1, 1, 0),
32
(1, 1, 1, 1)]
33
34
AUTHORS:
35
36
- David R. Kohel <[email protected]>: original version.
37
38
- John Cremona and Charlie Turner <[email protected]> (06-2010):
39
improvements to clarity and documentation.
40
"""
41
42
43
#*******************************************************************************
44
# Copyright (C) 2010 William Stein, David Kohel, John Cremona, Charlie Turner
45
# Distributed under the terms of the GNU General Public License (GPL)
46
# http://www.gnu.org/licenses/
47
#*******************************************************************************
48
49
50
from sage.rings.all import ZZ
51
from sage.misc.all import srange, cartesian_product_iterator
52
from sage.schemes.generic.scheme import is_Scheme
53
54
55
def enum_affine_rational_field(X,B):
56
"""
57
Enumerates affine rational points on scheme ``X`` (defined over `\QQ`) up
58
to bound ``B``.
59
60
INPUT:
61
62
- ``X`` - a scheme or set of abstract rational points of a scheme;
63
- ``B`` - a positive integer bound.
64
65
OUTPUT:
66
67
- a list containing the affine points of ``X`` of height up to ``B``,
68
sorted.
69
70
EXAMPLES::
71
72
sage: A.<x,y,z> = AffineSpace(3,QQ)
73
sage: from sage.schemes.affine.affine_rational_point import enum_affine_rational_field
74
sage: enum_affine_rational_field(A(QQ),1)
75
[(-1, -1, -1), (-1, -1, 0), (-1, -1, 1), (-1, 0, -1), (-1, 0, 0), (-1, 0, 1),
76
(-1, 1, -1), (-1, 1, 0), (-1, 1, 1), (0, -1, -1), (0, -1, 0), (0, -1, 1),
77
(0, 0, -1), (0, 0, 0), (0, 0, 1), (0, 1, -1), (0, 1, 0), (0, 1, 1), (1, -1, -1),
78
(1, -1, 0), (1, -1, 1), (1, 0, -1), (1, 0, 0), (1, 0, 1), (1, 1, -1), (1, 1, 0),
79
(1, 1, 1)]
80
81
::
82
83
sage: A.<w,x,y,z> = AffineSpace(4,QQ)
84
sage: S = A.subscheme([x^2-y*z+3,w^3+z+y^2])
85
sage: enum_affine_rational_field(S(QQ),2)
86
[]
87
sage: enum_affine_rational_field(S(QQ),3)
88
[(-2, 0, -3, -1)]
89
90
::
91
92
sage: A.<x,y> = AffineSpace(2,QQ)
93
sage: C = Curve(x^2+y-x)
94
sage: enum_affine_rational_field(C,10)
95
[(-2, -6), (-1, -2), (0, 0), (1, 0), (2, -2), (3, -6)]
96
97
98
AUTHORS:
99
100
- David R. Kohel <[email protected]>: original version.
101
102
- Charlie Turner (06-2010): small adjustments.
103
"""
104
if is_Scheme(X):
105
X = X(X.base_ring())
106
n = X.codomain().ambient_space().ngens()
107
if X.value_ring() is ZZ:
108
Q = [ 1 ]
109
else: # rational field
110
Q = range(1, B + 1)
111
R = [ 0 ] + [ s*k for k in range(1, B+1) for s in [1, -1] ]
112
pts = []
113
P = [0] * n
114
m = ZZ(0)
115
try:
116
pts.append(X(P))
117
except TypeError:
118
pass
119
iters = [ iter(R) for _ in range(n) ]
120
for it in iters:
121
it.next()
122
i = 0
123
while i < n:
124
try:
125
a = ZZ(iters[i].next())
126
except StopIteration:
127
iters[i] = iter(R) # reset
128
P[i] = iters[i].next() # reset P[i] to 0 and increment
129
i += 1
130
continue
131
m = m.gcd(a)
132
P[i] = a
133
for b in Q:
134
if m.gcd(b) == 1:
135
try:
136
pts.append(X([ num/b for num in P ]))
137
except TypeError:
138
pass
139
i = 0
140
m = ZZ(0)
141
pts.sort()
142
return pts
143
144
145
146
def enum_affine_finite_field(X):
147
r"""
148
Enumerates affine points on scheme ``X`` defined over a finite field.
149
150
INPUT:
151
152
- ``X`` - a scheme defined over a finite field or a set of abstract
153
rational points of such a scheme.
154
155
OUTPUT:
156
157
- a list containing the affine points of ``X`` over the finite field,
158
sorted.
159
160
EXAMPLES::
161
162
sage: F = GF(7)
163
sage: A.<w,x,y,z> = AffineSpace(4,F)
164
sage: C = A.subscheme([w^2+x+4,y*z*x-6,z*y+w*x])
165
sage: from sage.schemes.affine.affine_rational_point import enum_affine_finite_field
166
sage: enum_affine_finite_field(C(F))
167
[]
168
sage: C = A.subscheme([w^2+x+4,y*z*x-6])
169
sage: enum_affine_finite_field(C(F))
170
[(0, 3, 1, 2), (0, 3, 2, 1), (0, 3, 3, 3), (0, 3, 4, 4), (0, 3, 5, 6),
171
(0, 3, 6, 5), (1, 2, 1, 3), (1, 2, 2, 5), (1, 2, 3, 1), (1, 2, 4, 6),
172
(1, 2, 5, 2), (1, 2, 6, 4), (2, 6, 1, 1), (2, 6, 2, 4), (2, 6, 3, 5),
173
(2, 6, 4, 2), (2, 6, 5, 3), (2, 6, 6, 6), (3, 1, 1, 6), (3, 1, 2, 3),
174
(3, 1, 3, 2), (3, 1, 4, 5), (3, 1, 5, 4), (3, 1, 6, 1), (4, 1, 1, 6),
175
(4, 1, 2, 3), (4, 1, 3, 2), (4, 1, 4, 5), (4, 1, 5, 4), (4, 1, 6, 1),
176
(5, 6, 1, 1), (5, 6, 2, 4), (5, 6, 3, 5), (5, 6, 4, 2), (5, 6, 5, 3),
177
(5, 6, 6, 6), (6, 2, 1, 3), (6, 2, 2, 5), (6, 2, 3, 1), (6, 2, 4, 6),
178
(6, 2, 5, 2), (6, 2, 6, 4)]
179
180
::
181
182
sage: A.<x,y,z> = AffineSpace(3,GF(3))
183
sage: S = A.subscheme(x+y)
184
sage: enum_affine_finite_field(S)
185
[(0, 0, 0), (0, 0, 1), (0, 0, 2), (1, 2, 0), (1, 2, 1), (1, 2, 2),
186
(2, 1, 0), (2, 1, 1), (2, 1, 2)]
187
188
ALGORITHM:
189
190
Checks all points in affine space to see if they lie on X.
191
192
.. WARNING::
193
194
If ``X`` is defined over an infinite field, this code will not finish!
195
196
AUTHORS:
197
198
- John Cremona and Charlie Turner (06-2010)
199
"""
200
if is_Scheme(X):
201
X = X(X.base_ring())
202
n = X.codomain().ambient_space().ngens()
203
F = X.value_ring()
204
pts = []
205
for c in cartesian_product_iterator([F]*n):
206
try:
207
pts.append(X(c))
208
except Exception:
209
pass
210
pts.sort()
211
return pts
212
213