Morphisms on affine schemes#

This module implements morphisms from affine schemes. A morphism from an affine scheme to an affine scheme is determined by rational functions that define what the morphism does on points in the ambient affine space. A morphism from an affine scheme to a projective scheme is determined by homogeneous polynomials.

EXAMPLES:

sage: A2.<x,y> = AffineSpace(QQ, 2)
sage: P2.<x0,x1,x2> = ProjectiveSpace(QQ, 2)
sage: A2.hom([x, x + y], A2)
Scheme endomorphism of Affine Space of dimension 2 over Rational Field
  Defn: Defined on coordinates by sending (x, y) to (x, x + y)
sage: A2.hom([1, x, x + y], P2)
Scheme morphism:
  From: Affine Space of dimension 2 over Rational Field
  To:   Projective Space of dimension 2 over Rational Field
  Defn: Defined on coordinates by sending (x, y) to (1 : x : x + y)

AUTHORS:

  • David Kohel, William Stein: initial version

  • Volker Braun (2011-08-08): renamed classes, more documentation, misc cleanups

  • Ben Hutz (2013-03): iteration functionality and new directory structure for affine/projective

  • Kwankyu Lee (2020-02): added indeterminacy_locus() and image()

class sage.schemes.affine.affine_morphism.SchemeMorphism_polynomial_affine_space(parent, polys, check=True)#

Bases: SchemeMorphism_polynomial

A morphism of schemes determined by rational functions.

EXAMPLES:

sage: RA.<x,y> = QQ[]
sage: A2 = AffineSpace(RA)
sage: RP.<u,v,w> = QQ[]
sage: P2 = ProjectiveSpace(RP)
sage: H = A2.Hom(P2)
sage: f = H([x, y, 1])
sage: f
Scheme morphism:
  From: Affine Space of dimension 2 over Rational Field
  To:   Projective Space of dimension 2 over Rational Field
  Defn: Defined on coordinates by sending (x, y) to (x : y : 1)
as_dynamical_system()#

Return this endomorphism as a DynamicalSystem_affine.

OUTPUT:

EXAMPLES:

sage: A.<x,y,z> = AffineSpace(ZZ, 3)
sage: H = End(A)
sage: f = H([x^2, y^2, z^2])
sage: type(f.as_dynamical_system())                                         # needs sage.schemes
<class 'sage.dynamics.arithmetic_dynamics.affine_ds.DynamicalSystem_affine'>
sage: A.<x,y> = AffineSpace(ZZ, 2)
sage: H = End(A)
sage: f = H([x^2 - y^2, y^2])
sage: type(f.as_dynamical_system())                                         # needs sage.schemes
<class 'sage.dynamics.arithmetic_dynamics.affine_ds.DynamicalSystem_affine'>
sage: A.<x> = AffineSpace(GF(5), 1)
sage: H = End(A)
sage: f = H([x^2])
sage: type(f.as_dynamical_system())                                         # needs sage.schemes
<class 'sage.dynamics.arithmetic_dynamics.affine_ds.DynamicalSystem_affine_finite_field'>
sage: P.<x,y> = AffineSpace(RR, 2)
sage: f = DynamicalSystem([x^2 + y^2, y^2], P)                              # needs sage.schemes
sage: g = f.as_dynamical_system()                                           # needs sage.schemes
sage: g is f                                                                # needs sage.schemes
True
degree()#

Return the degree of the affine morphism.

EXAMPLES:

sage: R.<x> = AffineSpace(QQ, 1)
sage: H = Hom(R, R)
sage: f = H([x^7])
sage: f.degree()
7
sage: R.<x,y,z> = AffineSpace(QQ, 3)
sage: H = Hom(R, R)
sage: f = H([x^3, y^2 + 5, z^4 + y])
sage: f.degree()
4
global_height(prec=None)#

Take the height of the homogenization, and return the global height of the coefficients as a projective point.

INPUT:

  • prec – desired floating point precision (default: default RealField precision).

OUTPUT: A real number.

