]> gitweb.michael.orlitzky.com - sage.d.git/blobdiff - mjo/eja/eja_subalgebra.py
eja: drop custom gram_schmidt() routine that isn't noticeably faster.
[sage.d.git] / mjo / eja / eja_subalgebra.py
index 8f6e56b55f309d10967ee5aae2bb8b2fe7261566..c372e50072c73a50281dd45d07eec62a62848277 100644 (file)
@@ -3,7 +3,6 @@ from sage.matrix.constructor import matrix
 from mjo.eja.eja_algebra import FiniteDimensionalEuclideanJordanAlgebra
 from mjo.eja.eja_element import FiniteDimensionalEuclideanJordanAlgebraElement
 
-
 class FiniteDimensionalEuclideanJordanElementSubalgebraElement(FiniteDimensionalEuclideanJordanAlgebraElement):
     """
     SETUP::
@@ -23,6 +22,17 @@ class FiniteDimensionalEuclideanJordanElementSubalgebraElement(FiniteDimensional
         sage: actual == expected
         True
 
+    The left-multiplication-by operator for elements in the subalgebra
+    works like it does in the superalgebra, even if we orthonormalize
+    our basis::
+
+        sage: set_random_seed()
+        sage: x = random_eja(AA).random_element()
+        sage: A = x.subalgebra_generated_by(orthonormalize_basis=True)
+        sage: y = A.random_element()
+        sage: y.operator()(A.one()) == y
+        True
+
     """
 
     def superalgebra_element(self):
