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