X-Git-Url: http://gitweb.michael.orlitzky.com/?a=blobdiff_plain;f=mjo%2Feja%2Feja_algebra.py;h=d90d3f2dcb5f4adc5cd6135ac4151ef63c3166c8;hb=6c983c4d14a02b4eff37fd2a07ae6b32b93e611c;hp=f7983c857d8fdceffac3dbdfeaa6bf70fd7e811e;hpb=97d4dcfae8b0c5ab1577da2ba3629e0a4169c789;p=sage.d.git diff --git a/mjo/eja/eja_algebra.py b/mjo/eja/eja_algebra.py index f7983c8..d90d3f2 100644 --- a/mjo/eja/eja_algebra.py +++ b/mjo/eja/eja_algebra.py @@ -15,7 +15,7 @@ from sage.misc.prandom import choice from sage.misc.table import table from sage.modules.free_module import FreeModule, VectorSpace from sage.rings.integer_ring import ZZ -from sage.rings.number_field.number_field import NumberField +from sage.rings.number_field.number_field import NumberField, QuadraticField from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing from sage.rings.rational_field import QQ from sage.rings.real_lazy import CLF, RLF @@ -60,6 +60,9 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule): self._rank = rank self._natural_basis = natural_basis + # TODO: HACK for the charpoly.. needs redesign badly. + self._basis_normalizers = None + if category is None: category = MagmaticAlgebras(field).FiniteDimensional() category = category.WithBasis().Unital() @@ -224,6 +227,7 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule): return V.span_of_basis(b) + @cached_method def _charpoly_coeff(self, i): """ @@ -234,6 +238,19 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule): store the trace/determinant (a_{r-1} and a_{0} respectively) separate from the entire characteristic polynomial. """ + if self._basis_normalizers is not None: + # Must be a matrix class? + # WARNING/TODO: this whole mess is mis-designed. + n = self.natural_basis_space().nrows() + field = self.base_ring().base_ring() # yeeeeaaaahhh + J = self.__class__(n, field, False) + (_,x,_,_) = J._charpoly_matrix_system() + p = J._charpoly_coeff(i) + # p might be missing some vars, have to substitute "optionally" + pairs = zip(x.base_ring().gens(), self._basis_normalizers) + substitutions = { v: v*c for (v,c) in pairs } + return p.subs(substitutions) + (A_of_x, x, xr, detA) = self._charpoly_matrix_system() R = A_of_x.base_ring() if i >= self.rank(): @@ -416,9 +433,9 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule): True """ - if (not x in self) or (not y in self): - raise TypeError("arguments must live in this algebra") - return x.trace_inner_product(y) + X = x.natural_representation() + Y = y.natural_representation() + return self.natural_inner_product(X,Y) def is_trivial(self): @@ -537,6 +554,20 @@ class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule): return self._natural_basis[0].matrix_space() + @staticmethod + def natural_inner_product(X,Y): + """ + Compute the inner product of two naturally-represented elements. + + For example in the real symmetric matrix EJA, this will compute + the trace inner-product of two n-by-n symmetric matrices. The + default should work for the real cartesian product EJA, the + Jordan spin EJA, and the real symmetric matrices. The others + will have to be overridden. + """ + return (X.conjugate_transpose()*Y).trace() + + @cached_method def one(self): """ @@ -748,7 +779,30 @@ class RealCartesianProductEJA(FiniteDimensionalEuclideanJordanAlgebra): return fdeja.__init__(field, mult_table, rank=n, **kwargs) def inner_product(self, x, y): - return _usual_ip(x,y) + """ + Faster to reimplement than to use natural representations. + + SETUP:: + + sage: from mjo.eja.eja_algebra import RealCartesianProductEJA + + TESTS: + + Ensure that this is the usual inner product for the algebras + over `R^n`:: + + 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: X = x.natural_representation() + sage: Y = y.natural_representation() + sage: x.inner_product(y) == J.natural_inner_product(X,Y) + True + + """ + return x.to_vector().inner_product(y.to_vector()) def random_eja(): @@ -814,7 +868,7 @@ def _real_symmetric_basis(n, field): sage: set_random_seed() sage: n = ZZ.random_element(1,5) - sage: B = _real_symmetric_basis(n, QQbar) + sage: B = _real_symmetric_basis(n, QQ) sage: all( M.is_symmetric() for M in B) True @@ -829,8 +883,6 @@ def _real_symmetric_basis(n, field): Sij = Eij else: Sij = Eij + Eij.transpose() - # Now normalize it. - Sij = Sij / _real_symmetric_matrix_ip(Sij,Sij).sqrt() S.append(Sij) return tuple(S) @@ -853,9 +905,7 @@ def _complex_hermitian_basis(n, field): sage: set_random_seed() sage: n = ZZ.random_element(1,5) - sage: R = PolynomialRing(QQ, 'z') - sage: z = R.gen() - sage: field = NumberField(z**2 - 2, 'sqrt2', embedding=RLF(2).sqrt()) + sage: field = QuadraticField(2, 'sqrt2') sage: B = _complex_hermitian_basis(n, field) sage: all( M.is_symmetric() for M in B) True @@ -879,18 +929,15 @@ def _complex_hermitian_basis(n, field): Sij = _embed_complex_matrix(Eij) S.append(Sij) else: - # Beware, orthogonal but not normalized! The second one - # has a minus because it's conjugated. + # The second one has a minus because it's conjugated. Sij_real = _embed_complex_matrix(Eij + Eij.transpose()) S.append(Sij_real) Sij_imag = _embed_complex_matrix(I*Eij - I*Eij.transpose()) S.append(Sij_imag) - # Normalize these with our inner product before handing them back. - # And since we embedded them, we can drop back to the "field" that - # we started with instead of the complex extension "F". - return tuple( (s / _complex_hermitian_matrix_ip(s,s).sqrt()).change_ring(field) - for s in S ) + # Since we embedded these, we can drop back to the "field" that we + # started with instead of the complex extension "F". + return tuple( s.change_ring(field) for s in S ) @@ -989,9 +1036,7 @@ def _embed_complex_matrix(M): EXAMPLES:: - sage: R = PolynomialRing(QQ, 'z') - sage: z = R.gen() - sage: F = NumberField(z**2 + 1, 'i', embedding=CLF(-1).sqrt()) + sage: F = QuadraticField(-1, 'i') sage: x1 = F(4 - 2*i) sage: x2 = F(1 + 2*i) sage: x3 = F(-i) @@ -1010,9 +1055,7 @@ def _embed_complex_matrix(M): sage: set_random_seed() sage: n = ZZ.random_element(5) - sage: R = PolynomialRing(QQ, 'z') - sage: z = R.gen() - sage: F = NumberField(z**2 + 1, 'i', embedding=CLF(-1).sqrt()) + sage: F = QuadraticField(-1, 'i') sage: X = random_matrix(F, n) sage: Y = random_matrix(F, n) sage: actual = _embed_complex_matrix(X) * _embed_complex_matrix(Y) @@ -1059,9 +1102,7 @@ def _unembed_complex_matrix(M): Unembedding is the inverse of embedding:: sage: set_random_seed() - sage: R = PolynomialRing(QQ, 'z') - sage: z = R.gen() - sage: F = NumberField(z**2 + 1, 'i', embedding=CLF(-1).sqrt()) + sage: F = QuadraticField(-1, 'i') sage: M = random_matrix(F, 3) sage: _unembed_complex_matrix(_embed_complex_matrix(M)) == M True @@ -1137,9 +1178,7 @@ def _embed_quaternion_matrix(M): if M.ncols() != n: raise ValueError("the matrix 'M' must be square") - R = PolynomialRing(QQ, 'z') - z = R.gen() - F = NumberField(z**2 + 1, 'i', embedding=CLF(-1).sqrt()) + F = QuadraticField(-1, 'i') i = F.gen() blocks = [] @@ -1193,7 +1232,10 @@ def _unembed_quaternion_matrix(M): if not n.mod(4).is_zero(): raise ValueError("the matrix 'M' must be a complex embedding") - Q = QuaternionAlgebra(QQ,-1,-1) + # Use the base ring of the matrix to ensure that its entries can be + # multiplied by elements of the quaternion algebra. + field = M.base_ring() + Q = QuaternionAlgebra(field,-1,-1) i,j,k = Q.gens() # Go top-left to bottom-right (reading order), converting every @@ -1207,17 +1249,15 @@ def _unembed_quaternion_matrix(M): raise ValueError('bad on-diagonal submatrix') if submat[0,1] != -submat[1,0].conjugate(): raise ValueError('bad off-diagonal submatrix') - z = submat[0,0].real() + submat[0,0].imag()*i - z += submat[0,1].real()*j + submat[0,1].imag()*k + 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 elements.append(z) return matrix(Q, n/4, elements) -# The usual inner product on R^n. -def _usual_ip(x,y): - return x.to_vector().inner_product(y.to_vector()) - # The inner product used for the real symmetric simple EJA. # We keep it as a separate function because e.g. the complex # algebra uses the same inner product, except divided by 2. @@ -1226,15 +1266,6 @@ 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() - -def _complex_hermitian_matrix_ip(X,Y): - # This takes EMBEDDED matrices. - Xu = _unembed_complex_matrix(X) - Yu = _unembed_complex_matrix(Y) - # The trace need not be real; consider Xu = (i*I) and Yu = I. - return ((Xu*Yu).trace()).vector()[0] # real part, I guess class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra): """ @@ -1299,7 +1330,8 @@ class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra): sage: (x*y).inner_product(z) == y.inner_product(x*z) True - Our basis is normalized with respect to the natural inner product:: + Our natural basis is normalized with respect to the natural inner + product unless we specify otherwise:: sage: set_random_seed() sage: n = ZZ.random_element(1,5) @@ -1307,8 +1339,11 @@ class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra): sage: all( b.norm() == 1 for b in J.gens() ) True - Left-multiplication operators are symmetric because they satisfy - the Jordan axiom:: + Since our natural basis is normalized with respect to the natural + inner product, and since we know that this algebra is an EJA, any + left-multiplication operator's matrix will be symmetric because + natural->EJA basis representation is an isometry and within the EJA + the operator is self-adjoint by the Jordan axiom:: sage: set_random_seed() sage: n = ZZ.random_element(1,5) @@ -1317,8 +1352,10 @@ class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra): True """ - def __init__(self, n, field=QQ, **kwargs): - if n > 1: + def __init__(self, n, field=QQ, normalize_basis=True, **kwargs): + S = _real_symmetric_basis(n, field) + + if n > 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. @@ -1327,8 +1364,11 @@ class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra): p = z**2 - 2 if p.is_irreducible(): field = NumberField(p, 'sqrt2', embedding=RLF(2).sqrt()) + S = [ s.change_ring(field) for s in S ] + self._basis_normalizers = tuple( + ~(self.natural_inner_product(s,s).sqrt()) for s in S ) + S = tuple( s*c for (s,c) in zip(S,self._basis_normalizers) ) - S = _real_symmetric_basis(n, field) Qs = _multiplication_table_from_matrix_basis(S) fdeja = super(RealSymmetricEJA, self) @@ -1338,10 +1378,6 @@ class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra): natural_basis=S, **kwargs) - def inner_product(self, x, y): - X = x.natural_representation() - Y = y.natural_representation() - return _real_symmetric_matrix_ip(X,Y) class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): @@ -1397,7 +1433,8 @@ class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): sage: (x*y).inner_product(z) == y.inner_product(x*z) True - Our basis is normalized with respect to the natural inner product:: + Our natural basis is normalized with respect to the natural inner + product unless we specify otherwise:: sage: set_random_seed() sage: n = ZZ.random_element(1,4) @@ -1405,8 +1442,11 @@ class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): sage: all( b.norm() == 1 for b in J.gens() ) True - Left-multiplication operators are symmetric because they satisfy - the Jordan axiom:: + Since our natural basis is normalized with respect to the natural + inner product, and since we know that this algebra is an EJA, any + left-multiplication operator's matrix will be symmetric because + natural->EJA basis representation is an isometry and within the EJA + the operator is self-adjoint by the Jordan axiom:: sage: set_random_seed() sage: n = ZZ.random_element(1,5) @@ -1415,8 +1455,10 @@ class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): True """ - def __init__(self, n, field=QQ, **kwargs): - if n > 1: + def __init__(self, n, field=QQ, normalize_basis=True, **kwargs): + S = _complex_hermitian_basis(n, field) + + if n > 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. @@ -1425,8 +1467,11 @@ class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): p = z**2 - 2 if p.is_irreducible(): field = NumberField(p, 'sqrt2', embedding=RLF(2).sqrt()) + S = [ s.change_ring(field) for s in S ] + self._basis_normalizers = tuple( + ~(self.natural_inner_product(s,s).sqrt()) for s in S ) + S = tuple( s*c for (s,c) in zip(S,self._basis_normalizers) ) - S = _complex_hermitian_basis(n, field) Qs = _multiplication_table_from_matrix_basis(S) fdeja = super(ComplexHermitianEJA, self) @@ -1437,10 +1482,13 @@ class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): **kwargs) - def inner_product(self, x, y): - X = x.natural_representation() - Y = y.natural_representation() - return _complex_hermitian_matrix_ip(X,Y) + @staticmethod + def natural_inner_product(X,Y): + Xu = _unembed_complex_matrix(X) + Yu = _unembed_complex_matrix(Y) + # The trace need not be real; consider Xu = (i*I) and Yu = I. + return ((Xu*Yu).trace()).vector()[0] # real part, I guess + class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): @@ -1459,7 +1507,7 @@ class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): The dimension of this algebra is `n^2`:: sage: set_random_seed() - sage: n = ZZ.random_element(1,5) + sage: n = ZZ.random_element(1,4) sage: J = QuaternionHermitianEJA(n) sage: J.dimension() == 2*(n^2) - n True @@ -1467,7 +1515,7 @@ class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): The Jordan multiplication is what we think it is:: sage: set_random_seed() - sage: n = ZZ.random_element(1,5) + sage: n = ZZ.random_element(1,4) sage: J = QuaternionHermitianEJA(n) sage: x = J.random_element() sage: y = J.random_element() @@ -1488,7 +1536,7 @@ class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): Our inner product satisfies the Jordan axiom:: sage: set_random_seed() - sage: n = ZZ.random_element(1,5) + sage: n = ZZ.random_element(1,4) sage: J = QuaternionHermitianEJA(n) sage: x = J.random_element() sage: y = J.random_element() @@ -1496,9 +1544,45 @@ class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): sage: (x*y).inner_product(z) == y.inner_product(x*z) True + Our natural basis is normalized with respect to the natural inner + product unless we specify otherwise:: + + sage: set_random_seed() + sage: n = ZZ.random_element(1,4) + sage: J = QuaternionHermitianEJA(n) + sage: all( b.norm() == 1 for b in J.gens() ) + True + + Since our natural basis is normalized with respect to the natural + inner product, and since we know that this algebra is an EJA, any + left-multiplication operator's matrix will be symmetric because + natural->EJA basis representation is an isometry and within the EJA + the operator is self-adjoint by the Jordan axiom:: + + sage: set_random_seed() + sage: n = ZZ.random_element(1,5) + sage: x = QuaternionHermitianEJA(n).random_element() + sage: x.operator().matrix().is_symmetric() + True + """ - def __init__(self, n, field=QQ, **kwargs): + def __init__(self, n, field=QQ, normalize_basis=True, **kwargs): S = _quaternion_hermitian_basis(n, field) + + if n > 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. + R = PolynomialRing(field, 'z') + z = R.gen() + p = z**2 - 2 + if p.is_irreducible(): + field = NumberField(p, 'sqrt2', embedding=RLF(2).sqrt()) + S = [ s.change_ring(field) for s in S ] + self._basis_normalizers = tuple( + ~(self.natural_inner_product(s,s).sqrt()) for s in S ) + S = tuple( s*c for (s,c) in zip(S,self._basis_normalizers) ) + Qs = _multiplication_table_from_matrix_basis(S) fdeja = super(QuaternionHermitianEJA, self) @@ -1508,17 +1592,16 @@ class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra): natural_basis=S, **kwargs) - def inner_product(self, x, y): - # Since a+bi+cj+dk on the diagonal is represented as - # - # a + bi +cj + dk = [ a b c d] - # [ -b a -d c] - # [ -c d a -b] - # [ -d -c b a], - # - # we'll quadruple-count the "a" entries if we take the trace of - # the embedding. - return _matrix_ip(x,y)/4 + @staticmethod + def natural_inner_product(X,Y): + Xu = _unembed_quaternion_matrix(X) + Yu = _unembed_quaternion_matrix(Y) + # The trace need not be real; consider Xu = (i*I) and Yu = I. + # The result will be a quaternion algebra element, which doesn't + # have a "vector" method, but does have coefficient_tuple() method + # that returns the coefficients of 1, i, j, and k -- in that order. + return ((Xu*Yu).trace()).coefficient_tuple()[0] + class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra): @@ -1594,4 +1677,27 @@ class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra): return fdeja.__init__(field, mult_table, rank=min(n,2), **kwargs) def inner_product(self, x, y): - return _usual_ip(x,y) + """ + Faster to reimplement than to use natural representations. + + SETUP:: + + sage: from mjo.eja.eja_algebra import JordanSpinEJA + + TESTS: + + Ensure that this is the usual inner product for the algebras + over `R^n`:: + + 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: X = x.natural_representation() + sage: Y = y.natural_representation() + sage: x.inner_product(y) == J.natural_inner_product(X,Y) + True + + """ + return x.to_vector().inner_product(y.to_vector())