Affine curves

Affine curves in Sage are curves in an affine space or an affine plane.

EXAMPLES:

We can construct curves in either an affine plane:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve([y - x^2], A); C
Affine Plane Curve over Rational Field defined by -x^2 + y
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([y - x**Integer(2)], A); C
Affine Plane Curve over Rational Field defined by -x^2 + y

or in higher dimensional affine space:

sage: A.<x,y,z,w> = AffineSpace(QQ, 4)
sage: C = Curve([y - x^2, z - w^3, w - y^4], A); C
Affine Curve over Rational Field defined by -x^2 + y, -w^3 + z, -y^4 + w
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(4), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = A._first_ngens(4)
>>> C = Curve([y - x**Integer(2), z - w**Integer(3), w - y**Integer(4)], A); C
Affine Curve over Rational Field defined by -x^2 + y, -w^3 + z, -y^4 + w

Integral affine curves over finite fields

If the curve is defined over a finite field and integral, that is reduced and irreducible, its function field is tightly coupled with the curve so that advanced computations based on Sage’s global function field machinery are available.

EXAMPLES:

sage: k.<a> = GF(2)
sage: A.<x,y,z> = AffineSpace(k, 3)
sage: C = Curve([x^2 + x - y^3, y^4 - y - z^3], A)
sage: C.genus()
10
sage: C.function_field()
Function field in z defined by z^9 + x^8 + x^6 + x^5 + x^4 + x^3 + x
>>> from sage.all import *
>>> k = GF(Integer(2), names=('a',)); (a,) = k._first_ngens(1)
>>> A = AffineSpace(k, Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = A._first_ngens(3)
>>> C = Curve([x**Integer(2) + x - y**Integer(3), y**Integer(4) - y - z**Integer(3)], A)
>>> C.genus()
10
>>> C.function_field()
Function field in z defined by z^9 + x^8 + x^6 + x^5 + x^4 + x^3 + x

Closed points of arbitrary degree can be computed:

sage: # long time
sage: C.closed_points()
[Point (x, y, z), Point (x + 1, y, z)]
sage: C.closed_points(2)
[Point (x^2 + x + 1, y + 1, z),
 Point (y^2 + y + 1, x + y, z),
 Point (y^2 + y + 1, x + y + 1, z)]
sage: p = _[0]
sage: p.places()
[Place (x^2 + x + 1, (1/(x^4 + x^2 + 1))*z^7 + (1/(x^4 + x^2 + 1))*z^6 + 1)]
>>> from sage.all import *
>>> # long time
>>> C.closed_points()
[Point (x, y, z), Point (x + 1, y, z)]
>>> C.closed_points(Integer(2))
[Point (x^2 + x + 1, y + 1, z),
 Point (y^2 + y + 1, x + y, z),
 Point (y^2 + y + 1, x + y + 1, z)]
>>> p = _[Integer(0)]
>>> p.places()
[Place (x^2 + x + 1, (1/(x^4 + x^2 + 1))*z^7 + (1/(x^4 + x^2 + 1))*z^6 + 1)]

The places at infinity correspond to the extra closed points of the curve’s projective closure:

sage: C.places_at_infinity()                # long time
[Place (1/x, 1/x*z)]
>>> from sage.all import *
>>> C.places_at_infinity()                # long time
[Place (1/x, 1/x*z)]

It is easy to transit to and from the function field of the curve:

sage: fx = C(x)
sage: fy = C(y)
sage: fx^2 + fx - fy^3
0
sage: fx.divisor()
-9*Place (1/x, 1/x*z)
 + 9*Place (x, z)
sage: p, = fx.zeros()
sage: C.place_to_closed_point(p)
Point (x, y, z)
sage: _.rational_point()
(0, 0, 0)
sage: _.closed_point()
Point (x, y, z)
sage: _.place()
Place (x, z)
>>> from sage.all import *
>>> fx = C(x)
>>> fy = C(y)
>>> fx**Integer(2) + fx - fy**Integer(3)
0
>>> fx.divisor()
-9*Place (1/x, 1/x*z)
 + 9*Place (x, z)
>>> p, = fx.zeros()
>>> C.place_to_closed_point(p)
Point (x, y, z)
>>> _.rational_point()
(0, 0, 0)
>>> _.closed_point()
Point (x, y, z)
>>> _.place()
Place (x, z)

Integral affine curves over \(\QQ\)

An integral curve over \(\QQ\) is equipped also with the function field. Unlike over finite fields, it is not possible to enumerate closed points.

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve(x^2 + y^2 -1)
sage: p = C(0,1)
sage: p
(0, 1)
sage: p.closed_point()
Point (x, y - 1)
sage: pl = _.place()
sage: C.parametric_representation(pl)
(s + ..., 1 - 1/2*s^2 - 1/8*s^4 - 1/16*s^6 + ...)
sage: sx, sy = _
sage: sx = sx.polynomial(10); sx
s
sage: sy = sy.polynomial(10); sy
-7/256*s^10 - 5/128*s^8 - 1/16*s^6 - 1/8*s^4 - 1/2*s^2 + 1
sage: s = var('s')                                                                  # needs sage.symbolic
sage: P1 = parametric_plot([sx, sy], (s, -1, 1), color='red')                       # needs sage.plot sage.symbolic
sage: P2 = C.plot((x, -1, 1), (y, 0, 2))  # half circle                             # needs sage.plot sage.symbolic
sage: P1 + P2                                                                       # needs sage.plot sage.symbolic
Graphics object consisting of 2 graphics primitives
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(2) + y**Integer(2) -Integer(1))
>>> p = C(Integer(0),Integer(1))
>>> p
(0, 1)
>>> p.closed_point()
Point (x, y - 1)
>>> pl = _.place()
>>> C.parametric_representation(pl)
(s + ..., 1 - 1/2*s^2 - 1/8*s^4 - 1/16*s^6 + ...)
>>> sx, sy = _
>>> sx = sx.polynomial(Integer(10)); sx
s
>>> sy = sy.polynomial(Integer(10)); sy
-7/256*s^10 - 5/128*s^8 - 1/16*s^6 - 1/8*s^4 - 1/2*s^2 + 1
>>> s = var('s')                                                                  # needs sage.symbolic
>>> P1 = parametric_plot([sx, sy], (s, -Integer(1), Integer(1)), color='red')                       # needs sage.plot sage.symbolic
>>> P2 = C.plot((x, -Integer(1), Integer(1)), (y, Integer(0), Integer(2)))  # half circle                             # needs sage.plot sage.symbolic
>>> P1 + P2                                                                       # needs sage.plot sage.symbolic
Graphics object consisting of 2 graphics primitives

AUTHORS:

  • William Stein (2005-11-13)

  • David Joyner (2005-11-13)

  • David Kohel (2006-01)

  • Grayson Jorgenson (2016-08)

  • Kwankyu Lee (2019-05): added integral affine curves

class sage.schemes.curves.affine_curve.AffineCurve(A, X)[source]

Bases: Curve_generic, AlgebraicScheme_subscheme_affine

Affine curves.

EXAMPLES:

