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