Path: blob/master/src/sage/algebras/letterplace/letterplace_ideal.pyx
8822 views
###############################################################################1#2# Copyright (C) 2011 Simon King <[email protected]>3# Distributed under the terms of the GNU General Public License (GPL),4# version 2 or any later version. The full text of the GPL is available at:5# http://www.gnu.org/licenses/6#7###############################################################################89"""10Homogeneous ideals of free algebras.1112For twosided ideals and when the base ring is a field, this13implementation also provides Groebner bases and ideal containment14tests.1516EXAMPLES::1718sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')19sage: F20Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field21sage: I = F*[x*y+y*z,x^2+x*y-y*x-y^2]*F22sage: I23Twosided Ideal (x*y + y*z, x*x + x*y - y*x - y*y) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field2425One can compute Groebner bases out to a finite degree, can compute normal26forms and can test containment in the ideal::2728sage: I.groebner_basis(degbound=3)29Twosided Ideal (y*y*y - y*y*z + y*z*y - y*z*z, y*y*x + y*y*z + y*z*x + y*z*z, x*y + y*z, x*x - y*x - y*y - y*z) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field30sage: (x*y*z*y*x).normal_form(I)31y*z*z*y*z + y*z*z*z*x + y*z*z*z*z32sage: x*y*z*y*x - (x*y*z*y*x).normal_form(I) in I33True3435AUTHOR:3637- Simon King (2011-03-22): See :trac:`7797`.3839"""4041from sage.rings.noncommutative_ideals import Ideal_nc42from sage.libs.singular.function import lib, singular_function43from sage.algebras.letterplace.free_algebra_letterplace cimport FreeAlgebra_letterplace44from sage.algebras.letterplace.free_algebra_element_letterplace cimport FreeAlgebraElement_letterplace45from sage.all import Infinity4647#####################48# Define some singular functions49lib("freegb.lib")50singular_system=singular_function("system")51poly_reduce=singular_function("NF")5253class LetterplaceIdeal(Ideal_nc):54"""55Graded homogeneous ideals in free algebras.5657In the two-sided case over a field, one can compute Groebner bases58up to a degree bound, normal forms of graded homogeneous elements59of the free algebra, and ideal containment.6061EXAMPLES::6263sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')64sage: I = F*[x*y+y*z,x^2+x*y-y*x-y^2]*F65sage: I66Twosided Ideal (x*y + y*z, x*x + x*y - y*x - y*y) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field67sage: I.groebner_basis(2)68Twosided Ideal (x*y + y*z, x*x - y*x - y*y - y*z) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field69sage: I.groebner_basis(4)70Twosided Ideal (y*z*y*y - y*z*y*z + y*z*z*y - y*z*z*z, y*z*y*x + y*z*y*z + y*z*z*x + y*z*z*z, y*y*z*y - y*y*z*z + y*z*z*y - y*z*z*z, y*y*z*x + y*y*z*z + y*z*z*x + y*z*z*z, y*y*y - y*y*z + y*z*y - y*z*z, y*y*x + y*y*z + y*z*x + y*z*z, x*y + y*z, x*x - y*x - y*y - y*z) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field7172Groebner bases are cached. If one has computed a Groebner basis73out to a high degree then it will also be returned if a Groebner74basis with a lower degree bound is requested::7576sage: I.groebner_basis(2)77Twosided Ideal (y*z*y*y - y*z*y*z + y*z*z*y - y*z*z*z, y*z*y*x + y*z*y*z + y*z*z*x + y*z*z*z, y*y*z*y - y*y*z*z + y*z*z*y - y*z*z*z, y*y*z*x + y*y*z*z + y*z*z*x + y*z*z*z, y*y*y - y*y*z + y*z*y - y*z*z, y*y*x + y*y*z + y*z*x + y*z*z, x*y + y*z, x*x - y*x - y*y - y*z) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field7879Of course, the normal form of any element has to satisfy the following::8081sage: x*y*z*y*x - (x*y*z*y*x).normal_form(I) in I82True8384Left and right ideals can be constructed, but only twosided ideals provide85Groebner bases::8687sage: JL = F*[x*y+y*z,x^2+x*y-y*x-y^2]; JL88Left Ideal (x*y + y*z, x*x + x*y - y*x - y*y) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field89sage: JR = [x*y+y*z,x^2+x*y-y*x-y^2]*F; JR90Right Ideal (x*y + y*z, x*x + x*y - y*x - y*y) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field91sage: JR.groebner_basis(2)92Traceback (most recent call last):93...94TypeError: This ideal is not two-sided. We can only compute two-sided Groebner bases95sage: JL.groebner_basis(2)96Traceback (most recent call last):97...98TypeError: This ideal is not two-sided. We can only compute two-sided Groebner bases99100Also, it is currently not possible to compute a Groebner basis when the base101ring is not a field::102103sage: FZ.<a,b,c> = FreeAlgebra(ZZ, implementation='letterplace')104sage: J = FZ*[a^3-b^3]*FZ105sage: J.groebner_basis(2)106Traceback (most recent call last):107...108TypeError: Currently, we can only compute Groebner bases if the ring of coefficients is a field109110The letterplace implementation of free algebras also provides integral degree weights111for the generators, and we can compute Groebner bases for twosided graded homogeneous112ideals::113114sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace',degrees=[1,2,3])115sage: I = F*[x*y+z-y*x,x*y*z-x^6+y^3]*F116sage: I.groebner_basis(Infinity)117Twosided Ideal (x*z*z - y*x*x*z - y*x*y*y + y*x*z*x + y*y*y*x + z*x*z + z*y*y - z*z*x,118x*y - y*x + z,119x*x*x*x*z*y*y + x*x*x*z*y*y*x - x*x*x*z*y*z - x*x*z*y*x*z + x*x*z*y*y*x*x +120x*x*z*y*y*y - x*x*z*y*z*x - x*z*y*x*x*z - x*z*y*x*z*x +121x*z*y*y*x*x*x + 2*x*z*y*y*y*x - 2*x*z*y*y*z - x*z*y*z*x*x -122x*z*y*z*y + y*x*z*x*x*x*x*x - 4*y*x*z*x*x*z - 4*y*x*z*x*z*x +1234*y*x*z*y*x*x*x + 3*y*x*z*y*y*x - 4*y*x*z*y*z + y*y*x*x*x*x*z +124y*y*x*x*x*z*x - 3*y*y*x*x*z*x*x - y*y*x*x*z*y +1255*y*y*x*z*x*x*x + 4*y*y*x*z*y*x - 4*y*y*y*x*x*z +1264*y*y*y*x*z*x + 3*y*y*y*y*z + 4*y*y*y*z*x*x + 6*y*y*y*z*y +127y*y*z*x*x*x*x + y*y*z*x*z + 7*y*y*z*y*x*x + 7*y*y*z*y*y -1287*y*y*z*z*x - y*z*x*x*x*z - y*z*x*x*z*x + 3*y*z*x*z*x*x +129y*z*x*z*y + y*z*y*x*x*x*x - 3*y*z*y*x*z + 7*y*z*y*y*x*x +1303*y*z*y*y*y - 3*y*z*y*z*x - 5*y*z*z*x*x*x - 4*y*z*z*y*x +1314*y*z*z*z - z*y*x*x*x*z - z*y*x*x*z*x - z*y*x*z*x*x -132z*y*x*z*y + z*y*y*x*x*x*x - 3*z*y*y*x*z + 3*z*y*y*y*x*x +133z*y*y*y*y - 3*z*y*y*z*x - z*y*z*x*x*x - 2*z*y*z*y*x +1342*z*y*z*z - z*z*x*x*x*x*x + 4*z*z*x*x*z + 4*z*z*x*z*x -1354*z*z*y*x*x*x - 3*z*z*y*y*x + 4*z*z*y*z + 4*z*z*z*x*x +1362*z*z*z*y,137x*x*x*x*x*z + x*x*x*x*z*x + x*x*x*z*x*x + x*x*z*x*x*x + x*z*x*x*x*x +138y*x*z*y - y*y*x*z + y*z*z + z*x*x*x*x*x - z*z*y,139x*x*x*x*x*x - y*x*z - y*y*y + z*z)140of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field141142Again, we can compute normal forms::143144sage: (z*I.0-I.1).normal_form(I)1450146sage: (z*I.0-x*y*z).normal_form(I)147-y*x*z + z*z148149"""150def __init__(self, ring, gens, coerce=True, side = "twosided"):151"""152INPUT:153154- ``ring``: A free algebra in letterplace implementation.155- ``gens``: List, tuple or sequence of generators.156- ``coerce`` (optional bool, default ``True``):157Shall ``gens`` be coerced first?158- ``side``: optional string, one of ``"twosided"`` (default),159``"left"`` or ``"right"``. Determines whether the ideal160is a left, right or twosided ideal. Groebner bases or161only supported in the twosided case.162163TEST::164165sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')166sage: from sage.algebras.letterplace.letterplace_ideal import LetterplaceIdeal167sage: LetterplaceIdeal(F,x)168Twosided Ideal (x) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field169sage: LetterplaceIdeal(F,[x,y],side='left')170Left Ideal (x, y) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field171172It is not correctly detected that this class inherits from an173extension class. Therefore, we have to skip one item of the174test suite::175176sage: I = F*[x*y+y*z,x^2+x*y-y*x-y^2]*F177sage: TestSuite(I).run(skip=['_test_category'],verbose=True)178running ._test_eq() . . . pass179running ._test_not_implemented_methods() . . . pass180running ._test_pickling() . . . pass181182"""183Ideal_nc.__init__(self, ring, gens, coerce=coerce, side=side)184self.__GB = self185self.__uptodeg = 0186def groebner_basis(self, degbound=None):187"""188Twosided Groebner basis with degree bound.189190INPUT:191192- ``degbound`` (optional integer, or Infinity): If it is provided,193a Groebner basis at least out to that degree is returned. By194default, the current degree bound of the underlying ring is used.195196ASSUMPTIONS:197198Currently, we can only compute Groebner bases for twosided199ideals, and the ring of coefficients must be a field. A200`TypeError` is raised if one of these conditions is violated.201202NOTES:203204- The result is cached. The same Groebner basis is returned205if a smaller degree bound than the known one is requested.206- If the degree bound Infinity is requested, it is attempted to207compute a complete Groebner basis. But we can not guarantee208that the computation will terminate, since not all twosided209homogeneous ideals of a free algebra have a finite Groebner210basis.211212EXAMPLES::213214sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')215sage: I = F*[x*y+y*z,x^2+x*y-y*x-y^2]*F216217Since `F` was cached and since its degree bound can not be218decreased, it may happen that, as a side effect of other tests,219it already has a degree bound bigger than 3. So, we can not220test against the output of ``I.groebner_basis()``::221222sage: F.set_degbound(3)223sage: I.groebner_basis() # not tested224Twosided Ideal (y*y*y - y*y*z + y*z*y - y*z*z, y*y*x + y*y*z + y*z*x + y*z*z, x*y + y*z, x*x - y*x - y*y - y*z) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field225sage: I.groebner_basis(4)226Twosided Ideal (y*z*y*y - y*z*y*z + y*z*z*y - y*z*z*z, y*z*y*x + y*z*y*z + y*z*z*x + y*z*z*z, y*y*z*y - y*y*z*z + y*z*z*y - y*z*z*z, y*y*z*x + y*y*z*z + y*z*z*x + y*z*z*z, y*y*y - y*y*z + y*z*y - y*z*z, y*y*x + y*y*z + y*z*x + y*z*z, x*y + y*z, x*x - y*x - y*y - y*z) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field227sage: I.groebner_basis(2) is I.groebner_basis(4)228True229sage: G = I.groebner_basis(4)230sage: G.groebner_basis(3) is G231True232233If a finite complete Groebner basis exists, we can compute234it as follows::235236sage: I = F*[x*y-y*x,x*z-z*x,y*z-z*y,x^2*y-z^3,x*y^2+z*x^2]*F237sage: I.groebner_basis(Infinity)238Twosided Ideal (z*z*z*y*y + z*z*z*z*x, z*x*x*x + z*z*z*y, y*z - z*y, y*y*x + z*x*x, y*x*x - z*z*z, x*z - z*x, x*y - y*x) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field239240Since the commutators of the generators are contained in the ideal,241we can verify the above result by a computation in a polynomial ring242in negative lexicographic order::243244sage: P.<c,b,a> = PolynomialRing(QQ,order='neglex')245sage: J = P*[a^2*b-c^3,a*b^2+c*a^2]246sage: J.groebner_basis()247[b*a^2 - c^3, b^2*a + c*a^2, c*a^3 + c^3*b, c^3*b^2 + c^4*a]248249Aparently, the results are compatible, by sending `a` to `x`, `b`250to `y` and `c` to `z`.251252"""253cdef FreeAlgebra_letterplace A = self.ring()254cdef FreeAlgebraElement_letterplace x255if degbound is None:256degbound = A.degbound()257if self.__uptodeg >= degbound:258return self.__GB259if not A.base().is_field():260raise TypeError, "Currently, we can only compute Groebner bases if the ring of coefficients is a field"261if self.side()!='twosided':262raise TypeError, "This ideal is not two-sided. We can only compute two-sided Groebner bases"263if degbound == Infinity:264while self.__uptodeg<Infinity:265test_bound = 2*max([x._poly.degree() for x in self.__GB.gens()])266self.groebner_basis(test_bound)267return self.__GB268# Set the options required by letterplace269from sage.libs.singular.option import LibSingularOptions270libsingular_options = LibSingularOptions()271bck = (libsingular_options['redTail'],libsingular_options['redSB'])272libsingular_options['redTail'] = True273libsingular_options['redSB'] = True274A.set_degbound(degbound)275P = A._current_ring276out = [FreeAlgebraElement_letterplace(A,X,check=False) for X in277singular_system("freegb",P.ideal([x._poly for x in self.__GB.gens()]),278degbound,A.__ngens, ring = P)]279libsingular_options['redTail'] = bck[0]280libsingular_options['redSB'] = bck[1]281self.__GB = A.ideal(out,side='twosided',coerce=False)282if degbound >= 2*max([x._poly.degree() for x in out]):283degbound = Infinity284self.__uptodeg = degbound285self.__GB.__uptodeg = degbound286return self.__GB287288def __contains__(self,x):289"""290The containment test is based on a normal form computation.291292EXAMPLES::293294sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')295sage: I = F*[x*y+y*z,x^2+x*y-y*x-y^2]*F296sage: x*I.0-I.1*y+I.0*y in I # indirect doctest297True298sage: 1 in I299False300301"""302R = self.ring()303return (x in R) and R(x).normal_form(self).is_zero()304305def reduce(self, G):306"""307Reduction of this ideal by another ideal,308or normal form of an algebra element with respect to this ideal.309310INPUT:311312- ``G``: A list or tuple of elements, an ideal,313the ambient algebra, or a single element.314315OUTPUT:316317- The normal form of ``G`` with respect to this ideal, if318``G`` is an element of the algebra.319- The reduction of this ideal by the elements resp. generators320of ``G``, if ``G`` is a list, tuple or ideal.321- The zero ideal, if ``G`` is the algebra containing this ideal.322323EXAMPLES::324325sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')326sage: I = F*[x*y+y*z,x^2+x*y-y*x-y^2]*F327sage: I.reduce(F)328Twosided Ideal (0) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field329sage: I.reduce(x^3)330-y*z*x - y*z*y - y*z*z331sage: I.reduce([x*y])332Twosided Ideal (y*z, x*x - y*x - y*y) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field333sage: I.reduce(F*[x^2+x*y,y^2+y*z]*F)334Twosided Ideal (x*y + y*z, -y*x + y*z) of Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field335336"""337P = self.ring()338if not isinstance(G,(list,tuple)):339if G==P:340return P.ideal([P.zero_element()])341if G in P:342return G.normal_form(self)343G = G.gens()344C = P.current_ring()345sI = C.ideal([C(X.letterplace_polynomial()) for X in self.gens()], coerce=False)346selfdeg = max([x.degree() for x in sI.gens()])347gI = P._reductor_(G, selfdeg)348from sage.libs.singular.option import LibSingularOptions349libsingular_options = LibSingularOptions()350bck = (libsingular_options['redTail'],libsingular_options['redSB'])351libsingular_options['redTail'] = True352libsingular_options['redSB'] = True353sI = poly_reduce(sI,gI, ring=C, attributes={gI:{"isSB":1}})354libsingular_options['redTail'] = bck[0]355libsingular_options['redSB'] = bck[1]356return P.ideal([FreeAlgebraElement_letterplace(P,x,check=False) for x in sI], coerce=False)357358359