EXAMPLES:

sage: A.<x> = AffineSpace(QQ, 1)
sage: H = Hom(A, A)
sage: f = H([1/1331*x^2 + 4000])
sage: f.global_height()                                                     # needs sage.symbolic
15.4877354584971
sage: # needs sage.rings.number_field
sage: R.<x> = PolynomialRing(QQ)
sage: k.<w> = NumberField(x^2 + 5)
sage: A.<x,y> = AffineSpace(k, 2)
sage: H = Hom(A, A)
sage: f = H([13*w*x^2 + 4*y, 1/w*y^2])
sage: f.global_height(prec=2)
4.0
sage: A.<x> = AffineSpace(ZZ, 1)
sage: H = Hom(A, A)
sage: f = H([7*x^2 + 1513])
sage: f.global_height()                                                     # needs sage.symbolic
7.32184971378836
sage: A.<x> = AffineSpace(QQ, 1)
sage: B.<y,z> = AffineSpace(QQ, 2)
sage: H = Hom(A, B)
sage: f = H([1/3*x^2 + 10, 7*x^3])
sage: f.global_height()                                                     # needs sage.symbolic
3.40119738166216
sage: P.<x,y> = AffineSpace(QQ, 2)
sage: A.<z> = AffineSpace(QQ, 1)
sage: H = Hom(P, A)
sage: f = H([1/1331*x^2 + 4000*y])
sage: f.global_height()                                                     # needs sage.symbolic
15.4877354584971
homogenize(n)#

Return the homogenization of this map.

If it’s domain is a subscheme, the domain of the homogenized map is the projective embedding of the domain. The domain and codomain can be homogenized at different coordinates: n[0] for the domain and n[1] for the codomain.

INPUT:

  • n – a tuple of nonnegative integers. If n is an integer, then the two values of the tuple are assumed to be the same

OUTPUT: a morphism from the projective embedding of the domain of this map

EXAMPLES:

sage: A.<x,y> = AffineSpace(ZZ, 2)
sage: H = Hom(A, A)
sage: f = H([(x^2-2)/x^5, y^2])
sage: f.homogenize(2)
Scheme endomorphism of Projective Space of dimension 2 over Integer Ring
  Defn: Defined on coordinates by sending (x0 : x1 : x2) to
        (x0^2*x2^5 - 2*x2^7 : x0^5*x1^2 : x0^5*x2^2)
sage: # needs sage.rings.real_mpfr
sage: A.<x,y> = AffineSpace(CC, 2)
sage: H = Hom(A, A)
sage: f = H([(x^2-2)/(x*y), y^2 - x])
sage: f.homogenize((2, 0))
Scheme endomorphism of Projective Space of dimension 2
 over Complex Field with 53 bits of precision
  Defn: Defined on coordinates by sending (x0 : x1 : x2) to
        (x0*x1*x2^2 : x0^2*x2^2 + (-2.00000000000000)*x2^4 : x0*x1^3 - x0^2*x1*x2)
sage: A.<x,y> = AffineSpace(ZZ, 2)
sage: X = A.subscheme([x - y^2])
sage: H = Hom(X, X)
sage: f = H([9*y^2, 3*y])
sage: f.homogenize(2)                                                       # needs sage.libs.singular
Scheme endomorphism of Closed subscheme of Projective Space
 of dimension 2 over Integer Ring defined by: x1^2 - x0*x2
  Defn: Defined on coordinates by sending (x0 : x1 : x2) to
        (9*x1^2 : 3*x1*x2 : x2^2)
sage: R.<t> = PolynomialRing(ZZ)
sage: A.<x,y> = AffineSpace(R, 2)
sage: H = Hom(A, A)
sage: f = H([(x^2-2)/y, y^2 - x])
sage: f.homogenize((2, 0))
Scheme endomorphism of Projective Space of dimension 2
 over Univariate Polynomial Ring in t over Integer Ring
  Defn: Defined on coordinates by sending (x0 : x1 : x2) to
        (x1*x2^2 : x0^2*x2 + (-2)*x2^3 : x1^3 - x0*x1*x2)
