Orders in Number Fields#
AUTHORS:
William Stein and Robert Bradshaw (2007-09): initial version
EXAMPLES:
We define an absolute order:
sage: K.<a> = NumberField(x^2 + 1); O = K.order(2*a)
sage: O.basis()
[1, 2*a]
We compute a basis for an order in a relative extension that is generated by 2 elements:
sage: K.<a,b> = NumberField([x^2 + 1, x^2 - 3]); O = K.order([3*a,2*b])
sage: O.basis()
[1, 3*a - 2*b, -6*b*a + 6, 3*a]
We compute a maximal order of a degree 10 field:
sage: K.<a> = NumberField((x+1)^10 + 17)
sage: K.maximal_order()
Maximal Order in Number Field in a with defining polynomial x^10 + 10*x^9 + 45*x^8 + 120*x^7 + 210*x^6 + 252*x^5 + 210*x^4 + 120*x^3 + 45*x^2 + 10*x + 18
We compute a suborder, which has index a power of 17 in the maximal order:
sage: O = K.order(17*a); O
Order in Number Field in a with defining polynomial x^10 + 10*x^9 + 45*x^8 + 120*x^7 + 210*x^6 + 252*x^5 + 210*x^4 + 120*x^3 + 45*x^2 + 10*x + 18
sage: m = O.index_in(K.maximal_order()); m
23453165165327788911665591944416226304630809183732482257
sage: factor(m)
17^45
- class sage.rings.number_field.order.AbsoluteOrderFactory#
Bases:
OrderFactory
An order in an (absolute) number field.
EXAMPLES:
sage: K.<i> = NumberField(x^2 + 1) sage: K.order(i) Order in Number Field in i with defining polynomial x^2 + 1
- create_key_and_extra_args(K, module_rep, is_maximal=None, check=True, is_maximal_at=())#
Return normalized arguments to create an absolute order.
- create_object(version, key, is_maximal=None, is_maximal_at=())#
Create an absolute order.
- reduce_data(order)#
Return the data that can be used to pickle an order created by this factory.
This overrides the default implementation to update the latest knowledge about primes at which the order is maximal.
EXAMPLES:
This also works for relative orders since they are wrapping absolute orders:
sage: L.<a, b> = NumberField([x^2 - 1000003, x^2 - 5*1000099^2]) sage: O = L.maximal_order([5], assume_maximal=None) sage: s = dumps(O) sage: loads(s) is O True sage: N = L.maximal_order([7], assume_maximal=None) sage: dumps(N) == s False sage: loads(dumps(N)) is O True
- sage.rings.number_field.order.EisensteinIntegers(names='omega')#
Return the ring of Eisenstein integers.
This is the ring of all complex numbers of the form \(a + b \omega\) with \(a\) and \(b\) integers and \(omega = (-1 + \sqrt{-3})/2\).
EXAMPLES:
sage: R.<omega> = EisensteinIntegers() sage: R Eisenstein Integers in Number Field in omega with defining polynomial x^2 + x + 1 with omega = -0.50000000000000000? + 0.866025403784439?*I sage: factor(3 + omega) (-1) * (-omega - 3) sage: CC(omega) -0.500000000000000 + 0.866025403784439*I sage: omega.minpoly() x^2 + x + 1 sage: EisensteinIntegers().basis() [1, omega]
- sage.rings.number_field.order.EquationOrder(f, names, **kwds)#
Return the equation order generated by a root of the irreducible polynomial f or list of polynomials \(f\) (to construct a relative equation order).
IMPORTANT: Note that the generators of the returned order need not be roots of \(f\), since the generators of an order are – in Sage – module generators.
EXAMPLES:
sage: O.<a,b> = EquationOrder([x^2+1, x^2+2]) sage: O Relative Order in Number Field in a with defining polynomial x^2 + 1 over its base field sage: O.0 -b*a - 1 sage: O.1 -3*a + 2*b
Of course the input polynomial must be integral:
sage: R = EquationOrder(x^3 + x + 1/3, 'alpha'); R Traceback (most recent call last): ... ValueError: each generator must be integral sage: R = EquationOrder( [x^3 + x + 1, x^2 + 1/2], 'alpha'); R Traceback (most recent call last): ... ValueError: each generator must be integral
- sage.rings.number_field.order.GaussianIntegers(names='I', latex_name='i')#
Return the ring of Gaussian integers.
This is the ring of all complex numbers of the form \(a + b I\) with \(a\) and \(b\) integers and \(I = \sqrt{-1}\).
EXAMPLES:
sage: ZZI.<I> = GaussianIntegers() sage: ZZI Gaussian Integers in Number Field in I with defining polynomial x^2 + 1 with I = 1*I sage: factor(3 + I) (-I) * (I + 1) * (2*I + 1) sage: CC(I) 1.00000000000000*I sage: I.minpoly() x^2 + 1 sage: GaussianIntegers().basis() [1, I]
- class sage.rings.number_field.order.Order(K)#
Bases:
IntegralDomain
,Order
An order in a number field.
An order is a subring of the number field that has \(\ZZ\)-rank equal to the degree of the number field over \(\QQ\).
EXAMPLES:
sage: K.<theta> = NumberField(x^4 + x + 17) sage: K.maximal_order() Maximal Order in Number Field in theta with defining polynomial x^4 + x + 17 sage: R = K.order(17*theta); R Order in Number Field in theta with defining polynomial x^4 + x + 17 sage: R.basis() [1, 17*theta, 289*theta^2, 4913*theta^3] sage: R = K.order(17*theta, 13*theta); R Maximal Order in Number Field in theta with defining polynomial x^4 + x + 17 sage: R.basis() [1, theta, theta^2, theta^3] sage: R = K.order([34*theta, 17*theta + 17]); R Order in Number Field in theta with defining polynomial x^4 + x + 17 sage: K.<b> = NumberField(x^4 + x^2 + 2) sage: (b^2).charpoly().factor() (x^2 + x + 2)^2 sage: K.order(b^2) Traceback (most recent call last): ... ValueError: the rank of the span of gens is wrong
- absolute_degree()#
Return the absolute degree of this order, ie the degree of this order over \(\ZZ\).
EXAMPLES:
sage: K.<a> = NumberField(x^3 + 2) sage: O = K.maximal_order() sage: O.absolute_degree() 3
- ambient()#
Return the ambient number field that contains self.
This is the same as
self.number_field()
andself.fraction_field()
EXAMPLES:
sage: k.<z> = NumberField(x^2 - 389) sage: o = k.order(389*z + 1) sage: o Order in Number Field in z with defining polynomial x^2 - 389 sage: o.basis() [1, 389*z] sage: o.ambient() Number Field in z with defining polynomial x^2 - 389
- basis()#
Return a basis over \(\ZZ\) of this order.
EXAMPLES:
sage: K.<a> = NumberField(x^3 + x^2 - 16*x + 16) sage: O = K.maximal_order(); O Maximal Order in Number Field in a with defining polynomial x^3 + x^2 - 16*x + 16 sage: O.basis() [1, 1/4*a^2 + 1/4*a, a^2]
- class_group(proof=None, names='c')#
Return the class group of this order.
(Currently only implemented for the maximal order.)
EXAMPLES:
sage: k.<a> = NumberField(x^2 + 5077) sage: O = k.maximal_order(); O Maximal Order in Number Field in a with defining polynomial x^2 + 5077 sage: O.class_group() Class group of order 22 with structure C22 of Number Field in a with defining polynomial x^2 + 5077
- class_number(proof=None)#
Return the class number of this order.
EXAMPLES:
sage: ZZ[2^(1/3)].class_number() 1 sage: QQ[sqrt(-23)].maximal_order().class_number() 3 sage: ZZ[120*sqrt(-23)].class_number() 288
Note that non-maximal orders are only supported in quadratic fields:
sage: ZZ[120*sqrt(-23)].class_number() 288 sage: ZZ[100*sqrt(3)].class_number() 4 sage: ZZ[11*2^(1/3)].class_number() Traceback (most recent call last): ... NotImplementedError: computation of class numbers of non-maximal orders not in quadratic fields is not implemented
- coordinates(x)#
Return the coordinate vector of \(x\) with respect to this order.
INPUT:
x
– an element of the number field of this order.
OUTPUT:
A vector of length \(n\) (the degree of the field) giving the coordinates of \(x\) with respect to the integral basis of the order. In general this will be a vector of rationals; it will consist of integers if and only if \(x\) is in the order.
AUTHOR: John Cremona 2008-11-15
ALGORITHM:
Uses linear algebra. The change-of-basis matrix is cached. Provides simpler implementations for
_contains_()
,is_integral()
andsmallest_integer()
.EXAMPLES:
sage: K.<i> = QuadraticField(-1) sage: OK = K.ring_of_integers() sage: OK_basis = OK.basis(); OK_basis [1, i] sage: a = 23-14*i sage: acoords = OK.coordinates(a); acoords (23, -14) sage: sum([OK_basis[j]*acoords[j] for j in range(2)]) == a True sage: OK.coordinates((120+340*i)/8) (15, 85/2) sage: O = K.order(3*i) sage: O.is_maximal() False sage: O.index_in(OK) 3 sage: acoords = O.coordinates(a); acoords (23, -14/3) sage: sum([O.basis()[j]*acoords[j] for j in range(2)]) == a True
- degree()#
Return the degree of this order, which is the rank of this order as a \(\ZZ\)-module.
EXAMPLES:
sage: k.<c> = NumberField(x^3 + x^2 - 2*x+8) sage: o = k.maximal_order() sage: o.degree() 3 sage: o.rank() 3
- fraction_field()#
Return the fraction field of this order, which is the ambient number field.
EXAMPLES:
sage: K.<b> = NumberField(x^4 + 17*x^2 + 17) sage: O = K.order(17*b); O Order in Number Field in b with defining polynomial x^4 + 17*x^2 + 17 sage: O.fraction_field() Number Field in b with defining polynomial x^4 + 17*x^2 + 17
- fractional_ideal(*args, **kwds)#
Return the fractional ideal of the maximal order with given generators.
EXAMPLES:
sage: K.<a> = NumberField(x^2 + 2) sage: R = K.maximal_order() sage: R.fractional_ideal(2/3 + 7*a, a) Fractional ideal (1/3*a)
- free_module()#
Return the free \(\ZZ\)-module contained in the vector space associated to the ambient number field, that corresponds to this order.
EXAMPLES:
sage: K.<a> = NumberField(x^3 + x^2 - 2*x + 8) sage: O = K.maximal_order(); O.basis() [1, 1/2*a^2 + 1/2*a, a^2] sage: O.free_module() Free module of degree 3 and rank 3 over Integer Ring User basis matrix: [ 1 0 0] [ 0 1/2 1/2] [ 0 0 1]
An example in a relative extension. Notice that the module is a \(\ZZ\)-module in the absolute_field associated to the relative field:
sage: K.<a,b> = NumberField([x^2 + 1, x^2 + 2]) sage: O = K.maximal_order(); O.basis() [(-3/2*b - 5)*a + 7/2*b - 2, -3*a + 2*b, -2*b*a - 3, -7*a + 5*b] sage: O.free_module() Free module of degree 4 and rank 4 over Integer Ring User basis matrix: [1/4 1/4 3/4 3/4] [ 0 1/2 0 1/2] [ 0 0 1 0] [ 0 0 0 1]
- gen(i)#
Return \(i\)’th module generator of this order.
EXAMPLES:
sage: K.<c> = NumberField(x^3 + 2*x + 17) sage: O = K.maximal_order(); O Maximal Order in Number Field in c with defining polynomial x^3 + 2*x + 17 sage: O.basis() [1, c, c^2] sage: O.gen(1) c sage: O.gen(2) c^2 sage: O.gen(5) Traceback (most recent call last): ... IndexError: no 5th generator sage: O.gen(-1) Traceback (most recent call last): ... IndexError: no -1th generator
- ideal(*args, **kwds)#
Return the integral ideal with given generators.
EXAMPLES:
sage: K.<a> = NumberField(x^2 + 7) sage: R = K.maximal_order() sage: R.ideal(2/3 + 7*a, a) Traceback (most recent call last): ... ValueError: ideal must be integral; use fractional_ideal to create a non-integral ideal. sage: R.ideal(7*a, 77 + 28*a) Fractional ideal (7) sage: R = K.order(4*a) sage: R.ideal(8) Traceback (most recent call last): ... NotImplementedError: ideals of non-maximal orders not yet supported.
This function is called implicitly below:
sage: R = EquationOrder(x^2 + 2, 'a'); R Maximal Order in Number Field in a with defining polynomial x^2 + 2 sage: (3,15)*R Fractional ideal (3)
The zero ideal is handled properly:
sage: R.ideal(0) Ideal (0) of Number Field in a with defining polynomial x^2 + 2
- integral_closure()#
Return the integral closure of this order.
EXAMPLES:
sage: K.<a> = QuadraticField(5) sage: O2 = K.order(2*a); O2 Order in Number Field in a with defining polynomial x^2 - 5 with a = 2.236067977499790? sage: O2.integral_closure() Maximal Order in Number Field in a with defining polynomial x^2 - 5 with a = 2.236067977499790? sage: OK = K.maximal_order() sage: OK is OK.integral_closure() True
- is_field(proof=True)#
Return
False
(because an order is never a field).EXAMPLES:
sage: L.<alpha> = NumberField(x**4 - x**2 + 7) sage: O = L.maximal_order() ; O.is_field() False sage: CyclotomicField(12).ring_of_integers().is_field() False
- is_integrally_closed()#
Return
True
if this ring is integrally closed, i.e., is equal to the maximal order.EXAMPLES:
sage: K.<a> = NumberField(x^2 + 189*x + 394) sage: R = K.order(2*a) sage: R.is_integrally_closed() False sage: R Order in Number Field in a with defining polynomial x^2 + 189*x + 394 sage: S = K.maximal_order(); S Maximal Order in Number Field in a with defining polynomial x^2 + 189*x + 394 sage: S.is_integrally_closed() True
- is_noetherian()#
Return
True
(because orders are always Noetherian)EXAMPLES:
sage: L.<alpha> = NumberField(x**4 - x**2 + 7) sage: O = L.maximal_order() ; O.is_noetherian() True sage: E.<w> = NumberField(x^2 - x + 2) sage: OE = E.ring_of_integers(); OE.is_noetherian() True
- is_suborder(other)#
Return True if self and other are both orders in the same ambient number field and self is a subset of other.
EXAMPLES:
sage: W.<i> = NumberField(x^2 + 1) sage: O5 = W.order(5*i) sage: O10 = W.order(10*i) sage: O15 = W.order(15*i) sage: O15.is_suborder(O5) True sage: O5.is_suborder(O15) False sage: O10.is_suborder(O15) False
We create another isomorphic but different field:
sage: W2.<j> = NumberField(x^2 + 1) sage: P5 = W2.order(5*j)
This is False because the ambient number fields are not equal.:
sage: O5.is_suborder(P5) False
We create a field that contains (in no natural way!) W, and of course again is_suborder returns False:
sage: K.<z> = NumberField(x^4 + 1) sage: M = K.order(5*z) sage: O5.is_suborder(M) False
- krull_dimension()#
Return the Krull dimension of this order, which is 1.
EXAMPLES:
sage: K.<a> = QuadraticField(5) sage: OK = K.maximal_order() sage: OK.krull_dimension() 1 sage: O2 = K.order(2*a) sage: O2.krull_dimension() 1
- ngens()#
Return the number of module generators of this order.
EXAMPLES:
sage: K.<a> = NumberField(x^3 + x^2 - 2*x + 8) sage: O = K.maximal_order() sage: O.ngens() 3
- number_field()#
Return the number field of this order, which is the ambient number field that this order is embedded in.
EXAMPLES:
sage: K.<b> = NumberField(x^4 + x^2 + 2) sage: O = K.order(2*b); O Order in Number Field in b with defining polynomial x^4 + x^2 + 2 sage: O.basis() [1, 2*b, 4*b^2, 8*b^3] sage: O.number_field() Number Field in b with defining polynomial x^4 + x^2 + 2 sage: O.number_field() is K True
- random_element(*args, **kwds)#
Return a random element of this order.
INPUT:
args
,kwds
– parameters passed to the random integer function. See the documentation forZZ.random_element()
for details.
OUTPUT:
A random element of this order, computed as a random \(\ZZ\)-linear combination of the basis.
EXAMPLES:
sage: K.<a> = NumberField(x^3 + 2) sage: OK = K.ring_of_integers() sage: OK.random_element() # random output -2*a^2 - a - 2 sage: OK.random_element(distribution="uniform") # random output -a^2 - 1 sage: OK.random_element(-10,10) # random output -10*a^2 - 9*a - 2 sage: K.order(a).random_element() # random output a^2 - a - 3
sage: K.<z> = CyclotomicField(17) sage: OK = K.ring_of_integers() sage: OK.random_element() # random output z^15 - z^11 - z^10 - 4*z^9 + z^8 + 2*z^7 + z^6 - 2*z^5 - z^4 - 445*z^3 - 2*z^2 - 15*z - 2 sage: OK.random_element().is_integral() True sage: OK.random_element().parent() is OK True
A relative example:
sage: K.<a, b> = NumberField([x^2 + 2, x^2 + 1000*x + 1]) sage: OK = K.ring_of_integers() sage: OK.random_element() # random output (42221/2*b + 61/2)*a + 7037384*b + 7041 sage: OK.random_element().is_integral() # random output True sage: OK.random_element().parent() is OK # random output True
An example in a non-maximal order:
sage: K.<a> = QuadraticField(-3) sage: R = K.ring_of_integers() sage: A = K.order(a) sage: A.index_in(R) 2 sage: R.random_element() # random output -39/2*a - 1/2 sage: A.random_element() # random output 2*a - 1 sage: A.random_element().is_integral() True sage: A.random_element().parent() is A True
- rank()#
Return the rank of this order, which is the rank of the underlying \(\ZZ\)-module, or the degree of the ambient number field that contains this order.
This is a synonym for
degree()
.EXAMPLES:
sage: k.<c> = NumberField(x^5 + x^2 + 1) sage: o = k.maximal_order(); o Maximal Order in Number Field in c with defining polynomial x^5 + x^2 + 1 sage: o.rank() 5
- residue_field(prime, names=None, check=False)#
Return the residue field of this order at a given prime, ie \(O/pO\).
INPUT:
prime
– a prime ideal of the maximal order in this number field.names
– the name of the variable in the residue fieldcheck
– whether or not to check the primality of prime.
OUTPUT:
The residue field at this prime.
EXAMPLES:
sage: R.<x> = QQ[] sage: K.<a> = NumberField(x^4+3*x^2-17) sage: P = K.ideal(61).factor()[0][0] sage: OK = K.maximal_order() sage: OK.residue_field(P) Residue field in abar of Fractional ideal (61, a^2 + 30) sage: Fp.<b> = OK.residue_field(P) sage: Fp Residue field in b of Fractional ideal (61, a^2 + 30)
- ring_generators()#
Return generators for self as a ring.
EXAMPLES:
sage: K.<i> = NumberField(x^2 + 1) sage: O = K.maximal_order(); O Gaussian Integers in Number Field in i with defining polynomial x^2 + 1 sage: O.ring_generators() [i]
This is an example where 2 generators are required (because 2 is an essential discriminant divisor).:
sage: K.<a> = NumberField(x^3 + x^2 - 2*x + 8) sage: O = K.maximal_order(); O.basis() [1, 1/2*a^2 + 1/2*a, a^2] sage: O.ring_generators() [1/2*a^2 + 1/2*a, a^2]
An example in a relative number field:
sage: K.<a, b> = NumberField([x^2 + x + 1, x^3 - 3]) sage: O = K.maximal_order() sage: O.ring_generators() [(-5/3*b^2 + 3*b - 2)*a - 7/3*b^2 + b + 3, (-5*b^2 - 9)*a - 5*b^2 - b, (-6*b^2 - 11)*a - 6*b^2 - b]
- some_elements()#
Return a list of elements of the given order.
EXAMPLES:
sage: G = GaussianIntegers(); G Gaussian Integers in Number Field in I with defining polynomial x^2 + 1 with I = 1*I sage: G.some_elements() [1, I, 2*I, -1, 0, -I, 2, 4*I, -2, -2*I, -4] sage: R.<t> = QQ[] sage: K.<a> = QQ.extension(t^3 - 2); K Number Field in a with defining polynomial t^3 - 2 sage: Z = K.ring_of_integers(); Z Maximal Order in Number Field in a with defining polynomial t^3 - 2 sage: Z.some_elements() [1, a, a^2, 2*a, 0, 2, a^2 + 2*a + 1, ..., a^2 + 1, 2*a^2 + 2, a^2 + 2*a, 4*a^2 + 4]
- valuation(p)#
Return the
p
-adic valuation on this order.EXAMPLES:
The valuation can be specified with an integer
prime
that is completely ramified or unramified:sage: K.<a> = NumberField(x^2 + 1) sage: O = K.order(2*a) sage: valuations.pAdicValuation(O, 2) 2-adic valuation sage: GaussianIntegers().valuation(2) 2-adic valuation
sage: GaussianIntegers().valuation(3) 3-adic valuation
A
prime
that factors into pairwise distinct factors, results in an error:sage: GaussianIntegers().valuation(5) Traceback (most recent call last): ... ValueError: The valuation Gauss valuation induced by 5-adic valuation does not approximate a unique extension of 5-adic valuation with respect to x^2 + 1
The valuation can also be selected by giving a valuation on the base ring that extends uniquely:
sage: CyclotomicField(5).ring_of_integers().valuation(ZZ.valuation(5)) 5-adic valuation
When the extension is not unique, this does not work:
sage: GaussianIntegers().valuation(ZZ.valuation(5)) Traceback (most recent call last): ... ValueError: The valuation Gauss valuation induced by 5-adic valuation does not approximate a unique extension of 5-adic valuation with respect to x^2 + 1
If the fraction field is of the form \(K[x]/(G)\), you can specify a valuation by providing a discrete pseudo-valuation on \(K[x]\) which sends \(G\) to infinity:
sage: R.<x> = QQ[] sage: v = GaussianIntegers().valuation(GaussValuation(R, QQ.valuation(5)).augmentation(x + 2, infinity)) sage: w = GaussianIntegers().valuation(GaussValuation(R, QQ.valuation(5)).augmentation(x + 1/2, infinity)) sage: v == w False
- zeta(n=2, all=False)#
Return a primitive n-th root of unity in this order, if it contains one. If all is True, return all of them.
EXAMPLES:
sage: F.<alpha> = NumberField(x**2+3) sage: F.ring_of_integers().zeta(6) -1/2*alpha + 1/2 sage: O = F.order([3*alpha]) sage: O.zeta(3) Traceback (most recent call last): ... ArithmeticError: there are no 3rd roots of unity in self
- class sage.rings.number_field.order.OrderFactory#
Bases:
UniqueFactory
Abstract base class for factories creating orders, such as
AbsoluteOrderFactory
andRelativeOrderFactory
.- get_object(version, key, extra_args)#
Create the order identified by
key
.This overrides the default implementation to update the maximality of the order if it was explicitly specified.
EXAMPLES:
Even though orders are unique parents, this lets us update their internal state when they are recreated with more additional information available about them:
sage: L.<a, b> = NumberField([x^2 - 1000003, x^2 - 5*1000099^2]) sage: O = L.maximal_order([2], assume_maximal=None) sage: O._is_maximal_at(2) True sage: O._is_maximal_at(3) is None True sage: N = L.maximal_order([3], assume_maximal=None) sage: N is O True sage: N._is_maximal_at(2) True sage: N._is_maximal_at(3) True
- class sage.rings.number_field.order.Order_absolute(K, module_rep)#
Bases:
Order
EXAMPLES:
sage: from sage.rings.number_field.order import * sage: x = polygen(QQ) sage: K.<a> = NumberField(x^3+2) sage: V, from_v, to_v = K.vector_space() sage: M = span([to_v(a^2), to_v(a), to_v(1)],ZZ) sage: O = AbsoluteOrder(K, M); O Maximal Order in Number Field in a with defining polynomial x^3 + 2 sage: M = span([to_v(a^2), to_v(a), to_v(2)],ZZ) sage: O = AbsoluteOrder(K, M); O Traceback (most recent call last): ... ValueError: 1 is not in the span of the module, hence not an order
- absolute_discriminant()#
Return the discriminant of this order.
EXAMPLES:
sage: K.<a> = NumberField(x^8 + x^3 - 13*x + 26) sage: O = K.maximal_order() sage: factor(O.discriminant()) 3 * 11 * 13^2 * 613 * 1575917857 sage: L = K.order(13*a^2) sage: factor(L.discriminant()) 3^3 * 5^2 * 11 * 13^60 * 613 * 733^2 * 1575917857 sage: factor(L.index_in(O)) 3 * 5 * 13^29 * 733 sage: L.discriminant() / O.discriminant() == L.index_in(O)^2 True
- absolute_order()#
Return the absolute order associated to this order, which is just this order again since this is an absolute order.
EXAMPLES:
sage: K.<a> = NumberField(x^3 + 2) sage: O1 = K.order(a); O1 Maximal Order in Number Field in a with defining polynomial x^3 + 2 sage: O1.absolute_order() is O1 True
- basis()#
Return the basis over \(\ZZ\) for this order.
EXAMPLES:
sage: k.<c> = NumberField(x^3 + x^2 + 1) sage: O = k.maximal_order(); O Maximal Order in Number Field in c with defining polynomial x^3 + x^2 + 1 sage: O.basis() [1, c, c^2]
The basis is an immutable sequence:
sage: type(O.basis()) <class 'sage.structure.sequence.Sequence_generic'>
The generator functionality uses the basis method:
sage: O.0 1 sage: O.1 c sage: O.basis() [1, c, c^2] sage: O.ngens() 3
- change_names(names)#
Return a new order isomorphic to this one in the number field with given variable names.
EXAMPLES:
sage: R = EquationOrder(x^3 + x + 1, 'alpha'); R Order in Number Field in alpha with defining polynomial x^3 + x + 1 sage: R.basis() [1, alpha, alpha^2] sage: S = R.change_names('gamma'); S Order in Number Field in gamma with defining polynomial x^3 + x + 1 sage: S.basis() [1, gamma, gamma^2]
- discriminant()#
Return the discriminant of this order.
EXAMPLES:
sage: K.<a> = NumberField(x^8 + x^3 - 13*x + 26) sage: O = K.maximal_order() sage: factor(O.discriminant()) 3 * 11 * 13^2 * 613 * 1575917857 sage: L = K.order(13*a^2) sage: factor(L.discriminant()) 3^3 * 5^2 * 11 * 13^60 * 613 * 733^2 * 1575917857 sage: factor(L.index_in(O)) 3 * 5 * 13^29 * 733 sage: L.discriminant() / O.discriminant() == L.index_in(O)^2 True
- index_in(other)#
Return the index of self in other.
This is a lattice index, so it is a rational number if self is not contained in other.
INPUT:
other
– another absolute order with the same ambient number field.
OUTPUT:
a rational number
EXAMPLES:
sage: k.<i> = NumberField(x^2 + 1) sage: O1 = k.order(i) sage: O5 = k.order(5*i) sage: O5.index_in(O1) 5 sage: k.<a> = NumberField(x^3 + x^2 - 2*x+8) sage: o = k.maximal_order() sage: o Maximal Order in Number Field in a with defining polynomial x^3 + x^2 - 2*x + 8 sage: O1 = k.order(a); O1 Order in Number Field in a with defining polynomial x^3 + x^2 - 2*x + 8 sage: O1.index_in(o) 2 sage: O2 = k.order(1+2*a); O2 Order in Number Field in a with defining polynomial x^3 + x^2 - 2*x + 8 sage: O1.basis() [1, a, a^2] sage: O2.basis() [1, 2*a, 4*a^2] sage: o.index_in(O2) 1/16
- intersection(other)#
Return the intersection of this order with another order.
EXAMPLES:
sage: k.<i> = NumberField(x^2 + 1) sage: O6 = k.order(6*i) sage: O9 = k.order(9*i) sage: O6.basis() [1, 6*i] sage: O9.basis() [1, 9*i] sage: O6.intersection(O9).basis() [1, 18*i] sage: (O6 & O9).basis() [1, 18*i] sage: (O6 + O9).basis() [1, 3*i]
- is_maximal(p=None)#
Return whether this is the maximal order.
INPUT:
p
– an integer prime orNone
(default:None
); if set, return whether this order is maximal at the primep
.
EXAMPLES:
sage: K.<i> = NumberField(x^2 + 1) sage: K.order(3*i).is_maximal() False sage: K.order(5*i).is_maximal() False sage: (K.order(3*i) + K.order(5*i)).is_maximal() True sage: K.maximal_order().is_maximal() True
Maximality can be checked at primes when the order is maximal at that prime by construction:
sage: K.maximal_order().is_maximal(p=3) True
And also at other primes:
sage: K.order(3*i).is_maximal(p=3) False An example involving a relative order:: sage: K.<a, b> = NumberField([x^2 + 1, x^2 - 3]) sage: O = K.order([3*a,2*b]) sage: O.is_maximal() False
- module()#
Return the underlying free module corresponding to this order, embedded in the vector space corresponding to the ambient number field.
EXAMPLES:
sage: k.<a> = NumberField(x^3 + x + 3) sage: m = k.order(3*a); m Order in Number Field in a with defining polynomial x^3 + x + 3 sage: m.module() Free module of degree 3 and rank 3 over Integer Ring Echelon basis matrix: [1 0 0] [0 3 0] [0 0 9]
- class sage.rings.number_field.order.Order_relative(K, absolute_order)#
Bases:
Order
A relative order in a number field.
A relative order is an order in some relative number field
Invariants of this order may be computed with respect to the contained order.
- absolute_discriminant()#
Return the absolute discriminant of self, which is the discriminant of the absolute order associated to self.
OUTPUT:
an integer
EXAMPLES:
sage: R = EquationOrder([x^2 + 1, x^3 + 2], 'a,b') sage: d = R.absolute_discriminant(); d -746496 sage: d is R.absolute_discriminant() True sage: factor(d) -1 * 2^10 * 3^6
- absolute_order(names='z')#
Return underlying absolute order associated to this relative order.
INPUT:
names
– string (default: ‘z’); name of generator of absolute extension.
Note
There is a default variable name, since this absolute order is frequently used for internal algorithms.
EXAMPLES:
sage: R = EquationOrder([x^2 + 1, x^2 - 5], 'i,g'); R Relative Order in Number Field in i with defining polynomial x^2 + 1 over its base field sage: R.basis() [1, 6*i - g, -g*i + 2, 7*i - g] sage: S = R.absolute_order(); S Order in Number Field in z with defining polynomial x^4 - 8*x^2 + 36 sage: S.basis() [1, 5/12*z^3 + 1/6*z, 1/2*z^2, 1/2*z^3]
We compute a relative order in alpha0, alpha1, then make the number field that contains the absolute order be called gamma.:
sage: R = EquationOrder( [x^2 + 2, x^2 - 3], 'alpha'); R Relative Order in Number Field in alpha0 with defining polynomial x^2 + 2 over its base field sage: R.absolute_order('gamma') Order in Number Field in gamma with defining polynomial x^4 - 2*x^2 + 25 sage: R.absolute_order('gamma').basis() [1/2*gamma^2 + 1/2, 7/10*gamma^3 + 1/10*gamma, gamma^2, gamma^3]
- basis()#
Return a basis for this order as \(\ZZ\)-module.
EXAMPLES:
sage: K.<a,b> = NumberField([x^2+1, x^2+3]) sage: O = K.order([a,b]) sage: O.basis() [1, -2*a + b, -b*a - 2, -5*a + 3*b] sage: z = O.1; z -2*a + b sage: z.absolute_minpoly() x^4 + 14*x^2 + 1
- index_in(other)#
Return the index of self in other.
This is a lattice index, so it is a rational number if self is not contained in other.
INPUT:
other
– another order with the same ambient absolute number field.
OUTPUT:
a rational number
EXAMPLES:
sage: K.<a,b> = NumberField([x^3 + x + 3, x^2 + 1]) sage: R1 = K.order([3*a, 2*b]) sage: R2 = K.order([a, 4*b]) sage: R1.index_in(R2) 729/8 sage: R2.index_in(R1) 8/729
- is_maximal(p=None)#
Return whether this is the maximal order.
INPUT:
p
– an integer prime orNone
(default:None
); if set, return whether this order is maximal at the primep
.
EXAMPLES:
sage: K.<a, b> = NumberField([x^2 + 1, x^2 - 5]) sage: K.order(3*a, b).is_maximal() False sage: K.order(5*a, b/2 + 1/2).is_maximal() False sage: (K.order(3*a, b) + K.order(5*a, b/2 + 1/2)).is_maximal() True sage: K.maximal_order().is_maximal() True
Maximality can be checked at primes when the order is maximal at that prime by construction:
sage: K.maximal_order().is_maximal(p=3) True
And at other primes:
sage: K.order(3*a, b).is_maximal(p=3) False
- is_suborder(other)#
Return True if self is a subset of the order other.
EXAMPLES:
sage: K.<a,b> = NumberField([x^2 + 1, x^3 + 2]) sage: R1 = K.order([a,b]) sage: R2 = K.order([2*a,b]) sage: R3 = K.order([a + b, b + 2*a]) sage: R1.is_suborder(R2) False sage: R2.is_suborder(R1) True sage: R3.is_suborder(R1) True sage: R1.is_suborder(R3) True sage: R1 == R3 True
- class sage.rings.number_field.order.RelativeOrderFactory#
Bases:
OrderFactory
An order in a relative number field extension.
EXAMPLES:
sage: K.<i> = NumberField(x^2 + 1) sage: R.<j> = K[] sage: L.<j> = K.extension(j^2 - 2) sage: L.order([i, j]) Relative Order in Number Field in j with defining polynomial j^2 - 2 over its base field
- create_key_and_extra_args(K, absolute_order, is_maximal=None, check=True, is_maximal_at=())#
Return normalized arguments to create a relative order.
- create_object(version, key, is_maximal=None, is_maximal_at=())#
Create a relative order.
- sage.rings.number_field.order.absolute_order_from_module_generators(gens, check_integral=True, check_rank=True, check_is_ring=True, is_maximal=None, allow_subfield=False, is_maximal_at=())#
INPUT:
gens
– list of elements of an absolute number field that generates an order in that number field as a ZZ module.check_integral
– check that each gen is integralcheck_rank
– check that the gens span a module of the correct rankcheck_is_ring
– check that the module is closed under multiplication (this is very expensive)is_maximal
– bool (or None); set if maximality of the generated order is knownis_maximal_at
– a tuple of primes where this order is known to be maximal
OUTPUT:
an absolute order
EXAMPLES:
We have to explicitly import the function, since it is not meant for regular usage:
sage: from sage.rings.number_field.order import absolute_order_from_module_generators sage: K.<a> = NumberField(x^4 - 5) sage: O = K.maximal_order(); O Maximal Order in Number Field in a with defining polynomial x^4 - 5 sage: O.basis() [1/2*a^2 + 1/2, 1/2*a^3 + 1/2*a, a^2, a^3] sage: O.module() Free module of degree 4 and rank 4 over Integer Ring Echelon basis matrix: [1/2 0 1/2 0] [ 0 1/2 0 1/2] [ 0 0 1 0] [ 0 0 0 1] sage: g = O.basis(); g [1/2*a^2 + 1/2, 1/2*a^3 + 1/2*a, a^2, a^3] sage: absolute_order_from_module_generators(g) Maximal Order in Number Field in a with defining polynomial x^4 - 5
We illustrate each check flag – the output is the same but in case the function would run ever so slightly faster:
sage: absolute_order_from_module_generators(g, check_is_ring=False) Maximal Order in Number Field in a with defining polynomial x^4 - 5 sage: absolute_order_from_module_generators(g, check_rank=False) Maximal Order in Number Field in a with defining polynomial x^4 - 5 sage: absolute_order_from_module_generators(g, check_integral=False) Maximal Order in Number Field in a with defining polynomial x^4 - 5
Next we illustrate constructing “fake” orders to illustrate turning off various check flags:
sage: k.<i> = NumberField(x^2 + 1) sage: R = absolute_order_from_module_generators([2, 2*i], check_is_ring=False); R Order in Number Field in i with defining polynomial x^2 + 1 sage: R.basis() [2, 2*i] sage: R = absolute_order_from_module_generators([k(1)], check_rank=False); R Order in Number Field in i with defining polynomial x^2 + 1 sage: R.basis() [1]
If the order contains a non-integral element, even if we do not check that, we will find that the rank is wrong or that the order is not closed under multiplication:
sage: absolute_order_from_module_generators([1/2, i], check_integral=False) Traceback (most recent call last): ... ValueError: the module span of the gens is not closed under multiplication. sage: R = absolute_order_from_module_generators([1/2, i], check_is_ring=False, check_integral=False); R Order in Number Field in i with defining polynomial x^2 + 1 sage: R.basis() [1/2, i]
We turn off all check flags and make a really messed up order:
sage: R = absolute_order_from_module_generators([1/2, i], check_is_ring=False, check_integral=False, check_rank=False); R Order in Number Field in i with defining polynomial x^2 + 1 sage: R.basis() [1/2, i]
An order that lives in a subfield:
sage: F.<alpha> = NumberField(x**4+3) sage: F.order([alpha**2], allow_subfield=True) Order in Number Field in beta with defining polynomial ... with beta = ...
- sage.rings.number_field.order.absolute_order_from_ring_generators(gens, check_is_integral=True, check_rank=True, is_maximal=None, allow_subfield=False)#
INPUT:
gens
– list of integral elements of an absolute order.check_is_integral
– bool (default: True), whether to check that each generator is integral.check_rank
– bool (default: True), whether to check that the ring generated by gens is of full rank.is_maximal
– bool (or None); set if maximality of the generated order is knownallow_subfield
– bool (default: False), if True and the generators do not generate an order, i.e., they generate a subring of smaller rank, instead of raising an error, return an order in a smaller number field.
EXAMPLES:
sage: K.<a> = NumberField(x^4 - 5) sage: K.order(a) Order in Number Field in a with defining polynomial x^4 - 5
We have to explicitly import this function, since typically it is called with
K.order
as above.:sage: from sage.rings.number_field.order import absolute_order_from_ring_generators sage: absolute_order_from_ring_generators([a]) Order in Number Field in a with defining polynomial x^4 - 5 sage: absolute_order_from_ring_generators([3*a, 2, 6*a+1]) Order in Number Field in a with defining polynomial x^4 - 5
If one of the inputs is non-integral, it is an error.:
sage: absolute_order_from_ring_generators([a/2]) Traceback (most recent call last): ... ValueError: each generator must be integral
If the gens do not generate an order, i.e., generate a ring of full rank, then it is an error.:
sage: absolute_order_from_ring_generators([a^2]) Traceback (most recent call last): ... ValueError: the rank of the span of gens is wrong
Both checking for integrality and checking for full rank can be turned off in order to save time, though one can get nonsense as illustrated below.:
sage: absolute_order_from_ring_generators([a/2], check_is_integral=False) Order in Number Field in a with defining polynomial x^4 - 5 sage: absolute_order_from_ring_generators([a^2], check_rank=False) Order in Number Field in a with defining polynomial x^4 - 5
- sage.rings.number_field.order.each_is_integral(v)#
Return whether every element of the list
v
of elements of a number field is integral.EXAMPLES:
sage: W.<sqrt5> = NumberField(x^2 - 5) sage: from sage.rings.number_field.order import each_is_integral sage: each_is_integral([sqrt5, 2, (1+sqrt5)/2]) True sage: each_is_integral([sqrt5, (1+sqrt5)/3]) False
- sage.rings.number_field.order.is_NumberFieldOrder(R)#
Return True if R is either an order in a number field or is the ring \(\ZZ\) of integers.
EXAMPLES:
sage: from sage.rings.number_field.order import is_NumberFieldOrder sage: is_NumberFieldOrder(NumberField(x^2+1,'a').maximal_order()) True sage: is_NumberFieldOrder(ZZ) True sage: is_NumberFieldOrder(QQ) False sage: is_NumberFieldOrder(45) False
- sage.rings.number_field.order.relative_order_from_ring_generators(gens, check_is_integral=True, check_rank=True, is_maximal=None, allow_subfield=False, is_maximal_at=())#
INPUT:
gens
– list of integral elements of an absolute order.check_is_integral
– bool (default: True), whether to check that each generator is integral.check_rank
– bool (default: True), whether to check that the ring generated by gens is of full rank.is_maximal
– bool (or None); set if maximality of the generated order is known
EXAMPLES:
We have to explicitly import this function, since it is not meant for regular usage:
sage: from sage.rings.number_field.order import relative_order_from_ring_generators sage: K.<i, a> = NumberField([x^2 + 1, x^2 - 17]) sage: R = K.base_field().maximal_order() sage: S = relative_order_from_ring_generators([i,a]); S Relative Order in Number Field in i with defining polynomial x^2 + 1 over its base field
Basis for the relative order, which is obtained by computing the algebra generated by i and a:
sage: S.basis() [1, 7*i - 2*a, -a*i + 8, 25*i - 7*a]