]> gitweb.michael.orlitzky.com - sage.d.git/commitdiff
eja: make AA the default field because everything cool requires it.
authorMichael Orlitzky <michael@orlitzky.com>
Sun, 11 Oct 2020 21:24:17 +0000 (17:24 -0400)
committerMichael Orlitzky <michael@orlitzky.com>
Sun, 11 Oct 2020 21:24:17 +0000 (17:24 -0400)
mjo/eja/eja_algebra.py
mjo/eja/eja_element.py
mjo/eja/eja_element_subalgebra.py
mjo/eja/eja_operator.py
mjo/eja/eja_subalgebra.py

index 1f48f3c574454a1b59df75061f53b123c7ca986f..0a260653a7c5480acea6ffc00a9158169e4e7a18 100644 (file)
@@ -17,7 +17,7 @@ from sage.misc.lazy_import import lazy_import
 from sage.misc.prandom import choice
 from sage.misc.table import table
 from sage.modules.free_module import FreeModule, VectorSpace
-from sage.rings.all import (ZZ, QQ, RR, RLF, CLF,
+from sage.rings.all import (ZZ, QQ, AA, QQbar, RR, RLF, CLF,
                             PolynomialRing,
                             QuadraticField)
 from mjo.eja.eja_element import FiniteDimensionalEuclideanJordanAlgebraElement
@@ -207,8 +207,8 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule):
 
         Ensure that it says what we think it says::
 
-            sage: JordanSpinEJA(2, field=QQ)
-            Euclidean Jordan algebra of dimension 2 over Rational Field
+            sage: JordanSpinEJA(2, field=AA)
+            Euclidean Jordan algebra of dimension 2 over Algebraic Real Field
             sage: JordanSpinEJA(3, field=RDF)
             Euclidean Jordan algebra of dimension 3 over Real Double Field
 
@@ -551,8 +551,8 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule):
             Finite family {0: e0, 1: e1, 2: e2}
             sage: J.natural_basis()
             (
-            [1 0]  [        0 1/2*sqrt2]  [0 0]
-            [0 0], [1/2*sqrt2         0], [0 1]
+            [1 0]  [                  0 0.7071067811865475?]  [0 0]
+            [0 0], [0.7071067811865475?                   0], [0 1]
             )
 
         ::
@@ -1013,7 +1013,7 @@ class KnownRankEJA(object):
         return 5
 
     @classmethod
-    def random_instance(cls, field=QQ, **kwargs):
+    def random_instance(cls, field=AA, **kwargs):
         """
         Return a random instance of this type of algebra.
 
@@ -1069,7 +1069,7 @@ class HadamardEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
         (r0, r1, r2)
 
     """
-    def __init__(self, n, field=QQ, **kwargs):
+    def __init__(self, n, field=AA, **kwargs):
         V = VectorSpace(field, n)
         mult_table = [ [ V.gen(i)*(i == j) for j in range(n) ]
                        for i in range(n) ]
@@ -1102,7 +1102,7 @@ class HadamardEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
         return x.to_vector().inner_product(y.to_vector())
 
 
