Function Fields#

A function field (of one variable) is a finitely generated field extension of transcendence degree one. In Sage, a function field can be a rational function field or a finite extension of a function field.

EXAMPLES:

We create a rational function field:

sage: K.<x> = FunctionField(GF(5^2,'a')); K
Rational function field in x over Finite Field in a of size 5^2
sage: K.genus()
0
sage: f = (x^2 + x + 1) / (x^3 + 1)
sage: f
(x^2 + x + 1)/(x^3 + 1)
sage: f^3
(x^6 + 3*x^5 + x^4 + 2*x^3 + x^2 + 3*x + 1)/(x^9 + 3*x^6 + 3*x^3 + 1)

Then we create an extension of the rational function field, and do some simple arithmetic in it:

sage: R.<y> = K[]
sage: L.<y> = K.extension(y^3 - (x^3 + 2*x*y + 1/x)); L
Function field in y defined by y^3 + 3*x*y + (4*x^4 + 4)/x
sage: y^2
y^2
sage: y^3
2*x*y + (x^4 + 1)/x
sage: a = 1/y; a
(x/(x^4 + 1))*y^2 + 3*x^2/(x^4 + 1)
sage: a * y
1

We next make an extension of the above function field, illustrating that arithmetic with a tower of three fields is fully supported:

sage: S.<t> = L[]
sage: M.<t> = L.extension(t^2 - x*y)
sage: M
Function field in t defined by t^2 + 4*x*y
sage: t^2
x*y
sage: 1/t
((1/(x^4 + 1))*y^2 + 3*x/(x^4 + 1))*t
sage: M.base_field()
Function field in y defined by y^3 + 3*x*y + (4*x^4 + 4)/x
sage: M.base_field().base_field()
Rational function field in x over Finite Field in a of size 5^2

It is also possible to construct function fields over an imperfect base field:

sage: N.<u> = FunctionField(K)

and inseparable extension function fields:

sage: J.<x> = FunctionField(GF(5)); J
Rational function field in x over Finite Field of size 5
sage: T.<v> = J[]
sage: O.<v> = J.extension(v^5 - x); O
Function field in v defined by v^5 + 4*x

Function fields over the rational field are supported:

sage: F.<x> = FunctionField(QQ)
sage: R.<Y> = F[]
sage: L.<y> = F.extension(Y^2 - x^8 - 1)
sage: O = L.maximal_order()
sage: I = O.ideal(x, y - 1)
sage: P = I.place()
sage: D = P.divisor()
sage: D.basis_function_space()
[1]
sage: (2*D).basis_function_space()
[1]
sage: (3*D).basis_function_space()
[1]
sage: (4*D).basis_function_space()
[1, 1/x^4*y + 1/x^4]

sage: K.<x> = FunctionField(QQ); _.<Y> = K[]
sage: F.<y> = K.extension(Y^3 - x^2*(x^2 + x + 1)^2)
sage: O = F.maximal_order()
sage: I = O.ideal(y)
sage: I.divisor()
2*Place (x, y, (1/(x^3 + x^2 + x))*y^2)
 + 2*Place (x^2 + x + 1, y, (1/(x^3 + x^2 + x))*y^2)

sage: K.<x> = FunctionField(QQ); _.<Y> = K[]
sage: L.<y> = K.extension(Y^2 + Y + x + 1/x)
sage: O = L.maximal_order()
sage: I = O.ideal(y)
sage: I.divisor()
- Place (x, x*y)
 + Place (x^2 + 1, x*y)

Function fields over the algebraic field are supported:

sage: K.<x> = FunctionField(QQbar); _.<Y> = K[]
sage: L.<y> = K.extension(Y^2 + Y + x + 1/x)
sage: O = L.maximal_order()
sage: I = O.ideal(y)
sage: I.divisor()
Place (x - I, x*y)
 - Place (x, x*y)
 + Place (x + I, x*y)
sage: pl = I.divisor().support()[0]
sage: m = L.completion(pl, prec=5)
sage: m(x)
I + s + O(s^5)
sage: m(y)                             # long time (4s)
-2*s + (-4 - I)*s^2 + (-15 - 4*I)*s^3 + (-75 - 23*I)*s^4 + (-413 - 154*I)*s^5 + O(s^6)
sage: m(y)^2 + m(y) + m(x) + 1/m(x)    # long time (8s)
O(s^5)

Global function fields#

A global function field in Sage is an extension field of a rational function field over a finite constant field by an irreducible separable polynomial over the rational function field.

EXAMPLES:

A fundamental computation for a global or any function field is to get a basis of its maximal order and maximal infinite order, and then do arithmetic with ideals of those maximal orders:

sage: K.<x> = FunctionField(GF(3)); _.<t> = K[]
sage: L.<y> = K.extension(t^4 + t - x^5)
sage: O = L.maximal_order()
sage: O.basis()
(1, y, 1/x*y^2 + 1/x*y, 1/x^3*y^3 + 2/x^3*y^2 + 1/x^3*y)
sage: I = O.ideal(x,y); I
Ideal (x, y) of Maximal order of Function field in y defined by y^4 + y + 2*x^5
sage: J = I^-1
sage: J.basis_matrix()
[  1   0   0   0]
[1/x 1/x   0   0]
[  0   0   1   0]
[  0   0   0   1]
sage: L.maximal_order_infinite().basis()
(1, 1/x^2*y, 1/x^3*y^2, 1/x^4*y^3)

As an example of the most sophisticated computations that Sage can do with a global function field, we compute all the Weierstrass places of the Klein quartic over \(\GF{2}\) and gap numbers for ordinary places:

sage: K.<x> = FunctionField(GF(2)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^3 + x^3*Y + x)
sage: L.genus()
3
sage: L.weierstrass_places()
[Place (1/x, 1/x^3*y^2 + 1/x),
 Place (1/x, 1/x^3*y^2 + 1/x^2*y + 1),
 Place (x, y),
 Place (x + 1, (x^3 + 1)*y + x + 1),
 Place (x^3 + x + 1, y + 1),
 Place (x^3 + x + 1, y + x^2),
 Place (x^3 + x + 1, y + x^2 + 1),
 Place (x^3 + x^2 + 1, y + x),
 Place (x^3 + x^2 + 1, y + x^2 + 1),
 Place (x^3 + x^2 + 1, y + x^2 + x + 1)]
sage: L.gaps()
[1, 2, 3]

The gap numbers for Weierstrass places are of course not ordinary:

sage: p1,p2,p3 = L.weierstrass_places()[:3]
sage: p1.gaps()
[1, 2, 4]
sage: p2.gaps()
[1, 2, 4]
sage: p3.gaps()
[1, 2, 4]

AUTHORS:

  • William Stein (2010): initial version

  • Robert Bradshaw (2010-05-30): added is_finite()

  • Julian Rüth (2011-06-08, 2011-09-14, 2014-06-23, 2014-06-24, 2016-11-13): fixed hom(), extension(); use @cached_method; added derivation(); added support for relative vector spaces; fixed conversion to base fields

  • Maarten Derickx (2011-09-11): added doctests

  • Syed Ahmad Lavasani (2011-12-16): added genus(), is_RationalFunctionField()

  • Simon King (2014-10-29): Use the same generator names for a function field extension and the underlying polynomial ring.

  • Kwankyu Lee (2017-04-30): added global function fields

  • Brent Baccala (2019-12-20): added function fields over number fields and QQbar

class sage.rings.function_field.function_field.FunctionField(base_field, names, category=Category of function fields)#

Bases: Field

Abstract base class for all function fields.

INPUT:

  • base_field – field; the base of this function field

  • names – string that gives the name of the generator

EXAMPLES:

sage: K.<x> = FunctionField(QQ)
sage: K
Rational function field in x over Rational Field
basis_of_differentials_of_first_kind()#

Return a basis of the space of holomorphic differentials of this function field.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.basis_of_holomorphic_differentials()
[]

sage: K.<x> = FunctionField(GF(5)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^3 - (x^3 - 1)/(x^3 - 2))
sage: L.basis_of_holomorphic_differentials()
[((x/(x^3 + 4))*y) d(x),
 ((1/(x^3 + 4))*y) d(x),
 ((x/(x^3 + 4))*y^2) d(x),
 ((1/(x^3 + 4))*y^2) d(x)]
basis_of_holomorphic_differentials()#

Return a basis of the space of holomorphic differentials of this function field.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.basis_of_holomorphic_differentials()
[]

sage: K.<x> = FunctionField(GF(5)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^3 - (x^3 - 1)/(x^3 - 2))
sage: L.basis_of_holomorphic_differentials()
[((x/(x^3 + 4))*y) d(x),
 ((1/(x^3 + 4))*y) d(x),
 ((x/(x^3 + 4))*y^2) d(x),
 ((1/(x^3 + 4))*y^2) d(x)]
characteristic()#

Return the characteristic of the function field.

EXAMPLES:

sage: K.<x> = FunctionField(QQ)
sage: K.characteristic()
0
sage: K.<x> = FunctionField(QQbar)
sage: K.characteristic()
0
sage: K.<x> = FunctionField(GF(7))
sage: K.characteristic()
7
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^2 - x)
sage: L.characteristic()
7
completion(place, name=None, prec=None, gen_name=None)#

Return the completion of the function field at the place.

INPUT:

  • place – place

  • name – string; name of the series variable

  • prec – positive integer; default precision

  • gen_name – string; name of the generator of the residue field; used only when the place is non-rational

EXAMPLES:

sage: K.<x> = FunctionField(GF(2)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^2 + Y + x + 1/x)
sage: p = L.places_finite()[0]
sage: m = L.completion(p); m
Completion map:
  From: Function field in y defined by y^2 + y + (x^2 + 1)/x
  To:   Laurent Series Ring in s over Finite Field of size 2
sage: m(x,10)
s^2 + s^3 + s^4 + s^5 + s^7 + s^8 + s^9 + s^10 + O(s^12)
sage: m(y,10)
s^-1 + 1 + s^3 + s^5 + s^7 + O(s^9)