sage: A.<x> = AffineSpace(QQ, 1)
sage: H = End(A)
sage: f = H([x^2 - 1])
sage: f.homogenize((1, 0))
Scheme endomorphism of Projective Space of dimension 1 over Rational Field
  Defn: Defined on coordinates by sending (x0 : x1) to
        (x1^2 : x0^2 - x1^2)
sage: # needs sage.rings.number_field
sage: R.<a> = PolynomialRing(QQbar)
sage: A.<x,y> = AffineSpace(R, 2)
sage: H = End(A)
sage: f = H([QQbar(sqrt(2))*x*y, a*x^2])                                    # needs sage.symbolic
sage: f.homogenize(2)                                                       # needs sage.libs.singular sage.symbolic
Scheme endomorphism of Projective Space of dimension 2
 over Univariate Polynomial Ring in a over Algebraic Field
  Defn: Defined on coordinates by sending (x0 : x1 : x2) to
        (1.414213562373095?*x0*x1 : a*x0^2 : x2^2)
sage: P.<x,y,z> = AffineSpace(QQ, 3)
sage: H = End(P)
sage: f = H([x^2 - 2*x*y + z*x, z^2 -y^2 , 5*z*y])
sage: f.homogenize(2).dehomogenize(2) == f
True
sage: K.<c> = FunctionField(QQ)
sage: A.<x> = AffineSpace(K, 1)
sage: f = Hom(A, A)([x^2 + c])
sage: f.homogenize(1)
Scheme endomorphism of Projective Space of dimension 1
 over Rational function field in c over Rational Field
  Defn: Defined on coordinates by sending (x0 : x1) to
        (x0^2 + c*x1^2 : x1^2)
sage: # needs sage.rings.number_field
sage: A.<z> = AffineSpace(QQbar, 1)
sage: H = End(A)
sage: f = H([2*z / (z^2 + 2*z + 3)])
sage: f.homogenize(1)
Scheme endomorphism of Projective Space of dimension 1
 over Algebraic Field
  Defn: Defined on coordinates by sending (x0 : x1) to
        (x0*x1 : 1/2*x0^2 + x0*x1 + 3/2*x1^2)
sage: # needs sage.rings.number_field
sage: R.<c,d> = QQbar[]
sage: A.<x> = AffineSpace(R, 1)
sage: H = Hom(A, A)
sage: F = H([d*x^2 + c])
sage: F.homogenize(1)
Scheme endomorphism of Projective Space of dimension 1
 over Multivariate Polynomial Ring in c, d over Algebraic Field
  Defn: Defined on coordinates by sending (x0 : x1) to
        (d*x0^2 + c*x1^2 : x1^2)
jacobian()#

Return the Jacobian matrix of partial derivative of this map.

The \((i, j)\) entry of the Jacobian matrix is the partial derivative diff(functions[i], variables[j]).

OUTPUT:

  • matrix with coordinates in the coordinate ring of the map.

EXAMPLES:

sage: A.<z> = AffineSpace(QQ, 1)
sage: H = End(A)
sage: f = H([z^2 - 3/4])
sage: f.jacobian()                                                          # needs sage.modules
[2*z]
sage: A.<x,y> = AffineSpace(QQ, 2)
sage: H = End(A)
sage: f = H([x^3 - 25*x + 12*y, 5*y^2*x - 53*y + 24])
sage: f.jacobian()                                                          # needs sage.modules
[ 3*x^2 - 25          12]
[      5*y^2 10*x*y - 53]
sage: A.<x,y> = AffineSpace(ZZ, 2)
sage: H = End(A)
sage: f = H([(x^2 - x*y)/(1+y), (5+y)/(2+x)])
sage: f.jacobian()                                                          # needs sage.modules
[         (2*x - y)/(y + 1) (-x^2 - x)/(y^2 + 2*y + 1)]
[  (-y - 5)/(x^2 + 4*x + 4)                  1/(x + 2)]
local_height(v, prec=None)#

