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