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