sage: K.<x> = FunctionField(GF(2)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^2 + Y + x + 1/x)
sage: p = L.places_finite()[0]
sage: m = L.completion(p); m
Completion map:
  From: Function field in y defined by y^2 + y + (x^2 + 1)/x
  To:   Laurent Series Ring in s over Finite Field of size 2
sage: m(x,10)
s^2 + s^3 + s^4 + s^5 + s^7 + s^8 + s^9 + s^10 + O(s^12)
sage: m(y,10)
s^-1 + 1 + s^3 + s^5 + s^7 + O(s^9)

sage: K.<x> = FunctionField(GF(2))
sage: p = K.places_finite()[0]; p
Place (x)
sage: m = K.completion(p); m
Completion map:
  From: Rational function field in x over Finite Field of size 2
  To:   Laurent Series Ring in s over Finite Field of size 2
sage: m(1/(x+1))
1 + s + s^2 + s^3 + s^4 + s^5 + s^6 + s^7 + s^8 + s^9 + s^10 + s^11 + s^12
+ s^13 + s^14 + s^15 + s^16 + s^17 + s^18 + s^19 + O(s^20)

sage: p = K.place_infinite(); p
Place (1/x)
sage: m = K.completion(p); m
Completion map:
  From: Rational function field in x over Finite Field of size 2
  To:   Laurent Series Ring in s over Finite Field of size 2
sage: m(x)
s^-1 + O(s^19)

sage: m = K.completion(p, prec=infinity); m
Completion map:
  From: Rational function field in x over Finite Field of size 2
  To:   Lazy Laurent Series Ring in s over Finite Field of size 2
sage: f = m(x); f
s^-1 + ...
sage: f.coefficient(100)
0

sage: K.<x> = FunctionField(QQ); _.<Y> = K[]
sage: L.<y> = K.extension(Y^2 - x)
sage: O = L.maximal_order()
sage: decomp = O.decomposition(K.maximal_order().ideal(x - 1))
sage: pls = (decomp[0][0].place(), decomp[1][0].place())
sage: m = L.completion(pls[0]); m
Completion map:
  From: Function field in y defined by y^2 - x
  To:   Laurent Series Ring in s over Rational Field
sage: xe = m(x)
sage: ye = m(y)
sage: ye^2 - xe == 0
True

sage: decomp2 = O.decomposition(K.maximal_order().ideal(x^2 + 1))
sage: pls2 = decomp2[0][0].place()
sage: m = L.completion(pls2); m
Completion map:
  From: Function field in y defined by y^2 - x
  To:   Laurent Series Ring in s over Number Field in a with defining polynomial x^4 + 2*x^2 + 4*x + 2
sage: xe = m(x)
sage: ye = m(y)
sage: ye^2 - xe == 0
True
divisor_group()#

Return the group of divisors attached to the function field.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.divisor_group()
Divisor group of Rational function field in t over Rational Field

sage: _.<Y> = K[]
sage: L.<y> = K.extension(Y^3 - (t^3 - 1)/(t^3 - 2))
sage: L.divisor_group()
Divisor group of Function field in y defined by y^3 + (-t^3 + 1)/(t^3 - 2)

sage: K.<x> = FunctionField(GF(5)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^3 - (x^3 - 1)/(x^3 - 2))
sage: L.divisor_group()
Divisor group of Function field in y defined by y^3 + (4*x^3 + 1)/(x^3 + 3)
extension(f, names=None)#

Create an extension \(K(y)\) of this function field \(K\) extended with a root \(y\) of the univariate polynomial \(f\) over \(K\).

INPUT:

  • f – univariate polynomial over \(K\)

  • names – string or tuple of length 1 that names the variable \(y\)

OUTPUT:

  • a function field

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: K.extension(y^5 - x^3 - 3*x + x*y)
Function field in y defined by y^5 + x*y - x^3 - 3*x

A nonintegral defining polynomial:

sage: K.<t> = FunctionField(QQ); R.<y> = K[]
sage: K.extension(y^3 + (1/t)*y + t^3/(t+1), 'z')
Function field in z defined by z^3 + 1/t*z + t^3/(t + 1)

The defining polynomial need not be monic or integral:

sage: K.extension(t*y^3 + (1/t)*y + t^3/(t+1))
Function field in y defined by t*y^3 + 1/t*y + t^3/(t + 1)
extension_constant_field(k)#

Return the constant field extension with constant field \(k\).

INPUT:

  • k – an extension field of the constant field of this function field

EXAMPLES:

sage: K.<x> = FunctionField(GF(2)); _.<Y> = K[]
sage: F.<y> = K.extension(Y^2 + Y + x + 1/x)
sage: E = F.extension_constant_field(GF(2^4))
sage: E
Function field in y defined by y^2 + y + (x^2 + 1)/x over its base
sage: E.constant_base_field()
Finite Field in z4 of size 2^4
is_finite()#

Return whether the function field is finite, which is false.

EXAMPLES:

sage: R.<t> = FunctionField(QQ)
sage: R.is_finite()
False
sage: R.<t> = FunctionField(GF(7))
sage: R.is_finite()
False
is_global()#

Return whether the function field is global, that is, whether the constant field is finite.

EXAMPLES:

sage: R.<t> = FunctionField(QQ)
sage: R.is_global()
False
sage: R.<t> = FunctionField(QQbar)
sage: R.is_global()
False
sage: R.<t> = FunctionField(GF(7))
sage: R.is_global()
True
is_perfect()#

Return whether the field is perfect, i.e., its characteristic \(p\) is zero or every element has a \(p\)-th root.

EXAMPLES:

sage: FunctionField(QQ, 'x').is_perfect()
True
sage: FunctionField(GF(2), 'x').is_perfect()
False
order(x, check=True)#

Return the order generated by x over the base maximal order.

INPUT:

  • x – element or list of elements of the function field

  • check – boolean; if True, check that x really generates an order

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^3 + x^3 + 4*x + 1)
sage: O = L.order(y); O
Order in Function field in y defined by y^3 + x^3 + 4*x + 1
sage: O.basis()
(1, y, y^2)

sage: Z = K.order(x); Z
Order in Rational function field in x over Rational Field
sage: Z.basis()
(1,)

Orders with multiple generators are not yet supported:

sage: Z = K.order([x,x^2]); Z
Traceback (most recent call last):
...
NotImplementedError
order_infinite(x, check=True)#

Return the order generated by x over the maximal infinite order.

INPUT:

  • x – element or a list of elements of the function field

  • check – boolean; if True, check that x really generates an order

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^3 + x^3 + 4*x + 1)
sage: L.order_infinite(y)  # todo: not implemented

sage: Z = K.order(x); Z
Order in Rational function field in x over Rational Field
sage: Z.basis()
(1,)

Orders with multiple generators, not yet supported:

sage: Z = K.order_infinite([x,x^2]); Z
Traceback (most recent call last):
...
NotImplementedError
order_infinite_with_basis(basis, check=True)#

Return the order with given basis over the maximal infinite order of the base field.

INPUT:

  • basis – list of elements of the function field

  • check – boolean (default: True); if True, check that the basis is really linearly independent and that the module it spans is closed under multiplication, and contains the identity element.

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^3 + x^3 + 4*x + 1)
sage: O = L.order_infinite_with_basis([1, 1/x*y, 1/x^2*y^2]); O
Infinite order in Function field in y defined by y^3 + x^3 + 4*x + 1
sage: O.basis()
(1, 1/x*y, 1/x^2*y^2)

Note that 1 does not need to be an element of the basis, as long it is in the module spanned by it:

sage: O = L.order_infinite_with_basis([1+1/x*y,1/x*y, 1/x^2*y^2]); O
Infinite order in Function field in y defined by y^3 + x^3 + 4*x + 1
sage: O.basis()
(1/x*y + 1, 1/x*y, 1/x^2*y^2)

The following error is raised when the module spanned by the basis is not closed under multiplication:

sage: O = L.order_infinite_with_basis([1,y, 1/x^2*y^2]); O
Traceback (most recent call last):
...
ValueError: the module generated by basis (1, y, 1/x^2*y^2) must be closed under multiplication

and this happens when the identity is not in the module spanned by the basis:

sage: O = L.order_infinite_with_basis([1/x,1/x*y, 1/x^2*y^2])
Traceback (most recent call last):
...
ValueError: the identity element must be in the module spanned by basis (1/x, 1/x*y, 1/x^2*y^2)
order_with_basis(basis, check=True)#

Return the order with given basis over the maximal order of the base field.

INPUT:

  • basis – list of elements of this function field

  • check – boolean (default: True); if True, check that the basis is really linearly independent and that the module it spans is closed under multiplication, and contains the identity element.

OUTPUT:

  • an order in the function field

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]; L.<y> = K.extension(y^3 + x^3 + 4*x + 1)
sage: O = L.order_with_basis([1, y, y^2]); O
Order in Function field in y defined by y^3 + x^3 + 4*x + 1
sage: O.basis()
(1, y, y^2)

Note that 1 does not need to be an element of the basis, as long it is in the module spanned by it:

sage: O = L.order_with_basis([1+y, y, y^2]); O
Order in Function field in y defined by y^3 + x^3 + 4*x + 1
sage: O.basis()
(y + 1, y, y^2)

The following error is raised when the module spanned by the basis is not closed under multiplication:

sage: O = L.order_with_basis([1, x^2 + x*y, (2/3)*y^2]); O
Traceback (most recent call last):
...
ValueError: the module generated by basis (1, x*y + x^2, 2/3*y^2) must be closed under multiplication

and this happens when the identity is not in the module spanned by the basis:

sage: O = L.order_with_basis([x, x^2 + x*y, (2/3)*y^2])
Traceback (most recent call last):
...
ValueError: the identity element must be in the module spanned by basis (x, x*y + x^2, 2/3*y^2)
place_set()#

Return the set of all places of the function field.

EXAMPLES:

sage: K.<t> = FunctionField(GF(7))
sage: K.place_set()
Set of places of Rational function field in t over Finite Field of size 7

sage: K.<t> = FunctionField(QQ)
sage: K.place_set()
Set of places of Rational function field in t over Rational Field

