Morphisms on affine varieties

A morphism of schemes determined by rational functions that define what the morphism does on points in the ambient affine space.

AUTHORS:

  • David Kohel, William Stein
  • Volker Braun (2011-08-08): Renamed classes, more documentation, misc cleanups.
  • Ben Hutz (2013-03) iteration functionality and new directory structure for affine/projective
class sage.schemes.affine.affine_morphism.SchemeMorphism_polynomial_affine_space(parent, polys, check=True)

Bases: sage.schemes.generic.morphism.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())
<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())
<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())
<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)
sage: g = f.as_dynamical_system()
sage: g is f
True
dynatomic_polynomial(period)

Return the dynatomic polynomial.

EXAMPLES:

sage: A.<x> = AffineSpace(QQ, 1)
sage: H = End(A)
sage: f = H([x^2-10/9])
sage: f.dynatomic_polynomial([2, 1])
doctest:warning
...
531441*x^4 - 649539*x^2 - 524880
global_height(prec=None)

Returns the maximum of the heights of the coefficients in any of the coordinate functions of the affine morphism.

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()
8.29404964010203
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=100)
3.3696683136785869233538671082
sage: A.<x> = AffineSpace(ZZ, 1)
sage: H = Hom(A, A)
sage: f = H([7*x^2 + 1513]);
sage: f.global_height()
7.32184971378836
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:

  • SchemeMorphism_polynomial_projective_space.

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: 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)
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: 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])
sage: f.homogenize(2)
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: 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: 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: 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()
[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()
[ 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()
[         (2*x - y)/(y + 1) (-x^2 - x)/(y^2 + 2*y + 1)]
[  (-y - 5)/(x^2 + 4*x + 4)                  1/(x + 2)]
multiplier(P, n, check=True)

Return the multiplier of the point.

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: H = End(A)
sage: f = H([x^2, y^2])
sage: f.multiplier(A([1, 1]), 1)
doctest:warning
...
[2 0]
[0 2]
nth_iterate(P, n)

Return the nth iterate of the point.

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: H = End(A)
sage: f = H([(x-2*y^2)/x, 3*x*y])
sage: f.nth_iterate(A(9, 3), 3)
doctest:warning
...
(-104975/13123, -9566667)
nth_iterate_map(n)

Return the symbolic nth iterate.

EXAMPLES:

sage: A.<x,y> = AffineSpace(ZZ, 2)
sage: H = End(A)
sage: f = H([(x^2-2)/(2*y), y^2-3*x])
sage: f.nth_iterate_map(2)
doctest:warning
...
Dynamical System of Affine Space of dimension 2 over Integer Ring
  Defn: Defined on coordinates by sending (x, y) to
        ((x^4 - 4*x^2 - 8*y^2 + 4)/(8*y^4 - 24*x*y^2), (2*y^5 - 12*x*y^3
+ 18*x^2*y - 3*x^2 + 6)/(2*y))
orbit(P, n)

Return the orbit of the point.

EXAMPLES:

sage: A.<x,y> = AffineSpace(QQ, 2)
sage: H = End(A)
sage: f = H([(x-2*y^2)/x, 3*x*y])
sage: f.orbit(A(9, 3), 3)
doctest:warning
...
[(9, 3), (-1, 81), (13123, -243), (-104975/13123, -9566667)]
class sage.schemes.affine.affine_morphism.SchemeMorphism_polynomial_affine_space_field(parent, polys, check=True)

Bases: sage.schemes.affine.affine_morphism.SchemeMorphism_polynomial_affine_space

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: 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()
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: 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)
True
class sage.schemes.affine.affine_morphism.SchemeMorphism_polynomial_affine_space_finite_field(parent, polys, check=True)

Bases: sage.schemes.affine.affine_morphism.SchemeMorphism_polynomial_affine_space_field

cyclegraph()

Return the directed graph of the map.

EXAMPLES:

sage: A.<x,y> = AffineSpace(GF(5), 2)
sage: H = End(A)
sage: f = H([x^2-y, x*y+1])
sage: f.cyclegraph()
doctest:warning
...
Looped digraph on 25 vertices
orbit_structure(P)

Return the tail and period of the point.

EXAMPLES:

sage: A.<x,y> = AffineSpace(GF(13), 2)
sage: H = End(A)
sage: f = H([x^2 - 1, y^2])
sage: f.orbit_structure(A(2, 3))
doctest:warning
...
[1, 6]