-def random_eja(field=QQ, nontrivial=False):
+def random_eja(field=AA, nontrivial=False):
     """
     Return a "random" finite-dimensional Euclidean Jordan Algebra.
 
@@ -1322,8 +1322,8 @@ class RealSymmetricEJA(RealMatrixEuclideanJordanAlgebra, KnownRankEJA):
 
     In theory, our "field" can be any subfield of the reals::
 
-        sage: RealSymmetricEJA(2, AA)
-        Euclidean Jordan algebra of dimension 3 over Algebraic Real Field
+        sage: RealSymmetricEJA(2, RDF)
+        Euclidean Jordan algebra of dimension 3 over Real Double Field
         sage: RealSymmetricEJA(2, RR)
         Euclidean Jordan algebra of dimension 3 over Real Field with
         53 bits of precision
@@ -1415,7 +1415,7 @@ class RealSymmetricEJA(RealMatrixEuclideanJordanAlgebra, KnownRankEJA):
         return 4 # Dimension 10
 
 
-    def __init__(self, n, field=QQ, **kwargs):
+    def __init__(self, n, field=AA, **kwargs):
         basis = self._denormalized_basis(n, field)
         super(RealSymmetricEJA, self).__init__(field, basis, n, **kwargs)
 
@@ -1435,7 +1435,7 @@ class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
 
         EXAMPLES::
 
-            sage: F = QuadraticField(-1, 'i')
+            sage: F = QuadraticField(-1, 'I')
             sage: x1 = F(4 - 2*i)
             sage: x2 = F(1 + 2*i)
             sage: x3 = F(-i)
@@ -1455,7 +1455,7 @@ class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
             sage: set_random_seed()
             sage: n_max = ComplexMatrixEuclideanJordanAlgebra._max_test_case_size()
             sage: n = ZZ.random_element(n_max)
-            sage: F = QuadraticField(-1, 'i')
+            sage: F = QuadraticField(-1, 'I')
             sage: X = random_matrix(F, n)
             sage: Y = random_matrix(F, n)
             sage: Xe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X)
@@ -1498,15 +1498,15 @@ class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
             ....:                 [ 9,  10, 11, 12],
             ....:                 [-10, 9, -12, 11] ])
             sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(A)
-            [  2*i + 1   4*i + 3]
-            [ 10*i + 9 12*i + 11]
+            [  2*I + 1   4*I + 3]
+            [ 10*I + 9 12*I + 11]
 
         TESTS:
 
         Unembedding is the inverse of embedding::
 
             sage: set_random_seed()
-            sage: F = QuadraticField(-1, 'i')
+            sage: F = QuadraticField(-1, 'I')
             sage: M = random_matrix(F, 3)
             sage: Me = ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
             sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
@@ -1524,7 +1524,12 @@ class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
         field = M.base_ring()
         R = PolynomialRing(field, 'z')
         z = R.gen()
-        F = field.extension(z**2 + 1, 'i', embedding=CLF(-1).sqrt())
+        if field is AA:
+            # Sage doesn't know how to embed AA into QQbar, i.e. how
+            # to adjoin sqrt(-1) to AA.
+            F = QQbar
+        else:
+            F = field.extension(z**2 + 1, 'I', embedding=CLF(-1).sqrt())
         i = F.gen()
 
         # Go top-left to bottom-right (reading order), converting every
@@ -1565,7 +1570,7 @@ class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
             sage: Ye = y.natural_representation()
             sage: X = ComplexHermitianEJA.real_unembed(Xe)
             sage: Y = ComplexHermitianEJA.real_unembed(Ye)
-            sage: expected = (X*Y).trace().vector()[0]
+            sage: expected = (X*Y).trace().real()
             sage: actual = ComplexHermitianEJA.natural_inner_product(Xe,Ye)
             sage: actual == expected
             True
@@ -1589,8 +1594,8 @@ class ComplexHermitianEJA(ComplexMatrixEuclideanJordanAlgebra, KnownRankEJA):
 
     In theory, our "field" can be any subfield of the reals::
 
-        sage: ComplexHermitianEJA(2, AA)
-        Euclidean Jordan algebra of dimension 4 over Algebraic Real Field
+        sage: ComplexHermitianEJA(2, RDF)
+        Euclidean Jordan algebra of dimension 4 over Real Double Field
         sage: ComplexHermitianEJA(2, RR)
         Euclidean Jordan algebra of dimension 4 over Real Field with
         53 bits of precision
@@ -1700,7 +1705,7 @@ class ComplexHermitianEJA(ComplexMatrixEuclideanJordanAlgebra, KnownRankEJA):
         return ( s.change_ring(field) for s in S )
 
 
-    def __init__(self, n, field=QQ, **kwargs):
+    def __init__(self, n, field=AA, **kwargs):
         basis = self._denormalized_basis(n,field)
         super(ComplexHermitianEJA,self).__init__(field, basis, n, **kwargs)
 
@@ -1752,7 +1757,7 @@ class QuaternionMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
         if M.ncols() != n:
             raise ValueError("the matrix 'M' must be square")
 
-        F = QuadraticField(-1, 'i')
+        F = QuadraticField(-1, 'I')
         i = F.gen()
 
         blocks = []
@@ -1828,10 +1833,10 @@ class QuaternionMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
                     raise ValueError('bad on-diagonal submatrix')
                 if submat[0,1] != -submat[1,0].conjugate():
                     raise ValueError('bad off-diagonal submatrix')
-                z  = submat[0,0].vector()[0]   # real part
-                z += submat[0,0].vector()[1]*i # imag part
-                z += submat[0,1].vector()[0]*j # real part
-                z += submat[0,1].vector()[1]*k # imag part
+                z  = submat[0,0].real()
+                z += submat[0,0].imag()*i
+                z += submat[0,1].real()*j
+                z += submat[0,1].imag()*k
                 elements.append(z)
 
         return matrix(Q, n/4, elements)
@@ -1884,8 +1889,8 @@ class QuaternionHermitianEJA(QuaternionMatrixEuclideanJordanAlgebra,
 
     In theory, our "field" can be any subfield of the reals::
 
-        sage: QuaternionHermitianEJA(2, AA)
-        Euclidean Jordan algebra of dimension 6 over Algebraic Real Field
+        sage: QuaternionHermitianEJA(2, RDF)
+        Euclidean Jordan algebra of dimension 6 over Real Double Field
         sage: QuaternionHermitianEJA(2, RR)
         Euclidean Jordan algebra of dimension 6 over Real Field with
         53 bits of precision
@@ -1996,7 +2001,7 @@ class QuaternionHermitianEJA(QuaternionMatrixEuclideanJordanAlgebra,
         return ( s.change_ring(field) for s in S )
 
 
-    def __init__(self, n, field=QQ, **kwargs):
+    def __init__(self, n, field=AA, **kwargs):
         basis = self._denormalized_basis(n,field)
         super(QuaternionHermitianEJA,self).__init__(field, basis, n, **kwargs)
 
@@ -2056,7 +2061,7 @@ class BilinearFormEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
         sage: actual == expected
         True
     """