sage: K.<x> = FunctionField(GF(2)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^2 + Y + x + 1/x)
sage: L.place_set()
Set of places of Function field in y defined by y^2 + y + (x^2 + 1)/x
rational_function_field()#

Return the rational function field from which this field has been created as an extension.

EXAMPLES:

sage: K.<x> = FunctionField(QQ)
sage: K.rational_function_field()
Rational function field in x over Rational Field

sage: R.<y> = K[]
sage: L.<y> = K.extension(y^2-x)
sage: L.rational_function_field()
Rational function field in x over Rational Field

sage: R.<z> = L[]
sage: M.<z> = L.extension(z^2-y)
sage: M.rational_function_field()
Rational function field in x over Rational Field
some_elements()#

Return some elements in this function field.

EXAMPLES:

sage: K.<x> = FunctionField(QQ)
sage: K.some_elements()
[1,
 x,
 2*x,
 x/(x^2 + 2*x + 1),
 1/x^2,
 x/(x^2 - 1),
 x/(x^2 + 1),
 1/2*x/(x^2 + 1),
 0,
 1/x,
 ...]
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^2 - x)
sage: L.some_elements()
[1,
 y,
 1/x*y,
 ((x + 1)/(x^2 - 2*x + 1))*y - 2*x/(x^2 - 2*x + 1),
 1/x,
 (1/(x - 1))*y,
 (1/(x + 1))*y,
 (1/2/(x + 1))*y,
 0,
 ...]
space_of_differentials()#

Return the space of differentials attached to the function field.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.space_of_differentials()
Space of differentials of Rational function field in t over Rational Field

sage: K.<x> = FunctionField(GF(5)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^3 - (x^3 - 1)/(x^3 - 2))
sage: L.space_of_differentials()
Space of differentials of Function field in y defined by y^3 + (4*x^3 + 1)/(x^3 + 3)
space_of_differentials_of_first_kind()#

Return the space of holomorphic differentials of this function field.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.space_of_holomorphic_differentials()
(Vector space of dimension 0 over Rational Field,
 Linear map:
   From: Vector space of dimension 0 over Rational Field
   To:   Space of differentials of Rational function field in t over Rational Field,
 Section of linear map:
   From: Space of differentials of Rational function field in t over Rational Field
   To:   Vector space of dimension 0 over Rational Field)

sage: K.<x> = FunctionField(GF(5)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^3 - (x^3 - 1)/(x^3 - 2))
sage: L.space_of_holomorphic_differentials()
(Vector space of dimension 4 over Finite Field of size 5,
 Linear map:
   From: Vector space of dimension 4 over Finite Field of size 5
   To:   Space of differentials of Function field in y defined by y^3 + (4*x^3 + 1)/(x^3 + 3),
 Section of linear map:
   From: Space of differentials of Function field in y defined by y^3 + (4*x^3 + 1)/(x^3 + 3)
   To:   Vector space of dimension 4 over Finite Field of size 5)
space_of_holomorphic_differentials()#

Return the space of holomorphic differentials of this function field.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.space_of_holomorphic_differentials()
(Vector space of dimension 0 over Rational Field,
 Linear map:
   From: Vector space of dimension 0 over Rational Field
   To:   Space of differentials of Rational function field in t over Rational Field,
 Section of linear map:
   From: Space of differentials of Rational function field in t over Rational Field
   To:   Vector space of dimension 0 over Rational Field)

sage: K.<x> = FunctionField(GF(5)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^3 - (x^3 - 1)/(x^3 - 2))
sage: L.space_of_holomorphic_differentials()
(Vector space of dimension 4 over Finite Field of size 5,
 Linear map:
   From: Vector space of dimension 4 over Finite Field of size 5
   To:   Space of differentials of Function field in y defined by y^3 + (4*x^3 + 1)/(x^3 + 3),
 Section of linear map:
   From: Space of differentials of Function field in y defined by y^3 + (4*x^3 + 1)/(x^3 + 3)
   To:   Vector space of dimension 4 over Finite Field of size 5)
valuation(prime)#

Return the discrete valuation on this function field defined by prime.

INPUT:

  • prime – a place of the function field, a valuation on a subring, or a valuation on another function field together with information for isomorphisms to and from that function field

EXAMPLES:

We create valuations that correspond to finite rational places of a function field:

sage: K.<x> = FunctionField(QQ)
sage: v = K.valuation(1); v
(x - 1)-adic valuation
sage: v(x)
0
sage: v(x - 1)
1

A place can also be specified with an irreducible polynomial:

sage: v = K.valuation(x - 1); v
(x - 1)-adic valuation

Similarly, for a finite non-rational place:

sage: v = K.valuation(x^2 + 1); v
(x^2 + 1)-adic valuation
sage: v(x^2 + 1)
1
sage: v(x)
0

Or for the infinite place:

sage: v = K.valuation(1/x); v
Valuation at the infinite place
sage: v(x)
-1

Instead of specifying a generator of a place, we can define a valuation on a rational function field by giving a discrete valuation on the underlying polynomial ring:

sage: R.<x> = QQ[]
sage: w = valuations.GaussValuation(R, valuations.TrivialValuation(QQ)).augmentation(x - 1, 1)
sage: v = K.valuation(w); v
(x - 1)-adic valuation

Note that this allows us to specify valuations which do not correspond to a place of the function field:

sage: w = valuations.GaussValuation(R, QQ.valuation(2))
sage: v = K.valuation(w); v
2-adic valuation

The same is possible for valuations with \(v(1/x) > 0\) by passing in an extra pair of parameters, an isomorphism between this function field and an isomorphic function field. That way you can, for example, indicate that the valuation is to be understood as a valuation on \(K[1/x]\), i.e., after applying the substitution \(x \mapsto 1/x\) (here, the inverse map is also \(x \mapsto 1/x\)):

sage: w = valuations.GaussValuation(R, QQ.valuation(2)).augmentation(x, 1)
sage: w = K.valuation(w)
sage: v = K.valuation((w, K.hom([~K.gen()]), K.hom([~K.gen()]))); v
Valuation on rational function field induced by [ Gauss valuation induced by 2-adic valuation, v(x) = 1 ] (in Rational function field in x over Rational Field after x |--> 1/x)

Note that classical valuations at finite places or the infinite place are always normalized such that the uniformizing element has valuation 1:

sage: K.<t> = FunctionField(GF(3))
sage: M.<x> = FunctionField(K)
sage: v = M.valuation(x^3 - t)
sage: v(x^3 - t)
1

However, if such a valuation comes out of a base change of the ground field, this is not the case anymore. In the example below, the unique extension of v to L still has valuation 1 on \(x^3 - t\) but it has valuation 1/3 on its uniformizing element \(x - w\):

sage: R.<w> = K[]
sage: L.<w> = K.extension(w^3 - t)
sage: N.<x> = FunctionField(L)
sage: w = v.extension(N) # missing factorization, :trac:`16572`
Traceback (most recent call last):
...
NotImplementedError
sage: w(x^3 - t) # not tested
1
sage: w(x - w) # not tested
1/3

There are several ways to create valuations on extensions of rational function fields:

sage: K.<x> = FunctionField(QQ)
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^2 - x); L
Function field in y defined by y^2 - x

A place that has a unique extension can just be defined downstairs:

sage: v = L.valuation(x); v
(x)-adic valuation
class sage.rings.function_field.function_field.FunctionField_char_zero(polynomial, names, category=None)#

Bases: FunctionField_simple

Function fields of characteristic zero.

EXAMPLES:

sage: K.<x> = FunctionField(QQ); _.<Y> = K[]
sage: L.<y> = K.extension(Y^3 - (x^3 - 1)/(x^3 - 2))
sage: L
Function field in y defined by y^3 + (-x^3 + 1)/(x^3 - 2)
sage: L.characteristic()
0
higher_derivation()#

Return the higher derivation (also called the Hasse-Schmidt derivation) for the function field.

The higher derivation of the function field is uniquely determined with respect to the separating element \(x\) of the base rational function field \(k(x)\).

EXAMPLES:

sage: K.<x> = FunctionField(QQ); _.<Y> = K[]
sage: L.<y> = K.extension(Y^3 - (x^3 - 1)/(x^3 - 2))
sage: L.higher_derivation()
Higher derivation map:
  From: Function field in y defined by y^3 + (-x^3 + 1)/(x^3 - 2)
  To:   Function field in y defined by y^3 + (-x^3 + 1)/(x^3 - 2)
class sage.rings.function_field.function_field.FunctionField_char_zero_integral(polynomial, names, category=None)#

Bases: FunctionField_char_zero, FunctionField_integral

Function fields of characteristic zero, defined by an irreducible and separable polynomial, integral over the maximal order of the base rational function field with a finite constant field.

class sage.rings.function_field.function_field.FunctionField_global(polynomial, names)#

Bases: FunctionField_simple

Global function fields.

INPUT:

  • polynomial – monic irreducible and separable polynomial

  • names – name of the generator of the function field

EXAMPLES:

sage: K.<x> = FunctionField(GF(5)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^3 - (x^3 - 1)/(x^3 - 2))
sage: L
Function field in y defined by y^3 + (4*x^3 + 1)/(x^3 + 3)

The defining equation needs not be monic:

sage: K.<x> = FunctionField(GF(4)); _.<Y> = K[]
sage: L.<y> = K.extension((1 - x)*Y^7 - x^3)
sage: L.gaps()                         # long time (6s)
[1, 2, 3]

or may define a trivial extension:

sage: K.<x> = FunctionField(GF(5)); _.<Y> = K[]
sage: L.<y> = K.extension(Y-1)
sage: L.genus()
0
L_polynomial(name='t')#

Return the L-polynomial of the function field.

INPUT:

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

EXAMPLES:

sage: K.<x> = FunctionField(GF(2)); _.<Y> = K[]
sage: F.<y> = K.extension(Y^2 + Y + x + 1/x)
sage: F.L_polynomial()
2*t^2 + t + 1
gaps()#

Return the gaps of the function field.

