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