Projective curves#

Projective curves in Sage are curves in a projective space or a projective plane.

EXAMPLES:

We can construct curves in either a projective plane:

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

or in higher dimensional projective spaces:

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

Integral projective 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 = GF(2)
sage: P.<x,y,z> = ProjectiveSpace(k, 2)
sage: C = Curve(x^2*z - y^3, P)
sage: C.genus()
0
sage: C.function_field()
Function field in z defined by z + y^3
>>> from sage.all import *
>>> k = GF(Integer(2))
>>> P = ProjectiveSpace(k, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(x**Integer(2)*z - y**Integer(3), P)
>>> C.genus()
0
>>> C.function_field()
Function field in z defined by z + y^3

Closed points of arbitrary degree can be computed:

sage: C.closed_points()
[Point (x, y), Point (y, z), Point (x + z, y + z)]
sage: C.closed_points(2)
[Point (y^2 + y*z + z^2, x + z)]
sage: C.closed_points(3)
[Point (y^3 + y^2*z + z^3, x + y + z),
 Point (x^2 + y*z + z^2, x*y + x*z + y*z, y^2 + x*z + y*z + z^2)]
>>> from sage.all import *
>>> C.closed_points()
[Point (x, y), Point (y, z), Point (x + z, y + z)]
>>> C.closed_points(Integer(2))
[Point (y^2 + y*z + z^2, x + z)]
>>> C.closed_points(Integer(3))
[Point (y^3 + y^2*z + z^3, x + y + z),
 Point (x^2 + y*z + z^2, x*y + x*z + y*z, y^2 + x*z + y*z + z^2)]

All singular closed points can be found:

sage: C.singular_closed_points()
[Point (x, y)]
sage: p = _[0]
sage: p.places()  # a unibranch singularity, that is, a cusp
[Place (1/y)]
sage: pls = _[0]
sage: C.place_to_closed_point(pls)
Point (x, y)
>>> from sage.all import *
>>> C.singular_closed_points()
[Point (x, y)]
>>> p = _[Integer(0)]
>>> p.places()  # a unibranch singularity, that is, a cusp
[Place (1/y)]
>>> pls = _[Integer(0)]
>>> C.place_to_closed_point(pls)
Point (x, y)

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

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

Integral projective curves over \(\QQ\)#

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

EXAMPLES:

sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve(x^2*z^2 - x^4 - y^4, P)
sage: C.singular_closed_points()
[Point (x, y)]
sage: p, = _
sage: p.places()
[Place (1/y, 1/y^2*z - 1), Place (1/y, 1/y^2*z + 1)]
sage: fy = C.function(y/z)
sage: fy.divisor()
Place (1/y, 1/y^2*z - 1)
 + Place (1/y, 1/y^2*z + 1)
 + Place (y, z - 1)
 + Place (y, z + 1)
 - Place (y^4 + 1, z)
sage: supp = _.support()
sage: pl = supp[0]
sage: C.place_to_closed_point(pl)
Point (x, y)
sage: pl = supp[1]
sage: C.place_to_closed_point(pl)
Point (x, y)
sage: _.rational_point()
(0 : 0 : 1)
sage: _ in C
True
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(x**Integer(2)*z**Integer(2) - x**Integer(4) - y**Integer(4), P)
>>> C.singular_closed_points()
[Point (x, y)]
>>> p, = _
>>> p.places()
[Place (1/y, 1/y^2*z - 1), Place (1/y, 1/y^2*z + 1)]
>>> fy = C.function(y/z)
>>> fy.divisor()
Place (1/y, 1/y^2*z - 1)
 + Place (1/y, 1/y^2*z + 1)
 + Place (y, z - 1)
 + Place (y, z + 1)
 - Place (y^4 + 1, z)
>>> supp = _.support()
>>> pl = supp[Integer(0)]
>>> C.place_to_closed_point(pl)
Point (x, y)
>>> pl = supp[Integer(1)]
>>> C.place_to_closed_point(pl)
Point (x, y)
>>> _.rational_point()
(0 : 0 : 1)
>>> _ in C
True

AUTHORS:

  • William Stein (2005-11-13)

  • David Joyner (2005-11-13)

  • David Kohel (2006-01)

  • Moritz Minzlaff (2010-11)

  • Grayson Jorgenson (2016-08)

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

sage.schemes.curves.projective_curve.Hasse_bounds(q, genus=1)[source]#

Return the Hasse-Weil bounds for the cardinality of a nonsingular curve defined over \(\GF{q}\) of given genus.

INPUT:

  • q (int) – a prime power

  • genus (int, default 1) – a non-negative integer,

OUTPUT: A tuple. The Hasse bounds (lb,ub) for the cardinality of a curve of genus genus defined over \(\GF{q}\).

EXAMPLES:

sage: Hasse_bounds(2)
(1, 5)
sage: Hasse_bounds(next_prime(10^30))                                           # needs sage.libs.pari
(999999999999998000000000000058, 1000000000000002000000000000058)
>>> from sage.all import *
>>> Hasse_bounds(Integer(2))
(1, 5)
>>> Hasse_bounds(next_prime(Integer(10)**Integer(30)))                                           # needs sage.libs.pari
(999999999999998000000000000058, 1000000000000002000000000000058)
class sage.schemes.curves.projective_curve.IntegralProjectiveCurve(A, f)[source]#

Bases: ProjectiveCurve_field

Integral projective curve.

coordinate_functions(i=None)[source]#

Return the coordinate functions for the i-th affine patch.

If i is None, return the homogeneous coordinate functions.

EXAMPLES:

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

Return the function field element corresponding to f.

INPUT:

  • f – a fraction of homogeneous polynomials of the coordinate ring of the ambient space of the curve

OUTPUT: An element of the function field.

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: P.<x,y,z> = ProjectiveSpace(GF(4), 2)
sage: C = Curve(x^5 + y^5 + x*y*z^3 + z^5)
sage: f = C.function(x/y); f
1/y
sage: f.divisor()
Place (1/y, 1/y^2*z^2 + z2/y*z + 1)
 + Place (1/y, 1/y^2*z^2 + ((z2 + 1)/y)*z + 1)
 + Place (1/y, 1/y*z + 1)
 - Place (y, z^2 + z2*z + 1)
 - Place (y, z^2 + (z2 + 1)*z + 1)
 - Place (y, z + 1)
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> P = ProjectiveSpace(GF(Integer(4)), Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(x**Integer(5) + y**Integer(5) + x*y*z**Integer(3) + z**Integer(5))
>>> f = C.function(x/y); f
1/y
>>> f.divisor()
Place (1/y, 1/y^2*z^2 + z2/y*z + 1)
 + Place (1/y, 1/y^2*z^2 + ((z2 + 1)/y)*z + 1)
 + Place (1/y, 1/y*z + 1)
 - Place (y, z^2 + z2*z + 1)
 - Place (y, z^2 + (z2 + 1)*z + 1)
 - Place (y, z + 1)
function_field()[source]#

Return the function field of this curve.

EXAMPLES:

sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve(x^2 + y^2 + z^2, P)
sage: C.function_field()
Function field in z defined by z^2 + y^2 + 1
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(x**Integer(2) + y**Integer(2) + z**Integer(2), P)
>>> C.function_field()
Function field in z defined by z^2 + y^2 + 1
sage: # needs sage.rings.finite_rings
sage: P.<x,y,z> = ProjectiveSpace(GF(4), 2)
sage: C = Curve(x^5 + y^5 + x*y*z^3 + z^5)
sage: C.function_field()
Function field in z defined by z^5 + y*z^3 + y^5 + 1
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> P = ProjectiveSpace(GF(Integer(4)), Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(x**Integer(5) + y**Integer(5) + x*y*z**Integer(3) + z**Integer(5))
>>> C.function_field()
Function field in z defined by z^5 + y*z^3 + y^5 + 1
jacobian(model, base_div=None)[source]#

Return the Jacobian of this curve.

INPUT:

  • model – model to use for arithmetic

  • base_div – an effective divisor for the model

The degree of the base divisor should satisfy certain degree condition corresponding to the model used. The following table lists these conditions. Let \(g\) be the geometric genus of the curve.

  • hess: ideal-based arithmetic; requires base divisor of degree \(g\)

  • km_large: Khuri-Makdisi’s large model; requires base divisor of degree at least \(2g + 1\)

  • km_medium: Khuri-Makdisi’s medium model; requires base divisor of degree at least \(2g + 1\)

  • km_small: Khuri-Makdisi’s small model requires base divisor of degree at least \(g + 1\)

We assume the curve (or its function field) has a rational place. If a base divisor is not given, one is chosen using a rational place.

EXAMPLES:

sage: A.<x,y> = AffineSpace(GF(5), 2)
sage: C = Curve(y^2*(x^3 - 1) - (x^3 - 2)).projective_closure()
sage: J = C.jacobian(model='hess'); J
Jacobian of Projective Plane Curve over Finite Field of size 5
 defined by 2*x0^5 - x0^2*x1^3 - x0^3*x2^2 + x1^3*x2^2 (Hess model)
sage: J.base_divisor().degree() == C.genus()
True
>>> from sage.all import *
>>> A = AffineSpace(GF(Integer(5)), Integer(2), names=('x', 'y',)); (x, y,) = A._first_ngens(2)
>>> C = Curve(y**Integer(2)*(x**Integer(3) - Integer(1)) - (x**Integer(3) - Integer(2))).projective_closure()
>>> J = C.jacobian(model='hess'); J
Jacobian of Projective Plane Curve over Finite Field of size 5
 defined by 2*x0^5 - x0^2*x1^3 - x0^3*x2^2 + x1^3*x2^2 (Hess model)
>>> J.base_divisor().degree() == C.genus()
True
place_to_closed_point(place)[source]#

Return the closed point at the place.

INPUT:

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

EXAMPLES:

sage: P.<x,y,z> = ProjectiveSpace(GF(5), 2)
sage: C = Curve(y^2*z^7 - x^9 - x*z^8)
sage: pls = C.places()
sage: C.place_to_closed_point(pls[-1])
Point (x - 2*z, y - 2*z)
sage: pls2 = C.places(2)
sage: C.place_to_closed_point(pls2[0])
Point (y^2 + y*z + z^2, x + y)
>>> from sage.all import *
>>> P = ProjectiveSpace(GF(Integer(5)), Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(y**Integer(2)*z**Integer(7) - x**Integer(9) - x*z**Integer(8))
>>> pls = C.places()
>>> C.place_to_closed_point(pls[-Integer(1)])
Point (x - 2*z, y - 2*z)
>>> pls2 = C.places(Integer(2))
>>> C.place_to_closed_point(pls2[Integer(0)])
Point (y^2 + y*z + z^2, x + y)
places_on(point)[source]#

Return the places on the closed point.

INPUT:

  • point – a closed point of the curve

EXAMPLES:

sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve(x*y*z^4 - x^6 - y^6)
sage: C.singular_closed_points()
[Point (x, y)]
sage: p, = _
sage: C.places_on(p)
[Place (1/y, 1/y^2*z, 1/y^3*z^2, 1/y^4*z^3),
 Place (y, y*z, y*z^2, y*z^3)]
sage: pl1, pl2 =_
sage: C.place_to_closed_point(pl1)
Point (x, y)
sage: C.place_to_closed_point(pl2)
Point (x, y)
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(x*y*z**Integer(4) - x**Integer(6) - y**Integer(6))
>>> C.singular_closed_points()
[Point (x, y)]
>>> p, = _
>>> C.places_on(p)
[Place (1/y, 1/y^2*z, 1/y^3*z^2, 1/y^4*z^3),
 Place (y, y*z, y*z^2, y*z^3)]
>>> pl1, pl2 =_
>>> C.place_to_closed_point(pl1)
Point (x, y)
>>> C.place_to_closed_point(pl2)
Point (x, y)
sage: P.<x,y,z> = ProjectiveSpace(GF(5), 2)
sage: C = Curve(x^2*z - y^3)
sage: [C.places_on(p) for p in C.closed_points()]
[[Place (1/y)],
 [Place (y)],
 [Place (y + 1)],
 [Place (y + 2)],
 [Place (y + 3)],
 [Place (y + 4)]]
>>> from sage.all import *
>>> P = ProjectiveSpace(GF(Integer(5)), Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(x**Integer(2)*z - y**Integer(3))
>>> [C.places_on(p) for p in C.closed_points()]
[[Place (1/y)],
 [Place (y)],
 [Place (y + 1)],
 [Place (y + 2)],
 [Place (y + 3)],
 [Place (y + 4)]]
pull_from_function_field(f)[source]#

Return the fraction corresponding to f.

INPUT:

  • f – an element of the function field

OUTPUT:

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

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: P.<x,y,z> = ProjectiveSpace(GF(4), 2)
sage: C = Curve(x^5 + y^5 + x*y*z^3 + z^5)
sage: F = C.function_field()
sage: C.pull_from_function_field(F.gen())
z/x
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
>>> P = ProjectiveSpace(GF(Integer(4)), Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(x**Integer(5) + y**Integer(5) + x*y*z**Integer(3) + z**Integer(5))
>>> F = C.function_field()
>>> C.pull_from_function_field(F.gen())
z/x
>>> 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: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve(y^2*z - x^3, P)
sage: C.singular_closed_points()
[Point (x, y)]
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(y**Integer(2)*z - x**Integer(3), P)
>>> C.singular_closed_points()
[Point (x, y)]
sage: P.<x,y,z> = ProjectiveSpace(GF(5), 2)
sage: C = Curve(y^2*z^7 - x^9 - x*z^8)
sage: C.singular_closed_points()
[Point (x, z)]
>>> from sage.all import *
>>> P = ProjectiveSpace(GF(Integer(5)), Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(y**Integer(2)*z**Integer(7) - x**Integer(9) - x*z**Integer(8))
>>> C.singular_closed_points()
[Point (x, z)]
class sage.schemes.curves.projective_curve.IntegralProjectiveCurve_finite_field(A, f)[source]#

Bases: IntegralProjectiveCurve

Integral projective curve over a finite field.

INPUT:

  • A – an ambient projective space

  • f – homogeneous polynomials defining the curve

EXAMPLES:

sage: P.<x,y,z> = ProjectiveSpace(GF(5), 2)
sage: C = Curve(y^2*z^7 - x^9 - x*z^8)
sage: C.function_field()
Function field in z defined by z^8 + 4*y^2*z^7 + 1
sage: C.closed_points()
[Point (x, z),
 Point (x, y),
 Point (x - 2*z, y + 2*z),
 Point (x + 2*z, y + z),
 Point (x + 2*z, y - z),
 Point (x - 2*z, y - 2*z)]
>>> from sage.all import *
>>> P = ProjectiveSpace(GF(Integer(5)), Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(y**Integer(2)*z**Integer(7) - x**Integer(9) - x*z**Integer(8))
>>> C.function_field()
Function field in z defined by z^8 + 4*y^2*z^7 + 1
>>> C.closed_points()
[Point (x, z),
 Point (x, y),
 Point (x - 2*z, y + 2*z),
 Point (x + 2*z, y + z),
 Point (x + 2*z, y - z),
 Point (x - 2*z, y - 2*z)]
L_polynomial(name='t')[source]#

Return the L-polynomial of this possibly singular curve.

INPUT:

  • name – (default: t) name of the variable of the polynomial

EXAMPLES:

sage: A.<x,y> = AffineSpace(GF(3), 2)
sage: C = Curve(y^2 - x^5 - x^4 - 2*x^3 - 2*x - 2)
sage: Cbar = C.projective_closure()
sage: Cbar.L_polynomial()
9*t^4 - 3*t^3 + t^2 - t + 1
>>> from sage.all import *
>>> A = AffineSpace(GF(Integer(3)), 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))
>>> Cbar = C.projective_closure()
>>> Cbar.L_polynomial()
9*t^4 - 3*t^3 + t^2 - t + 1
closed_points(degree=1)[source]#

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

INPUT:

  • degree – a positive integer

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: A.<x,y> = AffineSpace(GF(9),2)
sage: C = Curve(y^2 - x^5 - x^4 - 2*x^3 - 2*x-2)
sage: Cp = C.projective_closure()
sage: Cp.closed_points()
[Point (x0, x1),
 Point (x0 + (-z2 - 1)*x2, x1),
 Point (x0 + (z2 + 1)*x2, x1),
 Point (x0 + z2*x2, x1 + (z2 - 1)*x2),
 Point (x0 + (-z2)*x2, x1 + (-z2 + 1)*x2),
 Point (x0 + (-z2 - 1)*x2, x1 + (-z2 - 1)*x2),
 Point (x0 + (z2 + 1)*x2, x1 + (z2 + 1)*x2),
 Point (x0 + (z2 - 1)*x2, x1 + z2*x2),
 Point (x0 + (-z2 + 1)*x2, x1 + (-z2)*x2),
 Point (x0 + x2, x1 - x2),
 Point (x0 - x2, x1 + x2)]
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> A = AffineSpace(GF(Integer(9)),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))
>>> Cp = C.projective_closure()
>>> Cp.closed_points()
[Point (x0, x1),
 Point (x0 + (-z2 - 1)*x2, x1),
 Point (x0 + (z2 + 1)*x2, x1),
 Point (x0 + z2*x2, x1 + (z2 - 1)*x2),
 Point (x0 + (-z2)*x2, x1 + (-z2 + 1)*x2),
 Point (x0 + (-z2 - 1)*x2, x1 + (-z2 - 1)*x2),
 Point (x0 + (z2 + 1)*x2, x1 + (z2 + 1)*x2),
 Point (x0 + (z2 - 1)*x2, x1 + z2*x2),
 Point (x0 + (-z2 + 1)*x2, x1 + (-z2)*x2),
 Point (x0 + x2, x1 - x2),
 Point (x0 - x2, x1 + x2)]
number_of_rational_points(r=1)[source]#

Return the number of rational points of the curve with constant field extended by degree r.

INPUT:

  • r – positive integer (default: \(1\))

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: A.<x,y> = AffineSpace(GF(3), 2)
sage: C = Curve(y^2 - x^5 - x^4 - 2*x^3 - 2*x - 2)
sage: Cbar = C.projective_closure()
sage: Cbar.number_of_rational_points(3)
21
sage: D = Cbar.change_ring(Cbar.base_ring().extension(3))
sage: D.base_ring()
Finite Field in z3 of size 3^3
sage: len(D.closed_points())
21
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> A = AffineSpace(GF(Integer(3)), 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))
>>> Cbar = C.projective_closure()
>>> Cbar.number_of_rational_points(Integer(3))
21
>>> D = Cbar.change_ring(Cbar.base_ring().extension(Integer(3)))
>>> D.base_ring()
Finite Field in z3 of size 3^3
>>> len(D.closed_points())
21
places(degree=1)[source]#

Return all places on the curve of the degree.

INPUT:

  • degree – positive integer

EXAMPLES:

sage: P.<x,y,z> = ProjectiveSpace(GF(5), 2)
sage: C = Curve(x^2*z - y^3)
sage: C.places()
[Place (1/y),
 Place (y),
 Place (y + 1),
 Place (y + 2),
 Place (y + 3),
 Place (y + 4)]
sage: C.places(2)
[Place (y^2 + 2),
 Place (y^2 + 3),
 Place (y^2 + y + 1),
 Place (y^2 + y + 2),
 Place (y^2 + 2*y + 3),
 Place (y^2 + 2*y + 4),
 Place (y^2 + 3*y + 3),
 Place (y^2 + 3*y + 4),
 Place (y^2 + 4*y + 1),
 Place (y^2 + 4*y + 2)]
>>> from sage.all import *
>>> P = ProjectiveSpace(GF(Integer(5)), Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(x**Integer(2)*z - y**Integer(3))
>>> C.places()
[Place (1/y),
 Place (y),
 Place (y + 1),
 Place (y + 2),
 Place (y + 3),
 Place (y + 4)]
>>> C.places(Integer(2))
[Place (y^2 + 2),
 Place (y^2 + 3),
 Place (y^2 + y + 1),
 Place (y^2 + y + 2),
 Place (y^2 + 2*y + 3),
 Place (y^2 + 2*y + 4),
 Place (y^2 + 3*y + 3),
 Place (y^2 + 3*y + 4),
 Place (y^2 + 4*y + 1),
 Place (y^2 + 4*y + 2)]
class sage.schemes.curves.projective_curve.IntegralProjectivePlaneCurve(A, f)[source]#

Bases: IntegralProjectiveCurve, ProjectivePlaneCurve_field

class sage.schemes.curves.projective_curve.IntegralProjectivePlaneCurve_finite_field(A, f)[source]#

Bases: IntegralProjectiveCurve_finite_field, ProjectivePlaneCurve_finite_field

Integral projective plane curve over a finite field.

INPUT:

  • A – ambient projective plane

  • f – a homogeneous equation that defines the curve

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: A.<x,y> = AffineSpace(GF(9), 2)
sage: C = Curve(y^2 - x^5 - x^4 - 2*x^3 - 2*x - 2)
sage: Cb = C.projective_closure()
sage: Cb.singular_closed_points()
[Point (x0, x1)]
sage: Cb.function_field()
Function field in y defined by y^2 + 2*x^5 + 2*x^4 + x^3 + x + 1
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> A = AffineSpace(GF(Integer(9)), 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))
>>> Cb = C.projective_closure()
>>> Cb.singular_closed_points()
[Point (x0, x1)]
>>> Cb.function_field()
Function field in y defined by y^2 + 2*x^5 + 2*x^4 + x^3 + x + 1
class sage.schemes.curves.projective_curve.ProjectiveCurve(A, X, category=None)[source]#

Bases: Curve_generic, AlgebraicScheme_subscheme_projective

Curves in projective spaces.

INPUT:

  • A – ambient projective space

  • X – list of multivariate polynomials; defining equations of the curve

EXAMPLES:

sage: P.<x,y,z,w,u> = ProjectiveSpace(GF(7), 4)
sage: C = Curve([y*u^2 - x^3, z*u^2 - x^3, w*u^2 - x^3, y^3 - x^3], P); C
Projective Curve over Finite Field of size 7 defined
 by -x^3 + y*u^2, -x^3 + z*u^2, -x^3 + w*u^2, -x^3 + y^3
>>> from sage.all import *
>>> P = ProjectiveSpace(GF(Integer(7)), Integer(4), names=('x', 'y', 'z', 'w', 'u',)); (x, y, z, w, u,) = P._first_ngens(5)
>>> C = Curve([y*u**Integer(2) - x**Integer(3), z*u**Integer(2) - x**Integer(3), w*u**Integer(2) - x**Integer(3), y**Integer(3) - x**Integer(3)], P); C
Projective Curve over Finite Field of size 7 defined
 by -x^3 + y*u^2, -x^3 + z*u^2, -x^3 + w*u^2, -x^3 + y^3
sage: # needs sage.rings.number_field
sage: K.<u> = CyclotomicField(11)
sage: P.<x,y,z,w> = ProjectiveSpace(K, 3)
sage: C = Curve([y*w - u*z^2 - x^2, x*w - 3*u^2*z*w], P); C
Projective Curve over Cyclotomic Field of order 11 and degree 10 defined
 by -x^2 + (-u)*z^2 + y*w, x*w + (-3*u^2)*z*w
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> K = CyclotomicField(Integer(11), names=('u',)); (u,) = K._first_ngens(1)
>>> P = ProjectiveSpace(K, Integer(3), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = P._first_ngens(4)
>>> C = Curve([y*w - u*z**Integer(2) - x**Integer(2), x*w - Integer(3)*u**Integer(2)*z*w], P); C
Projective Curve over Cyclotomic Field of order 11 and degree 10 defined
 by -x^2 + (-u)*z^2 + y*w, x*w + (-3*u^2)*z*w
affine_patch(i, AA=None)[source]#

Return the \(i\)-th affine patch of this projective curve.

INPUT:

  • i – affine coordinate chart of the projective ambient space of this curve to compute affine patch with respect to

  • AA – (default: None) ambient affine space, this is constructed if it is not given

OUTPUT: A curve in affine space.

EXAMPLES:

sage: P.<x,y,z,w> = ProjectiveSpace(CC, 3)
sage: C = Curve([y*z - x^2, w^2 - x*y], P)
sage: C.affine_patch(0)
Affine Curve over Complex Field with 53 bits of precision defined
 by y*z - 1.00000000000000, w^2 - y
>>> from sage.all import *
>>> P = ProjectiveSpace(CC, Integer(3), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = P._first_ngens(4)
>>> C = Curve([y*z - x**Integer(2), w**Integer(2) - x*y], P)
>>> C.affine_patch(Integer(0))
Affine Curve over Complex Field with 53 bits of precision defined
 by y*z - 1.00000000000000, w^2 - y
sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve(x^3 - x^2*y + y^3 - x^2*z, P)
sage: C.affine_patch(1)
Affine Plane Curve over Rational Field defined by x^3 - x^2*z - x^2 + 1
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(x**Integer(3) - x**Integer(2)*y + y**Integer(3) - x**Integer(2)*z, P)
>>> C.affine_patch(Integer(1))
Affine Plane Curve over Rational Field defined by x^3 - x^2*z - x^2 + 1
sage: A.<x,y> = AffineSpace(QQ, 2)
sage: P.<u,v,w> = ProjectiveSpace(QQ, 2)
sage: C = Curve([u^2 - v^2], P)
sage: C.affine_patch(1, A).ambient_space() == A
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([u**Integer(2) - v**Integer(2)], P)
>>> C.affine_patch(Integer(1), A).ambient_space() == A
True
plane_projection(PP=None)[source]#

Return a projection of this curve into a projective plane.

INPUT:

  • PP – (default: None) the projective plane the projected curve will be defined in. This space must be defined over the same base field as this curve, and must have dimension two. This space is constructed if not specified.

OUTPUT: A tuple of

  • a scheme morphism from this curve into a projective plane

  • the projective curve that is the image of that morphism

EXAMPLES:

sage: P.<x,y,z,w,u,v> = ProjectiveSpace(QQ, 5)
sage: C = P.curve([x*u - z*v, w - y, w*y - x^2, y^3*u*2*z - w^4*w])
sage: L.<a,b,c> = ProjectiveSpace(QQ, 2)
sage: proj1 = C.plane_projection(PP=L)
sage: proj1
(Scheme morphism:
   From: Projective Curve over Rational Field
         defined by x*u - z*v, -y + w, -x^2 + y*w, -w^5 + 2*y^3*z*u
   To:   Projective Space of dimension 2 over Rational Field
   Defn: Defined on coordinates by sending (x : y : z : w : u : v) to
         (x : -z + u : -z + v),
 Projective Plane Curve over Rational Field defined by a^8 + 6*a^7*b +
  4*a^5*b^3 - 4*a^7*c - 2*a^6*b*c - 4*a^5*b^2*c + 2*a^6*c^2)
sage: proj1[1].ambient_space() is L
True
sage: proj2 = C.projection()
sage: proj2[1].ambient_space() is L
False
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(5), names=('x', 'y', 'z', 'w', 'u', 'v',)); (x, y, z, w, u, v,) = P._first_ngens(6)
>>> C = P.curve([x*u - z*v, w - y, w*y - x**Integer(2), y**Integer(3)*u*Integer(2)*z - w**Integer(4)*w])
>>> L = ProjectiveSpace(QQ, Integer(2), names=('a', 'b', 'c',)); (a, b, c,) = L._first_ngens(3)
>>> proj1 = C.plane_projection(PP=L)
>>> proj1
(Scheme morphism:
   From: Projective Curve over Rational Field
         defined by x*u - z*v, -y + w, -x^2 + y*w, -w^5 + 2*y^3*z*u
   To:   Projective Space of dimension 2 over Rational Field
   Defn: Defined on coordinates by sending (x : y : z : w : u : v) to
         (x : -z + u : -z + v),
 Projective Plane Curve over Rational Field defined by a^8 + 6*a^7*b +
  4*a^5*b^3 - 4*a^7*c - 2*a^6*b*c - 4*a^5*b^2*c + 2*a^6*c^2)
>>> proj1[Integer(1)].ambient_space() is L
True
>>> proj2 = C.projection()
>>> proj2[Integer(1)].ambient_space() is L
False
sage: P.<x,y,z,w,u> = ProjectiveSpace(GF(7), 4)
sage: C = P.curve([x^2 - 6*y^2, w*z*u - y^3 + 4*y^2*z, u^2 - x^2])
sage: C.plane_projection()
(Scheme morphism:
   From: Projective Curve over Finite Field of size 7
         defined by x^2 + y^2, -y^3 - 3*y^2*z + z*w*u, -x^2 + u^2
   To:   Projective Space of dimension 2 over Finite Field of size 7
   Defn: Defined on coordinates by sending (x : y : z : w : u) to
         (x : z : -y + w),
 Projective Plane Curve over Finite Field of size 7
  defined by x0^10 + 2*x0^8*x1^2 + 2*x0^6*x1^4 - 3*x0^6*x1^3*x2
             + 2*x0^6*x1^2*x2^2 - 2*x0^4*x1^4*x2^2 + x0^2*x1^4*x2^4)
>>> from sage.all import *
>>> P = ProjectiveSpace(GF(Integer(7)), Integer(4), names=('x', 'y', 'z', 'w', 'u',)); (x, y, z, w, u,) = P._first_ngens(5)
>>> C = P.curve([x**Integer(2) - Integer(6)*y**Integer(2), w*z*u - y**Integer(3) + Integer(4)*y**Integer(2)*z, u**Integer(2) - x**Integer(2)])
>>> C.plane_projection()
(Scheme morphism:
   From: Projective Curve over Finite Field of size 7
         defined by x^2 + y^2, -y^3 - 3*y^2*z + z*w*u, -x^2 + u^2
   To:   Projective Space of dimension 2 over Finite Field of size 7
   Defn: Defined on coordinates by sending (x : y : z : w : u) to
         (x : z : -y + w),
 Projective Plane Curve over Finite Field of size 7
  defined by x0^10 + 2*x0^8*x1^2 + 2*x0^6*x1^4 - 3*x0^6*x1^3*x2
             + 2*x0^6*x1^2*x2^2 - 2*x0^4*x1^4*x2^2 + x0^2*x1^4*x2^4)
sage: P.<x,y,z> = ProjectiveSpace(GF(17), 2)
sage: C = P.curve(x^2 - y*z - z^2)
sage: C.plane_projection()
Traceback (most recent call last):
...
TypeError: this curve is already a plane curve
>>> from sage.all import *
>>> P = ProjectiveSpace(GF(Integer(17)), Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = P.curve(x**Integer(2) - y*z - z**Integer(2))
>>> C.plane_projection()
Traceback (most recent call last):
...
TypeError: this curve is already a plane curve
projection(P=None, PS=None)[source]#

Return a projection of this curve into projective space of dimension one less than the dimension of the ambient space of this curve.

This curve must not already be a plane curve. Over finite fields, if this curve contains all points in its ambient space, then an error will be returned.

INPUT:

  • P – (default: None) a point not on this curve that will be used to define the projection map; this is constructed if not specified.

  • PS – (default: None) the projective space the projected curve will be defined in. This space must be defined over the same base ring as this curve, and must have dimension one less than that of the ambient space of this curve. This space will be constructed if not specified.

OUTPUT: A tuple of

  • a scheme morphism from this curve into a projective space of dimension one less than that of the ambient space of this curve

  • the projective curve that is the image of that morphism

EXAMPLES:

sage: # needs sage.rings.number_field
sage: K.<a> = CyclotomicField(3)
sage: P.<x,y,z,w> = ProjectiveSpace(K, 3)
sage: C = Curve([y*w - x^2, z*w^2 - a*x^3], P)
sage: L.<a,b,c> = ProjectiveSpace(K, 2)
sage: proj1 = C.projection(PS=L)
sage: proj1
(Scheme morphism:
   From: Projective Curve over Cyclotomic Field of order 3 and degree 2
         defined by -x^2 + y*w, (-a)*x^3 + z*w^2
   To:   Projective Space of dimension 2
         over Cyclotomic Field of order 3 and degree 2
   Defn: Defined on coordinates by sending (x : y : z : w) to
         (x : y : -z + w),
 Projective Plane Curve over Cyclotomic Field of order 3 and degree 2
  defined by a^6 + (-a)*a^3*b^3 - a^4*b*c)
sage: proj1[1].ambient_space() is L
True
sage: proj2 = C.projection()
sage: proj2[1].ambient_space() is L
False
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> K = CyclotomicField(Integer(3), names=('a',)); (a,) = K._first_ngens(1)
>>> P = ProjectiveSpace(K, Integer(3), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = P._first_ngens(4)
>>> C = Curve([y*w - x**Integer(2), z*w**Integer(2) - a*x**Integer(3)], P)
>>> L = ProjectiveSpace(K, Integer(2), names=('a', 'b', 'c',)); (a, b, c,) = L._first_ngens(3)
>>> proj1 = C.projection(PS=L)
>>> proj1
(Scheme morphism:
   From: Projective Curve over Cyclotomic Field of order 3 and degree 2
         defined by -x^2 + y*w, (-a)*x^3 + z*w^2
   To:   Projective Space of dimension 2
         over Cyclotomic Field of order 3 and degree 2
   Defn: Defined on coordinates by sending (x : y : z : w) to
         (x : y : -z + w),
 Projective Plane Curve over Cyclotomic Field of order 3 and degree 2
  defined by a^6 + (-a)*a^3*b^3 - a^4*b*c)
>>> proj1[Integer(1)].ambient_space() is L
True
>>> proj2 = C.projection()
>>> proj2[Integer(1)].ambient_space() is L
False
sage: P.<x,y,z,w,a,b,c> = ProjectiveSpace(QQ, 6)
sage: C = Curve([y - x, z - a - b, w^2 - c^2, z - x - a, x^2 - w*z], P)
sage: C.projection()
(Scheme morphism:
   From: Projective Curve over Rational Field
         defined by -x + y, z - a - b, w^2 - c^2, -x + z - a, x^2 - z*w
   To:   Projective Space of dimension 5 over Rational Field
   Defn: Defined on coordinates by sending (x : y : z : w : a : b : c)
         to (x : y : -z + w : a : b : c),
 Projective Curve over Rational Field defined by x1 - x4, x0 - x4, x2*x3
  + x3^2 + x2*x4 + 2*x3*x4, x2^2 - x3^2 - 2*x3*x4 + x4^2 - x5^2, x2*x4^2 +
  x3*x4^2 + x4^3 - x3*x5^2 - x4*x5^2, x4^4 - x3^2*x5^2 - 2*x3*x4*x5^2 -
  x4^2*x5^2)
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(6), names=('x', 'y', 'z', 'w', 'a', 'b', 'c',)); (x, y, z, w, a, b, c,) = P._first_ngens(7)
>>> C = Curve([y - x, z - a - b, w**Integer(2) - c**Integer(2), z - x - a, x**Integer(2) - w*z], P)
>>> C.projection()
(Scheme morphism:
   From: Projective Curve over Rational Field
         defined by -x + y, z - a - b, w^2 - c^2, -x + z - a, x^2 - z*w
   To:   Projective Space of dimension 5 over Rational Field
   Defn: Defined on coordinates by sending (x : y : z : w : a : b : c)
         to (x : y : -z + w : a : b : c),
 Projective Curve over Rational Field defined by x1 - x4, x0 - x4, x2*x3
  + x3^2 + x2*x4 + 2*x3*x4, x2^2 - x3^2 - 2*x3*x4 + x4^2 - x5^2, x2*x4^2 +
  x3*x4^2 + x4^3 - x3*x5^2 - x4*x5^2, x4^4 - x3^2*x5^2 - 2*x3*x4*x5^2 -
  x4^2*x5^2)
sage: P.<x,y,z,w> = ProjectiveSpace(GF(2), 3)
sage: C = P.curve([(x - y)*(x - z)*(x - w)*(y - z)*(y - w),
....:              x*y*z*w*(x + y + z + w)])
sage: C.projection()
Traceback (most recent call last):
...
NotImplementedError: this curve contains all points of its ambient space
>>> from sage.all import *
>>> P = ProjectiveSpace(GF(Integer(2)), Integer(3), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = P._first_ngens(4)
>>> C = P.curve([(x - y)*(x - z)*(x - w)*(y - z)*(y - w),
...              x*y*z*w*(x + y + z + w)])
>>> C.projection()
Traceback (most recent call last):
...
NotImplementedError: this curve contains all points of its ambient space
sage: P.<x,y,z,w,u> = ProjectiveSpace(GF(7), 4)
sage: C = P.curve([x^3 - y*z*u, w^2 - u^2 + 2*x*z, 3*x*w - y^2])
sage: L.<a,b,c,d> = ProjectiveSpace(GF(7), 3)
sage: C.projection(PS=L)
(Scheme morphism:
   From: Projective Curve over Finite Field of size 7
         defined by x^3 - y*z*u, 2*x*z + w^2 - u^2, -y^2 + 3*x*w
   To:   Projective Space of dimension 3 over Finite Field of size 7
   Defn: Defined on coordinates by sending (x : y : z : w : u) to
         (x : y : z : w),
 Projective Curve over Finite Field of size 7 defined by b^2 - 3*a*d,
  a^5*b + a*b*c^3*d - 3*b*c^2*d^3, a^6 + a^2*c^3*d - 3*a*c^2*d^3)
sage: Q.<a,b,c> = ProjectiveSpace(GF(7), 2)
sage: C.projection(PS=Q)
Traceback (most recent call last):
...
TypeError: (=Projective Space of dimension 2 over Finite Field of
size 7) must have dimension (=3)
>>> from sage.all import *
>>> P = ProjectiveSpace(GF(Integer(7)), Integer(4), names=('x', 'y', 'z', 'w', 'u',)); (x, y, z, w, u,) = P._first_ngens(5)
>>> C = P.curve([x**Integer(3) - y*z*u, w**Integer(2) - u**Integer(2) + Integer(2)*x*z, Integer(3)*x*w - y**Integer(2)])
>>> L = ProjectiveSpace(GF(Integer(7)), Integer(3), names=('a', 'b', 'c', 'd',)); (a, b, c, d,) = L._first_ngens(4)
>>> C.projection(PS=L)
(Scheme morphism:
   From: Projective Curve over Finite Field of size 7
         defined by x^3 - y*z*u, 2*x*z + w^2 - u^2, -y^2 + 3*x*w
   To:   Projective Space of dimension 3 over Finite Field of size 7
   Defn: Defined on coordinates by sending (x : y : z : w : u) to
         (x : y : z : w),
 Projective Curve over Finite Field of size 7 defined by b^2 - 3*a*d,
  a^5*b + a*b*c^3*d - 3*b*c^2*d^3, a^6 + a^2*c^3*d - 3*a*c^2*d^3)
>>> Q = ProjectiveSpace(GF(Integer(7)), Integer(2), names=('a', 'b', 'c',)); (a, b, c,) = Q._first_ngens(3)
>>> C.projection(PS=Q)
Traceback (most recent call last):
...
TypeError: (=Projective Space of dimension 2 over Finite Field of
size 7) must have dimension (=3)
sage: PP.<x,y,z,w> = ProjectiveSpace(QQ, 3)
sage: C = PP.curve([x^3 - z^2*y, w^2 - z*x])
sage: Q = PP([1,0,1,1])
sage: C.projection(P=Q)
(Scheme morphism:
   From: Projective Curve over Rational Field
         defined by x^3 - y*z^2, -x*z + w^2
   To:   Projective Space of dimension 2 over Rational Field
   Defn: Defined on coordinates by sending (x : y : z : w) to
         (y : -x + z : -x + w),
 Projective Plane Curve over Rational Field defined by x0*x1^5 -
  6*x0*x1^4*x2 + 14*x0*x1^3*x2^2 - 16*x0*x1^2*x2^3 + 9*x0*x1*x2^4 -
  2*x0*x2^5 - x2^6)
sage: LL.<a,b,c> = ProjectiveSpace(QQ, 2)
sage: Q = PP([0,0,0,1])
sage: C.projection(PS=LL, P=Q)
(Scheme morphism:
   From: Projective Curve over Rational Field
         defined by x^3 - y*z^2, -x*z + w^2
   To:   Projective Space of dimension 2 over Rational Field
   Defn: Defined on coordinates by sending (x : y : z : w) to
         (x : y : z),
 Projective Plane Curve over Rational Field defined by a^3 - b*c^2)
sage: Q = PP([0,0,1,0])
sage: C.projection(P=Q)
Traceback (most recent call last):
...
TypeError: (=(0 : 0 : 1 : 0)) must be a point not on this curve
>>> from sage.all import *
>>> PP = ProjectiveSpace(QQ, Integer(3), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = PP._first_ngens(4)
>>> C = PP.curve([x**Integer(3) - z**Integer(2)*y, w**Integer(2) - z*x])
>>> Q = PP([Integer(1),Integer(0),Integer(1),Integer(1)])
>>> C.projection(P=Q)
(Scheme morphism:
   From: Projective Curve over Rational Field
         defined by x^3 - y*z^2, -x*z + w^2
   To:   Projective Space of dimension 2 over Rational Field
   Defn: Defined on coordinates by sending (x : y : z : w) to
         (y : -x + z : -x + w),
 Projective Plane Curve over Rational Field defined by x0*x1^5 -
  6*x0*x1^4*x2 + 14*x0*x1^3*x2^2 - 16*x0*x1^2*x2^3 + 9*x0*x1*x2^4 -
  2*x0*x2^5 - x2^6)
>>> LL = ProjectiveSpace(QQ, Integer(2), names=('a', 'b', 'c',)); (a, b, c,) = LL._first_ngens(3)
>>> Q = PP([Integer(0),Integer(0),Integer(0),Integer(1)])
>>> C.projection(PS=LL, P=Q)
(Scheme morphism:
   From: Projective Curve over Rational Field
         defined by x^3 - y*z^2, -x*z + w^2
   To:   Projective Space of dimension 2 over Rational Field
   Defn: Defined on coordinates by sending (x : y : z : w) to
         (x : y : z),
 Projective Plane Curve over Rational Field defined by a^3 - b*c^2)
>>> Q = PP([Integer(0),Integer(0),Integer(1),Integer(0)])
>>> C.projection(P=Q)
Traceback (most recent call last):
...
TypeError: (=(0 : 0 : 1 : 0)) must be a point not on this curve
sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = P.curve(y^2 - x^2 + z^2)
sage: C.projection()
Traceback (most recent call last):
...
TypeError: this curve is already a plane curve
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = P.curve(y**Integer(2) - x**Integer(2) + z**Integer(2))
>>> C.projection()
Traceback (most recent call last):
...
TypeError: this curve is already a plane curve
class sage.schemes.curves.projective_curve.ProjectiveCurve_field(A, X, category=None)[source]#

Bases: ProjectiveCurve, AlgebraicScheme_subscheme_projective_field

Projective curves over fields.

arithmetic_genus()[source]#

Return the arithmetic genus of this projective curve.

This is the arithmetic genus \(p_a(C)\) as defined in [Har1977]. If \(P\) is the Hilbert polynomial of the defining ideal of this curve, then the arithmetic genus of this curve is \(1 - P(0)\).

EXAMPLES:

sage: P.<x,y,z,w> = ProjectiveSpace(QQ, 3)
sage: C = P.curve([w*z - x^2, w^2 + y^2 + z^2])
sage: C.arithmetic_genus()
1
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(3), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = P._first_ngens(4)
>>> C = P.curve([w*z - x**Integer(2), w**Integer(2) + y**Integer(2) + z**Integer(2)])
>>> C.arithmetic_genus()
1
sage: P.<x,y,z,w,t> = ProjectiveSpace(GF(7), 4)
sage: C = P.curve([t^3 - x*y*w, x^3 + y^3 + z^3, z - w])
sage: C.arithmetic_genus()
10
>>> from sage.all import *
>>> P = ProjectiveSpace(GF(Integer(7)), Integer(4), names=('x', 'y', 'z', 'w', 't',)); (x, y, z, w, t,) = P._first_ngens(5)
>>> C = P.curve([t**Integer(3) - x*y*w, x**Integer(3) + y**Integer(3) + z**Integer(3), z - w])
>>> C.arithmetic_genus()
10
is_complete_intersection()[source]#

Return whether this projective curve is a complete intersection.

EXAMPLES:

sage: P.<x,y,z,w> = ProjectiveSpace(QQ, 3)
sage: C = Curve([x*y - z*w, x^2 - y*w, y^2*w - x*z*w], P)
sage: C.is_complete_intersection()
False
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(3), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = P._first_ngens(4)
>>> C = Curve([x*y - z*w, x**Integer(2) - y*w, y**Integer(2)*w - x*z*w], P)
>>> C.is_complete_intersection()
False
sage: P.<x,y,z,w> = ProjectiveSpace(QQ, 3)
sage: C = Curve([y*w - x^2, z*w^2 - x^3], P)
sage: C.is_complete_intersection()
True

sage: P.<x,y,z,w> = ProjectiveSpace(QQ, 3)
sage: C = Curve([z^2 - y*w, y*z - x*w, y^2 - x*z], P)
sage: C.is_complete_intersection()
False
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(3), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = P._first_ngens(4)
>>> C = Curve([y*w - x**Integer(2), z*w**Integer(2) - x**Integer(3)], P)
>>> C.is_complete_intersection()
True

>>> P = ProjectiveSpace(QQ, Integer(3), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = P._first_ngens(4)
>>> C = Curve([z**Integer(2) - y*w, y*z - x*w, y**Integer(2) - x*z], P)
>>> C.is_complete_intersection()
False
tangent_line(p)[source]#

Return the tangent line at the point p.

INPUT:

  • p – a rational point of the curve

EXAMPLES:

sage: P.<x,y,z,w> = ProjectiveSpace(QQ, 3)
sage: C = Curve([x*y - z*w, x^2 - y*w, y^2*w - x*z*w], P)
sage: p = C(1,1,1,1)
sage: C.tangent_line(p)
Projective Curve over Rational Field
 defined by -2*x + y + w, -3*x + z + 2*w
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(3), names=('x', 'y', 'z', 'w',)); (x, y, z, w,) = P._first_ngens(4)
>>> C = Curve([x*y - z*w, x**Integer(2) - y*w, y**Integer(2)*w - x*z*w], P)
>>> p = C(Integer(1),Integer(1),Integer(1),Integer(1))
>>> C.tangent_line(p)
Projective Curve over Rational Field
 defined by -2*x + y + w, -3*x + z + 2*w
class sage.schemes.curves.projective_curve.ProjectivePlaneCurve(A, f, category=None)[source]#

Bases: ProjectiveCurve

Curves in projective planes.

INPUT:

  • A – projective plane

  • f – homogeneous polynomial in the homogeneous coordinate ring of the plane

EXAMPLES:

A projective plane curve defined over an algebraic closure of \(\QQ\):

sage: # needs sage.rings.number_field
sage: P.<x,y,z> = ProjectiveSpace(QQbar, 2)
sage: set_verbose(-1)  # suppress warnings for slow computation
sage: C = Curve([y*z - x^2 - QQbar.gen()*z^2], P); C
Projective Plane Curve over Algebraic Field
 defined by -x^2 + y*z + (-I)*z^2
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> P = ProjectiveSpace(QQbar, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> set_verbose(-Integer(1))  # suppress warnings for slow computation
>>> C = Curve([y*z - x**Integer(2) - QQbar.gen()*z**Integer(2)], P); C
Projective Plane Curve over Algebraic Field
 defined by -x^2 + y*z + (-I)*z^2

A projective plane curve defined over a finite field:

sage: # needs sage.rings.finite_rings
sage: P.<x,y,z> = ProjectiveSpace(GF(5^2, 'v'), 2)
sage: C = Curve([y^2*z - x*z^2 - z^3], P); C
Projective Plane Curve over Finite Field in v of size 5^2
 defined by y^2*z - x*z^2 - z^3
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> P = ProjectiveSpace(GF(Integer(5)**Integer(2), 'v'), Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([y**Integer(2)*z - x*z**Integer(2) - z**Integer(3)], P); C
Projective Plane Curve over Finite Field in v of size 5^2
 defined by y^2*z - x*z^2 - z^3
degree()[source]#

Return the degree of this projective curve.

For a plane curve, this is just the degree of its defining polynomial.

OUTPUT: An integer.

EXAMPLES:

sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = P.curve([y^7 - x^2*z^5 + 7*z^7])
sage: C.degree()
7
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = P.curve([y**Integer(7) - x**Integer(2)*z**Integer(5) + Integer(7)*z**Integer(7)])
>>> C.degree()
7
divisor_of_function(r)[source]#

Return the divisor of a function on a curve.

INPUT: r is a rational function on X

OUTPUT: A 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: FF = FiniteField(5)
sage: P2 = ProjectiveSpace(2, FF, names=['x','y','z'])
sage: R = P2.coordinate_ring()
sage: x, y, z = R.gens()
sage: f = y^2*z^7 - x^9 - x*z^8
sage: C = Curve(f)
sage: K = FractionField(R)
sage: r = 1/x
sage: C.divisor_of_function(r)     # todo: not implemented  !!!!
[[-1, (0, 0, 1)]]
sage: r = 1/x^3
sage: C.divisor_of_function(r)     # todo: not implemented  !!!!
[[-3, (0, 0, 1)]]
>>> from sage.all import *
>>> FF = FiniteField(Integer(5))
>>> P2 = ProjectiveSpace(Integer(2), FF, names=['x','y','z'])
>>> R = P2.coordinate_ring()
>>> x, y, z = R.gens()
>>> f = y**Integer(2)*z**Integer(7) - x**Integer(9) - x*z**Integer(8)
>>> C = Curve(f)
>>> K = FractionField(R)
>>> r = Integer(1)/x
>>> C.divisor_of_function(r)     # todo: not implemented  !!!!
[[-1, (0, 0, 1)]]
>>> r = Integer(1)/x**Integer(3)
>>> C.divisor_of_function(r)     # todo: not implemented  !!!!
[[-3, (0, 0, 1)]]
excellent_position(Q)[source]#

Return a transformation of this curve into one in excellent position with respect to the point Q.

Here excellent position is defined as in [Ful1989]. A curve \(C\) of degree \(d\) containing the point \((0 : 0 : 1)\) with multiplicity \(r\) is said to be in excellent position if none of the coordinate lines are tangent to \(C\) at any of the fundamental points \((1 : 0 : 0)\), \((0 : 1 : 0)\), and \((0 : 0 : 1)\), and if the two coordinate lines containing \((0 : 0 : 1)\) intersect \(C\) transversally in \(d - r\) distinct non-fundamental points, and if the other coordinate line intersects \(C\) transversally at \(d\) distinct, non-fundamental points.

INPUT:

  • Q – a point on this curve.

OUTPUT:

  • a scheme morphism from this curve to a curve in excellent position that is a restriction of a change of coordinates map of the projective plane.

EXAMPLES:

sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve([x*y - z^2], P)
sage: Q = P([1,1,1])
sage: C.excellent_position(Q)
Scheme morphism:
  From: Projective Plane Curve over Rational Field defined by x*y - z^2
  To:   Projective Plane Curve over Rational Field
        defined by -x^2 - 3*x*y - 4*y^2 - x*z - 3*y*z
  Defn: Defined on coordinates by sending (x : y : z) to
        (-x + 1/2*y + 1/2*z : -1/2*y + 1/2*z : x + 1/2*y - 1/2*z)
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([x*y - z**Integer(2)], P)
>>> Q = P([Integer(1),Integer(1),Integer(1)])
>>> C.excellent_position(Q)
Scheme morphism:
  From: Projective Plane Curve over Rational Field defined by x*y - z^2
  To:   Projective Plane Curve over Rational Field
        defined by -x^2 - 3*x*y - 4*y^2 - x*z - 3*y*z
  Defn: Defined on coordinates by sending (x : y : z) to
        (-x + 1/2*y + 1/2*z : -1/2*y + 1/2*z : x + 1/2*y - 1/2*z)
sage: # needs sage.rings.number_field
sage: R.<a> = QQ[]
sage: K.<b> = NumberField(a^2 - 3)
sage: P.<x,y,z> = ProjectiveSpace(K, 2)
sage: C = P.curve([z^2*y^3*x^4 - y^6*x^3 - 4*z^2*y^4*x^3 - 4*z^4*y^2*x^3
....:              + 3*y^7*x^2 + 10*z^2*y^5*x^2 + 9*z^4*y^3*x^2
....:              + 5*z^6*y*x^2 - 3*y^8*x - 9*z^2*y^6*x - 11*z^4*y^4*x
....:              - 7*z^6*y^2*x - 2*z^8*x + y^9 + 2*z^2*y^7 + 3*z^4*y^5
....:              + 4*z^6*y^3 + 2*z^8*y])
sage: Q = P([1,0,0])
sage: C.excellent_position(Q)
Scheme morphism:
  From: Projective Plane Curve over Number Field in b
        with defining polynomial a^2 - 3
        defined by -x^3*y^6 + 3*x^2*y^7 - 3*x*y^8 + y^9 + x^4*y^3*z^2
                   - 4*x^3*y^4*z^2 + 10*x^2*y^5*z^2 - 9*x*y^6*z^2
                   + 2*y^7*z^2 - 4*x^3*y^2*z^4 + 9*x^2*y^3*z^4
                   - 11*x*y^4*z^4 + 3*y^5*z^4 + 5*x^2*y*z^6
                   - 7*x*y^2*z^6 + 4*y^3*z^6 - 2*x*z^8 + 2*y*z^8
  To:   Projective Plane Curve over Number Field in b
        with defining polynomial a^2 - 3
        defined by 900*x^9 - 7410*x^8*y + 29282*x^7*y^2 - 69710*x^6*y^3
                   + 110818*x^5*y^4 - 123178*x^4*y^5 + 96550*x^3*y^6
                   - 52570*x^2*y^7 + 18194*x*y^8 - 3388*y^9 - 1550*x^8*z
                   + 9892*x^7*y*z - 30756*x^6*y^2*z + 58692*x^5*y^3*z
                   - 75600*x^4*y^4*z + 67916*x^3*y^5*z - 42364*x^2*y^6*z
                   + 16844*x*y^7*z - 3586*y^8*z + 786*x^7*z^2
                   - 3958*x^6*y*z^2 + 9746*x^5*y^2*z^2 - 14694*x^4*y^3*z^2
                   + 15174*x^3*y^4*z^2 - 10802*x^2*y^5*z^2
                   + 5014*x*y^6*z^2 - 1266*y^7*z^2 - 144*x^6*z^3
                   + 512*x^5*y*z^3 - 912*x^4*y^2*z^3 + 1024*x^3*y^3*z^3
                   - 816*x^2*y^4*z^3 + 512*x*y^5*z^3 - 176*y^6*z^3
                   + 8*x^5*z^4 - 8*x^4*y*z^4 - 16*x^3*y^2*z^4
                   + 16*x^2*y^3*z^4 + 8*x*y^4*z^4 - 8*y^5*z^4
  Defn: Defined on coordinates by sending (x : y : z) to
        (1/4*y + 1/2*z : -1/4*y + 1/2*z : x + 1/4*y - 1/2*z)
>>> 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)
>>> P = ProjectiveSpace(K, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = P.curve([z**Integer(2)*y**Integer(3)*x**Integer(4) - y**Integer(6)*x**Integer(3) - Integer(4)*z**Integer(2)*y**Integer(4)*x**Integer(3) - Integer(4)*z**Integer(4)*y**Integer(2)*x**Integer(3)
...              + Integer(3)*y**Integer(7)*x**Integer(2) + Integer(10)*z**Integer(2)*y**Integer(5)*x**Integer(2) + Integer(9)*z**Integer(4)*y**Integer(3)*x**Integer(2)
...              + Integer(5)*z**Integer(6)*y*x**Integer(2) - Integer(3)*y**Integer(8)*x - Integer(9)*z**Integer(2)*y**Integer(6)*x - Integer(11)*z**Integer(4)*y**Integer(4)*x
...              - Integer(7)*z**Integer(6)*y**Integer(2)*x - Integer(2)*z**Integer(8)*x + y**Integer(9) + Integer(2)*z**Integer(2)*y**Integer(7) + Integer(3)*z**Integer(4)*y**Integer(5)
...              + Integer(4)*z**Integer(6)*y**Integer(3) + Integer(2)*z**Integer(8)*y])
>>> Q = P([Integer(1),Integer(0),Integer(0)])
>>> C.excellent_position(Q)
Scheme morphism:
  From: Projective Plane Curve over Number Field in b
        with defining polynomial a^2 - 3
        defined by -x^3*y^6 + 3*x^2*y^7 - 3*x*y^8 + y^9 + x^4*y^3*z^2
                   - 4*x^3*y^4*z^2 + 10*x^2*y^5*z^2 - 9*x*y^6*z^2
                   + 2*y^7*z^2 - 4*x^3*y^2*z^4 + 9*x^2*y^3*z^4
                   - 11*x*y^4*z^4 + 3*y^5*z^4 + 5*x^2*y*z^6
                   - 7*x*y^2*z^6 + 4*y^3*z^6 - 2*x*z^8 + 2*y*z^8
  To:   Projective Plane Curve over Number Field in b
        with defining polynomial a^2 - 3
        defined by 900*x^9 - 7410*x^8*y + 29282*x^7*y^2 - 69710*x^6*y^3
                   + 110818*x^5*y^4 - 123178*x^4*y^5 + 96550*x^3*y^6
                   - 52570*x^2*y^7 + 18194*x*y^8 - 3388*y^9 - 1550*x^8*z
                   + 9892*x^7*y*z - 30756*x^6*y^2*z + 58692*x^5*y^3*z
                   - 75600*x^4*y^4*z + 67916*x^3*y^5*z - 42364*x^2*y^6*z
                   + 16844*x*y^7*z - 3586*y^8*z + 786*x^7*z^2
                   - 3958*x^6*y*z^2 + 9746*x^5*y^2*z^2 - 14694*x^4*y^3*z^2
                   + 15174*x^3*y^4*z^2 - 10802*x^2*y^5*z^2
                   + 5014*x*y^6*z^2 - 1266*y^7*z^2 - 144*x^6*z^3
                   + 512*x^5*y*z^3 - 912*x^4*y^2*z^3 + 1024*x^3*y^3*z^3
                   - 816*x^2*y^4*z^3 + 512*x*y^5*z^3 - 176*y^6*z^3
                   + 8*x^5*z^4 - 8*x^4*y*z^4 - 16*x^3*y^2*z^4
                   + 16*x^2*y^3*z^4 + 8*x*y^4*z^4 - 8*y^5*z^4
  Defn: Defined on coordinates by sending (x : y : z) to
        (1/4*y + 1/2*z : -1/4*y + 1/2*z : x + 1/4*y - 1/2*z)
sage: # needs sage.rings.number_field sage.symbolic
sage: set_verbose(-1)
sage: a = QQbar(sqrt(2))
sage: P.<x,y,z> = ProjectiveSpace(QQbar, 2)
sage: C = Curve([(-1/4*a)*x^3 + (-3/4*a)*x^2*y
....:            + (-3/4*a)*x*y^2 + (-1/4*a)*y^3 - 2*x*y*z], P)
sage: Q = P([0,0,1])
sage: C.excellent_position(Q)
Scheme morphism:
  From: Projective Plane Curve over Algebraic Field defined
        by (-0.3535533905932738?)*x^3 + (-1.060660171779822?)*x^2*y
           + (-1.060660171779822?)*x*y^2 + (-0.3535533905932738?)*y^3
           + (-2)*x*y*z
  To:   Projective Plane Curve over Algebraic Field defined
        by (-2.828427124746190?)*x^3 + (-2)*x^2*y + 2*y^3
           + (-2)*x^2*z + 2*y^2*z
  Defn: Defined on coordinates by sending (x : y : z) to
        (1/2*x + 1/2*y : (-1/2)*x + 1/2*y : 1/2*x + (-1/2)*y + z)
>>> from sage.all import *
>>> # needs sage.rings.number_field sage.symbolic
>>> set_verbose(-Integer(1))
>>> a = QQbar(sqrt(Integer(2)))
>>> P = ProjectiveSpace(QQbar, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([(-Integer(1)/Integer(4)*a)*x**Integer(3) + (-Integer(3)/Integer(4)*a)*x**Integer(2)*y
...            + (-Integer(3)/Integer(4)*a)*x*y**Integer(2) + (-Integer(1)/Integer(4)*a)*y**Integer(3) - Integer(2)*x*y*z], P)
>>> Q = P([Integer(0),Integer(0),Integer(1)])
>>> C.excellent_position(Q)
Scheme morphism:
  From: Projective Plane Curve over Algebraic Field defined
        by (-0.3535533905932738?)*x^3 + (-1.060660171779822?)*x^2*y
           + (-1.060660171779822?)*x*y^2 + (-0.3535533905932738?)*y^3
           + (-2)*x*y*z
  To:   Projective Plane Curve over Algebraic Field defined
        by (-2.828427124746190?)*x^3 + (-2)*x^2*y + 2*y^3
           + (-2)*x^2*z + 2*y^2*z
  Defn: Defined on coordinates by sending (x : y : z) to
        (1/2*x + 1/2*y : (-1/2)*x + 1/2*y : 1/2*x + (-1/2)*y + z)
is_ordinary_singularity(P)[source]#

Return whether the singular point P of this projective 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:

  • Boolean. 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: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve([y^2*z^3 - x^5], P)
sage: Q = P([0,0,1])
sage: C.is_ordinary_singularity(Q)
False
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([y**Integer(2)*z**Integer(3) - x**Integer(5)], P)
>>> Q = P([Integer(0),Integer(0),Integer(1)])
>>> C.is_ordinary_singularity(Q)
False
sage: # needs sage.rings.number_field
sage: R.<a> = QQ[]
sage: K.<b> = NumberField(a^2 - 3)
sage: P.<x,y,z> = ProjectiveSpace(K, 2)
sage: C = P.curve([x^2*y^3*z^4 - y^6*z^3 - 4*x^2*y^4*z^3 - 4*x^4*y^2*z^3
....:              + 3*y^7*z^2 + 10*x^2*y^5*z^2 + 9*x^4*y^3*z^2
....:              + 5*x^6*y*z^2 - 3*y^8*z - 9*x^2*y^6*z - 11*x^4*y^4*z
....:              - 7*x^6*y^2*z - 2*x^8*z + y^9 + 2*x^2*y^7 + 3*x^4*y^5
....:              + 4*x^6*y^3 + 2*x^8*y])
sage: Q = P([0,1,1])
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)
>>> P = ProjectiveSpace(K, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = P.curve([x**Integer(2)*y**Integer(3)*z**Integer(4) - y**Integer(6)*z**Integer(3) - Integer(4)*x**Integer(2)*y**Integer(4)*z**Integer(3) - Integer(4)*x**Integer(4)*y**Integer(2)*z**Integer(3)
...              + Integer(3)*y**Integer(7)*z**Integer(2) + Integer(10)*x**Integer(2)*y**Integer(5)*z**Integer(2) + Integer(9)*x**Integer(4)*y**Integer(3)*z**Integer(2)
...              + Integer(5)*x**Integer(6)*y*z**Integer(2) - Integer(3)*y**Integer(8)*z - Integer(9)*x**Integer(2)*y**Integer(6)*z - Integer(11)*x**Integer(4)*y**Integer(4)*z
...              - Integer(7)*x**Integer(6)*y**Integer(2)*z - Integer(2)*x**Integer(8)*z + y**Integer(9) + Integer(2)*x**Integer(2)*y**Integer(7) + Integer(3)*x**Integer(4)*y**Integer(5)
...              + Integer(4)*x**Integer(6)*y**Integer(3) + Integer(2)*x**Integer(8)*y])
>>> Q = P([Integer(0),Integer(1),Integer(1)])
>>> C.is_ordinary_singularity(Q)
True
sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = P.curve([z^5 - y^5 + x^5 + x*y^2*z^2])
sage: Q = P([0,1,1])
sage: C.is_ordinary_singularity(Q)
Traceback (most recent call last):
...
TypeError: (=(0 : 1 : 1)) is not a singular point of (=Projective Plane
Curve over Rational Field defined by x^5 - y^5 + x*y^2*z^2 + z^5)
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = P.curve([z**Integer(5) - y**Integer(5) + x**Integer(5) + x*y**Integer(2)*z**Integer(2)])
>>> Q = P([Integer(0),Integer(1),Integer(1)])
>>> C.is_ordinary_singularity(Q)
Traceback (most recent call last):
...
TypeError: (=(0 : 1 : 1)) is not a singular point of (=Projective Plane
Curve over Rational Field defined by x^5 - y^5 + x*y^2*z^2 + z^5)
is_singular(P=None)[source]#

Return whether this curve is singular or not, or if a point P is provided, whether P is a singular point of this curve.

INPUT:

  • P – (default: None) a point on this curve

OUTPUT:

If no point P is provided, return True or False depending on whether this curve is singular or not. If a point P is provided, return True or False depending on whether P is or is not a singular point of this curve.

EXAMPLES:

Over \(\QQ\):

sage: F = QQ
sage: P2.<X,Y,Z> = ProjectiveSpace(F, 2)
sage: C = Curve(X^3 - Y^2*Z)
sage: C.is_singular()
True
>>> from sage.all import *
>>> F = QQ
>>> P2 = ProjectiveSpace(F, Integer(2), names=('X', 'Y', 'Z',)); (X, Y, Z,) = P2._first_ngens(3)
>>> C = Curve(X**Integer(3) - Y**Integer(2)*Z)
>>> C.is_singular()
True

Over a finite field:

sage: F = GF(19)
sage: P2.<X,Y,Z> = ProjectiveSpace(F, 2)
sage: C = Curve(X^3 + Y^3 + Z^3)
sage: C.is_singular()
False
sage: D = Curve(X^4 - X*Z^3)
sage: D.is_singular()
True
sage: E = Curve(X^5 + 19*Y^5 + Z^5)
sage: E.is_singular()
True
sage: E = Curve(X^5 + 9*Y^5 + Z^5)
sage: E.is_singular()
False
>>> from sage.all import *
>>> F = GF(Integer(19))
>>> P2 = ProjectiveSpace(F, Integer(2), names=('X', 'Y', 'Z',)); (X, Y, Z,) = P2._first_ngens(3)
>>> C = Curve(X**Integer(3) + Y**Integer(3) + Z**Integer(3))
>>> C.is_singular()
False
>>> D = Curve(X**Integer(4) - X*Z**Integer(3))
>>> D.is_singular()
True
>>> E = Curve(X**Integer(5) + Integer(19)*Y**Integer(5) + Z**Integer(5))
>>> E.is_singular()
True
>>> E = Curve(X**Integer(5) + Integer(9)*Y**Integer(5) + Z**Integer(5))
>>> E.is_singular()
False

Over \(\CC\):

sage: # needs sage.rings.function_field
sage: F = CC
sage: P2.<X,Y,Z> = ProjectiveSpace(F, 2)
sage: C = Curve(X)
sage: C.is_singular()
False
sage: D = Curve(Y^2*Z - X^3)
sage: D.is_singular()
True
sage: E = Curve(Y^2*Z - X^3 + Z^3)
sage: E.is_singular()
False
>>> from sage.all import *
>>> # needs sage.rings.function_field
>>> F = CC
>>> P2 = ProjectiveSpace(F, Integer(2), names=('X', 'Y', 'Z',)); (X, Y, Z,) = P2._first_ngens(3)
>>> C = Curve(X)
>>> C.is_singular()
False
>>> D = Curve(Y**Integer(2)*Z - X**Integer(3))
>>> D.is_singular()
True
>>> E = Curve(Y**Integer(2)*Z - X**Integer(3) + Z**Integer(3))
>>> E.is_singular()
False

Showing that Issue #12187 is fixed:

sage: F.<X,Y,Z> = GF(2)[]
sage: G = Curve(X^2 + Y*Z)
sage: G.is_singular()
False
>>> from sage.all import *
>>> F = GF(Integer(2))['X, Y, Z']; (X, Y, Z,) = F._first_ngens(3)
>>> G = Curve(X**Integer(2) + Y*Z)
>>> G.is_singular()
False
sage: # needs sage.fings.function_field
sage: P.<x,y,z> = ProjectiveSpace(CC, 2)
sage: C = Curve([y^4 - x^3*z], P)
sage: Q = P([0,0,1])
sage: C.is_singular()
True
>>> from sage.all import *
>>> # needs sage.fings.function_field
>>> P = ProjectiveSpace(CC, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([y**Integer(4) - x**Integer(3)*z], P)
>>> Q = P([Integer(0),Integer(0),Integer(1)])
>>> C.is_singular()
True
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: A boolean.

EXAMPLES:

sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve([x^2 - y^2], P)
sage: D = Curve([x - y], P)
sage: Q = P([1,1,0])
sage: C.is_transverse(D, Q)
False
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([x**Integer(2) - y**Integer(2)], P)
>>> D = Curve([x - y], P)
>>> Q = P([Integer(1),Integer(1),Integer(0)])
>>> C.is_transverse(D, Q)
False
sage: # needs sage.rings.number_field
sage: K = QuadraticField(-1)
sage: P.<x,y,z> = ProjectiveSpace(K, 2)
sage: C = Curve([y^2*z - K.0*x^3], P)
sage: D = Curve([z*x + y^2], P)
sage: Q = P([0,0,1])
sage: C.is_transverse(D, Q)
False
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> K = QuadraticField(-Integer(1))
>>> P = ProjectiveSpace(K, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([y**Integer(2)*z - K.gen(0)*x**Integer(3)], P)
>>> D = Curve([z*x + y**Integer(2)], P)
>>> Q = P([Integer(0),Integer(0),Integer(1)])
>>> C.is_transverse(D, Q)
False
sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve([x^2 - 2*y^2 - 2*z^2], P)
sage: D = Curve([y - z], P)
sage: Q = P([2,1,1])
sage: C.is_transverse(D, Q)
True
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([x**Integer(2) - Integer(2)*y**Integer(2) - Integer(2)*z**Integer(2)], P)
>>> D = Curve([y - z], P)
>>> Q = P([Integer(2),Integer(1),Integer(1)])
>>> 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 worse than others.

INPUT:

  • pt – a rational point on X which is not a point of ramification

    for the projection \((x,y) \to x\).

  • n – the number of terms desired

OUTPUT: \(x = x0 + t\), \(y = y0\) + power series in \(t\)

EXAMPLES:

sage: FF = FiniteField(5)
sage: P2 = ProjectiveSpace(2, FF, names=['x','y','z'])
sage: x, y, z = P2.coordinate_ring().gens()
sage: C = Curve(y^2*z^7 - x^9 - x*z^8)
sage: pt = C([2,3,1])
sage: C.local_coordinates(pt,9)     # todo: not implemented  !!!!
[2 + t,
 3 + 3*t^2 + t^3 + 3*t^4 + 3*t^6 + 3*t^7 + t^8 + 2*t^9 + 3*t^11 + 3*t^12]
>>> from sage.all import *
>>> FF = FiniteField(Integer(5))
>>> P2 = ProjectiveSpace(Integer(2), FF, names=['x','y','z'])
>>> x, y, z = P2.coordinate_ring().gens()
>>> C = Curve(y**Integer(2)*z**Integer(7) - x**Integer(9) - x*z**Integer(8))
>>> pt = C([Integer(2),Integer(3),Integer(1)])
>>> C.local_coordinates(pt,Integer(9))     # todo: not implemented  !!!!
[2 + t,
 3 + 3*t^2 + t^3 + 3*t^4 + 3*t^6 + 3*t^7 + t^8 + 2*t^9 + 3*t^11 + 3*t^12]
ordinary_model()[source]#

Return a birational map from this curve to a plane curve with only ordinary singularities.

Currently only implemented over number fields. If not all of the coordinates of the non-ordinary singularities of this curve are contained in its base field, then the domain and codomain of the map returned will be defined over an extension. This curve must be irreducible.

OUTPUT:

  • a scheme morphism from this curve to a curve with only ordinary singularities that defines a birational map between the two curves.

EXAMPLES:

sage: # needs sage.rings.number_field
sage: set_verbose(-1)
sage: K = QuadraticField(3)
sage: P.<x,y,z> = ProjectiveSpace(K, 2)
sage: C = Curve([x^5 - K.0*y*z^4], P)
sage: C.ordinary_model()
Scheme morphism:
  From: Projective Plane Curve over Number Field in a
        with defining polynomial x^2 - 3 with a = 1.732050807568878?
        defined by x^5 + (-a)*y*z^4
  To:   Projective Plane Curve over Number Field in a
        with defining polynomial x^2 - 3 with a = 1.732050807568878?
        defined by (-a)*x^5*y + (-4*a)*x^4*y^2 + (-6*a)*x^3*y^3
                   + (-4*a)*x^2*y^4 + (-a)*x*y^5 + (-a - 1)*x^5*z
                   + (-4*a + 5)*x^4*y*z + (-6*a - 10)*x^3*y^2*z
                   + (-4*a + 10)*x^2*y^3*z + (-a - 5)*x*y^4*z + y^5*z
  Defn: Defined on coordinates by sending (x : y : z) to
        (-1/4*x^2 - 1/2*x*y + 1/2*x*z + 1/2*y*z - 1/4*z^2 :
         1/4*x^2 + 1/2*x*y + 1/2*y*z - 1/4*z^2 :
         -1/4*x^2 + 1/4*z^2)
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> set_verbose(-Integer(1))
>>> K = QuadraticField(Integer(3))
>>> P = ProjectiveSpace(K, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([x**Integer(5) - K.gen(0)*y*z**Integer(4)], P)
>>> C.ordinary_model()
Scheme morphism:
  From: Projective Plane Curve over Number Field in a
        with defining polynomial x^2 - 3 with a = 1.732050807568878?
        defined by x^5 + (-a)*y*z^4
  To:   Projective Plane Curve over Number Field in a
        with defining polynomial x^2 - 3 with a = 1.732050807568878?
        defined by (-a)*x^5*y + (-4*a)*x^4*y^2 + (-6*a)*x^3*y^3
                   + (-4*a)*x^2*y^4 + (-a)*x*y^5 + (-a - 1)*x^5*z
                   + (-4*a + 5)*x^4*y*z + (-6*a - 10)*x^3*y^2*z
                   + (-4*a + 10)*x^2*y^3*z + (-a - 5)*x*y^4*z + y^5*z
  Defn: Defined on coordinates by sending (x : y : z) to
        (-1/4*x^2 - 1/2*x*y + 1/2*x*z + 1/2*y*z - 1/4*z^2 :
         1/4*x^2 + 1/2*x*y + 1/2*y*z - 1/4*z^2 :
         -1/4*x^2 + 1/4*z^2)
sage: set_verbose(-1)
sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve([y^2*z^2 - x^4 - x^3*z], P)
sage: D = C.ordinary_model(); D  # long time (2 seconds)
Scheme morphism:
  From: Projective Plane Curve over Rational Field defined
        by -x^4 - x^3*z + y^2*z^2
  To:   Projective Plane Curve over Rational Field defined
        by 4*x^6*y^3 - 24*x^5*y^4 + 36*x^4*y^5 + 8*x^6*y^2*z
           - 40*x^5*y^3*z + 24*x^4*y^4*z + 72*x^3*y^5*z - 4*x^6*y*z^2
           + 8*x^5*y^2*z^2 - 56*x^4*y^3*z^2 + 104*x^3*y^4*z^2
           + 44*x^2*y^5*z^2 + 8*x^6*z^3 - 16*x^5*y*z^3
           - 24*x^4*y^2*z^3 + 40*x^3*y^3*z^3 + 48*x^2*y^4*z^3
           + 8*x*y^5*z^3 - 8*x^5*z^4 + 36*x^4*y*z^4 - 56*x^3*y^2*z^4
           + 20*x^2*y^3*z^4 + 40*x*y^4*z^4 - 16*y^5*z^4
  Defn: Defined on coordinates by sending (x : y : z) to
        (-3/64*x^4 + 9/64*x^2*y^2 - 3/32*x*y^3 - 1/16*x^3*z
          - 1/8*x^2*y*z + 1/4*x*y^2*z - 1/16*y^3*z - 1/8*x*y*z^2
          + 1/16*y^2*z^2 :
         -1/64*x^4 + 3/64*x^2*y^2 - 1/32*x*y^3 + 1/16*x*y^2*z
          - 1/16*y^3*z + 1/16*y^2*z^2 :
         3/64*x^4 - 3/32*x^3*y + 3/64*x^2*y^2 + 1/16*x^3*z
          - 3/16*x^2*y*z + 1/8*x*y^2*z - 1/8*x*y*z^2 + 1/16*y^2*z^2)
sage: all(D.codomain().is_ordinary_singularity(Q)  # long time
....:     for Q in D.codomain().singular_points())
True
>>> from sage.all import *
>>> set_verbose(-Integer(1))
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([y**Integer(2)*z**Integer(2) - x**Integer(4) - x**Integer(3)*z], P)
>>> D = C.ordinary_model(); D  # long time (2 seconds)
Scheme morphism:
  From: Projective Plane Curve over Rational Field defined
        by -x^4 - x^3*z + y^2*z^2
  To:   Projective Plane Curve over Rational Field defined
        by 4*x^6*y^3 - 24*x^5*y^4 + 36*x^4*y^5 + 8*x^6*y^2*z
           - 40*x^5*y^3*z + 24*x^4*y^4*z + 72*x^3*y^5*z - 4*x^6*y*z^2
           + 8*x^5*y^2*z^2 - 56*x^4*y^3*z^2 + 104*x^3*y^4*z^2
           + 44*x^2*y^5*z^2 + 8*x^6*z^3 - 16*x^5*y*z^3
           - 24*x^4*y^2*z^3 + 40*x^3*y^3*z^3 + 48*x^2*y^4*z^3
           + 8*x*y^5*z^3 - 8*x^5*z^4 + 36*x^4*y*z^4 - 56*x^3*y^2*z^4
           + 20*x^2*y^3*z^4 + 40*x*y^4*z^4 - 16*y^5*z^4
  Defn: Defined on coordinates by sending (x : y : z) to
        (-3/64*x^4 + 9/64*x^2*y^2 - 3/32*x*y^3 - 1/16*x^3*z
          - 1/8*x^2*y*z + 1/4*x*y^2*z - 1/16*y^3*z - 1/8*x*y*z^2
          + 1/16*y^2*z^2 :
         -1/64*x^4 + 3/64*x^2*y^2 - 1/32*x*y^3 + 1/16*x*y^2*z
          - 1/16*y^3*z + 1/16*y^2*z^2 :
         3/64*x^4 - 3/32*x^3*y + 3/64*x^2*y^2 + 1/16*x^3*z
          - 3/16*x^2*y*z + 1/8*x*y^2*z - 1/8*x*y*z^2 + 1/16*y^2*z^2)
>>> all(D.codomain().is_ordinary_singularity(Q)  # long time
...     for Q in D.codomain().singular_points())
True
sage: set_verbose(-1)
sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve([(x^2 + y^2 - y*z - 2*z^2)*(y*z - x^2 + 2*z^2)*z + y^5], P)
sage: C.ordinary_model() # long time (5 seconds)
Scheme morphism:
  From: Projective Plane Curve over Number Field in a
        with defining polynomial y^2 - 2 defined
        by y^5 - x^4*z - x^2*y^2*z + 2*x^2*y*z^2 + y^3*z^2
           + 4*x^2*z^3 + y^2*z^3 - 4*y*z^4 - 4*z^5
  To:   Projective Plane Curve over Number Field in a
        with defining polynomial y^2 - 2 defined
        by (-29*a + 1)*x^8*y^6 + (10*a + 158)*x^7*y^7
           + (-109*a - 31)*x^6*y^8 + (-80*a - 198)*x^8*y^5*z
           + (531*a + 272)*x^7*y^6*z + (170*a - 718)*x^6*y^7*z
           + (19*a - 636)*x^5*y^8*z + (-200*a - 628)*x^8*y^4*z^2
           + (1557*a - 114)*x^7*y^5*z^2 + (2197*a - 2449)*x^6*y^6*z^2
           + (1223*a - 3800)*x^5*y^7*z^2 + (343*a - 1329)*x^4*y^8*z^2
           + (-323*a - 809)*x^8*y^3*z^3 + (1630*a - 631)*x^7*y^4*z^3
           + (4190*a - 3126)*x^6*y^5*z^3 + (3904*a - 7110)*x^5*y^6*z^3
           + (1789*a - 5161)*x^4*y^7*z^3 + (330*a - 1083)*x^3*y^8*z^3
           + (-259*a - 524)*x^8*y^2*z^4 + (720*a - 605)*x^7*y^3*z^4
           + (3082*a - 2011)*x^6*y^4*z^4 + (4548*a - 5462)*x^5*y^5*z^4
           + (2958*a - 6611)*x^4*y^6*z^4 + (994*a - 2931)*x^3*y^7*z^4
           + (117*a - 416)*x^2*y^8*z^4 + (-108*a - 184)*x^8*y*z^5
           + (169*a - 168)*x^7*y^2*z^5 + (831*a - 835)*x^6*y^3*z^5
           + (2225*a - 1725)*x^5*y^4*z^5 + (1970*a - 3316)*x^4*y^5*z^5
           + (952*a - 2442)*x^3*y^6*z^5 + (217*a - 725)*x^2*y^7*z^5
           + (16*a - 77)*x*y^8*z^5 + (-23*a - 35)*x^8*z^6
           + (43*a + 24)*x^7*y*z^6 + (21*a - 198)*x^6*y^2*z^6
           + (377*a - 179)*x^5*y^3*z^6 + (458*a - 537)*x^4*y^4*z^6
           + (288*a - 624)*x^3*y^5*z^6 + (100*a - 299)*x^2*y^6*z^6
           + (16*a - 67)*x*y^7*z^6 - 5*y^8*z^6
  Defn: Defined on coordinates by sending (x : y : z) to
        ((-5/128*a - 5/128)*x^4 + (-5/32*a + 5/32)*x^3*y
          + (-1/16*a + 3/32)*x^2*y^2 + (1/16*a - 1/16)*x*y^3
          + (1/32*a - 1/32)*y^4 - 1/32*x^3*z + (3/16*a - 5/8)*x^2*y*z
          + (1/8*a - 5/16)*x*y^2*z + (1/8*a + 5/32)*x^2*z^2
          + (-3/16*a + 5/16)*x*y*z^2 + (-3/16*a - 1/16)*y^2*z^2
          + 1/16*x*z^3 + (1/4*a + 1/4)*y*z^3 + (-3/32*a - 5/32)*z^4 :
         (-5/128*a - 5/128)*x^4 + (5/32*a)*x^3*y
          + (3/32*a + 3/32)*x^2*y^2 + (-1/16*a)*x*y^3
          + (-1/32*a - 1/32)*y^4 - 1/32*x^3*z + (-11/32*a)*x^2*y*z
          + (1/8*a + 5/16)*x*y^2*z + (3/16*a + 1/4)*y^3*z
          + (1/8*a + 5/32)*x^2*z^2 + (-1/16*a - 3/8)*x*y*z^2
          + (-3/8*a - 9/16)*y^2*z^2 + 1/16*x*z^3 + (5/16*a + 1/2)*y*z^3
          + (-3/32*a - 5/32)*z^4 :
         (1/64*a + 3/128)*x^4 + (-1/32*a - 1/32)*x^3*y
          + (3/32*a - 9/32)*x^2*y^2 + (1/16*a - 3/16)*x*y^3 - 1/32*y^4
          + (3/32*a + 1/8)*x^2*y*z + (-1/8*a + 1/8)*x*y^2*z
          + (-1/16*a)*y^3*z + (-1/16*a - 3/32)*x^2*z^2
          + (1/16*a + 1/16)*x*y*z^2 + (3/16*a + 3/16)*y^2*z^2
          + (-3/16*a - 1/4)*y*z^3 + (1/16*a + 3/32)*z^4)
>>> from sage.all import *
>>> set_verbose(-Integer(1))
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([(x**Integer(2) + y**Integer(2) - y*z - Integer(2)*z**Integer(2))*(y*z - x**Integer(2) + Integer(2)*z**Integer(2))*z + y**Integer(5)], P)
>>> C.ordinary_model() # long time (5 seconds)
Scheme morphism:
  From: Projective Plane Curve over Number Field in a
        with defining polynomial y^2 - 2 defined
        by y^5 - x^4*z - x^2*y^2*z + 2*x^2*y*z^2 + y^3*z^2
           + 4*x^2*z^3 + y^2*z^3 - 4*y*z^4 - 4*z^5
  To:   Projective Plane Curve over Number Field in a
        with defining polynomial y^2 - 2 defined
        by (-29*a + 1)*x^8*y^6 + (10*a + 158)*x^7*y^7
           + (-109*a - 31)*x^6*y^8 + (-80*a - 198)*x^8*y^5*z
           + (531*a + 272)*x^7*y^6*z + (170*a - 718)*x^6*y^7*z
           + (19*a - 636)*x^5*y^8*z + (-200*a - 628)*x^8*y^4*z^2
           + (1557*a - 114)*x^7*y^5*z^2 + (2197*a - 2449)*x^6*y^6*z^2
           + (1223*a - 3800)*x^5*y^7*z^2 + (343*a - 1329)*x^4*y^8*z^2
           + (-323*a - 809)*x^8*y^3*z^3 + (1630*a - 631)*x^7*y^4*z^3
           + (4190*a - 3126)*x^6*y^5*z^3 + (3904*a - 7110)*x^5*y^6*z^3
           + (1789*a - 5161)*x^4*y^7*z^3 + (330*a - 1083)*x^3*y^8*z^3
           + (-259*a - 524)*x^8*y^2*z^4 + (720*a - 605)*x^7*y^3*z^4
           + (3082*a - 2011)*x^6*y^4*z^4 + (4548*a - 5462)*x^5*y^5*z^4
           + (2958*a - 6611)*x^4*y^6*z^4 + (994*a - 2931)*x^3*y^7*z^4
           + (117*a - 416)*x^2*y^8*z^4 + (-108*a - 184)*x^8*y*z^5
           + (169*a - 168)*x^7*y^2*z^5 + (831*a - 835)*x^6*y^3*z^5
           + (2225*a - 1725)*x^5*y^4*z^5 + (1970*a - 3316)*x^4*y^5*z^5
           + (952*a - 2442)*x^3*y^6*z^5 + (217*a - 725)*x^2*y^7*z^5
           + (16*a - 77)*x*y^8*z^5 + (-23*a - 35)*x^8*z^6
           + (43*a + 24)*x^7*y*z^6 + (21*a - 198)*x^6*y^2*z^6
           + (377*a - 179)*x^5*y^3*z^6 + (458*a - 537)*x^4*y^4*z^6
           + (288*a - 624)*x^3*y^5*z^6 + (100*a - 299)*x^2*y^6*z^6
           + (16*a - 67)*x*y^7*z^6 - 5*y^8*z^6
  Defn: Defined on coordinates by sending (x : y : z) to
        ((-5/128*a - 5/128)*x^4 + (-5/32*a + 5/32)*x^3*y
          + (-1/16*a + 3/32)*x^2*y^2 + (1/16*a - 1/16)*x*y^3
          + (1/32*a - 1/32)*y^4 - 1/32*x^3*z + (3/16*a - 5/8)*x^2*y*z
          + (1/8*a - 5/16)*x*y^2*z + (1/8*a + 5/32)*x^2*z^2
          + (-3/16*a + 5/16)*x*y*z^2 + (-3/16*a - 1/16)*y^2*z^2
          + 1/16*x*z^3 + (1/4*a + 1/4)*y*z^3 + (-3/32*a - 5/32)*z^4 :
         (-5/128*a - 5/128)*x^4 + (5/32*a)*x^3*y
          + (3/32*a + 3/32)*x^2*y^2 + (-1/16*a)*x*y^3
          + (-1/32*a - 1/32)*y^4 - 1/32*x^3*z + (-11/32*a)*x^2*y*z
          + (1/8*a + 5/16)*x*y^2*z + (3/16*a + 1/4)*y^3*z
          + (1/8*a + 5/32)*x^2*z^2 + (-1/16*a - 3/8)*x*y*z^2
          + (-3/8*a - 9/16)*y^2*z^2 + 1/16*x*z^3 + (5/16*a + 1/2)*y*z^3
          + (-3/32*a - 5/32)*z^4 :
         (1/64*a + 3/128)*x^4 + (-1/32*a - 1/32)*x^3*y
          + (3/32*a - 9/32)*x^2*y^2 + (1/16*a - 3/16)*x*y^3 - 1/32*y^4
          + (3/32*a + 1/8)*x^2*y*z + (-1/8*a + 1/8)*x*y^2*z
          + (-1/16*a)*y^3*z + (-1/16*a - 3/32)*x^2*z^2
          + (1/16*a + 1/16)*x*y*z^2 + (3/16*a + 3/16)*y^2*z^2
          + (-3/16*a - 1/4)*y*z^3 + (1/16*a + 3/32)*z^4)
plot(*args, **kwds)[source]#

Plot the real points of an affine patch of this projective plane curve.

INPUT:

  • self – an affine plane curve

  • patch – (optional) the affine patch to be plotted; if not specified, the patch corresponding to the last projective coordinate being nonzero

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

The other affine patches of the same curve:

sage: # needs sage.plot
sage: C.plot(patch=0)
Graphics object consisting of 1 graphics primitive
sage: C.plot(patch=1)
Graphics object consisting of 1 graphics primitive
>>> from sage.all import *
>>> # needs sage.plot
>>> C.plot(patch=Integer(0))
Graphics object consisting of 1 graphics primitive
>>> C.plot(patch=Integer(1))
Graphics object consisting of 1 graphics primitive

An elliptic curve:

sage: # needs sage.plot
sage: E = EllipticCurve('101a')
sage: C = Curve(E)
sage: C.plot()
Graphics object consisting of 1 graphics primitive
sage: C.plot(patch=0)
Graphics object consisting of 1 graphics primitive
sage: C.plot(patch=1)
Graphics object consisting of 1 graphics primitive
>>> from sage.all import *
>>> # needs sage.plot
>>> E = EllipticCurve('101a')
>>> C = Curve(E)
>>> C.plot()
Graphics object consisting of 1 graphics primitive
>>> C.plot(patch=Integer(0))
Graphics object consisting of 1 graphics primitive
>>> C.plot(patch=Integer(1))
Graphics object consisting of 1 graphics primitive

A hyperelliptic curve:

sage: # needs sage.plot
sage: P.<x> = QQ[]
sage: f = 4*x^5 - 30*x^3 + 45*x - 22
sage: C = HyperellipticCurve(f)
sage: C.plot()
Graphics object consisting of 1 graphics primitive
sage: C.plot(patch=0)
Graphics object consisting of 1 graphics primitive
sage: C.plot(patch=1)
Graphics object consisting of 1 graphics primitive
>>> from sage.all import *
>>> # needs sage.plot
>>> P = QQ['x']; (x,) = P._first_ngens(1)
>>> f = Integer(4)*x**Integer(5) - Integer(30)*x**Integer(3) + Integer(45)*x - Integer(22)
>>> C = HyperellipticCurve(f)
>>> C.plot()
Graphics object consisting of 1 graphics primitive
>>> C.plot(patch=Integer(0))
Graphics object consisting of 1 graphics primitive
>>> C.plot(patch=Integer(1))
Graphics object consisting of 1 graphics primitive
quadratic_transform()[source]#

Return a birational map from this curve to the proper transform of this curve with respect to the standard Cremona transformation.

The standard Cremona transformation is the birational automorphism of \(\mathbb{P}^{2}\) defined \((x : y : z)\mapsto (yz : xz : xy)\).

OUTPUT:

  • a scheme morphism representing the restriction of the standard Cremona transformation from this curve to the proper transform.

EXAMPLES:

sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve(x^3*y - z^4 - z^2*x^2, P)
sage: C.quadratic_transform()
Scheme morphism:
  From: Projective Plane Curve over Rational Field
        defined by x^3*y - x^2*z^2 - z^4
  To:   Projective Plane Curve over Rational Field
        defined by -x^3*y - x*y*z^2 + z^4
  Defn: Defined on coordinates by sending (x : y : z) to
        (y*z : x*z : x*y)
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve(x**Integer(3)*y - z**Integer(4) - z**Integer(2)*x**Integer(2), P)
>>> C.quadratic_transform()
Scheme morphism:
  From: Projective Plane Curve over Rational Field
        defined by x^3*y - x^2*z^2 - z^4
  To:   Projective Plane Curve over Rational Field
        defined by -x^3*y - x*y*z^2 + z^4
  Defn: Defined on coordinates by sending (x : y : z) to
        (y*z : x*z : x*y)
sage: P.<x,y,z> = ProjectiveSpace(GF(17), 2)
sage: C = P.curve([y^7*z^2 - 16*x^9 + x*y*z^7 + 2*z^9])
sage: C.quadratic_transform()
Scheme morphism:
  From: Projective Plane Curve over Finite Field of size 17
        defined by x^9 + y^7*z^2 + x*y*z^7 + 2*z^9
  To:   Projective Plane Curve over Finite Field of size 17
        defined by 2*x^9*y^7 + x^8*y^6*z^2 + x^9*z^7 + y^7*z^9
  Defn: Defined on coordinates by sending (x : y : z) to
        (y*z : x*z : x*y)
>>> from sage.all import *
>>> P = ProjectiveSpace(GF(Integer(17)), Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = P.curve([y**Integer(7)*z**Integer(2) - Integer(16)*x**Integer(9) + x*y*z**Integer(7) + Integer(2)*z**Integer(9)])
>>> C.quadratic_transform()
Scheme morphism:
  From: Projective Plane Curve over Finite Field of size 17
        defined by x^9 + y^7*z^2 + x*y*z^7 + 2*z^9
  To:   Projective Plane Curve over Finite Field of size 17
        defined by 2*x^9*y^7 + x^8*y^6*z^2 + x^9*z^7 + y^7*z^9
  Defn: Defined on coordinates by sending (x : y : z) to
        (y*z : x*z : x*y)
tangents(P, factor=True)[source]#

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

These are found by homogenizing the tangents of an affine patch of this curve containing P. The point P must be a point on this curve.

INPUT:

  • P – a point on this curve.

  • factor – (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:

A list of polynomials in the coordinate ring of the ambient space of this curve.

EXAMPLES:

sage: # needs sage.rings.number_field
sage: set_verbose(-1)
sage: P.<x,y,z> = ProjectiveSpace(QQbar, 2)
sage: C = Curve([x^3*y + 2*x^2*y^2 + x*y^3 + x^3*z
....:            + 7*x^2*y*z + 14*x*y^2*z + 9*y^3*z], P)
sage: Q = P([0,0,1])
sage: C.tangents(Q)
[x + 4.147899035704788?*y,
 x + (1.426050482147607? + 0.3689894074818041?*I)*y,
 x + (1.426050482147607? - 0.3689894074818041?*I)*y]
sage: C.tangents(Q, factor=False)
[6*x^3 + 42*x^2*y + 84*x*y^2 + 54*y^3]
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> set_verbose(-Integer(1))
>>> P = ProjectiveSpace(QQbar, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([x**Integer(3)*y + Integer(2)*x**Integer(2)*y**Integer(2) + x*y**Integer(3) + x**Integer(3)*z
...            + Integer(7)*x**Integer(2)*y*z + Integer(14)*x*y**Integer(2)*z + Integer(9)*y**Integer(3)*z], P)
>>> Q = P([Integer(0),Integer(0),Integer(1)])
>>> C.tangents(Q)
[x + 4.147899035704788?*y,
 x + (1.426050482147607? + 0.3689894074818041?*I)*y,
 x + (1.426050482147607? - 0.3689894074818041?*I)*y]
>>> C.tangents(Q, factor=False)
[6*x^3 + 42*x^2*y + 84*x*y^2 + 54*y^3]
sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = P.curve([x^2*y^3*z^4 - y^6*z^3 - 4*x^2*y^4*z^3 - 4*x^4*y^2*z^3
....:              + 3*y^7*z^2 + 10*x^2*y^5*z^2 + 9*x^4*y^3*z^2 + 5*x^6*y*z^2
....:              - 3*y^8*z - 9*x^2*y^6*z - 11*x^4*y^4*z - 7*x^6*y^2*z
....:              - 2*x^8*z + y^9 + 2*x^2*y^7 + 3*x^4*y^5 + 4*x^6*y^3 + 2*x^8*y])
sage: Q = P([0,1,1])
sage: C.tangents(Q)
[-y + z, 3*x^2 - y^2 + 2*y*z - z^2]
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = P.curve([x**Integer(2)*y**Integer(3)*z**Integer(4) - y**Integer(6)*z**Integer(3) - Integer(4)*x**Integer(2)*y**Integer(4)*z**Integer(3) - Integer(4)*x**Integer(4)*y**Integer(2)*z**Integer(3)
...              + Integer(3)*y**Integer(7)*z**Integer(2) + Integer(10)*x**Integer(2)*y**Integer(5)*z**Integer(2) + Integer(9)*x**Integer(4)*y**Integer(3)*z**Integer(2) + Integer(5)*x**Integer(6)*y*z**Integer(2)
...              - Integer(3)*y**Integer(8)*z - Integer(9)*x**Integer(2)*y**Integer(6)*z - Integer(11)*x**Integer(4)*y**Integer(4)*z - Integer(7)*x**Integer(6)*y**Integer(2)*z
...              - Integer(2)*x**Integer(8)*z + y**Integer(9) + Integer(2)*x**Integer(2)*y**Integer(7) + Integer(3)*x**Integer(4)*y**Integer(5) + Integer(4)*x**Integer(6)*y**Integer(3) + Integer(2)*x**Integer(8)*y])
>>> Q = P([Integer(0),Integer(1),Integer(1)])
>>> C.tangents(Q)
[-y + z, 3*x^2 - y^2 + 2*y*z - z^2]
sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = P.curve([z^3*x + y^4 - x^2*z^2])
sage: Q = P([1,1,1])
sage: C.tangents(Q)
Traceback (most recent call last):
...
TypeError: (=(1 : 1 : 1)) is not a point on (=Projective Plane Curve
over Rational Field defined by y^4 - x^2*z^2 + x*z^3)
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = P.curve([z**Integer(3)*x + y**Integer(4) - x**Integer(2)*z**Integer(2)])
>>> Q = P([Integer(1),Integer(1),Integer(1)])
>>> C.tangents(Q)
Traceback (most recent call last):
...
TypeError: (=(1 : 1 : 1)) is not a point on (=Projective Plane Curve
over Rational Field defined by y^4 - x^2*z^2 + x*z^3)
class sage.schemes.curves.projective_curve.ProjectivePlaneCurve_field(A, f, category=None)[source]#

Bases: ProjectivePlaneCurve, ProjectiveCurve_field

Projective plane curves over fields.

arithmetic_genus()[source]#

Return the arithmetic genus of this projective curve.

This is the arithmetic genus \(p_a(C)\) as defined in [Har1977].

For an irreducible projective plane curve of degree \(d\), this is simply \((d - 1)(d - 2)/2\). It need not equal the geometric genus (the genus of the normalization of the curve).

EXAMPLES:

sage: x,y,z = PolynomialRing(GF(5), 3, 'xyz').gens()
sage: C = Curve(y^2*z^7 - x^9 - x*z^8); C
Projective Plane Curve over Finite Field of size 5
 defined by -x^9 + y^2*z^7 - x*z^8
sage: C.arithmetic_genus()
28
sage: C.genus()  # geometric
4
>>> from sage.all import *
>>> x,y,z = PolynomialRing(GF(Integer(5)), Integer(3), 'xyz').gens()
>>> C = Curve(y**Integer(2)*z**Integer(7) - x**Integer(9) - x*z**Integer(8)); C
Projective Plane Curve over Finite Field of size 5
 defined by -x^9 + y^2*z^7 - x*z^8
>>> C.arithmetic_genus()
28
>>> C.genus()  # geometric
4
sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve([y^3*x - x^2*y*z - 7*z^4])
sage: C.arithmetic_genus()
3
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([y**Integer(3)*x - x**Integer(2)*y*z - Integer(7)*z**Integer(4)])
>>> C.arithmetic_genus()
3
fundamental_group()[source]#

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

Note

The curve must be defined over the rationals or a number field with an embedding over \(\QQbar\).

Note

This functionality requires the sirocco package to be installed.

EXAMPLES:

sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = P.curve(x^2*z - y^3)
sage: C.fundamental_group()                                 # needs sirocco
Finitely presented group < x0 | x0^3 >
sage: P.curve(z*(x^2*z - y^3)).fundamental_group()          # needs sirocco
Finitely presented group < x0, x1 | x1*x0*x1*x0^-1*x1^-1*x0^-1 >
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = P.curve(x**Integer(2)*z - y**Integer(3))
>>> C.fundamental_group()                                 # needs sirocco
Finitely presented group < x0 | x0^3 >
>>> P.curve(z*(x**Integer(2)*z - y**Integer(3))).fundamental_group()          # needs sirocco
Finitely presented group < x0, x1 | x1*x0*x1*x0^-1*x1^-1*x0^-1 >

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

sage: # needs sage.rings.number_field
sage: a = QQ[x](x^2 + 5).roots(QQbar)[0][0]
sage: a
-2.236067977499790?*I
sage: F = NumberField(a.minpoly(), 'a', embedding=a)
sage: P.<x,y,z> = ProjectiveSpace(F, 2)
sage: F.inject_variables()
Defining a
sage: C = P.curve(x^2 + a * y^2)
sage: C.fundamental_group()                         # needs sirocco
Finitely presented group < x0 |  >
>>> from sage.all import *
>>> # needs sage.rings.number_field
>>> a = QQ[x](x**Integer(2) + Integer(5)).roots(QQbar)[Integer(0)][Integer(0)]
>>> a
-2.236067977499790?*I
>>> F = NumberField(a.minpoly(), 'a', embedding=a)
>>> P = ProjectiveSpace(F, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> F.inject_variables()
Defining a
>>> C = P.curve(x**Integer(2) + a * y**Integer(2))
>>> C.fundamental_group()                         # needs sirocco
Finitely presented group < x0 |  >
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{P}^{1}\) and this curve, given as a scheme morphism.

EXAMPLES:

sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve([y^2*z - x^3], P)
sage: C.rational_parameterization()
Scheme morphism:
  From: Projective Space of dimension 1 over Rational Field
  To:   Projective Plane Curve over Rational Field
        defined by -x^3 + y^2*z
  Defn: Defined on coordinates by sending (s : t) to
        (s^2*t : s^3 : t^3)
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([y**Integer(2)*z - x**Integer(3)], P)
>>> C.rational_parameterization()
Scheme morphism:
  From: Projective Space of dimension 1 over Rational Field
  To:   Projective Plane Curve over Rational Field
        defined by -x^3 + y^2*z
  Defn: Defined on coordinates by sending (s : t) to
        (s^2*t : s^3 : t^3)
sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve([x^3 - 4*y*z^2 + x*z^2 - x*y*z], P)
sage: C.rational_parameterization()
Scheme morphism:
  From: Projective Space of dimension 1 over Rational Field
  To:   Projective Plane Curve over Rational Field
        defined by x^3 - x*y*z + x*z^2 - 4*y*z^2
  Defn: Defined on coordinates by sending (s : t) to
        (4*s^2*t + s*t^2 : s^2*t + t^3 : 4*s^3 + s^2*t)
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([x**Integer(3) - Integer(4)*y*z**Integer(2) + x*z**Integer(2) - x*y*z], P)
>>> C.rational_parameterization()
Scheme morphism:
  From: Projective Space of dimension 1 over Rational Field
  To:   Projective Plane Curve over Rational Field
        defined by x^3 - x*y*z + x*z^2 - 4*y*z^2
  Defn: Defined on coordinates by sending (s : t) to
        (4*s^2*t + s*t^2 : s^2*t + t^3 : 4*s^3 + s^2*t)
sage: P.<x,y,z> = ProjectiveSpace(QQ, 2)
sage: C = Curve([x^2 + y^2 + z^2], P)
sage: C.rational_parameterization()
Scheme morphism:
  From: Projective Space of dimension 1 over Number Field in a
        with defining polynomial a^2 + 1
  To:   Projective Plane Curve over Number Field in a
        with defining polynomial a^2 + 1 defined by x^2 + y^2 + z^2
  Defn: Defined on coordinates by sending (s : t) to
        ((-a)*s^2 + (-a)*t^2 : s^2 - t^2 : 2*s*t)
>>> from sage.all import *
>>> P = ProjectiveSpace(QQ, Integer(2), names=('x', 'y', 'z',)); (x, y, z,) = P._first_ngens(3)
>>> C = Curve([x**Integer(2) + y**Integer(2) + z**Integer(2)], P)
>>> C.rational_parameterization()
Scheme morphism:
  From: Projective Space of dimension 1 over Number Field in a
        with defining polynomial a^2 + 1
  To:   Projective Plane Curve over Number Field in a
        with defining polynomial a^2 + 1 defined by x^2 + y^2 + z^2
  Defn: Defined on coordinates by sending (s : t) to
        ((-a)*s^2 + (-a)*t^2 : s^2 - t^2 : 2*s*t)
riemann_surface(**kwargs)[source]#

Return the complex Riemann surface determined by this curve

OUTPUT: A RiemannSurface object.

EXAMPLES:

sage: R.<x,y,z> = QQ[]
sage: C = Curve(x^3 + 3*y^3 + 5*z^3)
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, z']; (x, y, z,) = R._first_ngens(3)
>>> C = Curve(x**Integer(3) + Integer(3)*y**Integer(3) + Integer(5)*z**Integer(3))
>>> 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.projective_curve.ProjectivePlaneCurve_finite_field(A, f, category=None)[source]#

Bases: ProjectivePlaneCurve_field

Projective plane curves over finite fields

rational_points(algorithm='enum', sort=True)[source]#

Return the rational points on this curve.

INPUT:

  • algorithm – one of

    • 'enum' – straightforward enumeration

    • 'bn' – via Singular’s brnoeth package.

  • sort – boolean (default: True); whether the output points should be sorted. If False, the order of the output is non-deterministic.

OUTPUT: A list of all the rational points on the curve, possibly sorted.

Note

The Brill-Noether package does not always work (i.e., the ‘bn’ algorithm. When it fails a RuntimeError exception is raised.

EXAMPLES:

sage: x, y, z = PolynomialRing(GF(5), 3, 'xyz').gens()
sage: f = y^2*z^7 - x^9 - x*z^8
sage: C = Curve(f); C
Projective Plane Curve over Finite Field of size 5
 defined by -x^9 + y^2*z^7 - x*z^8
sage: C.rational_points()
[(0 : 0 : 1), (0 : 1 : 0), (2 : 2 : 1), (2 : 3 : 1),
 (3 : 1 : 1), (3 : 4 : 1)]
sage: C = Curve(x - y + z)
sage: C.rational_points()
[(0 : 1 : 1), (1 : 1 : 0), (1 : 2 : 1), (2 : 3 : 1),
 (3 : 4 : 1), (4 : 0 : 1)]
sage: C = Curve(x*z + z^2)
sage: C.rational_points('all')
[(0 : 1 : 0), (1 : 0 : 0), (1 : 1 : 0), (2 : 1 : 0),
 (3 : 1 : 0), (4 : 0 : 1), (4 : 1 : 0), (4 : 1 : 1),
 (4 : 2 : 1), (4 : 3 : 1), (4 : 4 : 1)]
>>> from sage.all import *
>>> x, y, z = PolynomialRing(GF(Integer(5)), Integer(3), 'xyz').gens()
>>> f = y**Integer(2)*z**Integer(7) - x**Integer(9) - x*z**Integer(8)
>>> C = Curve(f); C
Projective Plane Curve over Finite Field of size 5
 defined by -x^9 + y^2*z^7 - x*z^8
>>> C.rational_points()
[(0 : 0 : 1), (0 : 1 : 0), (2 : 2 : 1), (2 : 3 : 1),
 (3 : 1 : 1), (3 : 4 : 1)]
>>> C = Curve(x - y + z)
>>> C.rational_points()
[(0 : 1 : 1), (1 : 1 : 0), (1 : 2 : 1), (2 : 3 : 1),
 (3 : 4 : 1), (4 : 0 : 1)]
>>> C = Curve(x*z + z**Integer(2))
>>> C.rational_points('all')
[(0 : 1 : 0), (1 : 0 : 0), (1 : 1 : 0), (2 : 1 : 0),
 (3 : 1 : 0), (4 : 0 : 1), (4 : 1 : 0), (4 : 1 : 1),
 (4 : 2 : 1), (4 : 3 : 1), (4 : 4 : 1)]
sage: F = GF(7)
sage: P2.<X,Y,Z> = ProjectiveSpace(F, 2)
sage: C = Curve(X^3 + Y^3 - Z^3)
sage: C.rational_points()
[(0 : 1 : 1), (0 : 2 : 1), (0 : 4 : 1), (1 : 0 : 1), (2 : 0 : 1),
(3 : 1 : 0), (4 : 0 : 1), (5 : 1 : 0), (6 : 1 : 0)]
>>> from sage.all import *
>>> F = GF(Integer(7))
>>> P2 = ProjectiveSpace(F, Integer(2), names=('X', 'Y', 'Z',)); (X, Y, Z,) = P2._first_ngens(3)
>>> C = Curve(X**Integer(3) + Y**Integer(3) - Z**Integer(3))
>>> C.rational_points()
[(0 : 1 : 1), (0 : 2 : 1), (0 : 4 : 1), (1 : 0 : 1), (2 : 0 : 1),
(3 : 1 : 0), (4 : 0 : 1), (5 : 1 : 0), (6 : 1 : 0)]
sage: F = GF(1237)
sage: P2.<X,Y,Z> = ProjectiveSpace(F, 2)
sage: C = Curve(X^7 + 7*Y^6*Z + Z^4*X^2*Y*89)
sage: len(C.rational_points())
1237
>>> from sage.all import *
>>> F = GF(Integer(1237))
>>> P2 = ProjectiveSpace(F, Integer(2), names=('X', 'Y', 'Z',)); (X, Y, Z,) = P2._first_ngens(3)
>>> C = Curve(X**Integer(7) + Integer(7)*Y**Integer(6)*Z + Z**Integer(4)*X**Integer(2)*Y*Integer(89))
>>> len(C.rational_points())
1237
sage: # needs sage.rings.finite_rings
sage: F = GF(2^6,'a')
sage: P2.<X,Y,Z> = ProjectiveSpace(F, 2)
sage: C = Curve(X^5 + 11*X*Y*Z^3 + X^2*Y^3 - 13*Y^2*Z^3)
sage: len(C.rational_points())
104
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> F = GF(Integer(2)**Integer(6),'a')
>>> P2 = ProjectiveSpace(F, Integer(2), names=('X', 'Y', 'Z',)); (X, Y, Z,) = P2._first_ngens(3)
>>> C = Curve(X**Integer(5) + Integer(11)*X*Y*Z**Integer(3) + X**Integer(2)*Y**Integer(3) - Integer(13)*Y**Integer(2)*Z**Integer(3))
>>> len(C.rational_points())
104
sage: R.<x,y,z> = GF(2)[]
sage: f = x^3*y + y^3*z + x*z^3
sage: C = Curve(f); pts = C.rational_points()
sage: pts
[(0 : 0 : 1), (0 : 1 : 0), (1 : 0 : 0)]
>>> from sage.all import *
>>> R = GF(Integer(2))['x, y, z']; (x, y, z,) = R._first_ngens(3)
>>> f = x**Integer(3)*y + y**Integer(3)*z + x*z**Integer(3)
>>> C = Curve(f); pts = C.rational_points()
>>> pts
[(0 : 0 : 1), (0 : 1 : 0), (1 : 0 : 0)]
rational_points_iterator()[source]#

Return a generator object for the rational points on this curve.

INPUT:

  • self – a projective curve

OUTPUT:

A generator of all the rational points on the curve defined over its base field.

EXAMPLES:

sage: F = GF(37)
sage: P2.<X,Y,Z> = ProjectiveSpace(F, 2)
sage: C = Curve(X^7 + Y*X*Z^5*55 + Y^7*12)
sage: len(list(C.rational_points_iterator()))
37
>>> from sage.all import *
>>> F = GF(Integer(37))
>>> P2 = ProjectiveSpace(F, Integer(2), names=('X', 'Y', 'Z',)); (X, Y, Z,) = P2._first_ngens(3)
>>> C = Curve(X**Integer(7) + Y*X*Z**Integer(5)*Integer(55) + Y**Integer(7)*Integer(12))
>>> len(list(C.rational_points_iterator()))
37
sage: F = GF(2)
sage: P2.<X,Y,Z> = ProjectiveSpace(F, 2)
sage: C = Curve(X*Y*Z)
sage: a = C.rational_points_iterator()
sage: next(a)
(1 : 0 : 0)
sage: next(a)
(0 : 1 : 0)
sage: next(a)
(1 : 1 : 0)
sage: next(a)
(0 : 0 : 1)
sage: next(a)
(1 : 0 : 1)
sage: next(a)
(0 : 1 : 1)
sage: next(a)
Traceback (most recent call last):
...
StopIteration
>>> from sage.all import *
>>> F = GF(Integer(2))
>>> P2 = ProjectiveSpace(F, Integer(2), names=('X', 'Y', 'Z',)); (X, Y, Z,) = P2._first_ngens(3)
>>> C = Curve(X*Y*Z)
>>> a = C.rational_points_iterator()
>>> next(a)
(1 : 0 : 0)
>>> next(a)
(0 : 1 : 0)
>>> next(a)
(1 : 1 : 0)
>>> next(a)
(0 : 0 : 1)
>>> next(a)
(1 : 0 : 1)
>>> next(a)
(0 : 1 : 1)
>>> next(a)
Traceback (most recent call last):
...
StopIteration
sage: # needs sage.rings.finite_rings
sage: F = GF(3^2,'a')
sage: P2.<X,Y,Z> = ProjectiveSpace(F, 2)
sage: C = Curve(X^3 + 5*Y^2*Z - 33*X*Y*X)
sage: b = C.rational_points_iterator()
sage: next(b)
(0 : 1 : 0)
sage: next(b)
(0 : 0 : 1)
sage: next(b)
(2*a + 2 : a : 1)
sage: next(b)
(2 : a + 1 : 1)
sage: next(b)
(a + 1 : 2*a + 1 : 1)
sage: next(b)
(1 : 2 : 1)
sage: next(b)
(2*a + 2 : 2*a : 1)
sage: next(b)
(2 : 2*a + 2 : 1)
sage: next(b)
(a + 1 : a + 2 : 1)
sage: next(b)
(1 : 1 : 1)
sage: next(b)
Traceback (most recent call last):
...
StopIteration
>>> from sage.all import *
>>> # needs sage.rings.finite_rings
>>> F = GF(Integer(3)**Integer(2),'a')
>>> P2 = ProjectiveSpace(F, Integer(2), names=('X', 'Y', 'Z',)); (X, Y, Z,) = P2._first_ngens(3)
>>> C = Curve(X**Integer(3) + Integer(5)*Y**Integer(2)*Z - Integer(33)*X*Y*X)
>>> b = C.rational_points_iterator()
>>> next(b)
(0 : 1 : 0)
>>> next(b)
(0 : 0 : 1)
>>> next(b)
(2*a + 2 : a : 1)
>>> next(b)
(2 : a + 1 : 1)
>>> next(b)
(a + 1 : 2*a + 1 : 1)
>>> next(b)
(1 : 2 : 1)
>>> next(b)
(2*a + 2 : 2*a : 1)
>>> next(b)
(2 : 2*a + 2 : 1)
>>> next(b)
(a + 1 : a + 2 : 1)
>>> next(b)
(1 : 1 : 1)
>>> next(b)
Traceback (most recent call last):
...
StopIteration
riemann_roch_basis(D)[source]#

Return a basis for the Riemann-Roch space corresponding to \(D\).

This uses Singular’s Brill-Noether implementation.

INPUT:

  • D – a divisor

OUTPUT: A list of function field elements that form a basis of the Riemann-Roch space.

EXAMPLES:

sage: R.<x,y,z> = GF(2)[]
sage: f = x^3*y + y^3*z + x*z^3
sage: C = Curve(f); pts = C.rational_points()
sage: D = C.divisor([ (4, pts[0]), (4, pts[2]) ])
sage: C.riemann_roch_basis(D)
[x/y, 1, z/y, z^2/y^2, z/x, z^2/(x*y)]
>>> from sage.all import *
>>> R = GF(Integer(2))['x, y, z']; (x, y, z,) = R._first_ngens(3)
>>> f = x**Integer(3)*y + y**Integer(3)*z + x*z**Integer(3)
>>> C = Curve(f); pts = C.rational_points()
>>> D = C.divisor([ (Integer(4), pts[Integer(0)]), (Integer(4), pts[Integer(2)]) ])
>>> C.riemann_roch_basis(D)
[x/y, 1, z/y, z^2/y^2, z/x, z^2/(x*y)]
sage: R.<x,y,z> = GF(5)[]
sage: f = x^7 + y^7 + z^7
sage: C = Curve(f); pts = C.rational_points()
sage: D = C.divisor([ (3, pts[0]), (-1,pts[1]), (10, pts[5]) ])
sage: C.riemann_roch_basis(D)
[(-2*x + y)/(x + y), (-x + z)/(x + y)]
>>> from sage.all import *
>>> R = GF(Integer(5))['x, y, z']; (x, y, z,) = R._first_ngens(3)
>>> f = x**Integer(7) + y**Integer(7) + z**Integer(7)
>>> C = Curve(f); pts = C.rational_points()
>>> D = C.divisor([ (Integer(3), pts[Integer(0)]), (-Integer(1),pts[Integer(1)]), (Integer(10), pts[Integer(5)]) ])
>>> C.riemann_roch_basis(D)
[(-2*x + y)/(x + y), (-x + z)/(x + y)]

Note

Currently this only works over prime field and divisors supported on rational points.