-    def __init__(self, n, field=QQ, B=None, **kwargs):
+    def __init__(self, n, field=AA, B=None, **kwargs):
         if B is None:
             self._B = matrix.identity(field, max(0,n-1))
         else:
@@ -2169,7 +2174,7 @@ class JordanSpinEJA(BilinearFormEJA):
             True
 
     """
-    def __init__(self, n, field=QQ, **kwargs):
+    def __init__(self, n, field=AA, **kwargs):
         # This is a special case of the BilinearFormEJA with the identity
         # matrix as its bilinear form.
         return super(JordanSpinEJA, self).__init__(n, field, **kwargs)
@@ -2199,12 +2204,12 @@ class TrivialEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
         sage: J.one().norm()
         0
         sage: J.one().subalgebra_generated_by()
-        Euclidean Jordan algebra of dimension 0 over Rational Field
+        Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
         sage: J.rank()
         0
 
     """
-    def __init__(self, field=QQ, **kwargs):
+    def __init__(self, field=AA, **kwargs):
         mult_table = []
         fdeja = super(TrivialEJA, self)
         # The rank is zero using my definition, namely the dimension of the
index 7c4c79ddcd7315e654620a0be8f8bccf5ab9ac11..0f6a47cd4f10efbcb0298725c4ae26537eae6372 100644 (file)
@@ -113,7 +113,7 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement):
         We should always get back an element of the algebra::
 
             sage: set_random_seed()
-            sage: p = PolynomialRing(QQ, 't').random_element()
+            sage: p = PolynomialRing(AA, 't').random_element()
             sage: J = random_eja()
             sage: x = J.random_element()
             sage: x.apply_univariate_polynomial(p) in J
@@ -575,7 +575,7 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement):
         The spectral decomposition of a non-regular element should always
         contain at least one non-minimal idempotent::
 
-            sage: J = RealSymmetricEJA(3, AA)
+            sage: J = RealSymmetricEJA(3)
             sage: x = sum(J.gens())
             sage: x.is_regular()
             False
@@ -586,7 +586,7 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement):
         On the other hand, the spectral decomposition of a regular
         element should always be in terms of minimal idempotents::
 
