Surfaces in P3

Contents

Introduction

This section describes several packages of functionality developed for working with (hyper)surfaces in three-dimensional projective space P3 which can have arbitrary singularities (at least in the characteristic 0 case). It was designed and implemented by Tobias Beck and has at its core, routines to compute a formal desingularization of a hypersurface X, expressed via a collection of algebraic power series giving the formal completion of the components of some desingularization lying over the components of the singular subscheme of the hypersurface. The algorithm is based on the method of Jung. It is fully described in [Bec07] and will be partly sketched out in the following two subsections.

As of magma V2.21, the method of resolution of singularities by local blow-up is also available in the case that X has only point singularities. This can often be faster than the formal desingularisation in practice, although the method is less elegant, and provides additional information (see section Desingularisation by Blow Up). It is now the default method (given the point singularity restriction) when the general ResolveSingularSurface intrinsic is used. If desired, formal desingularisation can be used instead by calling that intrinsic with the UseFormalDesing parameter set to true, or the intrinsic FormallyResolveProjectiveHypersurface described below can be called. It should also be noted that formal resolution requires that the base field of the surface is characteristic 0. This is not required for blow-up resolution.

Major functionality to be described in this section, which is enabled by our resolutions, includes the computation of important birational invariants of any desingularization of X like the arithmetic and geometric genera and higher geometric plurigenera. Furthermore, we also provide an intrinsic for the computation of (i,j)-adjoint maps as rational maps on X. For formal desingularisation, the theoretical and algorithmic details may be found in [BS08]. Blow-up resolutions allow similar computations and the intrinsics for all of the above now work with either type of resolution data.

There are functions to determine whether X is of Kodaira dimension -∞, i.e., birationally ruled. For the important special case of rational surfaces, there is a suite of functions to determine whether a parametrization exists over the base field and to explicitly construct one in the affirmative case. This is based on the work of Josef Schicho described in [Sch98] and [Sch00]. As with the rest of the package, it was originally coded to make use of formal desingularisation data but now also works with blow-up resolution data.

A surface X is mapped to a standard model by applying an appropriate m-adjoint map. These are then parameterized by special case code. The main special cases are Del Pezzo surfaces (including some singular cases) and line and conic bundles. The functions can be called directly by the user. Apart from the previously existing Del Pezzo code (for degrees 6, 8 and 9) and the special singular code for degrees 3 and 4, these functions were implemented by Tobias Beck and Josef Schicho.

Embedded Formal Desingularization of Curves

A formal embedded desingularization of plane curves, as described below, is used in the Jung surface resolution process. The main function is available to the user and provides another alternative to the existing function field and resolution graph curve functionality.

Before describing the function, we introduce some terminology. Let C ⊂P be a plane algebraic curve (where P = AE2 or P = PE2 for some field E of characteristic zero) and π: Q to P an embedded desingularization, i.e., π is proper birational, Q is regular and D := π - 1(C) ⊂Q is a normal crossing divisor. Further let {p1, ..., pr} ∈Q be the generic points of the decomposition of D into irreducible components and {q1, ..., qs} ∈Q the closed points of the normal crossings of D.

From π we can construct morphisms Spec widehat(OOQ, pi) to P and Spec widehat(OOQ, qi) to P. The set of all these morphisms (up to isomorphism of the domain) is called a formal embedded desingularization of C ⊂P. Each of these morphisms has a centre on P which is defined to be the image of the closed point.

The two classes of morphisms are represented, respectively, by homomorphisms A to widehat(OOQ, pi) and A to widehat(OOQ, qi) (where A is either the normal polynomial ring E[x, y] or the graded polynomial ring E[x, y, w] and the inverse image of the maximal ideal of the completion ring is the prime ideal defining the centre), and we are free to choose an isomorphic representation of the codomain. We refer to the homomorphisms as μi and νi respectively.

ResolveAffineCurve(p) : RngMPolElt -> List, List, List, RngIntElt
    Factors: SeqEnum                    Default: []
    Ps: RngMPolElt                      Default: 0
    Focus: RngMPolElt                   Default: 0
    ExtName: MonStgElt                  Default: "alpha"
    ExtCount: RngIntElt                 Default: 0
    SetVerbose("Resolve", n):           Maximum: 1

Given the curve defined by p ∈E[x, y] (a non-zero bivariate polynomial over a number field), this intrinsic essentially computes a formal embedded resolution of the curve using a succession of point blow ups. Only morphisms whose centres vanish on the ideal generated by Focus are considered. Note that Focus may be a single polynomial or a sequence of polynomials.

The three returned lists contain elements of the form (b1, (y, m11), (p1, m12)), (b2, (y, m2)) and (b2, (p3, m3)) respectively. Here b1, b2 and b3 are homomorphisms E[x, y] to E'[x, y] to some bivariate polynomial ring over an algebraic field extension E' over E.

The first list gathers normal crossings. Precisely, the extended homomorphism b1: E[x, y] to E'll x, y rr corresponds to a νi from above. Moreover we have < b1(p) > = < ym11p1m12 > where y = 0 and p1 = 0 have a normal crossing.

The second list gathers exceptional divisors. The extended homomorphism b2: E[x, y] to E'(x)ll y rr corresponds to a μi from above. Moreover we have < b2(p) > = < ym2 > and y=0 corresponds to an exceptional divisor.

