]> gitweb.michael.orlitzky.com - sage.d.git/blobdiff - mjo/eja/eja_algebra.py
eja: fix tests and pre-cache ranks.
[sage.d.git] / mjo / eja / eja_algebra.py
index 0a260653a7c5480acea6ffc00a9158169e4e7a18..19db8b0ccefef297623d81721d4a774f99a9d3d1 100644 (file)
@@ -54,7 +54,6 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule):
     def __init__(self,
                  field,
                  mult_table,
-                 rank,
                  prefix='e',
                  category=None,
                  natural_basis=None,
@@ -91,7 +90,6 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule):
                 # a real embedding.
                 raise ValueError('field is not real')
 
-        self._rank = rank
         self._natural_basis = natural_basis
 
         if category is None:
@@ -773,108 +771,6 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule):
         return (J0, J5, J1)
 
 
-    def a_jordan_frame(self):
-        r"""
-        Generate a Jordan frame for this algebra.
-
-        This implementation is based on the so-called "central
-        orthogonal idempotents" implemented for (semisimple) centers
-        of SageMath ``FiniteDimensionalAlgebrasWithBasis``. Since all
-        Euclidean Jordan algebas are commutative (and thus equal to
-        their own centers) and semisimple, the method should work more
-        or less as implemented, if it ever worked in the first place.
-        (I don't know the justification for the original implementation.
-        yet).
-
-        How it works: we loop through the algebras generators, looking
-        for their eigenspaces. If there's more than one eigenspace,
-        and if they result in more than one subalgebra, then we split
-        those subalgebras recursively until we get to subalgebras of
-        dimension one (whose idempotent is the unit element). Why does
-        some generator have to produce at least two subalgebras? I
-        dunno. But it seems to work.
-
-        Beware that Koecher defines the "center" of a Jordan algebra to
-        be something else, because the usual definition is stupid in a
-        (necessarily commutative) Jordan algebra.
-
-        SETUP::
-
-            sage: from mjo.eja.eja_algebra import (random_eja,
-            ....:                                  JordanSpinEJA,
-            ....:                                  TrivialEJA)
-
-        EXAMPLES:
-
-        A Jordan frame for the trivial algebra has to be empty
-        (zero-length) since its rank is zero. More to the point, there
-        are no non-zero idempotents in the trivial EJA. This does not
-        cause any problems so long as we adopt the convention that the
-        empty sum is zero, since then the sole element of the trivial
-        EJA has an (empty) spectral decomposition::
-
-            sage: J = TrivialEJA()
-            sage: J.a_jordan_frame()
-            ()
-
-        A one-dimensional algebra has rank one (equal to its dimension),
-        and only one primitive idempotent, namely the algebra's unit
-        element::
-
-            sage: J = JordanSpinEJA(1)
-            sage: J.a_jordan_frame()
-            (e0,)
-
-        TESTS::
-
-            sage: J = random_eja()
-            sage: c = J.a_jordan_frame()
-            sage: all( x^2 == x for x in c )
-            True
-            sage: r = len(c)
-            sage: all( c[i]*c[j] == c[i]*(i==j) for i in range(r)
-            ....:                               for j in range(r) )
-            True
-
-        """
-        if self.dimension() == 0:
-            return ()
-        if self.dimension() == 1:
-            return (self.one(),)
-
-        for g in self.gens():
-            eigenpairs = g.operator().matrix().right_eigenspaces()
-            if len(eigenpairs) >= 2:
-                subalgebras = []
-                for eigval, eigspace in eigenpairs:
-                    # Make sub-EJAs from the matrix eigenspaces...
-                    sb = tuple( self.from_vector(b) for b in eigspace.basis() )
-                    try:
-                        # This will fail if e.g. the eigenspace basis
-                        # contains two elements and their product
-                        # isn't a linear combination of the two of
-                        # them (i.e. the generated EJA isn't actually
-                        # two dimensional).
-                        s = FiniteDimensionalEuclideanJordanSubalgebra(self, sb)
-                        subalgebras.append(s)
-                    except ArithmeticError as e:
-                        if str(e) == "vector is not in free module":
-                            # Ignore only the "not a sub-EJA" error
-                            pass
-
-                if len(subalgebras) >= 2:
-                    # apply this method recursively.
-                    return tuple( c.superalgebra_element()
-                                  for subalgebra in subalgebras
-                                  for c in subalgebra.a_jordan_frame() )
-
-        # If we got here, the algebra didn't decompose, at least not when we looked at
-        # the eigenspaces corresponding only to basis elements of the algebra. The
-        # implementation I stole says that this should work because of Schur's Lemma,
-        # so I personally blame Schur's Lemma if it does not.
-        raise Exception("Schur's Lemma didn't work!")
-
-
     def random_elements(self, count):
         """
         Return ``count`` random elements as a tuple.
@@ -893,23 +789,26 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule):
             True
 
         """
