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