Return the maximum of the local heights of the coefficients in any of the coordinate functions of this map.

INPUT:

  • v – a prime or prime ideal of the base ring.

  • prec – desired floating point precision (default: default RealField precision).

OUTPUT:

  • a real number.

EXAMPLES:

sage: P.<x,y> = AffineSpace(QQ, 2)
sage: H = Hom(P, P)
sage: f = H([1/1331*x^2 + 1/4000*y^2, 210*x*y])
sage: f.local_height(1331)                                                  # needs sage.rings.real_mpfr
7.19368581839511
sage: P.<x,y,z> = AffineSpace(QQ, 3)
sage: H = Hom(P, P)
sage: f = H([4*x^2 + 3/100*y^2, 8/210*x*y, 1/10000*z^2])
sage: f.local_height(2)                                                     # needs sage.rings.real_mpfr
2.77258872223978
sage: P.<x,y,z> = AffineSpace(QQ, 3)
sage: H = Hom(P, P)
sage: f = H([4*x^2 + 3/100*y^2, 8/210*x*y, 1/10000*z^2])
sage: f.local_height(2, prec=2)                                             # needs sage.rings.real_mpfr
3.0
sage: # needs sage.rings.number_field
sage: R.<z> = PolynomialRing(QQ)
sage: K.<w> = NumberField(z^2 - 2)
sage: P.<x,y> = AffineSpace(K, 2)
sage: H = Hom(P, P)
sage: f = H([2*x^2 + w/3*y^2, 1/w*y^2])
sage: f.local_height(K.ideal(3))
1.09861228866811
local_height_arch(i, prec=None)#

Return the maximum of the local height at the i-th infinite place of the coefficients in any of the coordinate functions of this map.

INPUT:

  • i – an integer.

  • prec – desired floating point precision (default: default RealField precision).

OUTPUT:

  • a real number.

EXAMPLES:

sage: P.<x,y> = AffineSpace(QQ, 2)
sage: H = Hom(P, P)
sage: f = H([1/1331*x^2 + 1/4000*y^2, 210*x*y]);
sage: f.local_height_arch(0)                                                # needs sage.rings.real_mpfr
5.34710753071747
sage: P.<x,y> = AffineSpace(QQ, 2)
sage: H = Hom(P, P)
sage: f = H([1/1331*x^2 + 1/4000*y^2, 210*x*y]);
sage: f.local_height_arch(0, prec=5)                                        # needs sage.rings.real_mpfr
5.2
sage: # needs sage.rings.number_field
sage: R.<z> = PolynomialRing(QQ)
sage: K.<w> = NumberField(z^2 - 2)
sage: P.<x,y> = AffineSpace(K, 2)
sage: H = Hom(P, P)
sage: f = H([2*x^2 + w/3*y^2, 1/w*y^2])
sage: f.local_height_arch(1)
0.6931471805599453094172321214582
class sage.schemes.affine.affine_morphism.SchemeMorphism_polynomial_affine_space_field(parent, polys, check=True)#

Bases: SchemeMorphism_polynomial_affine_space

image()#

Return the scheme-theoretic image of the morphism.

OUTPUT: a subscheme of the ambient space of the codomain

EXAMPLES:

sage: # needs sage.libs.singular
sage: A1.<w> = AffineSpace(QQ, 1)
sage: A2.<x,y> = AffineSpace(QQ, 2)
sage: f = A2.hom([x + y], A1)
sage: f.image()
Closed subscheme of Affine Space of dimension 1 over Rational Field defined by:
  (no polynomials)
sage: f = A2.hom([x, x], A2)
sage: f.image()
Closed subscheme of Affine Space of dimension 2 over Rational Field defined by:
  x - y
sage: f = A2.hom([x^2, x^3], A2)
sage: f.image()
Closed subscheme of Affine Space of dimension 2 over Rational Field defined by:
  x^3 - y^2
