Lie Algebras With Basis¶
AUTHORS:
- Travis Scrimshaw (07-15-2013): Initial implementation 
- class sage.categories.lie_algebras_with_basis.LieAlgebrasWithBasis(base_category)[source]¶
- Bases: - CategoryWithAxiom_over_base_ring- Category of Lie algebras with a basis. - class ElementMethods[source]¶
- Bases: - object- lift()[source]¶
- Lift - selfto the universal enveloping algebra.- EXAMPLES: - sage: # needs sage.combinat sage.groups sage: S = SymmetricGroup(3).algebra(QQ) sage: L = LieAlgebra(associative=S) sage: x = L.gen(3) sage: y = L.gen(1) sage: x.lift() b3 sage: y.lift() b1 sage: x * y b1*b3 + b4 - b5 - >>> from sage.all import * >>> # needs sage.combinat sage.groups >>> S = SymmetricGroup(Integer(3)).algebra(QQ) >>> L = LieAlgebra(associative=S) >>> x = L.gen(Integer(3)) >>> y = L.gen(Integer(1)) >>> x.lift() b3 >>> y.lift() b1 >>> x * y b1*b3 + b4 - b5 
 - to_vector(order=None)[source]¶
- Return the vector in - g.module()corresponding to the element- selfof- g(where- gis the parent of- self).- Implement this if you implement - g.module(). See- sage.categories.lie_algebras.LieAlgebras.module()for how this is to be done.- EXAMPLES: - sage: L = LieAlgebras(QQ).FiniteDimensional().WithBasis().example() # needs sage.modules sage: L.an_element().to_vector() # needs sage.modules (0, 0, 0) - >>> from sage.all import * >>> L = LieAlgebras(QQ).FiniteDimensional().WithBasis().example() # needs sage.modules >>> L.an_element().to_vector() # needs sage.modules (0, 0, 0) - Todo - Doctest this implementation on an example not overshadowed. 
 
 - Graded[source]¶
- alias of - GradedLieAlgebrasWithBasis
 - class ParentMethods[source]¶
- Bases: - object- bracket_on_basis(x, y)[source]¶
- Return the bracket of basis elements indexed by - xand- ywhere- x < y. If this is not implemented, then the method- _bracket_()for the elements must be overwritten.- EXAMPLES: - sage: L = LieAlgebras(QQ).WithBasis().example() # needs sage.combinat sage.modules sage: L.bracket_on_basis(Partition([3,1]), Partition([2,2,1,1])) # needs sage.combinat sage.modules 0 - >>> from sage.all import * >>> L = LieAlgebras(QQ).WithBasis().example() # needs sage.combinat sage.modules >>> L.bracket_on_basis(Partition([Integer(3),Integer(1)]), Partition([Integer(2),Integer(2),Integer(1),Integer(1)])) # needs sage.combinat sage.modules 0 
 - dimension()[source]¶
- Return the dimension of - self.- EXAMPLES: - sage: L = LieAlgebras(QQ).FiniteDimensional().WithBasis().example() # needs sage.modules sage: L.dimension() # needs sage.modules 3 - >>> from sage.all import * >>> L = LieAlgebras(QQ).FiniteDimensional().WithBasis().example() # needs sage.modules >>> L.dimension() # needs sage.modules 3 - sage: L = LieAlgebra(QQ, 'x,y', {('x','y'): {'x':1}}) # needs sage.combinat sage.modules sage: L.dimension() # needs sage.combinat sage.modules 2 - >>> from sage.all import * >>> L = LieAlgebra(QQ, 'x,y', {('x','y'): {'x':Integer(1)}}) # needs sage.combinat sage.modules >>> L.dimension() # needs sage.combinat sage.modules 2 
 - from_vector(v, order=None, coerce=False)[source]¶
- Return the element of - selfcorresponding to the vector- vin- self.module().- Implement this if you implement - module(); see the documentation of- sage.categories.lie_algebras.LieAlgebras.module()for how this is to be done.- EXAMPLES: - sage: L = LieAlgebras(QQ).FiniteDimensional().WithBasis().example() # needs sage.modules sage: u = L.from_vector(vector(QQ, (1, 0, 0))); u # needs sage.modules (1, 0, 0) sage: parent(u) is L # needs sage.modules True - >>> from sage.all import * >>> L = LieAlgebras(QQ).FiniteDimensional().WithBasis().example() # needs sage.modules >>> u = L.from_vector(vector(QQ, (Integer(1), Integer(0), Integer(0)))); u # needs sage.modules (1, 0, 0) >>> parent(u) is L # needs sage.modules True 
 - module()[source]¶
- Return an \(R\)-module which is isomorphic to the underlying \(R\)-module of - self.- See - sage.categories.lie_algebras.LieAlgebras.module()for an explanation.- EXAMPLES: - sage: L = LieAlgebras(QQ).WithBasis().example() # needs sage.combinat sage.modules sage: L.module() # needs sage.combinat sage.modules Free module generated by Partitions over Rational Field - >>> from sage.all import * >>> L = LieAlgebras(QQ).WithBasis().example() # needs sage.combinat sage.modules >>> L.module() # needs sage.combinat sage.modules Free module generated by Partitions over Rational Field 
 - pbw_basis(basis_key=None, **kwds)[source]¶
- Return the Poincare-Birkhoff-Witt basis of the universal enveloping algebra corresponding to - self.- EXAMPLES: - sage: L = lie_algebras.sl(QQ, 2) # needs sage.combinat sage.modules sage: PBW = L.pbw_basis() # needs sage.combinat sage.modules - >>> from sage.all import * >>> L = lie_algebras.sl(QQ, Integer(2)) # needs sage.combinat sage.modules >>> PBW = L.pbw_basis() # needs sage.combinat sage.modules 
 - poincare_birkhoff_witt_basis(basis_key=None, **kwds)[source]¶
- Return the Poincare-Birkhoff-Witt basis of the universal enveloping algebra corresponding to - self.- EXAMPLES: - sage: L = lie_algebras.sl(QQ, 2) # needs sage.combinat sage.modules sage: PBW = L.pbw_basis() # needs sage.combinat sage.modules - >>> from sage.all import * >>> L = lie_algebras.sl(QQ, Integer(2)) # needs sage.combinat sage.modules >>> PBW = L.pbw_basis() # needs sage.combinat sage.modules 
 
 - example(gens=None)[source]¶
- Return an example of a Lie algebra as per - Category.example.- EXAMPLES: - sage: LieAlgebras(QQ).WithBasis().example() # needs sage.combinat sage.modules An example of a Lie algebra: the abelian Lie algebra on the generators indexed by Partitions over Rational Field - >>> from sage.all import * >>> LieAlgebras(QQ).WithBasis().example() # needs sage.combinat sage.modules An example of a Lie algebra: the abelian Lie algebra on the generators indexed by Partitions over Rational Field - Another set of generators can be specified as an optional argument: - sage: LieAlgebras(QQ).WithBasis().example(Compositions()) # needs sage.combinat sage.modules An example of a Lie algebra: the abelian Lie algebra on the generators indexed by Compositions of nonnegative integers over Rational Field - >>> from sage.all import * >>> LieAlgebras(QQ).WithBasis().example(Compositions()) # needs sage.combinat sage.modules An example of a Lie algebra: the abelian Lie algebra on the generators indexed by Compositions of nonnegative integers over Rational Field