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