sage: P2.<x0,x1,x2> = ProjectiveSpace(QQ, 2)
sage: f = A2.hom([x, x^2, x^3], P2)
sage: f.image()
Closed subscheme of Projective Space of dimension 2 over Rational Field defined by:
  x1^2 - x0*x2
indeterminacy_locus()#

Return the indeterminacy locus of this map as a rational map on the domain.

The indeterminacy locus is the intersection of all the base indeterminacy locuses of maps that define the same rational map as by this map.

OUTPUT: a subscheme of the domain of the map

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: H = End(A)
sage: f = H([x - y, x^2 - y^2])
sage: f.indeterminacy_locus()                                               # needs sage.libs.singular
Closed subscheme of Affine Space of dimension 2 over Rational Field defined by:
  1
sage: A.<x,y> = AffineSpace(QQ, 2)
sage: f = A.hom([x, x/y], A)
sage: f.indeterminacy_locus()                                               # needs sage.libs.singular
Closed subscheme of Affine Space of dimension 2 over Rational Field defined by:
  y
indeterminacy_points(F=None)#

Return the points in the indeterminacy locus of this map.

If the dimension of the indeterminacy locus is not zero, an error is raised.

INPUT:

  • F – a field; if not given, the base ring of the domain is assumed

OUTPUT: indeterminacy points of the map defined over F

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: H = End(A)
sage: f = H([x - y, x^2 - y^2])
sage: f.indeterminacy_points()                                              # needs sage.libs.singular
[]
sage: A2.<x,y> = AffineSpace(QQ, 2)
sage: P2.<x0,x1,x2> = ProjectiveSpace(QQ, 2)
sage: f = A2.hom([x*y, y, x], P2)
sage: f.indeterminacy_points()                                              # needs sage.libs.singular
[(0, 0)]
reduce_base_field()#

Return this map defined over the field of definition of the coefficients.

The base field of the map could be strictly larger than the field where all of the coefficients are defined. This function reduces the base field to the minimal possible. This can be done when the base ring is a number field, QQbar, a finite field, or algebraic closure of a finite field.

OUTPUT: a scheme morphism

EXAMPLES:

sage: # needs sage.rings.finite_rings
sage: K.<t> = GF(5^4)
sage: A.<x> = AffineSpace(K, 1)
sage: A2.<a,b> = AffineSpace(K, 2)
sage: H = End(A)
sage: H2 = Hom(A, A2)
sage: H3 = Hom(A2, A)
sage: f = H([x^2 + 2*(t^3 + t^2 + t + 3)])
sage: f.reduce_base_field()
Scheme endomorphism of Affine Space of dimension 1
 over Finite Field in t2 of size 5^2
  Defn: Defined on coordinates by sending (x) to (x^2 + (2*t2))
sage: f2 = H2([x^2 + 4, 2*x])
sage: f2.reduce_base_field()
Scheme morphism:
  From: Affine Space of dimension 1 over Finite Field of size 5
  To:   Affine Space of dimension 2 over Finite Field of size 5
  Defn: Defined on coordinates by sending (x) to (x^2 - 1, 2*x)
sage: f3 = H3([a^2 + t*b])
sage: f3.reduce_base_field()
Scheme morphism:
  From: Affine Space of dimension 2 over Finite Field in t of size 5^4
  To:   Affine Space of dimension 1 over Finite Field in t of size 5^4
  Defn: Defined on coordinates by sending (a, b) to (a^2 + t*b)
sage: # needs sage.rings.number_field
sage: K.<v> = CyclotomicField(4)
sage: A.<x> = AffineSpace(K, 1)
sage: H = End(A)
sage: f = H([x^2 + v])
sage: g = f.reduce_base_field(); g
Scheme endomorphism of Affine Space of dimension 1
 over Cyclotomic Field of order 4 and degree 2
  Defn: Defined on coordinates by sending (x) to (x^2 + v)