@@ -99,7 +109,7 @@ class FiniteDimensionalEuclideanJordanElementSubalgebra(FiniteDimensionalEuclide
         1
 
     """
-    def __init__(self, elt):
+    def __init__(self, elt, orthonormalize_basis):
         self._superalgebra = elt.parent()
         category = self._superalgebra.category().Associative()
         V = self._superalgebra.vector_space()
@@ -117,44 +127,44 @@ class FiniteDimensionalEuclideanJordanElementSubalgebra(FiniteDimensionalEuclide
         except ValueError:
             prefix = prefixen[0]
 
-        if elt.is_zero():
-            # Short circuit because 0^0 == 1 is going to make us
-            # think we have a one-dimensional algebra otherwise.
-            natural_basis = tuple()
-            mult_table = tuple()
-            rank = 0
-            self._vector_space = V.zero_subspace()
-            self._superalgebra_basis = []
-            fdeja = super(FiniteDimensionalEuclideanJordanElementSubalgebra,
-                          self)
-            return fdeja.__init__(field,
-                                  mult_table,
-                                  rank,
-                                  prefix=prefix,
-                                  category=category,
-                                  natural_basis=natural_basis)
-
-
-        # First compute the vector subspace spanned by the powers of
-        # the given element.
+        # This list is guaranteed to contain all independent powers,
+        # because it's the maximal set of powers that could possibly
+        # be independent (by a dimension argument).
         powers = [ elt**k for k in range(V.dimension()) ]
         power_vectors = [ p.to_vector() for p in powers ]
+        P = matrix(field, power_vectors)
 
-        # Figure out which powers form a linearly-independent set.
-        ind_rows = matrix(field, power_vectors).pivot_rows()
+        if orthonormalize_basis == False:
+            # In this case, we just need to figure out which elements
+            # of the "powers" list are redundant... First compute the
+            # vector subspace spanned by the powers of the given
+            # element.
 
-        # Pick those out of the list of all powers.
-        superalgebra_basis = tuple(map(powers.__getitem__, ind_rows))
+            # Figure out which powers form a linearly-independent set.
+            ind_rows = P.pivot_rows()
 
-        # If our superalgebra is a subalgebra of something else, then
-        # these vectors won't have the right coordinates for
-        # V.span_of_basis() unless we use V.from_vector() on them.
-        basis_vectors = map(power_vectors.__getitem__, ind_rows)
-        W = V.span_of_basis( V.from_vector(v) for v in basis_vectors )
+            # Pick those out of the list of all powers.
+            superalgebra_basis = tuple(map(powers.__getitem__, ind_rows))
 
-        # Now figure out the entries of the right-multiplication
-        # matrix for the successive basis elements b0, b1,... of
-        # that subspace.
+            # If our superalgebra is a subalgebra of something else, then
+            # these vectors won't have the right coordinates for
+            # V.span_of_basis() unless we use V.from_vector() on them.
+            basis_vectors = map(power_vectors.__getitem__, ind_rows)
+        else:
+            # If we're going to orthonormalize the basis anyway, we
+            # might as well just do Gram-Schmidt on the whole list of
+            # powers. The redundant ones will get zero'd out. If this
+            # looks like a roundabout way to orthonormalize, it is.
+            # But converting everything from algebra elements to vectors
+            # to matrices and then back again turns out to be about
+            # as fast as reimplementing our own Gram-Schmidt that
+            # works in an EJA.
+            G,_ = P.gram_schmidt(orthonormal=True)
+            basis_vectors = [ g for g in G.rows() if not g.is_zero() ]
+            superalgebra_basis = [ self._superalgebra.from_vector(b)
+                                   for b in basis_vectors ]
+
+        W = V.span_of_basis( V.from_vector(v) for v in basis_vectors )
         n = len(superalgebra_basis)
         mult_table = [[W.zero() for i in range(n)] for j in range(n)]
         for i in range(n):
@@ -231,7 +241,7 @@ class FiniteDimensionalEuclideanJordanElementSubalgebra(FiniteDimensionalEuclide
 
             sage: J = RealSymmetricEJA(3)
             sage: x = sum( i*J.gens()[i] for i in range(6) )
-            sage: K = FiniteDimensionalEuclideanJordanElementSubalgebra(x)
+            sage: K = FiniteDimensionalEuclideanJordanElementSubalgebra(x,False)
             sage: [ K(x^k) for k in range(J.rank()) ]
             [f0, f1, f2]
 
@@ -249,22 +259,16 @@ class FiniteDimensionalEuclideanJordanElementSubalgebra(FiniteDimensionalEuclide
             return self.from_vector(coords)
 
 
-    def one_basis(self):
-        """
-        Return the basis-element-index of this algebra's unit element.
-        """
-        return 0
-
 
     def one(self):
         """
         Return the multiplicative identity element of this algebra.
 
         The superclass method computes the identity element, which is
-        beyond overkill in this case: the algebra identity should be our
-        first basis element. We implement this via :meth:`one_basis`
-        because that method can optionally be used by other parts of the
-        category framework.
+        beyond overkill in this case: the superalgebra identity
+        restricted to this algebra is its identity. Note that we can't
+        count on the first basis element being the identity -- it migth
+        have been scaled if we orthonormalized the basis.
 
         SETUP::
 
@@ -285,27 +289,54 @@ class FiniteDimensionalEuclideanJordanElementSubalgebra(FiniteDimensionalEuclide
 
         TESTS:
 
-        The identity element acts like the identity::
+        The identity element acts like the identity over the rationals::
 
             sage: set_random_seed()
-            sage: J = random_eja().random_element().subalgebra_generated_by()
-            sage: x = J.random_element()
-            sage: J.one()*x == x and x*J.one() == x
+            sage: x = random_eja().random_element()
+            sage: A = x.subalgebra_generated_by()
+            sage: x = A.random_element()
+            sage: A.one()*x == x and x*A.one() == x
             True
 
-        The matrix of the unit element's operator is the identity::
+        The identity element acts like the identity over the algebraic
+        reals with an orthonormal basis::
 
             sage: set_random_seed()
-            sage: J = random_eja().random_element().subalgebra_generated_by()
-            sage: actual = J.one().operator().matrix()
-            sage: expected = matrix.identity(J.base_ring(), J.dimension())
+            sage: x = random_eja(AA).random_element()
+            sage: A = x.subalgebra_generated_by(orthonormalize_basis=True)
+            sage: x = A.random_element()
+            sage: A.one()*x == x and x*A.one() == x
+            True
+
+        The matrix of the unit element's operator is the identity over
+        the rationals::
+
+            sage: set_random_seed()
+            sage: x = random_eja().random_element()
+            sage: A = x.subalgebra_generated_by()
+            sage: actual = A.one().operator().matrix()
+            sage: expected = matrix.identity(A.base_ring(), A.dimension())
             sage: actual == expected
             True
+
+        The matrix of the unit element's operator is the identity over
+        the algebraic reals with an orthonormal basis::
+
+            sage: set_random_seed()
+            sage: x = random_eja(AA).random_element()
+            sage: A = x.subalgebra_generated_by(orthonormalize_basis=True)
+            sage: actual = A.one().operator().matrix()
+            sage: expected = matrix.identity(A.base_ring(), A.dimension())
+            sage: actual == expected
+            True
+
         """
         if self.dimension() == 0:
             return self.zero()
         else:
-            return self.monomial(self.one_basis())
+            sa_one = self.superalgebra().one().to_vector()
+            sa_coords = self.vector_space().coordinate_vector(sa_one)
+            return self.from_vector(sa_coords)
 
 
     def natural_basis_space(self):
@@ -338,7 +369,7 @@ class FiniteDimensionalEuclideanJordanElementSubalgebra(FiniteDimensionalEuclide
 
             sage: J = RealSymmetricEJA(3)
             sage: x = J.monomial(0) + 2*J.monomial(2) + 5*J.monomial(5)
-            sage: K = FiniteDimensionalEuclideanJordanElementSubalgebra(x)
+            sage: K = FiniteDimensionalEuclideanJordanElementSubalgebra(x,False)
             sage: K.vector_space()
             Vector space of degree 6 and dimension 3 over...
             User basis matrix: