X-Git-Url: http://gitweb.michael.orlitzky.com/?a=blobdiff_plain;f=mjo%2Feja%2Feja_algebra.py;h=992174e45b763e95677588beee2787859c0e64b1;hb=d4abf92e1e275554019be8987c6e837dfdc40150;hp=3855f0ef41e0e95efcadb359d6c7b908960a13b2;hpb=13a49fd59d57cf34b026460ace004ddbd60d4b68;p=sage.d.git diff --git a/mjo/eja/eja_algebra.py b/mjo/eja/eja_algebra.py index 3855f0e..992174e 100644 --- a/mjo/eja/eja_algebra.py +++ b/mjo/eja/eja_algebra.py @@ -9,10 +9,11 @@ from sage.algebras.quatalg.quaternion_algebra import QuaternionAlgebra from sage.categories.magmatic_algebras import MagmaticAlgebras from sage.combinat.free_module import CombinatorialFreeModule from sage.matrix.constructor import matrix +from sage.matrix.matrix_space import MatrixSpace from sage.misc.cachefunc import cached_method from sage.misc.prandom import choice from sage.misc.table import table -from sage.modules.free_module import VectorSpace +from sage.modules.free_module import FreeModule, VectorSpace from sage.rings.integer_ring import ZZ from sage.rings.number_field.number_field import QuadraticField from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing @@ -134,13 +135,17 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule): return self.zero() natural_basis = self.natural_basis() - if elt not in natural_basis[0].matrix_space(): + basis_space = natural_basis[0].matrix_space() + if elt not in basis_space: raise ValueError("not a naturally-represented algebra element") - # Thanks for nothing! Matrix spaces aren't vector - # spaces in Sage, so we have to figure out its - # natural-basis coordinates ourselves. - V = VectorSpace(elt.base_ring(), elt.nrows()*elt.ncols()) + # Thanks for nothing! Matrix spaces aren't vector spaces in + # Sage, so we have to figure out its natural-basis coordinates + # ourselves. We use the basis space's ring instead of the + # element's ring because the basis space might be an algebraic + # closure whereas the base ring of the 3-by-3 identity matrix + # could be QQ instead of QQbar. + V = VectorSpace(basis_space.base_ring(), elt.nrows()*elt.ncols()) W = V.span_of_basis( _mat2vec(s) for s in natural_basis ) coords = W.coordinate_vector(_mat2vec(elt)) return self.from_vector(coords) @@ -207,8 +212,13 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule): determinant). """ z = self._a_regular_element() - V = self.vector_space() - V1 = V.span_of_basis( (z**k).to_vector() for k in range(self.rank()) ) + # Don't use the parent vector space directly here in case this + # happens to be a subalgebra. In that case, we would be e.g. + # two-dimensional but span_of_basis() would expect three + # coordinates. + V = VectorSpace(self.base_ring(), self.vector_space().dimension()) + basis = [ (z**k).to_vector() for k in range(self.rank()) ] + V1 = V.span_of_basis( basis ) b = (V1.basis() + V1.complement().basis()) return V.span_of_basis(b) @@ -263,7 +273,13 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule): # have multivatiate polynomial entries. names = tuple('X' + str(i) for i in range(1,n+1)) R = PolynomialRing(self.base_ring(), names) - V = self.vector_space().change_ring(R) + + # Using change_ring() on the parent's vector space doesn't work + # here because, in a subalgebra, that vector space has a basis + # and change_ring() tries to bring the basis along with it. And + # that doesn't work unless the new ring is a PID, which it usually + # won't be. + V = FreeModule(R,n) # Now let x = (X1,X2,...,Xn) be the vector whose entries are # indeterminates... @@ -404,6 +420,29 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule): return x.trace_inner_product(y) + def is_trivial(self): + """ + Return whether or not this algebra is trivial. + + A trivial algebra contains only the zero element. + + SETUP:: + + sage: from mjo.eja.eja_algebra import ComplexHermitianEJA + + EXAMPLES:: + + sage: J = ComplexHermitianEJA(3) + sage: J.is_trivial() + False + sage: A = J.zero().subalgebra_generated_by() + sage: A.is_trivial() + True + + """ + return self.dimension() == 0 + + def multiplication_table(self): """ Return a visual representation of this algebra's multiplication @@ -463,8 +502,8 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule): Finite family {0: e0, 1: e1, 2: e2} sage: J.natural_basis() ( - [1 0] [0 1] [0 0] - [0 0], [1 0], [0 1] + [1 0] [ 0 1/2*sqrt2] [0 0] + [0 0], [1/2*sqrt2 0], [0 1] ) :: @@ -480,11 +519,23 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule): """ if self._natural_basis is None: - return tuple( b.to_vector().column() for b in self.basis() ) + M = self.natural_basis_space() + return tuple( M(b.to_vector()) for b in self.basis() ) else: return self._natural_basis + def natural_basis_space(self): + """ + Return the matrix space in which this algebra's natural basis + elements live. + """ + if self._natural_basis is None or len(self._natural_basis) == 0: + return MatrixSpace(self.base_ring(), self.dimension(), 1) + else: + return self._natural_basis[0].matrix_space() + + @cached_method def one(self): """ @@ -546,6 +597,15 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule): return self.linear_combination(zip(self.gens(), coeffs)) + def random_element(self): + # Temporary workaround for https://trac.sagemath.org/ticket/28327 + if self.is_trivial(): + return self.zero() + else: + s = super(FiniteDimensionalEuclideanJordanAlgebra, self) + return s.random_element() + + def rank(self): """ Return the rank of this EJA. @@ -618,7 +678,7 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule): sage: J = RealSymmetricEJA(2) sage: J.vector_space() - Vector space of dimension 3 over Rational Field + Vector space of dimension 3 over... """ return self.zero().to_vector().parent().ambient_vector_space() @@ -659,14 +719,32 @@ class RealCartesianProductEJA(FiniteDimensionalEuclideanJordanAlgebra): sage: e2*e2 e2 + TESTS: + + We can change the generator prefix:: + + sage: RealCartesianProductEJA(3, prefix='r').gens() + (r0, r1, r2) + + Our inner product satisfies the Jordan axiom:: + + sage: set_random_seed() + sage: n = ZZ.random_element(1,5) + sage: J = RealCartesianProductEJA(n) + sage: x = J.random_element() + sage: y = J.random_element() + sage: z = J.random_element() + sage: (x*y).inner_product(z) == y.inner_product(x*z) + True + """ - def __init__(self, n, field=QQ): + def __init__(self, n, field=QQ, **kwargs): V = VectorSpace(field, n) mult_table = [ [ V.gen(i)*(i == j) for j in range(n) ] for i in range(n) ] fdeja = super(RealCartesianProductEJA, self) - return fdeja.__init__(field, mult_table, rank=n) + return fdeja.__init__(field, mult_table, rank=n, **kwargs) def inner_product(self, x, y): return _usual_ip(x,y) @@ -723,9 +801,22 @@ def random_eja(): -def _real_symmetric_basis(n, field=QQ): +def _real_symmetric_basis(n, field): """ Return a basis for the space of real symmetric n-by-n matrices. + + SETUP:: + + sage: from mjo.eja.eja_algebra import _real_symmetric_basis + + TESTS:: + + sage: set_random_seed() + sage: n = ZZ.random_element(1,5) + sage: B = _real_symmetric_basis(n, QQbar) + sage: all( M.is_symmetric() for M in B) + True + """ # The basis of symmetric matrices, as matrices, in their R^(n-by-n) # coordinates. @@ -736,13 +827,14 @@ def _real_symmetric_basis(n, field=QQ): if i == j: Sij = Eij else: - # Beware, orthogonal but not normalized! Sij = Eij + Eij.transpose() + # Now normalize it. + Sij = Sij / _real_symmetric_matrix_ip(Sij,Sij).sqrt() S.append(Sij) return tuple(S) -def _complex_hermitian_basis(n, field=QQ): +def _complex_hermitian_basis(n, field): """ Returns a basis for the space of complex Hermitian n-by-n matrices. @@ -754,7 +846,8 @@ def _complex_hermitian_basis(n, field=QQ): sage: set_random_seed() sage: n = ZZ.random_element(1,5) - sage: all( M.is_symmetric() for M in _complex_hermitian_basis(n) ) + sage: B = _complex_hermitian_basis(n, QQ) + sage: all( M.is_symmetric() for M in B) True """ @@ -783,7 +876,7 @@ def _complex_hermitian_basis(n, field=QQ): return tuple(S) -def _quaternion_hermitian_basis(n, field=QQ): +def _quaternion_hermitian_basis(n, field): """ Returns a basis for the space of quaternion Hermitian n-by-n matrices. @@ -795,7 +888,8 @@ def _quaternion_hermitian_basis(n, field=QQ): sage: set_random_seed() sage: n = ZZ.random_element(1,5) - sage: all( M.is_symmetric() for M in _quaternion_hermitian_basis(n) ) + sage: B = _quaternion_hermitian_basis(n, QQbar) + sage: all( M.is_symmetric() for M in B ) True """ @@ -1097,6 +1191,9 @@ def _matrix_ip(X,Y): Y_mat = Y.natural_representation() return (X_mat*Y_mat).trace() +def _real_symmetric_matrix_ip(X,Y): + return (X*Y).trace() + class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra): """ @@ -1115,7 +1212,7 @@ class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra): sage: e0*e0 e0 sage: e1*e1 - e0 + e2 + 1/2*e0 + 1/2*e2 sage: e2*e2 e2 @@ -1145,19 +1242,61 @@ class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra): sage: J(expected) == x*y True + We can change the generator prefix:: + + sage: RealSymmetricEJA(3, prefix='q').gens() + (q0, q1, q2, q3, q4, q5) + + Our inner product satisfies the Jordan axiom:: + + sage: set_random_seed() + sage: n = ZZ.random_element(1,5) + sage: J = RealSymmetricEJA(n) + sage: x = J.random_element() + sage: y = J.random_element() + sage: z = J.random_element() + sage: (x*y).inner_product(z) == y.inner_product(x*z) + True + + Our basis is normalized with respect to the natural inner product:: + + sage: set_random_seed() + sage: n = ZZ.random_element(1,5) + sage: J = RealSymmetricEJA(n) + sage: all( b.norm() == 1 for b in J.gens() ) + True + + Left-multiplication operators are symmetric because they satisfy + the Jordan axiom:: + + sage: set_random_seed() + sage: n = ZZ.random_element(1,5) + sage: x = RealSymmetricEJA(n).random_element() + sage: x.operator().matrix().is_symmetric() + True + """ - def __init__(self, n, field=QQ): - S = _real_symmetric_basis(n, field=field) + def __init__(self, n, field=QQ, **kwargs): + if n > 1 and field is QQ: + # We'll need sqrt(2) to normalize the basis, and this + # winds up in the multiplication table, so the whole + # algebra needs to be over the field extension. + field = QuadraticField(2, 'sqrt2') + + S = _real_symmetric_basis(n, field) Qs = _multiplication_table_from_matrix_basis(S) fdeja = super(RealSymmetricEJA, self) return fdeja.__init__(field, Qs, rank=n, - natural_basis=S) + natural_basis=S, + **kwargs) def inner_product(self, x, y): - return _matrix_ip(x,y) + X = x.natural_representation() + Y = y.natural_representation() + return _real_symmetric_matrix_ip(X,Y) class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): @@ -1197,16 +1336,33 @@ class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): sage: J(expected) == x*y True + We can change the generator prefix:: + + sage: ComplexHermitianEJA(2, prefix='z').gens() + (z0, z1, z2, z3) + + Our inner product satisfies the Jordan axiom:: + + sage: set_random_seed() + sage: n = ZZ.random_element(1,5) + sage: J = ComplexHermitianEJA(n) + sage: x = J.random_element() + sage: y = J.random_element() + sage: z = J.random_element() + sage: (x*y).inner_product(z) == y.inner_product(x*z) + True + """ - def __init__(self, n, field=QQ): - S = _complex_hermitian_basis(n) + def __init__(self, n, field=QQ, **kwargs): + S = _complex_hermitian_basis(n, field) Qs = _multiplication_table_from_matrix_basis(S) fdeja = super(ComplexHermitianEJA, self) return fdeja.__init__(field, Qs, rank=n, - natural_basis=S) + natural_basis=S, + **kwargs) def inner_product(self, x, y): @@ -1257,16 +1413,33 @@ class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): sage: J(expected) == x*y True + We can change the generator prefix:: + + sage: QuaternionHermitianEJA(2, prefix='a').gens() + (a0, a1, a2, a3, a4, a5) + + Our inner product satisfies the Jordan axiom:: + + sage: set_random_seed() + sage: n = ZZ.random_element(1,5) + sage: J = QuaternionHermitianEJA(n) + sage: x = J.random_element() + sage: y = J.random_element() + sage: z = J.random_element() + sage: (x*y).inner_product(z) == y.inner_product(x*z) + True + """ - def __init__(self, n, field=QQ): - S = _quaternion_hermitian_basis(n) + def __init__(self, n, field=QQ, **kwargs): + S = _quaternion_hermitian_basis(n, field) Qs = _multiplication_table_from_matrix_basis(S) fdeja = super(QuaternionHermitianEJA, self) return fdeja.__init__(field, Qs, rank=n, - natural_basis=S) + natural_basis=S, + **kwargs) def inner_product(self, x, y): # Since a+bi+cj+dk on the diagonal is represented as @@ -1313,8 +1486,24 @@ class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra): sage: e2*e3 0 + We can change the generator prefix:: + + sage: JordanSpinEJA(2, prefix='B').gens() + (B0, B1) + + Our inner product satisfies the Jordan axiom:: + + sage: set_random_seed() + sage: n = ZZ.random_element(1,5) + sage: J = JordanSpinEJA(n) + sage: x = J.random_element() + sage: y = J.random_element() + sage: z = J.random_element() + sage: (x*y).inner_product(z) == y.inner_product(x*z) + True + """ - def __init__(self, n, field=QQ): + def __init__(self, n, field=QQ, **kwargs): V = VectorSpace(field, n) mult_table = [[V.zero() for j in range(n)] for i in range(n)] for i in range(n): @@ -1335,7 +1524,7 @@ class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra): # one-dimensional ambient space (because the rank is bounded by # the ambient dimension). fdeja = super(JordanSpinEJA, self) - return fdeja.__init__(field, mult_table, rank=min(n,2)) + return fdeja.__init__(field, mult_table, rank=min(n,2), **kwargs) def inner_product(self, x, y): return _usual_ip(x,y)