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