Finally, the last list corresponds to the components of the original curve. The extended homomorphism b3: E[x, y] to widehat(E'[x, y]< p3 >) (where E'=E in this case) corresponds to another μi from above. Moreover we have that b3(p) has multiplicity m3 in widehat(E'[x, y]< p3 >) and corresponds to an original curve component.

If known, a factorization of p (as returned by the Factorization command) can be passed using the parameter Factors and the squarefree part of p (as returned by the SquarefreePart command) using Ps.

If the ground field has to be extended, the algebraic elements will be displayed as ExtName_i where i starts from ExtCount. The last return value is the value of ExtCount plus the number of field extensions that have been introduced, which can be useful for consecutive naming when making a series of resolution calls.

Example AlgSrf_aff_crv_res (H123E19)

We compute an embedded resolution of an affine plane curve.
> Q := Rationals();
> Qxy<x,y> := PolynomialRing(Q, 2, "glex");
> f := (y^2-x^3)*(x^2-y^2-y^3);
> NCs, EXs, DCs := ResolveAffineCurve(f : Factors := Factorization(f));
> #NCs, #EXs, #DCs;
7 4 2
> NCs[2]; EXs[3]; DCs[1];
[*
    Mapping from: RngMPol: Qxy to RngMPol: Qxy,
    <y, 4>,
    <-x^2 + 2*x + y, 1>
*]
[*
    Mapping from: RngMPol: Qxy to RngMPol: Qxy,
    <y, 10>
*]
[*
    Mapping from: RngMPol: Qxy to RngMPol: Qxy,
    <y^3 - x^2 + y^2, 1>
*]
> NCs[2][1](x), NCs[2][1](y);
x*y - y
y

Here we have passed the factorization of f only for illustrative purposes. The curve is the union of a cusp and a node at the origin. It has two singular points over Q, the origin and another intersection point of the two curves which has a residue field of degree 5 over Q.

We have computed the local information of a (not necessarily minimal) embedded resolution and find that it contains the 2 components of the strict transform, further 4 exceptional divisors and 7 normal crossings. For example, the pushforward of f under the chart map x |-> xy - y, y |-> y is equal to y4( - x2 + 2x + y) up to a local unit. The corresponding germ is isomorphic to a normal crossing in the embedded desingularization. We also see that one of the exceptional divisors has multiplicity 10.

If we were only interested in a local resolution, we would do the following:

> NCs, EXs, DCs := ResolveAffineCurve(f : Focus := [x,y]);
> #NCs, #EXs, #DCs;
5 3 0

We focus on the origin, hence, any curve components are not considered . We have 1 less exceptional divisor and 2 less normal crossings. This is because the second intersection point of the above two curve components was already a normal crossing, but our algorithm has nevertheless blown it up in the previous example.

ResolveProjectiveCurve(p) : RngMPolElt -> List, List, List, RngIntElt
    Focus: RngMPolElt                   Default: 0
    ExtName: MonStgElt                  Default: "alpha"
    ExtCount: RngIntElt                 Default: 0
    SetVerbose("Resolve", n):           Maximum: 1
Given the curve defined by p ∈E[x, y, z] (a non-zero trivariate homogeneous polynomial over a number field), this intrinsic essentially computes a formal embedded resolution of the curve using a succession of point blow ups. This is the same as ResolveAffineCurve above, but now p is a homogeneous polynomial in three variables that defines a projective curve. Accordingly, the bj map from the respective homogeneous coordinate ring to some E'[x, y].

Example AlgSrf_prj_crv_res (H123E20)

We can also desingularize the projectivisation of the above curve.
> Q := RationalField();
> QXYZ<X,Y,Z> := PolynomialRing(Q, 3);
> F := (Y^2*Z-X^3)*(X^2*Z-Y^2*Z-Y^3);
> NCs, EXs, DCs := ResolveProjectiveCurve(F); #NCs, #EXs, #DCs;
7 4 2
> NCs[3];
[*
    Mapping from: RngMPol: QXYZ to Polynomial ring of rank 2 over
    Rational Field ...,
    <y, 4>,
    <x^2 + 2*x - y, 1>
*]
> NCs[3][1](X);
x*y + y
> NCs[3][1](Y);
y
> NCs[3][1](Z);
1

The homomorphisms take a slightly different shape (because they have now Q[X, Y, Z] as domain), but otherwise they are the same. This is because the curve has no singularities at infinity.

Formal Desingularization of Surfaces

For curves we have described embedded formal desingularization. For surfaces instead we produce only formal desingularizations. Let S ⊂P be a hypersurface (where P = AE3 or P = PE3) and C ⊂S a closed subset (which typically contains the singular locus). Further let π: T to S be a desingularization, i.e., π is proper birational and T is regular. By {p1, ..., pr} ∈T we denote the generic points of the curve components of the decomposition of D := π - 1(C) into irreducibles.

From π we can construct morphisms Spec widehat(OOT, pi) to S. The set of all these morphisms (up to isomorphism of the domain) is called a formal desingularization of S over C ⊂S. Such a morphism has a centre on S which is defined as the image of the closed point (and actually is contained in C).

The morphisms are represented by homomorphisms A to widehat(OOT, pi) (where A is either the algebra E[x, y, z]/< p > or the graded algebra E[x, y, z, w]/< p > with p a defining polynomial), and we are free to choose an isomorphic representation of the codomain. We refer to such a homomorphisms as μi.

In the actual algorithm, C is the ramification locus of a finite projection, pr, to an affine or projective plane (C contains the singular subscheme of S). The underlying desingularization T (which is not computed explicitly) is a Jung resolution which is constructed in two stages. Firstly, an embedded resolution of the image of C in the plane is performed by blow-ups and T1 is taken as the normalization of the pullback of this by pr. So T1 then has only point singularities of a simple type (toric singularities), lying over the (normal-crossing) intersections of components of the embedded resolution. These are resolved by a finite succession of blow-ups on T1 to give T.

The algorithm computes the formal desingularization, as described above, corresponding to T, using the embedded formal desingularization for curves followed by algebraic power series operations for the normalization and final resolution of the toric singularities. This is described fully in [Bec07]. The μi homomorphisms are defined by algebraic power series images of the variables of P.

It is important to note that the Jung desingularization T is not a minimal desingularization and, in any case, the set of morphisms returned for the formal desingularization generally contain some elements whose centre on S is already non-singular (because, for example, components of C are often generically non-singular). However, there is a parameter option with the main function FormallyResolveProjectiveHyperSurface, which removes "non-singular" morphisms and possibly others that have no effect on the computation of birational invariants and m-adjoint maps.

ResolveAffineMonicSurface(s) : RngUPolElt -> List, RngIntElt
    Focus: RngMPolElt                   Default: 0
    ExtName: MonStgElt                  Default: "alpha"
    ExtCount: RngIntElt                 Default: 0
    SetVerbose("Resolve", n):           Maximum: 1
The main user resolution function FormallyResolveProjectiveHyperSurface is for projective hypersurfaces. This affine version, however, may be useful in some circumstances. The input is a monic, squarefree polynomial s ∈E[x, y][z] where E is a number field (i.e., s is univariate over a bivariate polynomial ring). Let S ⊂AE3 denote the surface defined by it and C ⊂S the closed subset defined by discz(s) (i.e., the intersection of S with the cylinder over the discriminant curve when considering the projection S to AE2 in z-direction). The function computes a formal desingularization of S over C (see above).

The first return value is a list of elements of the form ((X, Y, Z), o) where X, Y, Z ∈F ll t rr are univariate power series (over some field extension F of transcendence degree 1 over E) s.t. s(X, Y, Z)=0 and o is an integer. The induced homomorphism E[x, y][z]/(s) to F ll t rr corresponds to a μi from above and o is its adjoint order, i.e., the negation of the order of a special differential form (see Section Adjoint Systems and Birational Invariants).

One can specify a focus ideal FF ⊂E[x, y] by passing a single generator or sequence of generators in Focus (as for ResolveAffineCurve). In this case C is taken to be the intersection of S and the cylinder over the zero set of FF + < discz(s) >.

If the ground field has to be extended, the algebraic elements will be displayed as ExtName_i where i starts from ExtCount. The last return value is the value of ExtCount plus the number of field extensions that have been introduced, which can be useful for consecutive naming when making a series of resolution calls. A transcendental element will always be displayed as s.

Example AlgSrf_aff_res (H123E21)

We compute a formal desingularization for the affine surface z2 - xy = 0.
> Q := Rationals();
> Qxy<x,y> := PolynomialRing(Q, 2, "glex");
> Qxyz<z> := PolynomialRing(Qxy);
> f := z^2 - x*y;
> desing := ResolveAffineMonicSurface(f); #desing;
3

We have computed 3 morphisms. Two of them are centred over the coordinate axes x=0 and y=0. But they might not be of interest, because the surface is normal and has an isolated singularity over the origin.

> #ResolveAffineMonicSurface(f : Focus := [x,y]);
1

The only remaining morphism corresponds to the exceptional divisor obtained by blowing up the singularity.

Elements in the returned list which define the morphisms of the formal desingularization are examined more closely in the projective surface example below.

FormallyResolveProjectiveHyperSurface(S): Srfc -> List, RngIntElt
    AdjComp: BoolElt                    Default: false
    ExtName: MonStgElt                  Default: "gamma"
    ExtCount: RngIntElt                 Default: 0
    SetVerbose("Resolve", n):           Maximum: 1
The main intrinsic is now ResolveSingularSurface for resolution of singularities on surfaces and this can be forced to produce a formal resolution (if possible) by calling it with the parameter UseFormalDesing set to true. That is essentially equivalent to calling this intrinsic directly.

This gives a direct option to get the formal resolution and allows the use of the ExtName and ExtCount parameters.

It is similar in description to ResolveAffineMonicSurface above. The argument is a projective surface S ⊂PE3 over a field of characteristic 0. Computes and returns a formal desingularization (see above) of S. It will be a formal desingularization over an automatically chosen subset C ⊂S (using again the cylinder over the discriminant curve w.r.t. a nice projection onto some PE2). Accordingly the elements of the return list of formal desingularization data are now of the form ((X, Y, Z, W), o).

If AdjComp is true, then only a sublist is returned that is still sufficient for the computation of birational invariants and adjoint spaces (see Section Adjoint Systems and Birational Invariants). The parameters ExtName and ExtCount and the second return value have the same meaning as in the affine case.

As stated above, the algorithm is based on formally computing a Jung resolution and is described in [Bec07].

The returned resolution list is also cached on S in a surface attribute (alongside any blow-up desingularisation data that has been computed). Thus, it will not be recomputed on a further call to this intrinsic or on calls to any other intrinsics using formal desingularisation data. The only time it may be recomputed is when this (or ResolveSingularSurface) is first called with AdjComp true and then called again with AdjComp false. S can cache both an `adjcomp only' formal desingularisation list and a general one and the cached data is labelled correspondingly.

Example AlgSrf_prj_res (H123E22)

Computing a formal desingularization is easy. We also see that setting AdjComp to true produces fewer formal divisors needed to only compute birational invariants or adjoint spaces.
> P<x,y,z,w> := ProjectiveSpace(Rationals(), 3);
> F := w^3*y^2*z+(x*z+w^2)^3;
> S := Surface(P, F);
> #FormallyResolveProjectiveHypersurface(S : AdjComp := true);
18
> desing := FormallyResolveProjectiveHypersurface(S); #desing;
26
> // could have used ResolveSingularSurface(S : UseFormalDesing := true);

Hence, the formal desingularization of the projective surface defined by F contains 26 morphisms. They are represented by tuples of power series that vanish on F. We have a closer look at the first morphism.

> prm, ord := Explode(desing[1]);
> IsZero(AlgComb(F, prm)); ord;
true
4
> X, Y, Z, W := Explode(prm);
> Expand(X, 6); Expand(Y, 6); Expand(Z, 6); Expand(W, 6);
true 1
true -s*t^2
true -t^2
true -1/64*s^2*gamma_0*t^5 + 1/2*gamma_0^2*t^3 + gamma_0*t^2 + t
> Domain(W);
Polynomial ring of rank 1 over Algebraic function field defined
over Univariate rational function field over Rational Field
by $.1^3 - 1/8*s^2
Graded Lexicographical Order
Variables: t

One of the morphisms is of type Spec Q(s)[gamma0]ll t rr to Proj Q[x, y, z, w] / < F > where gamma03 - 1/8 s2 = 0. In particular, Q(s)[gamma0] is isomorphic to the residue field of the corresponding prime divisor on the desingularization. From this one can for example deduce that it is a rational curve. The morphism is given by the ring homomorphism x |-> 1, y |-> - st2, z |-> - t2 and w |-> t + gamma0t2 + 1/2gamma02t3 - 1/64s2gamma0t5 + ... .

The adjoint order for this morphism is 4. Consider the chart x != 0. The special differential form (see Section Adjoint Systems and Birational Invariants) in this chart obtained by dehomogenizing is

frac(x5)((∂F / ∂w)(x, y, z, w)) d y/x ^ d z/x.

Substituting the values X, Y, Z and W we see that it is mapped to

frac(X5)((∂F / ∂w)(X, Y, Z, W)) d Y/X ^ d Z/X

= frac(1)((∂F / ∂w)(X, Y, Z, W)) d ( - st2) ^ d ( - t2)

= frac(1)((∂F / ∂w)(X, Y, Z, W)) (2st d t + t2 d s) ^ 2t d t

= frac(1)((∂F / ∂w)(X, Y, Z, W)) 2t3 d s ^ d t

The adjoint order is minus the overall order of this expression, hence, -3 plus the order of (∂F / ∂w)(X, Y, Z, W). We check the computation.

> Order(AlgComb(Derivative(F,w), prm));
7

Adjoint Systems and Birational Invariants

In this section we describe computation of adjoint spaces. Let S ⊂PE3 be a surface defined by a homogeneous irreducible polynomial F ∈E[x0, x1, x2, x3] of degree d and ΩE(S) | E the vector space of rational differential forms of the function field (over the ground field E of characteristic zero). We can consider ΩE(S) | E a constant sheaf of OOS-modules. Let Ui ⊂S be the affine open subsets of the standard covering w.r.t. this choice of variables.

Let ωS0 ⊂ΩE(S) | E^ 2 be the subsheaf which is locally generated on Ui by

(frac(∂F/ ∂xj)(xid - 1)) - 1 bigwedge_(k ∈{0, ..., 3} - {i, j}) d frac(xk)(xi)

(for an arbitrary choice of j != i). By sending this generator to xid - 4 one finds that ωS0 isomorphic to OOS(d - 4). Further let FFS, m ⊂(ΩE(S) | E^ 2) tensor m be the subsheaf of those forms whose pullbacks are regular on some desingularization of S. It is called the sheaf of m-adjoints. It is in fact well-defined, i.e., doesn't depend on any specific desingularization, when m ≥0, and one can show FFS, m ⊆(ωS0) tensor m. For more details we refer to [BS08].

Now since FFS, m is a coherent sheaf on the projective scheme S ⊂PE3 it can be defined by its associated graded module MS, m and by the above discussion FFS, m is isomorphic to a subsheaf of OOS(m(d - 4)). The module MS, m is thus naturally a graded submodule of (E[x0, x1, x2, x3]/< F >)(m(d - 4)). The n-th graded piece of MS, m, a linear subsystem of the standard linear system of degree n + m(d - 4) homogeneous polynomials on S, corresponds to global sections of the Serre twist FFS, m(n). This, under pullback, corresponds to the space of global sections of the twisted m-adjoint sheaf (ωX) tensor m(n) for any desingularization X of S, where (n) now signifies twisting by the n-th tensor power of LL, the invertible sheaf on X which gives the map into projective space projecting X down onto S.

These adjoint linear systems immediately give the plurigenera of any desingularization X as well as an explicit representation of the important twisted m-adjoint maps into projective space as rational maps from S (defined by the sequence of homogeneous polynomials forming a basis of the adjoint system). These maps are used to take any rational hypersurface to a standard model, as described in the next section. The untwisted m-adjoint map is just the the m-pluri-canonical map when m > 0 and the (-m)-pluri-anticanonical map when m < 0.

Note that the computations are all valid for m < 0 but in that case, the adjoint and twisted adjoint linear systems do depend on the desingularisation, so the results are not birationally invariant but are specific to the particular formal or blow-up desingularisation that has been computed.

The intrinsics below have all been updated so that they can either use formal desingularisation data or the newer blow-up desingularisation data. Once a singularity resolution has been computed, the desingularisation data is cached with the hypersurface X and so will not need to be recomputed. Furthermore, as of Magma V2.25, most of them now work for surfaces in ordinary projective ambients of dimension ≥4 with blow-up desingularisation data.

HomAdjoints(m,n,S) : RngIntElt, RngIntElt, Srfc -> SeqEnum
    UseFormalDesing: BoolElt            Default: false
    SetVerbose("Classify", n):          Maximum: 1
Given a surface S of degree d in P3 defined over a number field E together with integers m and n, the intrinsic returns a basis for the vector space of the degree-n graded summand of the graded ring associated to FFS, m (i.e., Γ(S, OOS(n) tensor FFS, m)) as a subspace of the homogeneous forms in E[x0, x1, x2, x3] (the coordinate ring of the P3 ambient) of degree n + m(d - 4) (see above).

For S in Pr, r ≥4, with only point singularities this will also work, although the common degree d of the homogeneous polynomials in the sequence returned has no simple description. The sequence still defines the (m, n) adjoint-map on S and its length is equal to the dimension of global sections H0(X, ωX tensor m(n)) where X is the desingularisation of S.

The intrinsic needs desingularisation data to work with and can use either formal desingularisation data or blow-up desingularisation data. It will first try to use desingularisation data (of either type) cached on the surface. If there is no cached precomputed data, it will call ResolveSingularSurface at the beginning to get such data (which will then have been cached with X). Use/generation of blow-up desingularisation data is now the default when S has only point singularities. To force the use of formal desingularisation data, the UseFormalDesing parameter should be set to true.

When using formal data, the intrinsic computes the adjoint space as a linear subspace of homogeneous polynomials of the appropriate degree by using the formal divisor morphisms of the formal desingularization to give additional linear conditions at the singular places of S. This is explained fully in [BS08].

When using blow-up data, the same linear conditions are computed with the LinearSystemDivisorRestriction intrinsic and either the differential multiplicities intrinsic for hypersurfaces or the divisor of the base meromorphic 2-form computed in BasisOfHolomorphicTwoForms for surfaces in Pr, r ≥4. These are all described in section Desingularisation by Blow Up.

GeometricGenusOfDesingularization(S) : Srfc -> RngIntElt
    UseFormalDesing: BoolElt            Default: false
Given a hypersurface S in P3 or surface S with only point singularities in Pn, n ≥4, the intrinsic returns the geometric genus of (any) desingularization of S. The function just computes the dimension of the (1, 0) adjoint space. The parameter is the same as for HomAdjoints.
PlurigenusOfDesingularization(S,m) : Srfc, RngIntElt -> RngIntElt
    UseFormalDesing: BoolElt            Default: false
Given a hypersurface S in P3 or surface S with only point singularities in Pn, n ≥4, the intrinsic returns the m-th plurigenus of (any) desingularization, X, of S. This is the dimension of the global sections of the sheaf (ωX) tensor m and is just computed as the dimension of the (m, 0) adjoint space.The parameter is the same as for HomAdjoints.
ArithmeticGenusOfDesingularization(S) : Srfc -> RngIntElt
    UseFormalDesing: BoolElt            Default: false
Given a hypersurface S in P3, the intrinsic returns the arithmetic genus of (any) desingularization of S. This is computed from a simple formula coming from the Riemann-Roch theorem, involving the dimension of the (1,1)-adjoints and the degree of S. The parameter is the same as for HomAdjoints.

Example AlgSrf_adj_ex (H123E23)

We compute several adjoint spaces of a surface. We compute the resolution at the start, after which it is cached for the HomAdjoints calls. The resolution computes a formal desingularisation since the surface has a singular subscheme of dimension 1.
> P<x,y,z,w> := ProjectiveSpace(Rationals(), 3);
> F := w^3*y^2*z+(x*z+w^2)^3;
> S := Surface(P,F);
> desing := ResolveSingularSurface(S : AdjComp := true);
> #desing;
18
> HomAdjoints(1, 0, S);
[]
> HomAdjoints(1, 1, S);
[
    x*z*w + w^3
]
> HomAdjoints(1, 2, S);
[
    x^2*z^2 - w^4, x^2*z*w + x*w^3, x*y*z*w, x*z^2*w + z*w^3,
    x*z*w^2 + w^4, y*z*w^2, y*w^3
]
> HomAdjoints(1, 3, S);
[
    x^3*z^2 - x*w^4, x^2*y*z^2, x^2*z^3 - z*w^4,
    x^3*z*w + x^2*w^3, x^2*y*z*w, x*y^2*z*w, x^2*z^2*w - w^5,
    x*y*z^2*w, x*z^3*w + z^2*w^3, x^2*z*w^2 + x*w^4, x*y*z*w^2,
    y^2*z*w^2, x*z^2*w^2 + z*w^4, y*z^2*w^2, x*y*w^3, y^2*w^3,
    x*z*w^3 + w^5, y*z*w^3, y*w^4
]
>
> HomAdjoints(2, 0, S);
[]
> HomAdjoints(2, 1, S);
[]
> HomAdjoints(2, 2, S);
[
    x^2*z^2*w^2 + 2*x*z*w^4 + w^6
]
> HomAdjoints(2, 3, S);
[
    x^3*z^2*w^2 + 2*x^2*z*w^4 + x*w^6, x^2*y*z^2*w^2 - y*w^6,
    x^2*z^3*w^2 + 2*x*z^2*w^4 + z*w^6,
    x^2*z^2*w^3 + 2*x*z*w^5 + w^7, x*y*z^2*w^3 + y*z*w^5,
    x*y*z*w^4 + y*w^6, y^2*z*w^4
]

Classification and Parameterization of Rational Surfaces

This section contains intrinsics for the recognition of rational surfaces in Prj3, the classification and transformation to standard models using appropriate m-adjunction maps and, finally, special case code for these standard models, to determine a parametrization of the original hypersurface.

The intrinsic to test whether a surface is rational actually also applies to surfaces in higher-dimensional ordinary projective spaces as long as the surface model has no worse than simple singularities.

A non-singular surface in Prjr with r ≥4 may be transformed to a standard model using the intrinsic MinimalModelRationalSurface (see Section Minimal Models).

IsRational(X) : Srfc -> BoolElt
    UseFormalDesing: BoolElt            Default: false
    KnownADE: BoolElt                   Default: false
Returns true if the ordinary projective surface X is (geometrically) rational, i.e. birationally isomorphic to the projective plane over the algebraic closure of its base field. This simply uses the Castelnuevo criterion that X is rational if and only if both the arithmetic genus and second plurigenus of any desingularization are zero.

If the ambient of X is P3, there is no assumption about the singularity of X and a formal or blow-up desingularisation will be used to compute plurigenera, as long as either the characteristic of the base field is 0 or the singular subscheme is of dimension 0. If the ambient is Pr, r ≥4, then the intrinsic can use a blow-up desingularisation if X has singular subscheme of dimension 0.

If X has at worst simple (A-D-E) singularities, the generally faster sheaf algorithms of the first section may be used for the computations. The user can allow this method to be applied by setting the parameter KnownADE to true if they know this to be the case. Otherwise, simple singularities are not checked for and the sheaf method is not tried.

If a desingularisation first needs to be computed (if one is not already cached with X), as usual, blow-up desingularisation is the default (when it can be computed). To force the use of formal desingularisation in the hypersurface case, parameter UseFormalDesing should be set to true.

Reduction to Special Models

In this section, we describe the function for the birational transformation over the base field of a rational hypersurface in Prj3 into a special model of one of the types in the standard classification as listed by Josef Schicho in [Sch98]. He enumerates 5 basic cases [Sch98, p. 17 and Lem. 5.2-5.7] and splits the last case in two, choosing labels "1", "2", "3", "4", "5A" and "5B". The lemmas describing cases 3 and 5A involve a further case distinction. Also, a label "0" is useful for the non-rational case. From this, we get the label set

LL := {"0", "1", "2", "3a", "3b", "4", "5Aa", "5Ab", "5Ac", "5B"}

Let S be a surface in Prj3. In each of the above cases the author specifies a set of adjoint spaces defining interesting maps, either birationally to a special surface or to a rational normal curve (giving a pencil of rational curves on the surface). More precisely, the maps can be computed using Vn, m := HomAdjoints (m, n, S) for certain choices of m and n. Let μ to be the smallest integer s.t. V1, μ + 1 != [ ]. Then the important Vn, m for the different cases are as follows:

vbox{tabskip=0pt offinterlineskip halign {#& vrule#tabskip=1em plus2em& hfil#& vrule#& hfil#hfil& vrule#& hfil#& vrule#& hfil#hfil& vrule#& hfil#& vrule#& hfil#& vrule#tabskip=0pt
tablerule &&omit hidewidth 0hidewidth&&omit hidewidth 1hidewidth&& omit hidewidth 2hidewidth&&omit hidewidth 3a hidewidth&& omit hidewidth 3b hidewidth&&omit hidewidth 4hidewidth&
tablerule &&[ ]&&[V1, μ]&&[V1, μ]&&[V2, 2μ + 1, V1, μ]&& [V2, 2μ + 1]&&[V1, μ, V2, 2μ - 1]&
tablerule}}

vbox{tabskip=0pt offinterlineskip halign {#& vrule#tabskip=1em plus2em& hfil#& vrule#& hfil#hfil& vrule#& hfil#& vrule#& hfil#& vrule#tabskip=0pt
tablerule &&omit hidewidth 5Aa hidewidth&&omit hidewidth 5Ab hidewidth&& omit hidewidth 5Ac hidewidth&&omit hidewidth 5B hidewidth&
tablerule &&[V1, μ - 1]&&[V1, μ - 1, V2, 2μ - 2]&& [V1, μ - 1, V2, 2μ - 2, V3, 3μ - 3]&&[V2, 1]&
tablerule}}

The function to find a parameterization of a rational hypersurface, which uses the reduction function below as a first stage, is described in the next section.

ClassifyRationalSurface(S) : Srfc -> Srfc, List, MonStgElt
    UseFormalDesing: BoolElt            Default: false
    SetVerbose("Classify", n):          Maximum: 1
Given an ordinary projective surface S in Prj3 over a number field, the intrinsic returns the special rational surface type to which S is birationally equivalent and associated data. If S is not (geometrically) rational, the return values are S itself, a list containing only the identity map on S and the string "Not rational".

If S is rational and Schicho's algorithm reduces it to a standard surface Y over the base field k, Y is the first return value. The second return value is a list of one or two scheme maps. The first is always a birational map from S to Y. There is a second map if and only if Y is a rational scroll or conic bundle. Then S (and Y) have fibration maps to a rational normal curve such that the general fibre is a rational curve (and a line or conic for Y). The fibration map on S is the second return value. Note that, if the base field is Q, in these cases, S can be parameterized by calling ParametrizePencil with the fibration map as argument.

The third return value is a string describing the type of Y. It is "P2" (for the projective plane!), "Quadric surface" (for a degree 2 surface in Prj3), "Rational scroll", "Conic bundle" or "Del Pezzo of degree d" where 1 ≤d ≤9. The Del Pezzos might be degenerate (with simple singularities) and are anticanonically embedded in Prjd except for degrees 1 and 2 when they have their standard weighted-projective embeddings.

If formal or blow-up desingularisation data has not already been cached with S, it is computed at the start. Parameter UseFormalDesing is as for earlier intrinsics.

Example AlgSrf_class (H123E24)

Here are a few examples.
> P<x,y,z,w> := ProjectiveSpace(Rationals(),3);

The first surface:

> p1 := x^4 + y^4 - z^2*w^2;
> _,_,typ := ClassifyRationalSurface(Surface(P,p1));
> typ;
Not rational

The second surface:

> p2 := 2*x + y + 8*z + 5*w;
> Y,_,typ := ClassifyRationalSurface(Surface(P,p2));
> typ; Y;
P^2
Surface over Rational Field defined by
0

The third surface:

> p3 := x^2 - 4*x*z + 3*x*w + y*z - y*w + 2*z^2 - 3*z*w + w^2;
> _,_,typ := ClassifyRationalSurface(Surface(P,p3));
> typ;
Quadric surface

The fourth surface:

> p4 := (y^2 - w*z)*(w^2 - y*x) + (x*z - y*w)^2;
> S := Surface(P,p4);
> Y,mps,typ := ClassifyRationalSurface(S);
> typ;
Rational scroll
> mps[2]; // the fibration map
Mapping from: Srfc: S to Scheme over Rational Field defined by
$.1*$.2 - $.3^2
with equations :
x*y - w^2
y^2 - z*w
x*z - y*w

The fifth surface:

> p5 := x^3*y - 4*x^3*z - 6*x^3*w - 3*x^2*y^2 - 2*x^2*y*z
>     - 3*x^2*y*w + 50*x^2*z^2 + 146*x^2*z*w + 108*x^2*w^2
>     - 11*x*y^2*z + 2*x*y^2*w + 61*x*y*z^2 + 149*x*y*z*w
>     + 65*x*y*w^2 + 68*x*z^3 + 228*x*z^2*w + 260*x*z*w^2
>     + 112*x*w^3 + 4*y^4 - 13*y^3*z - 19*y^3*w + 20*y^2*z^2
>     + 77*y^2*z*w + 55*y^2*w^2 + 40*y*z^3 + 106*y*z^2*w
>     + 58*y*z*w^2 - 2*y*w^3 + 22*z^4 + 84*z^3*w + 130*z^2*w^2
>     + 108*z*w^3 + 38*w^4;
> S := Surface(P,p5);
> _,mps,typ := ClassifyRationalSurface(S);
> typ;
P2
> mps[1]; //birational map from Y to P2
Mapping from: Srfc: S to Surface over Rational Field defined by
0
with equations :
x^2-1114/45*x*z-232/15*y*z-241/15*z^2-1543/45*x*w-319/15*y*w-1327/45*z*w-
    457/45*w^2
x*y-182/45*x*z-11/15*y*z-38/15*z^2-284/45*x*w-17/15*y*w-266/45*z*w-146/45*w^2
y^2-16/45*x*z-28/15*y*z-4/15*z^2-22/45*x*w-61/15*y*w+32/45*z*w+92/45*w^2

The sixth surface:

> p6 := x^2*y^2 + 8*x^3*y + 4*x^4 + x*y*z^2 - x^2*z^2 - y^2*w^2
>     - 7*x*y*w^2 + 8*x^2*w^2;
> S := Surface(P,p6);
> Y,mps,typ := ClassifyRationalSurface(S);
> typ; Y;
Conic bundle
Surface over Rational Field defined by
$.1^2 + 2*$.1*$.2 + 1/4*$.1*$.3 - 1/4*$.4^2 + 1/4*$.4*$.5 +
    2*$.6^2 - 7/4*$.6*$.7 - 1/4*$.7^2,
$.2^2 - $.1*$.3,
$.2*$.4 - $.1*$.5,
$.3*$.4 - $.2*$.5,
$.2*$.6 - $.1*$.7,
$.3*$.6 - $.2*$.7,
$.5*$.6 - $.4*$.7
> mps[2]; //the fibration map
Mapping from: Srfc: S to Projective Space of dimension 1 over Rational Field
Variables: $.1, $.2
with equations :
x
y

The seventh surface:

> p7 := x^2*w^3 + y^3*w^2 + z^5;
> Y,_,typ := ClassifyRationalSurface(Surface(P,p7));
> typ; Y; Ambient(Y);
Del Pezzo degree 1
Surface over Rational Field defined by
$.1^5*$.2 + $.3^3 + $.4^2
Projective Space of dimension 3 over Rational Field
Variables: $.1, $.2, $.3, $.4
The grading is:
    1, 1, 2, 3

The eighth surface:

> p8 := w^3*y^2*z + (x*z + w^2)^3;
> Y,_,typ := ClassifyRationalSurface(Surface(P,p8));
> typ; Y;
Del Pezzo degree 6
Surface over Rational Field defined by
$.1^2 - 4*$.5^2 + $.3*$.6 - 3*$.6*$.7,
$.1*$.2 + 2*$.2*$.5 + $.3*$.7,
$.1*$.4 + 2*$.4*$.5 + $.6^2,
$.2*$.4 + $.5^2 + $.6*$.7,
$.3*$.4 - $.1*$.6 - $.4*$.7,
$.1*$.5 + 2*$.5^2 + $.6*$.7,
$.3*$.5 - $.1*$.7 - $.5*$.7,
$.2*$.6 - $.1*$.7 - $.5*$.7,
$.5*$.6 - $.4*$.7

Parametrization of Rational Surfaces

The package also includes functions to directly parametrize rational surfaces (over Q). These use the reduction to special type, described in the last section, followed by specialised algorithms for the special cases, which are described in the following sections and the section on Del Pezzo surfaces. There is a version for hypersurfaces in P3 and one for more general rational surfaces which are first birationally projected to a hypersurface. Note, however, that the projection method can be very inefficient because it often introduces nasty singularities that cause the resolution process to hang. If it is known that the surface S is non-singular, it is usually much better to use MinimalModelRationalSurface to get a birational map from S to a standard model Y and call the relevant special case parametrization routine for Y directly.

ParametrizeProjectiveHypersurface(X, P2) : Srfc, Prj -> BoolElt, MapSch
    UseFormalDesing: BoolElt            Default: false
    SetVerbose("Classify", n):          Maximum: 1
Given a surface X in PQ3 and a projective plane P2 over Q, the intrinsic returns false if the surface is not rational over Q, otherwise returns true and a birational parameterization P2 -> X.

The intrinsic begins by mapping X, as in ClassifyRationalSurface, to a surface of special type. Then specialised parametrization routines are called depending on the type.

It is assumed that X is defined over Q because some of the special type routines assume this, partly for simplicity. This will probably be generalised in future releases.

As usual, if formal or blow-up desingularisation data has not already been cached with S, it is computed at the start. Parameter UseFormalDesing is as for earlier intrinsics.

ParametrizeProjectiveSurface(X, P2) : Srfc, Prj -> BoolElt, MapSch
    SetVerbose("Classify", n):          Maximum: 1
Given an ordinary projective surface X in PrjQn for some n ≥2 and a projective plane P2 over Q, returns false if the surface is not rational over Q, otherwise return true and a birational parametrization P2 -> X.

The function finds a birational projection to a hypersurface in P3 and then calls ParametrizeProjectiveHypersurface.

It should be noted that the birational projection may produce a very singular hypersurface defined by a polynomial with large coefficients. Then, the desingularisation, classification and special parameterization routines may each be very slow. As noted above, it may be better to try to use MinimalModelRationalSurface followed by one of the specialised parametrization routines for a non-singular rational X for larger n.

Example AlgSrf_prm (H123E25)

We try to parameterize the hypersurfaces given by polynomials p1 - p8 from the previous example.

The surface defined by p1:

> P2<X,Y,W> := ProjectiveSpace(Rationals(), 2);
> ParametrizeProjectiveHypersurface(Surface(P, p1), P2);
false

The surface defined by p2:

> ParametrizeProjectiveHypersurface(Surface(P, p2), P2);
true Mapping from: Prj: P2 to Surface over Rational Field
defined by 2*x + y + 8*z + 5*w
with equations :
-1/2*X - 4*Y - 5/2*W
X
Y
W

The surface defined by p3:

> ParametrizeProjectiveHypersurface(Surface(P, p3), P2);
true Mapping from: Prj: P2 to Surface over Rational Field
defined by x^2 - 4*x*z + 3*x*w + y*z - y*w + 2*z^2 - 3*z*w + w^2
with equations :
X^2 - 2*X*W
X^2 + X*Y - 4*X*W - Y*W + 2*W^2
X^2 - 3*X*W + Y*W
X^2 - 4*X*W + Y*W + 2*W^2

The surface defined by p4:

> ParametrizeProjectiveHypersurface(Surface(P, p4), P2);
true Mapping from: Prj: P2 to Surface over Rational Field
defined by x^2*z^2 - x*y^3 - x*y*z*w + 2*y^2*w^2 - z*w^3
with equations :
2*X^2*Y^2*W^2 - Y*W^5
X^4*Y^2 - X^2*Y*W^3
X^3*Y*W^2
X^3*Y^2*W

The surface defined by p5:

> ParametrizeProjectiveHypersurface(Surface(P, p5), P2);
true Mapping from: Prj: P2 to Surface over Rational Field
defined by ...
with equations :
-1/4*X^2 + 9/2*X*Y + 1/2*X*W - 69/4*Y^2 - 87/8*Y*W + 1/2*W^2
1/2*X*Y + 11/4*X*W - 5/8*Y^2 - 131/8*Y*W - 63/8*W^2
-11/8*X*Y + 7/8*X*W + 13/4*Y^2 - 11/2*W^2
X*Y - 1/8*X*W - 23/8*Y^2 + 4*W^2

The surface defined by p6:

> ParametrizeProjectiveHypersurface(Surface(P, p6), P2);
true Mapping from: Prj: P2 to Surface over Rational Field
defined by 4*x^4 + 8*x^3*y + x^2*y^2 - x^2*z^2 + 8*x^2*w^2
            + x*y*z^2 - 7*x*y*w^2 - y^2*w^2
with equations :
-3/8*X^2*Y^2*W - 2*X^2*Y*W^2 - 8/3*X^2*W^3 - 59/24*X*Y^2*W^2
    - 38/3*X*Y*W^3 - 16*X*W^4 + 17/6*Y^2*W^3 + 44/3*Y*W^4
    + 56/3*W^5
-3/8*X^3*Y^2 - 2*X^3*Y*W - 8/3*X^3*W^2 - 59/24*X^2*Y^2*W
    - 38/3*X^2*Y*W^2 - 16*X^2*W^3 + 17/6*X*Y^2*W^2
    + 44/3*X*Y*W^3 + 56/3*X*W^4
-1/8*X^2*Y^2*W - 4/3*X^2*Y*W^2 - 8/3*X^2*W^3 - 1/12*X*Y^2*W^2
    - 16/3*X*Y*W^3 - 40/3*X*W^4 + 19/3*Y^2*W^3 + 104/3*Y*W^4
    + 48*W^5
-3/8*X^2*Y^2*W - 2*X^2*Y*W^2 - 8/3*X^2*W^3 - 8/3*X*Y^2*W^2
    - 14*X*Y*W^3 - 56/3*X*W^4 + Y^2*W^3 + 20/3*Y*W^4 + 32/3*W^5

The surface defined by p7:

> ParametrizeProjectiveHypersurface(Surface(P, p7), P2);
true Mapping from: Prj: P2 to Surface over Rational Field defined by
z^5 + y^3*w^2 + x^2*w^3
with equations :
-X^362*Y^144*W^79 - 9*X^363*Y^141*W^81 - 36*X^364*Y^138*W^83 -
    84*X^365*Y^135*W^85 - 126*X^366*Y^132*W^87 - 126*X^367*Y^129*W^89 -
    84*X^368*Y^126*W^91 - 36*X^369*Y^123*W^93 - 9*X^370*Y^120*W^95 -
    X^371*Y^117*W^97
X^359*Y^148*W^78 + 10*X^360*Y^145*W^80 + 45*X^361*Y^142*W^82 +
    120*X^362*Y^139*W^84 + 210*X^363*Y^136*W^86 + 252*X^364*Y^133*W^88 +
    210*X^365*Y^130*W^90 + 120*X^366*Y^127*W^92 + 45*X^367*Y^124*W^94 +
    10*X^368*Y^121*W^96 + X^369*Y^118*W^98
-X^357*Y^150*W^78 - 11*X^358*Y^147*W^80 - 55*X^359*Y^144*W^82 -
    165*X^360*Y^141*W^84 - 330*X^361*Y^138*W^86 - 462*X^362*Y^135*W^88 -
    462*X^363*Y^132*W^90 - 330*X^364*Y^129*W^92 - 165*X^365*Y^126*W^94 -
    55*X^366*Y^123*W^96 - 11*X^367*Y^120*W^98 - X^368*Y^117*W^100
X^354*Y^153*W^78 + 12*X^355*Y^150*W^80 + 66*X^356*Y^147*W^82 +
    220*X^357*Y^144*W^84 + 495*X^358*Y^141*W^86 + 792*X^359*Y^138*W^88 +
    924*X^360*Y^135*W^90 + 792*X^361*Y^132*W^92 + 495*X^362*Y^129*W^94 +
    220*X^363*Y^126*W^96 + 66*X^364*Y^123*W^98 + 12*X^365*Y^120*W^100 +
    X^366*Y^117*W^102
and inverse
z^4*w^8
y*z^2*w^9
x*z*w^10

The surface defined by p8:

> ParametrizeProjectiveHypersurface(Surface(P, p8), P2);
true Mapping from: Prj: P2 to Surface over Rational Field defined by
x^3*z^3 + 3*x^2*z^2*w^2 + 3*x*z*w^4 + y^2*z*w^3 + w^6
with equations :
2*X^20*Y^3*W - 4*X^19*Y^3*W^2 + 4*X^17*Y^3*W^4 - 2*X^16*Y^3*W^5
4*X^19*Y^4*W - 8*X^18*Y^4*W^2 + 4*X^17*Y^4*W^3
1/2*X^16*Y^7*W - X^15*Y^7*W^2 + 1/2*X^14*Y^7*W^3
-X^18*Y^5*W + 3*X^17*Y^5*W^2 - 3*X^16*Y^5*W^3 + X^15*Y^5*W^4
> IsRational(Surface(P, p7));
true

Here are two easy examples of parameterizing non-hypersurfaces.

> P4<u,v,w,x,y> := ProjectiveSpace(Rationals(),4);
> P2<X,Y,Z> := ProjectiveSpace(Rationals(), 2);
> S := Surface(P4,[u^2 + v^2 + w^2 - x^2, y - x]);
> ParametrizeProjectiveSurface(S, P2);
true Mapping from: Prj: P2 to Srfc: S
with equations :
-2*X*Z + 2*Z^2
-2*Y*Z
X^2 + Y^2 - 2*X*Z
-X^2 - Y^2 + 2*X*Z - 2*Z^2
-X^2 - Y^2 + 2*X*Z - 2*Z^2
and inverse
u + w + y
v
w + x

Here we parametrize a particularly easy surface -- P2 itself!

> S := Surface(P2,[]);
> ParametrizeProjectiveSurface(S, P2);
true Mapping from: Prj: P2 to Srfc: S
with equations :
X
Y
Z
and inverse
X
Y
Z
Solve(p, F) : RngMPolElt, FldFunRat -> SeqEnum
For convenience, this intrinsic provides a purely algebraic version of parameterization of an affine hypersurface.

Given a polynomial p ∈Q[x, y, z], the equation of a (not necessarily irreducible affine hypersurface S) and a two-variable rational function field F = Q(u, v), the function finds birational parameterizations of the irreducible components of S (that are parametrizable over Q).

A sequence of triples (X, Y, Z) ∈F3 is returned such that p(X, Y, Z) = 0 and each triple gives an isomorphism of F to the function field of a component of S.

The routine may again result in a runtime error if it involves parameterizations of special surface types that are not yet implemented, as for the preceding functions.

Example AlgSrf_prm2 (H123E26)

The following affine hypersurface has three irreducible factors: one not rational and two that are rational and parameterizable over Q.
> Q := RationalField();
> P<x,y,z> := PolynomialRing(Q, 3);
> F<s,t> := RationalFunctionField(Q, 2);
> p := (x^4+y^4-z^2)*(2*x + y + 8*z + 5)
>       *(x^2 - 4*x*z + 3*x + y*z - y + 2*z^2 - 3*z + 1);
> Solve(p, F);
[
    [ -1/2*s - 4*t - 5/2, s, t ],
    [
        (s^2 - 2*s)/(s^2 - 4*s + t + 2),
        (s^2 + s*t - 4*s - t + 2)/(s^2 - 4*s + t + 2),
        (s^2 - 3*s + t)/(s^2 - 4*s + t + 2)
    ]
]

Parametrization of Special Surfaces

In this section we describe routines for the explicit parameterization of the special classes of rational surfaces that arise from the reduction of the general case via m-adjoint maps. The algorithms are the work of Josef Schicho, in collaboration with others in some cases. The functions are used in the general parameterization routines but can also be called directly by the user.

ParametrizeQuadric(X,P2) : Sch, Prj -> BoolElt, MapSch
Suppose the scheme X ⊂PQ3 is a geometrically irreducible quadric (degree 2) projective hypersurface and P2 is a projective plane. The intrinsic returns false if X is not parametrizable over the rationals, otherwise it returns true together with a birational parameterization P2 -> X.

Given a rational point p on X, the intrinsic is based on a simple, well-known algorithm (see [Sch98, Sec. 3.1]). Finding the point p is equivalent to finding a non-trivial isotropic vector for F, the quadric form in four variables defining X. This is achieved by a reduction to the solution of two quadrics in three variables, which is performed using standard lattice methods. The solubility routines here assume that the quadric is defined over Q.

Example AlgSrf_solve_quad (H123E27)

We give a few simple examples.

> Q := Rationals();
> P2<X,Y,W> := ProjectiveSpace(Q, 2);
> P3<x,y,z,w> := ProjectiveSpace(Q, 3);
> X1 := Scheme(P3, x^2 + y^2 + z^2 + w^2);
> X2 := Scheme(P3, x^2 + y^2 + z^2 - w^2);
> X3 := Scheme(P3, x^2 + y^2 + z^2);
> X4 := Scheme(P3, x^2 + y^2 - z^2);
> X5 := Scheme(P3, x^2 - 4*x*z + 3*x*w + y*z - y*w + 2*z^2
>                   - 3*z*w + w^2);
> ParametrizeQuadric(X1, P2);
false
> ParametrizeQuadric(X2, P2);
true Mapping from: Prj: P2 to Sch: X2
with equations :
2*X*W
2*Y*W
-X^2 - Y^2 + W^2
X^2 + Y^2 + W^2
> ParametrizeQuadric(X3, P2);
false
> ParametrizeQuadric(X4, P2);
true Mapping from: Prj: P2 to Sch: X4
with equations :
X*Y
-1/2*X^2 + 1/2*Y^2
1/2*X^2 + 1/2*Y^2
Y*W
> ParametrizeQuadric(X5, P2);
true Mapping from: Prj: P2 to Sch: X5
with equations :
X^2 - 2*X*W
X^2 + X*Y - 4*X*W - Y*W + 2*W^2
X^2 - 3*X*W + Y*W
X^2 - 4*X*W + Y*W + 2*W^2
ParametrizePencil(phi, P2) : MapSch, Prj -> BoolElt, MapSch
Let X be an ordinary projective birationally ruled surface, given as the domain of a rational pencil φ defined over Q (i.e., a rational map X -> PQn for some n with image a rational normal curve) and P2 is a projective plane over Q. The intrinsic returns false if X is not parameterizable over the rationals. Otherwise, it returns true and a birational parameterization P2 -> X.

These intrinsics take care of rational scrolls and conic bundles in the classification of special surfaces. The algorithm is described in [Sch00]. The scheme X is not assumed to be non-singular.

Example AlgSrf_ruled_ex (H123E28)

We start with a (singular) degree 4 hypersurface X in P3 and construct a pencil from P2 to X.
> Q := Rationals();
> P3<x,y,z,w> := ProjectiveSpace(Q, 3);
> P2<X,Y,Z> := ProjectiveSpace(Q, 2);
> X := Scheme(P3, x^2*z^2 - x*y^3 - x*y*z*w + 2*y^2*w^2 - z*w^3);
> pencil := map<X -> P2 | [x*y - w^2, y^2 - z*w, x*z - y*w]>;
> DefiningPolynomial(Image(pencil));
X*Y - Z^2
> ParametrizePencil(pencil, P2);
true Mapping from: Prj: P2 to Sch: X
with equations :
-X^4*Y + 2*X^2*Z^3
-X^2*Y*Z^2 + Z^5
X*Y*Z^3
X*Z^4
ParametrizeDelPezzo(X, P2) : Sch, Prj -> BoolElt, MapSch
The argument X should be a (anticanonically-embedded) Del Pezzo surface (of type Sch or Srfc and P2 a projective plane both defined over Q. The function returns false if X is not parametrizable over the rationals and returns true and a birational parameterization P2 -> X otherwise.

This intrinsic is the main interface to a suite of functions parameterizing Del Pezzo surfaces (over Q). These include the anticanonically Del Pezzos of degrees 1 and 2, which lie in non-trivially weighted projective spaces. A degree d Del Pezzo surface refers to a rational surface that is embedded in projective space by its anti-canonical divisor. For degrees 1 and 2 this means an ample embedding into weighted projective space. For 3 <= d <= 9, this is a very-ample embedding giving X as a degree d surface in ordinary projective space of dimension d.

It should be noted that not only do the routines handle the usual non-singular cases, but they also deal with degenerate singular cases (arising in degrees d >= 3. In the latter case, rather than blowing up 9 - d distinct points in the plane, some of the blown-up points are "infinitely near" points: lying on the exceptional curves corresponding to already blown-up points). This is important as these degenerate cases can arise in the general parameterization of hypersurfaces in P3.

The package contains routines to find and blow down sets (defined over Q) of exceptional curves, reducing to a DelPezzo of degree d, 5 <= d <= 9. After this reduction, these cases are handled by the Magma routines described in Section Parametrization of Del Pezzo Surfaces (which now also handle singular Del Pezzos). There is an exception to the above rule. For singular degree 3 and degree 4 Del Pezzo surfaces, it is more efficient to apply special case code directly rather than trying to blow down lines to get to higher degree. Starting in V2.17, special functions are provided do this that can be called directly for a singular anti-canonical degree 3 or 4 surface and are described in the next section.

The routines are not yet fully documented. For the location of exceptional curves and the blowing-down, some details may be found in [Sch98, Sec. 3.5] while [Man86] contains the general theory. There are also implementation notes in the appendix of the software documentation report [Bec08] from which this documentation has been adapted. For the nonsingular degrees 6, 8 and 9 cases, which use the Lie algebra method and the degree 5 case, see the references in Section Parametrization of Del Pezzo Surfaces.

V2.28, 13 July 2023