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