These are the gaps at the ordinary places, that is, places which are not Weierstrass places.

EXAMPLES:

sage: K.<x> = FunctionField(GF(2)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^3 + x^3 * Y + x)
sage: L.gaps()
[1, 2, 3]
get_place(degree)#

Return a place of degree.

INPUT:

  • degree – a positive integer

OUTPUT: a place of degree if any exists; otherwise None

EXAMPLES:

sage: F.<a> = GF(2)
sage: K.<x> = FunctionField(F)
sage: R.<Y> = PolynomialRing(K)
sage: L.<y> = K.extension(Y^4 + Y - x^5)
sage: L.get_place(1)
Place (x, y)
sage: L.get_place(2)
Place (x, y^2 + y + 1)
sage: L.get_place(3)
Place (x^3 + x^2 + 1, y + x^2 + x)
sage: L.get_place(4)
Place (x + 1, x^5 + 1)
sage: L.get_place(5)
Place (x^5 + x^3 + x^2 + x + 1, y + x^4 + 1)
sage: L.get_place(6)
Place (x^3 + x^2 + 1, y^2 + y + x^2)
sage: L.get_place(7)
Place (x^7 + x + 1, y + x^6 + x^5 + x^4 + x^3 + x)
sage: L.get_place(8)
higher_derivation()#

Return the higher derivation (also called the Hasse-Schmidt derivation) for the function field.

The higher derivation of the function field is uniquely determined with respect to the separating element \(x\) of the base rational function field \(k(x)\).

EXAMPLES:

sage: K.<x>=FunctionField(GF(5)); _.<Y>=K[]
sage: L.<y>=K.extension(Y^3 - (x^3 - 1)/(x^3 - 2))
sage: L.higher_derivation()
Higher derivation map:
  From: Function field in y defined by y^3 + (4*x^3 + 1)/(x^3 + 3)
  To:   Function field in y defined by y^3 + (4*x^3 + 1)/(x^3 + 3)
maximal_order()#

Return the maximal order of the function field.

EXAMPLES:

sage: K.<x> = FunctionField(GF(2));
sage: R.<t> = PolynomialRing(K);
sage: F.<y> = K.extension(t^4 + x^12*t^2 + x^18*t + x^21 + x^18);
sage: O = F.maximal_order()
sage: O.basis()
(1, 1/x^4*y, 1/x^11*y^2 + 1/x^2, 1/x^15*y^3 + 1/x^6*y)
number_of_rational_places(r=1)#

Return the number of rational places of the function field whose constant field extended by degree r.

INPUT:

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

EXAMPLES:

sage: K.<x> = FunctionField(GF(2)); _.<Y> = K[]
sage: F.<y> = K.extension(Y^2 + Y + x + 1/x)
sage: F.number_of_rational_places()
4
sage: [F.number_of_rational_places(r) for r in [1..10]]
[4, 8, 4, 16, 44, 56, 116, 288, 508, 968]
places(degree=1)#

Return a list of the places with degree.

INPUT:

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

EXAMPLES:

sage: F.<a> = GF(2)
sage: K.<x> = FunctionField(F)
sage: R.<t> = PolynomialRing(K)
sage: L.<y> = K.extension(t^4 + t - x^5)
sage: L.places(1)
[Place (1/x, 1/x^4*y^3), Place (x, y), Place (x, y + 1)]
places_finite(degree=1)#

Return a list of the finite places with degree.

INPUT:

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

EXAMPLES:

sage: F.<a> = GF(2)
sage: K.<x> = FunctionField(F)
sage: R.<t> = PolynomialRing(K)
sage: L.<y> = K.extension(t^4+t-x^5)
sage: L.places_finite(1)
[Place (x, y), Place (x, y + 1)]
places_infinite(degree=1)#

Return a list of the infinite places with degree.

INPUT:

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

EXAMPLES:

sage: F.<a> = GF(2)
sage: K.<x> = FunctionField(F)
sage: R.<t> = PolynomialRing(K)
sage: L.<y> = K.extension(t^4+t-x^5)
sage: L.places_infinite(1)
[Place (1/x, 1/x^4*y^3)]
weierstrass_places()#

Return all Weierstrass places of the function field.

EXAMPLES:

sage: K.<x> = FunctionField(GF(2)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^3 + x^3 * Y + x)
sage: L.weierstrass_places()
[Place (1/x, 1/x^3*y^2 + 1/x),
 Place (1/x, 1/x^3*y^2 + 1/x^2*y + 1),
 Place (x, y),
 Place (x + 1, (x^3 + 1)*y + x + 1),
 Place (x^3 + x + 1, y + 1),
 Place (x^3 + x + 1, y + x^2),
 Place (x^3 + x + 1, y + x^2 + 1),
 Place (x^3 + x^2 + 1, y + x),
 Place (x^3 + x^2 + 1, y + x^2 + 1),
 Place (x^3 + x^2 + 1, y + x^2 + x + 1)]
class sage.rings.function_field.function_field.FunctionField_global_integral(polynomial, names)#

Bases: FunctionField_global, FunctionField_integral

Global function fields, defined by an irreducible and separable polynomial, integral over the maximal order of the base rational function field with a finite constant field.

class sage.rings.function_field.function_field.FunctionField_integral(polynomial, names, category=None)#

Bases: FunctionField_simple

Integral function fields.

A function field is integral if it is defined by an irreducible separable polynomial, which is integral over the maximal order of the base rational function field.

equation_order()#

Return the equation order of the function field.

EXAMPLES:

sage: K.<x> = FunctionField(GF(2)); R.<t> = PolynomialRing(K)
sage: F.<y> = K.extension(t^3-x^2*(x^2+x+1)^2)
sage: F.equation_order()
Order in Function field in y defined by y^3 + x^6 + x^4 + x^2

sage: K.<x> = FunctionField(QQ); R.<t> = PolynomialRing(K)
sage: F.<y> = K.extension(t^3-x^2*(x^2+x+1)^2)
sage: F.equation_order()
Order in Function field in y defined by y^3 - x^6 - 2*x^5 - 3*x^4 - 2*x^3 - x^2
equation_order_infinite()#

Return the infinite equation order of the function field.

This is by definition \(o[b]\) where \(b\) is the primitive integral element from primitive_integal_element_infinite() and \(o\) is the maximal infinite order of the base rational function field.

EXAMPLES:

sage: K.<x> = FunctionField(GF(2)); R.<t> = PolynomialRing(K)
sage: F.<y> = K.extension(t^3-x^2*(x^2+x+1)^2)
sage: F.equation_order_infinite()
Infinite order in Function field in y defined by y^3 + x^6 + x^4 + x^2

sage: K.<x> = FunctionField(QQ); R.<t> = PolynomialRing(K)
sage: F.<y> = K.extension(t^3-x^2*(x^2+x+1)^2)
sage: F.equation_order_infinite()
Infinite order in Function field in y defined by y^3 - x^6 - 2*x^5 - 3*x^4 - 2*x^3 - x^2
primitive_integal_element_infinite()#

Return a primitive integral element over the base maximal infinite order.

This element is integral over the maximal infinite order of the base rational function field and the function field is a simple extension by this element over the base order.

EXAMPLES:

sage: K.<x> = FunctionField(GF(2)); R.<t> = PolynomialRing(K)
sage: F.<y> = K.extension(t^3-x^2*(x^2+x+1)^2)
sage: b = F.primitive_integal_element_infinite(); b
1/x^2*y
sage: b.minimal_polynomial('t')
t^3 + (x^4 + x^2 + 1)/x^4
class sage.rings.function_field.function_field.FunctionField_polymod(polynomial, names, category=None)#

Bases: FunctionField

Function fields defined by a univariate polynomial, as an extension of the base field.

INPUT:

  • polynomial – univariate polynomial over a function field

  • names – tuple of length 1 or string; variable names

  • category – category (default: category of function fields)

EXAMPLES:

We make a function field defined by a degree 5 polynomial over the rational function field over the rational numbers:

sage: K.<x> = FunctionField(QQ)
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x)); L
Function field in y defined by y^5 - 2*x*y + (-x^4 - 1)/x

We next make a function field over the above nontrivial function field L:

sage: S.<z> = L[]
sage: M.<z> = L.extension(z^2 + y*z + y); M
Function field in z defined by z^2 + y*z + y
sage: 1/z
((-x/(x^4 + 1))*y^4 + 2*x^2/(x^4 + 1))*z - 1
sage: z * (1/z)
1

We drill down the tower of function fields:

sage: M.base_field()
Function field in y defined by y^5 - 2*x*y + (-x^4 - 1)/x
sage: M.base_field().base_field()
Rational function field in x over Rational Field
sage: M.base_field().base_field().constant_field()
Rational Field
sage: M.constant_base_field()
Rational Field

Warning

It is not checked if the polynomial used to define the function field is irreducible Hence it is not guaranteed that this object really is a field! This is illustrated below.

sage: K.<x> = FunctionField(QQ)
sage: R.<y> = K[]
sage: L.<y> = K.extension(x^2 - y^2)
sage: (y - x)*(y + x)
0
sage: 1/(y - x)
1
sage: y - x == 0; y + x == 0
False
False
Element#

alias of FunctionFieldElement_polymod

base_field()#

Return the base field of the function field. This function field is presented as \(L = K[y]/(f(y))\), and the base field is by definition the field \(K\).

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x))
sage: L.base_field()
Rational function field in x over Rational Field
change_variable_name(name)#

Return a field isomorphic to this field with variable(s) name.

INPUT:

  • name – a string or a tuple consisting of a strings, the names of the new variables starting with a generator of this field and going down to the rational function field.

OUTPUT:

A triple F,f,t where F is a function field, f is an isomorphism from F to this field, and t is the inverse of f.

EXAMPLES:

sage: K.<x> = FunctionField(QQ)
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^2 - x)
sage: R.<z> = L[]
sage: M.<z> = L.extension(z^2 - y)

