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