2 Euclidean Jordan Algebras. These are formally-real Jordan Algebras;
3 specifically those where u^2 + v^2 = 0 implies that u = v = 0. They
4 are used in optimization, and have some additional nice methods beyond
5 what can be supported in a general Jordan Algebra.
8 from itertools
import repeat
10 from sage
.algebras
.quatalg
.quaternion_algebra
import QuaternionAlgebra
11 from sage
.categories
.magmatic_algebras
import MagmaticAlgebras
12 from sage
.combinat
.free_module
import CombinatorialFreeModule
13 from sage
.matrix
.constructor
import matrix
14 from sage
.matrix
.matrix_space
import MatrixSpace
15 from sage
.misc
.cachefunc
import cached_method
16 from sage
.misc
.lazy_import
import lazy_import
17 from sage
.misc
.prandom
import choice
18 from sage
.misc
.table
import table
19 from sage
.modules
.free_module
import FreeModule
, VectorSpace
20 from sage
.rings
.all
import (ZZ
, QQ
, AA
, QQbar
, RR
, RLF
, CLF
,
23 from mjo
.eja
.eja_element
import FiniteDimensionalEuclideanJordanAlgebraElement
24 lazy_import('mjo.eja.eja_subalgebra',
25 'FiniteDimensionalEuclideanJordanSubalgebra')
26 from mjo
.eja
.eja_utils
import _mat2vec
28 class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule
):
30 def _coerce_map_from_base_ring(self
):
32 Disable the map from the base ring into the algebra.
34 Performing a nonsense conversion like this automatically
35 is counterpedagogical. The fallback is to try the usual
36 element constructor, which should also fail.
40 sage: from mjo.eja.eja_algebra import random_eja
44 sage: set_random_seed()
45 sage: J = random_eja()
47 Traceback (most recent call last):
49 ValueError: not a naturally-represented algebra element
65 sage: from mjo.eja.eja_algebra import (JordanSpinEJA, random_eja)
69 By definition, Jordan multiplication commutes::
71 sage: set_random_seed()
72 sage: J = random_eja()
73 sage: x,y = J.random_elements(2)
79 The ``field`` we're given must be real::
81 sage: JordanSpinEJA(2,QQbar)
82 Traceback (most recent call last):
84 ValueError: field is not real
88 if not field
.is_subring(RR
):
89 # Note: this does return true for the real algebraic
90 # field, and any quadratic field where we've specified
92 raise ValueError('field is not real')
95 self
._natural
_basis
= natural_basis
98 category
= MagmaticAlgebras(field
).FiniteDimensional()
99 category
= category
.WithBasis().Unital()
101 fda
= super(FiniteDimensionalEuclideanJordanAlgebra
, self
)
103 range(len(mult_table
)),
106 self
.print_options(bracket
='')
108 # The multiplication table we're given is necessarily in terms
109 # of vectors, because we don't have an algebra yet for
110 # anything to be an element of. However, it's faster in the
111 # long run to have the multiplication table be in terms of
112 # algebra elements. We do this after calling the superclass
113 # constructor so that from_vector() knows what to do.
114 self
._multiplication
_table
= [
115 list(map(lambda x
: self
.from_vector(x
), ls
))
120 def _element_constructor_(self
, elt
):
122 Construct an element of this algebra from its natural
125 This gets called only after the parent element _call_ method
126 fails to find a coercion for the argument.
130 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
132 ....: RealSymmetricEJA)
136 The identity in `S^n` is converted to the identity in the EJA::
138 sage: J = RealSymmetricEJA(3)
139 sage: I = matrix.identity(QQ,3)
140 sage: J(I) == J.one()
143 This skew-symmetric matrix can't be represented in the EJA::
145 sage: J = RealSymmetricEJA(3)
146 sage: A = matrix(QQ,3, lambda i,j: i-j)
148 Traceback (most recent call last):
150 ArithmeticError: vector is not in free module
154 Ensure that we can convert any element of the two non-matrix
155 simple algebras (whose natural representations are their usual
156 vector representations) back and forth faithfully::
158 sage: set_random_seed()
159 sage: J = HadamardEJA.random_instance()
160 sage: x = J.random_element()
161 sage: J(x.to_vector().column()) == x
163 sage: J = JordanSpinEJA.random_instance()
164 sage: x = J.random_element()
165 sage: J(x.to_vector().column()) == x
169 msg
= "not a naturally-represented algebra element"
171 # The superclass implementation of random_element()
172 # needs to be able to coerce "0" into the algebra.
174 elif elt
in self
.base_ring():
175 # Ensure that no base ring -> algebra coercion is performed
176 # by this method. There's some stupidity in sage that would
177 # otherwise propagate to this method; for example, sage thinks
178 # that the integer 3 belongs to the space of 2-by-2 matrices.
179 raise ValueError(msg
)
181 natural_basis
= self
.natural_basis()
182 basis_space
= natural_basis
[0].matrix_space()
183 if elt
not in basis_space
:
184 raise ValueError(msg
)
186 # Thanks for nothing! Matrix spaces aren't vector spaces in
187 # Sage, so we have to figure out its natural-basis coordinates
188 # ourselves. We use the basis space's ring instead of the
189 # element's ring because the basis space might be an algebraic
190 # closure whereas the base ring of the 3-by-3 identity matrix
191 # could be QQ instead of QQbar.
192 V
= VectorSpace(basis_space
.base_ring(), elt
.nrows()*elt
.ncols())
193 W
= V
.span_of_basis( _mat2vec(s
) for s
in natural_basis
)
194 coords
= W
.coordinate_vector(_mat2vec(elt
))
195 return self
.from_vector(coords
)
200 Return a string representation of ``self``.
204 sage: from mjo.eja.eja_algebra import JordanSpinEJA
208 Ensure that it says what we think it says::
210 sage: JordanSpinEJA(2, field=AA)
211 Euclidean Jordan algebra of dimension 2 over Algebraic Real Field
212 sage: JordanSpinEJA(3, field=RDF)
213 Euclidean Jordan algebra of dimension 3 over Real Double Field
216 fmt
= "Euclidean Jordan algebra of dimension {} over {}"
217 return fmt
.format(self
.dimension(), self
.base_ring())
219 def product_on_basis(self
, i
, j
):
220 return self
._multiplication
_table
[i
][j
]
222 def _a_regular_element(self
):
224 Guess a regular element. Needed to compute the basis for our
225 characteristic polynomial coefficients.
229 sage: from mjo.eja.eja_algebra import random_eja
233 Ensure that this hacky method succeeds for every algebra that we
234 know how to construct::
236 sage: set_random_seed()
237 sage: J = random_eja()
238 sage: J._a_regular_element().is_regular()
243 z
= self
.sum( (i
+1)*gs
[i
] for i
in range(len(gs
)) )
244 if not z
.is_regular():
245 raise ValueError("don't know a regular element")
250 def _charpoly_basis_space(self
):
252 Return the vector space spanned by the basis used in our
253 characteristic polynomial coefficients. This is used not only to
254 compute those coefficients, but also any time we need to
255 evaluate the coefficients (like when we compute the trace or
258 z
= self
._a
_regular
_element
()
259 # Don't use the parent vector space directly here in case this
260 # happens to be a subalgebra. In that case, we would be e.g.
261 # two-dimensional but span_of_basis() would expect three
263 V
= VectorSpace(self
.base_ring(), self
.vector_space().dimension())
264 basis
= [ (z
**k
).to_vector() for k
in range(self
.rank()) ]
265 V1
= V
.span_of_basis( basis
)
266 b
= (V1
.basis() + V1
.complement().basis())
267 return V
.span_of_basis(b
)
272 def _charpoly_coeff(self
, i
):
274 Return the coefficient polynomial "a_{i}" of this algebra's
275 general characteristic polynomial.
277 Having this be a separate cached method lets us compute and
278 store the trace/determinant (a_{r-1} and a_{0} respectively)
279 separate from the entire characteristic polynomial.
281 (A_of_x
, x
, xr
, detA
) = self
._charpoly
_matrix
_system
()
282 R
= A_of_x
.base_ring()
287 # Guaranteed by theory
290 # Danger: the in-place modification is done for performance
291 # reasons (reconstructing a matrix with huge polynomial
292 # entries is slow), but I don't know how cached_method works,
293 # so it's highly possible that we're modifying some global
294 # list variable by reference, here. In other words, you
295 # probably shouldn't call this method twice on the same
296 # algebra, at the same time, in two threads
297 Ai_orig
= A_of_x
.column(i
)
298 A_of_x
.set_column(i
,xr
)
299 numerator
= A_of_x
.det()
300 A_of_x
.set_column(i
,Ai_orig
)
302 # We're relying on the theory here to ensure that each a_i is
303 # indeed back in R, and the added negative signs are to make
304 # the whole charpoly expression sum to zero.
305 return R(-numerator
/detA
)
309 def _charpoly_matrix_system(self
):
311 Compute the matrix whose entries A_ij are polynomials in
312 X1,...,XN, the vector ``x`` of variables X1,...,XN, the vector
313 corresponding to `x^r` and the determinent of the matrix A =
314 [A_ij]. In other words, all of the fixed (cachable) data needed
315 to compute the coefficients of the characteristic polynomial.
320 # Turn my vector space into a module so that "vectors" can
321 # have multivatiate polynomial entries.
322 names
= tuple('X' + str(i
) for i
in range(1,n
+1))
323 R
= PolynomialRing(self
.base_ring(), names
)
325 # Using change_ring() on the parent's vector space doesn't work
326 # here because, in a subalgebra, that vector space has a basis
327 # and change_ring() tries to bring the basis along with it. And
328 # that doesn't work unless the new ring is a PID, which it usually
332 # Now let x = (X1,X2,...,Xn) be the vector whose entries are
336 # And figure out the "left multiplication by x" matrix in
339 monomial_matrices
= [ self
.monomial(i
).operator().matrix()
340 for i
in range(n
) ] # don't recompute these!
342 ek
= self
.monomial(k
).to_vector()
344 sum( x
[i
]*(monomial_matrices
[i
]*ek
)
345 for i
in range(n
) ) )
346 Lx
= matrix
.column(R
, lmbx_cols
)
348 # Now we can compute powers of x "symbolically"
349 x_powers
= [self
.one().to_vector(), x
]
350 for d
in range(2, r
+1):
351 x_powers
.append( Lx
*(x_powers
[-1]) )
353 idmat
= matrix
.identity(R
, n
)
355 W
= self
._charpoly
_basis
_space
()
356 W
= W
.change_ring(R
.fraction_field())
358 # Starting with the standard coordinates x = (X1,X2,...,Xn)
359 # and then converting the entries to W-coordinates allows us
360 # to pass in the standard coordinates to the charpoly and get
361 # back the right answer. Specifically, with x = (X1,X2,...,Xn),
364 # W.coordinates(x^2) eval'd at (standard z-coords)
368 # W-coords of (standard coords of x^2 eval'd at std-coords of z)
370 # We want the middle equivalent thing in our matrix, but use
371 # the first equivalent thing instead so that we can pass in
372 # standard coordinates.
373 x_powers
= [ W
.coordinate_vector(xp
) for xp
in x_powers
]
374 l2
= [idmat
.column(k
-1) for k
in range(r
+1, n
+1)]
375 A_of_x
= matrix
.column(R
, n
, (x_powers
[:r
] + l2
))
376 return (A_of_x
, x
, x_powers
[r
], A_of_x
.det())
380 def characteristic_polynomial(self
):
382 Return a characteristic polynomial that works for all elements
385 The resulting polynomial has `n+1` variables, where `n` is the
386 dimension of this algebra. The first `n` variables correspond to
387 the coordinates of an algebra element: when evaluated at the
388 coordinates of an algebra element with respect to a certain
389 basis, the result is a univariate polynomial (in the one
390 remaining variable ``t``), namely the characteristic polynomial
395 sage: from mjo.eja.eja_algebra import JordanSpinEJA, TrivialEJA
399 The characteristic polynomial in the spin algebra is given in
400 Alizadeh, Example 11.11::
402 sage: J = JordanSpinEJA(3)
403 sage: p = J.characteristic_polynomial(); p
404 X1^2 - X2^2 - X3^2 + (-2*t)*X1 + t^2
405 sage: xvec = J.one().to_vector()
409 By definition, the characteristic polynomial is a monic
410 degree-zero polynomial in a rank-zero algebra. Note that
411 Cayley-Hamilton is indeed satisfied since the polynomial
412 ``1`` evaluates to the identity element of the algebra on
415 sage: J = TrivialEJA()
416 sage: J.characteristic_polynomial()
423 # The list of coefficient polynomials a_0, a_1, a_2, ..., a_n.
424 a
= [ self
._charpoly
_coeff
(i
) for i
in range(r
+1) ]
426 # We go to a bit of trouble here to reorder the
427 # indeterminates, so that it's easier to evaluate the
428 # characteristic polynomial at x's coordinates and get back
429 # something in terms of t, which is what we want.
431 S
= PolynomialRing(self
.base_ring(),'t')
433 S
= PolynomialRing(S
, R
.variable_names())
436 return sum( a
[k
]*(t
**k
) for k
in range(len(a
)) )
439 def inner_product(self
, x
, y
):
441 The inner product associated with this Euclidean Jordan algebra.
443 Defaults to the trace inner product, but can be overridden by
444 subclasses if they are sure that the necessary properties are
449 sage: from mjo.eja.eja_algebra import random_eja
453 Our inner product is "associative," which means the following for
454 a symmetric bilinear form::
456 sage: set_random_seed()
457 sage: J = random_eja()
458 sage: x,y,z = J.random_elements(3)
459 sage: (x*y).inner_product(z) == y.inner_product(x*z)
463 X
= x
.natural_representation()
464 Y
= y
.natural_representation()
465 return self
.natural_inner_product(X
,Y
)
468 def is_trivial(self
):
470 Return whether or not this algebra is trivial.
472 A trivial algebra contains only the zero element.
476 sage: from mjo.eja.eja_algebra import (ComplexHermitianEJA,
481 sage: J = ComplexHermitianEJA(3)
487 sage: J = TrivialEJA()
492 return self
.dimension() == 0
495 def multiplication_table(self
):
497 Return a visual representation of this algebra's multiplication
498 table (on basis elements).
502 sage: from mjo.eja.eja_algebra import JordanSpinEJA
506 sage: J = JordanSpinEJA(4)
507 sage: J.multiplication_table()
508 +----++----+----+----+----+
509 | * || e0 | e1 | e2 | e3 |
510 +====++====+====+====+====+
511 | e0 || e0 | e1 | e2 | e3 |
512 +----++----+----+----+----+
513 | e1 || e1 | e0 | 0 | 0 |
514 +----++----+----+----+----+
515 | e2 || e2 | 0 | e0 | 0 |
516 +----++----+----+----+----+
517 | e3 || e3 | 0 | 0 | e0 |
518 +----++----+----+----+----+
521 M
= list(self
._multiplication
_table
) # copy
522 for i
in range(len(M
)):
523 # M had better be "square"
524 M
[i
] = [self
.monomial(i
)] + M
[i
]
525 M
= [["*"] + list(self
.gens())] + M
526 return table(M
, header_row
=True, header_column
=True, frame
=True)
529 def natural_basis(self
):
531 Return a more-natural representation of this algebra's basis.
533 Every finite-dimensional Euclidean Jordan Algebra is a direct
534 sum of five simple algebras, four of which comprise Hermitian
535 matrices. This method returns the original "natural" basis
536 for our underlying vector space. (Typically, the natural basis
537 is used to construct the multiplication table in the first place.)
539 Note that this will always return a matrix. The standard basis
540 in `R^n` will be returned as `n`-by-`1` column matrices.
544 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
545 ....: RealSymmetricEJA)
549 sage: J = RealSymmetricEJA(2)
551 Finite family {0: e0, 1: e1, 2: e2}
552 sage: J.natural_basis()
554 [1 0] [ 0 0.7071067811865475?] [0 0]
555 [0 0], [0.7071067811865475? 0], [0 1]
560 sage: J = JordanSpinEJA(2)
562 Finite family {0: e0, 1: e1}
563 sage: J.natural_basis()
570 if self
._natural
_basis
is None:
571 M
= self
.natural_basis_space()
572 return tuple( M(b
.to_vector()) for b
in self
.basis() )
574 return self
._natural
_basis
577 def natural_basis_space(self
):
579 Return the matrix space in which this algebra's natural basis
582 if self
._natural
_basis
is None or len(self
._natural
_basis
) == 0:
583 return MatrixSpace(self
.base_ring(), self
.dimension(), 1)
585 return self
._natural
_basis
[0].matrix_space()
589 def natural_inner_product(X
,Y
):
591 Compute the inner product of two naturally-represented elements.
593 For example in the real symmetric matrix EJA, this will compute
594 the trace inner-product of two n-by-n symmetric matrices. The
595 default should work for the real cartesian product EJA, the
596 Jordan spin EJA, and the real symmetric matrices. The others
597 will have to be overridden.
599 return (X
.conjugate_transpose()*Y
).trace()
605 Return the unit element of this algebra.
609 sage: from mjo.eja.eja_algebra import (HadamardEJA,
614 sage: J = HadamardEJA(5)
616 e0 + e1 + e2 + e3 + e4
620 The identity element acts like the identity::
622 sage: set_random_seed()
623 sage: J = random_eja()
624 sage: x = J.random_element()
625 sage: J.one()*x == x and x*J.one() == x
628 The matrix of the unit element's operator is the identity::
630 sage: set_random_seed()
631 sage: J = random_eja()
632 sage: actual = J.one().operator().matrix()
633 sage: expected = matrix.identity(J.base_ring(), J.dimension())
634 sage: actual == expected
638 # We can brute-force compute the matrices of the operators
639 # that correspond to the basis elements of this algebra.
640 # If some linear combination of those basis elements is the
641 # algebra identity, then the same linear combination of
642 # their matrices has to be the identity matrix.
644 # Of course, matrices aren't vectors in sage, so we have to
645 # appeal to the "long vectors" isometry.
646 oper_vecs
= [ _mat2vec(g
.operator().matrix()) for g
in self
.gens() ]
648 # Now we use basis linear algebra to find the coefficients,
649 # of the matrices-as-vectors-linear-combination, which should
650 # work for the original algebra basis too.
651 A
= matrix
.column(self
.base_ring(), oper_vecs
)
653 # We used the isometry on the left-hand side already, but we
654 # still need to do it for the right-hand side. Recall that we
655 # wanted something that summed to the identity matrix.
656 b
= _mat2vec( matrix
.identity(self
.base_ring(), self
.dimension()) )
658 # Now if there's an identity element in the algebra, this should work.
659 coeffs
= A
.solve_right(b
)
660 return self
.linear_combination(zip(self
.gens(), coeffs
))
663 def peirce_decomposition(self
, c
):
665 The Peirce decomposition of this algebra relative to the
668 In the future, this can be extended to a complete system of
669 orthogonal idempotents.
673 - ``c`` -- an idempotent of this algebra.
677 A triple (J0, J5, J1) containing two subalgebras and one subspace
680 - ``J0`` -- the algebra on the eigenspace of ``c.operator()``
681 corresponding to the eigenvalue zero.
683 - ``J5`` -- the eigenspace (NOT a subalgebra) of ``c.operator()``
684 corresponding to the eigenvalue one-half.
686 - ``J1`` -- the algebra on the eigenspace of ``c.operator()``
687 corresponding to the eigenvalue one.
689 These are the only possible eigenspaces for that operator, and this
690 algebra is a direct sum of them. The spaces ``J0`` and ``J1`` are
691 orthogonal, and are subalgebras of this algebra with the appropriate
696 sage: from mjo.eja.eja_algebra import random_eja, RealSymmetricEJA
700 The canonical example comes from the symmetric matrices, which
701 decompose into diagonal and off-diagonal parts::
703 sage: J = RealSymmetricEJA(3)
704 sage: C = matrix(QQ, [ [1,0,0],
708 sage: J0,J5,J1 = J.peirce_decomposition(c)
710 Euclidean Jordan algebra of dimension 1...
712 Vector space of degree 6 and dimension 2...
714 Euclidean Jordan algebra of dimension 3...
718 Every algebra decomposes trivially with respect to its identity
721 sage: set_random_seed()
722 sage: J = random_eja()
723 sage: J0,J5,J1 = J.peirce_decomposition(J.one())
724 sage: J0.dimension() == 0 and J5.dimension() == 0
726 sage: J1.superalgebra() == J and J1.dimension() == J.dimension()
729 The identity elements in the two subalgebras are the
730 projections onto their respective subspaces of the
731 superalgebra's identity element::
733 sage: set_random_seed()
734 sage: J = random_eja()
735 sage: x = J.random_element()
736 sage: if not J.is_trivial():
737 ....: while x.is_nilpotent():
738 ....: x = J.random_element()
739 sage: c = x.subalgebra_idempotent()
740 sage: J0,J5,J1 = J.peirce_decomposition(c)
741 sage: J1(c) == J1.one()
743 sage: J0(J.one() - c) == J0.one()
747 if not c
.is_idempotent():
748 raise ValueError("element is not idempotent: %s" % c
)
750 # Default these to what they should be if they turn out to be
751 # trivial, because eigenspaces_left() won't return eigenvalues
752 # corresponding to trivial spaces (e.g. it returns only the
753 # eigenspace corresponding to lambda=1 if you take the
754 # decomposition relative to the identity element).
755 trivial
= FiniteDimensionalEuclideanJordanSubalgebra(self
, ())
756 J0
= trivial
# eigenvalue zero
757 J5
= VectorSpace(self
.base_ring(), 0) # eigenvalue one-half
758 J1
= trivial
# eigenvalue one
760 for (eigval
, eigspace
) in c
.operator().matrix().right_eigenspaces():
761 if eigval
== ~
(self
.base_ring()(2)):
764 gens
= tuple( self
.from_vector(b
) for b
in eigspace
.basis() )
765 subalg
= FiniteDimensionalEuclideanJordanSubalgebra(self
, gens
)
771 raise ValueError("unexpected eigenvalue: %s" % eigval
)
776 def random_elements(self
, count
):
778 Return ``count`` random elements as a tuple.
782 sage: from mjo.eja.eja_algebra import JordanSpinEJA
786 sage: J = JordanSpinEJA(3)
787 sage: x,y,z = J.random_elements(3)
788 sage: all( [ x in J, y in J, z in J ])
790 sage: len( J.random_elements(10) ) == 10
794 return tuple( self
.random_element() for idx
in range(count
) )
797 def _rank_computation1(self
):
799 Compute the rank of this algebra using highly suspicious voodoo.
803 We first compute the basis representation of the operator L_x
804 using polynomial indeterminates are placeholders for the
805 coordinates of "x", which is arbitrary. We then use that
806 matrix to compute the (polynomial) entries of x^0, x^1, ...,
807 x^d,... for increasing values of "d", starting at zero. The
808 idea is that. If we also add "coefficient variables" a_0,
809 a_1,... to the ring, we can form the linear combination
810 a_0*x^0 + ... + a_d*x^d = 0, and ask what dimension the
811 solution space has as an affine variety. When "d" is smaller
812 than the rank, we expect that dimension to be the number of
813 coordinates of "x", since we can set *those* to whatever we
814 want, but linear independence forces the coefficients a_i to
815 be zero. Eventually, when "d" passes the rank, the dimension
816 of the solution space begins to grow, because we can *still*
817 set the coordinates of "x" arbitrarily, but now there are some
818 coefficients that make the sum zero as well. So, when the
819 dimension of the variety jumps, we return the corresponding
820 "d" as the rank of the algebra. This appears to work.
824 sage: from mjo.eja.eja_algebra import (HadamardEJA,
826 ....: RealSymmetricEJA,
827 ....: ComplexHermitianEJA,
828 ....: QuaternionHermitianEJA)
832 sage: J = HadamardEJA(5)
833 sage: J._rank_computation() == J.rank()
835 sage: J = JordanSpinEJA(5)
836 sage: J._rank_computation() == J.rank()
838 sage: J = RealSymmetricEJA(4)
839 sage: J._rank_computation() == J.rank()
841 sage: J = ComplexHermitianEJA(3)
842 sage: J._rank_computation() == J.rank()
844 sage: J = QuaternionHermitianEJA(2)
845 sage: J._rank_computation() == J.rank()
850 var_names
= [ "X" + str(z
) for z
in range(1,n
+1) ]
852 ideal_dim
= len(var_names
)
854 # From a result in my book, these are the entries of the
855 # basis representation of L_x.
856 return sum( vars[d
+k
]*self
.monomial(k
).operator().matrix()[i
,j
]
859 while ideal_dim
== len(var_names
):
860 coeff_names
= [ "a" + str(z
) for z
in range(d
) ]
861 R
= PolynomialRing(self
.base_ring(), coeff_names
+ var_names
)
863 L_x
= matrix(R
, n
, n
, L_x_i_j
)
864 x_powers
= [ vars[k
]*(L_x
**k
)*self
.one().to_vector()
866 eqs
= [ sum(x_powers
[k
][j
] for k
in range(d
)) for j
in range(n
) ]
867 ideal_dim
= R
.ideal(eqs
).dimension()
870 # Subtract one because we increment one too many times, and
871 # subtract another one because "d" is one greater than the
872 # answer anyway; when d=3, we go up to x^2.
875 def _rank_computation2(self
):
877 Instead of using the dimension of an ideal, find the rank of a
878 matrix containing polynomials.
881 var_names
= [ "X" + str(z
) for z
in range(1,n
+1) ]
882 R
= PolynomialRing(self
.base_ring(), var_names
)
886 # From a result in my book, these are the entries of the
887 # basis representation of L_x.
888 return sum( vars[k
]*self
.monomial(k
).operator().matrix()[i
,j
]
891 L_x
= matrix(R
, n
, n
, L_x_i_j
)
892 x_powers
= [ (vars[k
]*(L_x
**k
)*self
.one().to_vector()).row()
895 from sage
.matrix
.constructor
import block_matrix
896 M
= block_matrix(n
,1,x_powers
)
899 def _rank_computation3(self
):
901 Similar to :meth:`_rank_computation2`, but it stops echelonizing
902 as soon as it hits the first zero row.
910 var_names
= [ "X" + str(z
) for z
in range(1,n
+1) ]
911 R
= PolynomialRing(self
.base_ring(), var_names
)
915 # From a result in my book, these are the entries of the
916 # basis representation of L_x.
917 return sum( vars[k
]*self
.monomial(k
).operator().matrix()[i
,j
]
920 L_x
= matrix(R
, n
, n
, L_x_i_j
)
921 x_powers
= [ vars[k
]*(L_x
**k
)*self
.one().to_vector()
930 rows
= M
.rows() + [x_powers
[d
]]
934 if new_rank
== old_rank
:
941 Return the rank of this EJA.
945 The author knows of no algorithm to compute the rank of an EJA
946 where only the multiplication table is known. In lieu of one, we
947 require the rank to be specified when the algebra is created,
948 and simply pass along that number here.
952 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
953 ....: RealSymmetricEJA,
954 ....: ComplexHermitianEJA,
955 ....: QuaternionHermitianEJA,
960 The rank of the Jordan spin algebra is always two::
962 sage: JordanSpinEJA(2).rank()
964 sage: JordanSpinEJA(3).rank()
966 sage: JordanSpinEJA(4).rank()
969 The rank of the `n`-by-`n` Hermitian real, complex, or
970 quaternion matrices is `n`::
972 sage: RealSymmetricEJA(4).rank()
974 sage: ComplexHermitianEJA(3).rank()
976 sage: QuaternionHermitianEJA(2).rank()
981 Ensure that every EJA that we know how to construct has a
982 positive integer rank, unless the algebra is trivial in
983 which case its rank will be zero::
985 sage: set_random_seed()
986 sage: J = random_eja()
990 sage: r > 0 or (r == 0 and J.is_trivial())
997 def vector_space(self
):
999 Return the vector space that underlies this algebra.
1003 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1007 sage: J = RealSymmetricEJA(2)
1008 sage: J.vector_space()
1009 Vector space of dimension 3 over...
1012 return self
.zero().to_vector().parent().ambient_vector_space()
1015 Element
= FiniteDimensionalEuclideanJordanAlgebraElement
1018 class KnownRankEJA(object):
1020 A class for algebras that we actually know we can construct. The
1021 main issue is that, for most of our methods to make sense, we need
1022 to know the rank of our algebra. Thus we can't simply generate a
1023 "random" algebra, or even check that a given basis and product
1024 satisfy the axioms; because even if everything looks OK, we wouldn't
1025 know the rank we need to actuallty build the thing.
1027 Not really a subclass of FDEJA because doing that causes method
1028 resolution errors, e.g.
1030 TypeError: Error when calling the metaclass bases
1031 Cannot create a consistent method resolution
1032 order (MRO) for bases FiniteDimensionalEuclideanJordanAlgebra,
1037 def _max_test_case_size():
1039 Return an integer "size" that is an upper bound on the size of
1040 this algebra when it is used in a random test
1041 case. Unfortunately, the term "size" is quite vague -- when
1042 dealing with `R^n` under either the Hadamard or Jordan spin
1043 product, the "size" refers to the dimension `n`. When dealing
1044 with a matrix algebra (real symmetric or complex/quaternion
1045 Hermitian), it refers to the size of the matrix, which is
1046 far less than the dimension of the underlying vector space.
1048 We default to five in this class, which is safe in `R^n`. The
1049 matrix algebra subclasses (or any class where the "size" is
1050 interpreted to be far less than the dimension) should override
1051 with a smaller number.
1056 def random_instance(cls
, field
=AA
, **kwargs
):
1058 Return a random instance of this type of algebra.
1060 Beware, this will crash for "most instances" because the
1061 constructor below looks wrong.
1063 if cls
is TrivialEJA
:
1064 # The TrivialEJA class doesn't take an "n" argument because
1068 n
= ZZ
.random_element(cls
._max
_test
_case
_size
()) + 1
1069 return cls(n
, field
, **kwargs
)
1072 class HadamardEJA(FiniteDimensionalEuclideanJordanAlgebra
, KnownRankEJA
):
1074 Return the Euclidean Jordan Algebra corresponding to the set
1075 `R^n` under the Hadamard product.
1077 Note: this is nothing more than the Cartesian product of ``n``
1078 copies of the spin algebra. Once Cartesian product algebras
1079 are implemented, this can go.
1083 sage: from mjo.eja.eja_algebra import HadamardEJA
1087 This multiplication table can be verified by hand::
1089 sage: J = HadamardEJA(3)
1090 sage: e0,e1,e2 = J.gens()
1106 We can change the generator prefix::
1108 sage: HadamardEJA(3, prefix='r').gens()
1112 def __init__(self
, n
, field
=AA
, **kwargs
):
1113 V
= VectorSpace(field
, n
)
1114 mult_table
= [ [ V
.gen(i
)*(i
== j
) for j
in range(n
) ]
1117 fdeja
= super(HadamardEJA
, self
)
1118 return fdeja
.__init
__(field
, mult_table
, rank
=n
, **kwargs
)
1120 def inner_product(self
, x
, y
):
1122 Faster to reimplement than to use natural representations.
1126 sage: from mjo.eja.eja_algebra import HadamardEJA
1130 Ensure that this is the usual inner product for the algebras
1133 sage: set_random_seed()
1134 sage: J = HadamardEJA.random_instance()
1135 sage: x,y = J.random_elements(2)
1136 sage: X = x.natural_representation()
1137 sage: Y = y.natural_representation()
1138 sage: x.inner_product(y) == J.natural_inner_product(X,Y)
1142 return x
.to_vector().inner_product(y
.to_vector())
1145 def random_eja(field
=AA
, nontrivial
=False):
1147 Return a "random" finite-dimensional Euclidean Jordan Algebra.
1151 sage: from mjo.eja.eja_algebra import random_eja
1156 Euclidean Jordan algebra of dimension...
1159 eja_classes
= KnownRankEJA
.__subclasses
__()
1161 eja_classes
.remove(TrivialEJA
)
1162 classname
= choice(eja_classes
)
1163 return classname
.random_instance(field
=field
)
1170 class MatrixEuclideanJordanAlgebra(FiniteDimensionalEuclideanJordanAlgebra
):
1172 def _max_test_case_size():
1173 # Play it safe, since this will be squared and the underlying
1174 # field can have dimension 4 (quaternions) too.
1177 def __init__(self
, field
, basis
, rank
, normalize_basis
=True, **kwargs
):
1179 Compared to the superclass constructor, we take a basis instead of
1180 a multiplication table because the latter can be computed in terms
1181 of the former when the product is known (like it is here).
1183 # Used in this class's fast _charpoly_coeff() override.
1184 self
._basis
_normalizers
= None
1186 # We're going to loop through this a few times, so now's a good
1187 # time to ensure that it isn't a generator expression.
1188 basis
= tuple(basis
)
1190 if rank
> 1 and normalize_basis
:
1191 # We'll need sqrt(2) to normalize the basis, and this
1192 # winds up in the multiplication table, so the whole
1193 # algebra needs to be over the field extension.
1194 R
= PolynomialRing(field
, 'z')
1197 if p
.is_irreducible():
1198 field
= field
.extension(p
, 'sqrt2', embedding
=RLF(2).sqrt())
1199 basis
= tuple( s
.change_ring(field
) for s
in basis
)
1200 self
._basis
_normalizers
= tuple(
1201 ~
(self
.natural_inner_product(s
,s
).sqrt()) for s
in basis
)
1202 basis
= tuple(s
*c
for (s
,c
) in zip(basis
,self
._basis
_normalizers
))
1204 Qs
= self
.multiplication_table_from_matrix_basis(basis
)
1206 fdeja
= super(MatrixEuclideanJordanAlgebra
, self
)
1207 return fdeja
.__init
__(field
,
1210 natural_basis
=basis
,
1215 def _charpoly_coeff(self
, i
):
1217 Override the parent method with something that tries to compute
1218 over a faster (non-extension) field.
1220 if self
._basis
_normalizers
is None:
1221 # We didn't normalize, so assume that the basis we started
1222 # with had entries in a nice field.
1223 return super(MatrixEuclideanJordanAlgebra
, self
)._charpoly
_coeff
(i
)
1225 basis
= ( (b
/n
) for (b
,n
) in zip(self
.natural_basis(),
1226 self
._basis
_normalizers
) )
1228 # Do this over the rationals and convert back at the end.
1229 J
= MatrixEuclideanJordanAlgebra(QQ
,
1232 normalize_basis
=False)
1233 (_
,x
,_
,_
) = J
._charpoly
_matrix
_system
()
1234 p
= J
._charpoly
_coeff
(i
)
1235 # p might be missing some vars, have to substitute "optionally"
1236 pairs
= zip(x
.base_ring().gens(), self
._basis
_normalizers
)
1237 substitutions
= { v: v*c for (v,c) in pairs }
1238 result
= p
.subs(substitutions
)
1240 # The result of "subs" can be either a coefficient-ring
1241 # element or a polynomial. Gotta handle both cases.
1243 return self
.base_ring()(result
)
1245 return result
.change_ring(self
.base_ring())
1249 def multiplication_table_from_matrix_basis(basis
):
1251 At least three of the five simple Euclidean Jordan algebras have the
1252 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
1253 multiplication on the right is matrix multiplication. Given a basis
1254 for the underlying matrix space, this function returns a
1255 multiplication table (obtained by looping through the basis
1256 elements) for an algebra of those matrices.
1258 # In S^2, for example, we nominally have four coordinates even
1259 # though the space is of dimension three only. The vector space V
1260 # is supposed to hold the entire long vector, and the subspace W
1261 # of V will be spanned by the vectors that arise from symmetric
1262 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
1263 field
= basis
[0].base_ring()
1264 dimension
= basis
[0].nrows()
1266 V
= VectorSpace(field
, dimension
**2)
1267 W
= V
.span_of_basis( _mat2vec(s
) for s
in basis
)
1269 mult_table
= [[W
.zero() for j
in range(n
)] for i
in range(n
)]
1272 mat_entry
= (basis
[i
]*basis
[j
] + basis
[j
]*basis
[i
])/2
1273 mult_table
[i
][j
] = W
.coordinate_vector(_mat2vec(mat_entry
))
1281 Embed the matrix ``M`` into a space of real matrices.
1283 The matrix ``M`` can have entries in any field at the moment:
1284 the real numbers, complex numbers, or quaternions. And although
1285 they are not a field, we can probably support octonions at some
1286 point, too. This function returns a real matrix that "acts like"
1287 the original with respect to matrix multiplication; i.e.
1289 real_embed(M*N) = real_embed(M)*real_embed(N)
1292 raise NotImplementedError
1296 def real_unembed(M
):
1298 The inverse of :meth:`real_embed`.
1300 raise NotImplementedError
1304 def natural_inner_product(cls
,X
,Y
):
1305 Xu
= cls
.real_unembed(X
)
1306 Yu
= cls
.real_unembed(Y
)
1307 tr
= (Xu
*Yu
).trace()
1310 # It's real already.
1313 # Otherwise, try the thing that works for complex numbers; and
1314 # if that doesn't work, the thing that works for quaternions.
1316 return tr
.vector()[0] # real part, imag part is index 1
1317 except AttributeError:
1318 # A quaternions doesn't have a vector() method, but does
1319 # have coefficient_tuple() method that returns the
1320 # coefficients of 1, i, j, and k -- in that order.
1321 return tr
.coefficient_tuple()[0]
1324 class RealMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra
):
1328 The identity function, for embedding real matrices into real
1334 def real_unembed(M
):
1336 The identity function, for unembedding real matrices from real
1342 class RealSymmetricEJA(RealMatrixEuclideanJordanAlgebra
, KnownRankEJA
):
1344 The rank-n simple EJA consisting of real symmetric n-by-n
1345 matrices, the usual symmetric Jordan product, and the trace inner
1346 product. It has dimension `(n^2 + n)/2` over the reals.
1350 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1354 sage: J = RealSymmetricEJA(2)
1355 sage: e0, e1, e2 = J.gens()
1363 In theory, our "field" can be any subfield of the reals::
1365 sage: RealSymmetricEJA(2, RDF)
1366 Euclidean Jordan algebra of dimension 3 over Real Double Field
1367 sage: RealSymmetricEJA(2, RR)
1368 Euclidean Jordan algebra of dimension 3 over Real Field with
1369 53 bits of precision
1373 The dimension of this algebra is `(n^2 + n) / 2`::
1375 sage: set_random_seed()
1376 sage: n_max = RealSymmetricEJA._max_test_case_size()
1377 sage: n = ZZ.random_element(1, n_max)
1378 sage: J = RealSymmetricEJA(n)
1379 sage: J.dimension() == (n^2 + n)/2
1382 The Jordan multiplication is what we think it is::
1384 sage: set_random_seed()
1385 sage: J = RealSymmetricEJA.random_instance()
1386 sage: x,y = J.random_elements(2)
1387 sage: actual = (x*y).natural_representation()
1388 sage: X = x.natural_representation()
1389 sage: Y = y.natural_representation()
1390 sage: expected = (X*Y + Y*X)/2
1391 sage: actual == expected
1393 sage: J(expected) == x*y
1396 We can change the generator prefix::
1398 sage: RealSymmetricEJA(3, prefix='q').gens()
1399 (q0, q1, q2, q3, q4, q5)
1401 Our natural basis is normalized with respect to the natural inner
1402 product unless we specify otherwise::
1404 sage: set_random_seed()
1405 sage: J = RealSymmetricEJA.random_instance()
1406 sage: all( b.norm() == 1 for b in J.gens() )
1409 Since our natural basis is normalized with respect to the natural
1410 inner product, and since we know that this algebra is an EJA, any
1411 left-multiplication operator's matrix will be symmetric because
1412 natural->EJA basis representation is an isometry and within the EJA
1413 the operator is self-adjoint by the Jordan axiom::
1415 sage: set_random_seed()
1416 sage: x = RealSymmetricEJA.random_instance().random_element()
1417 sage: x.operator().matrix().is_symmetric()
1422 def _denormalized_basis(cls
, n
, field
):
1424 Return a basis for the space of real symmetric n-by-n matrices.
1428 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1432 sage: set_random_seed()
1433 sage: n = ZZ.random_element(1,5)
1434 sage: B = RealSymmetricEJA._denormalized_basis(n,QQ)
1435 sage: all( M.is_symmetric() for M in B)
1439 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1443 for j
in range(i
+1):
1444 Eij
= matrix(field
, n
, lambda k
,l
: k
==i
and l
==j
)
1448 Sij
= Eij
+ Eij
.transpose()
1454 def _max_test_case_size():
1455 return 4 # Dimension 10
1458 def __init__(self
, n
, field
=AA
, **kwargs
):
1459 basis
= self
._denormalized
_basis
(n
, field
)
1460 super(RealSymmetricEJA
, self
).__init
__(field
, basis
, n
, **kwargs
)
1463 class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra
):
1467 Embed the n-by-n complex matrix ``M`` into the space of real
1468 matrices of size 2n-by-2n via the map the sends each entry `z = a +
1469 bi` to the block matrix ``[[a,b],[-b,a]]``.
1473 sage: from mjo.eja.eja_algebra import \
1474 ....: ComplexMatrixEuclideanJordanAlgebra
1478 sage: F = QuadraticField(-1, 'I')
1479 sage: x1 = F(4 - 2*i)
1480 sage: x2 = F(1 + 2*i)
1483 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
1484 sage: ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1493 Embedding is a homomorphism (isomorphism, in fact)::
1495 sage: set_random_seed()
1496 sage: n_max = ComplexMatrixEuclideanJordanAlgebra._max_test_case_size()
1497 sage: n = ZZ.random_element(n_max)
1498 sage: F = QuadraticField(-1, 'I')
1499 sage: X = random_matrix(F, n)
1500 sage: Y = random_matrix(F, n)
1501 sage: Xe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X)
1502 sage: Ye = ComplexMatrixEuclideanJordanAlgebra.real_embed(Y)
1503 sage: XYe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1510 raise ValueError("the matrix 'M' must be square")
1512 # We don't need any adjoined elements...
1513 field
= M
.base_ring().base_ring()
1517 a
= z
.list()[0] # real part, I guess
1518 b
= z
.list()[1] # imag part, I guess
1519 blocks
.append(matrix(field
, 2, [[a
,b
],[-b
,a
]]))
1521 return matrix
.block(field
, n
, blocks
)
1525 def real_unembed(M
):
1527 The inverse of _embed_complex_matrix().
1531 sage: from mjo.eja.eja_algebra import \
1532 ....: ComplexMatrixEuclideanJordanAlgebra
1536 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1537 ....: [-2, 1, -4, 3],
1538 ....: [ 9, 10, 11, 12],
1539 ....: [-10, 9, -12, 11] ])
1540 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(A)
1542 [ 10*I + 9 12*I + 11]
1546 Unembedding is the inverse of embedding::
1548 sage: set_random_seed()
1549 sage: F = QuadraticField(-1, 'I')
1550 sage: M = random_matrix(F, 3)
1551 sage: Me = ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1552 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1558 raise ValueError("the matrix 'M' must be square")
1559 if not n
.mod(2).is_zero():
1560 raise ValueError("the matrix 'M' must be a complex embedding")
1562 # If "M" was normalized, its base ring might have roots
1563 # adjoined and they can stick around after unembedding.
1564 field
= M
.base_ring()
1565 R
= PolynomialRing(field
, 'z')
1568 # Sage doesn't know how to embed AA into QQbar, i.e. how
1569 # to adjoin sqrt(-1) to AA.
1572 F
= field
.extension(z
**2 + 1, 'I', embedding
=CLF(-1).sqrt())
1575 # Go top-left to bottom-right (reading order), converting every
1576 # 2-by-2 block we see to a single complex element.
1578 for k
in range(n
/2):
1579 for j
in range(n
/2):
1580 submat
= M
[2*k
:2*k
+2,2*j
:2*j
+2]
1581 if submat
[0,0] != submat
[1,1]:
1582 raise ValueError('bad on-diagonal submatrix')
1583 if submat
[0,1] != -submat
[1,0]:
1584 raise ValueError('bad off-diagonal submatrix')
1585 z
= submat
[0,0] + submat
[0,1]*i
1588 return matrix(F
, n
/2, elements
)
1592 def natural_inner_product(cls
,X
,Y
):
1594 Compute a natural inner product in this algebra directly from
1599 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1603 This gives the same answer as the slow, default method implemented
1604 in :class:`MatrixEuclideanJordanAlgebra`::
1606 sage: set_random_seed()
1607 sage: J = ComplexHermitianEJA.random_instance()
1608 sage: x,y = J.random_elements(2)
1609 sage: Xe = x.natural_representation()
1610 sage: Ye = y.natural_representation()
1611 sage: X = ComplexHermitianEJA.real_unembed(Xe)
1612 sage: Y = ComplexHermitianEJA.real_unembed(Ye)
1613 sage: expected = (X*Y).trace().real()
1614 sage: actual = ComplexHermitianEJA.natural_inner_product(Xe,Ye)
1615 sage: actual == expected
1619 return RealMatrixEuclideanJordanAlgebra
.natural_inner_product(X
,Y
)/2
1622 class ComplexHermitianEJA(ComplexMatrixEuclideanJordanAlgebra
, KnownRankEJA
):
1624 The rank-n simple EJA consisting of complex Hermitian n-by-n
1625 matrices over the real numbers, the usual symmetric Jordan product,
1626 and the real-part-of-trace inner product. It has dimension `n^2` over
1631 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1635 In theory, our "field" can be any subfield of the reals::
1637 sage: ComplexHermitianEJA(2, RDF)
1638 Euclidean Jordan algebra of dimension 4 over Real Double Field
1639 sage: ComplexHermitianEJA(2, RR)
1640 Euclidean Jordan algebra of dimension 4 over Real Field with
1641 53 bits of precision
1645 The dimension of this algebra is `n^2`::
1647 sage: set_random_seed()
1648 sage: n_max = ComplexHermitianEJA._max_test_case_size()
1649 sage: n = ZZ.random_element(1, n_max)
1650 sage: J = ComplexHermitianEJA(n)
1651 sage: J.dimension() == n^2
1654 The Jordan multiplication is what we think it is::
1656 sage: set_random_seed()
1657 sage: J = ComplexHermitianEJA.random_instance()
1658 sage: x,y = J.random_elements(2)
1659 sage: actual = (x*y).natural_representation()
1660 sage: X = x.natural_representation()
1661 sage: Y = y.natural_representation()
1662 sage: expected = (X*Y + Y*X)/2
1663 sage: actual == expected
1665 sage: J(expected) == x*y
1668 We can change the generator prefix::
1670 sage: ComplexHermitianEJA(2, prefix='z').gens()
1673 Our natural basis is normalized with respect to the natural inner
1674 product unless we specify otherwise::
1676 sage: set_random_seed()
1677 sage: J = ComplexHermitianEJA.random_instance()
1678 sage: all( b.norm() == 1 for b in J.gens() )
1681 Since our natural basis is normalized with respect to the natural
1682 inner product, and since we know that this algebra is an EJA, any
1683 left-multiplication operator's matrix will be symmetric because
1684 natural->EJA basis representation is an isometry and within the EJA
1685 the operator is self-adjoint by the Jordan axiom::
1687 sage: set_random_seed()
1688 sage: x = ComplexHermitianEJA.random_instance().random_element()
1689 sage: x.operator().matrix().is_symmetric()
1695 def _denormalized_basis(cls
, n
, field
):
1697 Returns a basis for the space of complex Hermitian n-by-n matrices.
1699 Why do we embed these? Basically, because all of numerical linear
1700 algebra assumes that you're working with vectors consisting of `n`
1701 entries from a field and scalars from the same field. There's no way
1702 to tell SageMath that (for example) the vectors contain complex
1703 numbers, while the scalar field is real.
1707 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1711 sage: set_random_seed()
1712 sage: n = ZZ.random_element(1,5)
1713 sage: field = QuadraticField(2, 'sqrt2')
1714 sage: B = ComplexHermitianEJA._denormalized_basis(n, field)
1715 sage: all( M.is_symmetric() for M in B)
1719 R
= PolynomialRing(field
, 'z')
1721 F
= field
.extension(z
**2 + 1, 'I')
1724 # This is like the symmetric case, but we need to be careful:
1726 # * We want conjugate-symmetry, not just symmetry.
1727 # * The diagonal will (as a result) be real.
1731 for j
in range(i
+1):
1732 Eij
= matrix(F
, n
, lambda k
,l
: k
==i
and l
==j
)
1734 Sij
= cls
.real_embed(Eij
)
1737 # The second one has a minus because it's conjugated.
1738 Sij_real
= cls
.real_embed(Eij
+ Eij
.transpose())
1740 Sij_imag
= cls
.real_embed(I
*Eij
- I
*Eij
.transpose())
1743 # Since we embedded these, we can drop back to the "field" that we
1744 # started with instead of the complex extension "F".
1745 return ( s
.change_ring(field
) for s
in S
)
1748 def __init__(self
, n
, field
=AA
, **kwargs
):
1749 basis
= self
._denormalized
_basis
(n
,field
)
1750 super(ComplexHermitianEJA
,self
).__init
__(field
, basis
, n
, **kwargs
)
1753 class QuaternionMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra
):
1757 Embed the n-by-n quaternion matrix ``M`` into the space of real
1758 matrices of size 4n-by-4n by first sending each quaternion entry `z
1759 = a + bi + cj + dk` to the block-complex matrix ``[[a + bi,
1760 c+di],[-c + di, a-bi]]`, and then embedding those into a real
1765 sage: from mjo.eja.eja_algebra import \
1766 ....: QuaternionMatrixEuclideanJordanAlgebra
1770 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1771 sage: i,j,k = Q.gens()
1772 sage: x = 1 + 2*i + 3*j + 4*k
1773 sage: M = matrix(Q, 1, [[x]])
1774 sage: QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1780 Embedding is a homomorphism (isomorphism, in fact)::
1782 sage: set_random_seed()
1783 sage: n_max = QuaternionMatrixEuclideanJordanAlgebra._max_test_case_size()
1784 sage: n = ZZ.random_element(n_max)
1785 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1786 sage: X = random_matrix(Q, n)
1787 sage: Y = random_matrix(Q, n)
1788 sage: Xe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X)
1789 sage: Ye = QuaternionMatrixEuclideanJordanAlgebra.real_embed(Y)
1790 sage: XYe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1795 quaternions
= M
.base_ring()
1798 raise ValueError("the matrix 'M' must be square")
1800 F
= QuadraticField(-1, 'I')
1805 t
= z
.coefficient_tuple()
1810 cplxM
= matrix(F
, 2, [[ a
+ b
*i
, c
+ d
*i
],
1811 [-c
+ d
*i
, a
- b
*i
]])
1812 realM
= ComplexMatrixEuclideanJordanAlgebra
.real_embed(cplxM
)
1813 blocks
.append(realM
)
1815 # We should have real entries by now, so use the realest field
1816 # we've got for the return value.
1817 return matrix
.block(quaternions
.base_ring(), n
, blocks
)
1822 def real_unembed(M
):
1824 The inverse of _embed_quaternion_matrix().
1828 sage: from mjo.eja.eja_algebra import \
1829 ....: QuaternionMatrixEuclideanJordanAlgebra
1833 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
1834 ....: [-2, 1, -4, 3],
1835 ....: [-3, 4, 1, -2],
1836 ....: [-4, -3, 2, 1]])
1837 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(M)
1838 [1 + 2*i + 3*j + 4*k]
1842 Unembedding is the inverse of embedding::
1844 sage: set_random_seed()
1845 sage: Q = QuaternionAlgebra(QQ, -1, -1)
1846 sage: M = random_matrix(Q, 3)
1847 sage: Me = QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1848 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1854 raise ValueError("the matrix 'M' must be square")
1855 if not n
.mod(4).is_zero():
1856 raise ValueError("the matrix 'M' must be a quaternion embedding")
1858 # Use the base ring of the matrix to ensure that its entries can be
1859 # multiplied by elements of the quaternion algebra.
1860 field
= M
.base_ring()
1861 Q
= QuaternionAlgebra(field
,-1,-1)
1864 # Go top-left to bottom-right (reading order), converting every
1865 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
1868 for l
in range(n
/4):
1869 for m
in range(n
/4):
1870 submat
= ComplexMatrixEuclideanJordanAlgebra
.real_unembed(
1871 M
[4*l
:4*l
+4,4*m
:4*m
+4] )
1872 if submat
[0,0] != submat
[1,1].conjugate():
1873 raise ValueError('bad on-diagonal submatrix')
1874 if submat
[0,1] != -submat
[1,0].conjugate():
1875 raise ValueError('bad off-diagonal submatrix')
1876 z
= submat
[0,0].real()
1877 z
+= submat
[0,0].imag()*i
1878 z
+= submat
[0,1].real()*j
1879 z
+= submat
[0,1].imag()*k
1882 return matrix(Q
, n
/4, elements
)
1886 def natural_inner_product(cls
,X
,Y
):
1888 Compute a natural inner product in this algebra directly from
1893 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1897 This gives the same answer as the slow, default method implemented
1898 in :class:`MatrixEuclideanJordanAlgebra`::
1900 sage: set_random_seed()
1901 sage: J = QuaternionHermitianEJA.random_instance()
1902 sage: x,y = J.random_elements(2)
1903 sage: Xe = x.natural_representation()
1904 sage: Ye = y.natural_representation()
1905 sage: X = QuaternionHermitianEJA.real_unembed(Xe)
1906 sage: Y = QuaternionHermitianEJA.real_unembed(Ye)
1907 sage: expected = (X*Y).trace().coefficient_tuple()[0]
1908 sage: actual = QuaternionHermitianEJA.natural_inner_product(Xe,Ye)
1909 sage: actual == expected
1913 return RealMatrixEuclideanJordanAlgebra
.natural_inner_product(X
,Y
)/4
1916 class QuaternionHermitianEJA(QuaternionMatrixEuclideanJordanAlgebra
,
1919 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
1920 matrices, the usual symmetric Jordan product, and the
1921 real-part-of-trace inner product. It has dimension `2n^2 - n` over
1926 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1930 In theory, our "field" can be any subfield of the reals::
1932 sage: QuaternionHermitianEJA(2, RDF)
1933 Euclidean Jordan algebra of dimension 6 over Real Double Field
1934 sage: QuaternionHermitianEJA(2, RR)
1935 Euclidean Jordan algebra of dimension 6 over Real Field with
1936 53 bits of precision
1940 The dimension of this algebra is `2*n^2 - n`::
1942 sage: set_random_seed()
1943 sage: n_max = QuaternionHermitianEJA._max_test_case_size()
1944 sage: n = ZZ.random_element(1, n_max)
1945 sage: J = QuaternionHermitianEJA(n)
1946 sage: J.dimension() == 2*(n^2) - n
1949 The Jordan multiplication is what we think it is::
1951 sage: set_random_seed()
1952 sage: J = QuaternionHermitianEJA.random_instance()
1953 sage: x,y = J.random_elements(2)
1954 sage: actual = (x*y).natural_representation()
1955 sage: X = x.natural_representation()
1956 sage: Y = y.natural_representation()
1957 sage: expected = (X*Y + Y*X)/2
1958 sage: actual == expected
1960 sage: J(expected) == x*y
1963 We can change the generator prefix::
1965 sage: QuaternionHermitianEJA(2, prefix='a').gens()
1966 (a0, a1, a2, a3, a4, a5)
1968 Our natural basis is normalized with respect to the natural inner
1969 product unless we specify otherwise::
1971 sage: set_random_seed()
1972 sage: J = QuaternionHermitianEJA.random_instance()
1973 sage: all( b.norm() == 1 for b in J.gens() )
1976 Since our natural basis is normalized with respect to the natural
1977 inner product, and since we know that this algebra is an EJA, any
1978 left-multiplication operator's matrix will be symmetric because
1979 natural->EJA basis representation is an isometry and within the EJA
1980 the operator is self-adjoint by the Jordan axiom::
1982 sage: set_random_seed()
1983 sage: x = QuaternionHermitianEJA.random_instance().random_element()
1984 sage: x.operator().matrix().is_symmetric()
1989 def _denormalized_basis(cls
, n
, field
):
1991 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
1993 Why do we embed these? Basically, because all of numerical
1994 linear algebra assumes that you're working with vectors consisting
1995 of `n` entries from a field and scalars from the same field. There's
1996 no way to tell SageMath that (for example) the vectors contain
1997 complex numbers, while the scalar field is real.
2001 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
2005 sage: set_random_seed()
2006 sage: n = ZZ.random_element(1,5)
2007 sage: B = QuaternionHermitianEJA._denormalized_basis(n,QQ)
2008 sage: all( M.is_symmetric() for M in B )
2012 Q
= QuaternionAlgebra(QQ
,-1,-1)
2015 # This is like the symmetric case, but we need to be careful:
2017 # * We want conjugate-symmetry, not just symmetry.
2018 # * The diagonal will (as a result) be real.
2022 for j
in range(i
+1):
2023 Eij
= matrix(Q
, n
, lambda k
,l
: k
==i
and l
==j
)
2025 Sij
= cls
.real_embed(Eij
)
2028 # The second, third, and fourth ones have a minus
2029 # because they're conjugated.
2030 Sij_real
= cls
.real_embed(Eij
+ Eij
.transpose())
2032 Sij_I
= cls
.real_embed(I
*Eij
- I
*Eij
.transpose())
2034 Sij_J
= cls
.real_embed(J
*Eij
- J
*Eij
.transpose())
2036 Sij_K
= cls
.real_embed(K
*Eij
- K
*Eij
.transpose())
2039 # Since we embedded these, we can drop back to the "field" that we
2040 # started with instead of the quaternion algebra "Q".
2041 return ( s
.change_ring(field
) for s
in S
)
2044 def __init__(self
, n
, field
=AA
, **kwargs
):
2045 basis
= self
._denormalized
_basis
(n
,field
)
2046 super(QuaternionHermitianEJA
,self
).__init
__(field
, basis
, n
, **kwargs
)
2049 class BilinearFormEJA(FiniteDimensionalEuclideanJordanAlgebra
, KnownRankEJA
):
2051 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2052 with the half-trace inner product and jordan product ``x*y =
2053 (x0*y0 + <B*x_bar,y_bar>, x0*y_bar + y0*x_bar)`` where ``B`` is a
2054 symmetric positive-definite "bilinear form" matrix. It has
2055 dimension `n` over the reals, and reduces to the ``JordanSpinEJA``
2056 when ``B`` is the identity matrix of order ``n-1``.
2060 sage: from mjo.eja.eja_algebra import (BilinearFormEJA,
2061 ....: JordanSpinEJA)
2065 When no bilinear form is specified, the identity matrix is used,
2066 and the resulting algebra is the Jordan spin algebra::
2068 sage: J0 = BilinearFormEJA(3)
2069 sage: J1 = JordanSpinEJA(3)
2070 sage: J0.multiplication_table() == J0.multiplication_table()
2075 We can create a zero-dimensional algebra::
2077 sage: J = BilinearFormEJA(0)
2081 We can check the multiplication condition given in the Jordan, von
2082 Neumann, and Wigner paper (and also discussed on my "On the
2083 symmetry..." paper). Note that this relies heavily on the standard
2084 choice of basis, as does anything utilizing the bilinear form matrix::
2086 sage: set_random_seed()
2087 sage: n = ZZ.random_element(5)
2088 sage: M = matrix.random(QQ, max(0,n-1), algorithm='unimodular')
2089 sage: B = M.transpose()*M
2090 sage: J = BilinearFormEJA(n, B=B)
2091 sage: eis = VectorSpace(M.base_ring(), M.ncols()).basis()
2092 sage: V = J.vector_space()
2093 sage: sis = [ J.from_vector(V([0] + (M.inverse()*ei).list()))
2094 ....: for ei in eis ]
2095 sage: actual = [ sis[i]*sis[j]
2096 ....: for i in range(n-1)
2097 ....: for j in range(n-1) ]
2098 sage: expected = [ J.one() if i == j else J.zero()
2099 ....: for i in range(n-1)
2100 ....: for j in range(n-1) ]
2101 sage: actual == expected
2104 def __init__(self
, n
, field
=AA
, B
=None, **kwargs
):
2106 self
._B
= matrix
.identity(field
, max(0,n
-1))
2110 V
= VectorSpace(field
, n
)
2111 mult_table
= [[V
.zero() for j
in range(n
)] for i
in range(n
)]
2120 z0
= x0
*y0
+ (self
._B
*xbar
).inner_product(ybar
)
2121 zbar
= y0
*xbar
+ x0
*ybar
2122 z
= V([z0
] + zbar
.list())
2123 mult_table
[i
][j
] = z
2125 # The rank of this algebra is two, unless we're in a
2126 # one-dimensional ambient space (because the rank is bounded
2127 # by the ambient dimension).
2128 fdeja
= super(BilinearFormEJA
, self
)
2129 return fdeja
.__init
__(field
, mult_table
, rank
=min(n
,2), **kwargs
)
2131 def inner_product(self
, x
, y
):
2133 Half of the trace inner product.
2135 This is defined so that the special case of the Jordan spin
2136 algebra gets the usual inner product.
2140 sage: from mjo.eja.eja_algebra import BilinearFormEJA
2144 Ensure that this is one-half of the trace inner-product when
2145 the algebra isn't just the reals (when ``n`` isn't one). This
2146 is in Faraut and Koranyi, and also my "On the symmetry..."
2149 sage: set_random_seed()
2150 sage: n = ZZ.random_element(2,5)
2151 sage: M = matrix.random(QQ, max(0,n-1), algorithm='unimodular')
2152 sage: B = M.transpose()*M
2153 sage: J = BilinearFormEJA(n, B=B)
2154 sage: x = J.random_element()
2155 sage: y = J.random_element()
2156 sage: x.inner_product(y) == (x*y).trace()/2
2160 xvec
= x
.to_vector()
2162 yvec
= y
.to_vector()
2164 return x
[0]*y
[0] + (self
._B
*xbar
).inner_product(ybar
)
2167 class JordanSpinEJA(BilinearFormEJA
):
2169 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2170 with the usual inner product and jordan product ``x*y =
2171 (<x,y>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
2176 sage: from mjo.eja.eja_algebra import JordanSpinEJA
2180 This multiplication table can be verified by hand::
2182 sage: J = JordanSpinEJA(4)
2183 sage: e0,e1,e2,e3 = J.gens()
2199 We can change the generator prefix::
2201 sage: JordanSpinEJA(2, prefix='B').gens()
2206 Ensure that we have the usual inner product on `R^n`::
2208 sage: set_random_seed()
2209 sage: J = JordanSpinEJA.random_instance()
2210 sage: x,y = J.random_elements(2)
2211 sage: X = x.natural_representation()
2212 sage: Y = y.natural_representation()
2213 sage: x.inner_product(y) == J.natural_inner_product(X,Y)
2217 def __init__(self
, n
, field
=AA
, **kwargs
):
2218 # This is a special case of the BilinearFormEJA with the identity
2219 # matrix as its bilinear form.
2220 return super(JordanSpinEJA
, self
).__init
__(n
, field
, **kwargs
)
2223 class TrivialEJA(FiniteDimensionalEuclideanJordanAlgebra
, KnownRankEJA
):
2225 The trivial Euclidean Jordan algebra consisting of only a zero element.
2229 sage: from mjo.eja.eja_algebra import TrivialEJA
2233 sage: J = TrivialEJA()
2240 sage: 7*J.one()*12*J.one()
2242 sage: J.one().inner_product(J.one())
2244 sage: J.one().norm()
2246 sage: J.one().subalgebra_generated_by()
2247 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
2252 def __init__(self
, field
=AA
, **kwargs
):
2254 fdeja
= super(TrivialEJA
, self
)
2255 # The rank is zero using my definition, namely the dimension of the
2256 # largest subalgebra generated by any element.
2257 return fdeja
.__init
__(field
, mult_table
, rank
=0, **kwargs
)