sage: g.base_ring() is K
True
sage: # needs sage.rings.number_field
sage: A.<x> = AffineSpace(QQbar, 1)
sage: H = End(A)
sage: f = H([(QQbar(sqrt(2))*x^2 + 1/QQbar(sqrt(3))) / (5*x)])              # needs sage.symbolic
sage: f.reduce_base_field()                                                 # needs sage.symbolic
Scheme endomorphism of Affine Space of dimension 1 over Number Field in a
 with defining polynomial y^4 - 4*y^2 + 1 with a = ...?
  Defn: Defined on coordinates by sending (x) to
        (((a^3 - 3*a)*x^2 + (-1/3*a^2 + 2/3))/(5*x))
sage: # needs sage.rings.number_field
sage: R.<x> = PolynomialRing(QQ)
sage: A.<x> = AffineSpace(QQbar, 1)
sage: H = End(A)
sage: f = H([QQbar(3^(1/3))*x^2 + QQbar(sqrt(-2))])                         # needs sage.symbolic
sage: f.reduce_base_field()                                                 # needs sage.symbolic
Scheme endomorphism of Affine Space of dimension 1 over Number
Field in a with defining polynomial y^6 + 6*y^4 - 6*y^3 + 12*y^2 + 36*y + 17
 with a = 1.442249570307409? + 1.414213562373095?*I
  Defn: Defined on coordinates by sending (x) to
        ((-48/269*a^5 + 27/269*a^4 - 320/269*a^3 + 468/269*a^2 - 772/269*a
        - 1092/269)*x^2 + (48/269*a^5 - 27/269*a^4 + 320/269*a^3 - 468/269*a^2
        + 1041/269*a + 1092/269))
sage: # needs sage.rings.number_field
sage: R.<x> = PolynomialRing(QQ)
sage: K.<a> = NumberField(x^3 - x + 1,
....:                     embedding=(x^3 + x + 1).roots(ring=CC)[0][0])
sage: A.<x> = AffineSpace(K, 1)
sage: A2.<u,v> = AffineSpace(K, 2)
sage: H = Hom(A, A2)
sage: f = H([x^2 + a*x + 3, 5*x])
sage: f.reduce_base_field()
Scheme morphism:
  From: Affine Space of dimension 1 over Number Field in a with
        defining polynomial x^3 - x + 1 with a = -1.324717957244746?
  To:   Affine Space of dimension 2 over Number Field in a with
        defining polynomial x^3 - x + 1 with a = -1.324717957244746?
  Defn: Defined on coordinates by sending (x) to (x^2 + a*x + 3, 5*x)
sage: # needs sage.rings.number_field
sage: K.<v> = QuadraticField(2)
sage: A.<x> = AffineSpace(K, 1)
sage: H = End(A)
sage: f = H([3*x^2 + x + 1])
sage: f.reduce_base_field()
Scheme endomorphism of Affine Space of dimension 1 over Rational Field
  Defn: Defined on coordinates by sending (x) to (3*x^2 + x + 1)
sage: # needs sage.rings.finite_rings
sage: K.<t> = GF(5^6)
sage: A.<x> = AffineSpace(K, 1)
sage: H = End(A)
sage: f = H([x^2 + x*(t^3 + 2*t^2 + 4*t) + (t^5 + 3*t^4 + t^2 + 4*t)])
sage: f.reduce_base_field()
Scheme endomorphism of Affine Space of dimension 1
 over Finite Field in t of size 5^6
  Defn: Defined on coordinates by sending (x) to
        (x^2 + (t^3 + 2*t^2 - t)*x + (t^5 - 2*t^4 + t^2 - t))
weil_restriction()#

Compute the Weil restriction of this morphism over some extension field.

If the field is a finite field, then this computes the Weil restriction to the prime subfield.

A Weil restriction of scalars - denoted \(Res_{L/k}\) - is a functor which, for any finite extension of fields \(L/k\) and any algebraic variety \(X\) over \(L\), produces another corresponding variety \(Res_{L/k}(X)\), defined over \(k\). It is useful for reducing questions about varieties over large fields to questions about more complicated varieties over smaller fields. Since it is a functor it also applied to morphisms. In particular, the functor applied to a morphism gives the equivalent morphism from the Weil restriction of the domain to the Weil restriction of the codomain.

