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