-            sage: J = JordanSpinEJA(4, AA)
+            sage: J = JordanSpinEJA(4)
             sage: x = sum( i*J.gens()[i] for i in range(len(J.gens())) )
             sage: x.is_regular()
             True
@@ -909,9 +909,9 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement):
             sage: set_random_seed()
             sage: n_max = RealSymmetricEJA._max_test_case_size()
             sage: n = ZZ.random_element(1, n_max)
-            sage: J1 = RealSymmetricEJA(n,QQ)
-            sage: J2 = RealSymmetricEJA(n,QQ,normalize_basis=False)
-            sage: X = random_matrix(QQ,n)
+            sage: J1 = RealSymmetricEJA(n)
+            sage: J2 = RealSymmetricEJA(n,normalize_basis=False)
+            sage: X = random_matrix(AA,n)
             sage: X = X*X.transpose()
             sage: x1 = J1(X)
             sage: x2 = J2(X)
@@ -1003,7 +1003,7 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement):
             sage: J = HadamardEJA(2)
             sage: x = sum(J.gens())
             sage: x.norm()
-            sqrt(2)
+            1.414213562373095?
 
         ::
 
@@ -1065,10 +1065,10 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement):
             sage: n = x_vec.degree()
             sage: x0 = x_vec[0]
             sage: x_bar = x_vec[1:]
-            sage: A = matrix(QQ, 1, [x_vec.inner_product(x_vec)])
+            sage: A = matrix(AA, 1, [x_vec.inner_product(x_vec)])
             sage: B = 2*x0*x_bar.row()
             sage: C = 2*x0*x_bar.column()
-            sage: D = matrix.identity(QQ, n-1)
+            sage: D = matrix.identity(AA, n-1)
             sage: D = (x0^2 - x_bar.inner_product(x_bar))*D
             sage: D = D + 2*x_bar.tensor_product(x_bar)
             sage: Q = matrix.block(2,2,[A,B,C,D])
@@ -1192,7 +1192,7 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement):
         The spectral decomposition of the identity is ``1`` times itself,
         and the spectral decomposition of zero is ``0`` times the identity::
 
-            sage: J = RealSymmetricEJA(3,AA)
+            sage: J = RealSymmetricEJA(3)
             sage: J.one()
             e0 + e2 + e5
             sage: J.one().spectral_decomposition()
@@ -1202,7 +1202,7 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement):
 
         TESTS::
 
-            sage: J = RealSymmetricEJA(4,AA)
+            sage: J = RealSymmetricEJA(4)
             sage: x = sum(J.gens())
             sage: sd = x.spectral_decomposition()
             sage: l0 = sd[0][0]
@@ -1453,14 +1453,14 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement):
             sage: J = HadamardEJA(2)
             sage: x = sum(J.gens())
             sage: x.trace_norm()
-            sqrt(2)
+            1.414213562373095?
 
         ::
 
             sage: J = JordanSpinEJA(4)
             sage: x = sum(J.gens())
             sage: x.trace_norm()
