Path: blob/master/sage/schemes/plane_curves/affine_curve.py
4158 views
"""1Affine plane curves over a general ring23AUTHORS:45- William Stein (2005-11-13)67- David Joyner (2005-11-13)89- David Kohel (2006-01)10"""1112#*****************************************************************************13# Copyright (C) 2005 William Stein <[email protected]>14#15# Distributed under the terms of the GNU General Public License (GPL)16#17# The full text of the GPL is available at:18#19# http://www.gnu.org/licenses/20#*****************************************************************************2122from sage.interfaces.all import singular2324from sage.misc.all import add2526from sage.rings.all import degree_lowest_rational_function2728from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing2930from sage.schemes.generic.affine_space import is_AffineSpace31from sage.schemes.generic.algebraic_scheme import AlgebraicScheme_subscheme_affine323334from curve import Curve_generic3536class AffineSpaceCurve_generic(Curve_generic, AlgebraicScheme_subscheme_affine):37def _repr_type(self):38return "Affine Space"3940def __init__(self, A, X):41if not is_AffineSpace(A):42raise TypeError, "A (=%s) must be an affine space"%A43Curve_generic.__init__(self, A, X)44d = self.dimension()45if d != 1:46raise ValueError, "defining equations (=%s) define a scheme of dimension %s != 1"%(X,d)4748class AffineCurve_generic(Curve_generic):49def __init__(self, A, f):50P = f.parent()51if not (is_AffineSpace(A) and A.dimension != 2):52raise TypeError, "Argument A (= %s) must be an affine plane."%A53Curve_generic.__init__(self, A, [f])5455def _repr_type(self):56return "Affine"5758def divisor_of_function(self, r):59"""60Return the divisor of a function on a curve.6162INPUT: r is a rational function on X6364OUTPUT:656667- ``list`` - The divisor of r represented as a list of68coefficients and points. (TODO: This will change to a more69structural output in the future.)707172EXAMPLES::7374sage: F = GF(5)75sage: P2 = AffineSpace(2, F, names = 'xy')76sage: R = P2.coordinate_ring()77sage: x, y = R.gens()78sage: f = y^2 - x^9 - x79sage: C = Curve(f)80sage: K = FractionField(R)81sage: r = 1/x82sage: C.divisor_of_function(r) # todo: not implemented (broken)83[[-1, (0, 0, 1)]]84sage: r = 1/x^385sage: C.divisor_of_function(r) # todo: not implemented (broken)86[[-3, (0, 0, 1)]]87"""88F = self.base_ring()89f = self.defining_polynomial()90pts = self.places_on_curve()91numpts = len(pts)92R = f.parent()93x,y = R.gens()94R0 = PolynomialRing(F,3,names = [str(x),str(y),"t"])95vars0 = R0.gens()96t = vars0[2]97divf = []98for pt0 in pts:99if pt0[2] != F(0):100lcs = self.local_coordinates(pt0,5)101yt = lcs[1]102xt = lcs[0]103ldg = degree_lowest_rational_function(r(xt,yt),t)104if ldg[0] != 0:105divf.append([ldg[0],pt0])106return divf107108def local_coordinates(self, pt, n):109r"""110Return local coordinates to precision n at the given point.111112Behaviour is flaky - some choices of `n` are worst that113others.114115116INPUT:117118119- ``pt`` - an F-rational point on X which is not a120point of ramification for the projection (x,y) - x.121122- ``n`` - the number of terms desired123124125OUTPUT: x = x0 + t y = y0 + power series in t126127EXAMPLES::128129sage: F = GF(5)130sage: pt = (2,3)131sage: R = PolynomialRing(F,2, names = ['x','y'])132sage: x,y = R.gens()133sage: f = y^2-x^9-x134sage: C = Curve(f)135sage: C.local_coordinates(pt, 9)136[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]137"""138f = self.defining_polynomial()139R = f.parent()140F = self.base_ring()141p = F.characteristic()142x0 = F(pt[0])143y0 = F(pt[1])144astr = ["a"+str(i) for i in range(1,2*n)]145x,y = R.gens()146R0 = PolynomialRing(F,2*n+2,names = [str(x),str(y),"t"]+astr)147vars0 = R0.gens()148t = vars0[2]149yt = y0*t**0+add([vars0[i]*t**(i-2) for i in range(3,2*n+2)])150xt = x0+t151ft = f(xt,yt)152S = singular153S.eval('ring s = '+str(p)+','+str(R0.gens())+',lp;')154S.eval('poly f = '+str(ft) + ';')155c = S('coeffs(%s, t)'%ft)156N = int(c.size())157b = ["%s[%s,1],"%(c.name(), i) for i in range(2,N/2-4)]158b = ''.join(b)159b = b[:len(b)-1] # to cut off the trailing comma160cmd = 'ideal I = '+b161S.eval(cmd)162S.eval('short=0') # print using *'s and ^'s.163c = S.eval('slimgb(I)')164d = c.split("=")165d = d[1:]166d[len(d)-1] += "\n"167e = [x[:x.index("\n")] for x in d]168vals = []169for x in e:170for y in vars0:171if str(y) in x:172if len(x.replace(str(y),"")) != 0:173i = x.find("-")174if i>0:175vals.append([eval(x[1:i]),x[:i],F(eval(x[i+1:]))])176i = x.find("+")177if i>0:178vals.append([eval(x[1:i]),x[:i],-F(eval(x[i+1:]))])179else:180vals.append([eval(str(y)[1:]),str(y),F(0)])181vals.sort()182k = len(vals)183v = [x0+t,y0+add([vals[i][2]*t**(i+1) for i in range(k)])]184return v185186def plot(self, *args, **kwds):187"""188Plot the real points on this affine plane curve.189190INPUT:191192193- ``self`` - an affine plane curve194195- ``*args`` - optional tuples (variable, minimum, maximum) for196plotting dimensions197198- ``**kwds`` - optional keyword arguments passed on to199``implicit_plot``200201202EXAMPLES:203204A cuspidal curve::205206sage: R.<x, y> = QQ[]207sage: C = Curve(x^3 - y^2)208sage: C.plot()209210A 5-nodal curve of degree 11. This example also illustrates211some of the optional arguments::212213sage: R.<x, y> = ZZ[]214sage: C = Curve(32*x^2 - 2097152*y^11 + 1441792*y^9 - 360448*y^7 + 39424*y^5 - 1760*y^3 + 22*y - 1)215sage: C.plot((x, -1, 1), (y, -1, 1), plot_points=400)216217A line over `\mathbf{RR}`::218219sage: R.<x, y> = RR[]220sage: C = Curve(R(y - sqrt(2)*x))221sage: C.plot()222"""223I = self.defining_ideal()224return I.plot(*args, **kwds)225226class AffineCurve_finite_field(AffineCurve_generic):227def rational_points(self, algorithm="enum"):228r"""229Return sorted list of all rational points on this curve.230231Use *very* naive point enumeration to find all rational points on232this curve over a finite field.233234EXAMPLE::235236sage: A, (x,y) = AffineSpace(2,GF(9,'a')).objgens()237sage: C = Curve(x^2 + y^2 - 1)238sage: C239Affine Curve over Finite Field in a of size 3^2 defined by x0^2 + x1^2 - 1240sage: C.rational_points()241[(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)]242"""243f = self.defining_polynomial()244R = f.parent()245K = R.base_ring()246points = []247for x in K:248for y in K:249if f(x,y) == 0:250points.append(self((x,y)))251points.sort()252return points253254255class AffineCurve_prime_finite_field(AffineCurve_finite_field):256# CHECK WHAT ASSUMPTIONS ARE MADE REGARDING AFFINE VS. PROJECTIVE MODELS!!!257# THIS IS VERY DIRTY STILL -- NO DATASTRUCTURES FOR DIVISORS.258259def riemann_roch_basis(self,D):260"""261Interfaces with Singular's BrillNoether command.262263INPUT:264265266- ``self`` - a plane curve defined by a polynomial eqn f(x,y)267= 0 over a prime finite field F = GF(p) in 2 variables x,y268representing a curve X: f(x,y) = 0 having n F-rational269points (see the Sage function places_on_curve)270271- ``D`` - an n-tuple of integers272`(d1, ..., dn)` representing the divisor273`Div = d1*P1+...+dn*Pn`, where274`X(F) = \{P1,...,Pn\}`.275*The ordering is that dictated by places_on_curve.*276277278OUTPUT: basis of L(Div)279280EXAMPLE::281282sage: R = PolynomialRing(GF(5),2,names = ["x","y"])283sage: x, y = R.gens()284sage: f = y^2 - x^9 - x285sage: C = Curve(f)286sage: D = [6,0,0,0,0,0]287sage: C.riemann_roch_basis(D)288[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]289"""290f = self.defining_polynomial()291R = f.parent()292F = self.base_ring()293p = F.characteristic()294Dstr = str(tuple(D))295G = singular(','.join([str(x) for x in D]), type='intvec')296297singular.LIB('brnoeth.lib')298299S = singular.ring(p, R.gens(), 'lp')300fsing = singular(str(f))301302X = fsing.Adj_div()303P = singular.NSplaces(1, X)304T = P[1][2]305T.set_ring()306307LG = G.BrillNoether(P)308309dim = len(LG)310basis = [(LG[i][1], LG[i][2]) for i in range(1,dim+1)]311x, y, z = PolynomialRing(F, 3, names = ["x","y","z"]).gens()312V = []313for g in basis:314T.set_ring() # necessary...315V.append(eval(g[0].sage_polystring())/eval(g[1].sage_polystring()))316return V317318319def rational_points(self, algorithm="enum"):320r"""321Return sorted list of all rational points on this curve.322323INPUT:324325326- ``algorithm`` - string:327328+ ``'enum'`` - straightforward enumeration329330+ ``'bn'`` - via Singular's Brill-Noether package.331332+ ``'all'`` - use all implemented algorithms and333verify that they give the same answer, then return it334335336.. note::337338The Brill-Noether package does not always work. When it339fails a RuntimeError exception is raised.340341EXAMPLE::342343sage: x, y = (GF(5)['x,y']).gens()344sage: f = y^2 - x^9 - x345sage: C = Curve(f); C346Affine Curve over Finite Field of size 5 defined by -x^9 + y^2 - x347sage: C.rational_points(algorithm='bn')348[(0, 0), (2, 2), (2, 3), (3, 1), (3, 4)]349sage: C = Curve(x - y + 1)350sage: C.rational_points()351[(0, 1), (1, 2), (2, 3), (3, 4), (4, 0)]352353We compare Brill-Noether and enumeration::354355sage: x, y = (GF(17)['x,y']).gens()356sage: C = Curve(x^2 + y^5 + x*y - 19)357sage: v = C.rational_points(algorithm='bn')358sage: w = C.rational_points(algorithm='enum')359sage: len(v)36020361sage: v == w362True363"""364if algorithm == "enum":365366return AffineCurve_finite_field.rational_points(self, algorithm="enum")367368elif algorithm == "bn":369f = self.defining_polynomial()._singular_()370singular = f.parent()371singular.lib('brnoeth')372try:373X1 = f.Adj_div()374except (TypeError, RuntimeError), s:375raise RuntimeError, str(s) + "\n\n ** Unable to use the Brill-Noether Singular package to compute all points (see above)."376377X2 = singular.NSplaces(1, X1)378R = X2[5][1][1]379singular.set_ring(R)380381# We use sage_flattened_str_list since iterating through382# the entire list through the sage/singular interface directly383# would involve hundreds of calls to singular, and timing issues384# with the expect interface could crop up. Also, this is vastly385# faster (and more robust).386v = singular('POINTS').sage_flattened_str_list()387pnts = [self(int(v[3*i]), int(v[3*i+1])) for i in range(len(v)/3) if int(v[3*i+2])!=0]388# remove multiple points389pnts = list(set(pnts))390pnts.sort()391return pnts392393elif algorithm == "all":394395S_enum = self.rational_points(algorithm = "enum")396S_bn = self.rational_points(algorithm = "bn")397if S_enum != S_bn:398raise RuntimeError, "Bug in rational_points -- different algorithms give different answers for curve %s!"%self399return S_enum400401else:402raise ValueError, "No algorithm '%s' known"%algorithm403404405