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