-            2*sqrt(2)
+            2.828427124746190?
 
         """
         return self.trace_inner_product(self).sqrt()
index 292871afb86415b845f1e1feee2a5560e75f3cf2..67342c6be3df71a63eb4cf82337f52ff7edf5fea 100644 (file)
@@ -130,7 +130,7 @@ class FiniteDimensionalEuclideanJordanElementSubalgebra(FiniteDimensionalEuclide
         reals with an orthonormal basis::
 
             sage: set_random_seed()
-            sage: x = random_eja(AA).random_element()
+            sage: x = random_eja().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
@@ -151,7 +151,7 @@ class FiniteDimensionalEuclideanJordanElementSubalgebra(FiniteDimensionalEuclide
         the algebraic reals with an orthonormal basis::
 
             sage: set_random_seed()
-            sage: x = random_eja(AA).random_element()
+            sage: x = random_eja().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())
index ee33dbf53b36fd9851d31169d5699041f460328e..667e3d5acba051e08bf461e12aafd9ae2437cc74 100644 (file)
@@ -125,9 +125,9 @@ class FiniteDimensionalEuclideanJordanAlgebraOperator(Map):
             sage: J1 = JordanSpinEJA(3)
             sage: J2 = HadamardEJA(2)
             sage: J3 = RealSymmetricEJA(1)
-            sage: mat1 = matrix(QQ, [[1,2,3],
+            sage: mat1 = matrix(AA, [[1,2,3],
             ....:                    [4,5,6]])
-            sage: mat2 = matrix(QQ, [[7,8]])
+            sage: mat2 = matrix(AA, [[7,8]])
             sage: g = FiniteDimensionalEuclideanJordanAlgebraOperator(J1,
             ....:                                                     J2,
             ....:                                                     mat1)
@@ -139,9 +139,9 @@ class FiniteDimensionalEuclideanJordanAlgebraOperator(Map):
             algebras represented by the matrix:
             [39 54 69]
             Domain: Euclidean Jordan algebra of dimension 3 over
-            Rational Field
+            Algebraic Real Field
             Codomain: Euclidean Jordan algebra of dimension 1 over
-            Rational Field
+            Algebraic Real Field
 
         """
         return FiniteDimensionalEuclideanJordanAlgebraOperator(
@@ -341,9 +341,9 @@ class FiniteDimensionalEuclideanJordanAlgebraOperator(Map):
             [1 0]
             [0 1]
             Domain: Euclidean Jordan algebra of dimension 2 over
-            Rational Field
+            Algebraic Real Field
             Codomain: Euclidean Jordan algebra of dimension 2 over
-            Rational Field
+            Algebraic Real Field
 
         """
         msg = ("Linear operator between finite-dimensional Euclidean Jordan "
@@ -542,7 +542,7 @@ class FiniteDimensionalEuclideanJordanAlgebraOperator(Map):
 
         EXAMPLES::
 
-            sage: J = RealSymmetricEJA(4,AA)
+            sage: J = RealSymmetricEJA(4)
             sage: x = sum(J.gens())
             sage: A = x.subalgebra_generated_by(orthonormalize_basis=True)
             sage: L0x = A(x).operator()
index 0be85616678d4f7f7ffa276ddcd13ab907c6e116..ac77f22a691e7cc04138e34cae46f148909de0f8 100644 (file)
@@ -95,9 +95,9 @@ class FiniteDimensionalEuclideanJordanSubalgebra(FiniteDimensionalEuclideanJorda
     matrices do not contain the superalgebra's identity element::
 
         sage: J = RealSymmetricEJA(2)
-        sage: E11 = matrix(QQ, [ [1,0],
+        sage: E11 = matrix(AA, [ [1,0],
         ....:                    [0,0] ])
-        sage: E22 = matrix(QQ, [ [0,0],
+        sage: E22 = matrix(AA, [ [0,0],
         ....:                    [0,1] ])
         sage: K1 = FiniteDimensionalEuclideanJordanSubalgebra(J, (J(E11),))
         sage: K1.one().natural_representation()
@@ -198,7 +198,7 @@ class FiniteDimensionalEuclideanJordanSubalgebra(FiniteDimensionalEuclideanJorda
         EXAMPLES::
 
             sage: J = RealSymmetricEJA(3)
-            sage: X = matrix(QQ, [ [0,0,1],
+            sage: X = matrix(AA, [ [0,0,1],
             ....:                  [0,1,0],
             ....:                  [1,0,0] ])
             sage: x = J(X)
@@ -249,10 +249,10 @@ class FiniteDimensionalEuclideanJordanSubalgebra(FiniteDimensionalEuclideanJorda
         EXAMPLES::
 
             sage: J = RealSymmetricEJA(3)
-            sage: E11 = matrix(QQ, [ [1,0,0],
+            sage: E11 = matrix(ZZ, [ [1,0,0],
             ....:                    [0,0,0],
             ....:                    [0,0,0] ])
-            sage: E22 = matrix(QQ, [ [0,0,0],
+            sage: E22 = matrix(ZZ, [ [0,0,0],
             ....:                    [0,1,0],
             ....:                    [0,0,0] ])
             sage: b1 = J(E11)