-        return  tuple( self.random_element() for idx in range(count) )
-
+        return tuple( self.random_element() for idx in range(count) )
 
+    @cached_method
     def rank(self):
         """
         Return the rank of this EJA.
 
         ALGORITHM:
 
-        The author knows of no algorithm to compute the rank of an EJA
-        where only the multiplication table is known. In lieu of one, we
-        require the rank to be specified when the algebra is created,
-        and simply pass along that number here.
+        We first compute the polynomial "column matrices" `p_{k}` that
+        evaluate to `x^k` on the coordinates of `x`. Then, we begin
+        adding them to a matrix one at a time, and trying to solve the
+        system that makes `p_{0}`,`p_{1}`,..., `p_{s-1}` add up to
+        `p_{s}`. This will succeed only when `s` is the rank of the
+        algebra, as proven in a recent draft paper of mine.
 
         SETUP::
 
-            sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
+            sage: from mjo.eja.eja_algebra import (HadamardEJA,
+            ....:                                  JordanSpinEJA,
             ....:                                  RealSymmetricEJA,
             ....:                                  ComplexHermitianEJA,
             ....:                                  QuaternionHermitianEJA,
@@ -950,8 +849,80 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule):
             sage: r > 0 or (r == 0 and J.is_trivial())
             True
 
+        Ensure that computing the rank actually works, since the ranks
+        of all simple algebras are known and will be cached by default::
+
+            sage: J = HadamardEJA(4)
+            sage: J.rank.clear_cache()
+            sage: J.rank()
+            4
+
+        ::
+
+            sage: J = JordanSpinEJA(4)
+            sage: J.rank.clear_cache()
+            sage: J.rank()
+            2
+
+        ::
+
+            sage: J = RealSymmetricEJA(3)
+            sage: J.rank.clear_cache()
+            sage: J.rank()
+            3
+
+        ::
+
+            sage: J = ComplexHermitianEJA(2)
+            sage: J.rank.clear_cache()
+            sage: J.rank()
+            2
+
+        ::
+
+            sage: J = QuaternionHermitianEJA(2)
+            sage: J.rank.clear_cache()
+            sage: J.rank()
+            2
+
         """
-        return self._rank
+        n = self.dimension()
+        if n == 0:
+            return 0
+        elif n == 1:
+            return 1
+
+        var_names = [ "X" + str(z) for z in range(1,n+1) ]
+        R = PolynomialRing(self.base_ring(), var_names)
+        vars = R.gens()
+
+        def L_x_i_j(i,j):
+            # From a result in my book, these are the entries of the
+            # basis representation of L_x.
+            return sum( vars[k]*self.monomial(k).operator().matrix()[i,j]
+                        for k in range(n) )
+
+        L_x = matrix(R, n, n, L_x_i_j)
+        x_powers = [ vars[k]*(L_x**k)*self.one().to_vector()
+                     for k in range(n) ]
+
+        # Can assume n >= 2
+        M = matrix([x_powers[0]])
+        old_rank = 1
+
+        for d in range(1,n):
+            M = matrix(M.rows() + [x_powers[d]])
+            M.echelonize()
+            # TODO: we've basically solved the system here.
+            # We should save the echelonized matrix somehow
+            # so that it can be reused in the charpoly method.
+            new_rank = M.rank()
+            if new_rank == old_rank:
+                return new_rank
+            else:
+                old_rank = new_rank
+
+        return n
 
 
     def vector_space(self):
@@ -1075,7 +1046,8 @@ class HadamardEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
                        for i in range(n) ]
 
         fdeja = super(HadamardEJA, self)
-        return fdeja.__init__(field, mult_table, rank=n, **kwargs)
+        fdeja.__init__(field, mult_table, **kwargs)
+        self.rank.set_cache(n)
 
     def inner_product(self, x, y):
         """
@@ -1134,7 +1106,7 @@ class MatrixEuclideanJordanAlgebra(FiniteDimensionalEuclideanJordanAlgebra):
         # field can have dimension 4 (quaternions) too.
         return 2
 
-    def __init__(self, field, basis, rank, normalize_basis=True, **kwargs):
+    def __init__(self, field, basis, normalize_basis=True, **kwargs):
         """
         Compared to the superclass constructor, we take a basis instead of
         a multiplication table because the latter can be computed in terms
