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