]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
ff2b5d7a9c52ff5c13df7d7aa3682899b9a30559
[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 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=AA)
211 Euclidean Jordan algebra of dimension 2 over Algebraic Real 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 0.7071067811865475?] [0 0]
555 [0 0], [0.7071067811865475? 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().right_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 random_elements(self, count):
777 """
778 Return ``count`` random elements as a tuple.
779
780 SETUP::
781
782 sage: from mjo.eja.eja_algebra import JordanSpinEJA
783
784 EXAMPLES::
785
786 sage: J = JordanSpinEJA(3)
787 sage: x,y,z = J.random_elements(3)
788 sage: all( [ x in J, y in J, z in J ])
789 True
790 sage: len( J.random_elements(10) ) == 10
791 True
792
793 """
794 return tuple( self.random_element() for idx in range(count) )
795
796 @cached_method
797 def rank(self):
798 """
799 Return the rank of this EJA.
800
801 ALGORITHM:
802
803 We first compute the polynomial "column matrices" `p_{k}` that
804 evaluate to `x^k` on the coordinates of `x`. Then, we begin
805 adding them to a matrix one at a time, and trying to solve the
806 system that makes `p_{0}`,`p_{1}`,..., `p_{s-1}` add up to
807 `p_{s}`. This will succeed only when `s` is the rank of the
808 algebra, as proven in a recent draft paper of mine.
809
810 SETUP::
811
812 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
813 ....: RealSymmetricEJA,
814 ....: ComplexHermitianEJA,
815 ....: QuaternionHermitianEJA,
816 ....: random_eja)
817
818 EXAMPLES:
819
820 The rank of the Jordan spin algebra is always two::
821
822 sage: JordanSpinEJA(2).rank()
823 2
824 sage: JordanSpinEJA(3).rank()
825 2
826 sage: JordanSpinEJA(4).rank()
827 2
828
829 The rank of the `n`-by-`n` Hermitian real, complex, or
830 quaternion matrices is `n`::
831
832 sage: RealSymmetricEJA(4).rank()
833 4
834 sage: ComplexHermitianEJA(3).rank()
835 3
836 sage: QuaternionHermitianEJA(2).rank()
837 2
838
839 TESTS:
840
841 Ensure that every EJA that we know how to construct has a
842 positive integer rank, unless the algebra is trivial in
843 which case its rank will be zero::
844
845 sage: set_random_seed()
846 sage: J = random_eja()
847 sage: r = J.rank()
848 sage: r in ZZ
849 True
850 sage: r > 0 or (r == 0 and J.is_trivial())
851 True
852
853 Ensure that computing the rank actually works, since the ranks
854 of all simple algebras are known and will be cached by default::
855
856 sage: J = HadamardEJA(4)
857 sage: J.rank.clear_cache()
858 sage: J.rank()
859 4
860
861 ::
862
863 sage: J = JordanSpinEJA(4)
864 sage: J.rank.clear_cache()
865 sage: J.rank()
866 2
867
868 ::
869
870 sage: J = RealSymmetricEJA(3)
871 sage: J.rank.clear_cache()
872 sage: J.rank()
873 3
874
875 ::
876
877 sage: J = ComplexHermitianEJA(2)
878 sage: J.rank.clear_cache()
879 sage: J.rank()
880 2
881
882 ::
883
884 sage: J = QuaternionHermitianEJA(2)
885 sage: J.rank.clear_cache()
886 sage: J.rank()
887 2
888
889 """
890 n = self.dimension()
891 if n == 0:
892 return 0
893 elif n == 1:
894 return 1
895
896 var_names = [ "X" + str(z) for z in range(1,n+1) ]
897 R = PolynomialRing(self.base_ring(), var_names)
898 vars = R.gens()
899
900 def L_x_i_j(i,j):
901 # From a result in my book, these are the entries of the
902 # basis representation of L_x.
903 return sum( vars[k]*self.monomial(k).operator().matrix()[i,j]
904 for k in range(n) )
905
906 L_x = matrix(R, n, n, L_x_i_j)
907 x_powers = [ vars[k]*(L_x**k)*self.one().to_vector()
908 for k in range(n) ]
909
910 # Can assume n >= 2
911 M = matrix([x_powers[0]])
912 old_rank = 1
913
914 for d in range(1,n):
915 M = matrix(M.rows() + [x_powers[d]])
916 M.echelonize()
917 # TODO: we've basically solved the system here.
918 # We should save the echelonized matrix somehow
919 # so that it can be reused in the charpoly method.
920 new_rank = M.rank()
921 if new_rank == old_rank:
922 return new_rank
923 else:
924 old_rank = new_rank
925
926 return n
927
928
929 def vector_space(self):
930 """
931 Return the vector space that underlies this algebra.
932
933 SETUP::
934
935 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
936
937 EXAMPLES::
938
939 sage: J = RealSymmetricEJA(2)
940 sage: J.vector_space()
941 Vector space of dimension 3 over...
942
943 """
944 return self.zero().to_vector().parent().ambient_vector_space()
945
946
947 Element = FiniteDimensionalEuclideanJordanAlgebraElement
948
949
950 class KnownRankEJA(object):
951 """
952 A class for algebras that we actually know we can construct. The
953 main issue is that, for most of our methods to make sense, we need
954 to know the rank of our algebra. Thus we can't simply generate a
955 "random" algebra, or even check that a given basis and product
956 satisfy the axioms; because even if everything looks OK, we wouldn't
957 know the rank we need to actuallty build the thing.
958
959 Not really a subclass of FDEJA because doing that causes method
960 resolution errors, e.g.
961
962 TypeError: Error when calling the metaclass bases
963 Cannot create a consistent method resolution
964 order (MRO) for bases FiniteDimensionalEuclideanJordanAlgebra,
965 KnownRankEJA
966
967 """
968 @staticmethod
969 def _max_test_case_size():
970 """
971 Return an integer "size" that is an upper bound on the size of
972 this algebra when it is used in a random test
973 case. Unfortunately, the term "size" is quite vague -- when
974 dealing with `R^n` under either the Hadamard or Jordan spin
975 product, the "size" refers to the dimension `n`. When dealing
976 with a matrix algebra (real symmetric or complex/quaternion
977 Hermitian), it refers to the size of the matrix, which is
978 far less than the dimension of the underlying vector space.
979
980 We default to five in this class, which is safe in `R^n`. The
981 matrix algebra subclasses (or any class where the "size" is
982 interpreted to be far less than the dimension) should override
983 with a smaller number.
984 """
985 return 5
986
987 @classmethod
988 def random_instance(cls, field=AA, **kwargs):
989 """
990 Return a random instance of this type of algebra.
991
992 Beware, this will crash for "most instances" because the
993 constructor below looks wrong.
994 """
995 if cls is TrivialEJA:
996 # The TrivialEJA class doesn't take an "n" argument because
997 # there's only one.
998 return cls(field)
999
1000 n = ZZ.random_element(cls._max_test_case_size()) + 1
1001 return cls(n, field, **kwargs)
1002
1003
1004 class HadamardEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
1005 """
1006 Return the Euclidean Jordan Algebra corresponding to the set
1007 `R^n` under the Hadamard product.
1008
1009 Note: this is nothing more than the Cartesian product of ``n``
1010 copies of the spin algebra. Once Cartesian product algebras
1011 are implemented, this can go.
1012
1013 SETUP::
1014
1015 sage: from mjo.eja.eja_algebra import HadamardEJA
1016
1017 EXAMPLES:
1018
1019 This multiplication table can be verified by hand::
1020
1021 sage: J = HadamardEJA(3)
1022 sage: e0,e1,e2 = J.gens()
1023 sage: e0*e0
1024 e0
1025 sage: e0*e1
1026 0
1027 sage: e0*e2
1028 0
1029 sage: e1*e1
1030 e1
1031 sage: e1*e2
1032 0
1033 sage: e2*e2
1034 e2
1035
1036 TESTS:
1037
1038 We can change the generator prefix::
1039
1040 sage: HadamardEJA(3, prefix='r').gens()
1041 (r0, r1, r2)
1042
1043 """
1044 def __init__(self, n, field=AA, **kwargs):
1045 V = VectorSpace(field, n)
1046 mult_table = [ [ V.gen(i)*(i == j) for j in range(n) ]
1047 for i in range(n) ]
1048
1049 fdeja = super(HadamardEJA, self)
1050 return fdeja.__init__(field, mult_table, rank=n, **kwargs)
1051
1052 def inner_product(self, x, y):
1053 """
1054 Faster to reimplement than to use natural representations.
1055
1056 SETUP::
1057
1058 sage: from mjo.eja.eja_algebra import HadamardEJA
1059
1060 TESTS:
1061
1062 Ensure that this is the usual inner product for the algebras
1063 over `R^n`::
1064
1065 sage: set_random_seed()
1066 sage: J = HadamardEJA.random_instance()
1067 sage: x,y = J.random_elements(2)
1068 sage: X = x.natural_representation()
1069 sage: Y = y.natural_representation()
1070 sage: x.inner_product(y) == J.natural_inner_product(X,Y)
1071 True
1072
1073 """
1074 return x.to_vector().inner_product(y.to_vector())
1075
1076
1077 def random_eja(field=AA, nontrivial=False):
1078 """
1079 Return a "random" finite-dimensional Euclidean Jordan Algebra.
1080
1081 SETUP::
1082
1083 sage: from mjo.eja.eja_algebra import random_eja
1084
1085 TESTS::
1086
1087 sage: random_eja()
1088 Euclidean Jordan algebra of dimension...
1089
1090 """
1091 eja_classes = KnownRankEJA.__subclasses__()
1092 if nontrivial:
1093 eja_classes.remove(TrivialEJA)
1094 classname = choice(eja_classes)
1095 return classname.random_instance(field=field)
1096
1097
1098
1099
1100
1101
1102 class MatrixEuclideanJordanAlgebra(FiniteDimensionalEuclideanJordanAlgebra):
1103 @staticmethod
1104 def _max_test_case_size():
1105 # Play it safe, since this will be squared and the underlying
1106 # field can have dimension 4 (quaternions) too.
1107 return 2
1108
1109 def __init__(self, field, basis, rank, normalize_basis=True, **kwargs):
1110 """
1111 Compared to the superclass constructor, we take a basis instead of
1112 a multiplication table because the latter can be computed in terms
1113 of the former when the product is known (like it is here).
1114 """
1115 # Used in this class's fast _charpoly_coeff() override.
1116 self._basis_normalizers = None
1117
1118 # We're going to loop through this a few times, so now's a good
1119 # time to ensure that it isn't a generator expression.
1120 basis = tuple(basis)
1121
1122 if rank > 1 and normalize_basis:
1123 # We'll need sqrt(2) to normalize the basis, and this
1124 # winds up in the multiplication table, so the whole
1125 # algebra needs to be over the field extension.
1126 R = PolynomialRing(field, 'z')
1127 z = R.gen()
1128 p = z**2 - 2
1129 if p.is_irreducible():
1130 field = field.extension(p, 'sqrt2', embedding=RLF(2).sqrt())
1131 basis = tuple( s.change_ring(field) for s in basis )
1132 self._basis_normalizers = tuple(
1133 ~(self.natural_inner_product(s,s).sqrt()) for s in basis )
1134 basis = tuple(s*c for (s,c) in zip(basis,self._basis_normalizers))
1135
1136 Qs = self.multiplication_table_from_matrix_basis(basis)
1137
1138 fdeja = super(MatrixEuclideanJordanAlgebra, self)
1139 return fdeja.__init__(field,
1140 Qs,
1141 rank=rank,
1142 natural_basis=basis,
1143 **kwargs)
1144
1145
1146 def _rank_computation(self):
1147 r"""
1148 Override the parent method with something that tries to compute
1149 over a faster (non-extension) field.
1150 """
1151 if self._basis_normalizers is None:
1152 # We didn't normalize, so assume that the basis we started
1153 # with had entries in a nice field.
1154 return super(MatrixEuclideanJordanAlgebra, self)._rank_computation()
1155 else:
1156 basis = ( (b/n) for (b,n) in zip(self.natural_basis(),
1157 self._basis_normalizers) )
1158
1159 # Do this over the rationals and convert back at the end.
1160 # Only works because we know the entries of the basis are
1161 # integers.
1162 J = MatrixEuclideanJordanAlgebra(QQ,
1163 basis,
1164 self.rank(),
1165 normalize_basis=False)
1166 return J._rank_computation()
1167
1168 @cached_method
1169 def _charpoly_coeff(self, i):
1170 """
1171 Override the parent method with something that tries to compute
1172 over a faster (non-extension) field.
1173 """
1174 if self._basis_normalizers is None:
1175 # We didn't normalize, so assume that the basis we started
1176 # with had entries in a nice field.
1177 return super(MatrixEuclideanJordanAlgebra, self)._charpoly_coeff(i)
1178 else:
1179 basis = ( (b/n) for (b,n) in zip(self.natural_basis(),
1180 self._basis_normalizers) )
1181
1182 # Do this over the rationals and convert back at the end.
1183 J = MatrixEuclideanJordanAlgebra(QQ,
1184 basis,
1185 self.rank(),
1186 normalize_basis=False)
1187 (_,x,_,_) = J._charpoly_matrix_system()
1188 p = J._charpoly_coeff(i)
1189 # p might be missing some vars, have to substitute "optionally"
1190 pairs = zip(x.base_ring().gens(), self._basis_normalizers)
1191 substitutions = { v: v*c for (v,c) in pairs }
1192 result = p.subs(substitutions)
1193
1194 # The result of "subs" can be either a coefficient-ring
1195 # element or a polynomial. Gotta handle both cases.
1196 if result in QQ:
1197 return self.base_ring()(result)
1198 else:
1199 return result.change_ring(self.base_ring())
1200
1201
1202 @staticmethod
1203 def multiplication_table_from_matrix_basis(basis):
1204 """
1205 At least three of the five simple Euclidean Jordan algebras have the
1206 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
1207 multiplication on the right is matrix multiplication. Given a basis
1208 for the underlying matrix space, this function returns a
1209 multiplication table (obtained by looping through the basis
1210 elements) for an algebra of those matrices.
1211 """
1212 # In S^2, for example, we nominally have four coordinates even
1213 # though the space is of dimension three only. The vector space V
1214 # is supposed to hold the entire long vector, and the subspace W
1215 # of V will be spanned by the vectors that arise from symmetric
1216 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
1217 field = basis[0].base_ring()
1218 dimension = basis[0].nrows()
1219
1220 V = VectorSpace(field, dimension**2)
1221 W = V.span_of_basis( _mat2vec(s) for s in basis )
1222 n = len(basis)
1223 mult_table = [[W.zero() for j in range(n)] for i in range(n)]
1224 for i in range(n):
1225 for j in range(n):
1226 mat_entry = (basis[i]*basis[j] + basis[j]*basis[i])/2
1227 mult_table[i][j] = W.coordinate_vector(_mat2vec(mat_entry))
1228
1229 return mult_table
1230
1231
1232 @staticmethod
1233 def real_embed(M):
1234 """
1235 Embed the matrix ``M`` into a space of real matrices.
1236
1237 The matrix ``M`` can have entries in any field at the moment:
1238 the real numbers, complex numbers, or quaternions. And although
1239 they are not a field, we can probably support octonions at some
1240 point, too. This function returns a real matrix that "acts like"
1241 the original with respect to matrix multiplication; i.e.
1242
1243 real_embed(M*N) = real_embed(M)*real_embed(N)
1244
1245 """
1246 raise NotImplementedError
1247
1248
1249 @staticmethod
1250 def real_unembed(M):
1251 """
1252 The inverse of :meth:`real_embed`.
1253 """
1254 raise NotImplementedError
1255
1256
1257 @classmethod
1258 def natural_inner_product(cls,X,Y):
1259 Xu = cls.real_unembed(X)
1260 Yu = cls.real_unembed(Y)
1261 tr = (Xu*Yu).trace()
1262
1263 if tr in RLF:
1264 # It's real already.
1265 return tr
1266
1267 # Otherwise, try the thing that works for complex numbers; and
1268 # if that doesn't work, the thing that works for quaternions.
1269 try:
1270 return tr.vector()[0] # real part, imag part is index 1
1271 except AttributeError:
1272 # A quaternions doesn't have a vector() method, but does
1273 # have coefficient_tuple() method that returns the
1274 # coefficients of 1, i, j, and k -- in that order.
1275 return tr.coefficient_tuple()[0]
1276
1277
1278 class RealMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1279 @staticmethod
1280 def real_embed(M):
1281 """
1282 The identity function, for embedding real matrices into real
1283 matrices.
1284 """
1285 return M
1286
1287 @staticmethod
1288 def real_unembed(M):
1289 """
1290 The identity function, for unembedding real matrices from real
1291 matrices.
1292 """
1293 return M
1294
1295
1296 class RealSymmetricEJA(RealMatrixEuclideanJordanAlgebra, KnownRankEJA):
1297 """
1298 The rank-n simple EJA consisting of real symmetric n-by-n
1299 matrices, the usual symmetric Jordan product, and the trace inner
1300 product. It has dimension `(n^2 + n)/2` over the reals.
1301
1302 SETUP::
1303
1304 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1305
1306 EXAMPLES::
1307
1308 sage: J = RealSymmetricEJA(2)
1309 sage: e0, e1, e2 = J.gens()
1310 sage: e0*e0
1311 e0
1312 sage: e1*e1
1313 1/2*e0 + 1/2*e2
1314 sage: e2*e2
1315 e2
1316
1317 In theory, our "field" can be any subfield of the reals::
1318
1319 sage: RealSymmetricEJA(2, RDF)
1320 Euclidean Jordan algebra of dimension 3 over Real Double Field
1321 sage: RealSymmetricEJA(2, RR)
1322 Euclidean Jordan algebra of dimension 3 over Real Field with
1323 53 bits of precision
1324
1325 TESTS:
1326
1327 The dimension of this algebra is `(n^2 + n) / 2`::
1328
1329 sage: set_random_seed()
1330 sage: n_max = RealSymmetricEJA._max_test_case_size()
1331 sage: n = ZZ.random_element(1, n_max)
1332 sage: J = RealSymmetricEJA(n)
1333 sage: J.dimension() == (n^2 + n)/2
1334 True
1335
1336 The Jordan multiplication is what we think it is::
1337
1338 sage: set_random_seed()
1339 sage: J = RealSymmetricEJA.random_instance()
1340 sage: x,y = J.random_elements(2)
1341 sage: actual = (x*y).natural_representation()
1342 sage: X = x.natural_representation()
1343 sage: Y = y.natural_representation()
1344 sage: expected = (X*Y + Y*X)/2
1345 sage: actual == expected
1346 True
1347 sage: J(expected) == x*y
1348 True
1349
1350 We can change the generator prefix::
1351
1352 sage: RealSymmetricEJA(3, prefix='q').gens()
1353 (q0, q1, q2, q3, q4, q5)
1354
1355 Our natural basis is normalized with respect to the natural inner
1356 product unless we specify otherwise::
1357
1358 sage: set_random_seed()
1359 sage: J = RealSymmetricEJA.random_instance()
1360 sage: all( b.norm() == 1 for b in J.gens() )
1361 True
1362
1363 Since our natural basis is normalized with respect to the natural
1364 inner product, and since we know that this algebra is an EJA, any
1365 left-multiplication operator's matrix will be symmetric because
1366 natural->EJA basis representation is an isometry and within the EJA
1367 the operator is self-adjoint by the Jordan axiom::
1368
1369 sage: set_random_seed()
1370 sage: x = RealSymmetricEJA.random_instance().random_element()
1371 sage: x.operator().matrix().is_symmetric()
1372 True
1373
1374 """
1375 @classmethod
1376 def _denormalized_basis(cls, n, field):
1377 """
1378 Return a basis for the space of real symmetric n-by-n matrices.
1379
1380 SETUP::
1381
1382 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1383
1384 TESTS::
1385
1386 sage: set_random_seed()
1387 sage: n = ZZ.random_element(1,5)
1388 sage: B = RealSymmetricEJA._denormalized_basis(n,QQ)
1389 sage: all( M.is_symmetric() for M in B)
1390 True
1391
1392 """
1393 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1394 # coordinates.
1395 S = []
1396 for i in range(n):
1397 for j in range(i+1):
1398 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1399 if i == j:
1400 Sij = Eij
1401 else:
1402 Sij = Eij + Eij.transpose()
1403 S.append(Sij)
1404 return S
1405
1406
1407 @staticmethod
1408 def _max_test_case_size():
1409 return 4 # Dimension 10
1410
1411
1412 def __init__(self, n, field=AA, **kwargs):
1413 basis = self._denormalized_basis(n, field)
1414 super(RealSymmetricEJA, self).__init__(field, basis, n, **kwargs)
1415
1416
1417 class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1418 @staticmethod
1419 def real_embed(M):
1420 """
1421 Embed the n-by-n complex matrix ``M`` into the space of real
1422 matrices of size 2n-by-2n via the map the sends each entry `z = a +
1423 bi` to the block matrix ``[[a,b],[-b,a]]``.
1424
1425 SETUP::
1426
1427 sage: from mjo.eja.eja_algebra import \
1428 ....: ComplexMatrixEuclideanJordanAlgebra
1429
1430 EXAMPLES::
1431
1432 sage: F = QuadraticField(-1, 'I')
1433 sage: x1 = F(4 - 2*i)
1434 sage: x2 = F(1 + 2*i)
1435 sage: x3 = F(-i)
1436 sage: x4 = F(6)
1437 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
1438 sage: ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1439 [ 4 -2| 1 2]
1440 [ 2 4|-2 1]
1441 [-----+-----]
1442 [ 0 -1| 6 0]
1443 [ 1 0| 0 6]
1444
1445 TESTS:
1446
1447 Embedding is a homomorphism (isomorphism, in fact)::
1448
1449 sage: set_random_seed()
1450 sage: n_max = ComplexMatrixEuclideanJordanAlgebra._max_test_case_size()
1451 sage: n = ZZ.random_element(n_max)
1452 sage: F = QuadraticField(-1, 'I')
1453 sage: X = random_matrix(F, n)
1454 sage: Y = random_matrix(F, n)
1455 sage: Xe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X)
1456 sage: Ye = ComplexMatrixEuclideanJordanAlgebra.real_embed(Y)
1457 sage: XYe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1458 sage: Xe*Ye == XYe
1459 True
1460
1461 """
1462 n = M.nrows()
1463 if M.ncols() != n:
1464 raise ValueError("the matrix 'M' must be square")
1465
1466 # We don't need any adjoined elements...
1467 field = M.base_ring().base_ring()
1468
1469 blocks = []
1470 for z in M.list():
1471 a = z.list()[0] # real part, I guess
1472 b = z.list()[1] # imag part, I guess
1473 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1474
1475 return matrix.block(field, n, blocks)
1476
1477
1478 @staticmethod
1479 def real_unembed(M):
1480 """
1481 The inverse of _embed_complex_matrix().
1482
1483 SETUP::
1484
1485 sage: from mjo.eja.eja_algebra import \
1486 ....: ComplexMatrixEuclideanJordanAlgebra
1487
1488 EXAMPLES::
1489
1490 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1491 ....: [-2, 1, -4, 3],
1492 ....: [ 9, 10, 11, 12],
1493 ....: [-10, 9, -12, 11] ])
1494 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(A)
1495 [ 2*I + 1 4*I + 3]
1496 [ 10*I + 9 12*I + 11]
1497
1498 TESTS:
1499
1500 Unembedding is the inverse of embedding::
1501
1502 sage: set_random_seed()
1503 sage: F = QuadraticField(-1, 'I')
1504 sage: M = random_matrix(F, 3)
1505 sage: Me = ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1506 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1507 True
1508
1509 """
1510 n = ZZ(M.nrows())
1511 if M.ncols() != n:
1512 raise ValueError("the matrix 'M' must be square")
1513 if not n.mod(2).is_zero():
1514 raise ValueError("the matrix 'M' must be a complex embedding")
1515
1516 # If "M" was normalized, its base ring might have roots
1517 # adjoined and they can stick around after unembedding.
1518 field = M.base_ring()
1519 R = PolynomialRing(field, 'z')
1520 z = R.gen()
1521 if field is AA:
1522 # Sage doesn't know how to embed AA into QQbar, i.e. how
1523 # to adjoin sqrt(-1) to AA.
1524 F = QQbar
1525 else:
1526 F = field.extension(z**2 + 1, 'I', embedding=CLF(-1).sqrt())
1527 i = F.gen()
1528
1529 # Go top-left to bottom-right (reading order), converting every
1530 # 2-by-2 block we see to a single complex element.
1531 elements = []
1532 for k in range(n/2):
1533 for j in range(n/2):
1534 submat = M[2*k:2*k+2,2*j:2*j+2]
1535 if submat[0,0] != submat[1,1]:
1536 raise ValueError('bad on-diagonal submatrix')
1537 if submat[0,1] != -submat[1,0]:
1538 raise ValueError('bad off-diagonal submatrix')
1539 z = submat[0,0] + submat[0,1]*i
1540 elements.append(z)
1541
1542 return matrix(F, n/2, elements)
1543
1544
1545 @classmethod
1546 def natural_inner_product(cls,X,Y):
1547 """
1548 Compute a natural inner product in this algebra directly from
1549 its real embedding.
1550
1551 SETUP::
1552
1553 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1554
1555 TESTS:
1556
1557 This gives the same answer as the slow, default method implemented
1558 in :class:`MatrixEuclideanJordanAlgebra`::
1559
1560 sage: set_random_seed()
1561 sage: J = ComplexHermitianEJA.random_instance()
1562 sage: x,y = J.random_elements(2)
1563 sage: Xe = x.natural_representation()
1564 sage: Ye = y.natural_representation()
1565 sage: X = ComplexHermitianEJA.real_unembed(Xe)
1566 sage: Y = ComplexHermitianEJA.real_unembed(Ye)
1567 sage: expected = (X*Y).trace().real()
1568 sage: actual = ComplexHermitianEJA.natural_inner_product(Xe,Ye)
1569 sage: actual == expected
1570 True
1571
1572 """
1573 return RealMatrixEuclideanJordanAlgebra.natural_inner_product(X,Y)/2
1574
1575
1576 class ComplexHermitianEJA(ComplexMatrixEuclideanJordanAlgebra, KnownRankEJA):
1577 """
1578 The rank-n simple EJA consisting of complex Hermitian n-by-n
1579 matrices over the real numbers, the usual symmetric Jordan product,
1580 and the real-part-of-trace inner product. It has dimension `n^2` over
1581 the reals.
1582
1583 SETUP::
1584
1585 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1586
1587 EXAMPLES:
1588
1589 In theory, our "field" can be any subfield of the reals::
1590
1591 sage: ComplexHermitianEJA(2, RDF)
1592 Euclidean Jordan algebra of dimension 4 over Real Double Field
1593 sage: ComplexHermitianEJA(2, RR)
1594 Euclidean Jordan algebra of dimension 4 over Real Field with
1595 53 bits of precision
1596
1597 TESTS:
1598
1599 The dimension of this algebra is `n^2`::
1600
1601 sage: set_random_seed()
1602 sage: n_max = ComplexHermitianEJA._max_test_case_size()
1603 sage: n = ZZ.random_element(1, n_max)
1604 sage: J = ComplexHermitianEJA(n)
1605 sage: J.dimension() == n^2
1606 True
1607
1608 The Jordan multiplication is what we think it is::
1609
1610 sage: set_random_seed()
1611 sage: J = ComplexHermitianEJA.random_instance()
1612 sage: x,y = J.random_elements(2)
1613 sage: actual = (x*y).natural_representation()
1614 sage: X = x.natural_representation()
1615 sage: Y = y.natural_representation()
1616 sage: expected = (X*Y + Y*X)/2
1617 sage: actual == expected
1618 True
1619 sage: J(expected) == x*y
1620 True
1621
1622 We can change the generator prefix::
1623
1624 sage: ComplexHermitianEJA(2, prefix='z').gens()
1625 (z0, z1, z2, z3)
1626
1627 Our natural basis is normalized with respect to the natural inner
1628 product unless we specify otherwise::
1629
1630 sage: set_random_seed()
1631 sage: J = ComplexHermitianEJA.random_instance()
1632 sage: all( b.norm() == 1 for b in J.gens() )
1633 True
1634
1635 Since our natural basis is normalized with respect to the natural
1636 inner product, and since we know that this algebra is an EJA, any
1637 left-multiplication operator's matrix will be symmetric because
1638 natural->EJA basis representation is an isometry and within the EJA
1639 the operator is self-adjoint by the Jordan axiom::
1640
1641 sage: set_random_seed()
1642 sage: x = ComplexHermitianEJA.random_instance().random_element()
1643 sage: x.operator().matrix().is_symmetric()
1644 True
1645
1646 """
1647
1648 @classmethod
1649 def _denormalized_basis(cls, n, field):
1650 """
1651 Returns a basis for the space of complex Hermitian n-by-n matrices.
1652
1653 Why do we embed these? Basically, because all of numerical linear
1654 algebra assumes that you're working with vectors consisting of `n`
1655 entries from a field and scalars from the same field. There's no way
1656 to tell SageMath that (for example) the vectors contain complex
1657 numbers, while the scalar field is real.
1658
1659 SETUP::
1660
1661 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1662
1663 TESTS::
1664
1665 sage: set_random_seed()
1666 sage: n = ZZ.random_element(1,5)
1667 sage: field = QuadraticField(2, 'sqrt2')
1668 sage: B = ComplexHermitianEJA._denormalized_basis(n, field)
1669 sage: all( M.is_symmetric() for M in B)
1670 True
1671
1672 """
1673 R = PolynomialRing(field, 'z')
1674 z = R.gen()
1675 F = field.extension(z**2 + 1, 'I')
1676 I = F.gen()
1677
1678 # This is like the symmetric case, but we need to be careful:
1679 #
1680 # * We want conjugate-symmetry, not just symmetry.
1681 # * The diagonal will (as a result) be real.
1682 #
1683 S = []
1684 for i in range(n):
1685 for j in range(i+1):
1686 Eij = matrix(F, n, lambda k,l: k==i and l==j)
1687 if i == j:
1688 Sij = cls.real_embed(Eij)
1689 S.append(Sij)
1690 else:
1691 # The second one has a minus because it's conjugated.
1692 Sij_real = cls.real_embed(Eij + Eij.transpose())
1693 S.append(Sij_real)
1694 Sij_imag = cls.real_embed(I*Eij - I*Eij.transpose())
1695 S.append(Sij_imag)
1696
1697 # Since we embedded these, we can drop back to the "field" that we
1698 # started with instead of the complex extension "F".
1699 return ( s.change_ring(field) for s in S )
1700
1701
1702 def __init__(self, n, field=AA, **kwargs):
1703 basis = self._denormalized_basis(n,field)
1704 super(ComplexHermitianEJA,self).__init__(field, basis, n, **kwargs)
1705
1706
1707 class QuaternionMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1708 @staticmethod
1709 def real_embed(M):
1710 """
1711 Embed the n-by-n quaternion matrix ``M`` into the space of real
1712 matrices of size 4n-by-4n by first sending each quaternion entry `z
1713 = a + bi + cj + dk` to the block-complex matrix ``[[a + bi,
1714 c+di],[-c + di, a-bi]]`, and then embedding those into a real
1715 matrix.
1716
1717 SETUP::
1718
1719 sage: from mjo.eja.eja_algebra import \
1720 ....: QuaternionMatrixEuclideanJordanAlgebra
1721
1722 EXAMPLES::
1723
1724 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1725 sage: i,j,k = Q.gens()
1726 sage: x = 1 + 2*i + 3*j + 4*k
1727 sage: M = matrix(Q, 1, [[x]])
1728 sage: QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1729 [ 1 2 3 4]
1730 [-2 1 -4 3]
1731 [-3 4 1 -2]
1732 [-4 -3 2 1]
1733
1734 Embedding is a homomorphism (isomorphism, in fact)::
1735
1736 sage: set_random_seed()
1737 sage: n_max = QuaternionMatrixEuclideanJordanAlgebra._max_test_case_size()
1738 sage: n = ZZ.random_element(n_max)
1739 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1740 sage: X = random_matrix(Q, n)
1741 sage: Y = random_matrix(Q, n)
1742 sage: Xe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X)
1743 sage: Ye = QuaternionMatrixEuclideanJordanAlgebra.real_embed(Y)
1744 sage: XYe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1745 sage: Xe*Ye == XYe
1746 True
1747
1748 """
1749 quaternions = M.base_ring()
1750 n = M.nrows()
1751 if M.ncols() != n:
1752 raise ValueError("the matrix 'M' must be square")
1753
1754 F = QuadraticField(-1, 'I')
1755 i = F.gen()
1756
1757 blocks = []
1758 for z in M.list():
1759 t = z.coefficient_tuple()
1760 a = t[0]
1761 b = t[1]
1762 c = t[2]
1763 d = t[3]
1764 cplxM = matrix(F, 2, [[ a + b*i, c + d*i],
1765 [-c + d*i, a - b*i]])
1766 realM = ComplexMatrixEuclideanJordanAlgebra.real_embed(cplxM)
1767 blocks.append(realM)
1768
1769 # We should have real entries by now, so use the realest field
1770 # we've got for the return value.
1771 return matrix.block(quaternions.base_ring(), n, blocks)
1772
1773
1774
1775 @staticmethod
1776 def real_unembed(M):
1777 """
1778 The inverse of _embed_quaternion_matrix().
1779
1780 SETUP::
1781
1782 sage: from mjo.eja.eja_algebra import \
1783 ....: QuaternionMatrixEuclideanJordanAlgebra
1784
1785 EXAMPLES::
1786
1787 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
1788 ....: [-2, 1, -4, 3],
1789 ....: [-3, 4, 1, -2],
1790 ....: [-4, -3, 2, 1]])
1791 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(M)
1792 [1 + 2*i + 3*j + 4*k]
1793
1794 TESTS:
1795
1796 Unembedding is the inverse of embedding::
1797
1798 sage: set_random_seed()
1799 sage: Q = QuaternionAlgebra(QQ, -1, -1)
1800 sage: M = random_matrix(Q, 3)
1801 sage: Me = QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1802 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1803 True
1804
1805 """
1806 n = ZZ(M.nrows())
1807 if M.ncols() != n:
1808 raise ValueError("the matrix 'M' must be square")
1809 if not n.mod(4).is_zero():
1810 raise ValueError("the matrix 'M' must be a quaternion embedding")
1811
1812 # Use the base ring of the matrix to ensure that its entries can be
1813 # multiplied by elements of the quaternion algebra.
1814 field = M.base_ring()
1815 Q = QuaternionAlgebra(field,-1,-1)
1816 i,j,k = Q.gens()
1817
1818 # Go top-left to bottom-right (reading order), converting every
1819 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
1820 # quaternion block.
1821 elements = []
1822 for l in range(n/4):
1823 for m in range(n/4):
1824 submat = ComplexMatrixEuclideanJordanAlgebra.real_unembed(
1825 M[4*l:4*l+4,4*m:4*m+4] )
1826 if submat[0,0] != submat[1,1].conjugate():
1827 raise ValueError('bad on-diagonal submatrix')
1828 if submat[0,1] != -submat[1,0].conjugate():
1829 raise ValueError('bad off-diagonal submatrix')
1830 z = submat[0,0].real()
1831 z += submat[0,0].imag()*i
1832 z += submat[0,1].real()*j
1833 z += submat[0,1].imag()*k
1834 elements.append(z)
1835
1836 return matrix(Q, n/4, elements)
1837
1838
1839 @classmethod
1840 def natural_inner_product(cls,X,Y):
1841 """
1842 Compute a natural inner product in this algebra directly from
1843 its real embedding.
1844
1845 SETUP::
1846
1847 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1848
1849 TESTS:
1850
1851 This gives the same answer as the slow, default method implemented
1852 in :class:`MatrixEuclideanJordanAlgebra`::
1853
1854 sage: set_random_seed()
1855 sage: J = QuaternionHermitianEJA.random_instance()
1856 sage: x,y = J.random_elements(2)
1857 sage: Xe = x.natural_representation()
1858 sage: Ye = y.natural_representation()
1859 sage: X = QuaternionHermitianEJA.real_unembed(Xe)
1860 sage: Y = QuaternionHermitianEJA.real_unembed(Ye)
1861 sage: expected = (X*Y).trace().coefficient_tuple()[0]
1862 sage: actual = QuaternionHermitianEJA.natural_inner_product(Xe,Ye)
1863 sage: actual == expected
1864 True
1865
1866 """
1867 return RealMatrixEuclideanJordanAlgebra.natural_inner_product(X,Y)/4
1868
1869
1870 class QuaternionHermitianEJA(QuaternionMatrixEuclideanJordanAlgebra,
1871 KnownRankEJA):
1872 """
1873 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
1874 matrices, the usual symmetric Jordan product, and the
1875 real-part-of-trace inner product. It has dimension `2n^2 - n` over
1876 the reals.
1877
1878 SETUP::
1879
1880 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1881
1882 EXAMPLES:
1883
1884 In theory, our "field" can be any subfield of the reals::
1885
1886 sage: QuaternionHermitianEJA(2, RDF)
1887 Euclidean Jordan algebra of dimension 6 over Real Double Field
1888 sage: QuaternionHermitianEJA(2, RR)
1889 Euclidean Jordan algebra of dimension 6 over Real Field with
1890 53 bits of precision
1891
1892 TESTS:
1893
1894 The dimension of this algebra is `2*n^2 - n`::
1895
1896 sage: set_random_seed()
1897 sage: n_max = QuaternionHermitianEJA._max_test_case_size()
1898 sage: n = ZZ.random_element(1, n_max)
1899 sage: J = QuaternionHermitianEJA(n)
1900 sage: J.dimension() == 2*(n^2) - n
1901 True
1902
1903 The Jordan multiplication is what we think it is::
1904
1905 sage: set_random_seed()
1906 sage: J = QuaternionHermitianEJA.random_instance()
1907 sage: x,y = J.random_elements(2)
1908 sage: actual = (x*y).natural_representation()
1909 sage: X = x.natural_representation()
1910 sage: Y = y.natural_representation()
1911 sage: expected = (X*Y + Y*X)/2
1912 sage: actual == expected
1913 True
1914 sage: J(expected) == x*y
1915 True
1916
1917 We can change the generator prefix::
1918
1919 sage: QuaternionHermitianEJA(2, prefix='a').gens()
1920 (a0, a1, a2, a3, a4, a5)
1921
1922 Our natural basis is normalized with respect to the natural inner
1923 product unless we specify otherwise::
1924
1925 sage: set_random_seed()
1926 sage: J = QuaternionHermitianEJA.random_instance()
1927 sage: all( b.norm() == 1 for b in J.gens() )
1928 True
1929
1930 Since our natural basis is normalized with respect to the natural
1931 inner product, and since we know that this algebra is an EJA, any
1932 left-multiplication operator's matrix will be symmetric because
1933 natural->EJA basis representation is an isometry and within the EJA
1934 the operator is self-adjoint by the Jordan axiom::
1935
1936 sage: set_random_seed()
1937 sage: x = QuaternionHermitianEJA.random_instance().random_element()
1938 sage: x.operator().matrix().is_symmetric()
1939 True
1940
1941 """
1942 @classmethod
1943 def _denormalized_basis(cls, n, field):
1944 """
1945 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
1946
1947 Why do we embed these? Basically, because all of numerical
1948 linear algebra assumes that you're working with vectors consisting
1949 of `n` entries from a field and scalars from the same field. There's
1950 no way to tell SageMath that (for example) the vectors contain
1951 complex numbers, while the scalar field is real.
1952
1953 SETUP::
1954
1955 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1956
1957 TESTS::
1958
1959 sage: set_random_seed()
1960 sage: n = ZZ.random_element(1,5)
1961 sage: B = QuaternionHermitianEJA._denormalized_basis(n,QQ)
1962 sage: all( M.is_symmetric() for M in B )
1963 True
1964
1965 """
1966 Q = QuaternionAlgebra(QQ,-1,-1)
1967 I,J,K = Q.gens()
1968
1969 # This is like the symmetric case, but we need to be careful:
1970 #
1971 # * We want conjugate-symmetry, not just symmetry.
1972 # * The diagonal will (as a result) be real.
1973 #
1974 S = []
1975 for i in range(n):
1976 for j in range(i+1):
1977 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
1978 if i == j:
1979 Sij = cls.real_embed(Eij)
1980 S.append(Sij)
1981 else:
1982 # The second, third, and fourth ones have a minus
1983 # because they're conjugated.
1984 Sij_real = cls.real_embed(Eij + Eij.transpose())
1985 S.append(Sij_real)
1986 Sij_I = cls.real_embed(I*Eij - I*Eij.transpose())
1987 S.append(Sij_I)
1988 Sij_J = cls.real_embed(J*Eij - J*Eij.transpose())
1989 S.append(Sij_J)
1990 Sij_K = cls.real_embed(K*Eij - K*Eij.transpose())
1991 S.append(Sij_K)
1992
1993 # Since we embedded these, we can drop back to the "field" that we
1994 # started with instead of the quaternion algebra "Q".
1995 return ( s.change_ring(field) for s in S )
1996
1997
1998 def __init__(self, n, field=AA, **kwargs):
1999 basis = self._denormalized_basis(n,field)
2000 super(QuaternionHermitianEJA,self).__init__(field, basis, n, **kwargs)
2001
2002
2003 class BilinearFormEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
2004 r"""
2005 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2006 with the half-trace inner product and jordan product ``x*y =
2007 (x0*y0 + <B*x_bar,y_bar>, x0*y_bar + y0*x_bar)`` where ``B`` is a
2008 symmetric positive-definite "bilinear form" matrix. It has
2009 dimension `n` over the reals, and reduces to the ``JordanSpinEJA``
2010 when ``B`` is the identity matrix of order ``n-1``.
2011
2012 SETUP::
2013
2014 sage: from mjo.eja.eja_algebra import (BilinearFormEJA,
2015 ....: JordanSpinEJA)
2016
2017 EXAMPLES:
2018
2019 When no bilinear form is specified, the identity matrix is used,
2020 and the resulting algebra is the Jordan spin algebra::
2021
2022 sage: J0 = BilinearFormEJA(3)
2023 sage: J1 = JordanSpinEJA(3)
2024 sage: J0.multiplication_table() == J0.multiplication_table()
2025 True
2026
2027 TESTS:
2028
2029 We can create a zero-dimensional algebra::
2030
2031 sage: J = BilinearFormEJA(0)
2032 sage: J.basis()
2033 Finite family {}
2034
2035 We can check the multiplication condition given in the Jordan, von
2036 Neumann, and Wigner paper (and also discussed on my "On the
2037 symmetry..." paper). Note that this relies heavily on the standard
2038 choice of basis, as does anything utilizing the bilinear form matrix::
2039
2040 sage: set_random_seed()
2041 sage: n = ZZ.random_element(5)
2042 sage: M = matrix.random(QQ, max(0,n-1), algorithm='unimodular')
2043 sage: B = M.transpose()*M
2044 sage: J = BilinearFormEJA(n, B=B)
2045 sage: eis = VectorSpace(M.base_ring(), M.ncols()).basis()
2046 sage: V = J.vector_space()
2047 sage: sis = [ J.from_vector(V([0] + (M.inverse()*ei).list()))
2048 ....: for ei in eis ]
2049 sage: actual = [ sis[i]*sis[j]
2050 ....: for i in range(n-1)
2051 ....: for j in range(n-1) ]
2052 sage: expected = [ J.one() if i == j else J.zero()
2053 ....: for i in range(n-1)
2054 ....: for j in range(n-1) ]
2055 sage: actual == expected
2056 True
2057 """
2058 def __init__(self, n, field=AA, B=None, **kwargs):
2059 if B is None:
2060 self._B = matrix.identity(field, max(0,n-1))
2061 else:
2062 self._B = B
2063
2064 V = VectorSpace(field, n)
2065 mult_table = [[V.zero() for j in range(n)] for i in range(n)]
2066 for i in range(n):
2067 for j in range(n):
2068 x = V.gen(i)
2069 y = V.gen(j)
2070 x0 = x[0]
2071 xbar = x[1:]
2072 y0 = y[0]
2073 ybar = y[1:]
2074 z0 = x0*y0 + (self._B*xbar).inner_product(ybar)
2075 zbar = y0*xbar + x0*ybar
2076 z = V([z0] + zbar.list())
2077 mult_table[i][j] = z
2078
2079 # The rank of this algebra is two, unless we're in a
2080 # one-dimensional ambient space (because the rank is bounded
2081 # by the ambient dimension).
2082 fdeja = super(BilinearFormEJA, self)
2083 return fdeja.__init__(field, mult_table, rank=min(n,2), **kwargs)
2084
2085 def inner_product(self, x, y):
2086 r"""
2087 Half of the trace inner product.
2088
2089 This is defined so that the special case of the Jordan spin
2090 algebra gets the usual inner product.
2091
2092 SETUP::
2093
2094 sage: from mjo.eja.eja_algebra import BilinearFormEJA
2095
2096 TESTS:
2097
2098 Ensure that this is one-half of the trace inner-product when
2099 the algebra isn't just the reals (when ``n`` isn't one). This
2100 is in Faraut and Koranyi, and also my "On the symmetry..."
2101 paper::
2102
2103 sage: set_random_seed()
2104 sage: n = ZZ.random_element(2,5)
2105 sage: M = matrix.random(QQ, max(0,n-1), algorithm='unimodular')
2106 sage: B = M.transpose()*M
2107 sage: J = BilinearFormEJA(n, B=B)
2108 sage: x = J.random_element()
2109 sage: y = J.random_element()
2110 sage: x.inner_product(y) == (x*y).trace()/2
2111 True
2112
2113 """
2114 xvec = x.to_vector()
2115 xbar = xvec[1:]
2116 yvec = y.to_vector()
2117 ybar = yvec[1:]
2118 return x[0]*y[0] + (self._B*xbar).inner_product(ybar)
2119
2120
2121 class JordanSpinEJA(BilinearFormEJA):
2122 """
2123 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2124 with the usual inner product and jordan product ``x*y =
2125 (<x,y>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
2126 the reals.
2127
2128 SETUP::
2129
2130 sage: from mjo.eja.eja_algebra import JordanSpinEJA
2131
2132 EXAMPLES:
2133
2134 This multiplication table can be verified by hand::
2135
2136 sage: J = JordanSpinEJA(4)
2137 sage: e0,e1,e2,e3 = J.gens()
2138 sage: e0*e0
2139 e0
2140 sage: e0*e1
2141 e1
2142 sage: e0*e2
2143 e2
2144 sage: e0*e3
2145 e3
2146 sage: e1*e2
2147 0
2148 sage: e1*e3
2149 0
2150 sage: e2*e3
2151 0
2152
2153 We can change the generator prefix::
2154
2155 sage: JordanSpinEJA(2, prefix='B').gens()
2156 (B0, B1)
2157
2158 TESTS:
2159
2160 Ensure that we have the usual inner product on `R^n`::
2161
2162 sage: set_random_seed()
2163 sage: J = JordanSpinEJA.random_instance()
2164 sage: x,y = J.random_elements(2)
2165 sage: X = x.natural_representation()
2166 sage: Y = y.natural_representation()
2167 sage: x.inner_product(y) == J.natural_inner_product(X,Y)
2168 True
2169
2170 """
2171 def __init__(self, n, field=AA, **kwargs):
2172 # This is a special case of the BilinearFormEJA with the identity
2173 # matrix as its bilinear form.
2174 return super(JordanSpinEJA, self).__init__(n, field, **kwargs)
2175
2176
2177 class TrivialEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
2178 """
2179 The trivial Euclidean Jordan algebra consisting of only a zero element.
2180
2181 SETUP::
2182
2183 sage: from mjo.eja.eja_algebra import TrivialEJA
2184
2185 EXAMPLES::
2186
2187 sage: J = TrivialEJA()
2188 sage: J.dimension()
2189 0
2190 sage: J.zero()
2191 0
2192 sage: J.one()
2193 0
2194 sage: 7*J.one()*12*J.one()
2195 0
2196 sage: J.one().inner_product(J.one())
2197 0
2198 sage: J.one().norm()
2199 0
2200 sage: J.one().subalgebra_generated_by()
2201 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
2202 sage: J.rank()
2203 0
2204
2205 """
2206 def __init__(self, field=AA, **kwargs):
2207 mult_table = []
2208 fdeja = super(TrivialEJA, self)
2209 # The rank is zero using my definition, namely the dimension of the
2210 # largest subalgebra generated by any element.
2211 return fdeja.__init__(field, mult_table, rank=0, **kwargs)