sage: # needs sage.rings.number_field
sage: R.<v> = QQ[]
sage: K.<u> = NumberField(v^2 + 3)
sage: A.<x,y,z> = AffineSpace(K, 3)
sage: C = Curve([z - u*x^2, y^2], A); C
Affine Curve over Number Field in u with defining polynomial v^2 + 3
defined by (-u)*x^2 + z, y^2
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> R = QQ['v']; (v,) = R._first_ngens(1)
>>> K = NumberField(v**Integer(2) + Integer(3), names=('u',)); (u,) = K._first_ngens(1)
>>> A = AffineSpace(K, Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = A._first_ngens(3)
>>> C = Curve([z - u*x**Integer(2), y**Integer(2)], A); C
Affine Curve over Number Field in u with defining polynomial v^2 + 3
defined by (-u)*x^2 + z, y^2

sage: A.<x,y,z> = AffineSpace(GF(7), 3)
sage: C = Curve([x^2 - z, z - 8*x], A); C
Affine Curve over Finite Field of size 7 defined by x^2 - z, -x + z
>>> from sage.all import *
>>> A = AffineSpace(GF(Integer(7)), Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = A._first_ngens(3)
>>> C = Curve([x**Integer(2) - z, z - Integer(8)*x], A); C
Affine Curve over Finite Field of size 7 defined by x^2 - z, -x + z
projective_closure(i=0, PP=None)[source]

Return the projective closure of this affine curve.

INPUT:

  • i – (default: 0) the index of the affine coordinate chart of the projective space that the affine ambient space of this curve embeds into

  • PP – (default: None) ambient projective space to compute the projective closure in; this is constructed if it is not given

OUTPUT: a curve in projective space

EXAMPLES:

sage: A.<x,y,z> = AffineSpace(QQ, 3)
sage: C = Curve([y-x^2,z-x^3], A)
sage: C.projective_closure()
Projective Curve over Rational Field defined by x1^2 - x0*x2,
x1*x2 - x0*x3, x2^2 - x1*x3
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = A._first_ngens(3)
>>> C = Curve([y-x**Integer(2),z-x**Integer(3)], A)
>>> C.projective_closure()
Projective Curve over Rational Field defined by x1^2 - x0*x2,
x1*x2 - x0*x3, x2^2 - x1*x3

sage: A.<x,y,z> = AffineSpace(QQ, 3)
sage: C = Curve([y - x^2, z - x^3], A)
sage: C.projective_closure()
Projective Curve over Rational Field defined by
x1^2 - x0*x2, x1*x2 - x0*x3, x2^2 - x1*x3
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = A._first_ngens(3)
>>> C = Curve([y - x**Integer(2), z - x**Integer(3)], A)
>>> C.projective_closure()
Projective Curve over Rational Field defined by
x1^2 - x0*x2, x1*x2 - x0*x3, x2^2 - x1*x3

sage: A.<x,y> = AffineSpace(CC, 2)
sage: C = Curve(y - x^3 + x - 1, A)
sage: C.projective_closure(1)
Projective Plane Curve over Complex Field with 53 bits of precision defined by
x0^3 - x0*x1^2 + x1^3 - x1^2*x2
>>> from sage.all import *
>>> A = AffineSpace(CC, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(y - x**Integer(3) + x - Integer(1), A)
>>> C.projective_closure(Integer(1))
Projective Plane Curve over Complex Field with 53 bits of precision defined by
x0^3 - x0*x1^2 + x1^3 - x1^2*x2

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: P.<u,v,w> = ProjectiveSpace(QQ, 2)
sage: C = Curve([y - x^2], A)
sage: C.projective_closure(1, P).ambient_space() == P
True
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> P = ProjectiveSpace(QQ, Integer(2), names=('u', 'v', 'w',)); (u, v, w,) = P._first_ngens(3)
>>> C = Curve([y - x**Integer(2)], A)
>>> C.projective_closure(Integer(1), P).ambient_space() == P
True
class sage.schemes.curves.affine_curve.AffineCurve_field(A, X)[source]

Bases: AffineCurve, AlgebraicScheme_subscheme_affine_field

Affine curves over fields.

blowup(P=None)[source]

Return the blow up of this affine curve at the point P.

The blow up is described by affine charts. This curve must be irreducible.

INPUT:

  • P – (default: None) a point on this curve at which to blow up; if None, then P is taken to be the origin

OUTPUT: a tuple of

  • a tuple of curves in affine space of the same dimension as the ambient space of this curve, which define the blow up in each affine chart.

  • a tuple of tuples such that the j-th element of the i-th tuple is the transition map from the i-th affine patch to the j-th affine patch.

  • a tuple consisting of the restrictions of the projection map from the blow up back to the original curve, restricted to each affine patch. There the i-th element will be the projection from the i-th affine patch.

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve([y^2 - x^3], A)
sage: C.blowup()
((Affine Plane Curve over Rational Field defined by s1^2 - x,
  Affine Plane Curve over Rational Field defined by y*s0^3 - 1),
([Scheme endomorphism of Affine Plane Curve over Rational Field
   defined by s1^2 - x
    Defn: Defined on coordinates by sending (x, s1) to (x, s1),
  Scheme morphism:
    From: Affine Plane Curve over Rational Field defined by s1^2 - x
    To:   Affine Plane Curve over Rational Field defined by y*s0^3 - 1
    Defn: Defined on coordinates by sending (x, s1) to (x*s1, 1/s1)],
 [Scheme morphism:
    From: Affine Plane Curve over Rational Field defined by y*s0^3 - 1
    To:   Affine Plane Curve over Rational Field defined by s1^2 - x
    Defn: Defined on coordinates by sending (y, s0) to (y*s0, 1/s0),
  Scheme endomorphism of Affine Plane Curve over Rational Field
   defined by y*s0^3 - 1
    Defn: Defined on coordinates by sending (y, s0) to (y, s0)]),
(Scheme morphism:
   From: Affine Plane Curve over Rational Field defined by s1^2 - x
   To:   Affine Plane Curve over Rational Field defined by -x^3 + y^2
   Defn: Defined on coordinates by sending (x, s1) to (x, x*s1),
 Scheme morphism:
   From: Affine Plane Curve over Rational Field defined by y*s0^3 - 1
   To:   Affine Plane Curve over Rational Field defined by -x^3 + y^2
   Defn: Defined on coordinates by sending (y, s0) to (y*s0, y)))
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([y**Integer(2) - x**Integer(3)], A)
>>> C.blowup()
((Affine Plane Curve over Rational Field defined by s1^2 - x,
  Affine Plane Curve over Rational Field defined by y*s0^3 - 1),
([Scheme endomorphism of Affine Plane Curve over Rational Field
   defined by s1^2 - x
    Defn: Defined on coordinates by sending (x, s1) to (x, s1),
  Scheme morphism:
    From: Affine Plane Curve over Rational Field defined by s1^2 - x
    To:   Affine Plane Curve over Rational Field defined by y*s0^3 - 1
    Defn: Defined on coordinates by sending (x, s1) to (x*s1, 1/s1)],
 [Scheme morphism:
    From: Affine Plane Curve over Rational Field defined by y*s0^3 - 1
    To:   Affine Plane Curve over Rational Field defined by s1^2 - x
    Defn: Defined on coordinates by sending (y, s0) to (y*s0, 1/s0),
  Scheme endomorphism of Affine Plane Curve over Rational Field
   defined by y*s0^3 - 1
    Defn: Defined on coordinates by sending (y, s0) to (y, s0)]),
(Scheme morphism:
   From: Affine Plane Curve over Rational Field defined by s1^2 - x
   To:   Affine Plane Curve over Rational Field defined by -x^3 + y^2
   Defn: Defined on coordinates by sending (x, s1) to (x, x*s1),
 Scheme morphism:
   From: Affine Plane Curve over Rational Field defined by y*s0^3 - 1
   To:   Affine Plane Curve over Rational Field defined by -x^3 + y^2
   Defn: Defined on coordinates by sending (y, s0) to (y*s0, y)))

sage: # needs sage.rings.number_field
sage: K.<a> = QuadraticField(2)
sage: A.<x,y,z> = AffineSpace(K, 3)
sage: C = Curve([y^2 - a*x^5, x - z], A)
sage: B = C.blowup()
sage: B[0]
(Affine Curve over Number Field in a with defining polynomial x^2 - 2
  with a = 1.414213562373095? defined by s2 - 1, 2*x^3 + (-a)*s1^2,
 Affine Curve over Number Field in a with defining polynomial x^2 - 2
  with a = 1.414213562373095? defined by s0 - s2, 2*y^3*s2^5 + (-a),
 Affine Curve over Number Field in a with defining polynomial x^2 - 2
  with a = 1.414213562373095? defined by s0 - 1, 2*z^3 + (-a)*s1^2)
sage: B[1][0][2]
Scheme morphism:
  From: Affine Curve over Number Field in a
        with defining polynomial x^2 - 2 with a = 1.414213562373095?
        defined by s2 - 1, 2*x^3 + (-a)*s1^2
  To:   Affine Curve over Number Field in a
        with defining polynomial x^2 - 2 with a = 1.414213562373095?
        defined by s0 - 1, 2*z^3 + (-a)*s1^2
  Defn: Defined on coordinates by sending (x, s1, s2) to
        (x*s2, 1/s2, s1/s2)
sage: B[1][2][0]
Scheme morphism:
  From: Affine Curve over Number Field in a
        with defining polynomial x^2 - 2 with a = 1.414213562373095?
        defined by s0 - 1, 2*z^3 + (-a)*s1^2
  To:   Affine Curve over Number Field in a
        with defining polynomial x^2 - 2 with a = 1.414213562373095?
        defined by s2 - 1, 2*x^3 + (-a)*s1^2
  Defn: Defined on coordinates by sending (z, s0, s1) to
        (z*s0, s1/s0, 1/s0)
sage: B[2]
(Scheme morphism:
   From: Affine Curve over Number Field in a
         with defining polynomial x^2 - 2 with a = 1.414213562373095?
         defined by s2 - 1, 2*x^3 + (-a)*s1^2
   To:   Affine Curve over Number Field in a
         with defining polynomial x^2 - 2 with a = 1.414213562373095?
         defined by (-a)*x^5 + y^2, x - z
   Defn: Defined on coordinates by sending (x, s1, s2) to
         (x, x*s1, x*s2),
 Scheme morphism:
   From: Affine Curve over Number Field in a
         with defining polynomial x^2 - 2 with a = 1.414213562373095?
         defined by s0 - s2, 2*y^3*s2^5 + (-a)
   To:   Affine Curve over Number Field in a
         with defining polynomial x^2 - 2 with a = 1.414213562373095?
         defined by (-a)*x^5 + y^2, x - z
   Defn: Defined on coordinates by sending (y, s0, s2) to
         (y*s0, y, y*s2),
 Scheme morphism:
   From: Affine Curve over Number Field in a
         with defining polynomial x^2 - 2 with a = 1.414213562373095?
         defined by s0 - 1, 2*z^3 + (-a)*s1^2
   To:   Affine Curve over Number Field in a
         with defining polynomial x^2 - 2 with a = 1.414213562373095?
         defined by (-a)*x^5 + y^2, x - z
   Defn: Defined on coordinates by sending (z, s0, s1) to
         (z*s0, z*s1, z))
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> K = QuadraticField(Integer(2), names=('a',)); (a,) = K._first_ngens(1)
>>> A = AffineSpace(K, Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = A._first_ngens(3)
>>> C = Curve([y**Integer(2) - a*x**Integer(5), x - z], A)
>>> B = C.blowup()
>>> B[Integer(0)]
(Affine Curve over Number Field in a with defining polynomial x^2 - 2
  with a = 1.414213562373095? defined by s2 - 1, 2*x^3 + (-a)*s1^2,
 Affine Curve over Number Field in a with defining polynomial x^2 - 2
  with a = 1.414213562373095? defined by s0 - s2, 2*y^3*s2^5 + (-a),
 Affine Curve over Number Field in a with defining polynomial x^2 - 2
  with a = 1.414213562373095? defined by s0 - 1, 2*z^3 + (-a)*s1^2)
>>> B[Integer(1)][Integer(0)][Integer(2)]
Scheme morphism:
  From: Affine Curve over Number Field in a
        with defining polynomial x^2 - 2 with a = 1.414213562373095?
        defined by s2 - 1, 2*x^3 + (-a)*s1^2
  To:   Affine Curve over Number Field in a
        with defining polynomial x^2 - 2 with a = 1.414213562373095?
        defined by s0 - 1, 2*z^3 + (-a)*s1^2
  Defn: Defined on coordinates by sending (x, s1, s2) to
        (x*s2, 1/s2, s1/s2)
>>> B[Integer(1)][Integer(2)][Integer(0)]
Scheme morphism:
  From: Affine Curve over Number Field in a
        with defining polynomial x^2 - 2 with a = 1.414213562373095?
        defined by s0 - 1, 2*z^3 + (-a)*s1^2
  To:   Affine Curve over Number Field in a
        with defining polynomial x^2 - 2 with a = 1.414213562373095?
        defined by s2 - 1, 2*x^3 + (-a)*s1^2
  Defn: Defined on coordinates by sending (z, s0, s1) to
        (z*s0, s1/s0, 1/s0)
>>> B[Integer(2)]
(Scheme morphism:
   From: Affine Curve over Number Field in a
         with defining polynomial x^2 - 2 with a = 1.414213562373095?
         defined by s2 - 1, 2*x^3 + (-a)*s1^2
   To:   Affine Curve over Number Field in a
         with defining polynomial x^2 - 2 with a = 1.414213562373095?
         defined by (-a)*x^5 + y^2, x - z
   Defn: Defined on coordinates by sending (x, s1, s2) to
         (x, x*s1, x*s2),
 Scheme morphism:
   From: Affine Curve over Number Field in a
         with defining polynomial x^2 - 2 with a = 1.414213562373095?
         defined by s0 - s2, 2*y^3*s2^5 + (-a)
   To:   Affine Curve over Number Field in a
         with defining polynomial x^2 - 2 with a = 1.414213562373095?
         defined by (-a)*x^5 + y^2, x - z
   Defn: Defined on coordinates by sending (y, s0, s2) to
         (y*s0, y, y*s2),
 Scheme morphism:
   From: Affine Curve over Number Field in a
         with defining polynomial x^2 - 2 with a = 1.414213562373095?
         defined by s0 - 1, 2*z^3 + (-a)*s1^2
   To:   Affine Curve over Number Field in a
         with defining polynomial x^2 - 2 with a = 1.414213562373095?
         defined by (-a)*x^5 + y^2, x - z
   Defn: Defined on coordinates by sending (z, s0, s1) to
         (z*s0, z*s1, z))

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = A.curve((y - 3/2)^3 - (x + 2)^5 - (x + 2)^6)
sage: Q = A([-2,3/2])
sage: C.blowup(Q)
((Affine Plane Curve over Rational Field
   defined by x^3 - s1^3 + 7*x^2 + 16*x + 12,
  Affine Plane Curve over Rational Field
   defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5
              + 54*y*s0^6 - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8),
 ([Scheme endomorphism of Affine Plane Curve over Rational Field
    defined by x^3 - s1^3 + 7*x^2 + 16*x + 12
     Defn: Defined on coordinates by sending (x, s1) to (x, s1),
   Scheme morphism:
     From: Affine Plane Curve over Rational Field
           defined by x^3 - s1^3 + 7*x^2 + 16*x + 12
     To:   Affine Plane Curve over Rational Field
           defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5
                      + 54*y*s0^6 - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8
     Defn: Defined on coordinates by sending (x, s1) to
           (x*s1 + 2*s1 + 3/2, 1/s1)],
  [Scheme morphism:
     From: Affine Plane Curve over Rational Field
           defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5
                      + 54*y*s0^6 - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8
     To:   Affine Plane Curve over Rational Field
           defined by x^3 - s1^3 + 7*x^2 + 16*x + 12
     Defn: Defined on coordinates by sending (y, s0) to
           (y*s0 - 3/2*s0 - 2, 1/s0),
   Scheme endomorphism of Affine Plane Curve over Rational Field
    defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5 + 54*y*s0^6
               - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8
     Defn: Defined on coordinates by sending (y, s0) to (y, s0)]),
 (Scheme morphism:
    From: Affine Plane Curve over Rational Field
          defined by x^3 - s1^3 + 7*x^2 + 16*x + 12
    To:   Affine Plane Curve over Rational Field
          defined by -x^6 - 13*x^5 - 70*x^4 - 200*x^3 + y^3
                     - 320*x^2 - 9/2*y^2 - 272*x + 27/4*y - 795/8
    Defn: Defined on coordinates by sending (x, s1) to
          (x, x*s1 + 2*s1 + 3/2),
  Scheme morphism:
    From: Affine Plane Curve over Rational Field
          defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5
                     + 54*y*s0^6 - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8
    To:   Affine Plane Curve over Rational Field
          defined by -x^6 - 13*x^5 - 70*x^4 - 200*x^3 + y^3
                     - 320*x^2 - 9/2*y^2 - 272*x + 27/4*y - 795/8
    Defn: Defined on coordinates by sending (y, s0) to
          (y*s0 - 3/2*s0 - 2, y)))
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = A.curve((y - Integer(3)/Integer(2))**Integer(3) - (x + Integer(2))**Integer(5) - (x + Integer(2))**Integer(6))
>>> Q = A([-Integer(2),Integer(3)/Integer(2)])
>>> C.blowup(Q)
((Affine Plane Curve over Rational Field
   defined by x^3 - s1^3 + 7*x^2 + 16*x + 12,
  Affine Plane Curve over Rational Field
   defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5
              + 54*y*s0^6 - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8),
 ([Scheme endomorphism of Affine Plane Curve over Rational Field
    defined by x^3 - s1^3 + 7*x^2 + 16*x + 12
     Defn: Defined on coordinates by sending (x, s1) to (x, s1),
   Scheme morphism:
     From: Affine Plane Curve over Rational Field
           defined by x^3 - s1^3 + 7*x^2 + 16*x + 12
     To:   Affine Plane Curve over Rational Field
           defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5
                      + 54*y*s0^6 - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8
     Defn: Defined on coordinates by sending (x, s1) to
           (x*s1 + 2*s1 + 3/2, 1/s1)],
  [Scheme morphism:
     From: Affine Plane Curve over Rational Field
           defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5
                      + 54*y*s0^6 - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8
     To:   Affine Plane Curve over Rational Field
           defined by x^3 - s1^3 + 7*x^2 + 16*x + 12
     Defn: Defined on coordinates by sending (y, s0) to
           (y*s0 - 3/2*s0 - 2, 1/s0),
   Scheme endomorphism of Affine Plane Curve over Rational Field
    defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5 + 54*y*s0^6
               - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8
     Defn: Defined on coordinates by sending (y, s0) to (y, s0)]),
 (Scheme morphism:
    From: Affine Plane Curve over Rational Field
          defined by x^3 - s1^3 + 7*x^2 + 16*x + 12
    To:   Affine Plane Curve over Rational Field
          defined by -x^6 - 13*x^5 - 70*x^4 - 200*x^3 + y^3
                     - 320*x^2 - 9/2*y^2 - 272*x + 27/4*y - 795/8
    Defn: Defined on coordinates by sending (x, s1) to
          (x, x*s1 + 2*s1 + 3/2),
  Scheme morphism:
    From: Affine Plane Curve over Rational Field
          defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5
                     + 54*y*s0^6 - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8
    To:   Affine Plane Curve over Rational Field
          defined by -x^6 - 13*x^5 - 70*x^4 - 200*x^3 + y^3
                     - 320*x^2 - 9/2*y^2 - 272*x + 27/4*y - 795/8
    Defn: Defined on coordinates by sending (y, s0) to
          (y*s0 - 3/2*s0 - 2, y)))

sage: A.<x,y,z,w> = AffineSpace(QQ, 4)
sage: C = A.curve([((x + 1)^2 + y^2)^3 - 4*(x + 1)^2*y^2, y - z, w - 4])
sage: Q = C([-1,0,0,4])
sage: B = C.blowup(Q)
sage: B[0]
(Affine Curve over Rational Field defined by s3, s1 - s2,
  x^2*s2^6 + 2*x*s2^6 + 3*x^2*s2^4 + s2^6 + 6*x*s2^4
  + 3*x^2*s2^2 + 3*s2^4 + 6*x*s2^2 + x^2 - s2^2 + 2*x + 1,
 Affine Curve over Rational Field defined by s3, s2 - 1,
  y^2*s0^6 + 3*y^2*s0^4 + 3*y^2*s0^2 + y^2 - 4*s0^2,
 Affine Curve over Rational Field defined by s3, s1 - 1,
  z^2*s0^6 + 3*z^2*s0^4 + 3*z^2*s0^2 + z^2 - 4*s0^2,
 Closed subscheme of Affine Space of dimension 4 over Rational Field
  defined by: 1)
sage: Q = A([6,2,3,1])
sage: B = C.blowup(Q)
Traceback (most recent call last):
...
TypeError: (=(6, 2, 3, 1)) must be a point on this curve
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(4), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = A._first_ngens(4)
>>> C = A.curve([((x + Integer(1))**Integer(2) + y**Integer(2))**Integer(3) - Integer(4)*(x + Integer(1))**Integer(2)*y**Integer(2), y - z, w - Integer(4)])
>>> Q = C([-Integer(1),Integer(0),Integer(0),Integer(4)])
>>> B = C.blowup(Q)
>>> B[Integer(0)]
(Affine Curve over Rational Field defined by s3, s1 - s2,
  x^2*s2^6 + 2*x*s2^6 + 3*x^2*s2^4 + s2^6 + 6*x*s2^4
  + 3*x^2*s2^2 + 3*s2^4 + 6*x*s2^2 + x^2 - s2^2 + 2*x + 1,
 Affine Curve over Rational Field defined by s3, s2 - 1,
  y^2*s0^6 + 3*y^2*s0^4 + 3*y^2*s0^2 + y^2 - 4*s0^2,
 Affine Curve over Rational Field defined by s3, s1 - 1,
  z^2*s0^6 + 3*z^2*s0^4 + 3*z^2*s0^2 + z^2 - 4*s0^2,
 Closed subscheme of Affine Space of dimension 4 over Rational Field
  defined by: 1)
>>> Q = A([Integer(6),Integer(2),Integer(3),Integer(1)])
>>> B = C.blowup(Q)
Traceback (most recent call last):
...
TypeError: (=(6, 2, 3, 1)) must be a point on this curve

sage: # needs sage.rings.number_field
sage: A.<x,y> = AffineSpace(QuadraticField(-1), 2)
sage: C = A.curve([y^2 + x^2])
sage: C.blowup()
Traceback (most recent call last):
...
TypeError: this curve must be irreducible
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> A = AffineSpace(QuadraticField(-Integer(1)), Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = A.curve([y**Integer(2) + x**Integer(2)])
>>> C.blowup()
Traceback (most recent call last):
...
TypeError: this curve must be irreducible
plane_projection(AP=None)[source]

Return a projection of this curve into an affine plane so that the image of the projection is a plane curve.

INPUT:

  • AP – (default: None) the affine plane to project this curve into. This space must be defined over the same base field as this curve, and must have dimension two. This space will be constructed if not specified.

OUTPUT: a tuple of

  • a scheme morphism from this curve into an affine plane

  • the plane curve that defines the image of that morphism

EXAMPLES:

sage: A.<x,y,z,w> = AffineSpace(QQ, 4)
sage: C = Curve([x^2 - y*z*w, z^3 - w, w + x*y - 3*z^3], A)
sage: C.plane_projection()
(Scheme morphism:
  From: Affine Curve over Rational Field defined by
        -y*z*w + x^2, z^3 - w, -3*z^3 + x*y + w
  To:   Affine Space of dimension 2 over Rational Field
  Defn: Defined on coordinates by sending (x, y, z, w) to (x, y),
 Affine Plane Curve over Rational Field defined by
 x0^2*x1^7 - 16*x0^4)
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(4), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = A._first_ngens(4)
>>> C = Curve([x**Integer(2) - y*z*w, z**Integer(3) - w, w + x*y - Integer(3)*z**Integer(3)], A)
>>> C.plane_projection()
(Scheme morphism:
  From: Affine Curve over Rational Field defined by
        -y*z*w + x^2, z^3 - w, -3*z^3 + x*y + w
  To:   Affine Space of dimension 2 over Rational Field
  Defn: Defined on coordinates by sending (x, y, z, w) to (x, y),
 Affine Plane Curve over Rational Field defined by
 x0^2*x1^7 - 16*x0^4)

sage: # needs sage.rings.number_field
sage: R.<a> = QQ[]
sage: K.<b> = NumberField(a^2 + 2)
sage: A.<x,y,z> = AffineSpace(K, 3)
sage: C = A.curve([x - b, y - 2])
sage: B.<a,b> = AffineSpace(K, 2)
sage: proj1 = C.plane_projection(AP=B)
sage: proj1
(Scheme morphism:
   From: Affine Curve over Number Field in b
         with defining polynomial a^2 + 2 defined by x + (-b), y - 2
   To:   Affine Space of dimension 2 over Number Field in b
         with defining polynomial a^2 + 2
   Defn: Defined on coordinates by sending (x, y, z) to
         (x, z),
 Affine Plane Curve over Number Field in b
 with defining polynomial a^2 + 2 defined by a + (-b))
sage: proj1[1].ambient_space() is B
True
sage: proj2 = C.plane_projection()
sage: proj2[1].ambient_space() is B
False
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> R = QQ['a']; (a,) = R._first_ngens(1)
>>> K = NumberField(a**Integer(2) + Integer(2), names=('b',)); (b,) = K._first_ngens(1)
>>> A = AffineSpace(K, Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = A._first_ngens(3)
>>> C = A.curve([x - b, y - Integer(2)])
>>> B = AffineSpace(K, Integer(2), names=('a', 'b',)); (a, b,) = B._first_ngens(2)
>>> proj1 = C.plane_projection(AP=B)
>>> proj1
(Scheme morphism:
   From: Affine Curve over Number Field in b
         with defining polynomial a^2 + 2 defined by x + (-b), y - 2
   To:   Affine Space of dimension 2 over Number Field in b
         with defining polynomial a^2 + 2
   Defn: Defined on coordinates by sending (x, y, z) to
         (x, z),
 Affine Plane Curve over Number Field in b
 with defining polynomial a^2 + 2 defined by a + (-b))
>>> proj1[Integer(1)].ambient_space() is B
True
>>> proj2 = C.plane_projection()
>>> proj2[Integer(1)].ambient_space() is B
False
projection(indices, AS=None)[source]

Return the projection of this curve onto the coordinates specified by indices.

INPUT:

  • indices – list or tuple of distinct integers specifying the indices of the coordinates to use in the projection. Can also be a list or tuple consisting of variables of the coordinate ring of the ambient space of this curve. If integers are used to specify the coordinates, 0 denotes the first coordinate. The length of indices must be between two and one less than the dimension of the ambient space of this curve, inclusive.

  • AS – (default: None) the affine space the projected curve will be defined in. This space must be defined over the same base field as this curve, and must have dimension equal to the length of indices. This space is constructed if not specified.

OUTPUT: a tuple of

  • a scheme morphism from this curve to affine space of dimension equal to the number of coordinates specified in indices

  • the affine subscheme that is the image of that morphism. If the image is a curve, the second element of the tuple will be a curve.

EXAMPLES:

sage: A.<x,y,z> = AffineSpace(QQ, 3)
sage: C = Curve([y^7 - x^2 + x^3 - 2*z, z^2 - x^7 - y^2], A)
sage: C.projection([0,1])
(Scheme morphism:
   From: Affine Curve over Rational Field
         defined by y^7 + x^3 - x^2 - 2*z, -x^7 - y^2 + z^2
   To:   Affine Space of dimension 2 over Rational Field
   Defn: Defined on coordinates by sending (x, y, z) to
         (x, y),
 Affine Plane Curve over Rational Field defined by x1^14 + 2*x0^3*x1^7 -
2*x0^2*x1^7 - 4*x0^7 + x0^6 - 2*x0^5 + x0^4 - 4*x1^2)
sage: C.projection([0,1,3,4])
Traceback (most recent call last):
...
ValueError: (=[0, 1, 3, 4]) must be a list or tuple of length between 2
and (=2), inclusive
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = A._first_ngens(3)
>>> C = Curve([y**Integer(7) - x**Integer(2) + x**Integer(3) - Integer(2)*z, z**Integer(2) - x**Integer(7) - y**Integer(2)], A)
>>> C.projection([Integer(0),Integer(1)])
(Scheme morphism:
   From: Affine Curve over Rational Field
         defined by y^7 + x^3 - x^2 - 2*z, -x^7 - y^2 + z^2
   To:   Affine Space of dimension 2 over Rational Field
   Defn: Defined on coordinates by sending (x, y, z) to
         (x, y),
 Affine Plane Curve over Rational Field defined by x1^14 + 2*x0^3*x1^7 -
2*x0^2*x1^7 - 4*x0^7 + x0^6 - 2*x0^5 + x0^4 - 4*x1^2)
>>> C.projection([Integer(0),Integer(1),Integer(3),Integer(4)])
Traceback (most recent call last):
...
ValueError: (=[0, 1, 3, 4]) must be a list or tuple of length between 2
and (=2), inclusive

sage: A.<x,y,z,w> = AffineSpace(QQ, 4)
sage: C = Curve([x - 2, y - 3, z - 1], A)
sage: B.<a,b,c> = AffineSpace(QQ, 3)
sage: C.projection([0,1,2], AS=B)
(Scheme morphism:
   From: Affine Curve over Rational Field defined by x - 2, y - 3, z - 1
   To:   Affine Space of dimension 3 over Rational Field
   Defn: Defined on coordinates by sending (x, y, z, w) to
         (x, y, z),
 Closed subscheme of Affine Space of dimension 3 over Rational Field defined by:
   c - 1,
   b - 3,
   a - 2)
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(4), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = A._first_ngens(4)
>>> C = Curve([x - Integer(2), y - Integer(3), z - Integer(1)], A)
>>> B = AffineSpace(QQ, Integer(3), names=('a', 'b', 'c',)); (a, b, c,) = B._first_ngens(3)
>>> C.projection([Integer(0),Integer(1),Integer(2)], AS=B)
(Scheme morphism:
   From: Affine Curve over Rational Field defined by x - 2, y - 3, z - 1
   To:   Affine Space of dimension 3 over Rational Field
   Defn: Defined on coordinates by sending (x, y, z, w) to
         (x, y, z),
 Closed subscheme of Affine Space of dimension 3 over Rational Field defined by:
   c - 1,
   b - 3,
   a - 2)

sage: A.<x,y,z,w,u> = AffineSpace(GF(11), 5)
sage: C = Curve([x^3 - 5*y*z + u^2, x - y^2 + 3*z^2,
....:            w^2 + 2*u^3*y, y - u^2 + z*x], A)
sage: B.<a,b,c> = AffineSpace(GF(11), 3)
sage: proj1 = C.projection([1,2,4], AS=B); proj1
(Scheme morphism:
   From: Affine Curve over Finite Field of size 11 defined by x^3 -
         5*y*z + u^2, -y^2 + 3*z^2 + x, 2*y*u^3 + w^2, x*z - u^2 + y
   To:   Affine Space of dimension 3 over Finite Field of size 11
   Defn: Defined on coordinates by sending (x, y, z, w, u) to
         (y, z, u),
 Affine Curve over Finite Field of size 11 defined by a^2*b - 3*b^3 -
c^2 + a, c^6 - 5*a*b^4 + b^3*c^2 - 3*a*c^4 + 3*a^2*c^2 - a^3, a^2*c^4 -
3*b^2*c^4 - 2*a^3*c^2 - 5*a*b^2*c^2 + a^4 - 5*a*b^3 + 2*b^4 + b^2*c^2 -
3*b*c^2 + 3*a*b, a^4*c^2 + 2*b^4*c^2 - a^5 - 2*a*b^4 + 5*b*c^4 + a*b*c^2
- 5*a*b^2 + 4*b^3 + b*c^2 + 5*c^2 - 5*a, a^6 - 5*b^6 - 5*b^3*c^2 +
5*a*b^3 + 2*c^4 - 4*a*c^2 + 2*a^2 - 5*a*b + c^2)
sage: proj1[1].ambient_space() is B
True
sage: proj2 = C.projection([1,2,4])
sage: proj2[1].ambient_space() is B
False
sage: C.projection([1,2,3,5], AS=B)
Traceback (most recent call last):
...
TypeError: (=Affine Space of dimension 3 over Finite Field of size 11)
must have dimension (=4)
>>> from sage.all import *
>>> A = AffineSpace(GF(Integer(11)), Integer(5), names=('x', 'y', 'z', 'w', 'u',)); (x, y, z, w, u,) = A._first_ngens(5)
>>> C = Curve([x**Integer(3) - Integer(5)*y*z + u**Integer(2), x - y**Integer(2) + Integer(3)*z**Integer(2),
...            w**Integer(2) + Integer(2)*u**Integer(3)*y, y - u**Integer(2) + z*x], A)
>>> B = AffineSpace(GF(Integer(11)), Integer(3), names=('a', 'b', 'c',)); (a, b, c,) = B._first_ngens(3)
>>> proj1 = C.projection([Integer(1),Integer(2),Integer(4)], AS=B); proj1
(Scheme morphism:
   From: Affine Curve over Finite Field of size 11 defined by x^3 -
         5*y*z + u^2, -y^2 + 3*z^2 + x, 2*y*u^3 + w^2, x*z - u^2 + y
   To:   Affine Space of dimension 3 over Finite Field of size 11
   Defn: Defined on coordinates by sending (x, y, z, w, u) to
         (y, z, u),
 Affine Curve over Finite Field of size 11 defined by a^2*b - 3*b^3 -
c^2 + a, c^6 - 5*a*b^4 + b^3*c^2 - 3*a*c^4 + 3*a^2*c^2 - a^3, a^2*c^4 -
3*b^2*c^4 - 2*a^3*c^2 - 5*a*b^2*c^2 + a^4 - 5*a*b^3 + 2*b^4 + b^2*c^2 -
3*b*c^2 + 3*a*b, a^4*c^2 + 2*b^4*c^2 - a^5 - 2*a*b^4 + 5*b*c^4 + a*b*c^2
- 5*a*b^2 + 4*b^3 + b*c^2 + 5*c^2 - 5*a, a^6 - 5*b^6 - 5*b^3*c^2 +
5*a*b^3 + 2*c^4 - 4*a*c^2 + 2*a^2 - 5*a*b + c^2)
>>> proj1[Integer(1)].ambient_space() is B
True
>>> proj2 = C.projection([Integer(1),Integer(2),Integer(4)])
>>> proj2[Integer(1)].ambient_space() is B
False
>>> C.projection([Integer(1),Integer(2),Integer(3),Integer(5)], AS=B)
Traceback (most recent call last):
...
TypeError: (=Affine Space of dimension 3 over Finite Field of size 11)
must have dimension (=4)

sage: A.<x,y,z,w> = AffineSpace(QQ, 4)
sage: C = A.curve([x*y - z^3, x*z - w^3, w^2 - x^3])
sage: C.projection([y,z])
(Scheme morphism:
   From: Affine Curve over Rational Field defined by
         -z^3 + x*y, -w^3 + x*z, -x^3 + w^2
   To:   Affine Space of dimension 2 over Rational Field
   Defn: Defined on coordinates by sending (x, y, z, w) to (y, z),
 Affine Plane Curve over Rational Field defined by x1^23 - x0^7*x1^4)
sage: B.<x,y,z> = AffineSpace(QQ, 3)
sage: C.projection([x,y,z], AS=B)
(Scheme morphism:
   From: Affine Curve over Rational Field defined by
         -z^3 + x*y, -w^3 + x*z, -x^3 + w^2
   To:   Affine Space of dimension 3 over Rational Field
   Defn: Defined on coordinates by sending (x, y, z, w) to
         (x, y, z),
 Affine Curve over Rational Field defined by
 z^3 - x*y, x^8 - x*z^2, x^7*z^2 - x*y*z)
sage: C.projection([y,z,z])
Traceback (most recent call last):
...
ValueError: (=[y, z, z]) must be a list or tuple of distinct indices or
variables
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(4), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = A._first_ngens(4)
>>> C = A.curve([x*y - z**Integer(3), x*z - w**Integer(3), w**Integer(2) - x**Integer(3)])
>>> C.projection([y,z])
(Scheme morphism:
   From: Affine Curve over Rational Field defined by
         -z^3 + x*y, -w^3 + x*z, -x^3 + w^2
   To:   Affine Space of dimension 2 over Rational Field
   Defn: Defined on coordinates by sending (x, y, z, w) to (y, z),
 Affine Plane Curve over Rational Field defined by x1^23 - x0^7*x1^4)
>>> B = AffineSpace(QQ, Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = B._first_ngens(3)
>>> C.projection([x,y,z], AS=B)
(Scheme morphism:
   From: Affine Curve over Rational Field defined by
         -z^3 + x*y, -w^3 + x*z, -x^3 + w^2
   To:   Affine Space of dimension 3 over Rational Field
   Defn: Defined on coordinates by sending (x, y, z, w) to
         (x, y, z),
 Affine Curve over Rational Field defined by
 z^3 - x*y, x^8 - x*z^2, x^7*z^2 - x*y*z)
>>> C.projection([y,z,z])
Traceback (most recent call last):
...
ValueError: (=[y, z, z]) must be a list or tuple of distinct indices or
variables
resolution_of_singularities(extend=False)[source]

Return a nonsingular model for this affine curve created by blowing up its singular points.

The nonsingular model is given as a collection of affine patches that cover it. If extend is False and if the base field is a number field, or if the base field is a finite field, the model returned may have singularities with coordinates not contained in the base field. An error is returned if this curve is already nonsingular, or if it has no singular points over its base field. This curve must be irreducible, and must be defined over a number field or finite field.

INPUT:

  • extend – boolean (default: False); specifies whether to extend the base field when necessary to find all singular points when this curve is defined over a number field. If extend is False, then only singularities with coordinates in the base field of this curve will be resolved. However, setting extend to True will slow down computations.

OUTPUT: a tuple of

  • a tuple of curves in affine space of the same dimension as the ambient space of this curve, which represent affine patches of the resolution of singularities.

  • a tuple of tuples such that the j-th element of the i-th tuple is the transition map from the i-th patch to the j-th patch.

  • a tuple consisting of birational maps from the patches back to the original curve that were created by composing the projection maps generated from the blow up computations. There the i-th element will be a map from the i-th patch.

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve([y^2 - x^3], A)
sage: C.resolution_of_singularities()
((Affine Plane Curve over Rational Field defined by s1^2 - x,
  Affine Plane Curve over Rational Field defined by y*s0^3 - 1),
 ((Scheme endomorphism of Affine Plane Curve over Rational Field
    defined by s1^2 - x
     Defn: Defined on coordinates by sending (x, s1) to (x, s1),
   Scheme morphism:
     From: Affine Plane Curve over Rational Field defined by s1^2 - x
     To:   Affine Plane Curve over Rational Field defined by y*s0^3 - 1
     Defn: Defined on coordinates by sending (x, s1) to (x*s1, 1/s1)),
  (Scheme morphism:
     From: Affine Plane Curve over Rational Field defined by y*s0^3 - 1
     To:   Affine Plane Curve over Rational Field defined by s1^2 - x
     Defn: Defined on coordinates by sending (y, s0) to (y*s0, 1/s0),
   Scheme endomorphism of Affine Plane Curve over Rational Field
    defined by y*s0^3 - 1
     Defn: Defined on coordinates by sending (y, s0) to (y, s0))),
 (Scheme morphism:
    From: Affine Plane Curve over Rational Field defined by s1^2 - x
    To:   Affine Plane Curve over Rational Field defined by -x^3 + y^2
    Defn: Defined on coordinates by sending (x, s1) to (x, x*s1),
  Scheme morphism:
    From: Affine Plane Curve over Rational Field defined by y*s0^3 - 1
    To:   Affine Plane Curve over Rational Field defined by -x^3 + y^2
    Defn: Defined on coordinates by sending (y, s0) to (y*s0, y)))
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([y**Integer(2) - x**Integer(3)], A)
>>> C.resolution_of_singularities()
((Affine Plane Curve over Rational Field defined by s1^2 - x,
  Affine Plane Curve over Rational Field defined by y*s0^3 - 1),
 ((Scheme endomorphism of Affine Plane Curve over Rational Field
    defined by s1^2 - x
     Defn: Defined on coordinates by sending (x, s1) to (x, s1),
   Scheme morphism:
     From: Affine Plane Curve over Rational Field defined by s1^2 - x
     To:   Affine Plane Curve over Rational Field defined by y*s0^3 - 1
     Defn: Defined on coordinates by sending (x, s1) to (x*s1, 1/s1)),
  (Scheme morphism:
     From: Affine Plane Curve over Rational Field defined by y*s0^3 - 1
     To:   Affine Plane Curve over Rational Field defined by s1^2 - x
     Defn: Defined on coordinates by sending (y, s0) to (y*s0, 1/s0),
   Scheme endomorphism of Affine Plane Curve over Rational Field
    defined by y*s0^3 - 1
     Defn: Defined on coordinates by sending (y, s0) to (y, s0))),
 (Scheme morphism:
    From: Affine Plane Curve over Rational Field defined by s1^2 - x
    To:   Affine Plane Curve over Rational Field defined by -x^3 + y^2
    Defn: Defined on coordinates by sending (x, s1) to (x, x*s1),
  Scheme morphism:
    From: Affine Plane Curve over Rational Field defined by y*s0^3 - 1
    To:   Affine Plane Curve over Rational Field defined by -x^3 + y^2
    Defn: Defined on coordinates by sending (y, s0) to (y*s0, y)))

sage: # needs sage.rings.number_field
sage: set_verbose(-1)
sage: K.<a> = QuadraticField(3)
sage: A.<x,y> = AffineSpace(K, 2)
sage: C = A.curve(x^4 + 2*x^2 + a*y^3 + 1)
sage: C.resolution_of_singularities(extend=True)[0]         # long time (2 s)
(Affine Plane Curve over Number Field in a0
  with defining polynomial y^4 - 4*y^2 + 16
  defined by 24*x^2*ss1^3 + 24*ss1^3 + (a0^3 - 8*a0),
 Affine Plane Curve over Number Field in a0
  with defining polynomial y^4 - 4*y^2 + 16
  defined by 24*s1^2*ss0 + (a0^3 - 8*a0)*ss0^2 + (-6*a0^3)*s1,
 Affine Plane Curve over Number Field in a0
  with defining polynomial y^4 - 4*y^2 + 16
  defined by 8*y^2*s0^4 + (4*a0^3)*y*s0^3 - 32*s0^2 + (a0^3 - 8*a0)*y)
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> set_verbose(-Integer(1))
>>> K = QuadraticField(Integer(3), names=('a',)); (a,) = K._first_ngens(1)
>>> A = AffineSpace(K, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = A.curve(x**Integer(4) + Integer(2)*x**Integer(2) + a*y**Integer(3) + Integer(1))
>>> C.resolution_of_singularities(extend=True)[Integer(0)]         # long time (2 s)
(Affine Plane Curve over Number Field in a0
  with defining polynomial y^4 - 4*y^2 + 16
  defined by 24*x^2*ss1^3 + 24*ss1^3 + (a0^3 - 8*a0),
 Affine Plane Curve over Number Field in a0
  with defining polynomial y^4 - 4*y^2 + 16
  defined by 24*s1^2*ss0 + (a0^3 - 8*a0)*ss0^2 + (-6*a0^3)*s1,
 Affine Plane Curve over Number Field in a0
  with defining polynomial y^4 - 4*y^2 + 16
  defined by 8*y^2*s0^4 + (4*a0^3)*y*s0^3 - 32*s0^2 + (a0^3 - 8*a0)*y)

sage: A.<x,y,z> = AffineSpace(GF(5), 3)
sage: C = Curve([y - x^3, (z - 2)^2 - y^3 - x^3], A)
sage: R = C.resolution_of_singularities()
sage: R[0]
(Affine Curve over Finite Field of size 5
  defined by x^2 - s1, s1^4 - x*s2^2 + s1, x*s1^3 - s2^2 + x,
 Affine Curve over Finite Field of size 5
  defined by y*s2^2 - y^2 - 1, s2^4 - s0^3 - y^2 - 2, y*s0^3 - s2^2 + y,
 Affine Curve over Finite Field of size 5
  defined by s0^3*s1 + z*s1^3 + s1^4 - 2*s1^3 - 1,
             z*s0^3 + z*s1^3 - 2*s0^3 - 2*s1^3 - 1,
             z^2*s1^3 + z*s1^3 - s1^3 - z + s1 + 2)
>>> from sage.all import *
>>> A = AffineSpace(GF(Integer(5)), Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = A._first_ngens(3)
>>> C = Curve([y - x**Integer(3), (z - Integer(2))**Integer(2) - y**Integer(3) - x**Integer(3)], A)
>>> R = C.resolution_of_singularities()
>>> R[Integer(0)]
(Affine Curve over Finite Field of size 5
  defined by x^2 - s1, s1^4 - x*s2^2 + s1, x*s1^3 - s2^2 + x,
 Affine Curve over Finite Field of size 5
  defined by y*s2^2 - y^2 - 1, s2^4 - s0^3 - y^2 - 2, y*s0^3 - s2^2 + y,
 Affine Curve over Finite Field of size 5
  defined by s0^3*s1 + z*s1^3 + s1^4 - 2*s1^3 - 1,
             z*s0^3 + z*s1^3 - 2*s0^3 - 2*s1^3 - 1,
             z^2*s1^3 + z*s1^3 - s1^3 - z + s1 + 2)

sage: A.<x,y,z,w> = AffineSpace(QQ, 4)
sage: C = A.curve([((x - 2)^2 + y^2)^2 - (x - 2)^2 - y^2 + (x - 2)^3,
....:              z - y - 7, w - 4])
sage: B = C.resolution_of_singularities()
sage: B[0]
(Affine Curve over Rational Field defined by s3, s1 - s2,
  x^2*s2^4 - 4*x*s2^4 + 2*x^2*s2^2 + 4*s2^4 - 8*x*s2^2
  + x^2 + 7*s2^2 - 3*x + 1,
 Affine Curve over Rational Field defined by s3, s2 - 1,
  y^2*s0^4 + 2*y^2*s0^2 + y*s0^3 + y^2 - s0^2 - 1,
 Affine Curve over Rational Field defined by s3, s1 - 1,
  z^2*s0^4 - 14*z*s0^4 + 2*z^2*s0^2 + z*s0^3 + 49*s0^4
  - 28*z*s0^2 - 7*s0^3 + z^2 + 97*s0^2 - 14*z + 48,
 Closed subscheme of Affine Space of dimension 4 over Rational Field
  defined by: 1)
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(4), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = A._first_ngens(4)
>>> C = A.curve([((x - Integer(2))**Integer(2) + y**Integer(2))**Integer(2) - (x - Integer(2))**Integer(2) - y**Integer(2) + (x - Integer(2))**Integer(3),
...              z - y - Integer(7), w - Integer(4)])
>>> B = C.resolution_of_singularities()
>>> B[Integer(0)]
(Affine Curve over Rational Field defined by s3, s1 - s2,
  x^2*s2^4 - 4*x*s2^4 + 2*x^2*s2^2 + 4*s2^4 - 8*x*s2^2
  + x^2 + 7*s2^2 - 3*x + 1,
 Affine Curve over Rational Field defined by s3, s2 - 1,
  y^2*s0^4 + 2*y^2*s0^2 + y*s0^3 + y^2 - s0^2 - 1,
 Affine Curve over Rational Field defined by s3, s1 - 1,
  z^2*s0^4 - 14*z*s0^4 + 2*z^2*s0^2 + z*s0^3 + 49*s0^4
  - 28*z*s0^2 - 7*s0^3 + z^2 + 97*s0^2 - 14*z + 48,
 Closed subscheme of Affine Space of dimension 4 over Rational Field
  defined by: 1)

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve([y - x^2 + 1], A)
sage: C.resolution_of_singularities()
Traceback (most recent call last):
...
TypeError: this curve is already nonsingular
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([y - x**Integer(2) + Integer(1)], A)
>>> C.resolution_of_singularities()
Traceback (most recent call last):
...
TypeError: this curve is already nonsingular

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = A.curve([(x^2 + y^2 - y - 2)*(y - x^2 + 2) + y^3])
sage: C.resolution_of_singularities()
Traceback (most recent call last):
...
TypeError: this curve has no singular points over its base field. If
working over a number field use extend=True
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = A.curve([(x**Integer(2) + y**Integer(2) - y - Integer(2))*(y - x**Integer(2) + Integer(2)) + y**Integer(3)])
>>> C.resolution_of_singularities()
Traceback (most recent call last):
...
TypeError: this curve has no singular points over its base field. If
working over a number field use extend=True
tangent_line(p)[source]

Return the tangent line at the point p.

INPUT:

  • p – a rational point of the curve

EXAMPLES:

sage: A3.<x,y,z> = AffineSpace(3, QQ)
sage: C = Curve([x + y + z, x^2 - y^2*z^2 + z^3])
sage: p = C(0,0,0)
sage: C.tangent_line(p)
Traceback (most recent call last):
...
ValueError: the curve is not smooth at (0, 0, 0)
sage: p = C(1,0,-1)
sage: C.tangent_line(p)
Affine Curve over Rational Field defined by x + y + z, 2*x + 3*z + 1
>>> from sage.all import *
>>> A3 = AffineSpace(Integer(3), QQ, names=('x', 'y', 'z',)); (x, y, z,) = A3._first_ngens(3)
>>> C = Curve([x + y + z, x**Integer(2) - y**Integer(2)*z**Integer(2) + z**Integer(3)])
>>> p = C(Integer(0),Integer(0),Integer(0))
>>> C.tangent_line(p)
Traceback (most recent call last):
...
ValueError: the curve is not smooth at (0, 0, 0)
>>> p = C(Integer(1),Integer(0),-Integer(1))
>>> C.tangent_line(p)
Affine Curve over Rational Field defined by x + y + z, 2*x + 3*z + 1

We check that the tangent line at p is the tangent space at p, translated to p.

sage: Tp = C.tangent_space(p)
sage: Tp
Closed subscheme of Affine Space of dimension 3 over Rational Field
 defined by: x + y + z, 2*x + 3*z
sage: phi = A3.translation(A3.origin(), p)
sage: T = phi * Tp.embedding_morphism()
sage: T.image()
Closed subscheme of Affine Space of dimension 3 over Rational Field
 defined by: -2*y + z + 1, x + y + z
sage: _ == C.tangent_line(p)
True
>>> from sage.all import *
>>> Tp = C.tangent_space(p)
>>> Tp
Closed subscheme of Affine Space of dimension 3 over Rational Field
 defined by: x + y + z, 2*x + 3*z
>>> phi = A3.translation(A3.origin(), p)
>>> T = phi * Tp.embedding_morphism()
>>> T.image()
Closed subscheme of Affine Space of dimension 3 over Rational Field
 defined by: -2*y + z + 1, x + y + z
>>> _ == C.tangent_line(p)
True
class sage.schemes.curves.affine_curve.AffinePlaneCurve(A, f)[source]

Bases: AffineCurve

Affine plane curves.

divisor_of_function(r)[source]

Return the divisor of a function on a curve.

INPUT:

  • r – a rational function on X

OUTPUT:

  • list – the divisor of r represented as a list of coefficients and points. (TODO: This will change to a more structural output in the future.)

EXAMPLES:

sage: F = GF(5)
sage: P2 = AffineSpace(2, F, names='xy')
sage: R = P2.coordinate_ring()
sage: x, y = R.gens()
sage: f = y^2 - x^9 - x
sage: C = Curve(f)
sage: K = FractionField(R)
sage: r = 1/x
sage: C.divisor_of_function(r)      # not implemented (broken)
      [[-1, (0, 0, 1)]]
sage: r = 1/x^3
sage: C.divisor_of_function(r)      # not implemented (broken)
      [[-3, (0, 0, 1)]]
>>> from sage.all import *
>>> F = GF(Integer(5))
>>> P2 = AffineSpace(Integer(2), F, names='xy')
>>> R = P2.coordinate_ring()
>>> x, y = R.gens()
>>> f = y**Integer(2) - x**Integer(9) - x
>>> C = Curve(f)
>>> K = FractionField(R)
>>> r = Integer(1)/x
>>> C.divisor_of_function(r)      # not implemented (broken)
      [[-1, (0, 0, 1)]]
>>> r = Integer(1)/x**Integer(3)
>>> C.divisor_of_function(r)      # not implemented (broken)
      [[-3, (0, 0, 1)]]
is_ordinary_singularity(P)[source]

Return whether the singular point P of this affine plane curve is an ordinary singularity.

The point P is an ordinary singularity of this curve if it is a singular point, and if the tangents of this curve at P are distinct.

INPUT:

  • P – a point on this curve

OUTPUT:

True or False depending on whether P is or is not an ordinary singularity of this curve, respectively. An error is raised if P is not a singular point of this curve.

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve([y^2 - x^3], A)
sage: Q = A([0,0])
sage: C.is_ordinary_singularity(Q)
False
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([y**Integer(2) - x**Integer(3)], A)
>>> Q = A([Integer(0),Integer(0)])
>>> C.is_ordinary_singularity(Q)
False

sage: # needs sage.rings.number_field
sage: R.<a> = QQ[]
sage: K.<b> = NumberField(a^2 - 3)
sage: A.<x,y> = AffineSpace(K, 2)
sage: C = Curve([(x^2 + y^2 - 2*x)^2 - x^2 - y^2], A)
sage: Q = A([0,0])
sage: C.is_ordinary_singularity(Q)
True
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> R = QQ['a']; (a,) = R._first_ngens(1)
>>> K = NumberField(a**Integer(2) - Integer(3), names=('b',)); (b,) = K._first_ngens(1)
>>> A = AffineSpace(K, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([(x**Integer(2) + y**Integer(2) - Integer(2)*x)**Integer(2) - x**Integer(2) - y**Integer(2)], A)
>>> Q = A([Integer(0),Integer(0)])
>>> C.is_ordinary_singularity(Q)
True

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = A.curve([x^2*y - y^2*x + y^2 + x^3])
sage: Q = A([-1,-1])
sage: C.is_ordinary_singularity(Q)
Traceback (most recent call last):
...
TypeError: (=(-1, -1)) is not a singular point of (=Affine Plane Curve
over Rational Field defined by x^3 + x^2*y - x*y^2 + y^2)
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = A.curve([x**Integer(2)*y - y**Integer(2)*x + y**Integer(2) + x**Integer(3)])
>>> Q = A([-Integer(1),-Integer(1)])
>>> C.is_ordinary_singularity(Q)
Traceback (most recent call last):
...
TypeError: (=(-1, -1)) is not a singular point of (=Affine Plane Curve
over Rational Field defined by x^3 + x^2*y - x*y^2 + y^2)
is_transverse(C, P)[source]

Return whether the intersection of this curve with the curve C at the point P is transverse.

The intersection at P is transverse if P is a nonsingular point of both curves, and if the tangents of the curves at P are distinct.

INPUT:

  • C – a curve in the ambient space of this curve

  • P – a point in the intersection of both curves

OUTPUT: boolean

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve([x^2 + y^2 - 1], A)
sage: D = Curve([x - 1], A)
sage: Q = A([1,0])
sage: C.is_transverse(D, Q)
False
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([x**Integer(2) + y**Integer(2) - Integer(1)], A)
>>> D = Curve([x - Integer(1)], A)
>>> Q = A([Integer(1),Integer(0)])
>>> C.is_transverse(D, Q)
False

sage: # needs sage.rings.number_field
sage: R.<a> = QQ[]
sage: K.<b> = NumberField(a^3 + 2)
sage: A.<x,y> = AffineSpace(K, 2)
sage: C = A.curve([x*y])
sage: D = A.curve([y - b*x])
sage: Q = A([0,0])
sage: C.is_transverse(D, Q)
False
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> R = QQ['a']; (a,) = R._first_ngens(1)
>>> K = NumberField(a**Integer(3) + Integer(2), names=('b',)); (b,) = K._first_ngens(1)
>>> A = AffineSpace(K, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = A.curve([x*y])
>>> D = A.curve([y - b*x])
>>> Q = A([Integer(0),Integer(0)])
>>> C.is_transverse(D, Q)
False

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve([y - x^3], A)
sage: D = Curve([y + x], A)
sage: Q = A([0,0])
sage: C.is_transverse(D, Q)
True
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([y - x**Integer(3)], A)
>>> D = Curve([y + x], A)
>>> Q = A([Integer(0),Integer(0)])
>>> C.is_transverse(D, Q)
True
local_coordinates(pt, n)[source]

Return local coordinates to precision n at the given point.

Behaviour is flaky - some choices of \(n\) are worst that others.

INPUT:

  • pt – an F-rational point on X which is not a point of ramification for the projection (x,y) - x

  • n – the number of terms desired

OUTPUT: x = x0 + t y = y0 + power series in t

EXAMPLES:

sage: F = GF(5)
sage: pt = (2,3)
sage: R = PolynomialRing(F, 2, names = ['x','y'])
sage: x,y = R.gens()
sage: f = y^2 - x^9 - x
sage: C = Curve(f)
sage: C.local_coordinates(pt, 9)
[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]
>>> from sage.all import *
>>> F = GF(Integer(5))
>>> pt = (Integer(2),Integer(3))
>>> R = PolynomialRing(F, Integer(2), names = ['x','y'])
>>> x,y = R.gens()
>>> f = y**Integer(2) - x**Integer(9) - x
>>> C = Curve(f)
>>> C.local_coordinates(pt, Integer(9))
[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]
multiplicity(P)[source]

Return the multiplicity of this affine plane curve at the point P.

In the special case of affine plane curves, the multiplicity of an affine plane curve at the point (0,0) can be computed as the minimum of the degrees of the homogeneous components of its defining polynomial. To compute the multiplicity of a different point, a linear change of coordinates is used.

This curve must be defined over a field. An error if raised if P is not a point on this curve.

INPUT:

  • P – a point in the ambient space of this curve

OUTPUT: integer

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve([y^2 - x^3], A)
sage: Q1 = A([1,1])
sage: C.multiplicity(Q1)
1
sage: Q2 = A([0,0])
sage: C.multiplicity(Q2)
2
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([y**Integer(2) - x**Integer(3)], A)
>>> Q1 = A([Integer(1),Integer(1)])
>>> C.multiplicity(Q1)
1
>>> Q2 = A([Integer(0),Integer(0)])
>>> C.multiplicity(Q2)
2

sage: # needs sage.rings.number_field
sage: A.<x,y> = AffineSpace(QQbar,2)
sage: C = Curve([-x^7 + (-7)*x^6 + y^6 + (-21)*x^5 + 12*y^5
....:            + (-35)*x^4 + 60*y^4 + (-35)*x^3 + 160*y^3
....:            + (-21)*x^2 + 240*y^2 + (-7)*x + 192*y + 63], A)
sage: Q = A([-1,-2])
sage: C.multiplicity(Q)
6
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> A = AffineSpace(QQbar,Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([-x**Integer(7) + (-Integer(7))*x**Integer(6) + y**Integer(6) + (-Integer(21))*x**Integer(5) + Integer(12)*y**Integer(5)
...            + (-Integer(35))*x**Integer(4) + Integer(60)*y**Integer(4) + (-Integer(35))*x**Integer(3) + Integer(160)*y**Integer(3)
...            + (-Integer(21))*x**Integer(2) + Integer(240)*y**Integer(2) + (-Integer(7))*x + Integer(192)*y + Integer(63)], A)
>>> Q = A([-Integer(1),-Integer(2)])
>>> C.multiplicity(Q)
6

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = A.curve([y^3 - x^3 + x^6])
sage: Q = A([1,1])
sage: C.multiplicity(Q)
Traceback (most recent call last):
...
TypeError: (=(1, 1)) is not a point on (=Affine Plane Curve over
Rational Field defined by x^6 - x^3 + y^3)
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = A.curve([y**Integer(3) - x**Integer(3) + x**Integer(6)])
>>> Q = A([Integer(1),Integer(1)])
>>> C.multiplicity(Q)
Traceback (most recent call last):
...
TypeError: (=(1, 1)) is not a point on (=Affine Plane Curve over
Rational Field defined by x^6 - x^3 + y^3)
plot(*args, **kwds)[source]

Plot the real points on this affine plane curve.

INPUT:

  • *args – (optional) tuples (variable, minimum, maximum) for plotting dimensions

  • **kwds – optional keyword arguments passed on to implicit_plot

EXAMPLES:

A cuspidal curve:

sage: R.<x, y> = QQ[]
sage: C = Curve(x^3 - y^2)
sage: C.plot()                                      # needs sage.plot
Graphics object consisting of 1 graphics primitive
>>> from sage.all import *
>>> R = QQ['x, y']; (x, y,) = R._first_ngens(2)
>>> C = Curve(x**Integer(3) - y**Integer(2))
>>> C.plot()                                      # needs sage.plot
Graphics object consisting of 1 graphics primitive

A 5-nodal curve of degree 11. This example also illustrates some of the optional arguments:

sage: # needs sage.plot
sage: R.<x, y> = ZZ[]
sage: C = Curve(32*x^2 - 2097152*y^11 + 1441792*y^9
....:            - 360448*y^7 + 39424*y^5 - 1760*y^3 + 22*y - 1)
sage: C.plot((x, -1, 1), (y, -1, 1), plot_points=400)
Graphics object consisting of 1 graphics primitive
>>> from sage.all import *
>>> # needs sage.plot
>>> R = ZZ['x, y']; (x, y,) = R._first_ngens(2)
>>> C = Curve(Integer(32)*x**Integer(2) - Integer(2097152)*y**Integer(11) + Integer(1441792)*y**Integer(9)
...            - Integer(360448)*y**Integer(7) + Integer(39424)*y**Integer(5) - Integer(1760)*y**Integer(3) + Integer(22)*y - Integer(1))
>>> C.plot((x, -Integer(1), Integer(1)), (y, -Integer(1), Integer(1)), plot_points=Integer(400))
Graphics object consisting of 1 graphics primitive

A line over \(\mathbf{RR}\):

sage: # needs sage.symbolic sage.plot
sage: R.<x, y> = RR[]
sage: C = Curve(R(y - sqrt(2)*x))
sage: C.plot()
Graphics object consisting of 1 graphics primitive
>>> from sage.all import *
>>> # needs sage.symbolic sage.plot
>>> R = RR['x, y']; (x, y,) = R._first_ngens(2)
>>> C = Curve(R(y - sqrt(Integer(2))*x))
>>> C.plot()
Graphics object consisting of 1 graphics primitive
rational_parameterization()[source]

Return a rational parameterization of this curve.

This curve must have rational coefficients and be absolutely irreducible (i.e. irreducible over the algebraic closure of the rational field). The curve must also be rational (have geometric genus zero).

The rational parameterization may have coefficients in a quadratic extension of the rational field.

OUTPUT: a birational map between \(\mathbb{A}^{1}\) and this curve, given as a scheme morphism

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve([y^2 - x], A)
sage: C.rational_parameterization()
Scheme morphism:
  From: Affine Space of dimension 1 over Rational Field
  To:   Affine Plane Curve over Rational Field defined by y^2 - x
  Defn: Defined on coordinates by sending (t) to
        (t^2, t)
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([y**Integer(2) - x], A)
>>> C.rational_parameterization()
Scheme morphism:
  From: Affine Space of dimension 1 over Rational Field
  To:   Affine Plane Curve over Rational Field defined by y^2 - x
  Defn: Defined on coordinates by sending (t) to
        (t^2, t)

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve([(x^2 + y^2 - 2*x)^2 - x^2 - y^2], A)
sage: C.rational_parameterization()
Scheme morphism:
  From: Affine Space of dimension 1 over Rational Field
  To:   Affine Plane Curve over Rational Field defined by x^4 +
2*x^2*y^2 + y^4 - 4*x^3 - 4*x*y^2 + 3*x^2 - y^2
  Defn: Defined on coordinates by sending (t) to
        ((-12*t^4 + 6*t^3 + 4*t^2 - 2*t)/(-25*t^4 + 40*t^3 - 26*t^2 +
8*t - 1), (-9*t^4 + 12*t^3 - 4*t + 1)/(-25*t^4 + 40*t^3 - 26*t^2 + 8*t - 1))
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([(x**Integer(2) + y**Integer(2) - Integer(2)*x)**Integer(2) - x**Integer(2) - y**Integer(2)], A)
>>> C.rational_parameterization()
Scheme morphism:
  From: Affine Space of dimension 1 over Rational Field
  To:   Affine Plane Curve over Rational Field defined by x^4 +
2*x^2*y^2 + y^4 - 4*x^3 - 4*x*y^2 + 3*x^2 - y^2
  Defn: Defined on coordinates by sending (t) to
        ((-12*t^4 + 6*t^3 + 4*t^2 - 2*t)/(-25*t^4 + 40*t^3 - 26*t^2 +
8*t - 1), (-9*t^4 + 12*t^3 - 4*t + 1)/(-25*t^4 + 40*t^3 - 26*t^2 + 8*t - 1))

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve([x^2 + y^2 + 7], A)
sage: C.rational_parameterization()
Scheme morphism:
  From: Affine Space of dimension 1 over Number Field in a with defining polynomial a^2 + 7
  To:   Affine Plane Curve over Number Field in a with defining
polynomial a^2 + 7 defined by x^2 + y^2 + 7
  Defn: Defined on coordinates by sending (t) to
        ((-7*t^2 + 7)/((-a)*t^2 + (-a)), 14*t/((-a)*t^2 + (-a)))
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([x**Integer(2) + y**Integer(2) + Integer(7)], A)
>>> C.rational_parameterization()
Scheme morphism:
  From: Affine Space of dimension 1 over Number Field in a with defining polynomial a^2 + 7
  To:   Affine Plane Curve over Number Field in a with defining
polynomial a^2 + 7 defined by x^2 + y^2 + 7
  Defn: Defined on coordinates by sending (t) to
        ((-7*t^2 + 7)/((-a)*t^2 + (-a)), 14*t/((-a)*t^2 + (-a)))
tangents(P, factor=True)[source]

Return the tangents of this affine plane curve at the point P.

The point P must be a point on this curve.

INPUT:

  • P – a point on this curve

  • factor – boolean (default: True); whether to attempt computing the polynomials of the individual tangent lines over the base field of this curve, or to just return the polynomial corresponding to the union of the tangent lines (which requires fewer computations)

OUTPUT: list of polynomials in the coordinate ring of the ambient space

EXAMPLES:

sage: # needs sage.rings.number_field
sage: set_verbose(-1)
sage: A.<x,y> = AffineSpace(QQbar, 2)
sage: C = Curve([x^5*y^3 + 2*x^4*y^4 + x^3*y^5 + 3*x^4*y^3
....:            + 6*x^3*y^4 + 3*x^2*y^5 + 3*x^3*y^3
....:            + 6*x^2*y^4 + 3*x*y^5 + x^5 + 10*x^4*y
....:            + 40*x^3*y^2 + 81*x^2*y^3 + 82*x*y^4 + 33*y^5], A)
sage: Q = A([0,0])
sage: C.tangents(Q)
[x + 3.425299577684700?*y,
 x + (1.949159013086856? + 1.179307909383728?*I)*y,
 x + (1.949159013086856? - 1.179307909383728?*I)*y,
 x + (1.338191198070795? + 0.2560234251008043?*I)*y,
 x + (1.338191198070795? - 0.2560234251008043?*I)*y]
sage: C.tangents(Q, factor=False)
[120*x^5 + 1200*x^4*y + 4800*x^3*y^2 + 9720*x^2*y^3 + 9840*x*y^4 + 3960*y^5]
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> set_verbose(-Integer(1))
>>> A = AffineSpace(QQbar, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([x**Integer(5)*y**Integer(3) + Integer(2)*x**Integer(4)*y**Integer(4) + x**Integer(3)*y**Integer(5) + Integer(3)*x**Integer(4)*y**Integer(3)
...            + Integer(6)*x**Integer(3)*y**Integer(4) + Integer(3)*x**Integer(2)*y**Integer(5) + Integer(3)*x**Integer(3)*y**Integer(3)
...            + Integer(6)*x**Integer(2)*y**Integer(4) + Integer(3)*x*y**Integer(5) + x**Integer(5) + Integer(10)*x**Integer(4)*y
...            + Integer(40)*x**Integer(3)*y**Integer(2) + Integer(81)*x**Integer(2)*y**Integer(3) + Integer(82)*x*y**Integer(4) + Integer(33)*y**Integer(5)], A)
>>> Q = A([Integer(0),Integer(0)])
>>> C.tangents(Q)
[x + 3.425299577684700?*y,
 x + (1.949159013086856? + 1.179307909383728?*I)*y,
 x + (1.949159013086856? - 1.179307909383728?*I)*y,
 x + (1.338191198070795? + 0.2560234251008043?*I)*y,
 x + (1.338191198070795? - 0.2560234251008043?*I)*y]
>>> C.tangents(Q, factor=False)
[120*x^5 + 1200*x^4*y + 4800*x^3*y^2 + 9720*x^2*y^3 + 9840*x*y^4 + 3960*y^5]

sage: # needs sage.rings.number_field
sage: R.<a> = QQ[]
sage: K.<b> = NumberField(a^2 - 3)
sage: A.<x,y> = AffineSpace(K, 2)
sage: C = Curve([(x^2 + y^2 - 2*x)^2 - x^2 - y^2], A)
sage: Q = A([0,0])
sage: C.tangents(Q)
[x + (-1/3*b)*y, x + (1/3*b)*y]
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> R = QQ['a']; (a,) = R._first_ngens(1)
>>> K = NumberField(a**Integer(2) - Integer(3), names=('b',)); (b,) = K._first_ngens(1)
>>> A = AffineSpace(K, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve([(x**Integer(2) + y**Integer(2) - Integer(2)*x)**Integer(2) - x**Integer(2) - y**Integer(2)], A)
>>> Q = A([Integer(0),Integer(0)])
>>> C.tangents(Q)
[x + (-1/3*b)*y, x + (1/3*b)*y]

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = A.curve([y^2 - x^3 - x^2])
sage: Q = A([0,0])
sage: C.tangents(Q)
[x - y, x + y]
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = A.curve([y**Integer(2) - x**Integer(3) - x**Integer(2)])
>>> Q = A([Integer(0),Integer(0)])
>>> C.tangents(Q)
[x - y, x + y]

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = A.curve([y*x - x^4 + 2*x^2])
sage: Q = A([1,1])
sage: C.tangents(Q)
Traceback (most recent call last):
...
TypeError: (=(1, 1)) is not a point on (=Affine Plane Curve over
Rational Field defined by -x^4 + 2*x^2 + x*y)
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = A.curve([y*x - x**Integer(4) + Integer(2)*x**Integer(2)])
>>> Q = A([Integer(1),Integer(1)])
>>> C.tangents(Q)
Traceback (most recent call last):
...
TypeError: (=(1, 1)) is not a point on (=Affine Plane Curve over
Rational Field defined by -x^4 + 2*x^2 + x*y)
class sage.schemes.curves.affine_curve.AffinePlaneCurve_field(A, f)[source]

Bases: AffinePlaneCurve, AffineCurve_field

Affine plane curves over fields.

braid_monodromy()[source]

Compute the braid monodromy of a projection of the curve.

OUTPUT:

A list of braids. The braids correspond to paths based in the same point; each of this paths is the conjugated of a loop around one of the points in the discriminant of the projection of self.

Note

The projection over the \(x\) axis is used if there are no vertical asymptotes. Otherwise, a linear change of variables is done to fall into the previous case.

Note

This functionality requires the sirocco package to be installed.

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = A.curve((x^2-y^3)*(x+3*y-5))
sage: C.braid_monodromy()                                   # needs sirocco
[s1*s0*(s1*s2)^2*s0*s2^2*s0^-1*(s2^-1*s1^-1)^2*s0^-1*s1^-1,
 s1*s0*(s1*s2)^2*(s0*s2^-1*s1*s2*s1*s2^-1)^2*(s2^-1*s1^-1)^2*s0^-1*s1^-1,
 s1*s0*(s1*s2)^2*s2*s1^-1*s2^-1*s1^-1*s0^-1*s1^-1,
 s1*s0*s2*s0^-1*s2*s1^-1]
sage: T.<t> = QQ[]
sage: K.<a> = NumberField(t^3 + 2, 'a')
sage: A.<x, y> = AffineSpace(K, 2)
sage: Curve(y^2 + a * x).braid_monodromy()
Traceback (most recent call last):
...
NotImplementedError: the base field must have an embedding to the algebraic field
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = A.curve((x**Integer(2)-y**Integer(3))*(x+Integer(3)*y-Integer(5)))
>>> C.braid_monodromy()                                   # needs sirocco
[s1*s0*(s1*s2)^2*s0*s2^2*s0^-1*(s2^-1*s1^-1)^2*s0^-1*s1^-1,
 s1*s0*(s1*s2)^2*(s0*s2^-1*s1*s2*s1*s2^-1)^2*(s2^-1*s1^-1)^2*s0^-1*s1^-1,
 s1*s0*(s1*s2)^2*s2*s1^-1*s2^-1*s1^-1*s0^-1*s1^-1,
 s1*s0*s2*s0^-1*s2*s1^-1]
>>> T = QQ['t']; (t,) = T._first_ngens(1)
>>> K = NumberField(t**Integer(3) + Integer(2), 'a', names=('a',)); (a,) = K._first_ngens(1)
>>> A = AffineSpace(K, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> Curve(y**Integer(2) + a * x).braid_monodromy()
Traceback (most recent call last):
...
NotImplementedError: the base field must have an embedding to the algebraic field
fundamental_group(simplified=True, puiseux=True)[source]

Return a presentation of the fundamental group of the complement of self.

INPUT:

  • simplified – boolean (default: True); to simplify the presentation

  • puiseux – boolean (default: True); to decide if the presentation is constructed in the classical way or using Puiseux shortcut

OUTPUT:

A presentation with generators \(x_1, \dots, x_d\) and relations. If puiseux is False the relations are \((x_j\cdot \tau)\cdot x_j^{-1}\) for \(1\leq j<d\) and \(tau\) a braid in the braid monodromy; finally the presentation is simplified. If puiseux is True, each \(tau\) is decomposed as \(\alpha^{-1}\cdot\beta\cdot\alpha\), where \(\beta\) is a positive braid; the relations are \(((x_j\cdot \beta)\cdot x_j^{-1})\cdot \alpha\) where \(j\) is an integer of the Tietze word of \(\beta\). This presentation is not simplified by default since it represents the homotopy type of the complement of the curve.

Note

The curve must be defined over the rationals or a number field with an embedding over \(\QQbar\). This functionality requires the sirocco package to be installed.

EXAMPLES:

sage: # needs sirocco
sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = A.curve(y^2 - x^3 - x^2)
sage: C.fundamental_group(puiseux=False)
Finitely presented group < x0 |  >
sage: bm = C.braid_monodromy()
sage: g = C.fundamental_group(simplified=False)
sage: g.sorted_presentation()
Finitely presented group < x0, x1 | x1^-1*x0^-1*x1*x0, x1^-1*x0 >
sage: g.simplified()
Finitely presented group < x0 |  >
>>> from sage.all import *
>>> # needs sirocco
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = A.curve(y**Integer(2) - x**Integer(3) - x**Integer(2))
>>> C.fundamental_group(puiseux=False)
Finitely presented group < x0 |  >
>>> bm = C.braid_monodromy()
>>> g = C.fundamental_group(simplified=False)
>>> g.sorted_presentation()
Finitely presented group < x0, x1 | x1^-1*x0^-1*x1*x0, x1^-1*x0 >
>>> g.simplified()
Finitely presented group < x0 |  >

In the case of number fields, they need to have an embedding to the algebraic field:

sage: # needs sage.rings.number_field
sage: x = polygen(ZZ)
sage: a = QQ[x](x^2 + 5).roots(QQbar)[0][0]
sage: F = NumberField(a.minpoly(), 'a', embedding=a)
sage: F.inject_variables()
Defining a
sage: A.<x,y> = AffineSpace(F, 2)
sage: C = A.curve(y^2 - a*x^3 - x^2)
sage: C.fundamental_group()                     # needs sirocco
Finitely presented group < x0 |  >
sage: C = A.curve(x * (x - 1))
sage: C.fundamental_group()                     # needs sirocco
Finitely presented group < x0, x1 |  >
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> x = polygen(ZZ)
>>> a = QQ[x](x**Integer(2) + Integer(5)).roots(QQbar)[Integer(0)][Integer(0)]
>>> F = NumberField(a.minpoly(), 'a', embedding=a)
>>> F.inject_variables()
Defining a
>>> A = AffineSpace(F, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = A.curve(y**Integer(2) - a*x**Integer(3) - x**Integer(2))
>>> C.fundamental_group()                     # needs sirocco
Finitely presented group < x0 |  >
>>> C = A.curve(x * (x - Integer(1)))
>>> C.fundamental_group()                     # needs sirocco
Finitely presented group < x0, x1 |  >
has_vertical_asymptote()[source]

Check if the curve is not a line and has vertical asymptotes.

EXAMPLES:

sage: A2.<x,y> = AffineSpace(2, QQ)
sage: Curve(x).has_vertical_asymptote()
False
sage: Curve(y^2 * x + x + y).has_vertical_asymptote()
True
>>> from sage.all import *
>>> A2 = AffineSpace(Integer(2), QQ, names=('x', 'y',)); (x, y,) = A2._first_ngens(2)
>>> Curve(x).has_vertical_asymptote()
False
>>> Curve(y**Integer(2) * x + x + y).has_vertical_asymptote()
True
is_vertical_line()[source]

Check if the curve is a vertical line.

EXAMPLES:

sage: A2.<x, y> = AffineSpace(2, QQ)
sage: Curve(x - 1).is_vertical_line()
True
sage: Curve(x - y).is_vertical_line()
False
sage: Curve(y^2 * x + x + y).is_vertical_line()
False
>>> from sage.all import *
>>> A2 = AffineSpace(Integer(2), QQ, names=('x', 'y',)); (x, y,) = A2._first_ngens(2)
>>> Curve(x - Integer(1)).is_vertical_line()
True
>>> Curve(x - y).is_vertical_line()
False
>>> Curve(y**Integer(2) * x + x + y).is_vertical_line()
False
riemann_surface(**kwargs)[source]

Return the complex Riemann surface determined by this curve.

OUTPUT: a RiemannSurface object

EXAMPLES:

sage: R.<x,y> = QQ[]
sage: C = Curve(x^3 + 3*y^3 + 5)
sage: C.riemann_surface()
Riemann surface defined by polynomial f = x^3 + 3*y^3 + 5 = 0,
 with 53 bits of precision
>>> from sage.all import *
>>> R = QQ['x, y']; (x, y,) = R._first_ngens(2)
>>> C = Curve(x**Integer(3) + Integer(3)*y**Integer(3) + Integer(5))
>>> C.riemann_surface()
Riemann surface defined by polynomial f = x^3 + 3*y^3 + 5 = 0,
 with 53 bits of precision
class sage.schemes.curves.affine_curve.AffinePlaneCurve_finite_field(A, f)[source]

Bases: AffinePlaneCurve_field

Affine plane curves over finite fields.

rational_points(algorithm='enum')[source]

Return sorted list of all rational points on this curve.

INPUT:

  • algorithm – possible choices:

    • 'enum' – use very naive point enumeration to find all rational points on this curve over a finite field.

    • 'bn' – via Singular’s Brill-Noether package.

    • 'all' – use all implemented algorithms and verify that they give the same answer, then return it

Note

The Brill-Noether package does not always work. When it fails, a RuntimeError exception is raised.

EXAMPLES:

sage: x, y = (GF(5)['x,y']).gens()
sage: f = y^2 - x^9 - x
sage: C = Curve(f); C
Affine Plane Curve over Finite Field of size 5 defined by -x^9 + y^2 - x
sage: C.rational_points(algorithm='bn')
[(0, 0), (2, 2), (2, 3), (3, 1), (3, 4)]
sage: C = Curve(x - y + 1)
sage: C.rational_points()
[(0, 1), (1, 2), (2, 3), (3, 4), (4, 0)]
>>> from sage.all import *
>>> x, y = (GF(Integer(5))['x,y']).gens()
>>> f = y**Integer(2) - x**Integer(9) - x
>>> C = Curve(f); C
Affine Plane Curve over Finite Field of size 5 defined by -x^9 + y^2 - x
>>> C.rational_points(algorithm='bn')
[(0, 0), (2, 2), (2, 3), (3, 1), (3, 4)]
>>> C = Curve(x - y + Integer(1))
>>> C.rational_points()
[(0, 1), (1, 2), (2, 3), (3, 4), (4, 0)]

We compare Brill-Noether and enumeration:

sage: x, y = (GF(17)['x,y']).gens()
sage: C = Curve(x^2 + y^5 + x*y - 19)
sage: v = C.rational_points(algorithm='bn')
sage: w = C.rational_points(algorithm='enum')
sage: len(v)
20
sage: v == w
True

sage: # needs sage.rings.finite_rings
sage: A.<x,y> = AffineSpace(2, GF(9,'a'))
sage: C = Curve(x^2 + y^2 - 1); C
Affine Plane Curve over Finite Field in a of size 3^2
 defined by x^2 + y^2 - 1
sage: C.rational_points()
[(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)]
>>> from sage.all import *
>>> x, y = (GF(Integer(17))['x,y']).gens()
>>> C = Curve(x**Integer(2) + y**Integer(5) + x*y - Integer(19))
>>> v = C.rational_points(algorithm='bn')
>>> w = C.rational_points(algorithm='enum')
>>> len(v)
20
>>> v == w
True

>>> # needs sage.rings.finite_rings
>>> A = AffineSpace(Integer(2), GF(Integer(9),'a'), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(2) + y**Integer(2) - Integer(1)); C
Affine Plane Curve over Finite Field in a of size 3^2
 defined by x^2 + y^2 - 1
>>> C.rational_points()
[(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)]
riemann_roch_basis(D)[source]

Return a basis of the Riemann-Roch space of the divisor D.

This interfaces with Singular’s Brill-Noether command.

This curve is assumed to be a plane curve defined by a polynomial equation \(f(x,y) = 0\) over a prime finite field \(F = GF(p)\) in 2 variables \(x,y\) representing a curve \(X: f(x,y) = 0\) having \(n\) \(F\)-rational points (see the Sage function places_on_curve)

INPUT:

  • D – an \(n\)-tuple of integers \((d_1, ..., d_n)\) representing the divisor \(Div = d_1P_1 + \dots + d_nP_n\), where \(X(F) = \{P_1, \dots, P_n\}\). The ordering is that dictated by places_on_curve.

OUTPUT: a basis of \(L(Div)\).

EXAMPLES:

sage: R = PolynomialRing(GF(5), 2, names=["x","y"])
sage: x, y = R.gens()
sage: f = y^2 - x^9 - x
sage: C = Curve(f)
sage: D = [6,0,0,0,0,0]
sage: C.riemann_roch_basis(D)
[1, (-x*z^5 + y^2*z^4)/x^6, (-x*z^6 + y^2*z^5)/x^7, (-x*z^7 + y^2*z^6)/x^8]
>>> from sage.all import *
>>> R = PolynomialRing(GF(Integer(5)), Integer(2), names=["x","y"])
>>> x, y = R.gens()
>>> f = y**Integer(2) - x**Integer(9) - x
>>> C = Curve(f)
>>> D = [Integer(6),Integer(0),Integer(0),Integer(0),Integer(0),Integer(0)]
>>> C.riemann_roch_basis(D)
[1, (-x*z^5 + y^2*z^4)/x^6, (-x*z^6 + y^2*z^5)/x^7, (-x*z^7 + y^2*z^6)/x^8]
class sage.schemes.curves.affine_curve.IntegralAffineCurve(A, X)[source]

Bases: AffineCurve_field

Base class for integral affine curves.

coordinate_functions()[source]

Return the coordinate functions.

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: A.<x,y> = AffineSpace(GF(8), 2)
sage: C = Curve(x^5 + y^5 + x*y + 1)
sage: x, y = C.coordinate_functions()
sage: x^5 + y^5 + x*y + 1
0
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> A = AffineSpace(GF(Integer(8)), Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(5) + y**Integer(5) + x*y + Integer(1))
>>> x, y = C.coordinate_functions()
>>> x**Integer(5) + y**Integer(5) + x*y + Integer(1)
0
function(f)[source]

Return the function field element coerced from f.

INPUT:

  • f – an element of the fraction field of the coordinate ring of the ambient space or the coordinate ring of the curve

OUTPUT: an element of the function field of this curve

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: A.<x,y> = AffineSpace(GF(8), 2)
sage: C = Curve(x^5 + y^5 + x*y + 1)
sage: f = C.function(x/y)
sage: f
(x/(x^5 + 1))*y^4 + x^2/(x^5 + 1)
sage: df = f.differential(); df
((1/(x^10 + 1))*y^4 + x^6/(x^10 + 1)) d(x)
sage: df.divisor()                  # long time
2*Place (1/x, 1/x^4*y^4 + 1/x^3*y^3 + 1/x^2*y^2 + 1/x*y + 1)
 + 2*Place (1/x, 1/x*y + 1)
 - 2*Place (x + 1, y)
 - 2*Place (x^4 + x^3 + x^2 + x + 1, y)
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> A = AffineSpace(GF(Integer(8)), Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(5) + y**Integer(5) + x*y + Integer(1))
>>> f = C.function(x/y)
>>> f
(x/(x^5 + 1))*y^4 + x^2/(x^5 + 1)
>>> df = f.differential(); df
((1/(x^10 + 1))*y^4 + x^6/(x^10 + 1)) d(x)
>>> df.divisor()                  # long time
2*Place (1/x, 1/x^4*y^4 + 1/x^3*y^3 + 1/x^2*y^2 + 1/x*y + 1)
 + 2*Place (1/x, 1/x*y + 1)
 - 2*Place (x + 1, y)
 - 2*Place (x^4 + x^3 + x^2 + x + 1, y)
function_field()[source]

Return the function field of the curve.

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve(x^3 - y^2 - x^4 - y^4)
sage: C.function_field()
Function field in y defined by y^4 + y^2 + x^4 - x^3
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(3) - y**Integer(2) - x**Integer(4) - y**Integer(4))
>>> C.function_field()
Function field in y defined by y^4 + y^2 + x^4 - x^3

sage: # needs sage.rings.finite_rings
sage: A.<x,y> = AffineSpace(GF(8), 2)
sage: C = Curve(x^5 + y^5 + x*y + 1)
sage: C.function_field()
Function field in y defined by y^5 + x*y + x^5 + 1
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> A = AffineSpace(GF(Integer(8)), Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(5) + y**Integer(5) + x*y + Integer(1))
>>> C.function_field()
Function field in y defined by y^5 + x*y + x^5 + 1
parametric_representation(place, name=None)[source]

Return a power series representation of the branch of the curve given by place.

INPUT:

  • place – a place on the curve

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve(x^2 + y^2 -1)
sage: p = C(0,1)
sage: p.closed_point()
Point (x, y - 1)
sage: pl = _.place()
sage: C.parametric_representation(pl)
(s + ..., 1 - 1/2*s^2 - 1/8*s^4 - 1/16*s^6 + ...)
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(2) + y**Integer(2) -Integer(1))
>>> p = C(Integer(0),Integer(1))
>>> p.closed_point()
Point (x, y - 1)
>>> pl = _.place()
>>> C.parametric_representation(pl)
(s + ..., 1 - 1/2*s^2 - 1/8*s^4 - 1/16*s^6 + ...)

sage: # needs sage.rings.finite_rings
sage: A.<x,y> = AffineSpace(GF(7^2), 2)
sage: C = Curve(x^2 - x^4 - y^4)
sage: p, = C.singular_closed_points()
sage: b1, b2 = p.places()
sage: xs, ys = C.parametric_representation(b1)
sage: f = xs^2 - xs^4 - ys^4
sage: [f.coefficient(i) for i in range(5)]
[0, 0, 0, 0, 0]
sage: xs, ys = C.parametric_representation(b2)
sage: f = xs^2 - xs^4 - ys^4
sage: [f.coefficient(i) for i in range(5)]
[0, 0, 0, 0, 0]
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> A = AffineSpace(GF(Integer(7)**Integer(2)), Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(2) - x**Integer(4) - y**Integer(4))
>>> p, = C.singular_closed_points()
>>> b1, b2 = p.places()
>>> xs, ys = C.parametric_representation(b1)
>>> f = xs**Integer(2) - xs**Integer(4) - ys**Integer(4)
>>> [f.coefficient(i) for i in range(Integer(5))]
[0, 0, 0, 0, 0]
>>> xs, ys = C.parametric_representation(b2)
>>> f = xs**Integer(2) - xs**Integer(4) - ys**Integer(4)
>>> [f.coefficient(i) for i in range(Integer(5))]
[0, 0, 0, 0, 0]
place_to_closed_point(place)[source]

Return the closed point on the place.

INPUT:

  • place – a place of the function field of the curve

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: A.<x,y> = AffineSpace(GF(4), 2)
sage: C = Curve(x^5 + y^5 + x*y + 1)
sage: F = C.function_field()
sage: pls = F.places(1)
sage: C.place_to_closed_point(pls[-1])
Point (x + 1, y + 1)
sage: C.place_to_closed_point(pls[-2])
Point (x + 1, y + 1)
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> A = AffineSpace(GF(Integer(4)), Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(5) + y**Integer(5) + x*y + Integer(1))
>>> F = C.function_field()
>>> pls = F.places(Integer(1))
>>> C.place_to_closed_point(pls[-Integer(1)])
Point (x + 1, y + 1)
>>> C.place_to_closed_point(pls[-Integer(2)])
Point (x + 1, y + 1)
places_at_infinity()[source]

Return the places of the curve at infinity.

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve(x^3 - y^2 - x^4 - y^4)
sage: C.places_at_infinity()
[Place (1/x, 1/x^2*y, 1/x^3*y^2, 1/x^4*y^3)]
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(3) - y**Integer(2) - x**Integer(4) - y**Integer(4))
>>> C.places_at_infinity()
[Place (1/x, 1/x^2*y, 1/x^3*y^2, 1/x^4*y^3)]

sage: # needs sage.rings.finite_rings
sage: F = GF(9)
sage: A2.<x,y> = AffineSpace(F, 2)
sage: C = A2.curve(y^3 + y - x^4)
sage: C.places_at_infinity()
[Place (1/x, 1/x^3*y^2)]
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> F = GF(Integer(9))
>>> A2 = AffineSpace(F, Integer(2), names=('x', 'y',)); (x, y,) = A2._first_ngens(2)
>>> C = A2.curve(y**Integer(3) + y - x**Integer(4))
>>> C.places_at_infinity()
[Place (1/x, 1/x^3*y^2)]

sage: A.<x,y,z> = AffineSpace(GF(11), 3)
sage: C = Curve([x*z - y^2, y - z^2, x - y*z], A)
sage: C.places_at_infinity()
[Place (1/x, 1/x*z^2)]
>>> from sage.all import *
>>> A = AffineSpace(GF(Integer(11)), Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = A._first_ngens(3)
>>> C = Curve([x*z - y**Integer(2), y - z**Integer(2), x - y*z], A)
>>> C.places_at_infinity()
[Place (1/x, 1/x*z^2)]
places_on(point)[source]

Return the places on the closed point.

INPUT:

  • point – a closed point of the curve

OUTPUT: list of the places of the function field of the curve

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: C = Curve(x^3 - y^2 - x^4 - y^4)
sage: C.singular_closed_points()
[Point (x, y)]
sage: p, = _
sage: C.places_on(p)
[Place (x, y, y^2, 1/x*y^3 + 1/x*y)]
>>> from sage.all import *
>>> A = AffineSpace(QQ, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(3) - y**Integer(2) - x**Integer(4) - y**Integer(4))
>>> C.singular_closed_points()
[Point (x, y)]
>>> p, = _
>>> C.places_on(p)
[Place (x, y, y^2, 1/x*y^3 + 1/x*y)]

sage: # needs sage.rings.finite_rings
sage: k.<a> = GF(9)
sage: A.<x,y> = AffineSpace(k, 2)
sage: C = Curve(y^2 - x^5 - x^4 - 2*x^3 - 2*x - 2)
sage: pts = C.closed_points()
sage: pts
[Point (x, y + (a + 1)),
 Point (x, y + (-a - 1)),
 Point (x + (a + 1), y + (a - 1)),
 Point (x + (a + 1), y + (-a + 1)),
 Point (x - 1, y + (a + 1)),
 Point (x - 1, y + (-a - 1)),
 Point (x + (-a - 1), y + a),
 Point (x + (-a - 1), y + (-a)),
 Point (x + 1, y + 1),
 Point (x + 1, y - 1)]
sage: p1, p2, p3 = pts[:3]
sage: C.places_on(p1)
[Place (x, y + a + 1)]
sage: C.places_on(p2)
[Place (x, y + 2*a + 2)]
sage: C.places_on(p3)
[Place (x + a + 1, y + a + 2)]
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> k = GF(Integer(9), names=('a',)); (a,) = k._first_ngens(1)
>>> A = AffineSpace(k, Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(y**Integer(2) - x**Integer(5) - x**Integer(4) - Integer(2)*x**Integer(3) - Integer(2)*x - Integer(2))
>>> pts = C.closed_points()
>>> pts
[Point (x, y + (a + 1)),
 Point (x, y + (-a - 1)),
 Point (x + (a + 1), y + (a - 1)),
 Point (x + (a + 1), y + (-a + 1)),
 Point (x - 1, y + (a + 1)),
 Point (x - 1, y + (-a - 1)),
 Point (x + (-a - 1), y + a),
 Point (x + (-a - 1), y + (-a)),
 Point (x + 1, y + 1),
 Point (x + 1, y - 1)]
>>> p1, p2, p3 = pts[:Integer(3)]
>>> C.places_on(p1)
[Place (x, y + a + 1)]
>>> C.places_on(p2)
[Place (x, y + 2*a + 2)]
>>> C.places_on(p3)
[Place (x + a + 1, y + a + 2)]

sage: # needs sage.rings.finite_rings
sage: F.<a> = GF(8)
sage: P.<x,y,z> = ProjectiveSpace(F, 2)
sage: Cp = Curve(x^3*y + y^3*z + x*z^3)
sage: C = Cp.affine_patch(0)
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> F = GF(Integer(8), names=('a',)); (a,) = F._first_ngens(1)
>>> P = ProjectiveSpace(F, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> Cp = Curve(x**Integer(3)*y + y**Integer(3)*z + x*z**Integer(3))
>>> C = Cp.affine_patch(Integer(0))
pull_from_function_field(f)[source]

Return the fraction corresponding to f.

INPUT:

  • f – an element of the function field

OUTPUT:

A fraction of polynomials in the coordinate ring of the ambient space of the curve.

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: A.<x,y> = AffineSpace(GF(8), 2)
sage: C = Curve(x^5 + y^5 + x*y + 1)
sage: F = C.function_field()
sage: C.pull_from_function_field(F.gen())
y
sage: C.pull_from_function_field(F.one())
1
sage: C.pull_from_function_field(F.zero())
0
sage: f1 = F.gen()
sage: f2 = F.base_ring().gen()
sage: C.function(C.pull_from_function_field(f1)) == f1
True
sage: C.function(C.pull_from_function_field(f2)) == f2
True
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> A = AffineSpace(GF(Integer(8)), Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(5) + y**Integer(5) + x*y + Integer(1))
>>> F = C.function_field()
>>> C.pull_from_function_field(F.gen())
y
>>> C.pull_from_function_field(F.one())
1
>>> C.pull_from_function_field(F.zero())
0
>>> f1 = F.gen()
>>> f2 = F.base_ring().gen()
>>> C.function(C.pull_from_function_field(f1)) == f1
True
>>> C.function(C.pull_from_function_field(f2)) == f2
True
singular_closed_points()[source]

Return the singular closed points of the curve.

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: A.<x,y> = AffineSpace(GF(7^2), 2)
sage: C = Curve(x^2 - x^4 - y^4)
sage: C.singular_closed_points()
[Point (x, y)]
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> A = AffineSpace(GF(Integer(7)**Integer(2)), Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(2) - x**Integer(4) - y**Integer(4))
>>> C.singular_closed_points()
[Point (x, y)]

sage: A.<x,y,z> = AffineSpace(GF(11), 3)
sage: C = Curve([x*z - y^2, y - z^2, x - y*z], A)
sage: C.singular_closed_points()
[]
>>> from sage.all import *
>>> A = AffineSpace(GF(Integer(11)), Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = A._first_ngens(3)
>>> C = Curve([x*z - y**Integer(2), y - z**Integer(2), x - y*z], A)
>>> C.singular_closed_points()
[]
class sage.schemes.curves.affine_curve.IntegralAffineCurve_finite_field(A, X)[source]

Bases: IntegralAffineCurve

Integral affine curves.

INPUT:

  • A – an ambient space in which the curve lives

  • X – list of polynomials that define the curve

EXAMPLES:

sage: A.<x,y,z> = AffineSpace(GF(11), 3)
sage: C = Curve([x*z - y^2, y - z^2, x - y*z], A); C
Affine Curve over Finite Field of size 11
 defined by -y^2 + x*z, -z^2 + y, -y*z + x
sage: C.function_field()
Function field in z defined by z^3 + 10*x
>>> from sage.all import *
>>> A = AffineSpace(GF(Integer(11)), Integer(3), names=('x', 'y', 'z',)); (x, y, z,) = A._first_ngens(3)
>>> C = Curve([x*z - y**Integer(2), y - z**Integer(2), x - y*z], A); C
Affine Curve over Finite Field of size 11
 defined by -y^2 + x*z, -z^2 + y, -y*z + x
>>> C.function_field()
Function field in z defined by z^3 + 10*x
closed_points(degree=1)[source]

Return a list of the closed points of degree of the curve.

INPUT:

  • degree – positive integer

EXAMPLES:

sage: A.<x,y> = AffineSpace(GF(7), 2)
sage: C = Curve(x^2 - x^4 - y^4)
sage: C.closed_points()
[Point (x, y),
 Point (x + 1, y),
 Point (x + 2, y + 2),
 Point (x + 2, y - 2),
 Point (x - 2, y + 2),
 Point (x - 2, y - 2),
 Point (x - 1, y)]
>>> from sage.all import *
>>> A = AffineSpace(GF(Integer(7)), Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(2) - x**Integer(4) - y**Integer(4))
>>> C.closed_points()
[Point (x, y),
 Point (x + 1, y),
 Point (x + 2, y + 2),
 Point (x + 2, y - 2),
 Point (x - 2, y + 2),
 Point (x - 2, y - 2),
 Point (x - 1, y)]
places(degree=1)[source]

Return all places on the curve of the degree.

INPUT:

  • degree – positive integer

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: F = GF(9)
sage: A2.<x,y> = AffineSpace(F, 2)
sage: C = A2.curve(y^3 + y - x^4)
sage: C.places()
[Place (1/x, 1/x^3*y^2),
 Place (x, y),
 Place (x, y + z2 + 1),
 Place (x, y + 2*z2 + 2),
 Place (x + z2, y + 2),
 Place (x + z2, y + z2),
 Place (x + z2, y + 2*z2 + 1),
 Place (x + z2 + 1, y + 1),
 Place (x + z2 + 1, y + z2 + 2),
 Place (x + z2 + 1, y + 2*z2),
 Place (x + 2*z2 + 1, y + 2),
 Place (x + 2*z2 + 1, y + z2),
 Place (x + 2*z2 + 1, y + 2*z2 + 1),
 Place (x + 2, y + 1),
 Place (x + 2, y + z2 + 2),
 Place (x + 2, y + 2*z2),
 Place (x + 2*z2, y + 2),
 Place (x + 2*z2, y + z2),
 Place (x + 2*z2, y + 2*z2 + 1),
 Place (x + 2*z2 + 2, y + 1),
 Place (x + 2*z2 + 2, y + z2 + 2),
 Place (x + 2*z2 + 2, y + 2*z2),
 Place (x + z2 + 2, y + 2),
 Place (x + z2 + 2, y + z2),
 Place (x + z2 + 2, y + 2*z2 + 1),
 Place (x + 1, y + 1),
 Place (x + 1, y + z2 + 2),
 Place (x + 1, y + 2*z2)]
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> F = GF(Integer(9))
>>> A2 = AffineSpace(F, Integer(2), names=('x', 'y',)); (x, y,) = A2._first_ngens(2)
>>> C = A2.curve(y**Integer(3) + y - x**Integer(4))
>>> C.places()
[Place (1/x, 1/x^3*y^2),
 Place (x, y),
 Place (x, y + z2 + 1),
 Place (x, y + 2*z2 + 2),
 Place (x + z2, y + 2),
 Place (x + z2, y + z2),
 Place (x + z2, y + 2*z2 + 1),
 Place (x + z2 + 1, y + 1),
 Place (x + z2 + 1, y + z2 + 2),
 Place (x + z2 + 1, y + 2*z2),
 Place (x + 2*z2 + 1, y + 2),
 Place (x + 2*z2 + 1, y + z2),
 Place (x + 2*z2 + 1, y + 2*z2 + 1),
 Place (x + 2, y + 1),
 Place (x + 2, y + z2 + 2),
 Place (x + 2, y + 2*z2),
 Place (x + 2*z2, y + 2),
 Place (x + 2*z2, y + z2),
 Place (x + 2*z2, y + 2*z2 + 1),
 Place (x + 2*z2 + 2, y + 1),
 Place (x + 2*z2 + 2, y + z2 + 2),
 Place (x + 2*z2 + 2, y + 2*z2),
 Place (x + z2 + 2, y + 2),
 Place (x + z2 + 2, y + z2),
 Place (x + z2 + 2, y + 2*z2 + 1),
 Place (x + 1, y + 1),
 Place (x + 1, y + z2 + 2),
 Place (x + 1, y + 2*z2)]
class sage.schemes.curves.affine_curve.IntegralAffinePlaneCurve(A, f)[source]

Bases: IntegralAffineCurve, AffinePlaneCurve_field

class sage.schemes.curves.affine_curve.IntegralAffinePlaneCurve_finite_field(A, f)[source]

Bases: AffinePlaneCurve_finite_field, IntegralAffineCurve_finite_field

Integral affine plane curve over a finite field.

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: A.<x,y> = AffineSpace(GF(8), 2)
sage: C = Curve(x^5 + y^5 + x*y + 1); C
Affine Plane Curve over Finite Field in z3 of size 2^3
 defined by x^5 + y^5 + x*y + 1
sage: C.function_field()
Function field in y defined by y^5 + x*y + x^5 + 1
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> A = AffineSpace(GF(Integer(8)), Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(x**Integer(5) + y**Integer(5) + x*y + Integer(1)); C
Affine Plane Curve over Finite Field in z3 of size 2^3
 defined by x^5 + y^5 + x*y + 1
>>> C.function_field()
Function field in y defined by y^5 + x*y + x^5 + 1