Localization¶
Localization is an important ring construction tool. Whenever you have to extend a given integral domain such that it contains the inverses of a finite set of elements but should allow non injective homomorphic images this construction will be needed. See the example on Ariki-Koike algebras below for such an application.
EXAMPLES:
sage: LZ = Localization(ZZ, (5,11))
sage: m = matrix(LZ, [[5, 7], [0,11]])
sage: m.parent()
Full MatrixSpace of 2 by 2 dense matrices over Integer Ring localized at (5, 11)
sage: ~m # parent of inverse is different: see documentation of m.__invert__
[ 1/5 -7/55]
[ 0 1/11]
sage: _.parent()
Full MatrixSpace of 2 by 2 dense matrices over Rational Field
sage: mi = matrix(LZ, ~m)
sage: mi.parent()
Full MatrixSpace of 2 by 2 dense matrices over Integer Ring localized at (5, 11)
sage: mi == ~m
True
The next example defines the most general ring containing the coefficients of the irreducible representations of the Ariki-Koike algebra corresponding to the three colored permutations on three elements:
sage: R.<u0, u1, u2, q> = ZZ[]
sage: u = [u0, u1, u2]
sage: S = Set(u)
sage: I = S.cartesian_product(S)
sage: add_units = u + [q, q+1] + [ui -uj for ui, uj in I if ui != uj]\
+ [q*ui -uj for ui, uj in I if ui != uj]
sage: L = R.localization(tuple(add_units)); L
Multivariate Polynomial Ring in u0, u1, u2, q over Integer Ring localized at
(q, q + 1, u2, u1, u1 - u2, u0, u0 - u2, u0 - u1, u2*q - u1, u2*q - u0,
u1*q - u2, u1*q - u0, u0*q - u2, u0*q - u1)
Define the representation matrices (of one of the three dimensional irreducible representations):
sage: m1 = matrix(L, [[u1, 0, 0],[0, u0, 0],[0, 0, u0]])
sage: m2 = matrix(L, [[(u0*q - u0)/(u0 - u1), (u0*q - u1)/(u0 - u1), 0],\
[(-u1*q + u0)/(u0 - u1), (-u1*q + u1)/(u0 - u1), 0],\
[0, 0, -1]])
sage: m3 = matrix(L, [[-1, 0, 0],\
[0, u0*(1 - q)/(u1*q - u0), q*(u1 - u0)/(u1*q - u0)],\
[0, (u1*q^2 - u0)/(u1*q - u0), (u1*q^ 2 - u1*q)/(u1*q - u0)]])
sage: m1.base_ring() == L
True
Check relations of the Ariki-Koike algebra:
sage: m1*m2*m1*m2 == m2*m1*m2*m1
True
sage: m2*m3*m2 == m3*m2*m3
True
sage: m1*m3 == m3*m1
True
sage: m1**3 -(u0+u1+u2)*m1**2 +(u0*u1+u0*u2+u1*u2)*m1 - u0*u1*u2 == 0
True
sage: m2**2 -(q-1)*m2 - q == 0
True
sage: m3**2 -(q-1)*m3 - q == 0
True
sage: ~m1 in m1.parent()
True
sage: ~m2 in m2.parent()
True
sage: ~m3 in m3.parent()
True
Obtain specializations in positive characteristic:
sage: Fp = GF(17)
sage: f = L.hom((3,5,7,11), codomain=Fp); f
Ring morphism:
From: Multivariate Polynomial Ring in u0, u1, u2, q over Integer Ring localized at
(q, q + 1, u2, u1, u1 - u2, u0, u0 - u2, u0 - u1, u2*q - u1, u2*q - u0,
u1*q - u2, u1*q - u0, u0*q - u2, u0*q - u1)
To: Finite Field of size 17
Defn: u0 |--> 3
u1 |--> 5
u2 |--> 7
q |--> 11
sage: mFp1 = matrix({k:f(v) for k, v in m1.dict().items()}); mFp1
[5 0 0]
[0 3 0]
[0 0 3]
sage: mFp1.base_ring()
Finite Field of size 17
sage: mFp2 = matrix({k:f(v) for k, v in m2.dict().items()}); mFp2
[ 2 3 0]
[ 9 8 0]
[ 0 0 16]
sage: mFp3 = matrix({k:f(v) for k, v in m3.dict().items()}); mFp3
[16 0 0]
[ 0 4 5]
[ 0 7 6]
Obtain specializations in characteristic 0:
sage: fQ = L.hom((3,5,7,11), codomain=QQ); fQ
Ring morphism:
From: Multivariate Polynomial Ring in u0, u1, u2, q over Integer Ring localized at
(q, q + 1, u2, u1, u1 - u2, u0, u0 - u2, u0 - u1, u2*q - u1, u2*q - u0,
u1*q - u2, u1*q - u0, u0*q - u2, u0*q - u1)
To: Rational Field
Defn: u0 |--> 3
u1 |--> 5
u2 |--> 7
q |--> 11
sage: mQ1 = matrix({k:fQ(v) for k, v in m1.dict().items()}); mQ1
[5 0 0]
[0 3 0]
[0 0 3]
sage: mQ1.base_ring()
Rational Field
sage: mQ2 = matrix({k:fQ(v) for k, v in m2.dict().items()}); mQ2
[-15 -14 0]
[ 26 25 0]
[ 0 0 -1]
sage: mQ3 = matrix({k:fQ(v) for k, v in m3.dict().items()}); mQ3
[ -1 0 0]
[ 0 -15/26 11/26]
[ 0 301/26 275/26]
sage: S.<x, y, z, t> = QQ[]
sage: T = S.quo(x+y+z)
sage: F = T.fraction_field()
sage: fF = L.hom((x, y, z, t), codomain=F); fF
Ring morphism:
From: Multivariate Polynomial Ring in u0, u1, u2, q over Integer Ring localized at
(q, q + 1, u2, u1, u1 - u2, u0, u0 - u2, u0 - u1, u2*q - u1, u2*q - u0,
u1*q - u2, u1*q - u0, u0*q - u2, u0*q - u1)
To: Fraction Field of Quotient of Multivariate Polynomial Ring in x, y, z, t over
Rational Field by the ideal (x + y + z)
Defn: u0 |--> -ybar - zbar
u1 |--> ybar
u2 |--> zbar
q |--> tbar
sage: mF1 = matrix({k:fF(v) for k, v in m1.dict().items()}); mF1
[ ybar 0 0]
[ 0 -ybar - zbar 0]
[ 0 0 -ybar - zbar]
sage: mF1.base_ring() == F
True
AUTHORS:
Sebastian Oehms 2019-12-09: initial version.
- class sage.rings.localization.Localization(base_ring, additional_units, names=None, normalize=True, category=None, warning=True)¶
Bases:
sage.rings.ring.IntegralDomain
,sage.structure.unique_representation.UniqueRepresentation
The localization generalizes the construction of the field of fractions of an integral domain to an arbitrary ring. Given a (not necessarily commutative) ring \(R\) and a subset \(S\) of \(R\), there exists a ring \(R[S^{-1}]\) together with the ring homomorphism \(R \longrightarrow R[S^{-1}]\) that “inverts” \(S\); that is, the homomorphism maps elements in \(S\) to unit elements in \(R[S^{-1}]\) and, moreover, any ring homomorphism from \(R\) that “inverts” \(S\) uniquely factors through \(R[S^{-1}]\).
The ring \(R[S^{-1}]\) is called the localization of \(R\) with respect to \(S\). For example, if \(R\) is a commutative ring and \(f\) an element in \(R\), then the localization consists of elements of the form \(r/f, r\in R, n \geq 0\) (to be precise, \(R[f^{-1}] = R[t]/(ft-1)\).
The above text is taken from \(Wikipedia\). The construction here used for this class relies on the construction of the field of fraction and is therefore restricted to integral domains.
Accordingly, this class is inherited from
IntegralDomain
and can only be used in that context. Furthermore, the base ring should supportsage.structure.element.CommutativeRingElement.divides()
and the exact division operator \(//\) (sage.structure.element.Element.__floordiv__()
) in order to guarantee an successful application.INPUT:
base_ring
– an instance ofRing
allowing the construction offraction_field()
(that is an integral domain)additional_units
– tuple of elements ofbase_ring
which should be turned into unitsnames
– passed toIntegralDomain
normalize
– (optional, default: True) passed toIntegralDomain
category
– (optional, default: None) passed toIntegralDomain
warning
– (optional, default: True) to suppress a warning which is thrown if self cannot be represented uniquely
REFERENCES:
EXAMPLES:
sage: L = Localization(ZZ, (3,5)) sage: 1/45 in L True sage: 1/43 in L False sage: Localization(L, (7,11)) Integer Ring localized at (3, 5, 7, 11) sage: _.is_subring(QQ) True sage: L(~7) Traceback (most recent call last): ... ValueError: factor 7 of denominator is not a unit sage: Localization(Zp(7), (3, 5)) Traceback (most recent call last): ... ValueError: all given elements are invertible in 7-adic Ring with capped relative precision 20 sage: R.<x> = ZZ[] sage: L = R.localization(x**2+1) sage: s = (x+5)/(x**2+1) sage: s in L True sage: t = (x+5)/(x**2+2) sage: t in L False sage: L(t) Traceback (most recent call last): ... TypeError: fraction must have unit denominator sage: L(s) in R False sage: y = L(x) sage: g = L(s) sage: g.parent() Univariate Polynomial Ring in x over Integer Ring localized at (x^2 + 1,) sage: f = (y+5)/(y**2+1); f (x + 5)/(x^2 + 1) sage: f == g True sage: (y+5)/(y**2+2) Traceback (most recent call last): ... ValueError: factor x^2 + 2 of denominator is not a unit
More examples will be shown typing
sage.rings.localization?
- Element¶
alias of
LocalizationElement
- characteristic()¶
Return the characteristic of
self
.EXAMPLES:
sage: R.<a> = GF(5)[] sage: L = R.localization((a**2-3, a)) sage: L.characteristic() 5
- fraction_field()¶
Return the fraction field of
self
.EXAMPLES:
sage: R.<a> = GF(5)[] sage: L = Localization(R, (a**2-3, a)) sage: L.fraction_field() Fraction Field of Univariate Polynomial Ring in a over Finite Field of size 5 sage: L.is_subring(_) True
- gen(i)¶
Return the
i
-th generator ofself
which is thei
-th generator of the base ring.EXAMPLES:
sage: R.<x, y> = ZZ[] sage: R.localization((x**2+1, y-1)).gen(0) x sage: ZZ.localization(2).gen(0) 1
- gens()¶
Return a tuple whose entries are the generators for this object, in order.
EXAMPLES:
sage: R.<x, y> = ZZ[] sage: Localization(R, (x**2+1, y-1)).gens() (x, y) sage: Localization(ZZ, 2).gens() (1,)
- is_field(proof=True)¶
Return
True
if this ring is a field.INPUT:
proof
– (default:True
) Determines what to do in unknown cases
ALGORITHM:
If the parameter
proof
is set toTrue
, the returned value is correct but the method might throw an error. Otherwise, if it is set toFalse
, the method returns True if it can establish that self is a field and False otherwise.EXAMPLES:
sage: R = ZZ.localization((2,3)) sage: R.is_field() False
- krull_dimension()¶
Return the Krull dimension of this localization.
Since the current implementation just allows integral domains as base ring and localization at a finite set of elements the spectrum of
self
is open in the irreducible spectrum of its base ring. Therefore, by density we may take the dimension from there.EXAMPLES:
sage: R = ZZ.localization((2,3)) sage: R.krull_dimension() 1
- ngens()¶
Return the number of generators of
self
according to the same method for the base ring.EXAMPLES:
sage: R.<x, y> = ZZ[] sage: Localization(R, (x**2+1, y-1)).ngens() 2 sage: Localization(ZZ, 2).ngens() 1
- class sage.rings.localization.LocalizationElement(parent, x)¶
Bases:
sage.structure.element.IntegralDomainElement
Element class for localizations of integral domains
INPUT:
parent
– instance ofLocalization
x
– instance ofFractionFieldElement
whose parent is the fractionfield of the parent’s base ring
EXAMPLES:
sage: from sage.rings.localization import LocalizationElement sage: P.<x,y,z> = GF(5)[] sage: L = P.localization((x, y*z-x)) sage: LocalizationElement(L, 4/(y*z-x)**2) (-1)/(y^2*z^2 - 2*x*y*z + x^2) sage: _.parent() Multivariate Polynomial Ring in x, y, z over Finite Field of size 5 localized at (x, y*z - x)
- denominator()¶
Return the denominator of
self
.EXAMPLES:
sage: L = Localization(ZZ, (3,5)) sage: L(7/15).denominator() 15
- inverse_of_unit()¶
Return the inverse of
self
.EXAMPLES:
sage: P.<x,y,z> = ZZ[] sage: L = Localization(P, x*y*z) sage: L(x*y*z).inverse_of_unit() 1/(x*y*z) sage: L(z).inverse_of_unit() 1/z
- is_unit()¶
Return
True
ifself
is a unit.EXAMPLES:
sage: P.<x,y,z> = QQ[] sage: L = P.localization((x, y*z)) sage: L(y*z).is_unit() True sage: L(z).is_unit() True sage: L(x*y*z).is_unit() True
- numerator()¶
Return the numerator of
self
.EXAMPLES:
sage: L = ZZ.localization((3,5)) sage: L(7/15).numerator() 7
- sage.rings.localization.normalize_additional_units(base_ring, add_units, warning=True)¶
Function to normalize input data.
The given list will be replaced by a list of the involved prime factors (if possible).
INPUT:
base_ring
– an instance ofIntegralDomain
add_units
– list of elements from base ringwarning
– (optional, default: True) to suppress a warning which is thrown if no normalization was possible
OUTPUT:
List of all prime factors of the elements of the given list.
EXAMPLES:
sage: from sage.rings.localization import normalize_additional_units sage: normalize_additional_units(ZZ, [3, -15, 45, 9, 2, 50]) [2, 3, 5] sage: P.<x,y,z> = ZZ[] sage: normalize_additional_units(P, [3*x, z*y**2, 2*z, 18*(x*y*z)**2, x*z, 6*x*z, 5]) [2, 3, 5, z, y, x] sage: P.<x,y,z> = QQ[] sage: normalize_additional_units(P, [3*x, z*y**2, 2*z, 18*(x*y*z)**2, x*z, 6*x*z, 5]) [z, y, x] sage: R.<x, y> = ZZ[] sage: Q.<a, b> = R.quo(x**2-5) sage: p = b**2-5 sage: p == (b-a)*(b+a) True sage: normalize_additional_units(Q, [p]) doctest:...: UserWarning: Localization may not be represented uniquely [b^2 - 5] sage: normalize_additional_units(Q, [p], warning=False) [b^2 - 5]