]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: use symmetry when constructing the inner product matrix.
[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 n = len(basis)
1169 vector_basis = basis
1170
1171 from sage.structure.element import is_Matrix
1172 basis_is_matrices = False
1173
1174 degree = 0
1175 if n > 0:
1176 if is_Matrix(basis[0]):
1177 basis_is_matrices = True
1178 from mjo.eja.eja_utils import _vec2mat
1179 vector_basis = tuple( map(_mat2vec,basis) )
1180 degree = basis[0].nrows()**2
1181 else:
1182 degree = basis[0].degree()
1183
1184 V = VectorSpace(field, degree)
1185
1186 # If we were asked to orthonormalize, and if the orthonormal
1187 # basis is different from the given one, then we also want to
1188 # compute multiplication and inner-product tables for the
1189 # deorthonormalized basis. These can be used later to
1190 # construct a deorthonormalized copy of this algebra over QQ
1191 # in which several operations are much faster.
1192 self._deortho_multiplication_table = None
1193 self._deortho_inner_product_table = None
1194
1195 if orthonormalize:
1196 # Compute the deorthonormalized tables before we orthonormalize
1197 # the given basis.
1198 W = V.span_of_basis( vector_basis )
1199
1200 # TODO: use symmetry
1201 self._deortho_multiplication_table = [ [0 for j in range(n)]
1202 for i in range(n) ]
1203 self._deortho_inner_product_table = [ [0 for j in range(n)]
1204 for i in range(n) ]
1205
1206 # Note: the Jordan and inner-products are defined in terms
1207 # of the ambient basis. It's important that their arguments
1208 # are in ambient coordinates as well.
1209 for i in range(n):
1210 for j in range(i+1):
1211 # given basis w.r.t. ambient coords
1212 q_i = vector_basis[i]
1213 q_j = vector_basis[j]
1214
1215 if basis_is_matrices:
1216 q_i = _vec2mat(q_i)
1217 q_j = _vec2mat(q_j)
1218
1219 elt = jordan_product(q_i, q_j)
1220 ip = inner_product(q_i, q_j)
1221
1222 if basis_is_matrices:
1223 # do another mat2vec because the multiplication
1224 # table is in terms of vectors
1225 elt = _mat2vec(elt)
1226
1227 # TODO: use symmetry
1228 elt = W.coordinate_vector(elt)
1229 self._deortho_multiplication_table[i][j] = elt
1230 self._deortho_multiplication_table[j][i] = elt
1231 self._deortho_inner_product_table[i][j] = ip
1232 self._deortho_inner_product_table[j][i] = ip
1233
1234 if self._deortho_multiplication_table is not None:
1235 self._deortho_multiplication_table = tuple(map(tuple, self._deortho_multiplication_table))
1236 if self._deortho_inner_product_table is not None:
1237 self._deortho_inner_product_table = tuple(map(tuple, self._deortho_inner_product_table))
1238
1239 # We overwrite the name "vector_basis" in a second, but never modify it
1240 # in place, to this effectively makes a copy of it.
1241 deortho_vector_basis = vector_basis
1242 self._deortho_matrix = None
1243
1244 if orthonormalize:
1245 from mjo.eja.eja_utils import gram_schmidt
1246 if basis_is_matrices:
1247 vector_ip = lambda x,y: inner_product(_vec2mat(x), _vec2mat(y))
1248 vector_basis = gram_schmidt(vector_basis, vector_ip)
1249 else:
1250 vector_basis = gram_schmidt(vector_basis, inner_product)
1251
1252 W = V.span_of_basis( vector_basis )
1253
1254 # Normalize the "matrix" basis, too!
1255 basis = vector_basis
1256
1257 if basis_is_matrices:
1258 basis = tuple( map(_vec2mat,basis) )
1259
1260 W = V.span_of_basis( vector_basis )
1261
1262 # Now "W" is the vector space of our algebra coordinates. The
1263 # variables "X1", "X2",... refer to the entries of vectors in
1264 # W. Thus to convert back and forth between the orthonormal
1265 # coordinates and the given ones, we need to stick the original
1266 # basis in W.
1267 U = V.span_of_basis( deortho_vector_basis )
1268 self._deortho_matrix = matrix( U.coordinate_vector(q)
1269 for q in vector_basis )
1270
1271 # TODO: use symmetry
1272 mult_table = [ [0 for j in range(n)] for i in range(n) ]
1273 ip_table = [ [0 for j in range(n)] for i in range(n) ]
1274
1275 # Note: the Jordan and inner-products are defined in terms
1276 # of the ambient basis. It's important that their arguments
1277 # are in ambient coordinates as well.
1278 for i in range(n):
1279 for j in range(i+1):
1280 # ortho basis w.r.t. ambient coords
1281 q_i = vector_basis[i]
1282 q_j = vector_basis[j]
1283
1284 if basis_is_matrices:
1285 q_i = _vec2mat(q_i)
1286 q_j = _vec2mat(q_j)
1287
1288 elt = jordan_product(q_i, q_j)
1289 ip = inner_product(q_i, q_j)
1290
1291 if basis_is_matrices:
1292 # do another mat2vec because the multiplication
1293 # table is in terms of vectors
1294 elt = _mat2vec(elt)
1295
1296 # TODO: use symmetry
1297 elt = W.coordinate_vector(elt)
1298 mult_table[i][j] = elt
1299 mult_table[j][i] = elt
1300 ip_table[i][j] = ip
1301 ip_table[j][i] = ip
1302
1303 if basis_is_matrices:
1304 for m in basis:
1305 m.set_immutable()
1306 else:
1307 basis = tuple( x.column() for x in basis )
1308
1309 super().__init__(field,
1310 mult_table,
1311 ip_table,
1312 prefix,
1313 category,
1314 basis, # matrix basis
1315 check_field,
1316 check_axioms)
1317
1318 @cached_method
1319 def _charpoly_coefficients(self):
1320 r"""
1321 SETUP::
1322
1323 sage: from mjo.eja.eja_algebra import (BilinearFormEJA,
1324 ....: JordanSpinEJA)
1325
1326 EXAMPLES:
1327
1328 The base ring of the resulting polynomial coefficients is what
1329 it should be, and not the rationals (unless the algebra was
1330 already over the rationals)::
1331
1332 sage: J = JordanSpinEJA(3)
1333 sage: J._charpoly_coefficients()
1334 (X1^2 - X2^2 - X3^2, -2*X1)
1335 sage: a0 = J._charpoly_coefficients()[0]
1336 sage: J.base_ring()
1337 Algebraic Real Field
1338 sage: a0.base_ring()
1339 Algebraic Real Field
1340
1341 """
1342 if self.base_ring() is QQ:
1343 # There's no need to construct *another* algebra over the
1344 # rationals if this one is already over the rationals.
1345 superclass = super(RationalBasisEuclideanJordanAlgebra, self)
1346 return superclass._charpoly_coefficients()
1347
1348 # Do the computation over the rationals. The answer will be
1349 # the same, because all we've done is a change of basis.
1350 J = FiniteDimensionalEuclideanJordanAlgebra(QQ,
1351 self._deortho_multiplication_table,
1352 self._deortho_inner_product_table)
1353
1354 # Change back from QQ to our real base ring
1355 a = ( a_i.change_ring(self.base_ring())
1356 for a_i in J._charpoly_coefficients() )
1357
1358 # Now convert the coordinate variables back to the
1359 # deorthonormalized ones.
1360 R = self.coordinate_polynomial_ring()
1361 from sage.modules.free_module_element import vector
1362 X = vector(R, R.gens())
1363 BX = self._deortho_matrix*X
1364
1365 subs_dict = { X[i]: BX[i] for i in range(len(X)) }
1366 return tuple( a_i.subs(subs_dict) for a_i in a )
1367
1368 class ConcreteEuclideanJordanAlgebra(RationalBasisEuclideanJordanAlgebra):
1369 r"""
1370 A class for the Euclidean Jordan algebras that we know by name.
1371
1372 These are the Jordan algebras whose basis, multiplication table,
1373 rank, and so on are known a priori. More to the point, they are
1374 the Euclidean Jordan algebras for which we are able to conjure up
1375 a "random instance."
1376
1377 SETUP::
1378
1379 sage: from mjo.eja.eja_algebra import ConcreteEuclideanJordanAlgebra
1380
1381 TESTS:
1382
1383 Our basis is normalized with respect to the algebra's inner
1384 product, unless we specify otherwise::
1385
1386 sage: set_random_seed()
1387 sage: J = ConcreteEuclideanJordanAlgebra.random_instance()
1388 sage: all( b.norm() == 1 for b in J.gens() )
1389 True
1390
1391 Since our basis is orthonormal with respect to the algebra's inner
1392 product, and since we know that this algebra is an EJA, any
1393 left-multiplication operator's matrix will be symmetric because
1394 natural->EJA basis representation is an isometry and within the
1395 EJA the operator is self-adjoint by the Jordan axiom::
1396
1397 sage: set_random_seed()
1398 sage: J = ConcreteEuclideanJordanAlgebra.random_instance()
1399 sage: x = J.random_element()
1400 sage: x.operator().is_self_adjoint()
1401 True
1402 """
1403
1404 @staticmethod
1405 def _max_random_instance_size():
1406 """
1407 Return an integer "size" that is an upper bound on the size of
1408 this algebra when it is used in a random test
1409 case. Unfortunately, the term "size" is ambiguous -- when
1410 dealing with `R^n` under either the Hadamard or Jordan spin
1411 product, the "size" refers to the dimension `n`. When dealing
1412 with a matrix algebra (real symmetric or complex/quaternion
1413 Hermitian), it refers to the size of the matrix, which is far
1414 less than the dimension of the underlying vector space.
1415
1416 This method must be implemented in each subclass.
1417 """
1418 raise NotImplementedError
1419
1420 @classmethod
1421 def random_instance(cls, *args, **kwargs):
1422 """
1423 Return a random instance of this type of algebra.
1424
1425 This method should be implemented in each subclass.
1426 """
1427 from sage.misc.prandom import choice
1428 eja_class = choice(cls.__subclasses__())
1429
1430 # These all bubble up to the RationalBasisEuclideanJordanAlgebra
1431 # superclass constructor, so any (kw)args valid there are also
1432 # valid here.
1433 return eja_class.random_instance(*args, **kwargs)
1434
1435
1436 class MatrixEuclideanJordanAlgebra:
1437 @staticmethod
1438 def real_embed(M):
1439 """
1440 Embed the matrix ``M`` into a space of real matrices.
1441
1442 The matrix ``M`` can have entries in any field at the moment:
1443 the real numbers, complex numbers, or quaternions. And although
1444 they are not a field, we can probably support octonions at some
1445 point, too. This function returns a real matrix that "acts like"
1446 the original with respect to matrix multiplication; i.e.
1447
1448 real_embed(M*N) = real_embed(M)*real_embed(N)
1449
1450 """
1451 raise NotImplementedError
1452
1453
1454 @staticmethod
1455 def real_unembed(M):
1456 """
1457 The inverse of :meth:`real_embed`.
1458 """
1459 raise NotImplementedError
1460
1461 @staticmethod
1462 def jordan_product(X,Y):
1463 return (X*Y + Y*X)/2
1464
1465 @classmethod
1466 def trace_inner_product(cls,X,Y):
1467 Xu = cls.real_unembed(X)
1468 Yu = cls.real_unembed(Y)
1469 tr = (Xu*Yu).trace()
1470
1471 try:
1472 # Works in QQ, AA, RDF, et cetera.
1473 return tr.real()
1474 except AttributeError:
1475 # A quaternion doesn't have a real() method, but does
1476 # have coefficient_tuple() method that returns the
1477 # coefficients of 1, i, j, and k -- in that order.
1478 return tr.coefficient_tuple()[0]
1479
1480
1481 class RealMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1482 @staticmethod
1483 def real_embed(M):
1484 """
1485 The identity function, for embedding real matrices into real
1486 matrices.
1487 """
1488 return M
1489
1490 @staticmethod
1491 def real_unembed(M):
1492 """
1493 The identity function, for unembedding real matrices from real
1494 matrices.
1495 """
1496 return M
1497
1498
1499 class RealSymmetricEJA(ConcreteEuclideanJordanAlgebra,
1500 RealMatrixEuclideanJordanAlgebra):
1501 """
1502 The rank-n simple EJA consisting of real symmetric n-by-n
1503 matrices, the usual symmetric Jordan product, and the trace inner
1504 product. It has dimension `(n^2 + n)/2` over the reals.
1505
1506 SETUP::
1507
1508 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1509
1510 EXAMPLES::
1511
1512 sage: J = RealSymmetricEJA(2)
1513 sage: e0, e1, e2 = J.gens()
1514 sage: e0*e0
1515 e0
1516 sage: e1*e1
1517 1/2*e0 + 1/2*e2
1518 sage: e2*e2
1519 e2
1520
1521 In theory, our "field" can be any subfield of the reals::
1522
1523 sage: RealSymmetricEJA(2, RDF)
1524 Euclidean Jordan algebra of dimension 3 over Real Double Field
1525 sage: RealSymmetricEJA(2, RR)
1526 Euclidean Jordan algebra of dimension 3 over Real Field with
1527 53 bits of precision
1528
1529 TESTS:
1530
1531 The dimension of this algebra is `(n^2 + n) / 2`::
1532
1533 sage: set_random_seed()
1534 sage: n_max = RealSymmetricEJA._max_random_instance_size()
1535 sage: n = ZZ.random_element(1, n_max)
1536 sage: J = RealSymmetricEJA(n)
1537 sage: J.dimension() == (n^2 + n)/2
1538 True
1539
1540 The Jordan multiplication is what we think it is::
1541
1542 sage: set_random_seed()
1543 sage: J = RealSymmetricEJA.random_instance()
1544 sage: x,y = J.random_elements(2)
1545 sage: actual = (x*y).to_matrix()
1546 sage: X = x.to_matrix()
1547 sage: Y = y.to_matrix()
1548 sage: expected = (X*Y + Y*X)/2
1549 sage: actual == expected
1550 True
1551 sage: J(expected) == x*y
1552 True
1553
1554 We can change the generator prefix::
1555
1556 sage: RealSymmetricEJA(3, prefix='q').gens()
1557 (q0, q1, q2, q3, q4, q5)
1558
1559 We can construct the (trivial) algebra of rank zero::
1560
1561 sage: RealSymmetricEJA(0)
1562 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
1563
1564 """
1565 @classmethod
1566 def _denormalized_basis(cls, n, field):
1567 """
1568 Return a basis for the space of real symmetric n-by-n matrices.
1569
1570 SETUP::
1571
1572 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1573
1574 TESTS::
1575
1576 sage: set_random_seed()
1577 sage: n = ZZ.random_element(1,5)
1578 sage: B = RealSymmetricEJA._denormalized_basis(n,QQ)
1579 sage: all( M.is_symmetric() for M in B)
1580 True
1581
1582 """
1583 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1584 # coordinates.
1585 S = []
1586 for i in range(n):
1587 for j in range(i+1):
1588 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1589 if i == j:
1590 Sij = Eij
1591 else:
1592 Sij = Eij + Eij.transpose()
1593 S.append(Sij)
1594 return S
1595
1596
1597 @staticmethod
1598 def _max_random_instance_size():
1599 return 4 # Dimension 10
1600
1601 @classmethod
1602 def random_instance(cls, field=AA, **kwargs):
1603 """
1604 Return a random instance of this type of algebra.
1605 """
1606 n = ZZ.random_element(cls._max_random_instance_size() + 1)
1607 return cls(n, field, **kwargs)
1608
1609 def __init__(self, n, field=AA, **kwargs):
1610 basis = self._denormalized_basis(n, field)
1611 super(RealSymmetricEJA, self).__init__(field,
1612 basis,
1613 self.jordan_product,
1614 self.trace_inner_product,
1615 **kwargs)
1616 self.rank.set_cache(n)
1617 self.one.set_cache(self(matrix.identity(field,n)))
1618
1619
1620 class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1621 @staticmethod
1622 def real_embed(M):
1623 """
1624 Embed the n-by-n complex matrix ``M`` into the space of real
1625 matrices of size 2n-by-2n via the map the sends each entry `z = a +
1626 bi` to the block matrix ``[[a,b],[-b,a]]``.
1627
1628 SETUP::
1629
1630 sage: from mjo.eja.eja_algebra import \
1631 ....: ComplexMatrixEuclideanJordanAlgebra
1632
1633 EXAMPLES::
1634
1635 sage: F = QuadraticField(-1, 'I')
1636 sage: x1 = F(4 - 2*i)
1637 sage: x2 = F(1 + 2*i)
1638 sage: x3 = F(-i)
1639 sage: x4 = F(6)
1640 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
1641 sage: ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1642 [ 4 -2| 1 2]
1643 [ 2 4|-2 1]
1644 [-----+-----]
1645 [ 0 -1| 6 0]
1646 [ 1 0| 0 6]
1647
1648 TESTS:
1649
1650 Embedding is a homomorphism (isomorphism, in fact)::
1651
1652 sage: set_random_seed()
1653 sage: n = ZZ.random_element(3)
1654 sage: F = QuadraticField(-1, 'I')
1655 sage: X = random_matrix(F, n)
1656 sage: Y = random_matrix(F, n)
1657 sage: Xe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X)
1658 sage: Ye = ComplexMatrixEuclideanJordanAlgebra.real_embed(Y)
1659 sage: XYe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1660 sage: Xe*Ye == XYe
1661 True
1662
1663 """
1664 n = M.nrows()
1665 if M.ncols() != n:
1666 raise ValueError("the matrix 'M' must be square")
1667
1668 # We don't need any adjoined elements...
1669 field = M.base_ring().base_ring()
1670
1671 blocks = []
1672 for z in M.list():
1673 a = z.list()[0] # real part, I guess
1674 b = z.list()[1] # imag part, I guess
1675 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1676
1677 return matrix.block(field, n, blocks)
1678
1679
1680 @staticmethod
1681 def real_unembed(M):
1682 """
1683 The inverse of _embed_complex_matrix().
1684
1685 SETUP::
1686
1687 sage: from mjo.eja.eja_algebra import \
1688 ....: ComplexMatrixEuclideanJordanAlgebra
1689
1690 EXAMPLES::
1691
1692 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1693 ....: [-2, 1, -4, 3],
1694 ....: [ 9, 10, 11, 12],
1695 ....: [-10, 9, -12, 11] ])
1696 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(A)
1697 [ 2*I + 1 4*I + 3]
1698 [ 10*I + 9 12*I + 11]
1699
1700 TESTS:
1701
1702 Unembedding is the inverse of embedding::
1703
1704 sage: set_random_seed()
1705 sage: F = QuadraticField(-1, 'I')
1706 sage: M = random_matrix(F, 3)
1707 sage: Me = ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1708 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1709 True
1710
1711 """
1712 n = ZZ(M.nrows())
1713 if M.ncols() != n:
1714 raise ValueError("the matrix 'M' must be square")
1715 if not n.mod(2).is_zero():
1716 raise ValueError("the matrix 'M' must be a complex embedding")
1717
1718 # If "M" was normalized, its base ring might have roots
1719 # adjoined and they can stick around after unembedding.
1720 field = M.base_ring()
1721 R = PolynomialRing(field, 'z')
1722 z = R.gen()
1723 if field is AA:
1724 # Sage doesn't know how to embed AA into QQbar, i.e. how
1725 # to adjoin sqrt(-1) to AA.
1726 F = QQbar
1727 else:
1728 F = field.extension(z**2 + 1, 'I', embedding=CLF(-1).sqrt())
1729 i = F.gen()
1730
1731 # Go top-left to bottom-right (reading order), converting every
1732 # 2-by-2 block we see to a single complex element.
1733 elements = []
1734 for k in range(n/2):
1735 for j in range(n/2):
1736 submat = M[2*k:2*k+2,2*j:2*j+2]
1737 if submat[0,0] != submat[1,1]:
1738 raise ValueError('bad on-diagonal submatrix')
1739 if submat[0,1] != -submat[1,0]:
1740 raise ValueError('bad off-diagonal submatrix')
1741 z = submat[0,0] + submat[0,1]*i
1742 elements.append(z)
1743
1744 return matrix(F, n/2, elements)
1745
1746
1747 @classmethod
1748 def trace_inner_product(cls,X,Y):
1749 """
1750 Compute a matrix inner product in this algebra directly from
1751 its real embedding.
1752
1753 SETUP::
1754
1755 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1756
1757 TESTS:
1758
1759 This gives the same answer as the slow, default method implemented
1760 in :class:`MatrixEuclideanJordanAlgebra`::
1761
1762 sage: set_random_seed()
1763 sage: J = ComplexHermitianEJA.random_instance()
1764 sage: x,y = J.random_elements(2)
1765 sage: Xe = x.to_matrix()
1766 sage: Ye = y.to_matrix()
1767 sage: X = ComplexHermitianEJA.real_unembed(Xe)
1768 sage: Y = ComplexHermitianEJA.real_unembed(Ye)
1769 sage: expected = (X*Y).trace().real()
1770 sage: actual = ComplexHermitianEJA.trace_inner_product(Xe,Ye)
1771 sage: actual == expected
1772 True
1773
1774 """
1775 return RealMatrixEuclideanJordanAlgebra.trace_inner_product(X,Y)/2
1776
1777
1778 class ComplexHermitianEJA(ConcreteEuclideanJordanAlgebra,
1779 ComplexMatrixEuclideanJordanAlgebra):
1780 """
1781 The rank-n simple EJA consisting of complex Hermitian n-by-n
1782 matrices over the real numbers, the usual symmetric Jordan product,
1783 and the real-part-of-trace inner product. It has dimension `n^2` over
1784 the reals.
1785
1786 SETUP::
1787
1788 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1789
1790 EXAMPLES:
1791
1792 In theory, our "field" can be any subfield of the reals::
1793
1794 sage: ComplexHermitianEJA(2, RDF)
1795 Euclidean Jordan algebra of dimension 4 over Real Double Field
1796 sage: ComplexHermitianEJA(2, RR)
1797 Euclidean Jordan algebra of dimension 4 over Real Field with
1798 53 bits of precision
1799
1800 TESTS:
1801
1802 The dimension of this algebra is `n^2`::
1803
1804 sage: set_random_seed()
1805 sage: n_max = ComplexHermitianEJA._max_random_instance_size()
1806 sage: n = ZZ.random_element(1, n_max)
1807 sage: J = ComplexHermitianEJA(n)
1808 sage: J.dimension() == n^2
1809 True
1810
1811 The Jordan multiplication is what we think it is::
1812
1813 sage: set_random_seed()
1814 sage: J = ComplexHermitianEJA.random_instance()
1815 sage: x,y = J.random_elements(2)
1816 sage: actual = (x*y).to_matrix()
1817 sage: X = x.to_matrix()
1818 sage: Y = y.to_matrix()
1819 sage: expected = (X*Y + Y*X)/2
1820 sage: actual == expected
1821 True
1822 sage: J(expected) == x*y
1823 True
1824
1825 We can change the generator prefix::
1826
1827 sage: ComplexHermitianEJA(2, prefix='z').gens()
1828 (z0, z1, z2, z3)
1829
1830 We can construct the (trivial) algebra of rank zero::
1831
1832 sage: ComplexHermitianEJA(0)
1833 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
1834
1835 """
1836
1837 @classmethod
1838 def _denormalized_basis(cls, n, field):
1839 """
1840 Returns a basis for the space of complex Hermitian n-by-n matrices.
1841
1842 Why do we embed these? Basically, because all of numerical linear
1843 algebra assumes that you're working with vectors consisting of `n`
1844 entries from a field and scalars from the same field. There's no way
1845 to tell SageMath that (for example) the vectors contain complex
1846 numbers, while the scalar field is real.
1847
1848 SETUP::
1849
1850 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1851
1852 TESTS::
1853
1854 sage: set_random_seed()
1855 sage: n = ZZ.random_element(1,5)
1856 sage: field = QuadraticField(2, 'sqrt2')
1857 sage: B = ComplexHermitianEJA._denormalized_basis(n, field)
1858 sage: all( M.is_symmetric() for M in B)
1859 True
1860
1861 """
1862 R = PolynomialRing(field, 'z')
1863 z = R.gen()
1864 F = field.extension(z**2 + 1, 'I')
1865 I = F.gen()
1866
1867 # This is like the symmetric case, but we need to be careful:
1868 #
1869 # * We want conjugate-symmetry, not just symmetry.
1870 # * The diagonal will (as a result) be real.
1871 #
1872 S = []
1873 for i in range(n):
1874 for j in range(i+1):
1875 Eij = matrix(F, n, lambda k,l: k==i and l==j)
1876 if i == j:
1877 Sij = cls.real_embed(Eij)
1878 S.append(Sij)
1879 else:
1880 # The second one has a minus because it's conjugated.
1881 Sij_real = cls.real_embed(Eij + Eij.transpose())
1882 S.append(Sij_real)
1883 Sij_imag = cls.real_embed(I*Eij - I*Eij.transpose())
1884 S.append(Sij_imag)
1885
1886 # Since we embedded these, we can drop back to the "field" that we
1887 # started with instead of the complex extension "F".
1888 return tuple( s.change_ring(field) for s in S )
1889
1890
1891 def __init__(self, n, field=AA, **kwargs):
1892 basis = self._denormalized_basis(n,field)
1893 super(ComplexHermitianEJA, self).__init__(field,
1894 basis,
1895 self.jordan_product,
1896 self.trace_inner_product,
1897 **kwargs)
1898 self.rank.set_cache(n)
1899 # TODO: pre-cache the identity!
1900
1901 @staticmethod
1902 def _max_random_instance_size():
1903 return 3 # Dimension 9
1904
1905 @classmethod
1906 def random_instance(cls, field=AA, **kwargs):
1907 """
1908 Return a random instance of this type of algebra.
1909 """
1910 n = ZZ.random_element(cls._max_random_instance_size() + 1)
1911 return cls(n, field, **kwargs)
1912
1913 class QuaternionMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1914 @staticmethod
1915 def real_embed(M):
1916 """
1917 Embed the n-by-n quaternion matrix ``M`` into the space of real
1918 matrices of size 4n-by-4n by first sending each quaternion entry `z
1919 = a + bi + cj + dk` to the block-complex matrix ``[[a + bi,
1920 c+di],[-c + di, a-bi]]`, and then embedding those into a real
1921 matrix.
1922
1923 SETUP::
1924
1925 sage: from mjo.eja.eja_algebra import \
1926 ....: QuaternionMatrixEuclideanJordanAlgebra
1927
1928 EXAMPLES::
1929
1930 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1931 sage: i,j,k = Q.gens()
1932 sage: x = 1 + 2*i + 3*j + 4*k
1933 sage: M = matrix(Q, 1, [[x]])
1934 sage: QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1935 [ 1 2 3 4]
1936 [-2 1 -4 3]
1937 [-3 4 1 -2]
1938 [-4 -3 2 1]
1939
1940 Embedding is a homomorphism (isomorphism, in fact)::
1941
1942 sage: set_random_seed()
1943 sage: n = ZZ.random_element(2)
1944 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1945 sage: X = random_matrix(Q, n)
1946 sage: Y = random_matrix(Q, n)
1947 sage: Xe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X)
1948 sage: Ye = QuaternionMatrixEuclideanJordanAlgebra.real_embed(Y)
1949 sage: XYe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1950 sage: Xe*Ye == XYe
1951 True
1952
1953 """
1954 quaternions = M.base_ring()
1955 n = M.nrows()
1956 if M.ncols() != n:
1957 raise ValueError("the matrix 'M' must be square")
1958
1959 F = QuadraticField(-1, 'I')
1960 i = F.gen()
1961
1962 blocks = []
1963 for z in M.list():
1964 t = z.coefficient_tuple()
1965 a = t[0]
1966 b = t[1]
1967 c = t[2]
1968 d = t[3]
1969 cplxM = matrix(F, 2, [[ a + b*i, c + d*i],
1970 [-c + d*i, a - b*i]])
1971 realM = ComplexMatrixEuclideanJordanAlgebra.real_embed(cplxM)
1972 blocks.append(realM)
1973
1974 # We should have real entries by now, so use the realest field
1975 # we've got for the return value.
1976 return matrix.block(quaternions.base_ring(), n, blocks)
1977
1978
1979
1980 @staticmethod
1981 def real_unembed(M):
1982 """
1983 The inverse of _embed_quaternion_matrix().
1984
1985 SETUP::
1986
1987 sage: from mjo.eja.eja_algebra import \
1988 ....: QuaternionMatrixEuclideanJordanAlgebra
1989
1990 EXAMPLES::
1991
1992 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
1993 ....: [-2, 1, -4, 3],
1994 ....: [-3, 4, 1, -2],
1995 ....: [-4, -3, 2, 1]])
1996 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(M)
1997 [1 + 2*i + 3*j + 4*k]
1998
1999 TESTS:
2000
2001 Unembedding is the inverse of embedding::
2002
2003 sage: set_random_seed()
2004 sage: Q = QuaternionAlgebra(QQ, -1, -1)
2005 sage: M = random_matrix(Q, 3)
2006 sage: Me = QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
2007 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
2008 True
2009
2010 """
2011 n = ZZ(M.nrows())
2012 if M.ncols() != n:
2013 raise ValueError("the matrix 'M' must be square")
2014 if not n.mod(4).is_zero():
2015 raise ValueError("the matrix 'M' must be a quaternion embedding")
2016
2017 # Use the base ring of the matrix to ensure that its entries can be
2018 # multiplied by elements of the quaternion algebra.
2019 field = M.base_ring()
2020 Q = QuaternionAlgebra(field,-1,-1)
2021 i,j,k = Q.gens()
2022
2023 # Go top-left to bottom-right (reading order), converting every
2024 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
2025 # quaternion block.
2026 elements = []
2027 for l in range(n/4):
2028 for m in range(n/4):
2029 submat = ComplexMatrixEuclideanJordanAlgebra.real_unembed(
2030 M[4*l:4*l+4,4*m:4*m+4] )
2031 if submat[0,0] != submat[1,1].conjugate():
2032 raise ValueError('bad on-diagonal submatrix')
2033 if submat[0,1] != -submat[1,0].conjugate():
2034 raise ValueError('bad off-diagonal submatrix')
2035 z = submat[0,0].real()
2036 z += submat[0,0].imag()*i
2037 z += submat[0,1].real()*j
2038 z += submat[0,1].imag()*k
2039 elements.append(z)
2040
2041 return matrix(Q, n/4, elements)
2042
2043
2044 @classmethod
2045 def trace_inner_product(cls,X,Y):
2046 """
2047 Compute a matrix inner product in this algebra directly from
2048 its real embedding.
2049
2050 SETUP::
2051
2052 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
2053
2054 TESTS:
2055
2056 This gives the same answer as the slow, default method implemented
2057 in :class:`MatrixEuclideanJordanAlgebra`::
2058
2059 sage: set_random_seed()
2060 sage: J = QuaternionHermitianEJA.random_instance()
2061 sage: x,y = J.random_elements(2)
2062 sage: Xe = x.to_matrix()
2063 sage: Ye = y.to_matrix()
2064 sage: X = QuaternionHermitianEJA.real_unembed(Xe)
2065 sage: Y = QuaternionHermitianEJA.real_unembed(Ye)
2066 sage: expected = (X*Y).trace().coefficient_tuple()[0]
2067 sage: actual = QuaternionHermitianEJA.trace_inner_product(Xe,Ye)
2068 sage: actual == expected
2069 True
2070
2071 """
2072 return RealMatrixEuclideanJordanAlgebra.trace_inner_product(X,Y)/4
2073
2074
2075 class QuaternionHermitianEJA(ConcreteEuclideanJordanAlgebra,
2076 QuaternionMatrixEuclideanJordanAlgebra):
2077 r"""
2078 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
2079 matrices, the usual symmetric Jordan product, and the
2080 real-part-of-trace inner product. It has dimension `2n^2 - n` over
2081 the reals.
2082
2083 SETUP::
2084
2085 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
2086
2087 EXAMPLES:
2088
2089 In theory, our "field" can be any subfield of the reals::
2090
2091 sage: QuaternionHermitianEJA(2, RDF)
2092 Euclidean Jordan algebra of dimension 6 over Real Double Field
2093 sage: QuaternionHermitianEJA(2, RR)
2094 Euclidean Jordan algebra of dimension 6 over Real Field with
2095 53 bits of precision
2096
2097 TESTS:
2098
2099 The dimension of this algebra is `2*n^2 - n`::
2100
2101 sage: set_random_seed()
2102 sage: n_max = QuaternionHermitianEJA._max_random_instance_size()
2103 sage: n = ZZ.random_element(1, n_max)
2104 sage: J = QuaternionHermitianEJA(n)
2105 sage: J.dimension() == 2*(n^2) - n
2106 True
2107
2108 The Jordan multiplication is what we think it is::
2109
2110 sage: set_random_seed()
2111 sage: J = QuaternionHermitianEJA.random_instance()
2112 sage: x,y = J.random_elements(2)
2113 sage: actual = (x*y).to_matrix()
2114 sage: X = x.to_matrix()
2115 sage: Y = y.to_matrix()
2116 sage: expected = (X*Y + Y*X)/2
2117 sage: actual == expected
2118 True
2119 sage: J(expected) == x*y
2120 True
2121
2122 We can change the generator prefix::
2123
2124 sage: QuaternionHermitianEJA(2, prefix='a').gens()
2125 (a0, a1, a2, a3, a4, a5)
2126
2127 We can construct the (trivial) algebra of rank zero::
2128
2129 sage: QuaternionHermitianEJA(0)
2130 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
2131
2132 """
2133 @classmethod
2134 def _denormalized_basis(cls, n, field):
2135 """
2136 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
2137
2138 Why do we embed these? Basically, because all of numerical
2139 linear algebra assumes that you're working with vectors consisting
2140 of `n` entries from a field and scalars from the same field. There's
2141 no way to tell SageMath that (for example) the vectors contain
2142 complex numbers, while the scalar field is real.
2143
2144 SETUP::
2145
2146 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
2147
2148 TESTS::
2149
2150 sage: set_random_seed()
2151 sage: n = ZZ.random_element(1,5)
2152 sage: B = QuaternionHermitianEJA._denormalized_basis(n,QQ)
2153 sage: all( M.is_symmetric() for M in B )
2154 True
2155
2156 """
2157 Q = QuaternionAlgebra(QQ,-1,-1)
2158 I,J,K = Q.gens()
2159
2160 # This is like the symmetric case, but we need to be careful:
2161 #
2162 # * We want conjugate-symmetry, not just symmetry.
2163 # * The diagonal will (as a result) be real.
2164 #
2165 S = []
2166 for i in range(n):
2167 for j in range(i+1):
2168 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
2169 if i == j:
2170 Sij = cls.real_embed(Eij)
2171 S.append(Sij)
2172 else:
2173 # The second, third, and fourth ones have a minus
2174 # because they're conjugated.
2175 Sij_real = cls.real_embed(Eij + Eij.transpose())
2176 S.append(Sij_real)
2177 Sij_I = cls.real_embed(I*Eij - I*Eij.transpose())
2178 S.append(Sij_I)
2179 Sij_J = cls.real_embed(J*Eij - J*Eij.transpose())
2180 S.append(Sij_J)
2181 Sij_K = cls.real_embed(K*Eij - K*Eij.transpose())
2182 S.append(Sij_K)
2183
2184 # Since we embedded these, we can drop back to the "field" that we
2185 # started with instead of the quaternion algebra "Q".
2186 return tuple( s.change_ring(field) for s in S )
2187
2188
2189 def __init__(self, n, field=AA, **kwargs):
2190 basis = self._denormalized_basis(n,field)
2191 super(QuaternionHermitianEJA, self).__init__(field,
2192 basis,
2193 self.jordan_product,
2194 self.trace_inner_product,
2195 **kwargs)
2196 self.rank.set_cache(n)
2197 # TODO: cache one()!
2198
2199 @staticmethod
2200 def _max_random_instance_size():
2201 r"""
2202 The maximum rank of a random QuaternionHermitianEJA.
2203 """
2204 return 2 # Dimension 6
2205
2206 @classmethod
2207 def random_instance(cls, field=AA, **kwargs):
2208 """
2209 Return a random instance of this type of algebra.
2210 """
2211 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2212 return cls(n, field, **kwargs)
2213
2214
2215 class HadamardEJA(ConcreteEuclideanJordanAlgebra):
2216 """
2217 Return the Euclidean Jordan Algebra corresponding to the set
2218 `R^n` under the Hadamard product.
2219
2220 Note: this is nothing more than the Cartesian product of ``n``
2221 copies of the spin algebra. Once Cartesian product algebras
2222 are implemented, this can go.
2223
2224 SETUP::
2225
2226 sage: from mjo.eja.eja_algebra import HadamardEJA
2227
2228 EXAMPLES:
2229
2230 This multiplication table can be verified by hand::
2231
2232 sage: J = HadamardEJA(3)
2233 sage: e0,e1,e2 = J.gens()
2234 sage: e0*e0
2235 e0
2236 sage: e0*e1
2237 0
2238 sage: e0*e2
2239 0
2240 sage: e1*e1
2241 e1
2242 sage: e1*e2
2243 0
2244 sage: e2*e2
2245 e2
2246
2247 TESTS:
2248
2249 We can change the generator prefix::
2250
2251 sage: HadamardEJA(3, prefix='r').gens()
2252 (r0, r1, r2)
2253
2254 """
2255 def __init__(self, n, field=AA, **kwargs):
2256 V = VectorSpace(field, n)
2257 basis = V.basis()
2258
2259 def jordan_product(x,y):
2260 return V([ xi*yi for (xi,yi) in zip(x,y) ])
2261 def inner_product(x,y):
2262 return x.inner_product(y)
2263
2264 super(HadamardEJA, self).__init__(field,
2265 basis,
2266 jordan_product,
2267 inner_product,
2268 **kwargs)
2269 self.rank.set_cache(n)
2270
2271 if n == 0:
2272 self.one.set_cache( self.zero() )
2273 else:
2274 self.one.set_cache( sum(self.gens()) )
2275
2276 @staticmethod
2277 def _max_random_instance_size():
2278 r"""
2279 The maximum dimension of a random HadamardEJA.
2280 """
2281 return 5
2282
2283 @classmethod
2284 def random_instance(cls, field=AA, **kwargs):
2285 """
2286 Return a random instance of this type of algebra.
2287 """
2288 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2289 return cls(n, field, **kwargs)
2290
2291
2292 class BilinearFormEJA(ConcreteEuclideanJordanAlgebra):
2293 r"""
2294 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2295 with the half-trace inner product and jordan product ``x*y =
2296 (<Bx,y>,y_bar>, x0*y_bar + y0*x_bar)`` where `B = 1 \times B22` is
2297 a symmetric positive-definite "bilinear form" matrix. Its
2298 dimension is the size of `B`, and it has rank two in dimensions
2299 larger than two. It reduces to the ``JordanSpinEJA`` when `B` is
2300 the identity matrix of order ``n``.
2301
2302 We insist that the one-by-one upper-left identity block of `B` be
2303 passed in as well so that we can be passed a matrix of size zero
2304 to construct a trivial algebra.
2305
2306 SETUP::
2307
2308 sage: from mjo.eja.eja_algebra import (BilinearFormEJA,
2309 ....: JordanSpinEJA)
2310
2311 EXAMPLES:
2312
2313 When no bilinear form is specified, the identity matrix is used,
2314 and the resulting algebra is the Jordan spin algebra::
2315
2316 sage: B = matrix.identity(AA,3)
2317 sage: J0 = BilinearFormEJA(B)
2318 sage: J1 = JordanSpinEJA(3)
2319 sage: J0.multiplication_table() == J0.multiplication_table()
2320 True
2321
2322 An error is raised if the matrix `B` does not correspond to a
2323 positive-definite bilinear form::
2324
2325 sage: B = matrix.random(QQ,2,3)
2326 sage: J = BilinearFormEJA(B)
2327 Traceback (most recent call last):
2328 ...
2329 ValueError: bilinear form is not positive-definite
2330 sage: B = matrix.zero(QQ,3)
2331 sage: J = BilinearFormEJA(B)
2332 Traceback (most recent call last):
2333 ...
2334 ValueError: bilinear form is not positive-definite
2335
2336 TESTS:
2337
2338 We can create a zero-dimensional algebra::
2339
2340 sage: B = matrix.identity(AA,0)
2341 sage: J = BilinearFormEJA(B)
2342 sage: J.basis()
2343 Finite family {}
2344
2345 We can check the multiplication condition given in the Jordan, von
2346 Neumann, and Wigner paper (and also discussed on my "On the
2347 symmetry..." paper). Note that this relies heavily on the standard
2348 choice of basis, as does anything utilizing the bilinear form
2349 matrix. We opt not to orthonormalize the basis, because if we
2350 did, we would have to normalize the `s_{i}` in a similar manner::
2351
2352 sage: set_random_seed()
2353 sage: n = ZZ.random_element(5)
2354 sage: M = matrix.random(QQ, max(0,n-1), algorithm='unimodular')
2355 sage: B11 = matrix.identity(QQ,1)
2356 sage: B22 = M.transpose()*M
2357 sage: B = block_matrix(2,2,[ [B11,0 ],
2358 ....: [0, B22 ] ])
2359 sage: J = BilinearFormEJA(B, orthonormalize=False)
2360 sage: eis = VectorSpace(M.base_ring(), M.ncols()).basis()
2361 sage: V = J.vector_space()
2362 sage: sis = [ J( V([0] + (M.inverse()*ei).list()).column() )
2363 ....: for ei in eis ]
2364 sage: actual = [ sis[i]*sis[j]
2365 ....: for i in range(n-1)
2366 ....: for j in range(n-1) ]
2367 sage: expected = [ J.one() if i == j else J.zero()
2368 ....: for i in range(n-1)
2369 ....: for j in range(n-1) ]
2370 sage: actual == expected
2371 True
2372 """
2373 def __init__(self, B, field=AA, **kwargs):
2374 if not B.is_positive_definite():
2375 raise ValueError("bilinear form is not positive-definite")
2376
2377 n = B.nrows()
2378 V = VectorSpace(field, n)
2379
2380 def inner_product(x,y):
2381 return (B*x).inner_product(y)
2382
2383 def jordan_product(x,y):
2384 x0 = x[0]
2385 xbar = x[1:]
2386 y0 = y[0]
2387 ybar = y[1:]
2388 z0 = inner_product(x,y)
2389 zbar = y0*xbar + x0*ybar
2390 return V([z0] + zbar.list())
2391
2392 super(BilinearFormEJA, self).__init__(field,
2393 V.basis(),
2394 jordan_product,
2395 inner_product,
2396 **kwargs)
2397
2398 # The rank of this algebra is two, unless we're in a
2399 # one-dimensional ambient space (because the rank is bounded
2400 # by the ambient dimension).
2401 self.rank.set_cache(min(n,2))
2402
2403 if n == 0:
2404 self.one.set_cache( self.zero() )
2405 else:
2406 self.one.set_cache( self.monomial(0) )
2407
2408 @staticmethod
2409 def _max_random_instance_size():
2410 r"""
2411 The maximum dimension of a random BilinearFormEJA.
2412 """
2413 return 5
2414
2415 @classmethod
2416 def random_instance(cls, field=AA, **kwargs):
2417 """
2418 Return a random instance of this algebra.
2419 """
2420 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2421 if n.is_zero():
2422 B = matrix.identity(field, n)
2423 return cls(B, field, **kwargs)
2424
2425 B11 = matrix.identity(field,1)
2426 M = matrix.random(field, n-1)
2427 I = matrix.identity(field, n-1)
2428 alpha = field.zero()
2429 while alpha.is_zero():
2430 alpha = field.random_element().abs()
2431 B22 = M.transpose()*M + alpha*I
2432
2433 from sage.matrix.special import block_matrix
2434 B = block_matrix(2,2, [ [B11, ZZ(0) ],
2435 [ZZ(0), B22 ] ])
2436
2437 return cls(B, field, **kwargs)
2438
2439
2440 class JordanSpinEJA(BilinearFormEJA):
2441 """
2442 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2443 with the usual inner product and jordan product ``x*y =
2444 (<x,y>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
2445 the reals.
2446
2447 SETUP::
2448
2449 sage: from mjo.eja.eja_algebra import JordanSpinEJA
2450
2451 EXAMPLES:
2452
2453 This multiplication table can be verified by hand::
2454
2455 sage: J = JordanSpinEJA(4)
2456 sage: e0,e1,e2,e3 = J.gens()
2457 sage: e0*e0
2458 e0
2459 sage: e0*e1
2460 e1
2461 sage: e0*e2
2462 e2
2463 sage: e0*e3
2464 e3
2465 sage: e1*e2
2466 0
2467 sage: e1*e3
2468 0
2469 sage: e2*e3
2470 0
2471
2472 We can change the generator prefix::
2473
2474 sage: JordanSpinEJA(2, prefix='B').gens()
2475 (B0, B1)
2476
2477 TESTS:
2478
2479 Ensure that we have the usual inner product on `R^n`::
2480
2481 sage: set_random_seed()
2482 sage: J = JordanSpinEJA.random_instance()
2483 sage: x,y = J.random_elements(2)
2484 sage: actual = x.inner_product(y)
2485 sage: expected = x.to_vector().inner_product(y.to_vector())
2486 sage: actual == expected
2487 True
2488
2489 """
2490 def __init__(self, n, field=AA, **kwargs):
2491 # This is a special case of the BilinearFormEJA with the identity
2492 # matrix as its bilinear form.
2493 B = matrix.identity(field, n)
2494 super(JordanSpinEJA, self).__init__(B, field, **kwargs)
2495
2496 @staticmethod
2497 def _max_random_instance_size():
2498 r"""
2499 The maximum dimension of a random JordanSpinEJA.
2500 """
2501 return 5
2502
2503 @classmethod
2504 def random_instance(cls, field=AA, **kwargs):
2505 """
2506 Return a random instance of this type of algebra.
2507
2508 Needed here to override the implementation for ``BilinearFormEJA``.
2509 """
2510 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2511 return cls(n, field, **kwargs)
2512
2513
2514 class TrivialEJA(ConcreteEuclideanJordanAlgebra):
2515 """
2516 The trivial Euclidean Jordan algebra consisting of only a zero element.
2517
2518 SETUP::
2519
2520 sage: from mjo.eja.eja_algebra import TrivialEJA
2521
2522 EXAMPLES::
2523
2524 sage: J = TrivialEJA()
2525 sage: J.dimension()
2526 0
2527 sage: J.zero()
2528 0
2529 sage: J.one()
2530 0
2531 sage: 7*J.one()*12*J.one()
2532 0
2533 sage: J.one().inner_product(J.one())
2534 0
2535 sage: J.one().norm()
2536 0
2537 sage: J.one().subalgebra_generated_by()
2538 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
2539 sage: J.rank()
2540 0
2541
2542 """
2543 def __init__(self, field=AA, **kwargs):
2544 jordan_product = lambda x,y: x
2545 inner_product = lambda x,y: field(0)
2546 basis = ()
2547 super(TrivialEJA, self).__init__(field,
2548 basis,
2549 jordan_product,
2550 inner_product,
2551 **kwargs)
2552 # The rank is zero using my definition, namely the dimension of the
2553 # largest subalgebra generated by any element.
2554 self.rank.set_cache(0)
2555 self.one.set_cache( self.zero() )
2556
2557 @classmethod
2558 def random_instance(cls, field=AA, **kwargs):
2559 # We don't take a "size" argument so the superclass method is
2560 # inappropriate for us.
2561 return cls(field, **kwargs)
2562
2563 class DirectSumEJA(FiniteDimensionalEuclideanJordanAlgebra):
2564 r"""
2565 The external (orthogonal) direct sum of two other Euclidean Jordan
2566 algebras. Essentially the Cartesian product of its two factors.
2567 Every Euclidean Jordan algebra decomposes into an orthogonal
2568 direct sum of simple Euclidean Jordan algebras, so no generality
2569 is lost by providing only this construction.
2570
2571 SETUP::
2572
2573 sage: from mjo.eja.eja_algebra import (random_eja,
2574 ....: HadamardEJA,
2575 ....: RealSymmetricEJA,
2576 ....: DirectSumEJA)
2577
2578 EXAMPLES::
2579
2580 sage: J1 = HadamardEJA(2)
2581 sage: J2 = RealSymmetricEJA(3)
2582 sage: J = DirectSumEJA(J1,J2)
2583 sage: J.dimension()
2584 8
2585 sage: J.rank()
2586 5
2587
2588 TESTS:
2589
2590 The external direct sum construction is only valid when the two factors
2591 have the same base ring; an error is raised otherwise::
2592
2593 sage: set_random_seed()
2594 sage: J1 = random_eja(AA)
2595 sage: J2 = random_eja(QQ,orthonormalize=False)
2596 sage: J = DirectSumEJA(J1,J2)
2597 Traceback (most recent call last):
2598 ...
2599 ValueError: algebras must share the same base field
2600
2601 """
2602 def __init__(self, J1, J2, **kwargs):
2603 if J1.base_ring() != J2.base_ring():
2604 raise ValueError("algebras must share the same base field")
2605 field = J1.base_ring()
2606
2607 self._factors = (J1, J2)
2608 n1 = J1.dimension()
2609 n2 = J2.dimension()
2610 n = n1+n2
2611 V = VectorSpace(field, n)
2612 mult_table = [ [ V.zero() for j in range(n) ]
2613 for i in range(n) ]
2614 for i in range(n1):
2615 for j in range(n1):
2616 p = (J1.monomial(i)*J1.monomial(j)).to_vector()
2617 mult_table[i][j] = V(p.list() + [field.zero()]*n2)
2618
2619 for i in range(n2):
2620 for j in range(n2):
2621 p = (J2.monomial(i)*J2.monomial(j)).to_vector()
2622 mult_table[n1+i][n1+j] = V([field.zero()]*n1 + p.list())
2623
2624 # TODO: build the IP table here from the two constituent IP
2625 # matrices (it'll be block diagonal, I think).
2626 ip_table = None
2627 super(DirectSumEJA, self).__init__(field,
2628 mult_table,
2629 ip_table,
2630 check_axioms=False,
2631 **kwargs)
2632 self.rank.set_cache(J1.rank() + J2.rank())
2633
2634
2635 def factors(self):
2636 r"""
2637 Return the pair of this algebra's factors.
2638
2639 SETUP::
2640
2641 sage: from mjo.eja.eja_algebra import (HadamardEJA,
2642 ....: JordanSpinEJA,
2643 ....: DirectSumEJA)
2644
2645 EXAMPLES::
2646
2647 sage: J1 = HadamardEJA(2,QQ)
2648 sage: J2 = JordanSpinEJA(3,QQ)
2649 sage: J = DirectSumEJA(J1,J2)
2650 sage: J.factors()
2651 (Euclidean Jordan algebra of dimension 2 over Rational Field,
2652 Euclidean Jordan algebra of dimension 3 over Rational Field)
2653
2654 """
2655 return self._factors
2656
2657 def projections(self):
2658 r"""
2659 Return a pair of projections onto this algebra's factors.
2660
2661 SETUP::
2662
2663 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
2664 ....: ComplexHermitianEJA,
2665 ....: DirectSumEJA)
2666
2667 EXAMPLES::
2668
2669 sage: J1 = JordanSpinEJA(2)
2670 sage: J2 = ComplexHermitianEJA(2)
2671 sage: J = DirectSumEJA(J1,J2)
2672 sage: (pi_left, pi_right) = J.projections()
2673 sage: J.one().to_vector()
2674 (1, 0, 1, 0, 0, 1)
2675 sage: pi_left(J.one()).to_vector()
2676 (1, 0)
2677 sage: pi_right(J.one()).to_vector()
2678 (1, 0, 0, 1)
2679
2680 """
2681 (J1,J2) = self.factors()
2682 m = J1.dimension()
2683 n = J2.dimension()
2684 V_basis = self.vector_space().basis()
2685 # Need to specify the dimensions explicitly so that we don't
2686 # wind up with a zero-by-zero matrix when we want e.g. a
2687 # zero-by-two matrix (important for composing things).
2688 P1 = matrix(self.base_ring(), m, m+n, V_basis[:m])
2689 P2 = matrix(self.base_ring(), n, m+n, V_basis[m:])
2690 pi_left = FiniteDimensionalEuclideanJordanAlgebraOperator(self,J1,P1)
2691 pi_right = FiniteDimensionalEuclideanJordanAlgebraOperator(self,J2,P2)
2692 return (pi_left, pi_right)
2693
2694 def inclusions(self):
2695 r"""
2696 Return the pair of inclusion maps from our factors into us.
2697
2698 SETUP::
2699
2700 sage: from mjo.eja.eja_algebra import (random_eja,
2701 ....: JordanSpinEJA,
2702 ....: RealSymmetricEJA,
2703 ....: DirectSumEJA)
2704
2705 EXAMPLES::
2706
2707 sage: J1 = JordanSpinEJA(3)
2708 sage: J2 = RealSymmetricEJA(2)
2709 sage: J = DirectSumEJA(J1,J2)
2710 sage: (iota_left, iota_right) = J.inclusions()
2711 sage: iota_left(J1.zero()) == J.zero()
2712 True
2713 sage: iota_right(J2.zero()) == J.zero()
2714 True
2715 sage: J1.one().to_vector()
2716 (1, 0, 0)
2717 sage: iota_left(J1.one()).to_vector()
2718 (1, 0, 0, 0, 0, 0)
2719 sage: J2.one().to_vector()
2720 (1, 0, 1)
2721 sage: iota_right(J2.one()).to_vector()
2722 (0, 0, 0, 1, 0, 1)
2723 sage: J.one().to_vector()
2724 (1, 0, 0, 1, 0, 1)
2725
2726 TESTS:
2727
2728 Composing a projection with the corresponding inclusion should
2729 produce the identity map, and mismatching them should produce
2730 the zero map::
2731
2732 sage: set_random_seed()
2733 sage: J1 = random_eja()
2734 sage: J2 = random_eja()
2735 sage: J = DirectSumEJA(J1,J2)
2736 sage: (iota_left, iota_right) = J.inclusions()
2737 sage: (pi_left, pi_right) = J.projections()
2738 sage: pi_left*iota_left == J1.one().operator()
2739 True
2740 sage: pi_right*iota_right == J2.one().operator()
2741 True
2742 sage: (pi_left*iota_right).is_zero()
2743 True
2744 sage: (pi_right*iota_left).is_zero()
2745 True
2746
2747 """
2748 (J1,J2) = self.factors()
2749 m = J1.dimension()
2750 n = J2.dimension()
2751 V_basis = self.vector_space().basis()
2752 # Need to specify the dimensions explicitly so that we don't
2753 # wind up with a zero-by-zero matrix when we want e.g. a
2754 # two-by-zero matrix (important for composing things).
2755 I1 = matrix.column(self.base_ring(), m, m+n, V_basis[:m])
2756 I2 = matrix.column(self.base_ring(), n, m+n, V_basis[m:])
2757 iota_left = FiniteDimensionalEuclideanJordanAlgebraOperator(J1,self,I1)
2758 iota_right = FiniteDimensionalEuclideanJordanAlgebraOperator(J2,self,I2)
2759 return (iota_left, iota_right)
2760
2761 def inner_product(self, x, y):
2762 r"""
2763 The standard Cartesian inner-product.
2764
2765 We project ``x`` and ``y`` onto our factors, and add up the
2766 inner-products from the subalgebras.
2767
2768 SETUP::
2769
2770
2771 sage: from mjo.eja.eja_algebra import (HadamardEJA,
2772 ....: QuaternionHermitianEJA,
2773 ....: DirectSumEJA)
2774
2775 EXAMPLE::
2776
2777 sage: J1 = HadamardEJA(3,QQ)
2778 sage: J2 = QuaternionHermitianEJA(2,QQ,orthonormalize=False)
2779 sage: J = DirectSumEJA(J1,J2)
2780 sage: x1 = J1.one()
2781 sage: x2 = x1
2782 sage: y1 = J2.one()
2783 sage: y2 = y1
2784 sage: x1.inner_product(x2)
2785 3
2786 sage: y1.inner_product(y2)
2787 2
2788 sage: J.one().inner_product(J.one())
2789 5
2790
2791 """
2792 (pi_left, pi_right) = self.projections()
2793 x1 = pi_left(x)
2794 x2 = pi_right(x)
2795 y1 = pi_left(y)
2796 y2 = pi_right(y)
2797
2798 return (x1.inner_product(y1) + x2.inner_product(y2))
2799
2800
2801
2802 random_eja = ConcreteEuclideanJordanAlgebra.random_instance