Universal cyclotomic field¶
The universal cyclotomic field is the smallest subfield of the complex field containing all roots of unity. It is also the maximal Galois Abelian extension of the rational numbers.
The implementation simply wraps GAP Cyclotomic. As mentioned in their documentation: arithmetical operations are quite expensive, so the use of internally represented cyclotomics is not recommended for doing arithmetic over number fields, such as calculations with matrices of cyclotomics.
Note
There used to be a native Sage version of the universal cyclotomic field written by Christian Stump (see trac ticket #8327). It was slower on most operations and it was decided to use a version based on libGAP instead (see trac ticket #18152). One main difference in the design choices is that GAP stores dense vectors whereas the native ones used Python dictionaries (storing only nonzero coefficients). Most operations are faster with libGAP except some operation on very sparse elements. All details can be found in trac ticket #18152.
REFERENCES:
[Bre97] | T. Breuer “Integral bases for subfields of cyclotomic fields” AAECC 8, 279–289 (1997). |
EXAMPLES:
sage: UCF = UniversalCyclotomicField(); UCF
Universal Cyclotomic Field
To generate cyclotomic elements:
sage: UCF.gen(5)
E(5)
sage: UCF.gen(5,2)
E(5)^2
sage: E = UCF.gen
Equality and inequality checks:
sage: E(6,2) == E(6)^2 == E(3)
True
sage: E(6)^2 != E(3)
False
Addition and multiplication:
sage: E(2) * E(3)
-E(3)
sage: f = E(2) + E(3); f
2*E(3) + E(3)^2
Inverses:
sage: f^-1
1/3*E(3) + 2/3*E(3)^2
sage: f.inverse()
1/3*E(3) + 2/3*E(3)^2
sage: f * f.inverse()
1
Conjugation and Galois conjugates:
sage: f.conjugate()
E(3) + 2*E(3)^2
sage: f.galois_conjugates()
[2*E(3) + E(3)^2, E(3) + 2*E(3)^2]
sage: f.norm_of_galois_extension()
3
One can create matrices and polynomials:
sage: m = matrix(2,[E(3),1,1,E(4)]); m
[E(3) 1]
[ 1 E(4)]
sage: m.parent()
Full MatrixSpace of 2 by 2 dense matrices over Universal Cyclotomic Field
sage: m**2
[ -E(3) E(12)^4 - E(12)^7 - E(12)^11]
[E(12)^4 - E(12)^7 - E(12)^11 0]
sage: m.charpoly()
x^2 + (-E(12)^4 + E(12)^7 + E(12)^11)*x + E(12)^4 + E(12)^7 + E(12)^8
sage: m.echelon_form()
[1 0]
[0 1]
sage: m.pivots()
(0, 1)
sage: m.rank()
2
sage: R.<x> = PolynomialRing(UniversalCyclotomicField(), 'x')
sage: E(3) * x - 1
E(3)*x - 1
AUTHORS:
- Christian Stump (2013): initial Sage version (see trac ticket #8327)
- Vincent Delecroix (2015): complete rewriting using libgap (see trac ticket #18152)
-
sage.rings.universal_cyclotomic_field.
E
(n, k=1)¶ Return the
n
-th root of unity as an element of the universal cyclotomic field.EXAMPLES:
sage: E(3) E(3) sage: E(3) + E(5) -E(15)^2 - 2*E(15)^8 - E(15)^11 - E(15)^13 - E(15)^14
-
class
sage.rings.universal_cyclotomic_field.
UCFtoQQbar
(UCF)¶ Bases:
sage.categories.morphism.Morphism
Conversion to
QQbar
.EXAMPLES:
sage: UCF = UniversalCyclotomicField() sage: QQbar(UCF.gen(3)) -0.500000000000000? + 0.866025403784439?*I sage: CC(UCF.gen(7,2) + UCF.gen(7,6)) 0.400968867902419 + 0.193096429713794*I sage: complex(E(7)+E(7,2)) (0.40096886790241915+1.7567593946498534j) sage: complex(UCF.one()/2) (0.5+0j)
-
class
sage.rings.universal_cyclotomic_field.
UniversalCyclotomicField
(names=None)¶ Bases:
sage.structure.unique_representation.UniqueRepresentation
,sage.rings.ring.Field
The universal cyclotomic field.
The universal cyclotomic field is the infinite algebraic extension of \(\QQ\) generated by the roots of unity. It is also the maximal Abelian extension of \(\QQ\) in the sense that any Abelian Galois extension of \(\QQ\) is also a subfield of the universal cyclotomic field.
-
Element
¶ alias of
UniversalCyclotomicFieldElement
-
algebraic_closure
()¶ The algebraic closure.
EXAMPLES:
sage: UniversalCyclotomicField().algebraic_closure() Algebraic Field
-
an_element
()¶ Return an element.
EXAMPLES:
sage: UniversalCyclotomicField().an_element() E(5) - 3*E(5)^2
-
characteristic
()¶ Return the characteristic.
EXAMPLES:
sage: UniversalCyclotomicField().characteristic() 0 sage: parent(_) Integer Ring
-
degree
()¶ Return the degree of
self
as a field extension over the Rationals.EXAMPLES:
sage: UCF = UniversalCyclotomicField() sage: UCF.degree() +Infinity
-
gen
(n, k=1)¶ Return the standard primitive
n
-th root of unity.If
k
is notNone
, return thek
-th power of it.EXAMPLES:
sage: UCF = UniversalCyclotomicField() sage: UCF.gen(15) E(15) sage: UCF.gen(7,3) E(7)^3 sage: UCF.gen(4,2) -1
There is an alias
zeta
also available:sage: UCF.zeta(6) -E(3)^2
-
is_exact
()¶ Return
True
as this is an exact ring (i.e. not numerical).EXAMPLES:
sage: UniversalCyclotomicField().is_exact() True
-
is_finite
()¶ Return
True
.EXAMPLES:
sage: UniversalCyclotomicField().is_finite() True
Todo
this method should be provided by the category.
-
one
()¶ Return one.
EXAMPLES:
sage: UCF = UniversalCyclotomicField() sage: UCF.one() 1 sage: parent(_) Universal Cyclotomic Field
-
some_elements
()¶ Return a tuple of some elements in the universal cyclotomic field.
EXAMPLES:
sage: UniversalCyclotomicField().some_elements() (0, 1, -1, E(3), E(7) - 2/3*E(7)^2) sage: all(parent(x) is UniversalCyclotomicField() for x in _) True
-
zero
()¶ Return zero.
EXAMPLES:
sage: UCF = UniversalCyclotomicField() sage: UCF.zero() 0 sage: parent(_) Universal Cyclotomic Field
-
zeta
(n, k=1)¶ Return the standard primitive
n
-th root of unity.If
k
is notNone
, return thek
-th power of it.EXAMPLES:
sage: UCF = UniversalCyclotomicField() sage: UCF.gen(15) E(15) sage: UCF.gen(7,3) E(7)^3 sage: UCF.gen(4,2) -1
There is an alias
zeta
also available:sage: UCF.zeta(6) -E(3)^2
-
-
class
sage.rings.universal_cyclotomic_field.
UniversalCyclotomicFieldElement
(parent, obj)¶ Bases:
sage.structure.element.FieldElement
INPUT:
parent
- a universal cyclotomic fieldobj
- a libgap element (either an integer, a rational or a cyclotomic)
-
additive_order
()¶ The additive order.
EXAMPLES:
sage: UCF = UniversalCyclotomicField() sage: UCF.zero().additive_order() 0 sage: UCF.one().additive_order() +Infinity sage: UCF.gen(3).additive_order() +Infinity
-
conductor
()¶ Return the conductor of
self
.EXAMPLES:
sage: E(3).conductor() 3 sage: (E(5) + E(3)).conductor() 15
-
conjugate
()¶ Return the complex conjugate.
EXAMPLES:
sage: (E(7) + 3*E(7,2) - 5 * E(7,3)).conjugate() -5*E(7)^4 + 3*E(7)^5 + E(7)^6
-
denominator
()¶ Return the denominator of this element.
EXAMPLES:
sage: a = E(5) + 1/2*E(5,2) + 1/3*E(5,3) sage: a E(5) + 1/2*E(5)^2 + 1/3*E(5)^3 sage: a.denominator() 6 sage: parent(_) Integer Ring
-
field_order
(*args, **kwds)¶ Deprecated: Use
conductor()
instead. See trac ticket #18152 for details.
-
galois_conjugates
(n=None)¶ Return the Galois conjugates of
self
.INPUT:
n
– an optional integer. If provided, return the orbit of the Galois group of then
-th cyclotomic field over \(\QQ\). Note thatn
must be such that this element belongs to then
-th cyclotomic field (in other words, it must be a multiple of the conductor).
EXAMPLES:
sage: E(6).galois_conjugates() [-E(3)^2, -E(3)] sage: E(6).galois_conjugates() [-E(3)^2, -E(3)] sage: (E(9,2) - E(9,4)).galois_conjugates() [E(9)^2 - E(9)^4, E(9)^2 + E(9)^4 + E(9)^5, -E(9)^2 - E(9)^5 - E(9)^7, -E(9)^2 - E(9)^4 - E(9)^7, E(9)^4 + E(9)^5 + E(9)^7, -E(9)^5 + E(9)^7] sage: zeta = E(5) sage: zeta.galois_conjugates(5) [E(5), E(5)^2, E(5)^3, E(5)^4] sage: zeta.galois_conjugates(10) [E(5), E(5)^3, E(5)^2, E(5)^4] sage: zeta.galois_conjugates(15) [E(5), E(5)^2, E(5)^4, E(5)^2, E(5)^3, E(5), E(5)^3, E(5)^4] sage: zeta.galois_conjugates(17) Traceback (most recent call last): ... ValueError: n = 17 must be a multiple of the conductor (5)
-
imag
()¶ Return the imaginary part of this element.
EXAMPLES:
sage: E(3).imag() -1/2*E(12)^7 + 1/2*E(12)^11 sage: E(5).imag() 1/2*E(20) - 1/2*E(20)^9 sage: a = E(5) - 2*E(3) sage: AA(a.imag()) == QQbar(a).imag() True
-
imag_part
()¶ Return the imaginary part of this element.
EXAMPLES:
sage: E(3).imag() -1/2*E(12)^7 + 1/2*E(12)^11 sage: E(5).imag() 1/2*E(20) - 1/2*E(20)^9 sage: a = E(5) - 2*E(3) sage: AA(a.imag()) == QQbar(a).imag() True
-
inverse
()¶
-
is_rational
()¶ Test whether this element is a rational number.
EXAMPLES:
sage: E(3).is_rational() False sage: (E(3) + E(3,2)).is_rational() True
-
is_real
()¶ Test whether this element is real.
EXAMPLES:
sage: E(3).is_real() False sage: (E(3) + E(3,2)).is_real() True sage: a = E(3) - 2*E(7) sage: a.real_part().is_real() True sage: a.imag_part().is_real() True
-
minpoly
(var='x')¶ The minimal polynomial of
self
element over \(\QQ\).INPUT:
var
– (optional, default ‘x’) the name of the variable to use.
EXAMPLES:
sage: UCF.<E> = UniversalCyclotomicField() sage: UCF(4).minpoly() x - 4 sage: UCF(4).minpoly(var='y') y - 4 sage: E(3).minpoly() x^2 + x + 1 sage: E(3).minpoly(var='y') y^2 + y + 1
Todo
Polynomials with libgap currently does not implement a
.sage()
method (see trac ticket #18266). It would be faster/safer to not use string to construct the polynomial.
-
multiplicative_order
()¶ The multiplicative order.
EXAMPLES:
sage: E(5).multiplicative_order() 5 sage: (E(5) + E(12)).multiplicative_order() +Infinity sage: UniversalCyclotomicField().zero().multiplicative_order() Traceback (most recent call last): ... ValueError: libGAP: Error, argument must be nonzero
-
norm_of_galois_extension
()¶ Return the norm as a Galois extension of \(\QQ\), which is given by the product of all galois_conjugates.
EXAMPLES:
sage: E(3).norm_of_galois_extension() 1 sage: E(6).norm_of_galois_extension() 1 sage: (E(2) + E(3)).norm_of_galois_extension() 3 sage: parent(_) Integer Ring
-
real
()¶ Return the real part of this element.
EXAMPLES:
sage: E(3).real() -1/2 sage: E(5).real() 1/2*E(5) + 1/2*E(5)^4 sage: a = E(5) - 2*E(3) sage: AA(a.real()) == QQbar(a).real() True
-
real_part
()¶ Return the real part of this element.
EXAMPLES:
sage: E(3).real() -1/2 sage: E(5).real() 1/2*E(5) + 1/2*E(5)^4 sage: a = E(5) - 2*E(3) sage: AA(a.real()) == QQbar(a).real() True
-
to_cyclotomic_field
(R=None)¶ Return this element as an element of a cyclotomic field.
EXAMPLES:
sage: UCF = UniversalCyclotomicField() sage: UCF.gen(3).to_cyclotomic_field() zeta3 sage: UCF.gen(3,2).to_cyclotomic_field() -zeta3 - 1 sage: CF = CyclotomicField(5) sage: CF(E(5)) # indirect doctest zeta5 sage: CF = CyclotomicField(7) sage: CF(E(5)) # indirect doctest Traceback (most recent call last): ... TypeError: Cannot coerce zeta5 into Cyclotomic Field of order 7 and degree 6 sage: CF = CyclotomicField(10) sage: CF(E(5)) # indirect doctest zeta10^2
Matrices are correctly dealt with:
sage: M = Matrix(UCF,2,[E(3),E(4),E(5),E(6)]); M [ E(3) E(4)] [ E(5) -E(3)^2] sage: Matrix(CyclotomicField(60),M) # indirect doctest [zeta60^10 - 1 zeta60^15] [ zeta60^12 zeta60^10]
Using a non-standard embedding:
sage: CF = CyclotomicField(5,embedding=CC(exp(4*pi*i/5))) sage: x = E(5) sage: CC(x) 0.309016994374947 + 0.951056516295154*I sage: CC(CF(x)) 0.309016994374947 + 0.951056516295154*I
Test that the bug reported in trac ticket #19912 has been fixed:
sage: a = 1+E(4); a 1 + E(4) sage: a.to_cyclotomic_field() zeta4 + 1
-
sage.rings.universal_cyclotomic_field.
late_import
()¶ This function avoids importing libgap on startup. It is called once through the constructor of
UniversalCyclotomicField
.EXAMPLES:
sage: import sage.rings.universal_cyclotomic_field as ucf sage: _ = UniversalCyclotomicField() # indirect doctest sage: ucf.libgap is None # indirect doctest False