sage: M.change_variable_name('zz')
(Function field in zz defined by zz^2 - y,
 Function Field morphism:
  From: Function field in zz defined by zz^2 - y
  To:   Function field in z defined by z^2 - y
  Defn: zz |--> z
        y |--> y
        x |--> x,
 Function Field morphism:
  From: Function field in z defined by z^2 - y
  To:   Function field in zz defined by zz^2 - y
  Defn: z |--> zz
        y |--> y
        x |--> x)
sage: M.change_variable_name(('zz','yy'))
(Function field in zz defined by zz^2 - yy, Function Field morphism:
  From: Function field in zz defined by zz^2 - yy
  To:   Function field in z defined by z^2 - y
  Defn: zz |--> z
        yy |--> y
        x |--> x, Function Field morphism:
  From: Function field in z defined by z^2 - y
  To:   Function field in zz defined by zz^2 - yy
  Defn: z |--> zz
        y |--> yy
        x |--> x)
sage: M.change_variable_name(('zz','yy','xx'))
(Function field in zz defined by zz^2 - yy,
 Function Field morphism:
  From: Function field in zz defined by zz^2 - yy
  To:   Function field in z defined by z^2 - y
  Defn: zz |--> z
        yy |--> y
        xx |--> x,
 Function Field morphism:
  From: Function field in z defined by z^2 - y
  To:   Function field in zz defined by zz^2 - yy
  Defn: z |--> zz
        y |--> yy
        x |--> xx)
constant_base_field()#

Return the base constant field of the function field.

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x)); L
Function field in y defined by y^5 - 2*x*y + (-x^4 - 1)/x
sage: L.constant_base_field()
Rational Field
sage: S.<z> = L[]
sage: M.<z> = L.extension(z^2 - y)
sage: M.constant_base_field()
Rational Field
constant_field()#

Return the algebraic closure of the constant field of the function field.

EXAMPLES:

sage: K.<x> = FunctionField(GF(5)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^5 - x)
sage: L.constant_field()
Traceback (most recent call last):
...
NotImplementedError
degree(base=None)#

Return the degree of the function field over the function field base.

INPUT:

  • base – a function field (default: None), a function field from which this field has been constructed as a finite extension.

EXAMPLES:

sage: K.<x> = FunctionField(QQ)
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x)); L
Function field in y defined by y^5 - 2*x*y + (-x^4 - 1)/x
sage: L.degree()
5
sage: L.degree(L)
1

sage: R.<z> = L[]
sage: M.<z> = L.extension(z^2 - y)
sage: M.degree(L)
2
sage: M.degree(K)
10
different()#

Return the different of the function field.

EXAMPLES:

sage: K.<x> = FunctionField(GF(2)); _.<Y> = K[]
sage: F.<y> = K.extension(Y^3 - x^2*(x^2 + x + 1)^2)
sage: F.different()
2*Place (x, (1/(x^3 + x^2 + x))*y^2)
 + 2*Place (x^2 + x + 1, (1/(x^3 + x^2 + x))*y^2)
equation_order()#

Return the equation order of the function field.

If we view the function field as being presented as \(K[y]/(f(y))\), then the order generated by the class of \(y\) is returned. If \(f\) is not monic, then _make_monic_integral() is called, and instead we get the order generated by some integral multiple of a root of \(f\).

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x))
sage: O = L.equation_order()
sage: O.basis()
(1, x*y, x^2*y^2, x^3*y^3, x^4*y^4)

We try an example, in which the defining polynomial is not monic and is not integral:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(x^2*y^5 - 1/x); L
Function field in y defined by x^2*y^5 - 1/x
sage: O = L.equation_order()
sage: O.basis()
(1, x^3*y, x^6*y^2, x^9*y^3, x^12*y^4)
free_module(base=None, basis=None, map=True)#

Return a vector space and isomorphisms from the field to and from the vector space.

This function allows us to identify the elements of this field with elements of a vector space over the base field, which is useful for representation and arithmetic with orders, ideals, etc.

INPUT:

  • base – a function field (default: None), the returned vector space is over this subfield \(R\), which defaults to the base field of this function field.

  • basis – a basis for this field over the base.

  • maps – boolean (default True), whether to return \(R\)-linear maps to and from \(V\).

OUTPUT:

  • a vector space over the base function field

  • an isomorphism from the vector space to the field (if requested)

  • an isomorphism from the field to the vector space (if requested)

EXAMPLES:

We define a function field:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x)); L
Function field in y defined by y^5 - 2*x*y + (-x^4 - 1)/x

We get the vector spaces, and maps back and forth:

sage: V, from_V, to_V = L.free_module()
sage: V
Vector space of dimension 5 over Rational function field in x over Rational Field
sage: from_V
Isomorphism:
  From: Vector space of dimension 5 over Rational function field in x over Rational Field
  To:   Function field in y defined by y^5 - 2*x*y + (-x^4 - 1)/x
sage: to_V
Isomorphism:
  From: Function field in y defined by y^5 - 2*x*y + (-x^4 - 1)/x
  To:   Vector space of dimension 5 over Rational function field in x over Rational Field

We convert an element of the vector space back to the function field:

sage: from_V(V.1)
y

We define an interesting element of the function field:

sage: a = 1/L.0; a
(x/(x^4 + 1))*y^4 - 2*x^2/(x^4 + 1)

We convert it to the vector space, and get a vector over the base field:

sage: to_V(a)
(-2*x^2/(x^4 + 1), 0, 0, 0, x/(x^4 + 1))

We convert to and back, and get the same element:

sage: from_V(to_V(a)) == a
True

In the other direction:

sage: v = x*V.0 + (1/x)*V.1
sage: to_V(from_V(v)) == v
True

And we show how it works over an extension of an extension field:

sage: R2.<z> = L[]; M.<z> = L.extension(z^2 -y)
sage: M.free_module()
(Vector space of dimension 2 over Function field in y defined by y^5 - 2*x*y + (-x^4 - 1)/x, Isomorphism:
  From: Vector space of dimension 2 over Function field in y defined by y^5 - 2*x*y + (-x^4 - 1)/x
  To:   Function field in z defined by z^2 - y, Isomorphism:
  From: Function field in z defined by z^2 - y
  To:   Vector space of dimension 2 over Function field in y defined by y^5 - 2*x*y + (-x^4 - 1)/x)

We can also get the vector space of M over K:

sage: M.free_module(K)
(Vector space of dimension 10 over Rational function field in x over Rational Field, Isomorphism:
  From: Vector space of dimension 10 over Rational function field in x over Rational Field
  To:   Function field in z defined by z^2 - y, Isomorphism:
  From: Function field in z defined by z^2 - y
  To:   Vector space of dimension 10 over Rational function field in x over Rational Field)
gen(n=0)#

Return the \(n\)-th generator of the function field. By default, \(n\) is 0; any other value of \(n\) leads to an error. The generator is the class of \(y\), if we view the function field as being presented as \(K[y]/(f(y))\).

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x))
sage: L.gen()
y
sage: L.gen(1)
Traceback (most recent call last):
...
IndexError: there is only one generator
genus()#

Return the genus of the function field.

For now, the genus is computed using Singular.

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^3 - (x^3 + 2*x*y + 1/x))
sage: L.genus()
3
hom(im_gens, base_morphism=None)#

Create a homomorphism from the function field to another function field.

INPUT:

  • im_gens – list of images of the generators of the function field and of successive base rings.

  • base_morphism – homomorphism of the base ring, after the im_gens are used. Thus if im_gens has length 2, then base_morphism should be a morphism from the base ring of the base ring of the function field.

EXAMPLES:

We create a rational function field, and a quadratic extension of it:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^2 - x^3 - 1)

We make the field automorphism that sends y to -y:

sage: f = L.hom(-y); f
Function Field endomorphism of Function field in y defined by y^2 - x^3 - 1
  Defn: y |--> -y

Evaluation works:

sage: f(y*x - 1/x)
-x*y - 1/x

We try to define an invalid morphism:

sage: f = L.hom(y+1)
Traceback (most recent call last):
...
ValueError: invalid morphism

We make a morphism of the base rational function field:

sage: phi = K.hom(x+1); phi
Function Field endomorphism of Rational function field in x over Rational Field
  Defn: x |--> x + 1
sage: phi(x^3 - 3)
x^3 + 3*x^2 + 3*x - 2
sage: (x+1)^3-3
x^3 + 3*x^2 + 3*x - 2

We make a morphism by specifying where the generators and the base generators go:

sage: L.hom([-y, x])
Function Field endomorphism of Function field in y defined by y^2 - x^3 - 1
  Defn: y |--> -y
        x |--> x

You can also specify a morphism on the base:

sage: R1.<q> = K[]
sage: L1.<q> = K.extension(q^2 - (x+1)^3 - 1)
sage: L.hom(q, base_morphism=phi)
Function Field morphism:
  From: Function field in y defined by y^2 - x^3 - 1
  To:   Function field in q defined by q^2 - x^3 - 3*x^2 - 3*x - 2
  Defn: y |--> q
        x |--> x + 1

We make another extension of a rational function field:

sage: K2.<t> = FunctionField(QQ); R2.<w> = K2[]
sage: L2.<w> = K2.extension((4*w)^2 - (t+1)^3 - 1)

We define a morphism, by giving the images of generators:

sage: f = L.hom([4*w, t+1]); f
Function Field morphism:
  From: Function field in y defined by y^2 - x^3 - 1
  To:   Function field in w defined by 16*w^2 - t^3 - 3*t^2 - 3*t - 2
  Defn: y |--> 4*w
        x |--> t + 1

Evaluation works, as expected:

sage: f(y+x)
4*w + t + 1
sage: f(x*y + x/(x^2+1))
(4*t + 4)*w + (t + 1)/(t^2 + 2*t + 2)

We make another extension of a rational function field:

sage: K3.<yy> = FunctionField(QQ); R3.<xx> = K3[]
sage: L3.<xx> = K3.extension(yy^2 - xx^3 - 1)

This is the function field L with the generators exchanged. We define a morphism to L:

sage: g = L3.hom([x,y]); g
Function Field morphism:
  From: Function field in xx defined by -xx^3 + yy^2 - 1
  To:   Function field in y defined by y^2 - x^3 - 1
  Defn: xx |--> x
        yy |--> y
