Path: blob/master/src/sage/algebras/letterplace/free_algebra_letterplace.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"""10Free associative unital algebras, implemented via Singular's letterplace rings1112AUTHOR:1314- Simon King (2011-03-21): :trac:`7797`1516With this implementation, Groebner bases out to a degree bound and17normal forms can be computed for twosided weighted homogeneous ideals18of free algebras. For now, all computations are restricted to weighted19homogeneous elements, i.e., other elements can not be created by20arithmetic operations.2122EXAMPLES::2324sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')25sage: F26Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field27sage: I = F*[x*y+y*z,x^2+x*y-y*x-y^2]*F28sage: I29Twosided 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 Field30sage: x*(x*I.0-I.1*y+I.0*y)-I.1*y*z31x*y*x*y + x*y*y*y - x*y*y*z + x*y*z*y + y*x*y*z + y*y*y*z32sage: x^2*I.0-x*I.1*y+x*I.0*y-I.1*y*z in I33True3435The preceding containment test is based on the computation of Groebner36bases with degree bound::3738sage: I.groebner_basis(degbound=4)39Twosided 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 Field4041When reducing an element by `I`, the original generators are chosen::4243sage: (y*z*y*y).reduce(I)44y*z*y*y4546However, there is a method for computing the normal form of an47element, which is the same as reduction by the Groebner basis out to48the degree of that element::4950sage: (y*z*y*y).normal_form(I)51y*z*y*z - y*z*z*y + y*z*z*z52sage: (y*z*y*y).reduce(I.groebner_basis(4))53y*z*y*z - y*z*z*y + y*z*z*z5455The default term order derives from the degree reverse lexicographic56order on the commutative version of the free algebra::5758sage: F.commutative_ring().term_order()59Degree reverse lexicographic term order6061A different term order can be chosen, and of course may yield a62different normal form::6364sage: L.<a,b,c> = FreeAlgebra(QQ, implementation='letterplace', order='lex')65sage: L.commutative_ring().term_order()66Lexicographic term order67sage: J = L*[a*b+b*c,a^2+a*b-b*c-c^2]*L68sage: J.groebner_basis(4)69Twosided Ideal (2*b*c*b - b*c*c + c*c*b, a*c*c - 2*b*c*a - 2*b*c*c - c*c*a, a*b + b*c, a*a - 2*b*c - c*c) of Free Associative Unital Algebra on 3 generators (a, b, c) over Rational Field70sage: (b*c*b*b).normal_form(J)711/2*b*c*c*b - 1/2*c*c*b*b7273Here is an example with degree weights::7475sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace', degrees=[1,2,3])76sage: (x*y+z).degree()7737879TEST::8081sage: TestSuite(F).run()82sage: TestSuite(L).run()83sage: loads(dumps(F)) is F84True8586TODO:8788The computation of Groebner bases only works for global term89orderings, and all elements must be weighted homogeneous with respect90to positive integral degree weights. It is ongoing work in Singular to91lift these restrictions.9293We support coercion from the letterplace wrapper to the corresponding94generic implementation of a free algebra95(:class:`~sage.algebras.free_algebra.FreeAlgebra_generic`), but there96is no coercion in the opposite direction, since the generic97implementation also comprises non-homogeneous elements.9899We also do not support coercion from a subalgebra, or between free100algebras with different term orderings, yet.101102"""103104from sage.all import PolynomialRing, prod105from sage.libs.singular.function import lib, singular_function106from sage.rings.polynomial.term_order import TermOrder107from sage.rings.polynomial.multi_polynomial_ring_generic import MPolynomialRing_generic108from sage.categories.algebras import Algebras109from sage.rings.noncommutative_ideals import IdealMonoid_nc110111#####################112# Define some singular functions113lib("freegb.lib")114poly_reduce = singular_function("NF")115singular_system=singular_function("system")116117# unfortunately we can not set Singular attributes for MPolynomialRing_libsingular118# Hence, we must constantly work around Letterplace's sanity checks,119# and can not use the following library functions:120#set_letterplace_attributes = singular_function("setLetterplaceAttributes")121#lpMult = singular_function("lpMult")122123#####################124# Auxiliar functions125126cdef MPolynomialRing_libsingular make_letterplace_ring(base_ring,blocks):127"""128Create a polynomial ring in block order.129130INPUT:131132- ``base_ring``: A multivariate polynomial ring.133- ``blocks``: The number of blocks to be formed.134135OUTPUT:136137A multivariate polynomial ring in block order, all blocks138isomorphic (as ordered rings) with the given ring, and the139variable names of the `n`-th block (`n>0`) ending with140``"_%d"%n``.141142TEST:143144Note that, since the algebras are cached, we need to choose145a different base ring, since other doctests could have a146side effect on the atteined degree bound::147148sage: F.<x,y,z> = FreeAlgebra(GF(17), implementation='letterplace')149sage: L.<a,b,c> = FreeAlgebra(GF(17), implementation='letterplace', order='lex')150sage: F.set_degbound(4)151sage: F.current_ring() # indirect doctest152Multivariate Polynomial Ring in x, y, z, x_1, y_1, z_1, x_2, y_2, z_2, x_3, y_3, z_3 over Finite Field of size 17153sage: F.current_ring().term_order()154Block term order with blocks:155(Degree reverse lexicographic term order of length 3,156Degree reverse lexicographic term order of length 3,157Degree reverse lexicographic term order of length 3,158Degree reverse lexicographic term order of length 3)159sage: L.set_degbound(2)160sage: L.current_ring().term_order()161Block term order with blocks:162(Lexicographic term order of length 3,163Lexicographic term order of length 3)164165"""166n = base_ring.ngens()167T0 = base_ring.term_order()168T = T0169cdef i170cdef tuple names0 = base_ring.variable_names()171cdef list names = list(names0)172for i from 1<=i<blocks:173T += T0174names.extend([x+'_'+str(i) for x in names0])175return PolynomialRing(base_ring.base_ring(),len(names),names,order=T)176177#####################178# The free algebra179180cdef class FreeAlgebra_letterplace(Algebra):181"""182Finitely generated free algebra, with arithmetic restricted to weighted homogeneous elements.183184NOTE:185186The restriction to weighted homogeneous elements should be lifted187as soon as the restriction to homogeneous elements is lifted in188Singular's "Letterplace algebras".189190EXAMPLE::191192sage: K.<z> = GF(25)193sage: F.<a,b,c> = FreeAlgebra(K, implementation='letterplace')194sage: F195Free Associative Unital Algebra on 3 generators (a, b, c) over Finite Field in z of size 5^2196sage: P = F.commutative_ring()197sage: P198Multivariate Polynomial Ring in a, b, c over Finite Field in z of size 5^2199200We can do arithmetic as usual, as long as we stay (weighted) homogeneous::201202sage: (z*a+(z+1)*b+2*c)^2203(z + 3)*a*a + (2*z + 3)*a*b + (2*z)*a*c + (2*z + 3)*b*a + (3*z + 4)*b*b + (2*z + 2)*b*c + (2*z)*c*a + (2*z + 2)*c*b - c*c204sage: a+1205Traceback (most recent call last):206...207ArithmeticError: Can only add elements of the same weighted degree208209"""210# It is not really a free algebra over the given generators. Rather,211# it is a free algebra over the commutative monoid generated by the given generators.212def __init__(self, R, degrees=None):213"""214INPUT:215216A multivariate polynomial ring of type :class:`~sage.rings.polynomial.multipolynomial_libsingular.MPolynomialRing_libsingular`.217218OUTPUT:219220The free associative version of the given commutative ring.221222NOTE:223224One is supposed to use the `FreeAlgebra` constructor, in order to use the cache.225226TEST::227228sage: from sage.algebras.letterplace.free_algebra_letterplace import FreeAlgebra_letterplace229sage: FreeAlgebra_letterplace(QQ['x','y'])230Free Associative Unital Algebra on 2 generators (x, y) over Rational Field231sage: FreeAlgebra_letterplace(QQ['x'])232Traceback (most recent call last):233...234TypeError: A letterplace algebra must be provided by a polynomial ring of type <type 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomialRing_libsingular'>235236::237238sage: K.<z> = GF(25)239sage: F.<a,b,c> = FreeAlgebra(K, implementation='letterplace')240sage: TestSuite(F).run(verbose=True)241running ._test_additive_associativity() . . . pass242running ._test_an_element() . . . pass243running ._test_associativity() . . . pass244running ._test_category() . . . pass245running ._test_characteristic() . . . pass246running ._test_distributivity() . . . pass247running ._test_elements() . . .248Running the test suite of self.an_element()249running ._test_category() . . . pass250running ._test_eq() . . . pass251running ._test_nonzero_equal() . . . pass252running ._test_not_implemented_methods() . . . pass253running ._test_pickling() . . . pass254pass255running ._test_elements_eq_reflexive() . . . pass256running ._test_elements_eq_symmetric() . . . pass257running ._test_elements_eq_transitive() . . . pass258running ._test_elements_neq() . . . pass259running ._test_eq() . . . pass260running ._test_not_implemented_methods() . . . pass261running ._test_one() . . . pass262running ._test_pickling() . . . pass263running ._test_prod() . . . pass264running ._test_some_elements() . . . pass265running ._test_zero() . . . pass266267"""268if not isinstance(R,MPolynomialRing_libsingular):269raise TypeError, "A letterplace algebra must be provided by a polynomial ring of type %s"%MPolynomialRing_libsingular270self.__ngens = R.ngens()271if degrees is None:272varnames = R.variable_names()273self._nb_slackvars = 0274else:275varnames = R.variable_names()[:-1]276self._nb_slackvars = 1277base_ring = R.base_ring()278Algebra.__init__(self, base_ring, varnames,279normalize=False, category=Algebras(base_ring))280self._commutative_ring = R281self._current_ring = make_letterplace_ring(R,1)282self._degbound = 1283if degrees is None:284self._degrees = tuple([int(1)]*self.__ngens)285else:286if (not isinstance(degrees,(tuple,list))) or len(degrees)!=self.__ngens-1 or any([i<=0 for i in degrees]):287raise TypeError, "The generator degrees must be given by a list or tuple of %d positive integers"%(self.__ngens-1)288self._degrees = tuple([int(i) for i in degrees])289self.set_degbound(max(self._degrees))290self._populate_coercion_lists_(coerce_list=[base_ring])291def __reduce__(self):292"""293TEST::294295sage: K.<z> = GF(25)296sage: F.<a,b,c> = FreeAlgebra(K, implementation='letterplace')297sage: loads(dumps(F)) is F # indirect doctest298True299300"""301from sage.algebras.free_algebra import FreeAlgebra302if self._nb_slackvars==0:303return FreeAlgebra,(self._commutative_ring,)304return FreeAlgebra,(self._commutative_ring,None,None,None,None,None,None,None,self._degrees)305# Small methods306def ngens(self):307"""308Return the number of generators.309310EXAMPLE::311312sage: F.<a,b,c> = FreeAlgebra(QQ, implementation='letterplace')313sage: F.ngens()3143315316"""317return self.__ngens-self._nb_slackvars318def gen(self,i):319"""320Return the `i`-th generator.321322INPUT:323324`i` -- an integer.325326OUTPUT:327328Generator number `i`.329330EXAMPLE::331332sage: F.<a,b,c> = FreeAlgebra(QQ, implementation='letterplace')333sage: F.1 is F.1 # indirect doctest334True335sage: F.gen(2)336c337338"""339if i>=self.__ngens-self._nb_slackvars:340raise ValueError, "This free algebra only has %d generators"%(self.__ngens-self._nb_slackvars)341if self._gens is not None:342return self._gens[i]343deg = self._degrees[i]344#self.set_degbound(deg)345p = self._current_ring.gen(i)346cdef int n347cdef int j = self.__ngens-1348for n from 1<=n<deg:349j += self.__ngens350p *= self._current_ring.gen(j)351return FreeAlgebraElement_letterplace(self, p)352def current_ring(self):353"""354Return the commutative ring that is used to emulate355the non-commutative multiplication out to the current degree.356357EXAMPLE::358359sage: F.<a,b,c> = FreeAlgebra(QQ, implementation='letterplace')360sage: F.current_ring()361Multivariate Polynomial Ring in a, b, c over Rational Field362sage: a*b363a*b364sage: F.current_ring()365Multivariate Polynomial Ring in a, b, c, a_1, b_1, c_1 over Rational Field366sage: F.set_degbound(3)367sage: F.current_ring()368Multivariate Polynomial Ring in a, b, c, a_1, b_1, c_1, a_2, b_2, c_2 over Rational Field369370"""371return self._current_ring372def commutative_ring(self):373"""374Return the commutative version of this free algebra.375376NOTE:377378This commutative ring is used as a unique key of the free algebra.379380EXAMPLE::381382sage: K.<z> = GF(25)383sage: F.<a,b,c> = FreeAlgebra(K, implementation='letterplace')384sage: F385Free Associative Unital Algebra on 3 generators (a, b, c) over Finite Field in z of size 5^2386sage: F.commutative_ring()387Multivariate Polynomial Ring in a, b, c over Finite Field in z of size 5^2388sage: FreeAlgebra(F.commutative_ring()) is F389True390391"""392return self._commutative_ring393def term_order_of_block(self):394"""395Return the term order that is used for the commutative version of this free algebra.396397EXAMPLE::398399sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')400sage: F.term_order_of_block()401Degree reverse lexicographic term order402sage: L.<a,b,c> = FreeAlgebra(QQ, implementation='letterplace',order='lex')403sage: L.term_order_of_block()404Lexicographic term order405406"""407return self._commutative_ring.term_order()408409def generator_degrees(self):410return self._degrees411412# Some basic properties of this ring413def is_commutative(self):414"""415Tell whether this algebra is commutative, i.e., whether the generator number is one.416417EXAMPLE::418419sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')420sage: F.is_commutative()421False422sage: FreeAlgebra(QQ, implementation='letterplace', names=['x']).is_commutative()423True424425"""426return self.__ngens-self._nb_slackvars <= 1427428def is_field(self):429"""430Tell whether this free algebra is a field.431432NOTE:433434This would only be the case in the degenerate case of no generators.435But such an example can not be constructed in this implementation.436437TEST::438439sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')440sage: F.is_field()441False442443"""444return (not (self.__ngens-self._nb_slackvars)) and self._base.is_field()445446def _repr_(self):447"""448EXAMPLE::449450sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')451sage: F # indirect doctest452Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field453454The degree weights are not part of the string representation::455456sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace', degrees=[2,1,3])457sage: F458Free Associative Unital Algebra on 3 generators (x, y, z) over Rational Field459460461"""462return "Free Associative Unital Algebra on %d generators %s over %s"%(self.__ngens-self._nb_slackvars,self.gens(),self._base)463464def _latex_(self):465"""466Representation of this free algebra in LaTeX.467468EXAMPLE::469470sage: F.<bla,alpha,z> = FreeAlgebra(QQ, implementation='letterplace', degrees=[1,2,3])471sage: latex(F)472\Bold{Q}\langle \mbox{bla}, \alpha, z\rangle473474"""475from sage.all import latex476return "%s\\langle %s\\rangle"%(latex(self.base_ring()),', '.join(self.latex_variable_names()))477478def degbound(self):479"""480Return the degree bound that is currently used.481482NOTE:483484When multiplying two elements of this free algebra, the degree485bound will be dynamically adapted. It can also be set by486:meth:`set_degbound`.487488EXAMPLE:489490In order to avoid we get a free algebras from the cache that491was created in another doctest and has a different degree492bound, we choose a base ring that does not appear in other tests::493494sage: F.<x,y,z> = FreeAlgebra(ZZ, implementation='letterplace')495sage: F.degbound()4961497sage: x*y498x*y499sage: F.degbound()5002501sage: F.set_degbound(4)502sage: F.degbound()5034504505"""506return self._degbound507def set_degbound(self,d):508"""509Increase the degree bound that is currently in place.510511NOTE:512513The degree bound can not be decreased.514515EXAMPLE:516517In order to avoid we get a free algebras from the cache that518was created in another doctest and has a different degree519bound, we choose a base ring that does not appear in other tests::520521sage: F.<x,y,z> = FreeAlgebra(GF(251), implementation='letterplace')522sage: F.degbound()5231524sage: x*y525x*y526sage: F.degbound()5272528sage: F.set_degbound(4)529sage: F.degbound()5304531sage: F.set_degbound(2)532sage: F.degbound()5334534535"""536if d<=self._degbound:537return538self._degbound = d539self._current_ring = make_letterplace_ring(self._commutative_ring,d)540541# def base_extend(self, R):542# if self._base.has_coerce_map_from(R):543# return self544545################################################546## Ideals547548def _ideal_class_(self, n=0):549"""550Return the class :class:`~sage.algebras.letterplace.letterplace_ideal.LetterplaceIdeal`.551552EXAMPLE::553554sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')555sage: I = [x*y+y*z,x^2+x*y-y*x-y^2]*F556sage: I557Right 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 Field558sage: type(I) is F._ideal_class_()559True560561"""562from sage.algebras.letterplace.letterplace_ideal import LetterplaceIdeal563return LetterplaceIdeal564565def ideal_monoid(self):566"""567Return the monoid of ideals of this free algebra.568569EXAMPLE::570571sage: F.<x,y> = FreeAlgebra(GF(2), implementation='letterplace')572sage: F.ideal_monoid()573Monoid of ideals of Free Associative Unital Algebra on 2 generators (x, y) over Finite Field of size 2574sage: F.ideal_monoid() is F.ideal_monoid()575True576577"""578if self.__monoid is None:579self.__monoid = IdealMonoid_nc(self)580return self.__monoid581582# Auxiliar methods583cdef str exponents_to_string(self, E):584"""585This auxiliary method is used for the string representation of elements of this free algebra.586587EXAMPLE::588589sage: F.<x,y,z> = FreeAlgebra(GF(2), implementation='letterplace')590sage: x*y*x*z # indirect doctest591x*y*x*z592593It should be possible to use the letterplace algebra to implement the594free algebra generated by the elements of a finitely generated free abelian595monoid. However, we can not use it, yet. So, for now, we raise an error::596597sage: from sage.algebras.letterplace.free_algebra_element_letterplace import FreeAlgebraElement_letterplace598sage: P = F.commutative_ring()599sage: FreeAlgebraElement_letterplace(F, P.0*P.1^2+P.1^3) # indirect doctest600Traceback (most recent call last):601...602NotImplementedError:603Apparently you tried to view the letterplace algebra with604shift-multiplication as the free algebra over a finitely605generated free abelian monoid.606In principle, this is correct, but it is not implemented, yet.607608"""609cdef int ngens = self.__ngens610cdef int nblocks = len(E)/ngens611cdef int i,j,base, exp, var_ind612cdef list out = []613cdef list tmp614for i from 0<=i<nblocks:615base = i*ngens616tmp = [(j,E[base+j]) for j in xrange(ngens) if E[base+j]]617if not tmp:618continue619var_ind, exp = tmp[0]620if len(tmp)>1 or exp>1:621raise NotImplementedError, "\n Apparently you tried to view the letterplace algebra with\n shift-multiplication as the free algebra over a finitely\n generated free abelian monoid.\n In principle, this is correct, but it is not implemented, yet."622623out.append(self._names[var_ind])624i += (self._degrees[var_ind]-1)625### This was the original implementation, with "monoid hack" but without generator degrees626#s = '.'.join([('%s^%d'%(x,e) if e>1 else x) for x,e in zip(self._names,E[i*ngens:(i+1)*ngens]) if e])627#if s:628# out.append(s)629return '*'.join(out)630631# Auxiliar methods632cdef str exponents_to_latex(self, E):633"""634This auxiliary method is used for the representation of elements of this free algebra as a latex string.635636EXAMPLE::637638sage: K.<z> = GF(25)639sage: F.<a,b,c> = FreeAlgebra(K, implementation='letterplace', degrees=[1,2,3])640sage: -(a*b*(z+1)-c)^2641(2*z + 1)*a*b*a*b + (z + 1)*a*b*c + (z + 1)*c*a*b - c*c642sage: latex(-(a*b*(z+1)-c)^2) # indirect doctest643\left(2 z + 1\right) a b a b + \left(z + 1\right) a b c + \left(z + 1\right) c a b - c c644645"""646cdef int ngens = self.__ngens647cdef int nblocks = len(E)/ngens648cdef int i,j,base, exp, var_ind649cdef list out = []650cdef list tmp651cdef list names = self.latex_variable_names()652for i from 0<=i<nblocks:653base = i*ngens654tmp = [(j,E[base+j]) for j in xrange(ngens) if E[base+j]]655if not tmp:656continue657var_ind, exp = tmp[0]658if len(tmp)>1 or exp>1:659raise NotImplementedError, "\n Apparently you tried to view the letterplace algebra with\n shift-multiplication as the free algebra over a finitely\n generated free abelian monoid.\n In principle, this is correct, but it is not implemented, yet."660661out.append(names[var_ind])662i += (self._degrees[var_ind]-1)663return ' '.join(out)664665def _reductor_(self, g, d):666"""667Return a commutative ideal that can be used to compute the normal668form of a free algebra element of a given degree.669670INPUT:671672``g`` - a list of elements of this free algebra.673``d`` - an integer.674675OUTPUT:676677An ideal such that reduction of a letterplace polynomial by that ideal corresponds678to reduction of an element of degree at most ``d`` by ``g``.679680EXAMPLE::681682sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')683sage: I = F*[x*y+y*z,x^2+x*y-y*x-y^2]*F684sage: p = y*x*y + y*y*y + y*z*y - y*z*z685sage: p.reduce(I)686y*y*y - y*y*z + y*z*y - y*z*z687sage: G = F._reductor_(I.gens(),3); G688Ideal (x*y_1 + y*z_1, x_1*y_2 + y_1*z_2, x*x_1 + x*y_1 - y*x_1 - y*y_1, x_1*x_2 + x_1*y_2 - y_1*x_2 - y_1*y_2) of Multivariate Polynomial Ring in x, y, z, x_1, y_1, z_1, x_2, y_2, z_2... over Rational Field689690We do not use the usual reduction method for polynomials in691Sage, since it does the reductions in a different order692compared to Singular. Therefore, we call the original Singular693reduction method, and prevent a warning message by asserting694that `G` is a Groebner basis.695696sage: from sage.libs.singular.function import singular_function697sage: poly_reduce = singular_function("NF")698sage: q = poly_reduce(p.letterplace_polynomial(), G, ring=F.current_ring(), attributes={G:{"isSB":1}}); q699y*y_1*y_2 - y*y_1*z_2 + y*z_1*y_2 - y*z_1*z_2700sage: p.reduce(I).letterplace_polynomial() == q701True702703"""704cdef list out = []705C = self.current_ring()706cdef FreeAlgebraElement_letterplace x707ngens = self.__ngens708degbound = self._degbound709cdef list G = [C(x._poly) for x in g]710for y in G:711out.extend([y]+[singular_system("stest",y,n+1,degbound,ngens,ring=C) for n in xrange(d-y.degree())])712return C.ideal(out)713714###########################715## Coercion716cpdef _coerce_map_from_(self,S):717"""718A ring ``R`` coerces into self, if719720- it coerces into the current polynomial ring, or721- it is a free graded algebra in letterplace implementation,722the generator names of ``R`` are a proper subset of the723generator names of self, the degrees of equally named724generators are equal, and the base ring of ``R`` coerces725into the base ring of self.726727TEST:728729Coercion from the base ring::730731sage: F.<x,y,z> = FreeAlgebra(GF(5), implementation='letterplace')732sage: 5 == F.zero() # indirect doctest733True734735Coercion from another free graded algebra::736737sage: F.<t,y,z> = FreeAlgebra(ZZ, implementation='letterplace', degrees=[4,2,3])738sage: G = FreeAlgebra(GF(5), implementation='letterplace', names=['x','y','z','t'], degrees=[1,2,3,4])739sage: t*G.0 # indirect doctest740t*x741742"""743if self==S or self._current_ring.has_coerce_map_from(S):744return True745cdef int i746# Do we have another letterplace algebra?747if not isinstance(S, FreeAlgebra_letterplace):748return False749# Do the base rings coerce?750if not self.base_ring().has_coerce_map_from(S.base_ring()):751return False752# Do the names match?753cdef tuple degs, Sdegs, names, Snames754names = self.variable_names()755Snames = S.variable_names()756if not set(names).issuperset(Snames):757return False758# Do the degrees match759degs = self._degrees760Sdegs = (<FreeAlgebra_letterplace>S)._degrees761for i from 0<=i<S.ngens():762if degs[names.index(Snames[i])] != Sdegs[i]:763return False764return True765766def _an_element_(self):767"""768Return an element.769770EXAMPLE::771772sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')773sage: F.an_element() # indirect doctest774x775776"""777return FreeAlgebraElement_letterplace(self, self._current_ring.an_element(), check=False)778779# def random_element(self, degree=2, terms=5):780# """781# Return a random element of a given degree and with a given number of terms.782#783# INPUT:784#785# - ``degree`` -- the maximal degree of the output (default 2).786# - ``terms`` -- the maximal number of terms of the output (default 5).787#788# NOTE:789#790# This method is currently not useful at all.791#792# Not tested.793# """794# self.set_degbound(degree)795# while(1):796# p = self._current_ring.random_element(degree=degree,terms=terms)797# if p.is_homogeneous():798# break799# return FreeAlgebraElement_letterplace(self, p, check=False)800801def _from_dict_(self, D, check=True):802"""803Create an element from a dictionary.804805INPUT:806807- A dictionary. Keys: tuples of exponents. Values:808The coefficients of the corresponding monomial809in the to-be-created element.810- ``check`` (optional bool, default ``True``):811This is forwarded to the initialisation of812:class:`~sage.algebas.letterplace.free_algebra_element_letterplace.FreeAlgebraElement_letterplace`.813814TEST:815816This method applied to the dictionary of any element must817return the same element. This must hold true even if the818underlying letterplace ring has been extended in the meantime.819::820821sage: F.<x,y,z> = FreeAlgebra(QQ, implementation='letterplace')822sage: p = 3*x*y+2*z^2823sage: F.set_degbound(10)824sage: p == F._from_dict_(dict(p))825True826827For the empty dictionary, zero is returned::828829sage: F._from_dict_({})8300831832"""833if not D:834return self.zero_element()835cdef int l836for e in D.iterkeys():837l = len(e)838break839cdef dict out = {}840self.set_degbound(l/self.__ngens)841cdef int n = self._current_ring.ngens()842for e,c in D.iteritems():843out[tuple(e)+(0,)*(n-l)] = c844return FreeAlgebraElement_letterplace(self,self._current_ring(out),845check=check)846847def _element_constructor_(self, x):848"""849Return an element of this free algebra.850851INPUT:852853An element of a free algebra with a proper subset of generator854names, or anything that can be interpreted in the polynomial855ring that is used to implement the letterplace algebra out to856the current degree bound, or a string that can be interpreted857as an expression in the algebra (provided that the858coefficients are numerical).859860EXAMPLE::861862sage: F.<t,y,z> = FreeAlgebra(ZZ, implementation='letterplace', degrees=[4,2,3])863864Conversion of a number::865866sage: F(3)8673868869Interpretation of a string as an algebra element::870871sage: F('t*y+3*z^2')872t*y + 3*z*z873874Conversion from the currently underlying polynomial ring::875876sage: F.set_degbound(3)877sage: P = F.current_ring()878sage: F(P.0*P.7*P.11*P.15*P.17*P.23 - 2*P.2*P.7*P.11*P.14*P.19*P.23)879t*y - 2*z*z880881Conversion from a graded sub-algebra::882883sage: G = FreeAlgebra(GF(5), implementation='letterplace', names=['x','y','z','t'], degrees=[1,2,3,4])884sage: G(t*y + 2*y^3 - 4*z^2) # indirect doctest885(2)*y*y*y + z*z + t*y886887"""888if isinstance(x, basestring):889from sage.all import sage_eval890return sage_eval(x,locals=self.gens_dict())891try:892P = x.parent()893except AttributeError:894P = None895if P is self:896(<FreeAlgebraElement_letterplace>x)._poly = self._current_ring((<FreeAlgebraElement_letterplace>x)._poly)897return x898if isinstance(P, FreeAlgebra_letterplace):899self.set_degbound(P.degbound())900Ppoly = (<FreeAlgebra_letterplace>P)._current_ring901Gens = self._current_ring.gens()902Names = self._current_ring.variable_names()903PNames = list(Ppoly.variable_names())904# translate the slack variables905PNames[P.ngens(): len(PNames): P.ngens()+1] = list(Names[self.ngens(): len(Names): self.ngens()+1])[:P.degbound()]906x = Ppoly.hom([Gens[Names.index(asdf)] for asdf in PNames])(x.letterplace_polynomial())907return FreeAlgebraElement_letterplace(self,self._current_ring(x))908909910