OUTPUT: Scheme morphism on the Weil restrictions of the domain

and codomain of the map.

EXAMPLES:

sage: # needs sage.rings.number_field
sage: K.<v> = QuadraticField(5)
sage: A.<x,y> = AffineSpace(K, 2)
sage: H = End(A)
sage: f = H([x^2 - y^2, y^2])
sage: f.weil_restriction()                                                  # needs sage.libs.singular
Scheme endomorphism of Affine Space of dimension 4 over Rational Field
  Defn: Defined on coordinates by sending (z0, z1, z2, z3) to
        (z0^2 + 5*z1^2 - z2^2 - 5*z3^2, 2*z0*z1 - 2*z2*z3, z2^2 + 5*z3^2, 2*z2*z3)
sage: # needs sage.rings.number_field
sage: K.<v> = QuadraticField(5)
sage: PS.<x,y> = AffineSpace(K, 2)
sage: H = Hom(PS, PS)
sage: f = H([x, y])
sage: F = f.weil_restriction()
sage: P = PS(2, 1)
sage: Q = P.weil_restriction()
sage: f(P).weil_restriction() == F(Q)                                       # needs sage.libs.singular
True
class sage.schemes.affine.affine_morphism.SchemeMorphism_polynomial_affine_space_finite_field(parent, polys, check=True)#

Bases: SchemeMorphism_polynomial_affine_space_field

class sage.schemes.affine.affine_morphism.SchemeMorphism_polynomial_affine_subscheme_field(parent, polys, check=True)#

Bases: SchemeMorphism_polynomial_affine_space_field

Morphisms from subschemes of affine spaces defined over fields.

image()#

Return the scheme-theoretic image of the morphism.

OUTPUT: a subscheme of the ambient space of the codomain

EXAMPLES:

sage: A1.<w> = AffineSpace(QQ, 1)
sage: A2.<x,y> = AffineSpace(QQ, 2)
sage: X = A2.subscheme(0)
sage: f = X.hom([x + y], A1)
sage: f.image()                                                             # needs sage.libs.singular
Closed subscheme of Affine Space of dimension 1 over Rational Field defined by:
  (no polynomials)
sage: A2.<x,y> = AffineSpace(QQ, 2)
sage: X = A2.subscheme([x*y^2 - y^3 - 1])
sage: f = X.hom([y, y/x], A2)
sage: f.image()                                                             # needs sage.libs.singular
Closed subscheme of Affine Space of dimension 2 over Rational Field defined by:
  -x^3*y + x^3 - y
indeterminacy_locus()#

Return the indeterminacy locus of this map.

The map defines a rational map on the domain. The output is the subscheme of the domain on which the rational map is not defined by any representative of the rational map. See representatives().

EXAMPLES:

sage: A2.<x1,x2> = AffineSpace(QQ, 2)
sage: X = A2.subscheme(0)
sage: A1.<x> = AffineSpace(QQ, 1)
sage: f = X.hom([x1/x2], A1)
sage: f.indeterminacy_locus()                                               # needs sage.libs.singular
Closed subscheme of Affine Space of dimension 2 over Rational Field defined by:
  x2
sage: A2.<x1,x2> = AffineSpace(QQ, 2)
sage: X = A2.subscheme(0)
sage: P1.<a,b> = ProjectiveSpace(QQ, 1)
sage: f = X.hom([x1,x2], P1)
sage: L = f.indeterminacy_locus()                                           # needs sage.libs.singular
sage: L.rational_points()                                                   # needs sage.libs.singular
[(0, 0)]
sage: A2.<x,y> = AffineSpace(QQ, 2)
sage: X = A2.subscheme([x^2 - y^2 - y])
sage: A1.<a> = AffineSpace(QQ, 1)
sage: f = X.hom([x/y], A1)
sage: f.indeterminacy_locus()                                               # needs sage.libs.singular
Closed subscheme of Affine Space of dimension 2 over Rational Field defined by:
  y,
  x
