]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: fix the Hadamard/JordanSpinEJA fast path for charpoly_coefficients.
[sage.d.git] / mjo / eja / eja_algebra.py
1 """
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.
6
7
8 SETUP::
9
10 sage: from mjo.eja.eja_algebra import random_eja
11
12 EXAMPLES::
13
14 sage: random_eja()
15 Euclidean Jordan algebra of dimension...
16
17 """
18
19 from itertools import repeat
20
21 from sage.algebras.quatalg.quaternion_algebra import QuaternionAlgebra
22 from sage.categories.magmatic_algebras import MagmaticAlgebras
23 from sage.combinat.free_module import CombinatorialFreeModule
24 from sage.matrix.constructor import matrix
25 from sage.matrix.matrix_space import MatrixSpace
26 from sage.misc.cachefunc import cached_method
27 from sage.misc.table import table
28 from sage.modules.free_module import FreeModule, VectorSpace
29 from sage.rings.all import (ZZ, QQ, AA, QQbar, RR, RLF, CLF,
30 PolynomialRing,
31 QuadraticField)
32 from mjo.eja.eja_element import FiniteDimensionalEuclideanJordanAlgebraElement
33 from mjo.eja.eja_operator import FiniteDimensionalEuclideanJordanAlgebraOperator
34 from mjo.eja.eja_utils import _mat2vec
35
36 class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule):
37 r"""
38 The lowest-level class for representing a Euclidean Jordan algebra.
39 """
40 def _coerce_map_from_base_ring(self):
41 """
42 Disable the map from the base ring into the algebra.
43
44 Performing a nonsense conversion like this automatically
45 is counterpedagogical. The fallback is to try the usual
46 element constructor, which should also fail.
47
48 SETUP::
49
50 sage: from mjo.eja.eja_algebra import random_eja
51
52 TESTS::
53
54 sage: set_random_seed()
55 sage: J = random_eja()
56 sage: J(1)
57 Traceback (most recent call last):
58 ...
59 ValueError: not an element of this algebra
60
61 """
62 return None
63
64 def __init__(self,
65 field,
66 multiplication_table,
67 inner_product_table,
68 prefix='e',
69 category=None,
70 matrix_basis=None,
71 check_field=True,
72 check_axioms=True):
73 """
74 INPUT:
75
76 * field -- the scalar field for this algebra (must be real)
77
78 * multiplication_table -- the multiplication table for this
79 algebra's implicit basis. Only the lower-triangular portion
80 of the table is used, since the multiplication is assumed
81 to be commutative.
82
83 SETUP::
84
85 sage: from mjo.eja.eja_algebra import (
86 ....: FiniteDimensionalEuclideanJordanAlgebra,
87 ....: JordanSpinEJA,
88 ....: random_eja)
89
90 EXAMPLES:
91
92 By definition, Jordan multiplication commutes::
93
94 sage: set_random_seed()
95 sage: J = random_eja()
96 sage: x,y = J.random_elements(2)
97 sage: x*y == y*x
98 True
99
100 An error is raised if the Jordan product is not commutative::
101
102 sage: JP = ((1,2),(0,0))
103 sage: IP = ((1,0),(0,1))
104 sage: FiniteDimensionalEuclideanJordanAlgebra(QQ,JP,IP)
105 Traceback (most recent call last):
106 ...
107 ValueError: Jordan product is not commutative
108
109 An error is raised if the inner-product is not commutative::
110
111 sage: JP = ((1,0),(0,1))
112 sage: IP = ((1,2),(0,0))
113 sage: FiniteDimensionalEuclideanJordanAlgebra(QQ,JP,IP)
114 Traceback (most recent call last):
115 ...
116 ValueError: inner-product is not commutative
117
118 TESTS:
119
120 The ``field`` we're given must be real with ``check_field=True``::
121
122 sage: JordanSpinEJA(2, field=QQbar)
123 Traceback (most recent call last):
124 ...
125 ValueError: scalar field is not real
126 sage: JordanSpinEJA(2, field=QQbar, check_field=False)
127 Euclidean Jordan algebra of dimension 2 over Algebraic Field
128
129 The multiplication table must be square with ``check_axioms=True``::
130
131 sage: FiniteDimensionalEuclideanJordanAlgebra(QQ,((),()),((1,),))
132 Traceback (most recent call last):
133 ...
134 ValueError: multiplication table is not square
135
136 The multiplication and inner-product tables must be the same
137 size (and in particular, the inner-product table must also be
138 square) with ``check_axioms=True``::
139
140 sage: FiniteDimensionalEuclideanJordanAlgebra(QQ,((1,),),(()))
141 Traceback (most recent call last):
142 ...
143 ValueError: multiplication and inner-product tables are
144 different sizes
145 sage: FiniteDimensionalEuclideanJordanAlgebra(QQ,((1,),),((1,2),))
146 Traceback (most recent call last):
147 ...
148 ValueError: multiplication and inner-product tables are
149 different sizes
150
151 """
152 if check_field:
153 if not field.is_subring(RR):
154 # Note: this does return true for the real algebraic
155 # field, the rationals, and any quadratic field where
156 # we've specified a real embedding.
157 raise ValueError("scalar field is not real")
158
159
160 # The multiplication and inner-product tables should be square
161 # if the user wants us to verify them. And we verify them as
162 # soon as possible, because we want to exploit their symmetry.
163 n = len(multiplication_table)
164 if check_axioms:
165 if not all( len(l) == n for l in multiplication_table ):
166 raise ValueError("multiplication table is not square")
167
168 # If the multiplication table is square, we can check if
169 # the inner-product table is square by comparing it to the
170 # multiplication table's dimensions.
171 msg = "multiplication and inner-product tables are different sizes"
172 if not len(inner_product_table) == n:
173 raise ValueError(msg)
174
175 if not all( len(l) == n for l in inner_product_table ):
176 raise ValueError(msg)
177
178 # Check commutativity of the Jordan product (symmetry of
179 # the multiplication table) and the commutativity of the
180 # inner-product (symmetry of the inner-product table)
181 # first if we're going to check them at all.. This has to
182 # be done before we define product_on_basis(), because
183 # that method assumes that self._multiplication_table is
184 # symmetric. And it has to be done before we build
185 # self._inner_product_matrix, because the process used to
186 # construct it assumes symmetry as well.
187 if not all( multiplication_table[j][i]
188 == multiplication_table[i][j]
189 for i in range(n)
190 for j in range(i+1) ):
191 raise ValueError("Jordan product is not commutative")
192
193 if not all( inner_product_table[j][i]
194 == inner_product_table[i][j]
195 for i in range(n)
196 for j in range(i+1) ):
197 raise ValueError("inner-product is not commutative")
198
199 self._matrix_basis = matrix_basis
200
201 if category is None:
202 category = MagmaticAlgebras(field).FiniteDimensional()
203 category = category.WithBasis().Unital()
204
205 fda = super(FiniteDimensionalEuclideanJordanAlgebra, self)
206 fda.__init__(field,
207 range(n),
208 prefix=prefix,
209 category=category)
210 self.print_options(bracket='')
211
212 # The multiplication table we're given is necessarily in terms
213 # of vectors, because we don't have an algebra yet for
214 # anything to be an element of. However, it's faster in the
215 # long run to have the multiplication table be in terms of
216 # algebra elements. We do this after calling the superclass
217 # constructor so that from_vector() knows what to do.
218 #
219 # Note: we take advantage of symmetry here, and only store
220 # the lower-triangular portion of the table.
221 self._multiplication_table = [ [ self.vector_space().zero()
222 for j in range(i+1) ]
223 for i in range(n) ]
224
225 for i in range(n):
226 for j in range(i+1):
227 elt = self.from_vector(multiplication_table[i][j])
228 self._multiplication_table[i][j] = elt
229
230 self._multiplication_table = tuple(map(tuple, self._multiplication_table))
231
232 # Save our inner product as a matrix, since the efficiency of
233 # matrix multiplication will usually outweigh the fact that we
234 # have to store a redundant upper- or lower-triangular part.
235 # Pre-cache the fact that these are Hermitian (real symmetric,
236 # in fact) in case some e.g. matrix multiplication routine can
237 # take advantage of it.
238 ip_matrix_constructor = lambda i,j: inner_product_table[i][j] if j <= i else inner_product_table[j][i]
239 self._inner_product_matrix = matrix(field, n, ip_matrix_constructor)
240 self._inner_product_matrix._cache = {'hermitian': True}
241 self._inner_product_matrix.set_immutable()
242
243 if check_axioms:
244 if not self._is_jordanian():
245 raise ValueError("Jordan identity does not hold")
246 if not self._inner_product_is_associative():
247 raise ValueError("inner product is not associative")
248
249 def _element_constructor_(self, elt):
250 """
251 Construct an element of this algebra from its vector or matrix
252 representation.
253
254 This gets called only after the parent element _call_ method
255 fails to find a coercion for the argument.
256
257 SETUP::
258
259 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
260 ....: HadamardEJA,
261 ....: RealSymmetricEJA)
262
263 EXAMPLES:
264
265 The identity in `S^n` is converted to the identity in the EJA::
266
267 sage: J = RealSymmetricEJA(3)
268 sage: I = matrix.identity(QQ,3)
269 sage: J(I) == J.one()
270 True
271
272 This skew-symmetric matrix can't be represented in the EJA::
273
274 sage: J = RealSymmetricEJA(3)
275 sage: A = matrix(QQ,3, lambda i,j: i-j)
276 sage: J(A)
277 Traceback (most recent call last):
278 ...
279 ValueError: not an element of this algebra
280
281 TESTS:
282
283 Ensure that we can convert any element of the two non-matrix
284 simple algebras (whose matrix representations are columns)
285 back and forth faithfully::
286
287 sage: set_random_seed()
288 sage: J = HadamardEJA.random_instance()
289 sage: x = J.random_element()
290 sage: J(x.to_vector().column()) == x
291 True
292 sage: J = JordanSpinEJA.random_instance()
293 sage: x = J.random_element()
294 sage: J(x.to_vector().column()) == x
295 True
296
297 """
298 msg = "not an element of this algebra"
299 if elt == 0:
300 # The superclass implementation of random_element()
301 # needs to be able to coerce "0" into the algebra.
302 return self.zero()
303 elif elt in self.base_ring():
304 # Ensure that no base ring -> algebra coercion is performed
305 # by this method. There's some stupidity in sage that would
306 # otherwise propagate to this method; for example, sage thinks
307 # that the integer 3 belongs to the space of 2-by-2 matrices.
308 raise ValueError(msg)
309
310 if elt not in self.matrix_space():
311 raise ValueError(msg)
312
313 # Thanks for nothing! Matrix spaces aren't vector spaces in
314 # Sage, so we have to figure out its matrix-basis coordinates
315 # ourselves. We use the basis space's ring instead of the
316 # element's ring because the basis space might be an algebraic
317 # closure whereas the base ring of the 3-by-3 identity matrix
318 # could be QQ instead of QQbar.
319 #
320 # We pass check=False because the matrix basis is "guaranteed"
321 # to be linearly independent... right? Ha ha.
322 V = VectorSpace(self.base_ring(), elt.nrows()*elt.ncols())
323 W = V.span_of_basis( (_mat2vec(s) for s in self.matrix_basis()),
324 check=False)
325
326 try:
327 coords = W.coordinate_vector(_mat2vec(elt))
328 except ArithmeticError: # vector is not in free module
329 raise ValueError(msg)
330
331 return self.from_vector(coords)
332
333 def _repr_(self):
334 """
335 Return a string representation of ``self``.
336
337 SETUP::
338
339 sage: from mjo.eja.eja_algebra import JordanSpinEJA
340
341 TESTS:
342
343 Ensure that it says what we think it says::
344
345 sage: JordanSpinEJA(2, field=AA)
346 Euclidean Jordan algebra of dimension 2 over Algebraic Real Field
347 sage: JordanSpinEJA(3, field=RDF)
348 Euclidean Jordan algebra of dimension 3 over Real Double Field
349
350 """
351 fmt = "Euclidean Jordan algebra of dimension {} over {}"
352 return fmt.format(self.dimension(), self.base_ring())
353
354 def product_on_basis(self, i, j):
355 # We only stored the lower-triangular portion of the
356 # multiplication table.
357 if j <= i:
358 return self._multiplication_table[i][j]
359 else:
360 return self._multiplication_table[j][i]
361
362 def _is_commutative(self):
363 r"""
364 Whether or not this algebra's multiplication table is commutative.
365
366 This method should of course always return ``True``, unless
367 this algebra was constructed with ``check_axioms=False`` and
368 passed an invalid multiplication table.
369 """
370 return all( self.product_on_basis(i,j) == self.product_on_basis(i,j)
371 for i in range(self.dimension())
372 for j in range(self.dimension()) )
373
374 def _is_jordanian(self):
375 r"""
376 Whether or not this algebra's multiplication table respects the
377 Jordan identity `(x^{2})(xy) = x(x^{2}y)`.
378
379 We only check one arrangement of `x` and `y`, so for a
380 ``True`` result to be truly true, you should also check
381 :meth:`_is_commutative`. This method should of course always
382 return ``True``, unless this algebra was constructed with
383 ``check_axioms=False`` and passed an invalid multiplication table.
384 """
385 return all( (self.monomial(i)**2)*(self.monomial(i)*self.monomial(j))
386 ==
387 (self.monomial(i))*((self.monomial(i)**2)*self.monomial(j))
388 for i in range(self.dimension())
389 for j in range(self.dimension()) )
390
391 def _inner_product_is_associative(self):
392 r"""
393 Return whether or not this algebra's inner product `B` is
394 associative; that is, whether or not `B(xy,z) = B(x,yz)`.
395
396 This method should of course always return ``True``, unless
397 this algebra was constructed with ``check_axioms=False`` and
398 passed an invalid multiplication table.
399 """
400
401 # Used to check whether or not something is zero in an inexact
402 # ring. This number is sufficient to allow the construction of
403 # QuaternionHermitianEJA(2, field=RDF) with check_axioms=True.
404 epsilon = 1e-16
405
406 for i in range(self.dimension()):
407 for j in range(self.dimension()):
408 for k in range(self.dimension()):
409 x = self.monomial(i)
410 y = self.monomial(j)
411 z = self.monomial(k)
412 diff = (x*y).inner_product(z) - x.inner_product(y*z)
413
414 if self.base_ring().is_exact():
415 if diff != 0:
416 return False
417 else:
418 if diff.abs() > epsilon:
419 return False
420
421 return True
422
423 @cached_method
424 def characteristic_polynomial_of(self):
425 """
426 Return the algebra's "characteristic polynomial of" function,
427 which is itself a multivariate polynomial that, when evaluated
428 at the coordinates of some algebra element, returns that
429 element's characteristic polynomial.
430
431 The resulting polynomial has `n+1` variables, where `n` is the
432 dimension of this algebra. The first `n` variables correspond to
433 the coordinates of an algebra element: when evaluated at the
434 coordinates of an algebra element with respect to a certain
435 basis, the result is a univariate polynomial (in the one
436 remaining variable ``t``), namely the characteristic polynomial
437 of that element.
438
439 SETUP::
440
441 sage: from mjo.eja.eja_algebra import JordanSpinEJA, TrivialEJA
442
443 EXAMPLES:
444
445 The characteristic polynomial in the spin algebra is given in
446 Alizadeh, Example 11.11::
447
448 sage: J = JordanSpinEJA(3)
449 sage: p = J.characteristic_polynomial_of(); p
450 X1^2 - X2^2 - X3^2 + (-2*t)*X1 + t^2
451 sage: xvec = J.one().to_vector()
452 sage: p(*xvec)
453 t^2 - 2*t + 1
454
455 By definition, the characteristic polynomial is a monic
456 degree-zero polynomial in a rank-zero algebra. Note that
457 Cayley-Hamilton is indeed satisfied since the polynomial
458 ``1`` evaluates to the identity element of the algebra on
459 any argument::
460
461 sage: J = TrivialEJA()
462 sage: J.characteristic_polynomial_of()
463 1
464
465 """
466 r = self.rank()
467 n = self.dimension()
468
469 # The list of coefficient polynomials a_0, a_1, a_2, ..., a_(r-1).
470 a = self._charpoly_coefficients()
471
472 # We go to a bit of trouble here to reorder the
473 # indeterminates, so that it's easier to evaluate the
474 # characteristic polynomial at x's coordinates and get back
475 # something in terms of t, which is what we want.
476 S = PolynomialRing(self.base_ring(),'t')
477 t = S.gen(0)
478 if r > 0:
479 R = a[0].parent()
480 S = PolynomialRing(S, R.variable_names())
481 t = S(t)
482
483 return (t**r + sum( a[k]*(t**k) for k in range(r) ))
484
485 def coordinate_polynomial_ring(self):
486 r"""
487 The multivariate polynomial ring in which this algebra's
488 :meth:`characteristic_polynomial_of` lives.
489
490 SETUP::
491
492 sage: from mjo.eja.eja_algebra import (HadamardEJA,
493 ....: RealSymmetricEJA)
494
495 EXAMPLES::
496
497 sage: J = HadamardEJA(2)
498 sage: J.coordinate_polynomial_ring()
499 Multivariate Polynomial Ring in X1, X2...
500 sage: J = RealSymmetricEJA(3,field=QQ,orthonormalize=False)
501 sage: J.coordinate_polynomial_ring()
502 Multivariate Polynomial Ring in X1, X2, X3, X4, X5, X6...
503
504 """
505 var_names = tuple( "X%d" % z for z in range(1, self.dimension()+1) )
506 return PolynomialRing(self.base_ring(), var_names)
507
508 def inner_product(self, x, y):
509 """
510 The inner product associated with this Euclidean Jordan algebra.
511
512 Defaults to the trace inner product, but can be overridden by
513 subclasses if they are sure that the necessary properties are
514 satisfied.
515
516 SETUP::
517
518 sage: from mjo.eja.eja_algebra import (random_eja,
519 ....: HadamardEJA,
520 ....: BilinearFormEJA)
521
522 EXAMPLES:
523
524 Our inner product is "associative," which means the following for
525 a symmetric bilinear form::
526
527 sage: set_random_seed()
528 sage: J = random_eja()
529 sage: x,y,z = J.random_elements(3)
530 sage: (x*y).inner_product(z) == y.inner_product(x*z)
531 True
532
533 TESTS:
534
535 Ensure that this is the usual inner product for the algebras
536 over `R^n`::
537
538 sage: set_random_seed()
539 sage: J = HadamardEJA.random_instance()
540 sage: x,y = J.random_elements(2)
541 sage: actual = x.inner_product(y)
542 sage: expected = x.to_vector().inner_product(y.to_vector())
543 sage: actual == expected
544 True
545
546 Ensure that this is one-half of the trace inner-product in a
547 BilinearFormEJA that isn't just the reals (when ``n`` isn't
548 one). This is in Faraut and Koranyi, and also my "On the
549 symmetry..." paper::
550
551 sage: set_random_seed()
552 sage: J = BilinearFormEJA.random_instance()
553 sage: n = J.dimension()
554 sage: x = J.random_element()
555 sage: y = J.random_element()
556 sage: (n == 1) or (x.inner_product(y) == (x*y).trace()/2)
557 True
558 """
559 B = self._inner_product_matrix
560 return (B*x.to_vector()).inner_product(y.to_vector())
561
562
563 def is_trivial(self):
564 """
565 Return whether or not this algebra is trivial.
566
567 A trivial algebra contains only the zero element.
568
569 SETUP::
570
571 sage: from mjo.eja.eja_algebra import (ComplexHermitianEJA,
572 ....: TrivialEJA)
573
574 EXAMPLES::
575
576 sage: J = ComplexHermitianEJA(3)
577 sage: J.is_trivial()
578 False
579
580 ::
581
582 sage: J = TrivialEJA()
583 sage: J.is_trivial()
584 True
585
586 """
587 return self.dimension() == 0
588
589
590 def multiplication_table(self):
591 """
592 Return a visual representation of this algebra's multiplication
593 table (on basis elements).
594
595 SETUP::
596
597 sage: from mjo.eja.eja_algebra import JordanSpinEJA
598
599 EXAMPLES::
600
601 sage: J = JordanSpinEJA(4)
602 sage: J.multiplication_table()
603 +----++----+----+----+----+
604 | * || e0 | e1 | e2 | e3 |
605 +====++====+====+====+====+
606 | e0 || e0 | e1 | e2 | e3 |
607 +----++----+----+----+----+
608 | e1 || e1 | e0 | 0 | 0 |
609 +----++----+----+----+----+
610 | e2 || e2 | 0 | e0 | 0 |
611 +----++----+----+----+----+
612 | e3 || e3 | 0 | 0 | e0 |
613 +----++----+----+----+----+
614
615 """
616 n = self.dimension()
617 M = [ [ self.zero() for j in range(n) ]
618 for i in range(n) ]
619 for i in range(n):
620 for j in range(i+1):
621 M[i][j] = self._multiplication_table[i][j]
622 M[j][i] = M[i][j]
623
624 for i in range(n):
625 # Prepend the left "header" column entry Can't do this in
626 # the loop because it messes up the symmetry.
627 M[i] = [self.monomial(i)] + M[i]
628
629 # Prepend the header row.
630 M = [["*"] + list(self.gens())] + M
631 return table(M, header_row=True, header_column=True, frame=True)
632
633
634 def matrix_basis(self):
635 """
636 Return an (often more natural) representation of this algebras
637 basis as an ordered tuple of matrices.
638
639 Every finite-dimensional Euclidean Jordan Algebra is a, up to
640 Jordan isomorphism, a direct sum of five simple
641 algebras---four of which comprise Hermitian matrices. And the
642 last type of algebra can of course be thought of as `n`-by-`1`
643 column matrices (ambiguusly called column vectors) to avoid
644 special cases. As a result, matrices (and column vectors) are
645 a natural representation format for Euclidean Jordan algebra
646 elements.
647
648 But, when we construct an algebra from a basis of matrices,
649 those matrix representations are lost in favor of coordinate
650 vectors *with respect to* that basis. We could eventually
651 convert back if we tried hard enough, but having the original
652 representations handy is valuable enough that we simply store
653 them and return them from this method.
654
655 Why implement this for non-matrix algebras? Avoiding special
656 cases for the :class:`BilinearFormEJA` pays with simplicity in
657 its own right. But mainly, we would like to be able to assume
658 that elements of a :class:`DirectSumEJA` can be displayed
659 nicely, without having to have special classes for direct sums
660 one of whose components was a matrix algebra.
661
662 SETUP::
663
664 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
665 ....: RealSymmetricEJA)
666
667 EXAMPLES::
668
669 sage: J = RealSymmetricEJA(2)
670 sage: J.basis()
671 Finite family {0: e0, 1: e1, 2: e2}
672 sage: J.matrix_basis()
673 (
674 [1 0] [ 0 0.7071067811865475?] [0 0]
675 [0 0], [0.7071067811865475? 0], [0 1]
676 )
677
678 ::
679
680 sage: J = JordanSpinEJA(2)
681 sage: J.basis()
682 Finite family {0: e0, 1: e1}
683 sage: J.matrix_basis()
684 (
685 [1] [0]
686 [0], [1]
687 )
688 """
689 if self._matrix_basis is None:
690 M = self.matrix_space()
691 return tuple( M(b.to_vector()) for b in self.basis() )
692 else:
693 return self._matrix_basis
694
695
696 def matrix_space(self):
697 """
698 Return the matrix space in which this algebra's elements live, if
699 we think of them as matrices (including column vectors of the
700 appropriate size).
701
702 Generally this will be an `n`-by-`1` column-vector space,
703 except when the algebra is trivial. There it's `n`-by-`n`
704 (where `n` is zero), to ensure that two elements of the matrix
705 space (empty matrices) can be multiplied.
706
707 Matrix algebras override this with something more useful.
708 """
709 if self.is_trivial():
710 return MatrixSpace(self.base_ring(), 0)
711 elif self._matrix_basis is None or len(self._matrix_basis) == 0:
712 return MatrixSpace(self.base_ring(), self.dimension(), 1)
713 else:
714 return self._matrix_basis[0].matrix_space()
715
716
717 @cached_method
718 def one(self):
719 """
720 Return the unit element of this algebra.
721
722 SETUP::
723
724 sage: from mjo.eja.eja_algebra import (HadamardEJA,
725 ....: random_eja)
726
727 EXAMPLES::
728
729 sage: J = HadamardEJA(5)
730 sage: J.one()
731 e0 + e1 + e2 + e3 + e4
732
733 TESTS:
734
735 The identity element acts like the identity::
736
737 sage: set_random_seed()
738 sage: J = random_eja()
739 sage: x = J.random_element()
740 sage: J.one()*x == x and x*J.one() == x
741 True
742
743 The matrix of the unit element's operator is the identity::
744
745 sage: set_random_seed()
746 sage: J = random_eja()
747 sage: actual = J.one().operator().matrix()
748 sage: expected = matrix.identity(J.base_ring(), J.dimension())
749 sage: actual == expected
750 True
751
752 Ensure that the cached unit element (often precomputed by
753 hand) agrees with the computed one::
754
755 sage: set_random_seed()
756 sage: J = random_eja()
757 sage: cached = J.one()
758 sage: J.one.clear_cache()
759 sage: J.one() == cached
760 True
761
762 """
763 # We can brute-force compute the matrices of the operators
764 # that correspond to the basis elements of this algebra.
765 # If some linear combination of those basis elements is the
766 # algebra identity, then the same linear combination of
767 # their matrices has to be the identity matrix.
768 #
769 # Of course, matrices aren't vectors in sage, so we have to
770 # appeal to the "long vectors" isometry.
771 oper_vecs = [ _mat2vec(g.operator().matrix()) for g in self.gens() ]
772
773 # Now we use basic linear algebra to find the coefficients,
774 # of the matrices-as-vectors-linear-combination, which should
775 # work for the original algebra basis too.
776 A = matrix(self.base_ring(), oper_vecs)
777
778 # We used the isometry on the left-hand side already, but we
779 # still need to do it for the right-hand side. Recall that we
780 # wanted something that summed to the identity matrix.
781 b = _mat2vec( matrix.identity(self.base_ring(), self.dimension()) )
782
783 # Now if there's an identity element in the algebra, this
784 # should work. We solve on the left to avoid having to
785 # transpose the matrix "A".
786 return self.from_vector(A.solve_left(b))
787
788
789 def peirce_decomposition(self, c):
790 """
791 The Peirce decomposition of this algebra relative to the
792 idempotent ``c``.
793
794 In the future, this can be extended to a complete system of
795 orthogonal idempotents.
796
797 INPUT:
798
799 - ``c`` -- an idempotent of this algebra.
800
801 OUTPUT:
802
803 A triple (J0, J5, J1) containing two subalgebras and one subspace
804 of this algebra,
805
806 - ``J0`` -- the algebra on the eigenspace of ``c.operator()``
807 corresponding to the eigenvalue zero.
808
809 - ``J5`` -- the eigenspace (NOT a subalgebra) of ``c.operator()``
810 corresponding to the eigenvalue one-half.
811
812 - ``J1`` -- the algebra on the eigenspace of ``c.operator()``
813 corresponding to the eigenvalue one.
814
815 These are the only possible eigenspaces for that operator, and this
816 algebra is a direct sum of them. The spaces ``J0`` and ``J1`` are
817 orthogonal, and are subalgebras of this algebra with the appropriate
818 restrictions.
819
820 SETUP::
821
822 sage: from mjo.eja.eja_algebra import random_eja, RealSymmetricEJA
823
824 EXAMPLES:
825
826 The canonical example comes from the symmetric matrices, which
827 decompose into diagonal and off-diagonal parts::
828
829 sage: J = RealSymmetricEJA(3)
830 sage: C = matrix(QQ, [ [1,0,0],
831 ....: [0,1,0],
832 ....: [0,0,0] ])
833 sage: c = J(C)
834 sage: J0,J5,J1 = J.peirce_decomposition(c)
835 sage: J0
836 Euclidean Jordan algebra of dimension 1...
837 sage: J5
838 Vector space of degree 6 and dimension 2...
839 sage: J1
840 Euclidean Jordan algebra of dimension 3...
841 sage: J0.one().to_matrix()
842 [0 0 0]
843 [0 0 0]
844 [0 0 1]
845 sage: orig_df = AA.options.display_format
846 sage: AA.options.display_format = 'radical'
847 sage: J.from_vector(J5.basis()[0]).to_matrix()
848 [ 0 0 1/2*sqrt(2)]
849 [ 0 0 0]
850 [1/2*sqrt(2) 0 0]
851 sage: J.from_vector(J5.basis()[1]).to_matrix()
852 [ 0 0 0]
853 [ 0 0 1/2*sqrt(2)]
854 [ 0 1/2*sqrt(2) 0]
855 sage: AA.options.display_format = orig_df
856 sage: J1.one().to_matrix()
857 [1 0 0]
858 [0 1 0]
859 [0 0 0]
860
861 TESTS:
862
863 Every algebra decomposes trivially with respect to its identity
864 element::
865
866 sage: set_random_seed()
867 sage: J = random_eja()
868 sage: J0,J5,J1 = J.peirce_decomposition(J.one())
869 sage: J0.dimension() == 0 and J5.dimension() == 0
870 True
871 sage: J1.superalgebra() == J and J1.dimension() == J.dimension()
872 True
873
874 The decomposition is into eigenspaces, and its components are
875 therefore necessarily orthogonal. Moreover, the identity
876 elements in the two subalgebras are the projections onto their
877 respective subspaces of the superalgebra's identity element::
878
879 sage: set_random_seed()
880 sage: J = random_eja()
881 sage: x = J.random_element()
882 sage: if not J.is_trivial():
883 ....: while x.is_nilpotent():
884 ....: x = J.random_element()
885 sage: c = x.subalgebra_idempotent()
886 sage: J0,J5,J1 = J.peirce_decomposition(c)
887 sage: ipsum = 0
888 sage: for (w,y,z) in zip(J0.basis(), J5.basis(), J1.basis()):
889 ....: w = w.superalgebra_element()
890 ....: y = J.from_vector(y)
891 ....: z = z.superalgebra_element()
892 ....: ipsum += w.inner_product(y).abs()
893 ....: ipsum += w.inner_product(z).abs()
894 ....: ipsum += y.inner_product(z).abs()
895 sage: ipsum
896 0
897 sage: J1(c) == J1.one()
898 True
899 sage: J0(J.one() - c) == J0.one()
900 True
901
902 """
903 if not c.is_idempotent():
904 raise ValueError("element is not idempotent: %s" % c)
905
906 from mjo.eja.eja_subalgebra import FiniteDimensionalEuclideanJordanSubalgebra
907
908 # Default these to what they should be if they turn out to be
909 # trivial, because eigenspaces_left() won't return eigenvalues
910 # corresponding to trivial spaces (e.g. it returns only the
911 # eigenspace corresponding to lambda=1 if you take the
912 # decomposition relative to the identity element).
913 trivial = FiniteDimensionalEuclideanJordanSubalgebra(self, ())
914 J0 = trivial # eigenvalue zero
915 J5 = VectorSpace(self.base_ring(), 0) # eigenvalue one-half
916 J1 = trivial # eigenvalue one
917
918 for (eigval, eigspace) in c.operator().matrix().right_eigenspaces():
919 if eigval == ~(self.base_ring()(2)):
920 J5 = eigspace
921 else:
922 gens = tuple( self.from_vector(b) for b in eigspace.basis() )
923 subalg = FiniteDimensionalEuclideanJordanSubalgebra(self,
924 gens,
925 check_axioms=False)
926 if eigval == 0:
927 J0 = subalg
928 elif eigval == 1:
929 J1 = subalg
930 else:
931 raise ValueError("unexpected eigenvalue: %s" % eigval)
932
933 return (J0, J5, J1)
934
935
936 def random_element(self, thorough=False):
937 r"""
938 Return a random element of this algebra.
939
940 Our algebra superclass method only returns a linear
941 combination of at most two basis elements. We instead
942 want the vector space "random element" method that
943 returns a more diverse selection.
944
945 INPUT:
946
947 - ``thorough`` -- (boolean; default False) whether or not we
948 should generate irrational coefficients for the random
949 element when our base ring is irrational; this slows the
950 algebra operations to a crawl, but any truly random method
951 should include them
952
953 """
954 # For a general base ring... maybe we can trust this to do the
955 # right thing? Unlikely, but.
956 V = self.vector_space()
957 v = V.random_element()
958
959 if self.base_ring() is AA:
960 # The "random element" method of the algebraic reals is
961 # stupid at the moment, and only returns integers between
962 # -2 and 2, inclusive:
963 #
964 # https://trac.sagemath.org/ticket/30875
965 #
966 # Instead, we implement our own "random vector" method,
967 # and then coerce that into the algebra. We use the vector
968 # space degree here instead of the dimension because a
969 # subalgebra could (for example) be spanned by only two
970 # vectors, each with five coordinates. We need to
971 # generate all five coordinates.
972 if thorough:
973 v *= QQbar.random_element().real()
974 else:
975 v *= QQ.random_element()
976
977 return self.from_vector(V.coordinate_vector(v))
978
979 def random_elements(self, count, thorough=False):
980 """
981 Return ``count`` random elements as a tuple.
982
983 INPUT:
984
985 - ``thorough`` -- (boolean; default False) whether or not we
986 should generate irrational coefficients for the random
987 elements when our base ring is irrational; this slows the
988 algebra operations to a crawl, but any truly random method
989 should include them
990
991 SETUP::
992
993 sage: from mjo.eja.eja_algebra import JordanSpinEJA
994
995 EXAMPLES::
996
997 sage: J = JordanSpinEJA(3)
998 sage: x,y,z = J.random_elements(3)
999 sage: all( [ x in J, y in J, z in J ])
1000 True
1001 sage: len( J.random_elements(10) ) == 10
1002 True
1003
1004 """
1005 return tuple( self.random_element(thorough)
1006 for idx in range(count) )
1007
1008
1009 @cached_method
1010 def _charpoly_coefficients(self):
1011 r"""
1012 The `r` polynomial coefficients of the "characteristic polynomial
1013 of" function.
1014 """
1015 n = self.dimension()
1016 R = self.coordinate_polynomial_ring()
1017 vars = R.gens()
1018 F = R.fraction_field()
1019
1020 def L_x_i_j(i,j):
1021 # From a result in my book, these are the entries of the
1022 # basis representation of L_x.
1023 return sum( vars[k]*self.monomial(k).operator().matrix()[i,j]
1024 for k in range(n) )
1025
1026 L_x = matrix(F, n, n, L_x_i_j)
1027
1028 r = None
1029 if self.rank.is_in_cache():
1030 r = self.rank()
1031 # There's no need to pad the system with redundant
1032 # columns if we *know* they'll be redundant.
1033 n = r
1034
1035 # Compute an extra power in case the rank is equal to
1036 # the dimension (otherwise, we would stop at x^(r-1)).
1037 x_powers = [ (L_x**k)*self.one().to_vector()
1038 for k in range(n+1) ]
1039 A = matrix.column(F, x_powers[:n])
1040 AE = A.extended_echelon_form()
1041 E = AE[:,n:]
1042 A_rref = AE[:,:n]
1043 if r is None:
1044 r = A_rref.rank()
1045 b = x_powers[r]
1046
1047 # The theory says that only the first "r" coefficients are
1048 # nonzero, and they actually live in the original polynomial
1049 # ring and not the fraction field. We negate them because
1050 # in the actual characteristic polynomial, they get moved
1051 # to the other side where x^r lives.
1052 return -A_rref.solve_right(E*b).change_ring(R)[:r]
1053
1054 @cached_method
1055 def rank(self):
1056 r"""
1057 Return the rank of this EJA.
1058
1059 This is a cached method because we know the rank a priori for
1060 all of the algebras we can construct. Thus we can avoid the
1061 expensive ``_charpoly_coefficients()`` call unless we truly
1062 need to compute the whole characteristic polynomial.
1063
1064 SETUP::
1065
1066 sage: from mjo.eja.eja_algebra import (HadamardEJA,
1067 ....: JordanSpinEJA,
1068 ....: RealSymmetricEJA,
1069 ....: ComplexHermitianEJA,
1070 ....: QuaternionHermitianEJA,
1071 ....: random_eja)
1072
1073 EXAMPLES:
1074
1075 The rank of the Jordan spin algebra is always two::
1076
1077 sage: JordanSpinEJA(2).rank()
1078 2
1079 sage: JordanSpinEJA(3).rank()
1080 2
1081 sage: JordanSpinEJA(4).rank()
1082 2
1083
1084 The rank of the `n`-by-`n` Hermitian real, complex, or
1085 quaternion matrices is `n`::
1086
1087 sage: RealSymmetricEJA(4).rank()
1088 4
1089 sage: ComplexHermitianEJA(3).rank()
1090 3
1091 sage: QuaternionHermitianEJA(2).rank()
1092 2
1093
1094 TESTS:
1095
1096 Ensure that every EJA that we know how to construct has a
1097 positive integer rank, unless the algebra is trivial in
1098 which case its rank will be zero::
1099
1100 sage: set_random_seed()
1101 sage: J = random_eja()
1102 sage: r = J.rank()
1103 sage: r in ZZ
1104 True
1105 sage: r > 0 or (r == 0 and J.is_trivial())
1106 True
1107
1108 Ensure that computing the rank actually works, since the ranks
1109 of all simple algebras are known and will be cached by default::
1110
1111 sage: set_random_seed() # long time
1112 sage: J = random_eja() # long time
1113 sage: caches = J.rank() # long time
1114 sage: J.rank.clear_cache() # long time
1115 sage: J.rank() == cached # long time
1116 True
1117
1118 """
1119 return len(self._charpoly_coefficients())
1120
1121
1122 def vector_space(self):
1123 """
1124 Return the vector space that underlies this algebra.
1125
1126 SETUP::
1127
1128 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1129
1130 EXAMPLES::
1131
1132 sage: J = RealSymmetricEJA(2)
1133 sage: J.vector_space()
1134 Vector space of dimension 3 over...
1135
1136 """
1137 return self.zero().to_vector().parent().ambient_vector_space()
1138
1139
1140 Element = FiniteDimensionalEuclideanJordanAlgebraElement
1141
1142 class RationalBasisEuclideanJordanAlgebra(FiniteDimensionalEuclideanJordanAlgebra):
1143 r"""
1144 New class for algebras whose supplied basis elements have all rational entries.
1145
1146 SETUP::
1147
1148 sage: from mjo.eja.eja_algebra import BilinearFormEJA
1149
1150 EXAMPLES:
1151
1152 The supplied basis is orthonormalized by default::
1153
1154 sage: B = matrix(QQ, [[1, 0, 0], [0, 25, -32], [0, -32, 41]])
1155 sage: J = BilinearFormEJA(B)
1156 sage: J.matrix_basis()
1157 (
1158 [1] [ 0] [ 0]
1159 [0] [1/5] [32/5]
1160 [0], [ 0], [ 5]
1161 )
1162
1163 """
1164 def __init__(self,
1165 basis,
1166 jordan_product,
1167 inner_product,
1168 field=AA,
1169 orthonormalize=True,
1170 prefix='e',
1171 category=None,
1172 check_field=True,
1173 check_axioms=True):
1174
1175 if check_field:
1176 # Abuse the check_field parameter to check that the entries of
1177 # out basis (in ambient coordinates) are in the field QQ.
1178 if not all( all(b_i in QQ for b_i in b.list()) for b in basis ):
1179 raise TypeError("basis not rational")
1180
1181 # Temporary(?) hack to ensure that the matrix and vector bases
1182 # are over the same ring.
1183 basis = tuple( b.change_ring(field) for b in basis )
1184
1185 n = len(basis)
1186 vector_basis = basis
1187
1188 from sage.structure.element import is_Matrix
1189 basis_is_matrices = False
1190
1191 degree = 0
1192 if n > 0:
1193 if is_Matrix(basis[0]):
1194 basis_is_matrices = True
1195 from mjo.eja.eja_utils import _vec2mat
1196 vector_basis = tuple( map(_mat2vec,basis) )
1197 degree = basis[0].nrows()**2
1198 else:
1199 degree = basis[0].degree()
1200
1201 V = VectorSpace(field, degree)
1202
1203 # Save a copy of an algebra with the original, rational basis
1204 # and over QQ where computations are fast.
1205 self._rational_algebra = None
1206
1207 if field is not QQ:
1208 # There's no point in constructing the extra algebra if this
1209 # one is already rational.
1210 #
1211 # Note: the same Jordan and inner-products work here,
1212 # because they are necessarily defined with respect to
1213 # ambient coordinates and not any particular basis.
1214 self._rational_algebra = RationalBasisEuclideanJordanAlgebra(
1215 basis,
1216 jordan_product,
1217 inner_product,
1218 field=QQ,
1219 orthonormalize=False,
1220 prefix=prefix,
1221 category=category,
1222 check_field=False,
1223 check_axioms=False)
1224
1225 if orthonormalize:
1226 # Compute the deorthonormalized tables before we orthonormalize
1227 # the given basis. The "check" parameter here guarantees that
1228 # the basis is linearly-independent.
1229 W = V.span_of_basis( vector_basis, check=check_axioms)
1230
1231 # Note: the Jordan and inner-products are defined in terms
1232 # of the ambient basis. It's important that their arguments
1233 # are in ambient coordinates as well.
1234 for i in range(n):
1235 for j in range(i+1):
1236 # given basis w.r.t. ambient coords
1237 q_i = vector_basis[i]
1238 q_j = vector_basis[j]
1239
1240 if basis_is_matrices:
1241 q_i = _vec2mat(q_i)
1242 q_j = _vec2mat(q_j)
1243
1244 elt = jordan_product(q_i, q_j)
1245 ip = inner_product(q_i, q_j)
1246
1247 if basis_is_matrices:
1248 # do another mat2vec because the multiplication
1249 # table is in terms of vectors
1250 elt = _mat2vec(elt)
1251
1252 # We overwrite the name "vector_basis" in a second, but never modify it
1253 # in place, to this effectively makes a copy of it.
1254 deortho_vector_basis = vector_basis
1255 self._deortho_matrix = None
1256
1257 if orthonormalize:
1258 from mjo.eja.eja_utils import gram_schmidt
1259 if basis_is_matrices:
1260 vector_ip = lambda x,y: inner_product(_vec2mat(x), _vec2mat(y))
1261 vector_basis = gram_schmidt(vector_basis, vector_ip)
1262 else:
1263 vector_basis = gram_schmidt(vector_basis, inner_product)
1264
1265 # Normalize the "matrix" basis, too!
1266 basis = vector_basis
1267
1268 if basis_is_matrices:
1269 basis = tuple( map(_vec2mat,basis) )
1270
1271 W = V.span_of_basis( vector_basis, check=check_axioms)
1272
1273 # Now "W" is the vector space of our algebra coordinates. The
1274 # variables "X1", "X2",... refer to the entries of vectors in
1275 # W. Thus to convert back and forth between the orthonormal
1276 # coordinates and the given ones, we need to stick the original
1277 # basis in W.
1278 U = V.span_of_basis( deortho_vector_basis, check=check_axioms)
1279 self._deortho_matrix = matrix( U.coordinate_vector(q)
1280 for q in vector_basis )
1281
1282 # If the superclass constructor is going to verify the
1283 # symmetry of this table, it has better at least be
1284 # square...
1285 if check_axioms:
1286 mult_table = [ [0 for j in range(n)] for i in range(n) ]
1287 ip_table = [ [0 for j in range(n)] for i in range(n) ]
1288 else:
1289 mult_table = [ [0 for j in range(i+1)] for i in range(n) ]
1290 ip_table = [ [0 for j in range(i+1)] for i in range(n) ]
1291
1292 # Note: the Jordan and inner-products are defined in terms
1293 # of the ambient basis. It's important that their arguments
1294 # are in ambient coordinates as well.
1295 for i in range(n):
1296 for j in range(i+1):
1297 # ortho basis w.r.t. ambient coords
1298 q_i = vector_basis[i]
1299 q_j = vector_basis[j]
1300
1301 if basis_is_matrices:
1302 q_i = _vec2mat(q_i)
1303 q_j = _vec2mat(q_j)
1304
1305 elt = jordan_product(q_i, q_j)
1306 ip = inner_product(q_i, q_j)
1307
1308 if basis_is_matrices:
1309 # do another mat2vec because the multiplication
1310 # table is in terms of vectors
1311 elt = _mat2vec(elt)
1312
1313 elt = W.coordinate_vector(elt)
1314 mult_table[i][j] = elt
1315 ip_table[i][j] = ip
1316 if check_axioms:
1317 # The tables are square if we're verifying that they
1318 # are commutative.
1319 mult_table[j][i] = elt
1320 ip_table[j][i] = ip
1321
1322 if basis_is_matrices:
1323 for m in basis:
1324 m.set_immutable()
1325 else:
1326 basis = tuple( x.column() for x in basis )
1327
1328 super().__init__(field,
1329 mult_table,
1330 ip_table,
1331 prefix,
1332 category,
1333 basis, # matrix basis
1334 check_field,
1335 check_axioms)
1336
1337 @cached_method
1338 def _charpoly_coefficients(self):
1339 r"""
1340 SETUP::
1341
1342 sage: from mjo.eja.eja_algebra import (BilinearFormEJA,
1343 ....: JordanSpinEJA)
1344
1345 EXAMPLES:
1346
1347 The base ring of the resulting polynomial coefficients is what
1348 it should be, and not the rationals (unless the algebra was
1349 already over the rationals)::
1350
1351 sage: J = JordanSpinEJA(3)
1352 sage: J._charpoly_coefficients()
1353 (X1^2 - X2^2 - X3^2, -2*X1)
1354 sage: a0 = J._charpoly_coefficients()[0]
1355 sage: J.base_ring()
1356 Algebraic Real Field
1357 sage: a0.base_ring()
1358 Algebraic Real Field
1359
1360 """
1361 if self._rational_algebra is None:
1362 # There's no need to construct *another* algebra over the
1363 # rationals if this one is already over the
1364 # rationals. Likewise, if we never orthonormalized our
1365 # basis, we might as well just use the given one.
1366 superclass = super(RationalBasisEuclideanJordanAlgebra, self)
1367 return superclass._charpoly_coefficients()
1368
1369 # Do the computation over the rationals. The answer will be
1370 # the same, because all we've done is a change of basis.
1371 # Then, change back from QQ to our real base ring
1372 a = ( a_i.change_ring(self.base_ring())
1373 for a_i in self._rational_algebra._charpoly_coefficients() )
1374
1375 # Now convert the coordinate variables back to the
1376 # deorthonormalized ones.
1377 R = self.coordinate_polynomial_ring()
1378 from sage.modules.free_module_element import vector
1379 X = vector(R, R.gens())
1380 BX = self._deortho_matrix*X
1381
1382 subs_dict = { X[i]: BX[i] for i in range(len(X)) }
1383 return tuple( a_i.subs(subs_dict) for a_i in a )
1384
1385 class ConcreteEuclideanJordanAlgebra(RationalBasisEuclideanJordanAlgebra):
1386 r"""
1387 A class for the Euclidean Jordan algebras that we know by name.
1388
1389 These are the Jordan algebras whose basis, multiplication table,
1390 rank, and so on are known a priori. More to the point, they are
1391 the Euclidean Jordan algebras for which we are able to conjure up
1392 a "random instance."
1393
1394 SETUP::
1395
1396 sage: from mjo.eja.eja_algebra import ConcreteEuclideanJordanAlgebra
1397
1398 TESTS:
1399
1400 Our basis is normalized with respect to the algebra's inner
1401 product, unless we specify otherwise::
1402
1403 sage: set_random_seed()
1404 sage: J = ConcreteEuclideanJordanAlgebra.random_instance()
1405 sage: all( b.norm() == 1 for b in J.gens() )
1406 True
1407
1408 Since our basis is orthonormal with respect to the algebra's inner
1409 product, and since we know that this algebra is an EJA, any
1410 left-multiplication operator's matrix will be symmetric because
1411 natural->EJA basis representation is an isometry and within the
1412 EJA the operator is self-adjoint by the Jordan axiom::
1413
1414 sage: set_random_seed()
1415 sage: J = ConcreteEuclideanJordanAlgebra.random_instance()
1416 sage: x = J.random_element()
1417 sage: x.operator().is_self_adjoint()
1418 True
1419 """
1420
1421 @staticmethod
1422 def _max_random_instance_size():
1423 """
1424 Return an integer "size" that is an upper bound on the size of
1425 this algebra when it is used in a random test
1426 case. Unfortunately, the term "size" is ambiguous -- when
1427 dealing with `R^n` under either the Hadamard or Jordan spin
1428 product, the "size" refers to the dimension `n`. When dealing
1429 with a matrix algebra (real symmetric or complex/quaternion
1430 Hermitian), it refers to the size of the matrix, which is far
1431 less than the dimension of the underlying vector space.
1432
1433 This method must be implemented in each subclass.
1434 """
1435 raise NotImplementedError
1436
1437 @classmethod
1438 def random_instance(cls, *args, **kwargs):
1439 """
1440 Return a random instance of this type of algebra.
1441
1442 This method should be implemented in each subclass.
1443 """
1444 from sage.misc.prandom import choice
1445 eja_class = choice(cls.__subclasses__())
1446
1447 # These all bubble up to the RationalBasisEuclideanJordanAlgebra
1448 # superclass constructor, so any (kw)args valid there are also
1449 # valid here.
1450 return eja_class.random_instance(*args, **kwargs)
1451
1452
1453 class MatrixEuclideanJordanAlgebra:
1454 @staticmethod
1455 def dimension_over_reals():
1456 r"""
1457 The dimension of this matrix's base ring over the reals.
1458
1459 The reals are dimension one over themselves, obviously; that's
1460 just `\mathbb{R}^{1}`. Likewise, the complex numbers `a + bi`
1461 have dimension two. Finally, the quaternions have dimension
1462 four over the reals.
1463
1464 This is used to determine the size of the matrix returned from
1465 :meth:`real_embed`, among other things.
1466 """
1467 raise NotImplementedError
1468
1469 @classmethod
1470 def real_embed(cls,M):
1471 """
1472 Embed the matrix ``M`` into a space of real matrices.
1473
1474 The matrix ``M`` can have entries in any field at the moment:
1475 the real numbers, complex numbers, or quaternions. And although
1476 they are not a field, we can probably support octonions at some
1477 point, too. This function returns a real matrix that "acts like"
1478 the original with respect to matrix multiplication; i.e.
1479
1480 real_embed(M*N) = real_embed(M)*real_embed(N)
1481
1482 """
1483 if M.ncols() != M.nrows():
1484 raise ValueError("the matrix 'M' must be square")
1485 return M
1486
1487
1488 @classmethod
1489 def real_unembed(cls,M):
1490 """
1491 The inverse of :meth:`real_embed`.
1492 """
1493 if M.ncols() != M.nrows():
1494 raise ValueError("the matrix 'M' must be square")
1495 if not ZZ(M.nrows()).mod(cls.dimension_over_reals()).is_zero():
1496 raise ValueError("the matrix 'M' must be a real embedding")
1497 return M
1498
1499 @staticmethod
1500 def jordan_product(X,Y):
1501 return (X*Y + Y*X)/2
1502
1503 @classmethod
1504 def trace_inner_product(cls,X,Y):
1505 r"""
1506 Compute the trace inner-product of two real-embeddings.
1507
1508 SETUP::
1509
1510 sage: from mjo.eja.eja_algebra import (RealSymmetricEJA,
1511 ....: ComplexHermitianEJA,
1512 ....: QuaternionHermitianEJA)
1513
1514 EXAMPLES::
1515
1516 This gives the same answer as it would if we computed the trace
1517 from the unembedded (original) matrices::
1518
1519 sage: set_random_seed()
1520 sage: J = RealSymmetricEJA.random_instance()
1521 sage: x,y = J.random_elements(2)
1522 sage: Xe = x.to_matrix()
1523 sage: Ye = y.to_matrix()
1524 sage: X = J.real_unembed(Xe)
1525 sage: Y = J.real_unembed(Ye)
1526 sage: expected = (X*Y).trace()
1527 sage: actual = J.trace_inner_product(Xe,Ye)
1528 sage: actual == expected
1529 True
1530
1531 ::
1532
1533 sage: set_random_seed()
1534 sage: J = ComplexHermitianEJA.random_instance()
1535 sage: x,y = J.random_elements(2)
1536 sage: Xe = x.to_matrix()
1537 sage: Ye = y.to_matrix()
1538 sage: X = J.real_unembed(Xe)
1539 sage: Y = J.real_unembed(Ye)
1540 sage: expected = (X*Y).trace().real()
1541 sage: actual = J.trace_inner_product(Xe,Ye)
1542 sage: actual == expected
1543 True
1544
1545 ::
1546
1547 sage: set_random_seed()
1548 sage: J = QuaternionHermitianEJA.random_instance()
1549 sage: x,y = J.random_elements(2)
1550 sage: Xe = x.to_matrix()
1551 sage: Ye = y.to_matrix()
1552 sage: X = J.real_unembed(Xe)
1553 sage: Y = J.real_unembed(Ye)
1554 sage: expected = (X*Y).trace().coefficient_tuple()[0]
1555 sage: actual = J.trace_inner_product(Xe,Ye)
1556 sage: actual == expected
1557 True
1558
1559 """
1560 Xu = cls.real_unembed(X)
1561 Yu = cls.real_unembed(Y)
1562 tr = (Xu*Yu).trace()
1563
1564 try:
1565 # Works in QQ, AA, RDF, et cetera.
1566 return tr.real()
1567 except AttributeError:
1568 # A quaternion doesn't have a real() method, but does
1569 # have coefficient_tuple() method that returns the
1570 # coefficients of 1, i, j, and k -- in that order.
1571 return tr.coefficient_tuple()[0]
1572
1573
1574 class RealMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1575 @staticmethod
1576 def dimension_over_reals():
1577 return 1
1578
1579
1580 class RealSymmetricEJA(ConcreteEuclideanJordanAlgebra,
1581 RealMatrixEuclideanJordanAlgebra):
1582 """
1583 The rank-n simple EJA consisting of real symmetric n-by-n
1584 matrices, the usual symmetric Jordan product, and the trace inner
1585 product. It has dimension `(n^2 + n)/2` over the reals.
1586
1587 SETUP::
1588
1589 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1590
1591 EXAMPLES::
1592
1593 sage: J = RealSymmetricEJA(2)
1594 sage: e0, e1, e2 = J.gens()
1595 sage: e0*e0
1596 e0
1597 sage: e1*e1
1598 1/2*e0 + 1/2*e2
1599 sage: e2*e2
1600 e2
1601
1602 In theory, our "field" can be any subfield of the reals::
1603
1604 sage: RealSymmetricEJA(2, field=RDF)
1605 Euclidean Jordan algebra of dimension 3 over Real Double Field
1606 sage: RealSymmetricEJA(2, field=RR)
1607 Euclidean Jordan algebra of dimension 3 over Real Field with
1608 53 bits of precision
1609
1610 TESTS:
1611
1612 The dimension of this algebra is `(n^2 + n) / 2`::
1613
1614 sage: set_random_seed()
1615 sage: n_max = RealSymmetricEJA._max_random_instance_size()
1616 sage: n = ZZ.random_element(1, n_max)
1617 sage: J = RealSymmetricEJA(n)
1618 sage: J.dimension() == (n^2 + n)/2
1619 True
1620
1621 The Jordan multiplication is what we think it is::
1622
1623 sage: set_random_seed()
1624 sage: J = RealSymmetricEJA.random_instance()
1625 sage: x,y = J.random_elements(2)
1626 sage: actual = (x*y).to_matrix()
1627 sage: X = x.to_matrix()
1628 sage: Y = y.to_matrix()
1629 sage: expected = (X*Y + Y*X)/2
1630 sage: actual == expected
1631 True
1632 sage: J(expected) == x*y
1633 True
1634
1635 We can change the generator prefix::
1636
1637 sage: RealSymmetricEJA(3, prefix='q').gens()
1638 (q0, q1, q2, q3, q4, q5)
1639
1640 We can construct the (trivial) algebra of rank zero::
1641
1642 sage: RealSymmetricEJA(0)
1643 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
1644
1645 """
1646 @classmethod
1647 def _denormalized_basis(cls, n):
1648 """
1649 Return a basis for the space of real symmetric n-by-n matrices.
1650
1651 SETUP::
1652
1653 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1654
1655 TESTS::
1656
1657 sage: set_random_seed()
1658 sage: n = ZZ.random_element(1,5)
1659 sage: B = RealSymmetricEJA._denormalized_basis(n)
1660 sage: all( M.is_symmetric() for M in B)
1661 True
1662
1663 """
1664 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1665 # coordinates.
1666 S = []
1667 for i in range(n):
1668 for j in range(i+1):
1669 Eij = matrix(ZZ, n, lambda k,l: k==i and l==j)
1670 if i == j:
1671 Sij = Eij
1672 else:
1673 Sij = Eij + Eij.transpose()
1674 S.append(Sij)
1675 return tuple(S)
1676
1677
1678 @staticmethod
1679 def _max_random_instance_size():
1680 return 4 # Dimension 10
1681
1682 @classmethod
1683 def random_instance(cls, **kwargs):
1684 """
1685 Return a random instance of this type of algebra.
1686 """
1687 n = ZZ.random_element(cls._max_random_instance_size() + 1)
1688 return cls(n, **kwargs)
1689
1690 def __init__(self, n, **kwargs):
1691 # We know this is a valid EJA, but will double-check
1692 # if the user passes check_axioms=True.
1693 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
1694
1695 super(RealSymmetricEJA, self).__init__(self._denormalized_basis(n),
1696 self.jordan_product,
1697 self.trace_inner_product,
1698 **kwargs)
1699
1700 # TODO: this could be factored out somehow, but is left here
1701 # because the MatrixEuclideanJordanAlgebra is not presently
1702 # a subclass of the FDEJA class that defines rank() and one().
1703 self.rank.set_cache(n)
1704 idV = matrix.identity(ZZ, self.dimension_over_reals()*n)
1705 self.one.set_cache(self(idV))
1706
1707
1708
1709 class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1710 @staticmethod
1711 def dimension_over_reals():
1712 return 2
1713
1714 @classmethod
1715 def real_embed(cls,M):
1716 """
1717 Embed the n-by-n complex matrix ``M`` into the space of real
1718 matrices of size 2n-by-2n via the map the sends each entry `z = a +
1719 bi` to the block matrix ``[[a,b],[-b,a]]``.
1720
1721 SETUP::
1722
1723 sage: from mjo.eja.eja_algebra import \
1724 ....: ComplexMatrixEuclideanJordanAlgebra
1725
1726 EXAMPLES::
1727
1728 sage: F = QuadraticField(-1, 'I')
1729 sage: x1 = F(4 - 2*i)
1730 sage: x2 = F(1 + 2*i)
1731 sage: x3 = F(-i)
1732 sage: x4 = F(6)
1733 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
1734 sage: ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1735 [ 4 -2| 1 2]
1736 [ 2 4|-2 1]
1737 [-----+-----]
1738 [ 0 -1| 6 0]
1739 [ 1 0| 0 6]
1740
1741 TESTS:
1742
1743 Embedding is a homomorphism (isomorphism, in fact)::
1744
1745 sage: set_random_seed()
1746 sage: n = ZZ.random_element(3)
1747 sage: F = QuadraticField(-1, 'I')
1748 sage: X = random_matrix(F, n)
1749 sage: Y = random_matrix(F, n)
1750 sage: Xe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X)
1751 sage: Ye = ComplexMatrixEuclideanJordanAlgebra.real_embed(Y)
1752 sage: XYe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1753 sage: Xe*Ye == XYe
1754 True
1755
1756 """
1757 super(ComplexMatrixEuclideanJordanAlgebra,cls).real_embed(M)
1758 n = M.nrows()
1759
1760 # We don't need any adjoined elements...
1761 field = M.base_ring().base_ring()
1762
1763 blocks = []
1764 for z in M.list():
1765 a = z.list()[0] # real part, I guess
1766 b = z.list()[1] # imag part, I guess
1767 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1768
1769 return matrix.block(field, n, blocks)
1770
1771
1772 @classmethod
1773 def real_unembed(cls,M):
1774 """
1775 The inverse of _embed_complex_matrix().
1776
1777 SETUP::
1778
1779 sage: from mjo.eja.eja_algebra import \
1780 ....: ComplexMatrixEuclideanJordanAlgebra
1781
1782 EXAMPLES::
1783
1784 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1785 ....: [-2, 1, -4, 3],
1786 ....: [ 9, 10, 11, 12],
1787 ....: [-10, 9, -12, 11] ])
1788 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(A)
1789 [ 2*I + 1 4*I + 3]
1790 [ 10*I + 9 12*I + 11]
1791
1792 TESTS:
1793
1794 Unembedding is the inverse of embedding::
1795
1796 sage: set_random_seed()
1797 sage: F = QuadraticField(-1, 'I')
1798 sage: M = random_matrix(F, 3)
1799 sage: Me = ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1800 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1801 True
1802
1803 """
1804 super(ComplexMatrixEuclideanJordanAlgebra,cls).real_unembed(M)
1805 n = ZZ(M.nrows())
1806 d = cls.dimension_over_reals()
1807
1808 # If "M" was normalized, its base ring might have roots
1809 # adjoined and they can stick around after unembedding.
1810 field = M.base_ring()
1811 R = PolynomialRing(field, 'z')
1812 z = R.gen()
1813
1814 # Sage doesn't know how to adjoin the complex "i" (the root of
1815 # x^2 + 1) to a field in a general way. Here, we just enumerate
1816 # all of the cases that I have cared to support so far.
1817 if field is AA:
1818 # Sage doesn't know how to embed AA into QQbar, i.e. how
1819 # to adjoin sqrt(-1) to AA.
1820 F = QQbar
1821 elif not field.is_exact():
1822 # RDF or RR
1823 F = field.complex_field()
1824 else:
1825 # Works for QQ and... maybe some other fields.
1826 F = field.extension(z**2 + 1, 'I', embedding=CLF(-1).sqrt())
1827 i = F.gen()
1828
1829 # Go top-left to bottom-right (reading order), converting every
1830 # 2-by-2 block we see to a single complex element.
1831 elements = []
1832 for k in range(n/d):
1833 for j in range(n/d):
1834 submat = M[d*k:d*k+d,d*j:d*j+d]
1835 if submat[0,0] != submat[1,1]:
1836 raise ValueError('bad on-diagonal submatrix')
1837 if submat[0,1] != -submat[1,0]:
1838 raise ValueError('bad off-diagonal submatrix')
1839 z = submat[0,0] + submat[0,1]*i
1840 elements.append(z)
1841
1842 return matrix(F, n/d, elements)
1843
1844
1845 class ComplexHermitianEJA(ConcreteEuclideanJordanAlgebra,
1846 ComplexMatrixEuclideanJordanAlgebra):
1847 """
1848 The rank-n simple EJA consisting of complex Hermitian n-by-n
1849 matrices over the real numbers, the usual symmetric Jordan product,
1850 and the real-part-of-trace inner product. It has dimension `n^2` over
1851 the reals.
1852
1853 SETUP::
1854
1855 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1856
1857 EXAMPLES:
1858
1859 In theory, our "field" can be any subfield of the reals::
1860
1861 sage: ComplexHermitianEJA(2, field=RDF)
1862 Euclidean Jordan algebra of dimension 4 over Real Double Field
1863 sage: ComplexHermitianEJA(2, field=RR)
1864 Euclidean Jordan algebra of dimension 4 over Real Field with
1865 53 bits of precision
1866
1867 TESTS:
1868
1869 The dimension of this algebra is `n^2`::
1870
1871 sage: set_random_seed()
1872 sage: n_max = ComplexHermitianEJA._max_random_instance_size()
1873 sage: n = ZZ.random_element(1, n_max)
1874 sage: J = ComplexHermitianEJA(n)
1875 sage: J.dimension() == n^2
1876 True
1877
1878 The Jordan multiplication is what we think it is::
1879
1880 sage: set_random_seed()
1881 sage: J = ComplexHermitianEJA.random_instance()
1882 sage: x,y = J.random_elements(2)
1883 sage: actual = (x*y).to_matrix()
1884 sage: X = x.to_matrix()
1885 sage: Y = y.to_matrix()
1886 sage: expected = (X*Y + Y*X)/2
1887 sage: actual == expected
1888 True
1889 sage: J(expected) == x*y
1890 True
1891
1892 We can change the generator prefix::
1893
1894 sage: ComplexHermitianEJA(2, prefix='z').gens()
1895 (z0, z1, z2, z3)
1896
1897 We can construct the (trivial) algebra of rank zero::
1898
1899 sage: ComplexHermitianEJA(0)
1900 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
1901
1902 """
1903
1904 @classmethod
1905 def _denormalized_basis(cls, n):
1906 """
1907 Returns a basis for the space of complex Hermitian n-by-n matrices.
1908
1909 Why do we embed these? Basically, because all of numerical linear
1910 algebra assumes that you're working with vectors consisting of `n`
1911 entries from a field and scalars from the same field. There's no way
1912 to tell SageMath that (for example) the vectors contain complex
1913 numbers, while the scalar field is real.
1914
1915 SETUP::
1916
1917 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1918
1919 TESTS::
1920
1921 sage: set_random_seed()
1922 sage: n = ZZ.random_element(1,5)
1923 sage: field = QuadraticField(2, 'sqrt2')
1924 sage: B = ComplexHermitianEJA._denormalized_basis(n)
1925 sage: all( M.is_symmetric() for M in B)
1926 True
1927
1928 """
1929 field = ZZ
1930 R = PolynomialRing(field, 'z')
1931 z = R.gen()
1932 F = field.extension(z**2 + 1, 'I')
1933 I = F.gen(1)
1934
1935 # This is like the symmetric case, but we need to be careful:
1936 #
1937 # * We want conjugate-symmetry, not just symmetry.
1938 # * The diagonal will (as a result) be real.
1939 #
1940 S = []
1941 for i in range(n):
1942 for j in range(i+1):
1943 Eij = matrix(F, n, lambda k,l: k==i and l==j)
1944 if i == j:
1945 Sij = cls.real_embed(Eij)
1946 S.append(Sij)
1947 else:
1948 # The second one has a minus because it's conjugated.
1949 Sij_real = cls.real_embed(Eij + Eij.transpose())
1950 S.append(Sij_real)
1951 Sij_imag = cls.real_embed(I*Eij - I*Eij.transpose())
1952 S.append(Sij_imag)
1953
1954 # Since we embedded these, we can drop back to the "field" that we
1955 # started with instead of the complex extension "F".
1956 return tuple( s.change_ring(field) for s in S )
1957
1958
1959 def __init__(self, n, **kwargs):
1960 # We know this is a valid EJA, but will double-check
1961 # if the user passes check_axioms=True.
1962 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
1963
1964 super(ComplexHermitianEJA, self).__init__(self._denormalized_basis(n),
1965 self.jordan_product,
1966 self.trace_inner_product,
1967 **kwargs)
1968 # TODO: this could be factored out somehow, but is left here
1969 # because the MatrixEuclideanJordanAlgebra is not presently
1970 # a subclass of the FDEJA class that defines rank() and one().
1971 self.rank.set_cache(n)
1972 idV = matrix.identity(ZZ, self.dimension_over_reals()*n)
1973 self.one.set_cache(self(idV))
1974
1975 @staticmethod
1976 def _max_random_instance_size():
1977 return 3 # Dimension 9
1978
1979 @classmethod
1980 def random_instance(cls, **kwargs):
1981 """
1982 Return a random instance of this type of algebra.
1983 """
1984 n = ZZ.random_element(cls._max_random_instance_size() + 1)
1985 return cls(n, **kwargs)
1986
1987 class QuaternionMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1988 @staticmethod
1989 def dimension_over_reals():
1990 return 4
1991
1992 @classmethod
1993 def real_embed(cls,M):
1994 """
1995 Embed the n-by-n quaternion matrix ``M`` into the space of real
1996 matrices of size 4n-by-4n by first sending each quaternion entry `z
1997 = a + bi + cj + dk` to the block-complex matrix ``[[a + bi,
1998 c+di],[-c + di, a-bi]]`, and then embedding those into a real
1999 matrix.
2000
2001 SETUP::
2002
2003 sage: from mjo.eja.eja_algebra import \
2004 ....: QuaternionMatrixEuclideanJordanAlgebra
2005
2006 EXAMPLES::
2007
2008 sage: Q = QuaternionAlgebra(QQ,-1,-1)
2009 sage: i,j,k = Q.gens()
2010 sage: x = 1 + 2*i + 3*j + 4*k
2011 sage: M = matrix(Q, 1, [[x]])
2012 sage: QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
2013 [ 1 2 3 4]
2014 [-2 1 -4 3]
2015 [-3 4 1 -2]
2016 [-4 -3 2 1]
2017
2018 Embedding is a homomorphism (isomorphism, in fact)::
2019
2020 sage: set_random_seed()
2021 sage: n = ZZ.random_element(2)
2022 sage: Q = QuaternionAlgebra(QQ,-1,-1)
2023 sage: X = random_matrix(Q, n)
2024 sage: Y = random_matrix(Q, n)
2025 sage: Xe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X)
2026 sage: Ye = QuaternionMatrixEuclideanJordanAlgebra.real_embed(Y)
2027 sage: XYe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X*Y)
2028 sage: Xe*Ye == XYe
2029 True
2030
2031 """
2032 super(QuaternionMatrixEuclideanJordanAlgebra,cls).real_embed(M)
2033 quaternions = M.base_ring()
2034 n = M.nrows()
2035
2036 F = QuadraticField(-1, 'I')
2037 i = F.gen()
2038
2039 blocks = []
2040 for z in M.list():
2041 t = z.coefficient_tuple()
2042 a = t[0]
2043 b = t[1]
2044 c = t[2]
2045 d = t[3]
2046 cplxM = matrix(F, 2, [[ a + b*i, c + d*i],
2047 [-c + d*i, a - b*i]])
2048 realM = ComplexMatrixEuclideanJordanAlgebra.real_embed(cplxM)
2049 blocks.append(realM)
2050
2051 # We should have real entries by now, so use the realest field
2052 # we've got for the return value.
2053 return matrix.block(quaternions.base_ring(), n, blocks)
2054
2055
2056
2057 @classmethod
2058 def real_unembed(cls,M):
2059 """
2060 The inverse of _embed_quaternion_matrix().
2061
2062 SETUP::
2063
2064 sage: from mjo.eja.eja_algebra import \
2065 ....: QuaternionMatrixEuclideanJordanAlgebra
2066
2067 EXAMPLES::
2068
2069 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
2070 ....: [-2, 1, -4, 3],
2071 ....: [-3, 4, 1, -2],
2072 ....: [-4, -3, 2, 1]])
2073 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(M)
2074 [1 + 2*i + 3*j + 4*k]
2075
2076 TESTS:
2077
2078 Unembedding is the inverse of embedding::
2079
2080 sage: set_random_seed()
2081 sage: Q = QuaternionAlgebra(QQ, -1, -1)
2082 sage: M = random_matrix(Q, 3)
2083 sage: Me = QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
2084 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
2085 True
2086
2087 """
2088 super(QuaternionMatrixEuclideanJordanAlgebra,cls).real_unembed(M)
2089 n = ZZ(M.nrows())
2090 d = cls.dimension_over_reals()
2091
2092 # Use the base ring of the matrix to ensure that its entries can be
2093 # multiplied by elements of the quaternion algebra.
2094 field = M.base_ring()
2095 Q = QuaternionAlgebra(field,-1,-1)
2096 i,j,k = Q.gens()
2097
2098 # Go top-left to bottom-right (reading order), converting every
2099 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
2100 # quaternion block.
2101 elements = []
2102 for l in range(n/d):
2103 for m in range(n/d):
2104 submat = ComplexMatrixEuclideanJordanAlgebra.real_unembed(
2105 M[d*l:d*l+d,d*m:d*m+d] )
2106 if submat[0,0] != submat[1,1].conjugate():
2107 raise ValueError('bad on-diagonal submatrix')
2108 if submat[0,1] != -submat[1,0].conjugate():
2109 raise ValueError('bad off-diagonal submatrix')
2110 z = submat[0,0].real()
2111 z += submat[0,0].imag()*i
2112 z += submat[0,1].real()*j
2113 z += submat[0,1].imag()*k
2114 elements.append(z)
2115
2116 return matrix(Q, n/d, elements)
2117
2118
2119 class QuaternionHermitianEJA(ConcreteEuclideanJordanAlgebra,
2120 QuaternionMatrixEuclideanJordanAlgebra):
2121 r"""
2122 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
2123 matrices, the usual symmetric Jordan product, and the
2124 real-part-of-trace inner product. It has dimension `2n^2 - n` over
2125 the reals.
2126
2127 SETUP::
2128
2129 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
2130
2131 EXAMPLES:
2132
2133 In theory, our "field" can be any subfield of the reals::
2134
2135 sage: QuaternionHermitianEJA(2, field=RDF)
2136 Euclidean Jordan algebra of dimension 6 over Real Double Field
2137 sage: QuaternionHermitianEJA(2, field=RR)
2138 Euclidean Jordan algebra of dimension 6 over Real Field with
2139 53 bits of precision
2140
2141 TESTS:
2142
2143 The dimension of this algebra is `2*n^2 - n`::
2144
2145 sage: set_random_seed()
2146 sage: n_max = QuaternionHermitianEJA._max_random_instance_size()
2147 sage: n = ZZ.random_element(1, n_max)
2148 sage: J = QuaternionHermitianEJA(n)
2149 sage: J.dimension() == 2*(n^2) - n
2150 True
2151
2152 The Jordan multiplication is what we think it is::
2153
2154 sage: set_random_seed()
2155 sage: J = QuaternionHermitianEJA.random_instance()
2156 sage: x,y = J.random_elements(2)
2157 sage: actual = (x*y).to_matrix()
2158 sage: X = x.to_matrix()
2159 sage: Y = y.to_matrix()
2160 sage: expected = (X*Y + Y*X)/2
2161 sage: actual == expected
2162 True
2163 sage: J(expected) == x*y
2164 True
2165
2166 We can change the generator prefix::
2167
2168 sage: QuaternionHermitianEJA(2, prefix='a').gens()
2169 (a0, a1, a2, a3, a4, a5)
2170
2171 We can construct the (trivial) algebra of rank zero::
2172
2173 sage: QuaternionHermitianEJA(0)
2174 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
2175
2176 """
2177 @classmethod
2178 def _denormalized_basis(cls, n):
2179 """
2180 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
2181
2182 Why do we embed these? Basically, because all of numerical
2183 linear algebra assumes that you're working with vectors consisting
2184 of `n` entries from a field and scalars from the same field. There's
2185 no way to tell SageMath that (for example) the vectors contain
2186 complex numbers, while the scalar field is real.
2187
2188 SETUP::
2189
2190 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
2191
2192 TESTS::
2193
2194 sage: set_random_seed()
2195 sage: n = ZZ.random_element(1,5)
2196 sage: B = QuaternionHermitianEJA._denormalized_basis(n)
2197 sage: all( M.is_symmetric() for M in B )
2198 True
2199
2200 """
2201 field = ZZ
2202 Q = QuaternionAlgebra(QQ,-1,-1)
2203 I,J,K = Q.gens()
2204
2205 # This is like the symmetric case, but we need to be careful:
2206 #
2207 # * We want conjugate-symmetry, not just symmetry.
2208 # * The diagonal will (as a result) be real.
2209 #
2210 S = []
2211 for i in range(n):
2212 for j in range(i+1):
2213 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
2214 if i == j:
2215 Sij = cls.real_embed(Eij)
2216 S.append(Sij)
2217 else:
2218 # The second, third, and fourth ones have a minus
2219 # because they're conjugated.
2220 Sij_real = cls.real_embed(Eij + Eij.transpose())
2221 S.append(Sij_real)
2222 Sij_I = cls.real_embed(I*Eij - I*Eij.transpose())
2223 S.append(Sij_I)
2224 Sij_J = cls.real_embed(J*Eij - J*Eij.transpose())
2225 S.append(Sij_J)
2226 Sij_K = cls.real_embed(K*Eij - K*Eij.transpose())
2227 S.append(Sij_K)
2228
2229 # Since we embedded these, we can drop back to the "field" that we
2230 # started with instead of the quaternion algebra "Q".
2231 return tuple( s.change_ring(field) for s in S )
2232
2233
2234 def __init__(self, n, **kwargs):
2235 # We know this is a valid EJA, but will double-check
2236 # if the user passes check_axioms=True.
2237 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
2238
2239 super(QuaternionHermitianEJA, self).__init__(self._denormalized_basis(n),
2240 self.jordan_product,
2241 self.trace_inner_product,
2242 **kwargs)
2243 # TODO: this could be factored out somehow, but is left here
2244 # because the MatrixEuclideanJordanAlgebra is not presently
2245 # a subclass of the FDEJA class that defines rank() and one().
2246 self.rank.set_cache(n)
2247 idV = matrix.identity(ZZ, self.dimension_over_reals()*n)
2248 self.one.set_cache(self(idV))
2249
2250
2251 @staticmethod
2252 def _max_random_instance_size():
2253 r"""
2254 The maximum rank of a random QuaternionHermitianEJA.
2255 """
2256 return 2 # Dimension 6
2257
2258 @classmethod
2259 def random_instance(cls, **kwargs):
2260 """
2261 Return a random instance of this type of algebra.
2262 """
2263 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2264 return cls(n, **kwargs)
2265
2266
2267 class HadamardEJA(ConcreteEuclideanJordanAlgebra):
2268 """
2269 Return the Euclidean Jordan Algebra corresponding to the set
2270 `R^n` under the Hadamard product.
2271
2272 Note: this is nothing more than the Cartesian product of ``n``
2273 copies of the spin algebra. Once Cartesian product algebras
2274 are implemented, this can go.
2275
2276 SETUP::
2277
2278 sage: from mjo.eja.eja_algebra import HadamardEJA
2279
2280 EXAMPLES:
2281
2282 This multiplication table can be verified by hand::
2283
2284 sage: J = HadamardEJA(3)
2285 sage: e0,e1,e2 = J.gens()
2286 sage: e0*e0
2287 e0
2288 sage: e0*e1
2289 0
2290 sage: e0*e2
2291 0
2292 sage: e1*e1
2293 e1
2294 sage: e1*e2
2295 0
2296 sage: e2*e2
2297 e2
2298
2299 TESTS:
2300
2301 We can change the generator prefix::
2302
2303 sage: HadamardEJA(3, prefix='r').gens()
2304 (r0, r1, r2)
2305
2306 """
2307 def __init__(self, n, **kwargs):
2308 def jordan_product(x,y):
2309 P = x.parent()
2310 return P(tuple( xi*yi for (xi,yi) in zip(x,y) ))
2311 def inner_product(x,y):
2312 return x.inner_product(y)
2313
2314 # New defaults for keyword arguments. Don't orthonormalize
2315 # because our basis is already orthonormal with respect to our
2316 # inner-product. Don't check the axioms, because we know this
2317 # is a valid EJA... but do double-check if the user passes
2318 # check_axioms=True. Note: we DON'T override the "check_field"
2319 # default here, because the user can pass in a field!
2320 if "orthonormalize" not in kwargs: kwargs["orthonormalize"] = False
2321 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
2322
2323
2324 standard_basis = FreeModule(ZZ, n).basis()
2325 super(HadamardEJA, self).__init__(standard_basis,
2326 jordan_product,
2327 inner_product,
2328 **kwargs)
2329 self.rank.set_cache(n)
2330
2331 if n == 0:
2332 self.one.set_cache( self.zero() )
2333 else:
2334 self.one.set_cache( sum(self.gens()) )
2335
2336 @staticmethod
2337 def _max_random_instance_size():
2338 r"""
2339 The maximum dimension of a random HadamardEJA.
2340 """
2341 return 5
2342
2343 @classmethod
2344 def random_instance(cls, **kwargs):
2345 """
2346 Return a random instance of this type of algebra.
2347 """
2348 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2349 return cls(n, **kwargs)
2350
2351
2352 class BilinearFormEJA(ConcreteEuclideanJordanAlgebra):
2353 r"""
2354 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2355 with the half-trace inner product and jordan product ``x*y =
2356 (<Bx,y>,y_bar>, x0*y_bar + y0*x_bar)`` where `B = 1 \times B22` is
2357 a symmetric positive-definite "bilinear form" matrix. Its
2358 dimension is the size of `B`, and it has rank two in dimensions
2359 larger than two. It reduces to the ``JordanSpinEJA`` when `B` is
2360 the identity matrix of order ``n``.
2361
2362 We insist that the one-by-one upper-left identity block of `B` be
2363 passed in as well so that we can be passed a matrix of size zero
2364 to construct a trivial algebra.
2365
2366 SETUP::
2367
2368 sage: from mjo.eja.eja_algebra import (BilinearFormEJA,
2369 ....: JordanSpinEJA)
2370
2371 EXAMPLES:
2372
2373 When no bilinear form is specified, the identity matrix is used,
2374 and the resulting algebra is the Jordan spin algebra::
2375
2376 sage: B = matrix.identity(AA,3)
2377 sage: J0 = BilinearFormEJA(B)
2378 sage: J1 = JordanSpinEJA(3)
2379 sage: J0.multiplication_table() == J0.multiplication_table()
2380 True
2381
2382 An error is raised if the matrix `B` does not correspond to a
2383 positive-definite bilinear form::
2384
2385 sage: B = matrix.random(QQ,2,3)
2386 sage: J = BilinearFormEJA(B)
2387 Traceback (most recent call last):
2388 ...
2389 ValueError: bilinear form is not positive-definite
2390 sage: B = matrix.zero(QQ,3)
2391 sage: J = BilinearFormEJA(B)
2392 Traceback (most recent call last):
2393 ...
2394 ValueError: bilinear form is not positive-definite
2395
2396 TESTS:
2397
2398 We can create a zero-dimensional algebra::
2399
2400 sage: B = matrix.identity(AA,0)
2401 sage: J = BilinearFormEJA(B)
2402 sage: J.basis()
2403 Finite family {}
2404
2405 We can check the multiplication condition given in the Jordan, von
2406 Neumann, and Wigner paper (and also discussed on my "On the
2407 symmetry..." paper). Note that this relies heavily on the standard
2408 choice of basis, as does anything utilizing the bilinear form
2409 matrix. We opt not to orthonormalize the basis, because if we
2410 did, we would have to normalize the `s_{i}` in a similar manner::
2411
2412 sage: set_random_seed()
2413 sage: n = ZZ.random_element(5)
2414 sage: M = matrix.random(QQ, max(0,n-1), algorithm='unimodular')
2415 sage: B11 = matrix.identity(QQ,1)
2416 sage: B22 = M.transpose()*M
2417 sage: B = block_matrix(2,2,[ [B11,0 ],
2418 ....: [0, B22 ] ])
2419 sage: J = BilinearFormEJA(B, orthonormalize=False)
2420 sage: eis = VectorSpace(M.base_ring(), M.ncols()).basis()
2421 sage: V = J.vector_space()
2422 sage: sis = [ J( V([0] + (M.inverse()*ei).list()).column() )
2423 ....: for ei in eis ]
2424 sage: actual = [ sis[i]*sis[j]
2425 ....: for i in range(n-1)
2426 ....: for j in range(n-1) ]
2427 sage: expected = [ J.one() if i == j else J.zero()
2428 ....: for i in range(n-1)
2429 ....: for j in range(n-1) ]
2430 sage: actual == expected
2431 True
2432 """
2433 def __init__(self, B, **kwargs):
2434 if not B.is_positive_definite():
2435 raise ValueError("bilinear form is not positive-definite")
2436
2437 def inner_product(x,y):
2438 return (B*x).inner_product(y)
2439
2440 def jordan_product(x,y):
2441 P = x.parent()
2442 x0 = x[0]
2443 xbar = x[1:]
2444 y0 = y[0]
2445 ybar = y[1:]
2446 z0 = inner_product(x,y)
2447 zbar = y0*xbar + x0*ybar
2448 return P((z0,) + tuple(zbar))
2449
2450 # We know this is a valid EJA, but will double-check
2451 # if the user passes check_axioms=True.
2452 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
2453
2454 n = B.nrows()
2455 standard_basis = FreeModule(ZZ, n).basis()
2456 super(BilinearFormEJA, self).__init__(standard_basis,
2457 jordan_product,
2458 inner_product,
2459 **kwargs)
2460
2461 # The rank of this algebra is two, unless we're in a
2462 # one-dimensional ambient space (because the rank is bounded
2463 # by the ambient dimension).
2464 self.rank.set_cache(min(n,2))
2465
2466 if n == 0:
2467 self.one.set_cache( self.zero() )
2468 else:
2469 self.one.set_cache( self.monomial(0) )
2470
2471 @staticmethod
2472 def _max_random_instance_size():
2473 r"""
2474 The maximum dimension of a random BilinearFormEJA.
2475 """
2476 return 5
2477
2478 @classmethod
2479 def random_instance(cls, **kwargs):
2480 """
2481 Return a random instance of this algebra.
2482 """
2483 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2484 if n.is_zero():
2485 B = matrix.identity(ZZ, n)
2486 return cls(B, **kwargs)
2487
2488 B11 = matrix.identity(ZZ, 1)
2489 M = matrix.random(ZZ, n-1)
2490 I = matrix.identity(ZZ, n-1)
2491 alpha = ZZ.zero()
2492 while alpha.is_zero():
2493 alpha = ZZ.random_element().abs()
2494 B22 = M.transpose()*M + alpha*I
2495
2496 from sage.matrix.special import block_matrix
2497 B = block_matrix(2,2, [ [B11, ZZ(0) ],
2498 [ZZ(0), B22 ] ])
2499
2500 return cls(B, **kwargs)
2501
2502
2503 class JordanSpinEJA(BilinearFormEJA):
2504 """
2505 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2506 with the usual inner product and jordan product ``x*y =
2507 (<x,y>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
2508 the reals.
2509
2510 SETUP::
2511
2512 sage: from mjo.eja.eja_algebra import JordanSpinEJA
2513
2514 EXAMPLES:
2515
2516 This multiplication table can be verified by hand::
2517
2518 sage: J = JordanSpinEJA(4)
2519 sage: e0,e1,e2,e3 = J.gens()
2520 sage: e0*e0
2521 e0
2522 sage: e0*e1
2523 e1
2524 sage: e0*e2
2525 e2
2526 sage: e0*e3
2527 e3
2528 sage: e1*e2
2529 0
2530 sage: e1*e3
2531 0
2532 sage: e2*e3
2533 0
2534
2535 We can change the generator prefix::
2536
2537 sage: JordanSpinEJA(2, prefix='B').gens()
2538 (B0, B1)
2539
2540 TESTS:
2541
2542 Ensure that we have the usual inner product on `R^n`::
2543
2544 sage: set_random_seed()
2545 sage: J = JordanSpinEJA.random_instance()
2546 sage: x,y = J.random_elements(2)
2547 sage: actual = x.inner_product(y)
2548 sage: expected = x.to_vector().inner_product(y.to_vector())
2549 sage: actual == expected
2550 True
2551
2552 """
2553 def __init__(self, n, **kwargs):
2554 # This is a special case of the BilinearFormEJA with the
2555 # identity matrix as its bilinear form.
2556 B = matrix.identity(ZZ, n)
2557
2558 # Don't orthonormalize because our basis is already
2559 # orthonormal with respect to our inner-product.
2560 if "orthonormalize" not in kwargs: kwargs["orthonormalize"] = False
2561
2562 # But also don't pass check_field=False here, because the user
2563 # can pass in a field!
2564 super(JordanSpinEJA, self).__init__(B, **kwargs)
2565
2566 @staticmethod
2567 def _max_random_instance_size():
2568 r"""
2569 The maximum dimension of a random JordanSpinEJA.
2570 """
2571 return 5
2572
2573 @classmethod
2574 def random_instance(cls, **kwargs):
2575 """
2576 Return a random instance of this type of algebra.
2577
2578 Needed here to override the implementation for ``BilinearFormEJA``.
2579 """
2580 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2581 return cls(n, **kwargs)
2582
2583
2584 class TrivialEJA(ConcreteEuclideanJordanAlgebra):
2585 """
2586 The trivial Euclidean Jordan algebra consisting of only a zero element.
2587
2588 SETUP::
2589
2590 sage: from mjo.eja.eja_algebra import TrivialEJA
2591
2592 EXAMPLES::
2593
2594 sage: J = TrivialEJA()
2595 sage: J.dimension()
2596 0
2597 sage: J.zero()
2598 0
2599 sage: J.one()
2600 0
2601 sage: 7*J.one()*12*J.one()
2602 0
2603 sage: J.one().inner_product(J.one())
2604 0
2605 sage: J.one().norm()
2606 0
2607 sage: J.one().subalgebra_generated_by()
2608 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
2609 sage: J.rank()
2610 0
2611
2612 """
2613 def __init__(self, **kwargs):
2614 jordan_product = lambda x,y: x
2615 inner_product = lambda x,y: 0
2616 basis = ()
2617
2618 # New defaults for keyword arguments
2619 if "orthonormalize" not in kwargs: kwargs["orthonormalize"] = False
2620 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
2621
2622 super(TrivialEJA, self).__init__(basis,
2623 jordan_product,
2624 inner_product,
2625 **kwargs)
2626 # The rank is zero using my definition, namely the dimension of the
2627 # largest subalgebra generated by any element.
2628 self.rank.set_cache(0)
2629 self.one.set_cache( self.zero() )
2630
2631 @classmethod
2632 def random_instance(cls, **kwargs):
2633 # We don't take a "size" argument so the superclass method is
2634 # inappropriate for us.
2635 return cls(**kwargs)
2636
2637 class DirectSumEJA(FiniteDimensionalEuclideanJordanAlgebra):
2638 r"""
2639 The external (orthogonal) direct sum of two other Euclidean Jordan
2640 algebras. Essentially the Cartesian product of its two factors.
2641 Every Euclidean Jordan algebra decomposes into an orthogonal
2642 direct sum of simple Euclidean Jordan algebras, so no generality
2643 is lost by providing only this construction.
2644
2645 SETUP::
2646
2647 sage: from mjo.eja.eja_algebra import (random_eja,
2648 ....: HadamardEJA,
2649 ....: RealSymmetricEJA,
2650 ....: DirectSumEJA)
2651
2652 EXAMPLES::
2653
2654 sage: J1 = HadamardEJA(2)
2655 sage: J2 = RealSymmetricEJA(3)
2656 sage: J = DirectSumEJA(J1,J2)
2657 sage: J.dimension()
2658 8
2659 sage: J.rank()
2660 5
2661
2662 TESTS:
2663
2664 The external direct sum construction is only valid when the two factors
2665 have the same base ring; an error is raised otherwise::
2666
2667 sage: set_random_seed()
2668 sage: J1 = random_eja(field=AA)
2669 sage: J2 = random_eja(field=QQ,orthonormalize=False)
2670 sage: J = DirectSumEJA(J1,J2)
2671 Traceback (most recent call last):
2672 ...
2673 ValueError: algebras must share the same base field
2674
2675 """
2676 def __init__(self, J1, J2, **kwargs):
2677 if J1.base_ring() != J2.base_ring():
2678 raise ValueError("algebras must share the same base field")
2679 field = J1.base_ring()
2680
2681 self._factors = (J1, J2)
2682 n1 = J1.dimension()
2683 n2 = J2.dimension()
2684 n = n1+n2
2685 V = VectorSpace(field, n)
2686 mult_table = [ [ V.zero() for j in range(i+1) ]
2687 for i in range(n) ]
2688 for i in range(n1):
2689 for j in range(i+1):
2690 p = (J1.monomial(i)*J1.monomial(j)).to_vector()
2691 mult_table[i][j] = V(p.list() + [field.zero()]*n2)
2692
2693 for i in range(n2):
2694 for j in range(i+1):
2695 p = (J2.monomial(i)*J2.monomial(j)).to_vector()
2696 mult_table[n1+i][n1+j] = V([field.zero()]*n1 + p.list())
2697
2698 # TODO: build the IP table here from the two constituent IP
2699 # matrices (it'll be block diagonal, I think).
2700 ip_table = [ [ field.zero() for j in range(i+1) ]
2701 for i in range(n) ]
2702 super(DirectSumEJA, self).__init__(field,
2703 mult_table,
2704 ip_table,
2705 check_axioms=False,
2706 **kwargs)
2707 self.rank.set_cache(J1.rank() + J2.rank())
2708
2709
2710 def factors(self):
2711 r"""
2712 Return the pair of this algebra's factors.
2713
2714 SETUP::
2715
2716 sage: from mjo.eja.eja_algebra import (HadamardEJA,
2717 ....: JordanSpinEJA,
2718 ....: DirectSumEJA)
2719
2720 EXAMPLES::
2721
2722 sage: J1 = HadamardEJA(2, field=QQ)
2723 sage: J2 = JordanSpinEJA(3, field=QQ)
2724 sage: J = DirectSumEJA(J1,J2)
2725 sage: J.factors()
2726 (Euclidean Jordan algebra of dimension 2 over Rational Field,
2727 Euclidean Jordan algebra of dimension 3 over Rational Field)
2728
2729 """
2730 return self._factors
2731
2732 def projections(self):
2733 r"""
2734 Return a pair of projections onto this algebra's factors.
2735
2736 SETUP::
2737
2738 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
2739 ....: ComplexHermitianEJA,
2740 ....: DirectSumEJA)
2741
2742 EXAMPLES::
2743
2744 sage: J1 = JordanSpinEJA(2)
2745 sage: J2 = ComplexHermitianEJA(2)
2746 sage: J = DirectSumEJA(J1,J2)
2747 sage: (pi_left, pi_right) = J.projections()
2748 sage: J.one().to_vector()
2749 (1, 0, 1, 0, 0, 1)
2750 sage: pi_left(J.one()).to_vector()
2751 (1, 0)
2752 sage: pi_right(J.one()).to_vector()
2753 (1, 0, 0, 1)
2754
2755 """
2756 (J1,J2) = self.factors()
2757 m = J1.dimension()
2758 n = J2.dimension()
2759 V_basis = self.vector_space().basis()
2760 # Need to specify the dimensions explicitly so that we don't
2761 # wind up with a zero-by-zero matrix when we want e.g. a
2762 # zero-by-two matrix (important for composing things).
2763 P1 = matrix(self.base_ring(), m, m+n, V_basis[:m])
2764 P2 = matrix(self.base_ring(), n, m+n, V_basis[m:])
2765 pi_left = FiniteDimensionalEuclideanJordanAlgebraOperator(self,J1,P1)
2766 pi_right = FiniteDimensionalEuclideanJordanAlgebraOperator(self,J2,P2)
2767 return (pi_left, pi_right)
2768
2769 def inclusions(self):
2770 r"""
2771 Return the pair of inclusion maps from our factors into us.
2772
2773 SETUP::
2774
2775 sage: from mjo.eja.eja_algebra import (random_eja,
2776 ....: JordanSpinEJA,
2777 ....: RealSymmetricEJA,
2778 ....: DirectSumEJA)
2779
2780 EXAMPLES::
2781
2782 sage: J1 = JordanSpinEJA(3)
2783 sage: J2 = RealSymmetricEJA(2)
2784 sage: J = DirectSumEJA(J1,J2)
2785 sage: (iota_left, iota_right) = J.inclusions()
2786 sage: iota_left(J1.zero()) == J.zero()
2787 True
2788 sage: iota_right(J2.zero()) == J.zero()
2789 True
2790 sage: J1.one().to_vector()
2791 (1, 0, 0)
2792 sage: iota_left(J1.one()).to_vector()
2793 (1, 0, 0, 0, 0, 0)
2794 sage: J2.one().to_vector()
2795 (1, 0, 1)
2796 sage: iota_right(J2.one()).to_vector()
2797 (0, 0, 0, 1, 0, 1)
2798 sage: J.one().to_vector()
2799 (1, 0, 0, 1, 0, 1)
2800
2801 TESTS:
2802
2803 Composing a projection with the corresponding inclusion should
2804 produce the identity map, and mismatching them should produce
2805 the zero map::
2806
2807 sage: set_random_seed()
2808 sage: J1 = random_eja()
2809 sage: J2 = random_eja()
2810 sage: J = DirectSumEJA(J1,J2)
2811 sage: (iota_left, iota_right) = J.inclusions()
2812 sage: (pi_left, pi_right) = J.projections()
2813 sage: pi_left*iota_left == J1.one().operator()
2814 True
2815 sage: pi_right*iota_right == J2.one().operator()
2816 True
2817 sage: (pi_left*iota_right).is_zero()
2818 True
2819 sage: (pi_right*iota_left).is_zero()
2820 True
2821
2822 """
2823 (J1,J2) = self.factors()
2824 m = J1.dimension()
2825 n = J2.dimension()
2826 V_basis = self.vector_space().basis()
2827 # Need to specify the dimensions explicitly so that we don't
2828 # wind up with a zero-by-zero matrix when we want e.g. a
2829 # two-by-zero matrix (important for composing things).
2830 I1 = matrix.column(self.base_ring(), m, m+n, V_basis[:m])
2831 I2 = matrix.column(self.base_ring(), n, m+n, V_basis[m:])
2832 iota_left = FiniteDimensionalEuclideanJordanAlgebraOperator(J1,self,I1)
2833 iota_right = FiniteDimensionalEuclideanJordanAlgebraOperator(J2,self,I2)
2834 return (iota_left, iota_right)
2835
2836 def inner_product(self, x, y):
2837 r"""
2838 The standard Cartesian inner-product.
2839
2840 We project ``x`` and ``y`` onto our factors, and add up the
2841 inner-products from the subalgebras.
2842
2843 SETUP::
2844
2845
2846 sage: from mjo.eja.eja_algebra import (HadamardEJA,
2847 ....: QuaternionHermitianEJA,
2848 ....: DirectSumEJA)
2849
2850 EXAMPLE::
2851
2852 sage: J1 = HadamardEJA(3,field=QQ)
2853 sage: J2 = QuaternionHermitianEJA(2,field=QQ,orthonormalize=False)
2854 sage: J = DirectSumEJA(J1,J2)
2855 sage: x1 = J1.one()
2856 sage: x2 = x1
2857 sage: y1 = J2.one()
2858 sage: y2 = y1
2859 sage: x1.inner_product(x2)
2860 3
2861 sage: y1.inner_product(y2)
2862 2
2863 sage: J.one().inner_product(J.one())
2864 5
2865
2866 """
2867 (pi_left, pi_right) = self.projections()
2868 x1 = pi_left(x)
2869 x2 = pi_right(x)
2870 y1 = pi_left(y)
2871 y2 = pi_right(y)
2872
2873 return (x1.inner_product(y1) + x2.inner_product(y2))
2874
2875
2876
2877 random_eja = ConcreteEuclideanJordanAlgebra.random_instance