is_separable(base=None)#

Return whether this is a separable extension of base.

INPUT:

  • base – a function field from which this field has been created as an extension or None (default: None); if None, then return whether this is a separable extension over its base field.

EXAMPLES:

sage: K.<x> = FunctionField(GF(2))
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^2 - x)
sage: L.is_separable()
False
sage: R.<z> = L[]
sage: M.<z> = L.extension(z^3 - y)
sage: M.is_separable()
True
sage: M.is_separable(K)
False

sage: K.<x> = FunctionField(GF(5))
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x))
sage: L.is_separable()
True

sage: K.<x> = FunctionField(GF(5))
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^5 - 1)
sage: L.is_separable()
False
maximal_order()#

Return the maximal order of the function field.

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x))
sage: L.maximal_order()
Maximal order of Function field in y defined by y^5 - 2*x*y + (-x^4 - 1)/x
maximal_order_infinite()#

Return the maximal infinite order of the function field.

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x))
sage: L.maximal_order_infinite()
Maximal infinite order of Function field in y defined by y^5 - 2*x*y + (-x^4 - 1)/x

sage: K.<x> = FunctionField(GF(2)); _.<t> = K[]
sage: F.<y> = K.extension(t^3 - x^2*(x^2 + x + 1)^2)
sage: F.maximal_order_infinite()
Maximal infinite order of Function field in y defined by y^3 + x^6 + x^4 + x^2

sage: K.<x> = FunctionField(GF(2)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^2 + Y + x + 1/x)
sage: L.maximal_order_infinite()
Maximal infinite order of Function field in y defined by y^2 + y + (x^2 + 1)/x
monic_integral_model(names=None)#

Return a function field isomorphic to this field but which is an extension of a rational function field with defining polynomial that is monic and integral over the constant base field.

INPUT:

  • names – a string or a tuple of up to two strings (default: None), the name of the generator of the field, and the name of the generator of the underlying rational function field (if a tuple); if not given, then the names are chosen automatically.

OUTPUT:

A triple (F,f,t) where F is a function field, f is an isomorphism from F to this field, and t is the inverse of f.

EXAMPLES:

sage: K.<x> = FunctionField(QQ)
sage: R.<y> = K[]
sage: L.<y> = K.extension(x^2*y^5 - 1/x); L
Function field in y defined by x^2*y^5 - 1/x
sage: A, from_A, to_A = L.monic_integral_model('z')
sage: A
Function field in z defined by z^5 - x^12
sage: from_A
Function Field morphism:
  From: Function field in z defined by z^5 - x^12
  To:   Function field in y defined by x^2*y^5 - 1/x
  Defn: z |--> x^3*y
        x |--> x
sage: to_A
Function Field morphism:
  From: Function field in y defined by x^2*y^5 - 1/x
  To:   Function field in z defined by z^5 - x^12
  Defn: y |--> 1/x^3*z
        x |--> x
sage: to_A(y)
1/x^3*z
sage: from_A(to_A(y))
y
sage: from_A(to_A(1/y))
x^3*y^4
sage: from_A(to_A(1/y)) == 1/y
True

This also works for towers of function fields:

sage: R.<z> = L[]
sage: M.<z> = L.extension(z^2*y - 1/x)
sage: M.monic_integral_model()
(Function field in z_ defined by z_^10 - x^18, Function Field morphism:
  From: Function field in z_ defined by z_^10 - x^18
  To:   Function field in z defined by y*z^2 - 1/x
  Defn: z_ |--> x^2*z
        x |--> x, Function Field morphism:
  From: Function field in z defined by y*z^2 - 1/x
  To:   Function field in z_ defined by z_^10 - x^18
  Defn: z |--> 1/x^2*z_
        y |--> 1/x^15*z_^8
        x |--> x)
ngens()#

Return the number of generators of the function field over its base field. This is by definition 1.

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x))
sage: L.ngens()
1
polynomial()#

Return the univariate polynomial that defines the function field, that is, the polynomial \(f(y)\) so that the function field is of the form \(K[y]/(f(y))\).

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x))
sage: L.polynomial()
y^5 - 2*x*y + (-x^4 - 1)/x
polynomial_ring()#

Return the polynomial ring used to represent elements of the function field. If we view the function field as being presented as \(K[y]/(f(y))\), then this function returns the ring \(K[y]\).

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x))
sage: L.polynomial_ring()
Univariate Polynomial Ring in y over Rational function field in x over Rational Field
primitive_element()#

Return a primitive element over the underlying rational function field.

If this is a finite extension of a rational function field \(K(x)\) with \(K\) perfect, then this is a simple extension of \(K(x)\), i.e., there is a primitive element \(y\) which generates this field over \(K(x)\). This method returns such an element \(y\).

EXAMPLES:

sage: K.<x> = FunctionField(QQ)
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^2-x)
sage: R.<z> = L[]
sage: M.<z> = L.extension(z^2-y)
sage: R.<z> = L[]
sage: N.<u> = L.extension(z^2-x-1)
sage: N.primitive_element()
u + y
sage: M.primitive_element()
z
sage: L.primitive_element()
y

This also works for inseparable extensions:

sage: K.<x> = FunctionField(GF(2))
sage: R.<Y> = K[]
sage: L.<y> = K.extension(Y^2-x)
sage: R.<Z> = L[]
sage: M.<z> = L.extension(Z^2-y)
sage: M.primitive_element()
z
random_element(*args, **kwds)#

Create a random element of the function field. Parameters are passed onto the random_element method of the base_field.

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: L.<y> = K.extension(y^2 - (x^2 + x))
sage: L.random_element() # random
((x^2 - x + 2/3)/(x^2 + 1/3*x - 1))*y^2 + ((-1/4*x^2 + 1/2*x - 1)/(-5/2*x + 2/3))*y
+ (-1/2*x^2 - 4)/(-12*x^2 + 1/2*x - 1/95)
separable_model(names=None)#

Return a function field isomorphic to this field which is a separable extension of a rational function field.

INPUT:

  • names – a tuple of two strings or None (default: None); the second entry will be used as the variable name of the rational function field, the first entry will be used as the variable name of its separable extension. If None, then the variable names will be chosen automatically.

OUTPUT:

A triple (F,f,t) where F is a function field, f is an isomorphism from F to this function field, and t is the inverse of f.

ALGORITHM:

Suppose that the constant base field is perfect. If this is a monic integral inseparable extension of a rational function field, then the defining polynomial is separable if we swap the variables (Proposition 4.8 in Chapter VIII of [Lan2002].) The algorithm reduces to this case with monic_integral_model().

EXAMPLES:

sage: K.<x> = FunctionField(GF(2))
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^2 - x^3)
sage: L.separable_model(('t','w'))
(Function field in t defined by t^3 + w^2,
 Function Field morphism:
   From: Function field in t defined by t^3 + w^2
   To:   Function field in y defined by y^2 + x^3
   Defn: t |--> x
         w |--> y,
 Function Field morphism:
   From: Function field in y defined by y^2 + x^3
   To:   Function field in t defined by t^3 + w^2
   Defn: y |--> w
         x |--> t)

This also works for non-integral polynomials:

sage: K.<x> = FunctionField(GF(2))
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^2/x - x^2)
sage: L.separable_model()
(Function field in y_ defined by y_^3 + x_^2,
 Function Field morphism:
   From: Function field in y_ defined by y_^3 + x_^2
   To:   Function field in y defined by 1/x*y^2 + x^2
   Defn: y_ |--> x
         x_ |--> y,
 Function Field morphism:
   From: Function field in y defined by 1/x*y^2 + x^2
   To:   Function field in y_ defined by y_^3 + x_^2
   Defn: y |--> x_
         x |--> y_)

If the base field is not perfect this is only implemented in trivial cases:

sage: k.<t> = FunctionField(GF(2))
sage: k.is_perfect()
False
sage: K.<x> = FunctionField(k)
sage: R.<y> = K[]
sage: L.<y> = K.extension(y^3 - t)
sage: L.separable_model()
(Function field in y defined by y^3 + t,
 Function Field endomorphism of Function field in y defined by y^3 + t
   Defn: y |--> y
         x |--> x,
 Function Field endomorphism of Function field in y defined by y^3 + t
   Defn: y |--> y
         x |--> x)

Some other cases for which a separable model could be constructed are not supported yet:

sage: R.<y> = K[]
sage: L.<y> = K.extension(y^2 - t)
sage: L.separable_model()
Traceback (most recent call last):
...
NotImplementedError: constructing a separable model is only implemented for function fields over a perfect constant base field
simple_model(name=None)#

Return a function field isomorphic to this field which is a simple extension of a rational function field.

INPUT:

  • name – a string (default: None), the name of generator of the simple extension. If None, then the name of the generator will be the same as the name of the generator of this function field.

OUTPUT:

A triple (F,f,t) where F is a field isomorphic to this field, f is an isomorphism from F to this function field and t is the inverse of f.

EXAMPLES:

A tower of four function fields:

sage: K.<x> = FunctionField(QQ); R.<z> = K[]
sage: L.<z> = K.extension(z^2-x); R.<u> = L[]
sage: M.<u> = L.extension(u^2-z); R.<v> = M[]
sage: N.<v> = M.extension(v^2-u)

The fields N and M as simple extensions of K:

sage: N.simple_model()
(Function field in v defined by v^8 - x,
 Function Field morphism:
  From: Function field in v defined by v^8 - x
  To:   Function field in v defined by v^2 - u
  Defn: v |--> v,
 Function Field morphism:
  From: Function field in v defined by v^2 - u
  To:   Function field in v defined by v^8 - x
  Defn: v |--> v
        u |--> v^2
        z |--> v^4
        x |--> x)
sage: M.simple_model()
(Function field in u defined by u^4 - x,
 Function Field morphism:
  From: Function field in u defined by u^4 - x
  To:   Function field in u defined by u^2 - z
  Defn: u |--> u,
 Function Field morphism:
  From: Function field in u defined by u^2 - z
  To:   Function field in u defined by u^4 - x
  Defn: u |--> u
        z |--> u^2
        x |--> x)

