]> gitweb.michael.orlitzky.com - sage.d.git/blobdiff - mjo/eja/euclidean_jordan_algebra.py
eja: fix alphabetical ordering of element methods.
[sage.d.git] / mjo / eja / euclidean_jordan_algebra.py
index 4713ff0859876b5832e6e5e9551f632444c6a138..8d9b27e974fff75f769579853d7ac60716d26298 100644 (file)
@@ -20,7 +20,8 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
                               names='e',
                               assume_associative=False,
                               category=None,
-                              rank=None):
+                              rank=None,
+                              natural_basis=None):
         n = len(mult_table)
         mult_table = [b.base_extend(field) for b in mult_table]
         for b in mult_table:
@@ -43,7 +44,8 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
                                  assume_associative=assume_associative,
                                  names=names,
                                  category=cat,
-                                 rank=rank)
+                                 rank=rank,
+                                 natural_basis=natural_basis)
 
 
     def __init__(self, field,
@@ -51,7 +53,8 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
                  names='e',
                  assume_associative=False,
                  category=None,
-                 rank=None):
+                 rank=None,
+                 natural_basis=None):
         """
         EXAMPLES:
 
@@ -66,6 +69,7 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
 
         """
         self._rank = rank
+        self._natural_basis = natural_basis
         fda = super(FiniteDimensionalEuclideanJordanAlgebra, self)
         fda.__init__(field,
                      mult_table,
@@ -80,6 +84,49 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
         fmt = "Euclidean Jordan algebra of degree {} over {}"
         return fmt.format(self.degree(), self.base_ring())
 
+
+    def natural_basis(self):
+        """
+        Return a more-natural representation of this algebra's basis.
+
+        Every finite-dimensional Euclidean Jordan Algebra is a direct
+        sum of five simple algebras, four of which comprise Hermitian
+        matrices. This method returns the original "natural" basis
+        for our underlying vector space. (Typically, the natural basis
+        is used to construct the multiplication table in the first place.)
+
+        Note that this will always return a matrix. The standard basis
+        in `R^n` will be returned as `n`-by-`1` column matrices.
+
+        EXAMPLES::
+
+            sage: J = RealSymmetricSimpleEJA(2)
+            sage: J.basis()
+            Family (e0, e1, e2)
+            sage: J.natural_basis()
+            (
+            [1 0]  [0 1]  [0 0]
+            [0 0], [1 0], [0 1]
+            )
+
+        ::
+
+            sage: J = JordanSpinSimpleEJA(2)
+            sage: J.basis()
+            Family (e0, e1)
+            sage: J.natural_basis()
+            (
+            [1]  [0]
+            [0], [1]
+            )
+
+        """
+        if self._natural_basis is None:
+            return tuple( b.vector().column() for b in self.basis() )
+        else:
+            return self._natural_basis
+
+
     def rank(self):
         """
         Return the rank of this EJA.
@@ -112,7 +159,7 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
 
                 sage: set_random_seed()
                 sage: x = random_eja().random_element()
-                sage: x.matrix()*x.vector() == (x^2).vector()
+                sage: x.operator_matrix()*x.vector() == (x^2).vector()
                 True
 
             A few examples of power-associativity::
@@ -131,8 +178,8 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
                 sage: x = random_eja().random_element()
                 sage: m = ZZ.random_element(0,10)
                 sage: n = ZZ.random_element(0,10)
-                sage: Lxm = (x^m).matrix()
-                sage: Lxn = (x^n).matrix()
+                sage: Lxm = (x^m).operator_matrix()
+                sage: Lxn = (x^n).operator_matrix()
                 sage: Lxm*Lxn == Lxn*Lxm
                 True
 
@@ -143,7 +190,7 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
             elif n == 1:
                 return self
             else:
-                return A.element_class(A, (self.matrix()**(n-1))*self.vector())
+                return A( (self.operator_matrix()**(n-1))*self.vector() )
 
 
         def characteristic_polynomial(self):
@@ -193,8 +240,8 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
             if not other in self.parent():
                 raise ArgumentError("'other' must live in the same algebra")
 
-            A = self.matrix()
-            B = other.matrix()
+            A = self.operator_matrix()
+            B = other.operator_matrix()
             return (A*B == B*A)
 
 
@@ -414,73 +461,6 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
             return self.span_of_powers().dimension()
 
 
-        def matrix(self):
-            """
-            Return the matrix that represents left- (or right-)
-            multiplication by this element in the parent algebra.
-
-            We have to override this because the superclass method
-            returns a matrix that acts on row vectors (that is, on
-            the right).
-
-            EXAMPLES:
-
-            Test the first polarization identity from my notes, Koecher Chapter
-            III, or from Baes (2.3)::
-
-                sage: set_random_seed()
-                sage: J = random_eja()
-                sage: x = J.random_element()
-                sage: y = J.random_element()
-                sage: Lx = x.matrix()
-                sage: Ly = y.matrix()
-                sage: Lxx = (x*x).matrix()
-                sage: Lxy = (x*y).matrix()
-                sage: bool(2*Lx*Lxy + Ly*Lxx == 2*Lxy*Lx + Lxx*Ly)
-                True
-
-            Test the second polarization identity from my notes or from
-            Baes (2.4)::
-
-                sage: set_random_seed()
-                sage: J = random_eja()
-                sage: x = J.random_element()
-                sage: y = J.random_element()
-                sage: z = J.random_element()
-                sage: Lx = x.matrix()
-                sage: Ly = y.matrix()
-                sage: Lz = z.matrix()
-                sage: Lzy = (z*y).matrix()
-                sage: Lxy = (x*y).matrix()
-                sage: Lxz = (x*z).matrix()
-                sage: bool(Lx*Lzy + Lz*Lxy + Ly*Lxz == Lzy*Lx + Lxy*Lz + Lxz*Ly)
-                True
-
-            Test the third polarization identity from my notes or from
-            Baes (2.5)::
-
-                sage: set_random_seed()
-                sage: J = random_eja()
-                sage: u = J.random_element()
-                sage: y = J.random_element()
-                sage: z = J.random_element()
-                sage: Lu = u.matrix()
-                sage: Ly = y.matrix()
-                sage: Lz = z.matrix()
-                sage: Lzy = (z*y).matrix()
-                sage: Luy = (u*y).matrix()
-                sage: Luz = (u*z).matrix()
-                sage: Luyz = (u*(y*z)).matrix()
-                sage: lhs = Lu*Lzy + Lz*Luy + Ly*Luz
-                sage: rhs = Luyz + Ly*Lu*Lz + Lz*Lu*Ly
-                sage: bool(lhs == rhs)
-                True
-
-            """
-            fda_elt = FiniteDimensionalAlgebraElement(self.parent(), self)
-            return fda_elt.matrix().transpose()
-
-
         def minimal_polynomial(self):
             """
             EXAMPLES::
@@ -538,6 +518,102 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
             return elt.minimal_polynomial()
 
 
+        def natural_representation(self):
+            """
+            Return a more-natural representation of this element.
+
+            Every finite-dimensional Euclidean Jordan Algebra is a
+            direct sum of five simple algebras, four of which comprise
+            Hermitian matrices. This method returns the original
+            "natural" representation of this element as a Hermitian
+            matrix, if it has one. If not, you get the usual representation.
+
+            EXAMPLES::
+
+                sage: J = ComplexHermitianSimpleEJA(3)
+                sage: J.one()
+                e0 + e5 + e8
+                sage: J.one().natural_representation()
+                [1 0 0 0 0 0]
+                [0 1 0 0 0 0]
+                [0 0 1 0 0 0]
+                [0 0 0 1 0 0]
+                [0 0 0 0 1 0]
+                [0 0 0 0 0 1]
+
+            """
+            B = self.parent().natural_basis()
+            W = B[0].matrix_space()
+            return W.linear_combination(zip(self.vector(), B))
+
+
+        def operator_matrix(self):
+            """
+            Return the matrix that represents left- (or right-)
+            multiplication by this element in the parent algebra.
+
+            We have to override this because the superclass method
+            returns a matrix that acts on row vectors (that is, on
+            the right).
+
+            EXAMPLES:
+
+            Test the first polarization identity from my notes, Koecher Chapter
+            III, or from Baes (2.3)::
+
+                sage: set_random_seed()
+                sage: J = random_eja()
+                sage: x = J.random_element()
+                sage: y = J.random_element()
+                sage: Lx = x.operator_matrix()
+                sage: Ly = y.operator_matrix()
+                sage: Lxx = (x*x).operator_matrix()
+                sage: Lxy = (x*y).operator_matrix()
+                sage: bool(2*Lx*Lxy + Ly*Lxx == 2*Lxy*Lx + Lxx*Ly)
+                True
+
+            Test the second polarization identity from my notes or from
+            Baes (2.4)::
+
+                sage: set_random_seed()
+                sage: J = random_eja()
+                sage: x = J.random_element()
+                sage: y = J.random_element()
+                sage: z = J.random_element()
+                sage: Lx = x.operator_matrix()
+                sage: Ly = y.operator_matrix()
+                sage: Lz = z.operator_matrix()
+                sage: Lzy = (z*y).operator_matrix()
+                sage: Lxy = (x*y).operator_matrix()
+                sage: Lxz = (x*z).operator_matrix()
+                sage: bool(Lx*Lzy + Lz*Lxy + Ly*Lxz == Lzy*Lx + Lxy*Lz + Lxz*Ly)
+                True
+
+            Test the third polarization identity from my notes or from
+            Baes (2.5)::
+
+                sage: set_random_seed()
+                sage: J = random_eja()
+                sage: u = J.random_element()
+                sage: y = J.random_element()
+                sage: z = J.random_element()
+                sage: Lu = u.operator_matrix()
+                sage: Ly = y.operator_matrix()
+                sage: Lz = z.operator_matrix()
+                sage: Lzy = (z*y).operator_matrix()
+                sage: Luy = (u*y).operator_matrix()
+                sage: Luz = (u*z).operator_matrix()
+                sage: Luyz = (u*(y*z)).operator_matrix()
+                sage: lhs = Lu*Lzy + Lz*Luy + Ly*Luz
+                sage: rhs = Luyz + Ly*Lu*Lz + Lz*Lu*Ly
+                sage: bool(lhs == rhs)
+                True
+
+            """
+            fda_elt = FiniteDimensionalAlgebraElement(self.parent(), self)
+            return fda_elt.matrix().transpose()
+
+
         def quadratic_representation(self, other=None):
             """
             Return the quadratic representation of this element.
@@ -610,9 +686,9 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
             elif not other in self.parent():
                 raise ArgumentError("'other' must live in the same algebra")
 
-            return ( self.matrix()*other.matrix()
-                       + other.matrix()*self.matrix()
-                       - (self*other).matrix() )
+            L = self.operator_matrix()
+            M = other.operator_matrix()
+            return ( L*M + M*L - (self*other).operator_matrix() )
 
 
         def span_of_powers(self):
@@ -645,7 +721,7 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
                 sage: set_random_seed()
                 sage: x = random_eja().random_element()
                 sage: u = x.subalgebra_generated_by().random_element()
-                sage: u.matrix()*u.vector() == (u**2).vector()
+                sage: u.operator_matrix()*u.vector() == (u**2).vector()
                 True
 
             """
@@ -717,7 +793,7 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
             s = 0
             minimal_dim = V.dimension()
             for i in xrange(1, V.dimension()):
-                this_dim = (u**i).matrix().image().dimension()
+                this_dim = (u**i).operator_matrix().image().dimension()
                 if this_dim < minimal_dim:
                     minimal_dim = this_dim
                     s = i
@@ -734,7 +810,7 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
             # Beware, solve_right() means that we're using COLUMN vectors.
             # Our FiniteDimensionalAlgebraElement superclass uses rows.
             u_next = u**(s+1)
-            A = u_next.matrix()
+            A = u_next.operator_matrix()
             c_coordinates = A.solve_right(u_next.vector())
 
             # Now c_coordinates is the idempotent we want, but it's in
@@ -864,7 +940,7 @@ def _real_symmetric_basis(n, field=QQ):
                 # Beware, orthogonal but not normalized!
                 Sij = Eij + Eij.transpose()
             S.append(Sij)
-    return S
+    return tuple(S)
 
 
 def _complex_hermitian_basis(n, field=QQ):
@@ -901,7 +977,7 @@ def _complex_hermitian_basis(n, field=QQ):
                 S.append(Sij_real)
                 Sij_imag = _embed_complex_matrix(I*Eij - I*Eij.transpose())
                 S.append(Sij_imag)
-    return S
+    return tuple(S)
 
 
 def _multiplication_table_from_matrix_basis(basis):
@@ -911,7 +987,10 @@ def _multiplication_table_from_matrix_basis(basis):
     multiplication on the right is matrix multiplication. Given a basis
     for the underlying matrix space, this function returns a
     multiplication table (obtained by looping through the basis
-    elements) for an algebra of those matrices.
+    elements) for an algebra of those matrices. A reordered copy
+    of the basis is also returned to work around the fact that
+    the ``span()`` in this function will change the order of the basis
+    from what we think it is, to... something else.
     """
     # In S^2, for example, we nominally have four coordinates even
     # though the space is of dimension three only. The vector space V
@@ -933,7 +1012,7 @@ def _multiplication_table_from_matrix_basis(basis):
     # Taking the span above reorders our basis (thanks, jerk!) so we
     # need to put our "matrix basis" in the same order as the
     # (reordered) vector basis.
-    S = [ vec2mat(b) for b in W.basis() ]
+    S = tuple( vec2mat(b) for b in W.basis() )
 
     Qs = []
     for s in S:
@@ -951,7 +1030,7 @@ def _multiplication_table_from_matrix_basis(basis):
         Q = matrix(field, W.dimension(), Q_rows)
         Qs.append(Q)
 
-    return Qs
+    return (Qs, S)
 
 
 def _embed_complex_matrix(M):
@@ -1058,9 +1137,12 @@ def RealSymmetricSimpleEJA(n, field=QQ):
 
     """
     S = _real_symmetric_basis(n, field=field)
-    Qs = _multiplication_table_from_matrix_basis(S)
+    (Qs, T) = _multiplication_table_from_matrix_basis(S)
 
-    return FiniteDimensionalEuclideanJordanAlgebra(field,Qs,rank=n)
+    return FiniteDimensionalEuclideanJordanAlgebra(field,
+                                                   Qs,
+                                                   rank=n,
+                                                   natural_basis=T)
 
 
 def ComplexHermitianSimpleEJA(n, field=QQ):
@@ -1082,8 +1164,11 @@ def ComplexHermitianSimpleEJA(n, field=QQ):
 
     """
     S = _complex_hermitian_basis(n)
-    Qs = _multiplication_table_from_matrix_basis(S)
-    return FiniteDimensionalEuclideanJordanAlgebra(field, Qs, rank=n)
+    (Qs, T) = _multiplication_table_from_matrix_basis(S)
+    return FiniteDimensionalEuclideanJordanAlgebra(field,
+                                                   Qs,
+                                                   rank=n,
+                                                   natural_basis=T)
 
 
 def QuaternionHermitianSimpleEJA(n):