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