]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: move the TrivialEJA random_instance() into an override method.
[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_random_instance_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 raise NotImplementedError
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 n = ZZ.random_element(cls._max_random_instance_size() + 1)
843 return cls(n, field, **kwargs)
844
845 @cached_method
846 def _charpoly_coefficients(self):
847 r"""
848 The `r` polynomial coefficients of the "characteristic polynomial
849 of" function.
850 """
851 n = self.dimension()
852 var_names = [ "X" + str(z) for z in range(1,n+1) ]
853 R = PolynomialRing(self.base_ring(), var_names)
854 vars = R.gens()
855 F = R.fraction_field()
856
857 def L_x_i_j(i,j):
858 # From a result in my book, these are the entries of the
859 # basis representation of L_x.
860 return sum( vars[k]*self.monomial(k).operator().matrix()[i,j]
861 for k in range(n) )
862
863 L_x = matrix(F, n, n, L_x_i_j)
864
865 r = None
866 if self.rank.is_in_cache():
867 r = self.rank()
868 # There's no need to pad the system with redundant
869 # columns if we *know* they'll be redundant.
870 n = r
871
872 # Compute an extra power in case the rank is equal to
873 # the dimension (otherwise, we would stop at x^(r-1)).
874 x_powers = [ (L_x**k)*self.one().to_vector()
875 for k in range(n+1) ]
876 A = matrix.column(F, x_powers[:n])
877 AE = A.extended_echelon_form()
878 E = AE[:,n:]
879 A_rref = AE[:,:n]
880 if r is None:
881 r = A_rref.rank()
882 b = x_powers[r]
883
884 # The theory says that only the first "r" coefficients are
885 # nonzero, and they actually live in the original polynomial
886 # ring and not the fraction field. We negate them because
887 # in the actual characteristic polynomial, they get moved
888 # to the other side where x^r lives.
889 return -A_rref.solve_right(E*b).change_ring(R)[:r]
890
891 @cached_method
892 def rank(self):
893 r"""
894 Return the rank of this EJA.
895
896 This is a cached method because we know the rank a priori for
897 all of the algebras we can construct. Thus we can avoid the
898 expensive ``_charpoly_coefficients()`` call unless we truly
899 need to compute the whole characteristic polynomial.
900
901 SETUP::
902
903 sage: from mjo.eja.eja_algebra import (HadamardEJA,
904 ....: JordanSpinEJA,
905 ....: RealSymmetricEJA,
906 ....: ComplexHermitianEJA,
907 ....: QuaternionHermitianEJA,
908 ....: random_eja)
909
910 EXAMPLES:
911
912 The rank of the Jordan spin algebra is always two::
913
914 sage: JordanSpinEJA(2).rank()
915 2
916 sage: JordanSpinEJA(3).rank()
917 2
918 sage: JordanSpinEJA(4).rank()
919 2
920
921 The rank of the `n`-by-`n` Hermitian real, complex, or
922 quaternion matrices is `n`::
923
924 sage: RealSymmetricEJA(4).rank()
925 4
926 sage: ComplexHermitianEJA(3).rank()
927 3
928 sage: QuaternionHermitianEJA(2).rank()
929 2
930
931 TESTS:
932
933 Ensure that every EJA that we know how to construct has a
934 positive integer rank, unless the algebra is trivial in
935 which case its rank will be zero::
936
937 sage: set_random_seed()
938 sage: J = random_eja()
939 sage: r = J.rank()
940 sage: r in ZZ
941 True
942 sage: r > 0 or (r == 0 and J.is_trivial())
943 True
944
945 Ensure that computing the rank actually works, since the ranks
946 of all simple algebras are known and will be cached by default::
947
948 sage: J = HadamardEJA(4)
949 sage: J.rank.clear_cache()
950 sage: J.rank()
951 4
952
953 ::
954
955 sage: J = JordanSpinEJA(4)
956 sage: J.rank.clear_cache()
957 sage: J.rank()
958 2
959
960 ::
961
962 sage: J = RealSymmetricEJA(3)
963 sage: J.rank.clear_cache()
964 sage: J.rank()
965 3
966
967 ::
968
969 sage: J = ComplexHermitianEJA(2)
970 sage: J.rank.clear_cache()
971 sage: J.rank()
972 2
973
974 ::
975
976 sage: J = QuaternionHermitianEJA(2)
977 sage: J.rank.clear_cache()
978 sage: J.rank()
979 2
980 """
981 return len(self._charpoly_coefficients())
982
983
984 def vector_space(self):
985 """
986 Return the vector space that underlies this algebra.
987
988 SETUP::
989
990 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
991
992 EXAMPLES::
993
994 sage: J = RealSymmetricEJA(2)
995 sage: J.vector_space()
996 Vector space of dimension 3 over...
997
998 """
999 return self.zero().to_vector().parent().ambient_vector_space()
1000
1001
1002 Element = FiniteDimensionalEuclideanJordanAlgebraElement
1003
1004
1005
1006 def random_eja(field=AA):
1007 """
1008 Return a "random" finite-dimensional Euclidean Jordan Algebra.
1009
1010 SETUP::
1011
1012 sage: from mjo.eja.eja_algebra import random_eja
1013
1014 TESTS::
1015
1016 sage: random_eja()
1017 Euclidean Jordan algebra of dimension...
1018
1019 """
1020 classname = choice([TrivialEJA,
1021 HadamardEJA,
1022 JordanSpinEJA,
1023 RealSymmetricEJA,
1024 ComplexHermitianEJA,
1025 QuaternionHermitianEJA])
1026 return classname.random_instance(field=field)
1027
1028
1029
1030
1031 class RationalBasisEuclideanJordanAlgebra(FiniteDimensionalEuclideanJordanAlgebra):
1032 r"""
1033 Algebras whose basis consists of vectors with rational
1034 entries. Equivalently, algebras whose multiplication tables
1035 contain only rational coefficients.
1036
1037 When an EJA has a basis that can be made rational, we can speed up
1038 the computation of its characteristic polynomial by doing it over
1039 ``QQ``. All of the named EJA constructors that we provide fall
1040 into this category.
1041 """
1042 @cached_method
1043 def _charpoly_coefficients(self):
1044 r"""
1045 Override the parent method with something that tries to compute
1046 over a faster (non-extension) field.
1047
1048 SETUP::
1049
1050 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1051
1052 EXAMPLES:
1053
1054 The base ring of the resulting polynomial coefficients is what
1055 it should be, and not the rationals (unless the algebra was
1056 already over the rationals)::
1057
1058 sage: J = JordanSpinEJA(3)
1059 sage: J._charpoly_coefficients()
1060 (X1^2 - X2^2 - X3^2, -2*X1)
1061 sage: a0 = J._charpoly_coefficients()[0]
1062 sage: J.base_ring()
1063 Algebraic Real Field
1064 sage: a0.base_ring()
1065 Algebraic Real Field
1066
1067 """
1068 if self.base_ring() is QQ:
1069 # There's no need to construct *another* algebra over the
1070 # rationals if this one is already over the rationals.
1071 superclass = super(RationalBasisEuclideanJordanAlgebra, self)
1072 return superclass._charpoly_coefficients()
1073
1074 mult_table = tuple(
1075 map(lambda x: x.to_vector(), ls)
1076 for ls in self._multiplication_table
1077 )
1078
1079 # Do the computation over the rationals. The answer will be
1080 # the same, because our basis coordinates are (essentially)
1081 # rational.
1082 J = FiniteDimensionalEuclideanJordanAlgebra(QQ,
1083 mult_table,
1084 check_field=False,
1085 check_axioms=False)
1086 a = J._charpoly_coefficients()
1087 return tuple(map(lambda x: x.change_ring(self.base_ring()), a))
1088
1089
1090 class MatrixEuclideanJordanAlgebra(FiniteDimensionalEuclideanJordanAlgebra):
1091 @staticmethod
1092 def _max_random_instance_size():
1093 # Play it safe, since this will be squared and the underlying
1094 # field can have dimension 4 (quaternions) too.
1095 return 2
1096
1097 def __init__(self, field, basis, normalize_basis=True, **kwargs):
1098 """
1099 Compared to the superclass constructor, we take a basis instead of
1100 a multiplication table because the latter can be computed in terms
1101 of the former when the product is known (like it is here).
1102 """
1103 # Used in this class's fast _charpoly_coefficients() override.
1104 self._basis_normalizers = None
1105
1106 # We're going to loop through this a few times, so now's a good
1107 # time to ensure that it isn't a generator expression.
1108 basis = tuple(basis)
1109
1110 if len(basis) > 1 and normalize_basis:
1111 # We'll need sqrt(2) to normalize the basis, and this
1112 # winds up in the multiplication table, so the whole
1113 # algebra needs to be over the field extension.
1114 R = PolynomialRing(field, 'z')
1115 z = R.gen()
1116 p = z**2 - 2
1117 if p.is_irreducible():
1118 field = field.extension(p, 'sqrt2', embedding=RLF(2).sqrt())
1119 basis = tuple( s.change_ring(field) for s in basis )
1120 self._basis_normalizers = tuple(
1121 ~(self.natural_inner_product(s,s).sqrt()) for s in basis )
1122 basis = tuple(s*c for (s,c) in zip(basis,self._basis_normalizers))
1123
1124 Qs = self.multiplication_table_from_matrix_basis(basis)
1125
1126 super(MatrixEuclideanJordanAlgebra, self).__init__(field,
1127 Qs,
1128 natural_basis=basis,
1129 **kwargs)
1130
1131
1132 @cached_method
1133 def _charpoly_coefficients(self):
1134 r"""
1135 Override the parent method with something that tries to compute
1136 over a faster (non-extension) field.
1137 """
1138 if self._basis_normalizers is None or self.base_ring() is QQ:
1139 # We didn't normalize, or the basis we started with had
1140 # entries in a nice field already. Just compute the thing.
1141 return super(MatrixEuclideanJordanAlgebra, self)._charpoly_coefficients()
1142
1143 basis = ( (b/n) for (b,n) in zip(self.natural_basis(),
1144 self._basis_normalizers) )
1145
1146 # Do this over the rationals and convert back at the end.
1147 # Only works because we know the entries of the basis are
1148 # integers. The argument ``check_axioms=False`` is required
1149 # because the trace inner-product method for this
1150 # class is a stub and can't actually be checked.
1151 J = MatrixEuclideanJordanAlgebra(QQ,
1152 basis,
1153 normalize_basis=False,
1154 check_field=False,
1155 check_axioms=False)
1156 a = J._charpoly_coefficients()
1157
1158 # Unfortunately, changing the basis does change the
1159 # coefficients of the characteristic polynomial, but since
1160 # these are really the coefficients of the "characteristic
1161 # polynomial of" function, everything is still nice and
1162 # unevaluated. It's therefore "obvious" how scaling the
1163 # basis affects the coordinate variables X1, X2, et
1164 # cetera. Scaling the first basis vector up by "n" adds a
1165 # factor of 1/n into every "X1" term, for example. So here
1166 # we simply undo the basis_normalizer scaling that we
1167 # performed earlier.
1168 #
1169 # The a[0] access here is safe because trivial algebras
1170 # won't have any basis normalizers and therefore won't
1171 # make it to this "else" branch.
1172 XS = a[0].parent().gens()
1173 subs_dict = { XS[i]: self._basis_normalizers[i]*XS[i]
1174 for i in range(len(XS)) }
1175 return tuple( a_i.subs(subs_dict) for a_i in a )
1176
1177
1178 @staticmethod
1179 def multiplication_table_from_matrix_basis(basis):
1180 """
1181 At least three of the five simple Euclidean Jordan algebras have the
1182 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
1183 multiplication on the right is matrix multiplication. Given a basis
1184 for the underlying matrix space, this function returns a
1185 multiplication table (obtained by looping through the basis
1186 elements) for an algebra of those matrices.
1187 """
1188 # In S^2, for example, we nominally have four coordinates even
1189 # though the space is of dimension three only. The vector space V
1190 # is supposed to hold the entire long vector, and the subspace W
1191 # of V will be spanned by the vectors that arise from symmetric
1192 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
1193 if len(basis) == 0:
1194 return []
1195
1196 field = basis[0].base_ring()
1197 dimension = basis[0].nrows()
1198
1199 V = VectorSpace(field, dimension**2)
1200 W = V.span_of_basis( _mat2vec(s) for s in basis )
1201 n = len(basis)
1202 mult_table = [[W.zero() for j in range(n)] for i in range(n)]
1203 for i in range(n):
1204 for j in range(n):
1205 mat_entry = (basis[i]*basis[j] + basis[j]*basis[i])/2
1206 mult_table[i][j] = W.coordinate_vector(_mat2vec(mat_entry))
1207
1208 return mult_table
1209
1210
1211 @staticmethod
1212 def real_embed(M):
1213 """
1214 Embed the matrix ``M`` into a space of real matrices.
1215
1216 The matrix ``M`` can have entries in any field at the moment:
1217 the real numbers, complex numbers, or quaternions. And although
1218 they are not a field, we can probably support octonions at some
1219 point, too. This function returns a real matrix that "acts like"
1220 the original with respect to matrix multiplication; i.e.
1221
1222 real_embed(M*N) = real_embed(M)*real_embed(N)
1223
1224 """
1225 raise NotImplementedError
1226
1227
1228 @staticmethod
1229 def real_unembed(M):
1230 """
1231 The inverse of :meth:`real_embed`.
1232 """
1233 raise NotImplementedError
1234
1235
1236 @classmethod
1237 def natural_inner_product(cls,X,Y):
1238 Xu = cls.real_unembed(X)
1239 Yu = cls.real_unembed(Y)
1240 tr = (Xu*Yu).trace()
1241
1242 try:
1243 # Works in QQ, AA, RDF, et cetera.
1244 return tr.real()
1245 except AttributeError:
1246 # A quaternion doesn't have a real() method, but does
1247 # have coefficient_tuple() method that returns the
1248 # coefficients of 1, i, j, and k -- in that order.
1249 return tr.coefficient_tuple()[0]
1250
1251
1252 class RealMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1253 @staticmethod
1254 def real_embed(M):
1255 """
1256 The identity function, for embedding real matrices into real
1257 matrices.
1258 """
1259 return M
1260
1261 @staticmethod
1262 def real_unembed(M):
1263 """
1264 The identity function, for unembedding real matrices from real
1265 matrices.
1266 """
1267 return M
1268
1269
1270 class RealSymmetricEJA(RealMatrixEuclideanJordanAlgebra):
1271 """
1272 The rank-n simple EJA consisting of real symmetric n-by-n
1273 matrices, the usual symmetric Jordan product, and the trace inner
1274 product. It has dimension `(n^2 + n)/2` over the reals.
1275
1276 SETUP::
1277
1278 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1279
1280 EXAMPLES::
1281
1282 sage: J = RealSymmetricEJA(2)
1283 sage: e0, e1, e2 = J.gens()
1284 sage: e0*e0
1285 e0
1286 sage: e1*e1
1287 1/2*e0 + 1/2*e2
1288 sage: e2*e2
1289 e2
1290
1291 In theory, our "field" can be any subfield of the reals::
1292
1293 sage: RealSymmetricEJA(2, RDF)
1294 Euclidean Jordan algebra of dimension 3 over Real Double Field
1295 sage: RealSymmetricEJA(2, RR)
1296 Euclidean Jordan algebra of dimension 3 over Real Field with
1297 53 bits of precision
1298
1299 TESTS:
1300
1301 The dimension of this algebra is `(n^2 + n) / 2`::
1302
1303 sage: set_random_seed()
1304 sage: n_max = RealSymmetricEJA._max_random_instance_size()
1305 sage: n = ZZ.random_element(1, n_max)
1306 sage: J = RealSymmetricEJA(n)
1307 sage: J.dimension() == (n^2 + n)/2
1308 True
1309
1310 The Jordan multiplication is what we think it is::
1311
1312 sage: set_random_seed()
1313 sage: J = RealSymmetricEJA.random_instance()
1314 sage: x,y = J.random_elements(2)
1315 sage: actual = (x*y).natural_representation()
1316 sage: X = x.natural_representation()
1317 sage: Y = y.natural_representation()
1318 sage: expected = (X*Y + Y*X)/2
1319 sage: actual == expected
1320 True
1321 sage: J(expected) == x*y
1322 True
1323
1324 We can change the generator prefix::
1325
1326 sage: RealSymmetricEJA(3, prefix='q').gens()
1327 (q0, q1, q2, q3, q4, q5)
1328
1329 Our natural basis is normalized with respect to the natural inner
1330 product unless we specify otherwise::
1331
1332 sage: set_random_seed()
1333 sage: J = RealSymmetricEJA.random_instance()
1334 sage: all( b.norm() == 1 for b in J.gens() )
1335 True
1336
1337 Since our natural basis is normalized with respect to the natural
1338 inner product, and since we know that this algebra is an EJA, any
1339 left-multiplication operator's matrix will be symmetric because
1340 natural->EJA basis representation is an isometry and within the EJA
1341 the operator is self-adjoint by the Jordan axiom::
1342
1343 sage: set_random_seed()
1344 sage: x = RealSymmetricEJA.random_instance().random_element()
1345 sage: x.operator().matrix().is_symmetric()
1346 True
1347
1348 We can construct the (trivial) algebra of rank zero::
1349
1350 sage: RealSymmetricEJA(0)
1351 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
1352
1353 """
1354 @classmethod
1355 def _denormalized_basis(cls, n, field):
1356 """
1357 Return a basis for the space of real symmetric n-by-n matrices.
1358
1359 SETUP::
1360
1361 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1362
1363 TESTS::
1364
1365 sage: set_random_seed()
1366 sage: n = ZZ.random_element(1,5)
1367 sage: B = RealSymmetricEJA._denormalized_basis(n,QQ)
1368 sage: all( M.is_symmetric() for M in B)
1369 True
1370
1371 """
1372 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1373 # coordinates.
1374 S = []
1375 for i in range(n):
1376 for j in range(i+1):
1377 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1378 if i == j:
1379 Sij = Eij
1380 else:
1381 Sij = Eij + Eij.transpose()
1382 S.append(Sij)
1383 return S
1384
1385
1386 @staticmethod
1387 def _max_random_instance_size():
1388 return 4 # Dimension 10
1389
1390
1391 def __init__(self, n, field=AA, **kwargs):
1392 basis = self._denormalized_basis(n, field)
1393 super(RealSymmetricEJA, self).__init__(field,
1394 basis,
1395 check_axioms=False,
1396 **kwargs)
1397 self.rank.set_cache(n)
1398
1399
1400 class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1401 @staticmethod
1402 def real_embed(M):
1403 """
1404 Embed the n-by-n complex matrix ``M`` into the space of real
1405 matrices of size 2n-by-2n via the map the sends each entry `z = a +
1406 bi` to the block matrix ``[[a,b],[-b,a]]``.
1407
1408 SETUP::
1409
1410 sage: from mjo.eja.eja_algebra import \
1411 ....: ComplexMatrixEuclideanJordanAlgebra
1412
1413 EXAMPLES::
1414
1415 sage: F = QuadraticField(-1, 'I')
1416 sage: x1 = F(4 - 2*i)
1417 sage: x2 = F(1 + 2*i)
1418 sage: x3 = F(-i)
1419 sage: x4 = F(6)
1420 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
1421 sage: ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1422 [ 4 -2| 1 2]
1423 [ 2 4|-2 1]
1424 [-----+-----]
1425 [ 0 -1| 6 0]
1426 [ 1 0| 0 6]
1427
1428 TESTS:
1429
1430 Embedding is a homomorphism (isomorphism, in fact)::
1431
1432 sage: set_random_seed()
1433 sage: n_max = ComplexMatrixEuclideanJordanAlgebra._max_random_instance_size()
1434 sage: n = ZZ.random_element(n_max)
1435 sage: F = QuadraticField(-1, 'I')
1436 sage: X = random_matrix(F, n)
1437 sage: Y = random_matrix(F, n)
1438 sage: Xe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X)
1439 sage: Ye = ComplexMatrixEuclideanJordanAlgebra.real_embed(Y)
1440 sage: XYe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1441 sage: Xe*Ye == XYe
1442 True
1443
1444 """
1445 n = M.nrows()
1446 if M.ncols() != n:
1447 raise ValueError("the matrix 'M' must be square")
1448
1449 # We don't need any adjoined elements...
1450 field = M.base_ring().base_ring()
1451
1452 blocks = []
1453 for z in M.list():
1454 a = z.list()[0] # real part, I guess
1455 b = z.list()[1] # imag part, I guess
1456 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1457
1458 return matrix.block(field, n, blocks)
1459
1460
1461 @staticmethod
1462 def real_unembed(M):
1463 """
1464 The inverse of _embed_complex_matrix().
1465
1466 SETUP::
1467
1468 sage: from mjo.eja.eja_algebra import \
1469 ....: ComplexMatrixEuclideanJordanAlgebra
1470
1471 EXAMPLES::
1472
1473 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1474 ....: [-2, 1, -4, 3],
1475 ....: [ 9, 10, 11, 12],
1476 ....: [-10, 9, -12, 11] ])
1477 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(A)
1478 [ 2*I + 1 4*I + 3]
1479 [ 10*I + 9 12*I + 11]
1480
1481 TESTS:
1482
1483 Unembedding is the inverse of embedding::
1484
1485 sage: set_random_seed()
1486 sage: F = QuadraticField(-1, 'I')
1487 sage: M = random_matrix(F, 3)
1488 sage: Me = ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1489 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1490 True
1491
1492 """
1493 n = ZZ(M.nrows())
1494 if M.ncols() != n:
1495 raise ValueError("the matrix 'M' must be square")
1496 if not n.mod(2).is_zero():
1497 raise ValueError("the matrix 'M' must be a complex embedding")
1498
1499 # If "M" was normalized, its base ring might have roots
1500 # adjoined and they can stick around after unembedding.
1501 field = M.base_ring()
1502 R = PolynomialRing(field, 'z')
1503 z = R.gen()
1504 if field is AA:
1505 # Sage doesn't know how to embed AA into QQbar, i.e. how
1506 # to adjoin sqrt(-1) to AA.
1507 F = QQbar
1508 else:
1509 F = field.extension(z**2 + 1, 'I', embedding=CLF(-1).sqrt())
1510 i = F.gen()
1511
1512 # Go top-left to bottom-right (reading order), converting every
1513 # 2-by-2 block we see to a single complex element.
1514 elements = []
1515 for k in range(n/2):
1516 for j in range(n/2):
1517 submat = M[2*k:2*k+2,2*j:2*j+2]
1518 if submat[0,0] != submat[1,1]:
1519 raise ValueError('bad on-diagonal submatrix')
1520 if submat[0,1] != -submat[1,0]:
1521 raise ValueError('bad off-diagonal submatrix')
1522 z = submat[0,0] + submat[0,1]*i
1523 elements.append(z)
1524
1525 return matrix(F, n/2, elements)
1526
1527
1528 @classmethod
1529 def natural_inner_product(cls,X,Y):
1530 """
1531 Compute a natural inner product in this algebra directly from
1532 its real embedding.
1533
1534 SETUP::
1535
1536 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1537
1538 TESTS:
1539
1540 This gives the same answer as the slow, default method implemented
1541 in :class:`MatrixEuclideanJordanAlgebra`::
1542
1543 sage: set_random_seed()
1544 sage: J = ComplexHermitianEJA.random_instance()
1545 sage: x,y = J.random_elements(2)
1546 sage: Xe = x.natural_representation()
1547 sage: Ye = y.natural_representation()
1548 sage: X = ComplexHermitianEJA.real_unembed(Xe)
1549 sage: Y = ComplexHermitianEJA.real_unembed(Ye)
1550 sage: expected = (X*Y).trace().real()
1551 sage: actual = ComplexHermitianEJA.natural_inner_product(Xe,Ye)
1552 sage: actual == expected
1553 True
1554
1555 """
1556 return RealMatrixEuclideanJordanAlgebra.natural_inner_product(X,Y)/2
1557
1558
1559 class ComplexHermitianEJA(ComplexMatrixEuclideanJordanAlgebra):
1560 """
1561 The rank-n simple EJA consisting of complex Hermitian n-by-n
1562 matrices over the real numbers, the usual symmetric Jordan product,
1563 and the real-part-of-trace inner product. It has dimension `n^2` over
1564 the reals.
1565
1566 SETUP::
1567
1568 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1569
1570 EXAMPLES:
1571
1572 In theory, our "field" can be any subfield of the reals::
1573
1574 sage: ComplexHermitianEJA(2, RDF)
1575 Euclidean Jordan algebra of dimension 4 over Real Double Field
1576 sage: ComplexHermitianEJA(2, RR)
1577 Euclidean Jordan algebra of dimension 4 over Real Field with
1578 53 bits of precision
1579
1580 TESTS:
1581
1582 The dimension of this algebra is `n^2`::
1583
1584 sage: set_random_seed()
1585 sage: n_max = ComplexHermitianEJA._max_random_instance_size()
1586 sage: n = ZZ.random_element(1, n_max)
1587 sage: J = ComplexHermitianEJA(n)
1588 sage: J.dimension() == n^2
1589 True
1590
1591 The Jordan multiplication is what we think it is::
1592
1593 sage: set_random_seed()
1594 sage: J = ComplexHermitianEJA.random_instance()
1595 sage: x,y = J.random_elements(2)
1596 sage: actual = (x*y).natural_representation()
1597 sage: X = x.natural_representation()
1598 sage: Y = y.natural_representation()
1599 sage: expected = (X*Y + Y*X)/2
1600 sage: actual == expected
1601 True
1602 sage: J(expected) == x*y
1603 True
1604
1605 We can change the generator prefix::
1606
1607 sage: ComplexHermitianEJA(2, prefix='z').gens()
1608 (z0, z1, z2, z3)
1609
1610 Our natural basis is normalized with respect to the natural inner
1611 product unless we specify otherwise::
1612
1613 sage: set_random_seed()
1614 sage: J = ComplexHermitianEJA.random_instance()
1615 sage: all( b.norm() == 1 for b in J.gens() )
1616 True
1617
1618 Since our natural basis is normalized with respect to the natural
1619 inner product, and since we know that this algebra is an EJA, any
1620 left-multiplication operator's matrix will be symmetric because
1621 natural->EJA basis representation is an isometry and within the EJA
1622 the operator is self-adjoint by the Jordan axiom::
1623
1624 sage: set_random_seed()
1625 sage: x = ComplexHermitianEJA.random_instance().random_element()
1626 sage: x.operator().matrix().is_symmetric()
1627 True
1628
1629 We can construct the (trivial) algebra of rank zero::
1630
1631 sage: ComplexHermitianEJA(0)
1632 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
1633
1634 """
1635
1636 @classmethod
1637 def _denormalized_basis(cls, n, field):
1638 """
1639 Returns a basis for the space of complex Hermitian n-by-n matrices.
1640
1641 Why do we embed these? Basically, because all of numerical linear
1642 algebra assumes that you're working with vectors consisting of `n`
1643 entries from a field and scalars from the same field. There's no way
1644 to tell SageMath that (for example) the vectors contain complex
1645 numbers, while the scalar field is real.
1646
1647 SETUP::
1648
1649 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1650
1651 TESTS::
1652
1653 sage: set_random_seed()
1654 sage: n = ZZ.random_element(1,5)
1655 sage: field = QuadraticField(2, 'sqrt2')
1656 sage: B = ComplexHermitianEJA._denormalized_basis(n, field)
1657 sage: all( M.is_symmetric() for M in B)
1658 True
1659
1660 """
1661 R = PolynomialRing(field, 'z')
1662 z = R.gen()
1663 F = field.extension(z**2 + 1, 'I')
1664 I = F.gen()
1665
1666 # This is like the symmetric case, but we need to be careful:
1667 #
1668 # * We want conjugate-symmetry, not just symmetry.
1669 # * The diagonal will (as a result) be real.
1670 #
1671 S = []
1672 for i in range(n):
1673 for j in range(i+1):
1674 Eij = matrix(F, n, lambda k,l: k==i and l==j)
1675 if i == j:
1676 Sij = cls.real_embed(Eij)
1677 S.append(Sij)
1678 else:
1679 # The second one has a minus because it's conjugated.
1680 Sij_real = cls.real_embed(Eij + Eij.transpose())
1681 S.append(Sij_real)
1682 Sij_imag = cls.real_embed(I*Eij - I*Eij.transpose())
1683 S.append(Sij_imag)
1684
1685 # Since we embedded these, we can drop back to the "field" that we
1686 # started with instead of the complex extension "F".
1687 return ( s.change_ring(field) for s in S )
1688
1689
1690 def __init__(self, n, field=AA, **kwargs):
1691 basis = self._denormalized_basis(n,field)
1692 super(ComplexHermitianEJA,self).__init__(field,
1693 basis,
1694 check_axioms=False,
1695 **kwargs)
1696 self.rank.set_cache(n)
1697
1698
1699 class QuaternionMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1700 @staticmethod
1701 def real_embed(M):
1702 """
1703 Embed the n-by-n quaternion matrix ``M`` into the space of real
1704 matrices of size 4n-by-4n by first sending each quaternion entry `z
1705 = a + bi + cj + dk` to the block-complex matrix ``[[a + bi,
1706 c+di],[-c + di, a-bi]]`, and then embedding those into a real
1707 matrix.
1708
1709 SETUP::
1710
1711 sage: from mjo.eja.eja_algebra import \
1712 ....: QuaternionMatrixEuclideanJordanAlgebra
1713
1714 EXAMPLES::
1715
1716 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1717 sage: i,j,k = Q.gens()
1718 sage: x = 1 + 2*i + 3*j + 4*k
1719 sage: M = matrix(Q, 1, [[x]])
1720 sage: QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1721 [ 1 2 3 4]
1722 [-2 1 -4 3]
1723 [-3 4 1 -2]
1724 [-4 -3 2 1]
1725
1726 Embedding is a homomorphism (isomorphism, in fact)::
1727
1728 sage: set_random_seed()
1729 sage: n_max = QuaternionMatrixEuclideanJordanAlgebra._max_random_instance_size()
1730 sage: n = ZZ.random_element(n_max)
1731 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1732 sage: X = random_matrix(Q, n)
1733 sage: Y = random_matrix(Q, n)
1734 sage: Xe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X)
1735 sage: Ye = QuaternionMatrixEuclideanJordanAlgebra.real_embed(Y)
1736 sage: XYe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1737 sage: Xe*Ye == XYe
1738 True
1739
1740 """
1741 quaternions = M.base_ring()
1742 n = M.nrows()
1743 if M.ncols() != n:
1744 raise ValueError("the matrix 'M' must be square")
1745
1746 F = QuadraticField(-1, 'I')
1747 i = F.gen()
1748
1749 blocks = []
1750 for z in M.list():
1751 t = z.coefficient_tuple()
1752 a = t[0]
1753 b = t[1]
1754 c = t[2]
1755 d = t[3]
1756 cplxM = matrix(F, 2, [[ a + b*i, c + d*i],
1757 [-c + d*i, a - b*i]])
1758 realM = ComplexMatrixEuclideanJordanAlgebra.real_embed(cplxM)
1759 blocks.append(realM)
1760
1761 # We should have real entries by now, so use the realest field
1762 # we've got for the return value.
1763 return matrix.block(quaternions.base_ring(), n, blocks)
1764
1765
1766
1767 @staticmethod
1768 def real_unembed(M):
1769 """
1770 The inverse of _embed_quaternion_matrix().
1771
1772 SETUP::
1773
1774 sage: from mjo.eja.eja_algebra import \
1775 ....: QuaternionMatrixEuclideanJordanAlgebra
1776
1777 EXAMPLES::
1778
1779 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
1780 ....: [-2, 1, -4, 3],
1781 ....: [-3, 4, 1, -2],
1782 ....: [-4, -3, 2, 1]])
1783 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(M)
1784 [1 + 2*i + 3*j + 4*k]
1785
1786 TESTS:
1787
1788 Unembedding is the inverse of embedding::
1789
1790 sage: set_random_seed()
1791 sage: Q = QuaternionAlgebra(QQ, -1, -1)
1792 sage: M = random_matrix(Q, 3)
1793 sage: Me = QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1794 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1795 True
1796
1797 """
1798 n = ZZ(M.nrows())
1799 if M.ncols() != n:
1800 raise ValueError("the matrix 'M' must be square")
1801 if not n.mod(4).is_zero():
1802 raise ValueError("the matrix 'M' must be a quaternion embedding")
1803
1804 # Use the base ring of the matrix to ensure that its entries can be
1805 # multiplied by elements of the quaternion algebra.
1806 field = M.base_ring()
1807 Q = QuaternionAlgebra(field,-1,-1)
1808 i,j,k = Q.gens()
1809
1810 # Go top-left to bottom-right (reading order), converting every
1811 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
1812 # quaternion block.
1813 elements = []
1814 for l in range(n/4):
1815 for m in range(n/4):
1816 submat = ComplexMatrixEuclideanJordanAlgebra.real_unembed(
1817 M[4*l:4*l+4,4*m:4*m+4] )
1818 if submat[0,0] != submat[1,1].conjugate():
1819 raise ValueError('bad on-diagonal submatrix')
1820 if submat[0,1] != -submat[1,0].conjugate():
1821 raise ValueError('bad off-diagonal submatrix')
1822 z = submat[0,0].real()
1823 z += submat[0,0].imag()*i
1824 z += submat[0,1].real()*j
1825 z += submat[0,1].imag()*k
1826 elements.append(z)
1827
1828 return matrix(Q, n/4, elements)
1829
1830
1831 @classmethod
1832 def natural_inner_product(cls,X,Y):
1833 """
1834 Compute a natural inner product in this algebra directly from
1835 its real embedding.
1836
1837 SETUP::
1838
1839 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1840
1841 TESTS:
1842
1843 This gives the same answer as the slow, default method implemented
1844 in :class:`MatrixEuclideanJordanAlgebra`::
1845
1846 sage: set_random_seed()
1847 sage: J = QuaternionHermitianEJA.random_instance()
1848 sage: x,y = J.random_elements(2)
1849 sage: Xe = x.natural_representation()
1850 sage: Ye = y.natural_representation()
1851 sage: X = QuaternionHermitianEJA.real_unembed(Xe)
1852 sage: Y = QuaternionHermitianEJA.real_unembed(Ye)
1853 sage: expected = (X*Y).trace().coefficient_tuple()[0]
1854 sage: actual = QuaternionHermitianEJA.natural_inner_product(Xe,Ye)
1855 sage: actual == expected
1856 True
1857
1858 """
1859 return RealMatrixEuclideanJordanAlgebra.natural_inner_product(X,Y)/4
1860
1861
1862 class QuaternionHermitianEJA(QuaternionMatrixEuclideanJordanAlgebra):
1863 """
1864 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
1865 matrices, the usual symmetric Jordan product, and the
1866 real-part-of-trace inner product. It has dimension `2n^2 - n` over
1867 the reals.
1868
1869 SETUP::
1870
1871 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1872
1873 EXAMPLES:
1874
1875 In theory, our "field" can be any subfield of the reals::
1876
1877 sage: QuaternionHermitianEJA(2, RDF)
1878 Euclidean Jordan algebra of dimension 6 over Real Double Field
1879 sage: QuaternionHermitianEJA(2, RR)
1880 Euclidean Jordan algebra of dimension 6 over Real Field with
1881 53 bits of precision
1882
1883 TESTS:
1884
1885 The dimension of this algebra is `2*n^2 - n`::
1886
1887 sage: set_random_seed()
1888 sage: n_max = QuaternionHermitianEJA._max_random_instance_size()
1889 sage: n = ZZ.random_element(1, n_max)
1890 sage: J = QuaternionHermitianEJA(n)
1891 sage: J.dimension() == 2*(n^2) - n
1892 True
1893
1894 The Jordan multiplication is what we think it is::
1895
1896 sage: set_random_seed()
1897 sage: J = QuaternionHermitianEJA.random_instance()
1898 sage: x,y = J.random_elements(2)
1899 sage: actual = (x*y).natural_representation()
1900 sage: X = x.natural_representation()
1901 sage: Y = y.natural_representation()
1902 sage: expected = (X*Y + Y*X)/2
1903 sage: actual == expected
1904 True
1905 sage: J(expected) == x*y
1906 True
1907
1908 We can change the generator prefix::
1909
1910 sage: QuaternionHermitianEJA(2, prefix='a').gens()
1911 (a0, a1, a2, a3, a4, a5)
1912
1913 Our natural basis is normalized with respect to the natural inner
1914 product unless we specify otherwise::
1915
1916 sage: set_random_seed()
1917 sage: J = QuaternionHermitianEJA.random_instance()
1918 sage: all( b.norm() == 1 for b in J.gens() )
1919 True
1920
1921 Since our natural basis is normalized with respect to the natural
1922 inner product, and since we know that this algebra is an EJA, any
1923 left-multiplication operator's matrix will be symmetric because
1924 natural->EJA basis representation is an isometry and within the EJA
1925 the operator is self-adjoint by the Jordan axiom::
1926
1927 sage: set_random_seed()
1928 sage: x = QuaternionHermitianEJA.random_instance().random_element()
1929 sage: x.operator().matrix().is_symmetric()
1930 True
1931
1932 We can construct the (trivial) algebra of rank zero::
1933
1934 sage: QuaternionHermitianEJA(0)
1935 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
1936
1937 """
1938 @classmethod
1939 def _denormalized_basis(cls, n, field):
1940 """
1941 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
1942
1943 Why do we embed these? Basically, because all of numerical
1944 linear algebra assumes that you're working with vectors consisting
1945 of `n` entries from a field and scalars from the same field. There's
1946 no way to tell SageMath that (for example) the vectors contain
1947 complex numbers, while the scalar field is real.
1948
1949 SETUP::
1950
1951 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1952
1953 TESTS::
1954
1955 sage: set_random_seed()
1956 sage: n = ZZ.random_element(1,5)
1957 sage: B = QuaternionHermitianEJA._denormalized_basis(n,QQ)
1958 sage: all( M.is_symmetric() for M in B )
1959 True
1960
1961 """
1962 Q = QuaternionAlgebra(QQ,-1,-1)
1963 I,J,K = Q.gens()
1964
1965 # This is like the symmetric case, but we need to be careful:
1966 #
1967 # * We want conjugate-symmetry, not just symmetry.
1968 # * The diagonal will (as a result) be real.
1969 #
1970 S = []
1971 for i in range(n):
1972 for j in range(i+1):
1973 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
1974 if i == j:
1975 Sij = cls.real_embed(Eij)
1976 S.append(Sij)
1977 else:
1978 # The second, third, and fourth ones have a minus
1979 # because they're conjugated.
1980 Sij_real = cls.real_embed(Eij + Eij.transpose())
1981 S.append(Sij_real)
1982 Sij_I = cls.real_embed(I*Eij - I*Eij.transpose())
1983 S.append(Sij_I)
1984 Sij_J = cls.real_embed(J*Eij - J*Eij.transpose())
1985 S.append(Sij_J)
1986 Sij_K = cls.real_embed(K*Eij - K*Eij.transpose())
1987 S.append(Sij_K)
1988
1989 # Since we embedded these, we can drop back to the "field" that we
1990 # started with instead of the quaternion algebra "Q".
1991 return ( s.change_ring(field) for s in S )
1992
1993
1994 def __init__(self, n, field=AA, **kwargs):
1995 basis = self._denormalized_basis(n,field)
1996 super(QuaternionHermitianEJA,self).__init__(field,
1997 basis,
1998 check_axioms=False,
1999 **kwargs)
2000 self.rank.set_cache(n)
2001
2002
2003 class HadamardEJA(RationalBasisEuclideanJordanAlgebra):
2004 """
2005 Return the Euclidean Jordan Algebra corresponding to the set
2006 `R^n` under the Hadamard product.
2007
2008 Note: this is nothing more than the Cartesian product of ``n``
2009 copies of the spin algebra. Once Cartesian product algebras
2010 are implemented, this can go.
2011
2012 SETUP::
2013
2014 sage: from mjo.eja.eja_algebra import HadamardEJA
2015
2016 EXAMPLES:
2017
2018 This multiplication table can be verified by hand::
2019
2020 sage: J = HadamardEJA(3)
2021 sage: e0,e1,e2 = J.gens()
2022 sage: e0*e0
2023 e0
2024 sage: e0*e1
2025 0
2026 sage: e0*e2
2027 0
2028 sage: e1*e1
2029 e1
2030 sage: e1*e2
2031 0
2032 sage: e2*e2
2033 e2
2034
2035 TESTS:
2036
2037 We can change the generator prefix::
2038
2039 sage: HadamardEJA(3, prefix='r').gens()
2040 (r0, r1, r2)
2041
2042 """
2043 def __init__(self, n, field=AA, **kwargs):
2044 V = VectorSpace(field, n)
2045 mult_table = [ [ V.gen(i)*(i == j) for j in range(n) ]
2046 for i in range(n) ]
2047
2048 super(HadamardEJA, self).__init__(field,
2049 mult_table,
2050 check_axioms=False,
2051 **kwargs)
2052 self.rank.set_cache(n)
2053
2054 @staticmethod
2055 def _max_random_instance_size():
2056 return 5
2057
2058 def inner_product(self, x, y):
2059 """
2060 Faster to reimplement than to use natural representations.
2061
2062 SETUP::
2063
2064 sage: from mjo.eja.eja_algebra import HadamardEJA
2065
2066 TESTS:
2067
2068 Ensure that this is the usual inner product for the algebras
2069 over `R^n`::
2070
2071 sage: set_random_seed()
2072 sage: J = HadamardEJA.random_instance()
2073 sage: x,y = J.random_elements(2)
2074 sage: X = x.natural_representation()
2075 sage: Y = y.natural_representation()
2076 sage: x.inner_product(y) == J.natural_inner_product(X,Y)
2077 True
2078
2079 """
2080 return x.to_vector().inner_product(y.to_vector())
2081
2082
2083 class BilinearFormEJA(RationalBasisEuclideanJordanAlgebra):
2084 r"""
2085 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2086 with the half-trace inner product and jordan product ``x*y =
2087 (x0*y0 + <B*x_bar,y_bar>, x0*y_bar + y0*x_bar)`` where ``B`` is a
2088 symmetric positive-definite "bilinear form" matrix. It has
2089 dimension `n` over the reals, and reduces to the ``JordanSpinEJA``
2090 when ``B`` is the identity matrix of order ``n-1``.
2091
2092 SETUP::
2093
2094 sage: from mjo.eja.eja_algebra import (BilinearFormEJA,
2095 ....: JordanSpinEJA)
2096
2097 EXAMPLES:
2098
2099 When no bilinear form is specified, the identity matrix is used,
2100 and the resulting algebra is the Jordan spin algebra::
2101
2102 sage: J0 = BilinearFormEJA(3)
2103 sage: J1 = JordanSpinEJA(3)
2104 sage: J0.multiplication_table() == J0.multiplication_table()
2105 True
2106
2107 TESTS:
2108
2109 We can create a zero-dimensional algebra::
2110
2111 sage: J = BilinearFormEJA(0)
2112 sage: J.basis()
2113 Finite family {}
2114
2115 We can check the multiplication condition given in the Jordan, von
2116 Neumann, and Wigner paper (and also discussed on my "On the
2117 symmetry..." paper). Note that this relies heavily on the standard
2118 choice of basis, as does anything utilizing the bilinear form matrix::
2119
2120 sage: set_random_seed()
2121 sage: n = ZZ.random_element(5)
2122 sage: M = matrix.random(QQ, max(0,n-1), algorithm='unimodular')
2123 sage: B = M.transpose()*M
2124 sage: J = BilinearFormEJA(n, B=B)
2125 sage: eis = VectorSpace(M.base_ring(), M.ncols()).basis()
2126 sage: V = J.vector_space()
2127 sage: sis = [ J.from_vector(V([0] + (M.inverse()*ei).list()))
2128 ....: for ei in eis ]
2129 sage: actual = [ sis[i]*sis[j]
2130 ....: for i in range(n-1)
2131 ....: for j in range(n-1) ]
2132 sage: expected = [ J.one() if i == j else J.zero()
2133 ....: for i in range(n-1)
2134 ....: for j in range(n-1) ]
2135 sage: actual == expected
2136 True
2137 """
2138 def __init__(self, n, field=AA, B=None, **kwargs):
2139 if B is None:
2140 self._B = matrix.identity(field, max(0,n-1))
2141 else:
2142 self._B = B
2143
2144 V = VectorSpace(field, n)
2145 mult_table = [[V.zero() for j in range(n)] for i in range(n)]
2146 for i in range(n):
2147 for j in range(n):
2148 x = V.gen(i)
2149 y = V.gen(j)
2150 x0 = x[0]
2151 xbar = x[1:]
2152 y0 = y[0]
2153 ybar = y[1:]
2154 z0 = x0*y0 + (self._B*xbar).inner_product(ybar)
2155 zbar = y0*xbar + x0*ybar
2156 z = V([z0] + zbar.list())
2157 mult_table[i][j] = z
2158
2159 # The rank of this algebra is two, unless we're in a
2160 # one-dimensional ambient space (because the rank is bounded
2161 # by the ambient dimension).
2162 super(BilinearFormEJA, self).__init__(field,
2163 mult_table,
2164 check_axioms=False,
2165 **kwargs)
2166 self.rank.set_cache(min(n,2))
2167
2168 @staticmethod
2169 def _max_random_instance_size():
2170 return 5
2171
2172 def inner_product(self, x, y):
2173 r"""
2174 Half of the trace inner product.
2175
2176 This is defined so that the special case of the Jordan spin
2177 algebra gets the usual inner product.
2178
2179 SETUP::
2180
2181 sage: from mjo.eja.eja_algebra import BilinearFormEJA
2182
2183 TESTS:
2184
2185 Ensure that this is one-half of the trace inner-product when
2186 the algebra isn't just the reals (when ``n`` isn't one). This
2187 is in Faraut and Koranyi, and also my "On the symmetry..."
2188 paper::
2189
2190 sage: set_random_seed()
2191 sage: n = ZZ.random_element(2,5)
2192 sage: M = matrix.random(QQ, max(0,n-1), algorithm='unimodular')
2193 sage: B = M.transpose()*M
2194 sage: J = BilinearFormEJA(n, B=B)
2195 sage: x = J.random_element()
2196 sage: y = J.random_element()
2197 sage: x.inner_product(y) == (x*y).trace()/2
2198 True
2199
2200 """
2201 xvec = x.to_vector()
2202 xbar = xvec[1:]
2203 yvec = y.to_vector()
2204 ybar = yvec[1:]
2205 return x[0]*y[0] + (self._B*xbar).inner_product(ybar)
2206
2207
2208 class JordanSpinEJA(BilinearFormEJA):
2209 """
2210 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2211 with the usual inner product and jordan product ``x*y =
2212 (<x,y>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
2213 the reals.
2214
2215 SETUP::
2216
2217 sage: from mjo.eja.eja_algebra import JordanSpinEJA
2218
2219 EXAMPLES:
2220
2221 This multiplication table can be verified by hand::
2222
2223 sage: J = JordanSpinEJA(4)
2224 sage: e0,e1,e2,e3 = J.gens()
2225 sage: e0*e0
2226 e0
2227 sage: e0*e1
2228 e1
2229 sage: e0*e2
2230 e2
2231 sage: e0*e3
2232 e3
2233 sage: e1*e2
2234 0
2235 sage: e1*e3
2236 0
2237 sage: e2*e3
2238 0
2239
2240 We can change the generator prefix::
2241
2242 sage: JordanSpinEJA(2, prefix='B').gens()
2243 (B0, B1)
2244
2245 TESTS:
2246
2247 Ensure that we have the usual inner product on `R^n`::
2248
2249 sage: set_random_seed()
2250 sage: J = JordanSpinEJA.random_instance()
2251 sage: x,y = J.random_elements(2)
2252 sage: X = x.natural_representation()
2253 sage: Y = y.natural_representation()
2254 sage: x.inner_product(y) == J.natural_inner_product(X,Y)
2255 True
2256
2257 """
2258 def __init__(self, n, field=AA, **kwargs):
2259 # This is a special case of the BilinearFormEJA with the identity
2260 # matrix as its bilinear form.
2261 super(JordanSpinEJA, self).__init__(n, field, **kwargs)
2262
2263
2264 class TrivialEJA(FiniteDimensionalEuclideanJordanAlgebra):
2265 """
2266 The trivial Euclidean Jordan algebra consisting of only a zero element.
2267
2268 SETUP::
2269
2270 sage: from mjo.eja.eja_algebra import TrivialEJA
2271
2272 EXAMPLES::
2273
2274 sage: J = TrivialEJA()
2275 sage: J.dimension()
2276 0
2277 sage: J.zero()
2278 0
2279 sage: J.one()
2280 0
2281 sage: 7*J.one()*12*J.one()
2282 0
2283 sage: J.one().inner_product(J.one())
2284 0
2285 sage: J.one().norm()
2286 0
2287 sage: J.one().subalgebra_generated_by()
2288 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
2289 sage: J.rank()
2290 0
2291
2292 """
2293 def __init__(self, field=AA, **kwargs):
2294 mult_table = []
2295 super(TrivialEJA, self).__init__(field,
2296 mult_table,
2297 check_axioms=False,
2298 **kwargs)
2299 # The rank is zero using my definition, namely the dimension of the
2300 # largest subalgebra generated by any element.
2301 self.rank.set_cache(0)
2302
2303 @classmethod
2304 def random_instance(cls, field=AA, **kwargs):
2305 # We don't take a "size" argument so the superclass method is
2306 # inappropriate for us.
2307 return cls(field, **kwargs)
2308
2309 class DirectSumEJA(FiniteDimensionalEuclideanJordanAlgebra):
2310 r"""
2311 The external (orthogonal) direct sum of two other Euclidean Jordan
2312 algebras. Essentially the Cartesian product of its two factors.
2313 Every Euclidean Jordan algebra decomposes into an orthogonal
2314 direct sum of simple Euclidean Jordan algebras, so no generality
2315 is lost by providing only this construction.
2316
2317 SETUP::
2318
2319 sage: from mjo.eja.eja_algebra import (HadamardEJA,
2320 ....: RealSymmetricEJA,
2321 ....: DirectSumEJA)
2322
2323 EXAMPLES::
2324
2325 sage: J1 = HadamardEJA(2)
2326 sage: J2 = RealSymmetricEJA(3)
2327 sage: J = DirectSumEJA(J1,J2)
2328 sage: J.dimension()
2329 8
2330 sage: J.rank()
2331 5
2332
2333 """
2334 def __init__(self, J1, J2, field=AA, **kwargs):
2335 self._factors = (J1, J2)
2336 n1 = J1.dimension()
2337 n2 = J2.dimension()
2338 n = n1+n2
2339 V = VectorSpace(field, n)
2340 mult_table = [ [ V.zero() for j in range(n) ]
2341 for i in range(n) ]
2342 for i in range(n1):
2343 for j in range(n1):
2344 p = (J1.monomial(i)*J1.monomial(j)).to_vector()
2345 mult_table[i][j] = V(p.list() + [field.zero()]*n2)
2346
2347 for i in range(n2):
2348 for j in range(n2):
2349 p = (J2.monomial(i)*J2.monomial(j)).to_vector()
2350 mult_table[n1+i][n1+j] = V([field.zero()]*n1 + p.list())
2351
2352 super(DirectSumEJA, self).__init__(field,
2353 mult_table,
2354 check_axioms=False,
2355 **kwargs)
2356 self.rank.set_cache(J1.rank() + J2.rank())
2357
2358
2359 def factors(self):
2360 r"""
2361 Return the pair of this algebra's factors.
2362
2363 SETUP::
2364
2365 sage: from mjo.eja.eja_algebra import (HadamardEJA,
2366 ....: JordanSpinEJA,
2367 ....: DirectSumEJA)
2368
2369 EXAMPLES::
2370
2371 sage: J1 = HadamardEJA(2,QQ)
2372 sage: J2 = JordanSpinEJA(3,QQ)
2373 sage: J = DirectSumEJA(J1,J2)
2374 sage: J.factors()
2375 (Euclidean Jordan algebra of dimension 2 over Rational Field,
2376 Euclidean Jordan algebra of dimension 3 over Rational Field)
2377
2378 """
2379 return self._factors
2380
2381 def projections(self):
2382 r"""
2383 Return a pair of projections onto this algebra's factors.
2384
2385 SETUP::
2386
2387 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
2388 ....: ComplexHermitianEJA,
2389 ....: DirectSumEJA)
2390
2391 EXAMPLES::
2392
2393 sage: J1 = JordanSpinEJA(2)
2394 sage: J2 = ComplexHermitianEJA(2)
2395 sage: J = DirectSumEJA(J1,J2)
2396 sage: (pi_left, pi_right) = J.projections()
2397 sage: J.one().to_vector()
2398 (1, 0, 1, 0, 0, 1)
2399 sage: pi_left(J.one()).to_vector()
2400 (1, 0)
2401 sage: pi_right(J.one()).to_vector()
2402 (1, 0, 0, 1)
2403
2404 """
2405 (J1,J2) = self.factors()
2406 n = J1.dimension()
2407 pi_left = lambda x: J1.from_vector(x.to_vector()[:n])
2408 pi_right = lambda x: J2.from_vector(x.to_vector()[n:])
2409 return (pi_left, pi_right)
2410
2411 def inclusions(self):
2412 r"""
2413 Return the pair of inclusion maps from our factors into us.
2414
2415 SETUP::
2416
2417 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
2418 ....: RealSymmetricEJA,
2419 ....: DirectSumEJA)
2420
2421 EXAMPLES::
2422
2423 sage: J1 = JordanSpinEJA(3)
2424 sage: J2 = RealSymmetricEJA(2)
2425 sage: J = DirectSumEJA(J1,J2)
2426 sage: (iota_left, iota_right) = J.inclusions()
2427 sage: iota_left(J1.zero()) == J.zero()
2428 True
2429 sage: iota_right(J2.zero()) == J.zero()
2430 True
2431 sage: J1.one().to_vector()
2432 (1, 0, 0)
2433 sage: iota_left(J1.one()).to_vector()
2434 (1, 0, 0, 0, 0, 0)
2435 sage: J2.one().to_vector()
2436 (1, 0, 1)
2437 sage: iota_right(J2.one()).to_vector()
2438 (0, 0, 0, 1, 0, 1)
2439 sage: J.one().to_vector()
2440 (1, 0, 0, 1, 0, 1)
2441
2442 """
2443 (J1,J2) = self.factors()
2444 n = J1.dimension()
2445 V_basis = self.vector_space().basis()
2446 I1 = matrix.column(self.base_ring(), V_basis[:n])
2447 I2 = matrix.column(self.base_ring(), V_basis[n:])
2448 iota_left = lambda x: self.from_vector(I1*x.to_vector())
2449 iota_right = lambda x: self.from_vector(I2*+x.to_vector())
2450 return (iota_left, iota_right)
2451
2452 def inner_product(self, x, y):
2453 r"""
2454 The standard Cartesian inner-product.
2455
2456 We project ``x`` and ``y`` onto our factors, and add up the
2457 inner-products from the subalgebras.
2458
2459 SETUP::
2460
2461
2462 sage: from mjo.eja.eja_algebra import (HadamardEJA,
2463 ....: QuaternionHermitianEJA,
2464 ....: DirectSumEJA)
2465
2466 EXAMPLE::
2467
2468 sage: J1 = HadamardEJA(3)
2469 sage: J2 = QuaternionHermitianEJA(2,QQ,normalize_basis=False)
2470 sage: J = DirectSumEJA(J1,J2)
2471 sage: x1 = J1.one()
2472 sage: x2 = x1
2473 sage: y1 = J2.one()
2474 sage: y2 = y1
2475 sage: x1.inner_product(x2)
2476 3
2477 sage: y1.inner_product(y2)
2478 2
2479 sage: J.one().inner_product(J.one())
2480 5
2481
2482 """
2483 (pi_left, pi_right) = self.projections()
2484 x1 = pi_left(x)
2485 x2 = pi_right(x)
2486 y1 = pi_left(y)
2487 y2 = pi_right(y)
2488
2489 return (x1.inner_product(y1) + x2.inner_product(y2))