]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: update some docstrings.
[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
9
10 from sage.algebras.finite_dimensional_algebras.finite_dimensional_algebra import FiniteDimensionalAlgebra
11 from sage.algebras.finite_dimensional_algebras.finite_dimensional_algebra_element import FiniteDimensionalAlgebraElement
12 from sage.algebras.quatalg.quaternion_algebra import QuaternionAlgebra
13 from sage.categories.finite_dimensional_algebras_with_basis import FiniteDimensionalAlgebrasWithBasis
14 from sage.functions.other import sqrt
15 from sage.matrix.constructor import matrix
16 from sage.misc.cachefunc import cached_method
17 from sage.misc.prandom import choice
18 from sage.modules.free_module import VectorSpace
19 from sage.modules.free_module_element import vector
20 from sage.rings.integer_ring import ZZ
21 from sage.rings.number_field.number_field import QuadraticField
22 from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing
23 from sage.rings.rational_field import QQ
24 from sage.structure.element import is_Matrix
25 from sage.structure.category_object import normalize_names
26
27 from mjo.eja.eja_operator import FiniteDimensionalEuclideanJordanAlgebraOperator
28
29
30 class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
31 @staticmethod
32 def __classcall_private__(cls,
33 field,
34 mult_table,
35 names='e',
36 assume_associative=False,
37 category=None,
38 rank=None,
39 natural_basis=None):
40 n = len(mult_table)
41 mult_table = [b.base_extend(field) for b in mult_table]
42 for b in mult_table:
43 b.set_immutable()
44 if not (is_Matrix(b) and b.dimensions() == (n, n)):
45 raise ValueError("input is not a multiplication table")
46 mult_table = tuple(mult_table)
47
48 cat = FiniteDimensionalAlgebrasWithBasis(field)
49 cat.or_subcategory(category)
50 if assume_associative:
51 cat = cat.Associative()
52
53 names = normalize_names(n, names)
54
55 fda = super(FiniteDimensionalEuclideanJordanAlgebra, cls)
56 return fda.__classcall__(cls,
57 field,
58 mult_table,
59 assume_associative=assume_associative,
60 names=names,
61 category=cat,
62 rank=rank,
63 natural_basis=natural_basis)
64
65
66 def __init__(self,
67 field,
68 mult_table,
69 names='e',
70 assume_associative=False,
71 category=None,
72 rank=None,
73 natural_basis=None):
74 """
75 SETUP::
76
77 sage: from mjo.eja.eja_algebra import random_eja
78
79 EXAMPLES:
80
81 By definition, Jordan multiplication commutes::
82
83 sage: set_random_seed()
84 sage: J = random_eja()
85 sage: x = J.random_element()
86 sage: y = J.random_element()
87 sage: x*y == y*x
88 True
89
90 """
91 self._rank = rank
92 self._natural_basis = natural_basis
93 self._multiplication_table = mult_table
94 fda = super(FiniteDimensionalEuclideanJordanAlgebra, self)
95 fda.__init__(field,
96 mult_table,
97 names=names,
98 category=category)
99
100
101 def _repr_(self):
102 """
103 Return a string representation of ``self``.
104
105 SETUP::
106
107 sage: from mjo.eja.eja_algebra import JordanSpinEJA
108
109 TESTS:
110
111 Ensure that it says what we think it says::
112
113 sage: JordanSpinEJA(2, field=QQ)
114 Euclidean Jordan algebra of degree 2 over Rational Field
115 sage: JordanSpinEJA(3, field=RDF)
116 Euclidean Jordan algebra of degree 3 over Real Double Field
117
118 """
119 fmt = "Euclidean Jordan algebra of degree {} over {}"
120 return fmt.format(self.degree(), self.base_ring())
121
122
123 def _a_regular_element(self):
124 """
125 Guess a regular element. Needed to compute the basis for our
126 characteristic polynomial coefficients.
127
128 SETUP::
129
130 sage: from mjo.eja.eja_algebra import random_eja
131
132 TESTS:
133
134 Ensure that this hacky method succeeds for every algebra that we
135 know how to construct::
136
137 sage: set_random_seed()
138 sage: J = random_eja()
139 sage: J._a_regular_element().is_regular()
140 True
141
142 """
143 gs = self.gens()
144 z = self.sum( (i+1)*gs[i] for i in range(len(gs)) )
145 if not z.is_regular():
146 raise ValueError("don't know a regular element")
147 return z
148
149
150 @cached_method
151 def _charpoly_basis_space(self):
152 """
153 Return the vector space spanned by the basis used in our
154 characteristic polynomial coefficients. This is used not only to
155 compute those coefficients, but also any time we need to
156 evaluate the coefficients (like when we compute the trace or
157 determinant).
158 """
159 z = self._a_regular_element()
160 V = self.vector_space()
161 V1 = V.span_of_basis( (z**k).vector() for k in range(self.rank()) )
162 b = (V1.basis() + V1.complement().basis())
163 return V.span_of_basis(b)
164
165
166 @cached_method
167 def _charpoly_coeff(self, i):
168 """
169 Return the coefficient polynomial "a_{i}" of this algebra's
170 general characteristic polynomial.
171
172 Having this be a separate cached method lets us compute and
173 store the trace/determinant (a_{r-1} and a_{0} respectively)
174 separate from the entire characteristic polynomial.
175 """
176 (A_of_x, x, xr, detA) = self._charpoly_matrix_system()
177 R = A_of_x.base_ring()
178 if i >= self.rank():
179 # Guaranteed by theory
180 return R.zero()
181
182 # Danger: the in-place modification is done for performance
183 # reasons (reconstructing a matrix with huge polynomial
184 # entries is slow), but I don't know how cached_method works,
185 # so it's highly possible that we're modifying some global
186 # list variable by reference, here. In other words, you
187 # probably shouldn't call this method twice on the same
188 # algebra, at the same time, in two threads
189 Ai_orig = A_of_x.column(i)
190 A_of_x.set_column(i,xr)
191 numerator = A_of_x.det()
192 A_of_x.set_column(i,Ai_orig)
193
194 # We're relying on the theory here to ensure that each a_i is
195 # indeed back in R, and the added negative signs are to make
196 # the whole charpoly expression sum to zero.
197 return R(-numerator/detA)
198
199
200 @cached_method
201 def _charpoly_matrix_system(self):
202 """
203 Compute the matrix whose entries A_ij are polynomials in
204 X1,...,XN, the vector ``x`` of variables X1,...,XN, the vector
205 corresponding to `x^r` and the determinent of the matrix A =
206 [A_ij]. In other words, all of the fixed (cachable) data needed
207 to compute the coefficients of the characteristic polynomial.
208 """
209 r = self.rank()
210 n = self.dimension()
211
212 # Construct a new algebra over a multivariate polynomial ring...
213 names = ['X' + str(i) for i in range(1,n+1)]
214 R = PolynomialRing(self.base_ring(), names)
215 J = FiniteDimensionalEuclideanJordanAlgebra(R,
216 self._multiplication_table,
217 rank=r)
218
219 idmat = matrix.identity(J.base_ring(), n)
220
221 W = self._charpoly_basis_space()
222 W = W.change_ring(R.fraction_field())
223
224 # Starting with the standard coordinates x = (X1,X2,...,Xn)
225 # and then converting the entries to W-coordinates allows us
226 # to pass in the standard coordinates to the charpoly and get
227 # back the right answer. Specifically, with x = (X1,X2,...,Xn),
228 # we have
229 #
230 # W.coordinates(x^2) eval'd at (standard z-coords)
231 # =
232 # W-coords of (z^2)
233 # =
234 # W-coords of (standard coords of x^2 eval'd at std-coords of z)
235 #
236 # We want the middle equivalent thing in our matrix, but use
237 # the first equivalent thing instead so that we can pass in
238 # standard coordinates.
239 x = J(W(R.gens()))
240 l1 = [matrix.column(W.coordinates((x**k).vector())) for k in range(r)]
241 l2 = [idmat.column(k-1).column() for k in range(r+1, n+1)]
242 A_of_x = matrix.block(R, 1, n, (l1 + l2))
243 xr = W.coordinates((x**r).vector())
244 return (A_of_x, x, xr, A_of_x.det())
245
246
247 @cached_method
248 def characteristic_polynomial(self):
249 """
250 Return a characteristic polynomial that works for all elements
251 of this algebra.
252
253 The resulting polynomial has `n+1` variables, where `n` is the
254 dimension of this algebra. The first `n` variables correspond to
255 the coordinates of an algebra element: when evaluated at the
256 coordinates of an algebra element with respect to a certain
257 basis, the result is a univariate polynomial (in the one
258 remaining variable ``t``), namely the characteristic polynomial
259 of that element.
260
261 SETUP::
262
263 sage: from mjo.eja.eja_algebra import JordanSpinEJA
264
265 EXAMPLES:
266
267 The characteristic polynomial in the spin algebra is given in
268 Alizadeh, Example 11.11::
269
270 sage: J = JordanSpinEJA(3)
271 sage: p = J.characteristic_polynomial(); p
272 X1^2 - X2^2 - X3^2 + (-2*t)*X1 + t^2
273 sage: xvec = J.one().vector()
274 sage: p(*xvec)
275 t^2 - 2*t + 1
276
277 """
278 r = self.rank()
279 n = self.dimension()
280
281 # The list of coefficient polynomials a_1, a_2, ..., a_n.
282 a = [ self._charpoly_coeff(i) for i in range(n) ]
283
284 # We go to a bit of trouble here to reorder the
285 # indeterminates, so that it's easier to evaluate the
286 # characteristic polynomial at x's coordinates and get back
287 # something in terms of t, which is what we want.
288 R = a[0].parent()
289 S = PolynomialRing(self.base_ring(),'t')
290 t = S.gen(0)
291 S = PolynomialRing(S, R.variable_names())
292 t = S(t)
293
294 # Note: all entries past the rth should be zero. The
295 # coefficient of the highest power (x^r) is 1, but it doesn't
296 # appear in the solution vector which contains coefficients
297 # for the other powers (to make them sum to x^r).
298 if (r < n):
299 a[r] = 1 # corresponds to x^r
300 else:
301 # When the rank is equal to the dimension, trying to
302 # assign a[r] goes out-of-bounds.
303 a.append(1) # corresponds to x^r
304
305 return sum( a[k]*(t**k) for k in range(len(a)) )
306
307
308 def inner_product(self, x, y):
309 """
310 The inner product associated with this Euclidean Jordan algebra.
311
312 Defaults to the trace inner product, but can be overridden by
313 subclasses if they are sure that the necessary properties are
314 satisfied.
315
316 SETUP::
317
318 sage: from mjo.eja.eja_algebra import random_eja
319
320 EXAMPLES:
321
322 The inner product must satisfy its axiom for this algebra to truly
323 be a Euclidean Jordan Algebra::
324
325 sage: set_random_seed()
326 sage: J = random_eja()
327 sage: x = J.random_element()
328 sage: y = J.random_element()
329 sage: z = J.random_element()
330 sage: (x*y).inner_product(z) == y.inner_product(x*z)
331 True
332
333 """
334 if (not x in self) or (not y in self):
335 raise TypeError("arguments must live in this algebra")
336 return x.trace_inner_product(y)
337
338
339 def natural_basis(self):
340 """
341 Return a more-natural representation of this algebra's basis.
342
343 Every finite-dimensional Euclidean Jordan Algebra is a direct
344 sum of five simple algebras, four of which comprise Hermitian
345 matrices. This method returns the original "natural" basis
346 for our underlying vector space. (Typically, the natural basis
347 is used to construct the multiplication table in the first place.)
348
349 Note that this will always return a matrix. The standard basis
350 in `R^n` will be returned as `n`-by-`1` column matrices.
351
352 SETUP::
353
354 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
355 ....: RealSymmetricEJA)
356
357 EXAMPLES::
358
359 sage: J = RealSymmetricEJA(2)
360 sage: J.basis()
361 Family (e0, e1, e2)
362 sage: J.natural_basis()
363 (
364 [1 0] [0 1] [0 0]
365 [0 0], [1 0], [0 1]
366 )
367
368 ::
369
370 sage: J = JordanSpinEJA(2)
371 sage: J.basis()
372 Family (e0, e1)
373 sage: J.natural_basis()
374 (
375 [1] [0]
376 [0], [1]
377 )
378
379 """
380 if self._natural_basis is None:
381 return tuple( b.vector().column() for b in self.basis() )
382 else:
383 return self._natural_basis
384
385
386 def rank(self):
387 """
388 Return the rank of this EJA.
389 """
390 if self._rank is None:
391 raise ValueError("no rank specified at genesis")
392 else:
393 return self._rank
394
395
396 def vector_space(self):
397 """
398 Return the vector space that underlies this algebra.
399
400 SETUP::
401
402 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
403
404 EXAMPLES::
405
406 sage: J = RealSymmetricEJA(2)
407 sage: J.vector_space()
408 Vector space of dimension 3 over Rational Field
409
410 """
411 return self.zero().vector().parent().ambient_vector_space()
412
413
414 class Element(FiniteDimensionalAlgebraElement):
415 """
416 An element of a Euclidean Jordan algebra.
417 """
418
419 def __dir__(self):
420 """
421 Oh man, I should not be doing this. This hides the "disabled"
422 methods ``left_matrix`` and ``matrix`` from introspection;
423 in particular it removes them from tab-completion.
424 """
425 return filter(lambda s: s not in ['left_matrix', 'matrix'],
426 dir(self.__class__) )
427
428
429 def __init__(self, A, elt=None):
430 """
431
432 SETUP::
433
434 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
435
436 EXAMPLES:
437
438 The identity in `S^n` is converted to the identity in the EJA::
439
440 sage: J = RealSymmetricEJA(3)
441 sage: I = matrix.identity(QQ,3)
442 sage: J(I) == J.one()
443 True
444
445 This skew-symmetric matrix can't be represented in the EJA::
446
447 sage: J = RealSymmetricEJA(3)
448 sage: A = matrix(QQ,3, lambda i,j: i-j)
449 sage: J(A)
450 Traceback (most recent call last):
451 ...
452 ArithmeticError: vector is not in free module
453
454 """
455 # Goal: if we're given a matrix, and if it lives in our
456 # parent algebra's "natural ambient space," convert it
457 # into an algebra element.
458 #
459 # The catch is, we make a recursive call after converting
460 # the given matrix into a vector that lives in the algebra.
461 # This we need to try the parent class initializer first,
462 # to avoid recursing forever if we're given something that
463 # already fits into the algebra, but also happens to live
464 # in the parent's "natural ambient space" (this happens with
465 # vectors in R^n).
466 try:
467 FiniteDimensionalAlgebraElement.__init__(self, A, elt)
468 except ValueError:
469 natural_basis = A.natural_basis()
470 if elt in natural_basis[0].matrix_space():
471 # Thanks for nothing! Matrix spaces aren't vector
472 # spaces in Sage, so we have to figure out its
473 # natural-basis coordinates ourselves.
474 V = VectorSpace(elt.base_ring(), elt.nrows()**2)
475 W = V.span( _mat2vec(s) for s in natural_basis )
476 coords = W.coordinates(_mat2vec(elt))
477 FiniteDimensionalAlgebraElement.__init__(self, A, coords)
478
479 def __pow__(self, n):
480 """
481 Return ``self`` raised to the power ``n``.
482
483 Jordan algebras are always power-associative; see for
484 example Faraut and Koranyi, Proposition II.1.2 (ii).
485
486 .. WARNING:
487
488 We have to override this because our superclass uses row vectors
489 instead of column vectors! We, on the other hand, assume column
490 vectors everywhere.
491
492 SETUP::
493
494 sage: from mjo.eja.eja_algebra import random_eja
495
496 EXAMPLES::
497
498 sage: set_random_seed()
499 sage: x = random_eja().random_element()
500 sage: x.operator()(x) == (x^2)
501 True
502
503 A few examples of power-associativity::
504
505 sage: set_random_seed()
506 sage: x = random_eja().random_element()
507 sage: x*(x*x)*(x*x) == x^5
508 True
509 sage: (x*x)*(x*x*x) == x^5
510 True
511
512 We also know that powers operator-commute (Koecher, Chapter
513 III, Corollary 1)::
514
515 sage: set_random_seed()
516 sage: x = random_eja().random_element()
517 sage: m = ZZ.random_element(0,10)
518 sage: n = ZZ.random_element(0,10)
519 sage: Lxm = (x^m).operator()
520 sage: Lxn = (x^n).operator()
521 sage: Lxm*Lxn == Lxn*Lxm
522 True
523
524 """
525 if n == 0:
526 return self.parent().one()
527 elif n == 1:
528 return self
529 else:
530 return (self.operator()**(n-1))(self)
531
532
533 def apply_univariate_polynomial(self, p):
534 """
535 Apply the univariate polynomial ``p`` to this element.
536
537 A priori, SageMath won't allow us to apply a univariate
538 polynomial to an element of an EJA, because we don't know
539 that EJAs are rings (they are usually not associative). Of
540 course, we know that EJAs are power-associative, so the
541 operation is ultimately kosher. This function sidesteps
542 the CAS to get the answer we want and expect.
543
544 SETUP::
545
546 sage: from mjo.eja.eja_algebra import (RealCartesianProductEJA,
547 ....: random_eja)
548
549 EXAMPLES::
550
551 sage: R = PolynomialRing(QQ, 't')
552 sage: t = R.gen(0)
553 sage: p = t^4 - t^3 + 5*t - 2
554 sage: J = RealCartesianProductEJA(5)
555 sage: J.one().apply_univariate_polynomial(p) == 3*J.one()
556 True
557
558 TESTS:
559
560 We should always get back an element of the algebra::
561
562 sage: set_random_seed()
563 sage: p = PolynomialRing(QQ, 't').random_element()
564 sage: J = random_eja()
565 sage: x = J.random_element()
566 sage: x.apply_univariate_polynomial(p) in J
567 True
568
569 """
570 if len(p.variables()) > 1:
571 raise ValueError("not a univariate polynomial")
572 P = self.parent()
573 R = P.base_ring()
574 # Convert the coeficcients to the parent's base ring,
575 # because a priori they might live in an (unnecessarily)
576 # larger ring for which P.sum() would fail below.
577 cs = [ R(c) for c in p.coefficients(sparse=False) ]
578 return P.sum( cs[k]*(self**k) for k in range(len(cs)) )
579
580
581 def characteristic_polynomial(self):
582 """
583 Return the characteristic polynomial of this element.
584
585 SETUP::
586
587 sage: from mjo.eja.eja_algebra import RealCartesianProductEJA
588
589 EXAMPLES:
590
591 The rank of `R^3` is three, and the minimal polynomial of
592 the identity element is `(t-1)` from which it follows that
593 the characteristic polynomial should be `(t-1)^3`::
594
595 sage: J = RealCartesianProductEJA(3)
596 sage: J.one().characteristic_polynomial()
597 t^3 - 3*t^2 + 3*t - 1
598
599 Likewise, the characteristic of the zero element in the
600 rank-three algebra `R^{n}` should be `t^{3}`::
601
602 sage: J = RealCartesianProductEJA(3)
603 sage: J.zero().characteristic_polynomial()
604 t^3
605
606 The characteristic polynomial of an element should evaluate
607 to zero on that element::
608
609 sage: set_random_seed()
610 sage: x = RealCartesianProductEJA(3).random_element()
611 sage: p = x.characteristic_polynomial()
612 sage: x.apply_univariate_polynomial(p)
613 0
614
615 """
616 p = self.parent().characteristic_polynomial()
617 return p(*self.vector())
618
619
620 def inner_product(self, other):
621 """
622 Return the parent algebra's inner product of myself and ``other``.
623
624 SETUP::
625
626 sage: from mjo.eja.eja_algebra import (
627 ....: ComplexHermitianEJA,
628 ....: JordanSpinEJA,
629 ....: QuaternionHermitianEJA,
630 ....: RealSymmetricEJA,
631 ....: random_eja)
632
633 EXAMPLES:
634
635 The inner product in the Jordan spin algebra is the usual
636 inner product on `R^n` (this example only works because the
637 basis for the Jordan algebra is the standard basis in `R^n`)::
638
639 sage: J = JordanSpinEJA(3)
640 sage: x = vector(QQ,[1,2,3])
641 sage: y = vector(QQ,[4,5,6])
642 sage: x.inner_product(y)
643 32
644 sage: J(x).inner_product(J(y))
645 32
646
647 The inner product on `S^n` is `<X,Y> = trace(X*Y)`, where
648 multiplication is the usual matrix multiplication in `S^n`,
649 so the inner product of the identity matrix with itself
650 should be the `n`::
651
652 sage: J = RealSymmetricEJA(3)
653 sage: J.one().inner_product(J.one())
654 3
655
656 Likewise, the inner product on `C^n` is `<X,Y> =
657 Re(trace(X*Y))`, where we must necessarily take the real
658 part because the product of Hermitian matrices may not be
659 Hermitian::
660
661 sage: J = ComplexHermitianEJA(3)
662 sage: J.one().inner_product(J.one())
663 3
664
665 Ditto for the quaternions::
666
667 sage: J = QuaternionHermitianEJA(3)
668 sage: J.one().inner_product(J.one())
669 3
670
671 TESTS:
672
673 Ensure that we can always compute an inner product, and that
674 it gives us back a real number::
675
676 sage: set_random_seed()
677 sage: J = random_eja()
678 sage: x = J.random_element()
679 sage: y = J.random_element()
680 sage: x.inner_product(y) in RR
681 True
682
683 """
684 P = self.parent()
685 if not other in P:
686 raise TypeError("'other' must live in the same algebra")
687
688 return P.inner_product(self, other)
689
690
691 def operator_commutes_with(self, other):
692 """
693 Return whether or not this element operator-commutes
694 with ``other``.
695
696 SETUP::
697
698 sage: from mjo.eja.eja_algebra import random_eja
699
700 EXAMPLES:
701
702 The definition of a Jordan algebra says that any element
703 operator-commutes with its square::
704
705 sage: set_random_seed()
706 sage: x = random_eja().random_element()
707 sage: x.operator_commutes_with(x^2)
708 True
709
710 TESTS:
711
712 Test Lemma 1 from Chapter III of Koecher::
713
714 sage: set_random_seed()
715 sage: J = random_eja()
716 sage: u = J.random_element()
717 sage: v = J.random_element()
718 sage: lhs = u.operator_commutes_with(u*v)
719 sage: rhs = v.operator_commutes_with(u^2)
720 sage: lhs == rhs
721 True
722
723 Test the first polarization identity from my notes, Koecher Chapter
724 III, or from Baes (2.3)::
725
726 sage: set_random_seed()
727 sage: J = random_eja()
728 sage: x = J.random_element()
729 sage: y = J.random_element()
730 sage: Lx = x.operator()
731 sage: Ly = y.operator()
732 sage: Lxx = (x*x).operator()
733 sage: Lxy = (x*y).operator()
734 sage: bool(2*Lx*Lxy + Ly*Lxx == 2*Lxy*Lx + Lxx*Ly)
735 True
736
737 Test the second polarization identity from my notes or from
738 Baes (2.4)::
739
740 sage: set_random_seed()
741 sage: J = random_eja()
742 sage: x = J.random_element()
743 sage: y = J.random_element()
744 sage: z = J.random_element()
745 sage: Lx = x.operator()
746 sage: Ly = y.operator()
747 sage: Lz = z.operator()
748 sage: Lzy = (z*y).operator()
749 sage: Lxy = (x*y).operator()
750 sage: Lxz = (x*z).operator()
751 sage: bool(Lx*Lzy + Lz*Lxy + Ly*Lxz == Lzy*Lx + Lxy*Lz + Lxz*Ly)
752 True
753
754 Test the third polarization identity from my notes or from
755 Baes (2.5)::
756
757 sage: set_random_seed()
758 sage: J = random_eja()
759 sage: u = J.random_element()
760 sage: y = J.random_element()
761 sage: z = J.random_element()
762 sage: Lu = u.operator()
763 sage: Ly = y.operator()
764 sage: Lz = z.operator()
765 sage: Lzy = (z*y).operator()
766 sage: Luy = (u*y).operator()
767 sage: Luz = (u*z).operator()
768 sage: Luyz = (u*(y*z)).operator()
769 sage: lhs = Lu*Lzy + Lz*Luy + Ly*Luz
770 sage: rhs = Luyz + Ly*Lu*Lz + Lz*Lu*Ly
771 sage: bool(lhs == rhs)
772 True
773
774 """
775 if not other in self.parent():
776 raise TypeError("'other' must live in the same algebra")
777
778 A = self.operator()
779 B = other.operator()
780 return (A*B == B*A)
781
782
783 def det(self):
784 """
785 Return my determinant, the product of my eigenvalues.
786
787 SETUP::
788
789 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
790 ....: random_eja)
791
792 EXAMPLES::
793
794 sage: J = JordanSpinEJA(2)
795 sage: e0,e1 = J.gens()
796 sage: x = sum( J.gens() )
797 sage: x.det()
798 0
799
800 ::
801
802 sage: J = JordanSpinEJA(3)
803 sage: e0,e1,e2 = J.gens()
804 sage: x = sum( J.gens() )
805 sage: x.det()
806 -1
807
808 TESTS:
809
810 An element is invertible if and only if its determinant is
811 non-zero::
812
813 sage: set_random_seed()
814 sage: x = random_eja().random_element()
815 sage: x.is_invertible() == (x.det() != 0)
816 True
817
818 """
819 P = self.parent()
820 r = P.rank()
821 p = P._charpoly_coeff(0)
822 # The _charpoly_coeff function already adds the factor of
823 # -1 to ensure that _charpoly_coeff(0) is really what
824 # appears in front of t^{0} in the charpoly. However,
825 # we want (-1)^r times THAT for the determinant.
826 return ((-1)**r)*p(*self.vector())
827
828
829 def inverse(self):
830 """
831 Return the Jordan-multiplicative inverse of this element.
832
833 ALGORITHM:
834
835 We appeal to the quadratic representation as in Koecher's
836 Theorem 12 in Chapter III, Section 5.
837
838 SETUP::
839
840 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
841 ....: random_eja)
842
843 EXAMPLES:
844
845 The inverse in the spin factor algebra is given in Alizadeh's
846 Example 11.11::
847
848 sage: set_random_seed()
849 sage: n = ZZ.random_element(1,10)
850 sage: J = JordanSpinEJA(n)
851 sage: x = J.random_element()
852 sage: while not x.is_invertible():
853 ....: x = J.random_element()
854 sage: x_vec = x.vector()
855 sage: x0 = x_vec[0]
856 sage: x_bar = x_vec[1:]
857 sage: coeff = ~(x0^2 - x_bar.inner_product(x_bar))
858 sage: inv_vec = x_vec.parent()([x0] + (-x_bar).list())
859 sage: x_inverse = coeff*inv_vec
860 sage: x.inverse() == J(x_inverse)
861 True
862
863 TESTS:
864
865 The identity element is its own inverse::
866
867 sage: set_random_seed()
868 sage: J = random_eja()
869 sage: J.one().inverse() == J.one()
870 True
871
872 If an element has an inverse, it acts like one::
873
874 sage: set_random_seed()
875 sage: J = random_eja()
876 sage: x = J.random_element()
877 sage: (not x.is_invertible()) or (x.inverse()*x == J.one())
878 True
879
880 The inverse of the inverse is what we started with::
881
882 sage: set_random_seed()
883 sage: J = random_eja()
884 sage: x = J.random_element()
885 sage: (not x.is_invertible()) or (x.inverse().inverse() == x)
886 True
887
888 The zero element is never invertible::
889
890 sage: set_random_seed()
891 sage: J = random_eja().zero().inverse()
892 Traceback (most recent call last):
893 ...
894 ValueError: element is not invertible
895
896 """
897 if not self.is_invertible():
898 raise ValueError("element is not invertible")
899
900 return (~self.quadratic_representation())(self)
901
902
903 def is_invertible(self):
904 """
905 Return whether or not this element is invertible.
906
907 We can't use the superclass method because it relies on
908 the algebra being associative.
909
910 ALGORITHM:
911
912 The usual way to do this is to check if the determinant is
913 zero, but we need the characteristic polynomial for the
914 determinant. The minimal polynomial is a lot easier to get,
915 so we use Corollary 2 in Chapter V of Koecher to check
916 whether or not the paren't algebra's zero element is a root
917 of this element's minimal polynomial.
918
919 SETUP::
920
921 sage: from mjo.eja.eja_algebra import random_eja
922
923 TESTS:
924
925 The identity element is always invertible::
926
927 sage: set_random_seed()
928 sage: J = random_eja()
929 sage: J.one().is_invertible()
930 True
931
932 The zero element is never invertible::
933
934 sage: set_random_seed()
935 sage: J = random_eja()
936 sage: J.zero().is_invertible()
937 False
938
939 """
940 zero = self.parent().zero()
941 p = self.minimal_polynomial()
942 return not (p(zero) == zero)
943
944
945 def is_nilpotent(self):
946 """
947 Return whether or not some power of this element is zero.
948
949 The superclass method won't work unless we're in an
950 associative algebra, and we aren't. However, we generate
951 an assocoative subalgebra and we're nilpotent there if and
952 only if we're nilpotent here (probably).
953
954 SETUP::
955
956 sage: from mjo.eja.eja_algebra import random_eja
957
958 TESTS:
959
960 The identity element is never nilpotent::
961
962 sage: set_random_seed()
963 sage: random_eja().one().is_nilpotent()
964 False
965
966 The additive identity is always nilpotent::
967
968 sage: set_random_seed()
969 sage: random_eja().zero().is_nilpotent()
970 True
971
972 """
973 # The element we're going to call "is_nilpotent()" on.
974 # Either myself, interpreted as an element of a finite-
975 # dimensional algebra, or an element of an associative
976 # subalgebra.
977 elt = None
978
979 if self.parent().is_associative():
980 elt = FiniteDimensionalAlgebraElement(self.parent(), self)
981 else:
982 V = self.span_of_powers()
983 assoc_subalg = self.subalgebra_generated_by()
984 # Mis-design warning: the basis used for span_of_powers()
985 # and subalgebra_generated_by() must be the same, and in
986 # the same order!
987 elt = assoc_subalg(V.coordinates(self.vector()))
988
989 # Recursive call, but should work since elt lives in an
990 # associative algebra.
991 return elt.is_nilpotent()
992
993
994 def is_regular(self):
995 """
996 Return whether or not this is a regular element.
997
998 SETUP::
999
1000 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1001
1002 EXAMPLES:
1003
1004 The identity element always has degree one, but any element
1005 linearly-independent from it is regular::
1006
1007 sage: J = JordanSpinEJA(5)
1008 sage: J.one().is_regular()
1009 False
1010 sage: e0, e1, e2, e3, e4 = J.gens() # e0 is the identity
1011 sage: for x in J.gens():
1012 ....: (J.one() + x).is_regular()
1013 False
1014 True
1015 True
1016 True
1017 True
1018
1019 """
1020 return self.degree() == self.parent().rank()
1021
1022
1023 def degree(self):
1024 """
1025 Compute the degree of this element the straightforward way
1026 according to the definition; by appending powers to a list
1027 and figuring out its dimension (that is, whether or not
1028 they're linearly dependent).
1029
1030 SETUP::
1031
1032 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1033
1034 EXAMPLES::
1035
1036 sage: J = JordanSpinEJA(4)
1037 sage: J.one().degree()
1038 1
1039 sage: e0,e1,e2,e3 = J.gens()
1040 sage: (e0 - e1).degree()
1041 2
1042
1043 In the spin factor algebra (of rank two), all elements that
1044 aren't multiples of the identity are regular::
1045
1046 sage: set_random_seed()
1047 sage: n = ZZ.random_element(1,10)
1048 sage: J = JordanSpinEJA(n)
1049 sage: x = J.random_element()
1050 sage: x == x.coefficient(0)*J.one() or x.degree() == 2
1051 True
1052
1053 """
1054 return self.span_of_powers().dimension()
1055
1056
1057 def left_matrix(self):
1058 """
1059 Our parent class defines ``left_matrix`` and ``matrix``
1060 methods whose names are misleading. We don't want them.
1061 """
1062 raise NotImplementedError("use operator().matrix() instead")
1063
1064 matrix = left_matrix
1065
1066
1067 def minimal_polynomial(self):
1068 """
1069 Return the minimal polynomial of this element,
1070 as a function of the variable `t`.
1071
1072 ALGORITHM:
1073
1074 We restrict ourselves to the associative subalgebra
1075 generated by this element, and then return the minimal
1076 polynomial of this element's operator matrix (in that
1077 subalgebra). This works by Baes Proposition 2.3.16.
1078
1079 SETUP::
1080
1081 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
1082 ....: random_eja)
1083
1084 TESTS:
1085
1086 The minimal polynomial of the identity and zero elements are
1087 always the same::
1088
1089 sage: set_random_seed()
1090 sage: J = random_eja()
1091 sage: J.one().minimal_polynomial()
1092 t - 1
1093 sage: J.zero().minimal_polynomial()
1094 t
1095
1096 The degree of an element is (by one definition) the degree
1097 of its minimal polynomial::
1098
1099 sage: set_random_seed()
1100 sage: x = random_eja().random_element()
1101 sage: x.degree() == x.minimal_polynomial().degree()
1102 True
1103
1104 The minimal polynomial and the characteristic polynomial coincide
1105 and are known (see Alizadeh, Example 11.11) for all elements of
1106 the spin factor algebra that aren't scalar multiples of the
1107 identity::
1108
1109 sage: set_random_seed()
1110 sage: n = ZZ.random_element(2,10)
1111 sage: J = JordanSpinEJA(n)
1112 sage: y = J.random_element()
1113 sage: while y == y.coefficient(0)*J.one():
1114 ....: y = J.random_element()
1115 sage: y0 = y.vector()[0]
1116 sage: y_bar = y.vector()[1:]
1117 sage: actual = y.minimal_polynomial()
1118 sage: t = PolynomialRing(J.base_ring(),'t').gen(0)
1119 sage: expected = t^2 - 2*y0*t + (y0^2 - norm(y_bar)^2)
1120 sage: bool(actual == expected)
1121 True
1122
1123 The minimal polynomial should always kill its element::
1124
1125 sage: set_random_seed()
1126 sage: x = random_eja().random_element()
1127 sage: p = x.minimal_polynomial()
1128 sage: x.apply_univariate_polynomial(p)
1129 0
1130
1131 """
1132 V = self.span_of_powers()
1133 assoc_subalg = self.subalgebra_generated_by()
1134 # Mis-design warning: the basis used for span_of_powers()
1135 # and subalgebra_generated_by() must be the same, and in
1136 # the same order!
1137 elt = assoc_subalg(V.coordinates(self.vector()))
1138 return elt.operator().minimal_polynomial()
1139
1140
1141
1142 def natural_representation(self):
1143 """
1144 Return a more-natural representation of this element.
1145
1146 Every finite-dimensional Euclidean Jordan Algebra is a
1147 direct sum of five simple algebras, four of which comprise
1148 Hermitian matrices. This method returns the original
1149 "natural" representation of this element as a Hermitian
1150 matrix, if it has one. If not, you get the usual representation.
1151
1152 SETUP::
1153
1154 sage: from mjo.eja.eja_algebra import (ComplexHermitianEJA,
1155 ....: QuaternionHermitianEJA)
1156
1157 EXAMPLES::
1158
1159 sage: J = ComplexHermitianEJA(3)
1160 sage: J.one()
1161 e0 + e5 + e8
1162 sage: J.one().natural_representation()
1163 [1 0 0 0 0 0]
1164 [0 1 0 0 0 0]
1165 [0 0 1 0 0 0]
1166 [0 0 0 1 0 0]
1167 [0 0 0 0 1 0]
1168 [0 0 0 0 0 1]
1169
1170 ::
1171
1172 sage: J = QuaternionHermitianEJA(3)
1173 sage: J.one()
1174 e0 + e9 + e14
1175 sage: J.one().natural_representation()
1176 [1 0 0 0 0 0 0 0 0 0 0 0]
1177 [0 1 0 0 0 0 0 0 0 0 0 0]
1178 [0 0 1 0 0 0 0 0 0 0 0 0]
1179 [0 0 0 1 0 0 0 0 0 0 0 0]
1180 [0 0 0 0 1 0 0 0 0 0 0 0]
1181 [0 0 0 0 0 1 0 0 0 0 0 0]
1182 [0 0 0 0 0 0 1 0 0 0 0 0]
1183 [0 0 0 0 0 0 0 1 0 0 0 0]
1184 [0 0 0 0 0 0 0 0 1 0 0 0]
1185 [0 0 0 0 0 0 0 0 0 1 0 0]
1186 [0 0 0 0 0 0 0 0 0 0 1 0]
1187 [0 0 0 0 0 0 0 0 0 0 0 1]
1188
1189 """
1190 B = self.parent().natural_basis()
1191 W = B[0].matrix_space()
1192 return W.linear_combination(zip(self.vector(), B))
1193
1194
1195 def operator(self):
1196 """
1197 Return the left-multiplication-by-this-element
1198 operator on the ambient algebra.
1199
1200 SETUP::
1201
1202 sage: from mjo.eja.eja_algebra import random_eja
1203
1204 TESTS::
1205
1206 sage: set_random_seed()
1207 sage: J = random_eja()
1208 sage: x = J.random_element()
1209 sage: y = J.random_element()
1210 sage: x.operator()(y) == x*y
1211 True
1212 sage: y.operator()(x) == x*y
1213 True
1214
1215 """
1216 P = self.parent()
1217 fda_elt = FiniteDimensionalAlgebraElement(P, self)
1218 return FiniteDimensionalEuclideanJordanAlgebraOperator(
1219 P,
1220 P,
1221 fda_elt.matrix().transpose() )
1222
1223
1224 def quadratic_representation(self, other=None):
1225 """
1226 Return the quadratic representation of this element.
1227
1228 SETUP::
1229
1230 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
1231 ....: random_eja)
1232
1233 EXAMPLES:
1234
1235 The explicit form in the spin factor algebra is given by
1236 Alizadeh's Example 11.12::
1237
1238 sage: set_random_seed()
1239 sage: n = ZZ.random_element(1,10)
1240 sage: J = JordanSpinEJA(n)
1241 sage: x = J.random_element()
1242 sage: x_vec = x.vector()
1243 sage: x0 = x_vec[0]
1244 sage: x_bar = x_vec[1:]
1245 sage: A = matrix(QQ, 1, [x_vec.inner_product(x_vec)])
1246 sage: B = 2*x0*x_bar.row()
1247 sage: C = 2*x0*x_bar.column()
1248 sage: D = matrix.identity(QQ, n-1)
1249 sage: D = (x0^2 - x_bar.inner_product(x_bar))*D
1250 sage: D = D + 2*x_bar.tensor_product(x_bar)
1251 sage: Q = matrix.block(2,2,[A,B,C,D])
1252 sage: Q == x.quadratic_representation().matrix()
1253 True
1254
1255 Test all of the properties from Theorem 11.2 in Alizadeh::
1256
1257 sage: set_random_seed()
1258 sage: J = random_eja()
1259 sage: x = J.random_element()
1260 sage: y = J.random_element()
1261 sage: Lx = x.operator()
1262 sage: Lxx = (x*x).operator()
1263 sage: Qx = x.quadratic_representation()
1264 sage: Qy = y.quadratic_representation()
1265 sage: Qxy = x.quadratic_representation(y)
1266 sage: Qex = J.one().quadratic_representation(x)
1267 sage: n = ZZ.random_element(10)
1268 sage: Qxn = (x^n).quadratic_representation()
1269
1270 Property 1:
1271
1272 sage: 2*Qxy == (x+y).quadratic_representation() - Qx - Qy
1273 True
1274
1275 Property 2 (multiply on the right for :trac:`28272`):
1276
1277 sage: alpha = QQ.random_element()
1278 sage: (alpha*x).quadratic_representation() == Qx*(alpha^2)
1279 True
1280
1281 Property 3:
1282
1283 sage: not x.is_invertible() or ( Qx(x.inverse()) == x )
1284 True
1285
1286 sage: not x.is_invertible() or (
1287 ....: ~Qx
1288 ....: ==
1289 ....: x.inverse().quadratic_representation() )
1290 True
1291
1292 sage: Qxy(J.one()) == x*y
1293 True
1294
1295 Property 4:
1296
1297 sage: not x.is_invertible() or (
1298 ....: x.quadratic_representation(x.inverse())*Qx
1299 ....: == Qx*x.quadratic_representation(x.inverse()) )
1300 True
1301
1302 sage: not x.is_invertible() or (
1303 ....: x.quadratic_representation(x.inverse())*Qx
1304 ....: ==
1305 ....: 2*x.operator()*Qex - Qx )
1306 True
1307
1308 sage: 2*x.operator()*Qex - Qx == Lxx
1309 True
1310
1311 Property 5:
1312
1313 sage: Qy(x).quadratic_representation() == Qy*Qx*Qy
1314 True
1315
1316 Property 6:
1317
1318 sage: Qxn == (Qx)^n
1319 True
1320
1321 Property 7:
1322
1323 sage: not x.is_invertible() or (
1324 ....: Qx*x.inverse().operator() == Lx )
1325 True
1326
1327 Property 8:
1328
1329 sage: not x.operator_commutes_with(y) or (
1330 ....: Qx(y)^n == Qxn(y^n) )
1331 True
1332
1333 """
1334 if other is None:
1335 other=self
1336 elif not other in self.parent():
1337 raise TypeError("'other' must live in the same algebra")
1338
1339 L = self.operator()
1340 M = other.operator()
1341 return ( L*M + M*L - (self*other).operator() )
1342
1343
1344 def span_of_powers(self):
1345 """
1346 Return the vector space spanned by successive powers of
1347 this element.
1348 """
1349 # The dimension of the subalgebra can't be greater than
1350 # the big algebra, so just put everything into a list
1351 # and let span() get rid of the excess.
1352 #
1353 # We do the extra ambient_vector_space() in case we're messing
1354 # with polynomials and the direct parent is a module.
1355 V = self.parent().vector_space()
1356 return V.span( (self**d).vector() for d in xrange(V.dimension()) )
1357
1358
1359 def subalgebra_generated_by(self):
1360 """
1361 Return the associative subalgebra of the parent EJA generated
1362 by this element.
1363
1364 SETUP::
1365
1366 sage: from mjo.eja.eja_algebra import random_eja
1367
1368 TESTS::
1369
1370 sage: set_random_seed()
1371 sage: x = random_eja().random_element()
1372 sage: x.subalgebra_generated_by().is_associative()
1373 True
1374
1375 Squaring in the subalgebra should work the same as in
1376 the superalgebra::
1377
1378 sage: set_random_seed()
1379 sage: x = random_eja().random_element()
1380 sage: u = x.subalgebra_generated_by().random_element()
1381 sage: u.operator()(u) == u^2
1382 True
1383
1384 """
1385 # First get the subspace spanned by the powers of myself...
1386 V = self.span_of_powers()
1387 F = self.base_ring()
1388
1389 # Now figure out the entries of the right-multiplication
1390 # matrix for the successive basis elements b0, b1,... of
1391 # that subspace.
1392 mats = []
1393 for b_right in V.basis():
1394 eja_b_right = self.parent()(b_right)
1395 b_right_rows = []
1396 # The first row of the right-multiplication matrix by
1397 # b1 is what we get if we apply that matrix to b1. The
1398 # second row of the right multiplication matrix by b1
1399 # is what we get when we apply that matrix to b2...
1400 #
1401 # IMPORTANT: this assumes that all vectors are COLUMN
1402 # vectors, unlike our superclass (which uses row vectors).
1403 for b_left in V.basis():
1404 eja_b_left = self.parent()(b_left)
1405 # Multiply in the original EJA, but then get the
1406 # coordinates from the subalgebra in terms of its
1407 # basis.
1408 this_row = V.coordinates((eja_b_left*eja_b_right).vector())
1409 b_right_rows.append(this_row)
1410 b_right_matrix = matrix(F, b_right_rows)
1411 mats.append(b_right_matrix)
1412
1413 # It's an algebra of polynomials in one element, and EJAs
1414 # are power-associative.
1415 #
1416 # TODO: choose generator names intelligently.
1417 return FiniteDimensionalEuclideanJordanAlgebra(F, mats, assume_associative=True, names='f')
1418
1419
1420 def subalgebra_idempotent(self):
1421 """
1422 Find an idempotent in the associative subalgebra I generate
1423 using Proposition 2.3.5 in Baes.
1424
1425 SETUP::
1426
1427 sage: from mjo.eja.eja_algebra import random_eja
1428
1429 TESTS::
1430
1431 sage: set_random_seed()
1432 sage: J = random_eja()
1433 sage: x = J.random_element()
1434 sage: while x.is_nilpotent():
1435 ....: x = J.random_element()
1436 sage: c = x.subalgebra_idempotent()
1437 sage: c^2 == c
1438 True
1439
1440 """
1441 if self.is_nilpotent():
1442 raise ValueError("this only works with non-nilpotent elements!")
1443
1444 V = self.span_of_powers()
1445 J = self.subalgebra_generated_by()
1446 # Mis-design warning: the basis used for span_of_powers()
1447 # and subalgebra_generated_by() must be the same, and in
1448 # the same order!
1449 u = J(V.coordinates(self.vector()))
1450
1451 # The image of the matrix of left-u^m-multiplication
1452 # will be minimal for some natural number s...
1453 s = 0
1454 minimal_dim = V.dimension()
1455 for i in xrange(1, V.dimension()):
1456 this_dim = (u**i).operator().matrix().image().dimension()
1457 if this_dim < minimal_dim:
1458 minimal_dim = this_dim
1459 s = i
1460
1461 # Now minimal_matrix should correspond to the smallest
1462 # non-zero subspace in Baes's (or really, Koecher's)
1463 # proposition.
1464 #
1465 # However, we need to restrict the matrix to work on the
1466 # subspace... or do we? Can't we just solve, knowing that
1467 # A(c) = u^(s+1) should have a solution in the big space,
1468 # too?
1469 #
1470 # Beware, solve_right() means that we're using COLUMN vectors.
1471 # Our FiniteDimensionalAlgebraElement superclass uses rows.
1472 u_next = u**(s+1)
1473 A = u_next.operator().matrix()
1474 c_coordinates = A.solve_right(u_next.vector())
1475
1476 # Now c_coordinates is the idempotent we want, but it's in
1477 # the coordinate system of the subalgebra.
1478 #
1479 # We need the basis for J, but as elements of the parent algebra.
1480 #
1481 basis = [self.parent(v) for v in V.basis()]
1482 return self.parent().linear_combination(zip(c_coordinates, basis))
1483
1484
1485 def trace(self):
1486 """
1487 Return my trace, the sum of my eigenvalues.
1488
1489 SETUP::
1490
1491 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
1492 ....: RealCartesianProductEJA,
1493 ....: random_eja)
1494
1495 EXAMPLES::
1496
1497 sage: J = JordanSpinEJA(3)
1498 sage: x = sum(J.gens())
1499 sage: x.trace()
1500 2
1501
1502 ::
1503
1504 sage: J = RealCartesianProductEJA(5)
1505 sage: J.one().trace()
1506 5
1507
1508 TESTS:
1509
1510 The trace of an element is a real number::
1511
1512 sage: set_random_seed()
1513 sage: J = random_eja()
1514 sage: J.random_element().trace() in J.base_ring()
1515 True
1516
1517 """
1518 P = self.parent()
1519 r = P.rank()
1520 p = P._charpoly_coeff(r-1)
1521 # The _charpoly_coeff function already adds the factor of
1522 # -1 to ensure that _charpoly_coeff(r-1) is really what
1523 # appears in front of t^{r-1} in the charpoly. However,
1524 # we want the negative of THAT for the trace.
1525 return -p(*self.vector())
1526
1527
1528 def trace_inner_product(self, other):
1529 """
1530 Return the trace inner product of myself and ``other``.
1531
1532 SETUP::
1533
1534 sage: from mjo.eja.eja_algebra import random_eja
1535
1536 TESTS:
1537
1538 The trace inner product is commutative::
1539
1540 sage: set_random_seed()
1541 sage: J = random_eja()
1542 sage: x = J.random_element(); y = J.random_element()
1543 sage: x.trace_inner_product(y) == y.trace_inner_product(x)
1544 True
1545
1546 The trace inner product is bilinear::
1547
1548 sage: set_random_seed()
1549 sage: J = random_eja()
1550 sage: x = J.random_element()
1551 sage: y = J.random_element()
1552 sage: z = J.random_element()
1553 sage: a = QQ.random_element();
1554 sage: actual = (a*(x+z)).trace_inner_product(y)
1555 sage: expected = ( a*x.trace_inner_product(y) +
1556 ....: a*z.trace_inner_product(y) )
1557 sage: actual == expected
1558 True
1559 sage: actual = x.trace_inner_product(a*(y+z))
1560 sage: expected = ( a*x.trace_inner_product(y) +
1561 ....: a*x.trace_inner_product(z) )
1562 sage: actual == expected
1563 True
1564
1565 The trace inner product satisfies the compatibility
1566 condition in the definition of a Euclidean Jordan algebra::
1567
1568 sage: set_random_seed()
1569 sage: J = random_eja()
1570 sage: x = J.random_element()
1571 sage: y = J.random_element()
1572 sage: z = J.random_element()
1573 sage: (x*y).trace_inner_product(z) == y.trace_inner_product(x*z)
1574 True
1575
1576 """
1577 if not other in self.parent():
1578 raise TypeError("'other' must live in the same algebra")
1579
1580 return (self*other).trace()
1581
1582
1583 class RealCartesianProductEJA(FiniteDimensionalEuclideanJordanAlgebra):
1584 """
1585 Return the Euclidean Jordan Algebra corresponding to the set
1586 `R^n` under the Hadamard product.
1587
1588 Note: this is nothing more than the Cartesian product of ``n``
1589 copies of the spin algebra. Once Cartesian product algebras
1590 are implemented, this can go.
1591
1592 SETUP::
1593
1594 sage: from mjo.eja.eja_algebra import RealCartesianProductEJA
1595
1596 EXAMPLES:
1597
1598 This multiplication table can be verified by hand::
1599
1600 sage: J = RealCartesianProductEJA(3)
1601 sage: e0,e1,e2 = J.gens()
1602 sage: e0*e0
1603 e0
1604 sage: e0*e1
1605 0
1606 sage: e0*e2
1607 0
1608 sage: e1*e1
1609 e1
1610 sage: e1*e2
1611 0
1612 sage: e2*e2
1613 e2
1614
1615 """
1616 @staticmethod
1617 def __classcall_private__(cls, n, field=QQ):
1618 # The FiniteDimensionalAlgebra constructor takes a list of
1619 # matrices, the ith representing right multiplication by the ith
1620 # basis element in the vector space. So if e_1 = (1,0,0), then
1621 # right (Hadamard) multiplication of x by e_1 picks out the first
1622 # component of x; and likewise for the ith basis element e_i.
1623 Qs = [ matrix(field, n, n, lambda k,j: 1*(k == j == i))
1624 for i in xrange(n) ]
1625
1626 fdeja = super(RealCartesianProductEJA, cls)
1627 return fdeja.__classcall_private__(cls, field, Qs, rank=n)
1628
1629 def inner_product(self, x, y):
1630 return _usual_ip(x,y)
1631
1632
1633 def random_eja():
1634 """
1635 Return a "random" finite-dimensional Euclidean Jordan Algebra.
1636
1637 ALGORITHM:
1638
1639 For now, we choose a random natural number ``n`` (greater than zero)
1640 and then give you back one of the following:
1641
1642 * The cartesian product of the rational numbers ``n`` times; this is
1643 ``QQ^n`` with the Hadamard product.
1644
1645 * The Jordan spin algebra on ``QQ^n``.
1646
1647 * The ``n``-by-``n`` rational symmetric matrices with the symmetric
1648 product.
1649
1650 * The ``n``-by-``n`` complex-rational Hermitian matrices embedded
1651 in the space of ``2n``-by-``2n`` real symmetric matrices.
1652
1653 * The ``n``-by-``n`` quaternion-rational Hermitian matrices embedded
1654 in the space of ``4n``-by-``4n`` real symmetric matrices.
1655
1656 Later this might be extended to return Cartesian products of the
1657 EJAs above.
1658
1659 SETUP::
1660
1661 sage: from mjo.eja.eja_algebra import random_eja
1662
1663 TESTS::
1664
1665 sage: random_eja()
1666 Euclidean Jordan algebra of degree...
1667
1668 """
1669
1670 # The max_n component lets us choose different upper bounds on the
1671 # value "n" that gets passed to the constructor. This is needed
1672 # because e.g. R^{10} is reasonable to test, while the Hermitian
1673 # 10-by-10 quaternion matrices are not.
1674 (constructor, max_n) = choice([(RealCartesianProductEJA, 6),
1675 (JordanSpinEJA, 6),
1676 (RealSymmetricEJA, 5),
1677 (ComplexHermitianEJA, 4),
1678 (QuaternionHermitianEJA, 3)])
1679 n = ZZ.random_element(1, max_n)
1680 return constructor(n, field=QQ)
1681
1682
1683
1684 def _real_symmetric_basis(n, field=QQ):
1685 """
1686 Return a basis for the space of real symmetric n-by-n matrices.
1687 """
1688 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1689 # coordinates.
1690 S = []
1691 for i in xrange(n):
1692 for j in xrange(i+1):
1693 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1694 if i == j:
1695 Sij = Eij
1696 else:
1697 # Beware, orthogonal but not normalized!
1698 Sij = Eij + Eij.transpose()
1699 S.append(Sij)
1700 return tuple(S)
1701
1702
1703 def _complex_hermitian_basis(n, field=QQ):
1704 """
1705 Returns a basis for the space of complex Hermitian n-by-n matrices.
1706
1707 SETUP::
1708
1709 sage: from mjo.eja.eja_algebra import _complex_hermitian_basis
1710
1711 TESTS::
1712
1713 sage: set_random_seed()
1714 sage: n = ZZ.random_element(1,5)
1715 sage: all( M.is_symmetric() for M in _complex_hermitian_basis(n) )
1716 True
1717
1718 """
1719 F = QuadraticField(-1, 'I')
1720 I = F.gen()
1721
1722 # This is like the symmetric case, but we need to be careful:
1723 #
1724 # * We want conjugate-symmetry, not just symmetry.
1725 # * The diagonal will (as a result) be real.
1726 #
1727 S = []
1728 for i in xrange(n):
1729 for j in xrange(i+1):
1730 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1731 if i == j:
1732 Sij = _embed_complex_matrix(Eij)
1733 S.append(Sij)
1734 else:
1735 # Beware, orthogonal but not normalized! The second one
1736 # has a minus because it's conjugated.
1737 Sij_real = _embed_complex_matrix(Eij + Eij.transpose())
1738 S.append(Sij_real)
1739 Sij_imag = _embed_complex_matrix(I*Eij - I*Eij.transpose())
1740 S.append(Sij_imag)
1741 return tuple(S)
1742
1743
1744 def _quaternion_hermitian_basis(n, field=QQ):
1745 """
1746 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
1747
1748 SETUP::
1749
1750 sage: from mjo.eja.eja_algebra import _quaternion_hermitian_basis
1751
1752 TESTS::
1753
1754 sage: set_random_seed()
1755 sage: n = ZZ.random_element(1,5)
1756 sage: all( M.is_symmetric() for M in _quaternion_hermitian_basis(n) )
1757 True
1758
1759 """
1760 Q = QuaternionAlgebra(QQ,-1,-1)
1761 I,J,K = Q.gens()
1762
1763 # This is like the symmetric case, but we need to be careful:
1764 #
1765 # * We want conjugate-symmetry, not just symmetry.
1766 # * The diagonal will (as a result) be real.
1767 #
1768 S = []
1769 for i in xrange(n):
1770 for j in xrange(i+1):
1771 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
1772 if i == j:
1773 Sij = _embed_quaternion_matrix(Eij)
1774 S.append(Sij)
1775 else:
1776 # Beware, orthogonal but not normalized! The second,
1777 # third, and fourth ones have a minus because they're
1778 # conjugated.
1779 Sij_real = _embed_quaternion_matrix(Eij + Eij.transpose())
1780 S.append(Sij_real)
1781 Sij_I = _embed_quaternion_matrix(I*Eij - I*Eij.transpose())
1782 S.append(Sij_I)
1783 Sij_J = _embed_quaternion_matrix(J*Eij - J*Eij.transpose())
1784 S.append(Sij_J)
1785 Sij_K = _embed_quaternion_matrix(K*Eij - K*Eij.transpose())
1786 S.append(Sij_K)
1787 return tuple(S)
1788
1789
1790 def _mat2vec(m):
1791 return vector(m.base_ring(), m.list())
1792
1793 def _vec2mat(v):
1794 return matrix(v.base_ring(), sqrt(v.degree()), v.list())
1795
1796 def _multiplication_table_from_matrix_basis(basis):
1797 """
1798 At least three of the five simple Euclidean Jordan algebras have the
1799 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
1800 multiplication on the right is matrix multiplication. Given a basis
1801 for the underlying matrix space, this function returns a
1802 multiplication table (obtained by looping through the basis
1803 elements) for an algebra of those matrices. A reordered copy
1804 of the basis is also returned to work around the fact that
1805 the ``span()`` in this function will change the order of the basis
1806 from what we think it is, to... something else.
1807 """
1808 # In S^2, for example, we nominally have four coordinates even
1809 # though the space is of dimension three only. The vector space V
1810 # is supposed to hold the entire long vector, and the subspace W
1811 # of V will be spanned by the vectors that arise from symmetric
1812 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
1813 field = basis[0].base_ring()
1814 dimension = basis[0].nrows()
1815
1816 V = VectorSpace(field, dimension**2)
1817 W = V.span( _mat2vec(s) for s in basis )
1818
1819 # Taking the span above reorders our basis (thanks, jerk!) so we
1820 # need to put our "matrix basis" in the same order as the
1821 # (reordered) vector basis.
1822 S = tuple( _vec2mat(b) for b in W.basis() )
1823
1824 Qs = []
1825 for s in S:
1826 # Brute force the multiplication-by-s matrix by looping
1827 # through all elements of the basis and doing the computation
1828 # to find out what the corresponding row should be. BEWARE:
1829 # these multiplication tables won't be symmetric! It therefore
1830 # becomes REALLY IMPORTANT that the underlying algebra
1831 # constructor uses ROW vectors and not COLUMN vectors. That's
1832 # why we're computing rows here and not columns.
1833 Q_rows = []
1834 for t in S:
1835 this_row = _mat2vec((s*t + t*s)/2)
1836 Q_rows.append(W.coordinates(this_row))
1837 Q = matrix(field, W.dimension(), Q_rows)
1838 Qs.append(Q)
1839
1840 return (Qs, S)
1841
1842
1843 def _embed_complex_matrix(M):
1844 """
1845 Embed the n-by-n complex matrix ``M`` into the space of real
1846 matrices of size 2n-by-2n via the map the sends each entry `z = a +
1847 bi` to the block matrix ``[[a,b],[-b,a]]``.
1848
1849 SETUP::
1850
1851 sage: from mjo.eja.eja_algebra import _embed_complex_matrix
1852
1853 EXAMPLES::
1854
1855 sage: F = QuadraticField(-1,'i')
1856 sage: x1 = F(4 - 2*i)
1857 sage: x2 = F(1 + 2*i)
1858 sage: x3 = F(-i)
1859 sage: x4 = F(6)
1860 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
1861 sage: _embed_complex_matrix(M)
1862 [ 4 -2| 1 2]
1863 [ 2 4|-2 1]
1864 [-----+-----]
1865 [ 0 -1| 6 0]
1866 [ 1 0| 0 6]
1867
1868 TESTS:
1869
1870 Embedding is a homomorphism (isomorphism, in fact)::
1871
1872 sage: set_random_seed()
1873 sage: n = ZZ.random_element(5)
1874 sage: F = QuadraticField(-1, 'i')
1875 sage: X = random_matrix(F, n)
1876 sage: Y = random_matrix(F, n)
1877 sage: actual = _embed_complex_matrix(X) * _embed_complex_matrix(Y)
1878 sage: expected = _embed_complex_matrix(X*Y)
1879 sage: actual == expected
1880 True
1881
1882 """
1883 n = M.nrows()
1884 if M.ncols() != n:
1885 raise ValueError("the matrix 'M' must be square")
1886 field = M.base_ring()
1887 blocks = []
1888 for z in M.list():
1889 a = z.real()
1890 b = z.imag()
1891 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1892
1893 # We can drop the imaginaries here.
1894 return matrix.block(field.base_ring(), n, blocks)
1895
1896
1897 def _unembed_complex_matrix(M):
1898 """
1899 The inverse of _embed_complex_matrix().
1900
1901 SETUP::
1902
1903 sage: from mjo.eja.eja_algebra import (_embed_complex_matrix,
1904 ....: _unembed_complex_matrix)
1905
1906 EXAMPLES::
1907
1908 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1909 ....: [-2, 1, -4, 3],
1910 ....: [ 9, 10, 11, 12],
1911 ....: [-10, 9, -12, 11] ])
1912 sage: _unembed_complex_matrix(A)
1913 [ 2*i + 1 4*i + 3]
1914 [ 10*i + 9 12*i + 11]
1915
1916 TESTS:
1917
1918 Unembedding is the inverse of embedding::
1919
1920 sage: set_random_seed()
1921 sage: F = QuadraticField(-1, 'i')
1922 sage: M = random_matrix(F, 3)
1923 sage: _unembed_complex_matrix(_embed_complex_matrix(M)) == M
1924 True
1925
1926 """
1927 n = ZZ(M.nrows())
1928 if M.ncols() != n:
1929 raise ValueError("the matrix 'M' must be square")
1930 if not n.mod(2).is_zero():
1931 raise ValueError("the matrix 'M' must be a complex embedding")
1932
1933 F = QuadraticField(-1, 'i')
1934 i = F.gen()
1935
1936 # Go top-left to bottom-right (reading order), converting every
1937 # 2-by-2 block we see to a single complex element.
1938 elements = []
1939 for k in xrange(n/2):
1940 for j in xrange(n/2):
1941 submat = M[2*k:2*k+2,2*j:2*j+2]
1942 if submat[0,0] != submat[1,1]:
1943 raise ValueError('bad on-diagonal submatrix')
1944 if submat[0,1] != -submat[1,0]:
1945 raise ValueError('bad off-diagonal submatrix')
1946 z = submat[0,0] + submat[0,1]*i
1947 elements.append(z)
1948
1949 return matrix(F, n/2, elements)
1950
1951
1952 def _embed_quaternion_matrix(M):
1953 """
1954 Embed the n-by-n quaternion matrix ``M`` into the space of real
1955 matrices of size 4n-by-4n by first sending each quaternion entry
1956 `z = a + bi + cj + dk` to the block-complex matrix
1957 ``[[a + bi, c+di],[-c + di, a-bi]]`, and then embedding those into
1958 a real matrix.
1959
1960 SETUP::
1961
1962 sage: from mjo.eja.eja_algebra import _embed_quaternion_matrix
1963
1964 EXAMPLES::
1965
1966 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1967 sage: i,j,k = Q.gens()
1968 sage: x = 1 + 2*i + 3*j + 4*k
1969 sage: M = matrix(Q, 1, [[x]])
1970 sage: _embed_quaternion_matrix(M)
1971 [ 1 2 3 4]
1972 [-2 1 -4 3]
1973 [-3 4 1 -2]
1974 [-4 -3 2 1]
1975
1976 Embedding is a homomorphism (isomorphism, in fact)::
1977
1978 sage: set_random_seed()
1979 sage: n = ZZ.random_element(5)
1980 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1981 sage: X = random_matrix(Q, n)
1982 sage: Y = random_matrix(Q, n)
1983 sage: actual = _embed_quaternion_matrix(X)*_embed_quaternion_matrix(Y)
1984 sage: expected = _embed_quaternion_matrix(X*Y)
1985 sage: actual == expected
1986 True
1987
1988 """
1989 quaternions = M.base_ring()
1990 n = M.nrows()
1991 if M.ncols() != n:
1992 raise ValueError("the matrix 'M' must be square")
1993
1994 F = QuadraticField(-1, 'i')
1995 i = F.gen()
1996
1997 blocks = []
1998 for z in M.list():
1999 t = z.coefficient_tuple()
2000 a = t[0]
2001 b = t[1]
2002 c = t[2]
2003 d = t[3]
2004 cplx_matrix = matrix(F, 2, [[ a + b*i, c + d*i],
2005 [-c + d*i, a - b*i]])
2006 blocks.append(_embed_complex_matrix(cplx_matrix))
2007
2008 # We should have real entries by now, so use the realest field
2009 # we've got for the return value.
2010 return matrix.block(quaternions.base_ring(), n, blocks)
2011
2012
2013 def _unembed_quaternion_matrix(M):
2014 """
2015 The inverse of _embed_quaternion_matrix().
2016
2017 SETUP::
2018
2019 sage: from mjo.eja.eja_algebra import (_embed_quaternion_matrix,
2020 ....: _unembed_quaternion_matrix)
2021
2022 EXAMPLES::
2023
2024 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
2025 ....: [-2, 1, -4, 3],
2026 ....: [-3, 4, 1, -2],
2027 ....: [-4, -3, 2, 1]])
2028 sage: _unembed_quaternion_matrix(M)
2029 [1 + 2*i + 3*j + 4*k]
2030
2031 TESTS:
2032
2033 Unembedding is the inverse of embedding::
2034
2035 sage: set_random_seed()
2036 sage: Q = QuaternionAlgebra(QQ, -1, -1)
2037 sage: M = random_matrix(Q, 3)
2038 sage: _unembed_quaternion_matrix(_embed_quaternion_matrix(M)) == M
2039 True
2040
2041 """
2042 n = ZZ(M.nrows())
2043 if M.ncols() != n:
2044 raise ValueError("the matrix 'M' must be square")
2045 if not n.mod(4).is_zero():
2046 raise ValueError("the matrix 'M' must be a complex embedding")
2047
2048 Q = QuaternionAlgebra(QQ,-1,-1)
2049 i,j,k = Q.gens()
2050
2051 # Go top-left to bottom-right (reading order), converting every
2052 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
2053 # quaternion block.
2054 elements = []
2055 for l in xrange(n/4):
2056 for m in xrange(n/4):
2057 submat = _unembed_complex_matrix(M[4*l:4*l+4,4*m:4*m+4])
2058 if submat[0,0] != submat[1,1].conjugate():
2059 raise ValueError('bad on-diagonal submatrix')
2060 if submat[0,1] != -submat[1,0].conjugate():
2061 raise ValueError('bad off-diagonal submatrix')
2062 z = submat[0,0].real() + submat[0,0].imag()*i
2063 z += submat[0,1].real()*j + submat[0,1].imag()*k
2064 elements.append(z)
2065
2066 return matrix(Q, n/4, elements)
2067
2068
2069 # The usual inner product on R^n.
2070 def _usual_ip(x,y):
2071 return x.vector().inner_product(y.vector())
2072
2073 # The inner product used for the real symmetric simple EJA.
2074 # We keep it as a separate function because e.g. the complex
2075 # algebra uses the same inner product, except divided by 2.
2076 def _matrix_ip(X,Y):
2077 X_mat = X.natural_representation()
2078 Y_mat = Y.natural_representation()
2079 return (X_mat*Y_mat).trace()
2080
2081
2082 class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra):
2083 """
2084 The rank-n simple EJA consisting of real symmetric n-by-n
2085 matrices, the usual symmetric Jordan product, and the trace inner
2086 product. It has dimension `(n^2 + n)/2` over the reals.
2087
2088 SETUP::
2089
2090 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
2091
2092 EXAMPLES::
2093
2094 sage: J = RealSymmetricEJA(2)
2095 sage: e0, e1, e2 = J.gens()
2096 sage: e0*e0
2097 e0
2098 sage: e1*e1
2099 e0 + e2
2100 sage: e2*e2
2101 e2
2102
2103 TESTS:
2104
2105 The degree of this algebra is `(n^2 + n) / 2`::
2106
2107 sage: set_random_seed()
2108 sage: n = ZZ.random_element(1,5)
2109 sage: J = RealSymmetricEJA(n)
2110 sage: J.degree() == (n^2 + n)/2
2111 True
2112
2113 The Jordan multiplication is what we think it is::
2114
2115 sage: set_random_seed()
2116 sage: n = ZZ.random_element(1,5)
2117 sage: J = RealSymmetricEJA(n)
2118 sage: x = J.random_element()
2119 sage: y = J.random_element()
2120 sage: actual = (x*y).natural_representation()
2121 sage: X = x.natural_representation()
2122 sage: Y = y.natural_representation()
2123 sage: expected = (X*Y + Y*X)/2
2124 sage: actual == expected
2125 True
2126 sage: J(expected) == x*y
2127 True
2128
2129 """
2130 @staticmethod
2131 def __classcall_private__(cls, n, field=QQ):
2132 S = _real_symmetric_basis(n, field=field)
2133 (Qs, T) = _multiplication_table_from_matrix_basis(S)
2134
2135 fdeja = super(RealSymmetricEJA, cls)
2136 return fdeja.__classcall_private__(cls,
2137 field,
2138 Qs,
2139 rank=n,
2140 natural_basis=T)
2141
2142 def inner_product(self, x, y):
2143 return _matrix_ip(x,y)
2144
2145
2146 class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
2147 """
2148 The rank-n simple EJA consisting of complex Hermitian n-by-n
2149 matrices over the real numbers, the usual symmetric Jordan product,
2150 and the real-part-of-trace inner product. It has dimension `n^2` over
2151 the reals.
2152
2153 SETUP::
2154
2155 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
2156
2157 TESTS:
2158
2159 The degree of this algebra is `n^2`::
2160
2161 sage: set_random_seed()
2162 sage: n = ZZ.random_element(1,5)
2163 sage: J = ComplexHermitianEJA(n)
2164 sage: J.degree() == n^2
2165 True
2166
2167 The Jordan multiplication is what we think it is::
2168
2169 sage: set_random_seed()
2170 sage: n = ZZ.random_element(1,5)
2171 sage: J = ComplexHermitianEJA(n)
2172 sage: x = J.random_element()
2173 sage: y = J.random_element()
2174 sage: actual = (x*y).natural_representation()
2175 sage: X = x.natural_representation()
2176 sage: Y = y.natural_representation()
2177 sage: expected = (X*Y + Y*X)/2
2178 sage: actual == expected
2179 True
2180 sage: J(expected) == x*y
2181 True
2182
2183 """
2184 @staticmethod
2185 def __classcall_private__(cls, n, field=QQ):
2186 S = _complex_hermitian_basis(n)
2187 (Qs, T) = _multiplication_table_from_matrix_basis(S)
2188
2189 fdeja = super(ComplexHermitianEJA, cls)
2190 return fdeja.__classcall_private__(cls,
2191 field,
2192 Qs,
2193 rank=n,
2194 natural_basis=T)
2195
2196 def inner_product(self, x, y):
2197 # Since a+bi on the diagonal is represented as
2198 #
2199 # a + bi = [ a b ]
2200 # [ -b a ],
2201 #
2202 # we'll double-count the "a" entries if we take the trace of
2203 # the embedding.
2204 return _matrix_ip(x,y)/2
2205
2206
2207 class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
2208 """
2209 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
2210 matrices, the usual symmetric Jordan product, and the
2211 real-part-of-trace inner product. It has dimension `2n^2 - n` over
2212 the reals.
2213
2214 SETUP::
2215
2216 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
2217
2218 TESTS:
2219
2220 The degree of this algebra is `n^2`::
2221
2222 sage: set_random_seed()
2223 sage: n = ZZ.random_element(1,5)
2224 sage: J = QuaternionHermitianEJA(n)
2225 sage: J.degree() == 2*(n^2) - n
2226 True
2227
2228 The Jordan multiplication is what we think it is::
2229
2230 sage: set_random_seed()
2231 sage: n = ZZ.random_element(1,5)
2232 sage: J = QuaternionHermitianEJA(n)
2233 sage: x = J.random_element()
2234 sage: y = J.random_element()
2235 sage: actual = (x*y).natural_representation()
2236 sage: X = x.natural_representation()
2237 sage: Y = y.natural_representation()
2238 sage: expected = (X*Y + Y*X)/2
2239 sage: actual == expected
2240 True
2241 sage: J(expected) == x*y
2242 True
2243
2244 """
2245 @staticmethod
2246 def __classcall_private__(cls, n, field=QQ):
2247 S = _quaternion_hermitian_basis(n)
2248 (Qs, T) = _multiplication_table_from_matrix_basis(S)
2249
2250 fdeja = super(QuaternionHermitianEJA, cls)
2251 return fdeja.__classcall_private__(cls,
2252 field,
2253 Qs,
2254 rank=n,
2255 natural_basis=T)
2256
2257 def inner_product(self, x, y):
2258 # Since a+bi+cj+dk on the diagonal is represented as
2259 #
2260 # a + bi +cj + dk = [ a b c d]
2261 # [ -b a -d c]
2262 # [ -c d a -b]
2263 # [ -d -c b a],
2264 #
2265 # we'll quadruple-count the "a" entries if we take the trace of
2266 # the embedding.
2267 return _matrix_ip(x,y)/4
2268
2269
2270 class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra):
2271 """
2272 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2273 with the usual inner product and jordan product ``x*y =
2274 (<x_bar,y_bar>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
2275 the reals.
2276
2277 SETUP::
2278
2279 sage: from mjo.eja.eja_algebra import JordanSpinEJA
2280
2281 EXAMPLES:
2282
2283 This multiplication table can be verified by hand::
2284
2285 sage: J = JordanSpinEJA(4)
2286 sage: e0,e1,e2,e3 = J.gens()
2287 sage: e0*e0
2288 e0
2289 sage: e0*e1
2290 e1
2291 sage: e0*e2
2292 e2
2293 sage: e0*e3
2294 e3
2295 sage: e1*e2
2296 0
2297 sage: e1*e3
2298 0
2299 sage: e2*e3
2300 0
2301
2302 """
2303 @staticmethod
2304 def __classcall_private__(cls, n, field=QQ):
2305 Qs = []
2306 id_matrix = matrix.identity(field, n)
2307 for i in xrange(n):
2308 ei = id_matrix.column(i)
2309 Qi = matrix.zero(field, n)
2310 Qi.set_row(0, ei)
2311 Qi.set_column(0, ei)
2312 Qi += matrix.diagonal(n, [ei[0]]*n)
2313 # The addition of the diagonal matrix adds an extra ei[0] in the
2314 # upper-left corner of the matrix.
2315 Qi[0,0] = Qi[0,0] * ~field(2)
2316 Qs.append(Qi)
2317
2318 # The rank of the spin algebra is two, unless we're in a
2319 # one-dimensional ambient space (because the rank is bounded by
2320 # the ambient dimension).
2321 fdeja = super(JordanSpinEJA, cls)
2322 return fdeja.__classcall_private__(cls, field, Qs, rank=min(n,2))
2323
2324 def inner_product(self, x, y):
2325 return _usual_ip(x,y)