sage: A3.<x,y,z> = AffineSpace(QQ, 3)
sage: X = A3.subscheme(x^2 - y*z - x)
sage: A2.<a,b> = AffineSpace(QQ, 2)
sage: f = X.hom([y, y/x], A2)
sage: L = f.indeterminacy_locus(); L                                        # needs sage.libs.singular
Closed subscheme of Affine Space of dimension 3 over Rational Field defined by:
  x,
  y*z
sage: L.dimension()                                                         # needs sage.libs.singular
1
is_morphism()#

Return True if the map is defined everywhere on the domain.

EXAMPLES:

sage: P2.<x,y,z> = ProjectiveSpace(QQ,2)
sage: P1.<a,b> = ProjectiveSpace(QQ,1)
sage: X = P2.subscheme([x^2 - y^2 - y*z])
sage: f = X.hom([x,y], P1)
sage: f.is_morphism()                                                       # needs sage.libs.singular
True
representatives()#

Return all maps representing the same rational map as by this map.

EXAMPLES:

sage: A2.<x,y> = AffineSpace(QQ, 2)
sage: X = A2.subscheme(0)
sage: f = X.hom([x, x/y], A2)
sage: f.representatives()                                                   # needs sage.libs.singular
[Scheme morphism:
   From: Closed subscheme of Affine Space of dimension 2 over Rational Field
         defined by: 0
   To:   Affine Space of dimension 2 over Rational Field
   Defn: Defined on coordinates by sending (x, y) to (x, x/y)]
sage: # needs sage.libs.singular
sage: A2.<x,y> = AffineSpace(QQ, 2)
sage: A1.<a> = AffineSpace(QQ, 1)
sage: X = A2.subscheme([x^2 - y^2 - y])
sage: f = X.hom([x/y], A1)
sage: f.representatives()
[Scheme morphism:
   From: Closed subscheme of Affine Space of dimension 2 over Rational Field
         defined by: x^2 - y^2 - y
   To:   Affine Space of dimension 1 over Rational Field
   Defn: Defined on coordinates by sending (x, y) to (x/y),
 Scheme morphism:
   From: Closed subscheme of Affine Space of dimension 2 over Rational Field
         defined by: x^2 - y^2 - y
   To:   Affine Space of dimension 1 over Rational Field
   Defn: Defined on coordinates by sending (x, y) to ((y + 1)/x)]
sage: g = _[1]
sage: g.representatives()
[Scheme morphism:
   From: Closed subscheme of Affine Space of dimension 2 over Rational Field
         defined by: x^2 - y^2 - y
   To:   Affine Space of dimension 1 over Rational Field
   Defn: Defined on coordinates by sending (x, y) to (x/y),
 Scheme morphism:
   From: Closed subscheme of Affine Space of dimension 2 over Rational Field
         defined by: x^2 - y^2 - y
   To:   Affine Space of dimension 1 over Rational Field
   Defn: Defined on coordinates by sending (x, y) to ((y + 1)/x)]
sage: A2.<x,y> = AffineSpace(QQ, 2)
sage: P1.<a,b> = ProjectiveSpace(QQ, 1)
sage: X = A2.subscheme([x^2 - y^2 - y])
sage: f = X.hom([x, y], P1)
sage: f.representatives()                                                   # needs sage.libs.singular
[Scheme morphism:
   From: Closed subscheme of Affine Space of dimension 2 over Rational Field
         defined by: x^2 - y^2 - y
   To:   Projective Space of dimension 1 over Rational Field
   Defn: Defined on coordinates by sending (x, y) to
         (x : y),
 Scheme morphism:
   From: Closed subscheme of Affine Space of dimension 2 over Rational Field
         defined by: x^2 - y^2 - y
   To:   Projective Space of dimension 1 over Rational Field
   Defn: Defined on coordinates by sending (x, y) to (y + 1 : x)]