An optional parameter name can be used to set the name of the generator of the simple extension:

sage: M.simple_model(name='t')
(Function field in t defined by t^4 - x, Function Field morphism:
  From: Function field in t defined by t^4 - x
  To:   Function field in u defined by u^2 - z
  Defn: t |--> u, Function Field morphism:
  From: Function field in u defined by u^2 - z
  To:   Function field in t defined by t^4 - x
  Defn: u |--> t
        z |--> t^2
        x |--> x)

An example with higher degrees:

sage: K.<x> = FunctionField(GF(3)); R.<y> = K[]
sage: L.<y> = K.extension(y^5-x); R.<z> = L[]
sage: M.<z> = L.extension(z^3-x)
sage: M.simple_model()
(Function field in z defined by z^15 + x*z^12 + x^2*z^9 + 2*x^3*z^6 + 2*x^4*z^3 + 2*x^5 + 2*x^3,
 Function Field morphism:
   From: Function field in z defined by z^15 + x*z^12 + x^2*z^9 + 2*x^3*z^6 + 2*x^4*z^3 + 2*x^5 + 2*x^3
   To:   Function field in z defined by z^3 + 2*x
   Defn: z |--> z + y,
 Function Field morphism:
   From: Function field in z defined by z^3 + 2*x
   To:   Function field in z defined by z^15 + x*z^12 + x^2*z^9 + 2*x^3*z^6 + 2*x^4*z^3 + 2*x^5 + 2*x^3
   Defn: z |--> 2/x*z^6 + 2*z^3 + z + 2*x
         y |--> 1/x*z^6 + z^3 + x
         x |--> x)

This also works for inseparable extensions:

sage: K.<x> = FunctionField(GF(2)); R.<y> = K[]
sage: L.<y> = K.extension(y^2-x); R.<z> = L[]
sage: M.<z> = L.extension(z^2-y)
sage: M.simple_model()
(Function field in z defined by z^4 + x, Function Field morphism:
   From: Function field in z defined by z^4 + x
   To:   Function field in z defined by z^2 + y
   Defn: z |--> z, Function Field morphism:
   From: Function field in z defined by z^2 + y
   To:   Function field in z defined by z^4 + x
   Defn: z |--> z
         y |--> z^2
         x |--> x)
class sage.rings.function_field.function_field.FunctionField_simple(polynomial, names, category=None)#

Bases: FunctionField_polymod

Function fields defined by irreducible and separable polynomials over rational function fields.

constant_field()#

Return the algebraic closure of the base constant field in the function field.

EXAMPLES:

sage: K.<x> = FunctionField(GF(3)); _.<y> = K[]
sage: L.<y> = K.extension(y^5 - (x^3 + 2*x*y + 1/x))
sage: L.constant_field()
Finite Field of size 3
exact_constant_field(name='t')#

Return the exact constant field and its embedding into the function field.

INPUT:

  • name – name (default: \(t\)) of the generator of the exact constant field

EXAMPLES:

sage: K.<x> = FunctionField(GF(3)); _.<Y> = K[]
sage: f = Y^2 - x*Y + x^2 + 1 # irreducible but not absolutely irreducible
sage: L.<y> = K.extension(f)
sage: L.genus()
0
sage: L.exact_constant_field()
(Finite Field in t of size 3^2, Ring morphism:
   From: Finite Field in t of size 3^2
   To:   Function field in y defined by y^2 + 2*x*y + x^2 + 1
   Defn: t |--> y + x)
sage: (y+x).divisor()
0
genus()#

Return the genus of the function field.

EXAMPLES:

sage: F.<a> = GF(16)
sage: K.<x> = FunctionField(F); K
Rational function field in x over Finite Field in a of size 2^4
sage: R.<t> = PolynomialRing(K)
sage: L.<y> = K.extension(t^4+t-x^5)
sage: L.genus()
6

The genus is computed by the Hurwitz genus formula.

places_above(p)#

Return places lying above p.

INPUT:

  • p – place of the base rational function field.

EXAMPLES:

sage: K.<x> = FunctionField(GF(2)); _.<Y> = K[]
sage: F.<y> = K.extension(Y^3 - x^2*(x^2 + x + 1)^2)
sage: all(q.place_below() == p for p in K.places() for q in F.places_above(p))
True

sage: K.<x> = FunctionField(QQ); _.<Y> = K[]
sage: F.<y> = K.extension(Y^3 - x^2*(x^2 + x + 1)^2)
sage: O = K.maximal_order()
sage: pls = [O.ideal(x-c).place() for c in [-2, -1, 0, 1, 2]]
sage: all(q.place_below() == p for p in pls for q in F.places_above(p))
True

sage: K.<x> = FunctionField(QQbar); _.<Y> = K[]
sage: F.<y> = K.extension(Y^3 - x^2*(x^2 + x + 1)^2)
sage: O = K.maximal_order()
sage: pls = [O.ideal(x-QQbar(sqrt(c))).place() for c in [-2, -1, 0, 1, 2]]
sage: all(q.place_below() == p         # long time (4s)
....:     for p in pls for q in F.places_above(p))
True
residue_field(place, name=None)#

Return the residue field associated with the place along with the maps from and to the residue field.

INPUT:

  • place – place of the function field

  • name – string; name of the generator of the residue field

The domain of the map to the residue field is the discrete valuation ring associated with the place.

The discrete valuation ring is defined as the ring of all elements of the function field with nonnegative valuation at the place. The maximal ideal is the set of elements of positive valuation. The residue field is then the quotient of the discrete valuation ring by its maximal ideal.

If an element not in the valuation ring is applied to the map, an exception TypeError is raised.

EXAMPLES:

sage: K.<x> = FunctionField(GF(2)); _.<Y> = K[]
sage: L.<y> = K.extension(Y^2 + Y + x + 1/x)
sage: p = L.places_finite()[0]
sage: R, fr_R, to_R = L.residue_field(p)
sage: R
Finite Field of size 2
sage: f = 1 + y
sage: f.valuation(p)
-1
sage: to_R(f)
Traceback (most recent call last):
...
TypeError: ...
sage: (1+1/f).valuation(p)
0
sage: to_R(1 + 1/f)
1
sage: [fr_R(e) for e in R]
[0, 1]
class sage.rings.function_field.function_field.RationalFunctionField(constant_field, names, category=None)#

Bases: FunctionField

Rational function field in one variable, over an arbitrary base field.

INPUT:

  • constant_field – arbitrary field

  • names – string or tuple of length 1

EXAMPLES:

sage: K.<t> = FunctionField(GF(3)); K
Rational function field in t over Finite Field of size 3
sage: K.gen()
t
sage: 1/t + t^3 + 5
(t^4 + 2*t + 1)/t

sage: K.<t> = FunctionField(QQ); K
Rational function field in t over Rational Field
sage: K.gen()
t
sage: 1/t + t^3 + 5
(t^4 + 5*t + 1)/t

There are various ways to get at the underlying fields and rings associated to a rational function field:

sage: K.<t> = FunctionField(GF(7))
sage: K.base_field()
Rational function field in t over Finite Field of size 7
sage: K.field()
Fraction Field of Univariate Polynomial Ring in t over Finite Field of size 7
sage: K.constant_field()
Finite Field of size 7
sage: K.maximal_order()
Maximal order of Rational function field in t over Finite Field of size 7

sage: K.<t> = FunctionField(QQ)
sage: K.base_field()
Rational function field in t over Rational Field
sage: K.field()
Fraction Field of Univariate Polynomial Ring in t over Rational Field
sage: K.constant_field()
Rational Field
sage: K.maximal_order()
Maximal order of Rational function field in t over Rational Field

We define a morphism:

sage: K.<t> = FunctionField(QQ)
sage: L = FunctionField(QQ, 'tbar') # give variable name as second input
sage: K.hom(L.gen())
Function Field morphism:
  From: Rational function field in t over Rational Field
  To:   Rational function field in tbar over Rational Field
  Defn: t |--> tbar

Here are some calculations over a number field:

sage: R.<x> = FunctionField(QQ)
sage: L.<y> = R[]
sage: F.<y> = R.extension(y^2 - (x^2+1))
sage: (y/x).divisor()
- Place (x, y - 1)
 - Place (x, y + 1)
 + Place (x^2 + 1, y)

sage: A.<z> = QQ[]
sage: NF.<i> = NumberField(z^2+1)
sage: R.<x> = FunctionField(NF)
sage: L.<y> = R[]
sage: F.<y> = R.extension(y^2 - (x^2+1))

sage: (x/y*x.differential()).divisor()
-2*Place (1/x, 1/x*y - 1)
 - 2*Place (1/x, 1/x*y + 1)
 + Place (x, y - 1)
 + Place (x, y + 1)

sage: (x/y).divisor()
- Place (x - i, y)
 + Place (x, y - 1)
 + Place (x, y + 1)
 - Place (x + i, y)
Element#

alias of FunctionFieldElement_rational

base_field()#

Return the base field of the rational function field, which is just the function field itself.

EXAMPLES:

sage: K.<t> = FunctionField(GF(7))
sage: K.base_field()
Rational function field in t over Finite Field of size 7
change_variable_name(name)#

Return a field isomorphic to this field with variable name.

INPUT:

  • name – a string or a tuple consisting of a single string, the name of the new variable

OUTPUT:

A triple F,f,t where F is a rational function field, f is an isomorphism from F to this field, and t is the inverse of f.

EXAMPLES:

sage: K.<x> = FunctionField(QQ)
sage: L,f,t = K.change_variable_name('y')
sage: L,f,t
(Rational function field in y over Rational Field,
 Function Field morphism:
  From: Rational function field in y over Rational Field
  To:   Rational function field in x over Rational Field
  Defn: y |--> x,
 Function Field morphism:
  From: Rational function field in x over Rational Field
  To:   Rational function field in y over Rational Field
  Defn: x |--> y)
sage: L.change_variable_name('x')[0] is K
True
constant_base_field()#