@@ -1147,7 +1119,7 @@ class MatrixEuclideanJordanAlgebra(FiniteDimensionalEuclideanJordanAlgebra):
         # time to ensure that it isn't a generator expression.
         basis = tuple(basis)
 
-        if rank > 1 and normalize_basis:
+        if len(basis) > 1 and normalize_basis:
             # 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.
@@ -1164,12 +1136,31 @@ class MatrixEuclideanJordanAlgebra(FiniteDimensionalEuclideanJordanAlgebra):
         Qs = self.multiplication_table_from_matrix_basis(basis)
 
         fdeja = super(MatrixEuclideanJordanAlgebra, self)
-        return fdeja.__init__(field,
-                              Qs,
-                              rank=rank,
-                              natural_basis=basis,
-                              **kwargs)
+        fdeja.__init__(field, Qs, natural_basis=basis, **kwargs)
+        return
+
 
+    @cached_method
+    def rank(self):
+        r"""
+        Override the parent method with something that tries to compute
+        over a faster (non-extension) field.
+        """
+        if self._basis_normalizers is None:
+            # We didn't normalize, so assume that the basis we started
+            # with had entries in a nice field.
+            return super(MatrixEuclideanJordanAlgebra, self).rank()
+        else:
+            basis = ( (b/n) for (b,n) in zip(self.natural_basis(),
+                                             self._basis_normalizers) )
+
+            # Do this over the rationals and convert back at the end.
+            # Only works because we know the entries of the basis are
+            # integers.
+            J = MatrixEuclideanJordanAlgebra(QQ,
+                                             basis,
+                                             normalize_basis=False)
+            return J.rank()
 
     @cached_method
     def _charpoly_coeff(self, i):
@@ -1188,7 +1179,6 @@ class MatrixEuclideanJordanAlgebra(FiniteDimensionalEuclideanJordanAlgebra):
             # Do this over the rationals and convert back at the end.
             J = MatrixEuclideanJordanAlgebra(QQ,
                                              basis,
-                                             self.rank(),
                                              normalize_basis=False)
             (_,x,_,_) = J._charpoly_matrix_system()
             p = J._charpoly_coeff(i)
@@ -1417,7 +1407,8 @@ class RealSymmetricEJA(RealMatrixEuclideanJordanAlgebra, KnownRankEJA):
 
     def __init__(self, n, field=AA, **kwargs):
         basis = self._denormalized_basis(n, field)
-        super(RealSymmetricEJA, self).__init__(field, basis, n, **kwargs)
+        super(RealSymmetricEJA, self).__init__(field, basis, **kwargs)
+        self.rank.set_cache(n)
 
 
 class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
@@ -1707,7 +1698,8 @@ class ComplexHermitianEJA(ComplexMatrixEuclideanJordanAlgebra, KnownRankEJA):
 
     def __init__(self, n, field=AA, **kwargs):
         basis = self._denormalized_basis(n,field)
-        super(ComplexHermitianEJA,self).__init__(field, basis, n, **kwargs)
+        super(ComplexHermitianEJA,self).__init__(field, basis, **kwargs)
+        self.rank.set_cache(n)
 
 
 class QuaternionMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
@@ -2003,7 +1995,8 @@ class QuaternionHermitianEJA(QuaternionMatrixEuclideanJordanAlgebra,
 
     def __init__(self, n, field=AA, **kwargs):
         basis = self._denormalized_basis(n,field)
-        super(QuaternionHermitianEJA,self).__init__(field, basis, n, **kwargs)
+        super(QuaternionHermitianEJA,self).__init__(field, basis, **kwargs)
+        self.rank.set_cache(n)
 
 
 class BilinearFormEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
@@ -2086,7 +2079,8 @@ class BilinearFormEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
         # one-dimensional ambient space (because the rank is bounded
         # by the ambient dimension).
         fdeja = super(BilinearFormEJA, self)
-        return fdeja.__init__(field, mult_table, rank=min(n,2), **kwargs)
+        fdeja.__init__(field, mult_table, **kwargs)
+        self.rank.set_cache(min(n,2))
 
     def inner_product(self, x, y):
         r"""
@@ -2214,4 +2208,5 @@ class TrivialEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
         fdeja = super(TrivialEJA, self)
         # The rank is zero using my definition, namely the dimension of the
         # largest subalgebra generated by any element.
-        return fdeja.__init__(field, mult_table, rank=0, **kwargs)
+        fdeja.__init__(field, mult_table, **kwargs)
+        self.rank.set_cache(0)