Return the field of which the rational function field is a transcendental extension.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.constant_base_field()
Rational Field
constant_field()#

Return the field of which the rational function field is a transcendental extension.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.constant_base_field()
Rational Field
degree(base=None)#

Return the degree over the base field of the rational function field. Since the base field is the rational function field itself, the degree is 1.

INPUT:

  • base – the base field of the vector space; must be the function field itself (the default)

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.degree()
1
different()#

Return the different of the rational function field.

For a rational function field, the different is simply the zero divisor.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.different()
0
equation_order()#

Return the maximal order of the function field.

Since this is a rational function field it is of the form \(K(t)\), and the maximal order is by definition \(K[t]\), where \(K\) is the constant field.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.maximal_order()
Maximal order of Rational function field in t over Rational Field
sage: K.equation_order()
Maximal order of Rational function field in t over Rational Field
equation_order_infinite()#

Return the maximal infinite order of the function field.

By definition, this is the valuation ring of the degree valuation of the rational function field.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.maximal_order_infinite()
Maximal infinite order of Rational function field in t over Rational Field
sage: K.equation_order_infinite()
Maximal infinite order of Rational function field in t over Rational Field
extension(f, names=None)#

Create an extension \(L = K[y]/(f(y))\) of the rational function field.

INPUT:

  • f – univariate polynomial over self

  • names – string or length-1 tuple

OUTPUT:

  • a function field

EXAMPLES:

sage: K.<x> = FunctionField(QQ); R.<y> = K[]
sage: K.extension(y^5 - x^3 - 3*x + x*y)
Function field in y defined by y^5 + x*y - x^3 - 3*x

A nonintegral defining polynomial:

sage: K.<t> = FunctionField(QQ); R.<y> = K[]
sage: K.extension(y^3 + (1/t)*y + t^3/(t+1))
Function field in y defined by y^3 + 1/t*y + t^3/(t + 1)

The defining polynomial need not be monic or integral:

sage: K.extension(t*y^3 + (1/t)*y + t^3/(t+1))
Function field in y defined by t*y^3 + 1/t*y + t^3/(t + 1)
field()#

Return the underlying field, forgetting the function field structure.

EXAMPLES:

sage: K.<t> = FunctionField(GF(7))
sage: K.field()
Fraction Field of Univariate Polynomial Ring in t over Finite Field of size 7
free_module(base=None, basis=None, map=True)#

Return a vector space \(V\) and isomorphisms from the field to \(V\) and from \(V\) to the field.

This function allows us to identify the elements of this field with elements of a one-dimensional vector space over the field itself. This method exists so that all function fields (rational or not) have the same interface.

INPUT:

  • base – the base field of the vector space; must be the function field itself (the default)

  • basis – (ignored) a basis for the vector space

  • map – (default True), whether to return maps to and from the vector space

OUTPUT:

  • a vector space \(V\) over base field

  • an isomorphism from \(V\) to the field

  • the inverse isomorphism from the field to \(V\)

EXAMPLES:

sage: K.<x> = FunctionField(QQ)
sage: K.free_module()
(Vector space of dimension 1 over Rational function field in x over Rational Field, Isomorphism:
  From: Vector space of dimension 1 over Rational function field in x over Rational Field
  To:   Rational function field in x over Rational Field, Isomorphism:
  From: Rational function field in x over Rational Field
  To:   Vector space of dimension 1 over Rational function field in x over Rational Field)
gen(n=0)#

Return the n-th generator of the function field. If n is not 0, then an IndexError is raised.

EXAMPLES:

sage: K.<t> = FunctionField(QQ); K.gen()
t
sage: K.gen().parent()
Rational function field in t over Rational Field
sage: K.gen(1)
Traceback (most recent call last):
...
IndexError: Only one generator.
genus()#

Return the genus of the function field, namely 0.

EXAMPLES:

sage: K.<x> = FunctionField(QQ)
sage: K.genus()
0
hom(im_gens, base_morphism=None)#

Create a homomorphism from self to another ring.

INPUT:

  • im_gens – exactly one element of some ring. It must be invertible and transcendental over the image of base_morphism; this is not checked.

  • base_morphism – a homomorphism from the base field into the other ring. If None, try to use a coercion map.

OUTPUT:

  • a map between function fields

EXAMPLES:

We make a map from a rational function field to itself:

sage: K.<x> = FunctionField(GF(7))
sage: K.hom( (x^4 + 2)/x)
Function Field endomorphism of Rational function field in x over Finite Field of size 7
  Defn: x |--> (x^4 + 2)/x

We construct a map from a rational function field into a non-rational extension field:

sage: K.<x> = FunctionField(GF(7)); R.<y> = K[]
sage: L.<y> = K.extension(y^3 + 6*x^3 + x)
sage: f = K.hom(y^2 + y  + 2); f
Function Field morphism:
  From: Rational function field in x over Finite Field of size 7
  To:   Function field in y defined by y^3 + 6*x^3 + x
  Defn: x |--> y^2 + y + 2
sage: f(x)
y^2 + y + 2
sage: f(x^2)
5*y^2 + (x^3 + 6*x + 4)*y + 2*x^3 + 5*x + 4
maximal_order()#

Return the maximal order of the function field.

Since this is a rational function field it is of the form \(K(t)\), and the maximal order is by definition \(K[t]\), where \(K\) is the constant field.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.maximal_order()
Maximal order of Rational function field in t over Rational Field
sage: K.equation_order()
Maximal order of Rational function field in t over Rational Field
maximal_order_infinite()#

Return the maximal infinite order of the function field.

By definition, this is the valuation ring of the degree valuation of the rational function field.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.maximal_order_infinite()
Maximal infinite order of Rational function field in t over Rational Field
sage: K.equation_order_infinite()
Maximal infinite order of Rational function field in t over Rational Field
ngens()#

Return the number of generators, which is 1.

EXAMPLES:

sage: K.<t> = FunctionField(QQ)
sage: K.ngens()
1
polynomial_ring(var='x')#

Return a polynomial ring in one variable over the rational function field.

INPUT:

  • var – string; name of the variable

EXAMPLES:

sage: K.<x> = FunctionField(QQ)
sage: K.polynomial_ring()
Univariate Polynomial Ring in x over Rational function field in x over Rational Field
sage: K.polynomial_ring('T')
Univariate Polynomial Ring in T over Rational function field in x over Rational Field
random_element(*args, **kwds)#

Create a random element of the rational function field.

Parameters are passed to the random_element method of the underlying fraction field.

EXAMPLES:

sage: FunctionField(QQ,'alpha').random_element()   # random
(-1/2*alpha^2 - 4)/(-12*alpha^2 + 1/2*alpha - 1/95)
residue_field(place, name=None)#

Return the residue field of the place along with the maps from and to it.

INPUT:

  • place – place of the function field

  • name – string; name of the generator of the residue field

EXAMPLES:

sage: F.<x> = FunctionField(GF(5))
sage: p = F.places_finite(2)[0]
sage: R, fr_R, to_R = F.residue_field(p)
sage: R
Finite Field in z2 of size 5^2
sage: to_R(x) in R
True
class sage.rings.function_field.function_field.RationalFunctionField_char_zero(constant_field, names, category=None)#

Bases: RationalFunctionField

Rational function fields of characteristic zero.

higher_derivation()#

Return the higher derivation for the function field.

This is also called the Hasse-Schmidt derivation.

EXAMPLES:

sage: F.<x> = FunctionField(QQ)
sage: d = F.higher_derivation()
sage: [d(x^5,i) for i in range(10)]
[x^5, 5*x^4, 10*x^3, 10*x^2, 5*x, 1, 0, 0, 0, 0]
sage: [d(x^9,i) for i in range(10)]
[x^9, 9*x^8, 36*x^7, 84*x^6, 126*x^5, 126*x^4, 84*x^3, 36*x^2, 9*x, 1]
class sage.rings.function_field.function_field.RationalFunctionField_global(constant_field, names, category=None)#

Bases: RationalFunctionField

Rational function field over finite fields.

get_place(degree)#

Return a place of degree.

INPUT:

  • degree – a positive integer

EXAMPLES:

sage: F.<a> = GF(2)
sage: K.<x> = FunctionField(F)
sage: K.get_place(1)
Place (x)
sage: K.get_place(2)
Place (x^2 + x + 1)
sage: K.get_place(3)
Place (x^3 + x + 1)
sage: K.get_place(4)
Place (x^4 + x + 1)
sage: K.get_place(5)
Place (x^5 + x^2 + 1)
higher_derivation()#

Return the higher derivation for the function field.

This is also called the Hasse-Schmidt derivation.

EXAMPLES:

sage: F.<x> = FunctionField(GF(5))
sage: d = F.higher_derivation()
sage: [d(x^5,i) for i in range(10)]
[x^5, 0, 0, 0, 0, 1, 0, 0, 0, 0]
sage: [d(x^7,i) for i in range(10)]
[x^7, 2*x^6, x^5, 0, 0, x^2, 2*x, 1, 0, 0]
place_infinite()#

Return the unique place at infinity.

EXAMPLES:

sage: F.<x> = FunctionField(GF(5))
sage: F.place_infinite()
Place (1/x)
places(degree=1)#

Return all places of the degree.

INPUT:

  • degree – (default: 1) a positive integer

EXAMPLES:

sage: F.<x> = FunctionField(GF(5))
sage: F.places()
[Place (1/x),
 Place (x),
 Place (x + 1),
 Place (x + 2),
 Place (x + 3),
 Place (x + 4)]
places_finite(degree=1)#

Return the finite places of the degree.

INPUT:

  • degree – (default: 1) a positive integer

EXAMPLES:

sage: F.<x> = FunctionField(GF(5))
sage: F.places_finite()
[Place (x), Place (x + 1), Place (x + 2), Place (x + 3), Place (x + 4)]
sage.rings.function_field.function_field.is_FunctionField(x)#

Return True if x is a function field.

EXAMPLES:

sage: from sage.rings.function_field.function_field import is_FunctionField
sage: is_FunctionField(QQ)
False
sage: is_FunctionField(FunctionField(QQ, 't'))
True