]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: simplify and justify the implementation of is_nilpotent() for elements.
[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
1063 EXAMPLES:
1064
1065 The identity element always has degree one, but any element
1066 linearly-independent from it is regular::
1067
1068 sage: J = JordanSpinEJA(5)
1069 sage: J.one().is_regular()
1070 False
1071 sage: e0, e1, e2, e3, e4 = J.gens() # e0 is the identity
1072 sage: for x in J.gens():
1073 ....: (J.one() + x).is_regular()
1074 False
1075 True
1076 True
1077 True
1078 True
1079
1080 """
1081 return self.degree() == self.parent().rank()
1082
1083
1084 def degree(self):
1085 """
1086 Compute the degree of this element the straightforward way
1087 according to the definition; by appending powers to a list
1088 and figuring out its dimension (that is, whether or not
1089 they're linearly dependent).
1090
1091 SETUP::
1092
1093 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1094
1095 EXAMPLES::
1096
1097 sage: J = JordanSpinEJA(4)
1098 sage: J.one().degree()
1099 1
1100 sage: e0,e1,e2,e3 = J.gens()
1101 sage: (e0 - e1).degree()
1102 2
1103
1104 In the spin factor algebra (of rank two), all elements that
1105 aren't multiples of the identity are regular::
1106
1107 sage: set_random_seed()
1108 sage: n = ZZ.random_element(1,10)
1109 sage: J = JordanSpinEJA(n)
1110 sage: x = J.random_element()
1111 sage: x == x.coefficient(0)*J.one() or x.degree() == 2
1112 True
1113
1114 """
1115 return self.span_of_powers().dimension()
1116
1117
1118 def left_matrix(self):
1119 """
1120 Our parent class defines ``left_matrix`` and ``matrix``
1121 methods whose names are misleading. We don't want them.
1122 """
1123 raise NotImplementedError("use operator().matrix() instead")
1124
1125 matrix = left_matrix
1126
1127
1128 def minimal_polynomial(self):
1129 """
1130 Return the minimal polynomial of this element,
1131 as a function of the variable `t`.
1132
1133 ALGORITHM:
1134
1135 We restrict ourselves to the associative subalgebra
1136 generated by this element, and then return the minimal
1137 polynomial of this element's operator matrix (in that
1138 subalgebra). This works by Baes Proposition 2.3.16.
1139
1140 SETUP::
1141
1142 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
1143 ....: random_eja)
1144
1145 TESTS:
1146
1147 The minimal polynomial of the identity and zero elements are
1148 always the same::
1149
1150 sage: set_random_seed()
1151 sage: J = random_eja()
1152 sage: J.one().minimal_polynomial()
1153 t - 1
1154 sage: J.zero().minimal_polynomial()
1155 t
1156
1157 The degree of an element is (by one definition) the degree
1158 of its minimal polynomial::
1159
1160 sage: set_random_seed()
1161 sage: x = random_eja().random_element()
1162 sage: x.degree() == x.minimal_polynomial().degree()
1163 True
1164
1165 The minimal polynomial and the characteristic polynomial coincide
1166 and are known (see Alizadeh, Example 11.11) for all elements of
1167 the spin factor algebra that aren't scalar multiples of the
1168 identity::
1169
1170 sage: set_random_seed()
1171 sage: n = ZZ.random_element(2,10)
1172 sage: J = JordanSpinEJA(n)
1173 sage: y = J.random_element()
1174 sage: while y == y.coefficient(0)*J.one():
1175 ....: y = J.random_element()
1176 sage: y0 = y.vector()[0]
1177 sage: y_bar = y.vector()[1:]
1178 sage: actual = y.minimal_polynomial()
1179 sage: t = PolynomialRing(J.base_ring(),'t').gen(0)
1180 sage: expected = t^2 - 2*y0*t + (y0^2 - norm(y_bar)^2)
1181 sage: bool(actual == expected)
1182 True
1183
1184 The minimal polynomial should always kill its element::
1185
1186 sage: set_random_seed()
1187 sage: x = random_eja().random_element()
1188 sage: p = x.minimal_polynomial()
1189 sage: x.apply_univariate_polynomial(p)
1190 0
1191
1192 """
1193 V = self.span_of_powers()
1194 assoc_subalg = self.subalgebra_generated_by()
1195 # Mis-design warning: the basis used for span_of_powers()
1196 # and subalgebra_generated_by() must be the same, and in
1197 # the same order!
1198 elt = assoc_subalg(V.coordinates(self.vector()))
1199 return elt.operator().minimal_polynomial()
1200
1201
1202
1203 def natural_representation(self):
1204 """
1205 Return a more-natural representation of this element.
1206
1207 Every finite-dimensional Euclidean Jordan Algebra is a
1208 direct sum of five simple algebras, four of which comprise
1209 Hermitian matrices. This method returns the original
1210 "natural" representation of this element as a Hermitian
1211 matrix, if it has one. If not, you get the usual representation.
1212
1213 SETUP::
1214
1215 sage: from mjo.eja.eja_algebra import (ComplexHermitianEJA,
1216 ....: QuaternionHermitianEJA)
1217
1218 EXAMPLES::
1219
1220 sage: J = ComplexHermitianEJA(3)
1221 sage: J.one()
1222 e0 + e5 + e8
1223 sage: J.one().natural_representation()
1224 [1 0 0 0 0 0]
1225 [0 1 0 0 0 0]
1226 [0 0 1 0 0 0]
1227 [0 0 0 1 0 0]
1228 [0 0 0 0 1 0]
1229 [0 0 0 0 0 1]
1230
1231 ::
1232
1233 sage: J = QuaternionHermitianEJA(3)
1234 sage: J.one()
1235 e0 + e9 + e14
1236 sage: J.one().natural_representation()
1237 [1 0 0 0 0 0 0 0 0 0 0 0]
1238 [0 1 0 0 0 0 0 0 0 0 0 0]
1239 [0 0 1 0 0 0 0 0 0 0 0 0]
1240 [0 0 0 1 0 0 0 0 0 0 0 0]
1241 [0 0 0 0 1 0 0 0 0 0 0 0]
1242 [0 0 0 0 0 1 0 0 0 0 0 0]
1243 [0 0 0 0 0 0 1 0 0 0 0 0]
1244 [0 0 0 0 0 0 0 1 0 0 0 0]
1245 [0 0 0 0 0 0 0 0 1 0 0 0]
1246 [0 0 0 0 0 0 0 0 0 1 0 0]
1247 [0 0 0 0 0 0 0 0 0 0 1 0]
1248 [0 0 0 0 0 0 0 0 0 0 0 1]
1249
1250 """
1251 B = self.parent().natural_basis()
1252 W = B[0].matrix_space()
1253 return W.linear_combination(zip(self.vector(), B))
1254
1255
1256 def operator(self):
1257 """
1258 Return the left-multiplication-by-this-element
1259 operator on the ambient algebra.
1260
1261 SETUP::
1262
1263 sage: from mjo.eja.eja_algebra import random_eja
1264
1265 TESTS::
1266
1267 sage: set_random_seed()
1268 sage: J = random_eja()
1269 sage: x = J.random_element()
1270 sage: y = J.random_element()
1271 sage: x.operator()(y) == x*y
1272 True
1273 sage: y.operator()(x) == x*y
1274 True
1275
1276 """
1277 P = self.parent()
1278 fda_elt = FiniteDimensionalAlgebraElement(P, self)
1279 return FiniteDimensionalEuclideanJordanAlgebraOperator(
1280 P,
1281 P,
1282 fda_elt.matrix().transpose() )
1283
1284
1285 def quadratic_representation(self, other=None):
1286 """
1287 Return the quadratic representation of this element.
1288
1289 SETUP::
1290
1291 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
1292 ....: random_eja)
1293
1294 EXAMPLES:
1295
1296 The explicit form in the spin factor algebra is given by
1297 Alizadeh's Example 11.12::
1298
1299 sage: set_random_seed()
1300 sage: n = ZZ.random_element(1,10)
1301 sage: J = JordanSpinEJA(n)
1302 sage: x = J.random_element()
1303 sage: x_vec = x.vector()
1304 sage: x0 = x_vec[0]
1305 sage: x_bar = x_vec[1:]
1306 sage: A = matrix(QQ, 1, [x_vec.inner_product(x_vec)])
1307 sage: B = 2*x0*x_bar.row()
1308 sage: C = 2*x0*x_bar.column()
1309 sage: D = matrix.identity(QQ, n-1)
1310 sage: D = (x0^2 - x_bar.inner_product(x_bar))*D
1311 sage: D = D + 2*x_bar.tensor_product(x_bar)
1312 sage: Q = matrix.block(2,2,[A,B,C,D])
1313 sage: Q == x.quadratic_representation().matrix()
1314 True
1315
1316 Test all of the properties from Theorem 11.2 in Alizadeh::
1317
1318 sage: set_random_seed()
1319 sage: J = random_eja()
1320 sage: x = J.random_element()
1321 sage: y = J.random_element()
1322 sage: Lx = x.operator()
1323 sage: Lxx = (x*x).operator()
1324 sage: Qx = x.quadratic_representation()
1325 sage: Qy = y.quadratic_representation()
1326 sage: Qxy = x.quadratic_representation(y)
1327 sage: Qex = J.one().quadratic_representation(x)
1328 sage: n = ZZ.random_element(10)
1329 sage: Qxn = (x^n).quadratic_representation()
1330
1331 Property 1:
1332
1333 sage: 2*Qxy == (x+y).quadratic_representation() - Qx - Qy
1334 True
1335
1336 Property 2 (multiply on the right for :trac:`28272`):
1337
1338 sage: alpha = QQ.random_element()
1339 sage: (alpha*x).quadratic_representation() == Qx*(alpha^2)
1340 True
1341
1342 Property 3:
1343
1344 sage: not x.is_invertible() or ( Qx(x.inverse()) == x )
1345 True
1346
1347 sage: not x.is_invertible() or (
1348 ....: ~Qx
1349 ....: ==
1350 ....: x.inverse().quadratic_representation() )
1351 True
1352
1353 sage: Qxy(J.one()) == x*y
1354 True
1355
1356 Property 4:
1357
1358 sage: not x.is_invertible() or (
1359 ....: x.quadratic_representation(x.inverse())*Qx
1360 ....: == Qx*x.quadratic_representation(x.inverse()) )
1361 True
1362
1363 sage: not x.is_invertible() or (
1364 ....: x.quadratic_representation(x.inverse())*Qx
1365 ....: ==
1366 ....: 2*x.operator()*Qex - Qx )
1367 True
1368
1369 sage: 2*x.operator()*Qex - Qx == Lxx
1370 True
1371
1372 Property 5:
1373
1374 sage: Qy(x).quadratic_representation() == Qy*Qx*Qy
1375 True
1376
1377 Property 6:
1378
1379 sage: Qxn == (Qx)^n
1380 True
1381
1382 Property 7:
1383
1384 sage: not x.is_invertible() or (
1385 ....: Qx*x.inverse().operator() == Lx )
1386 True
1387
1388 Property 8:
1389
1390 sage: not x.operator_commutes_with(y) or (
1391 ....: Qx(y)^n == Qxn(y^n) )
1392 True
1393
1394 """
1395 if other is None:
1396 other=self
1397 elif not other in self.parent():
1398 raise TypeError("'other' must live in the same algebra")
1399
1400 L = self.operator()
1401 M = other.operator()
1402 return ( L*M + M*L - (self*other).operator() )
1403
1404
1405 def span_of_powers(self):
1406 """
1407 Return the vector space spanned by successive powers of
1408 this element.
1409 """
1410 # The dimension of the subalgebra can't be greater than
1411 # the big algebra, so just put everything into a list
1412 # and let span() get rid of the excess.
1413 #
1414 # We do the extra ambient_vector_space() in case we're messing
1415 # with polynomials and the direct parent is a module.
1416 V = self.parent().vector_space()
1417 return V.span( (self**d).vector() for d in xrange(V.dimension()) )
1418
1419
1420 def subalgebra_generated_by(self):
1421 """
1422 Return the associative subalgebra of the parent EJA generated
1423 by this element.
1424
1425 SETUP::
1426
1427 sage: from mjo.eja.eja_algebra import random_eja
1428
1429 TESTS::
1430
1431 sage: set_random_seed()
1432 sage: x = random_eja().random_element()
1433 sage: x.subalgebra_generated_by().is_associative()
1434 True
1435
1436 Squaring in the subalgebra should work the same as in
1437 the superalgebra::
1438
1439 sage: set_random_seed()
1440 sage: x = random_eja().random_element()
1441 sage: u = x.subalgebra_generated_by().random_element()
1442 sage: u.operator()(u) == u^2
1443 True
1444
1445 """
1446 # First get the subspace spanned by the powers of myself...
1447 V = self.span_of_powers()
1448 F = self.base_ring()
1449
1450 # Now figure out the entries of the right-multiplication
1451 # matrix for the successive basis elements b0, b1,... of
1452 # that subspace.
1453 mats = []
1454 for b_right in V.basis():
1455 eja_b_right = self.parent()(b_right)
1456 b_right_rows = []
1457 # The first row of the right-multiplication matrix by
1458 # b1 is what we get if we apply that matrix to b1. The
1459 # second row of the right multiplication matrix by b1
1460 # is what we get when we apply that matrix to b2...
1461 #
1462 # IMPORTANT: this assumes that all vectors are COLUMN
1463 # vectors, unlike our superclass (which uses row vectors).
1464 for b_left in V.basis():
1465 eja_b_left = self.parent()(b_left)
1466 # Multiply in the original EJA, but then get the
1467 # coordinates from the subalgebra in terms of its
1468 # basis.
1469 this_row = V.coordinates((eja_b_left*eja_b_right).vector())
1470 b_right_rows.append(this_row)
1471 b_right_matrix = matrix(F, b_right_rows)
1472 mats.append(b_right_matrix)
1473
1474 # It's an algebra of polynomials in one element, and EJAs
1475 # are power-associative.
1476 #
1477 # TODO: choose generator names intelligently.
1478 #
1479 # The rank is the highest possible degree of a minimal polynomial,
1480 # and is bounded above by the dimension. We know in this case that
1481 # there's an element whose minimal polynomial has the same degree
1482 # as the space's dimension, so that must be its rank too.
1483 return FiniteDimensionalEuclideanJordanAlgebra(
1484 F,
1485 mats,
1486 V.dimension(),
1487 assume_associative=True,
1488 names='f')
1489
1490
1491 def subalgebra_idempotent(self):
1492 """
1493 Find an idempotent in the associative subalgebra I generate
1494 using Proposition 2.3.5 in Baes.
1495
1496 SETUP::
1497
1498 sage: from mjo.eja.eja_algebra import random_eja
1499
1500 TESTS::
1501
1502 sage: set_random_seed()
1503 sage: J = random_eja()
1504 sage: x = J.random_element()
1505 sage: while x.is_nilpotent():
1506 ....: x = J.random_element()
1507 sage: c = x.subalgebra_idempotent()
1508 sage: c^2 == c
1509 True
1510
1511 """
1512 if self.is_nilpotent():
1513 raise ValueError("this only works with non-nilpotent elements!")
1514
1515 V = self.span_of_powers()
1516 J = self.subalgebra_generated_by()
1517 # Mis-design warning: the basis used for span_of_powers()
1518 # and subalgebra_generated_by() must be the same, and in
1519 # the same order!
1520 u = J(V.coordinates(self.vector()))
1521
1522 # The image of the matrix of left-u^m-multiplication
1523 # will be minimal for some natural number s...
1524 s = 0
1525 minimal_dim = V.dimension()
1526 for i in xrange(1, V.dimension()):
1527 this_dim = (u**i).operator().matrix().image().dimension()
1528 if this_dim < minimal_dim:
1529 minimal_dim = this_dim
1530 s = i
1531
1532 # Now minimal_matrix should correspond to the smallest
1533 # non-zero subspace in Baes's (or really, Koecher's)
1534 # proposition.
1535 #
1536 # However, we need to restrict the matrix to work on the
1537 # subspace... or do we? Can't we just solve, knowing that
1538 # A(c) = u^(s+1) should have a solution in the big space,
1539 # too?
1540 #
1541 # Beware, solve_right() means that we're using COLUMN vectors.
1542 # Our FiniteDimensionalAlgebraElement superclass uses rows.
1543 u_next = u**(s+1)
1544 A = u_next.operator().matrix()
1545 c_coordinates = A.solve_right(u_next.vector())
1546
1547 # Now c_coordinates is the idempotent we want, but it's in
1548 # the coordinate system of the subalgebra.
1549 #
1550 # We need the basis for J, but as elements of the parent algebra.
1551 #
1552 basis = [self.parent(v) for v in V.basis()]
1553 return self.parent().linear_combination(zip(c_coordinates, basis))
1554
1555
1556 def trace(self):
1557 """
1558 Return my trace, the sum of my eigenvalues.
1559
1560 SETUP::
1561
1562 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
1563 ....: RealCartesianProductEJA,
1564 ....: random_eja)
1565
1566 EXAMPLES::
1567
1568 sage: J = JordanSpinEJA(3)
1569 sage: x = sum(J.gens())
1570 sage: x.trace()
1571 2
1572
1573 ::
1574
1575 sage: J = RealCartesianProductEJA(5)
1576 sage: J.one().trace()
1577 5
1578
1579 TESTS:
1580
1581 The trace of an element is a real number::
1582
1583 sage: set_random_seed()
1584 sage: J = random_eja()
1585 sage: J.random_element().trace() in J.base_ring()
1586 True
1587
1588 """
1589 P = self.parent()
1590 r = P.rank()
1591 p = P._charpoly_coeff(r-1)
1592 # The _charpoly_coeff function already adds the factor of
1593 # -1 to ensure that _charpoly_coeff(r-1) is really what
1594 # appears in front of t^{r-1} in the charpoly. However,
1595 # we want the negative of THAT for the trace.
1596 return -p(*self.vector())
1597
1598
1599 def trace_inner_product(self, other):
1600 """
1601 Return the trace inner product of myself and ``other``.
1602
1603 SETUP::
1604
1605 sage: from mjo.eja.eja_algebra import random_eja
1606
1607 TESTS:
1608
1609 The trace inner product is commutative::
1610
1611 sage: set_random_seed()
1612 sage: J = random_eja()
1613 sage: x = J.random_element(); y = J.random_element()
1614 sage: x.trace_inner_product(y) == y.trace_inner_product(x)
1615 True
1616
1617 The trace inner product is bilinear::
1618
1619 sage: set_random_seed()
1620 sage: J = random_eja()
1621 sage: x = J.random_element()
1622 sage: y = J.random_element()
1623 sage: z = J.random_element()
1624 sage: a = QQ.random_element();
1625 sage: actual = (a*(x+z)).trace_inner_product(y)
1626 sage: expected = ( a*x.trace_inner_product(y) +
1627 ....: a*z.trace_inner_product(y) )
1628 sage: actual == expected
1629 True
1630 sage: actual = x.trace_inner_product(a*(y+z))
1631 sage: expected = ( a*x.trace_inner_product(y) +
1632 ....: a*x.trace_inner_product(z) )
1633 sage: actual == expected
1634 True
1635
1636 The trace inner product satisfies the compatibility
1637 condition in the definition of a Euclidean Jordan algebra::
1638
1639 sage: set_random_seed()
1640 sage: J = random_eja()
1641 sage: x = J.random_element()
1642 sage: y = J.random_element()
1643 sage: z = J.random_element()
1644 sage: (x*y).trace_inner_product(z) == y.trace_inner_product(x*z)
1645 True
1646
1647 """
1648 if not other in self.parent():
1649 raise TypeError("'other' must live in the same algebra")
1650
1651 return (self*other).trace()
1652
1653
1654 class RealCartesianProductEJA(FiniteDimensionalEuclideanJordanAlgebra):
1655 """
1656 Return the Euclidean Jordan Algebra corresponding to the set
1657 `R^n` under the Hadamard product.
1658
1659 Note: this is nothing more than the Cartesian product of ``n``
1660 copies of the spin algebra. Once Cartesian product algebras
1661 are implemented, this can go.
1662
1663 SETUP::
1664
1665 sage: from mjo.eja.eja_algebra import RealCartesianProductEJA
1666
1667 EXAMPLES:
1668
1669 This multiplication table can be verified by hand::
1670
1671 sage: J = RealCartesianProductEJA(3)
1672 sage: e0,e1,e2 = J.gens()
1673 sage: e0*e0
1674 e0
1675 sage: e0*e1
1676 0
1677 sage: e0*e2
1678 0
1679 sage: e1*e1
1680 e1
1681 sage: e1*e2
1682 0
1683 sage: e2*e2
1684 e2
1685
1686 """
1687 @staticmethod
1688 def __classcall_private__(cls, n, field=QQ):
1689 # The FiniteDimensionalAlgebra constructor takes a list of
1690 # matrices, the ith representing right multiplication by the ith
1691 # basis element in the vector space. So if e_1 = (1,0,0), then
1692 # right (Hadamard) multiplication of x by e_1 picks out the first
1693 # component of x; and likewise for the ith basis element e_i.
1694 Qs = [ matrix(field, n, n, lambda k,j: 1*(k == j == i))
1695 for i in xrange(n) ]
1696
1697 fdeja = super(RealCartesianProductEJA, cls)
1698 return fdeja.__classcall_private__(cls, field, Qs, rank=n)
1699
1700 def inner_product(self, x, y):
1701 return _usual_ip(x,y)
1702
1703
1704 def random_eja():
1705 """
1706 Return a "random" finite-dimensional Euclidean Jordan Algebra.
1707
1708 ALGORITHM:
1709
1710 For now, we choose a random natural number ``n`` (greater than zero)
1711 and then give you back one of the following:
1712
1713 * The cartesian product of the rational numbers ``n`` times; this is
1714 ``QQ^n`` with the Hadamard product.
1715
1716 * The Jordan spin algebra on ``QQ^n``.
1717
1718 * The ``n``-by-``n`` rational symmetric matrices with the symmetric
1719 product.
1720
1721 * The ``n``-by-``n`` complex-rational Hermitian matrices embedded
1722 in the space of ``2n``-by-``2n`` real symmetric matrices.
1723
1724 * The ``n``-by-``n`` quaternion-rational Hermitian matrices embedded
1725 in the space of ``4n``-by-``4n`` real symmetric matrices.
1726
1727 Later this might be extended to return Cartesian products of the
1728 EJAs above.
1729
1730 SETUP::
1731
1732 sage: from mjo.eja.eja_algebra import random_eja
1733
1734 TESTS::
1735
1736 sage: random_eja()
1737 Euclidean Jordan algebra of degree...
1738
1739 """
1740
1741 # The max_n component lets us choose different upper bounds on the
1742 # value "n" that gets passed to the constructor. This is needed
1743 # because e.g. R^{10} is reasonable to test, while the Hermitian
1744 # 10-by-10 quaternion matrices are not.
1745 (constructor, max_n) = choice([(RealCartesianProductEJA, 6),
1746 (JordanSpinEJA, 6),
1747 (RealSymmetricEJA, 5),
1748 (ComplexHermitianEJA, 4),
1749 (QuaternionHermitianEJA, 3)])
1750 n = ZZ.random_element(1, max_n)
1751 return constructor(n, field=QQ)
1752
1753
1754
1755 def _real_symmetric_basis(n, field=QQ):
1756 """
1757 Return a basis for the space of real symmetric n-by-n matrices.
1758 """
1759 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1760 # coordinates.
1761 S = []
1762 for i in xrange(n):
1763 for j in xrange(i+1):
1764 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1765 if i == j:
1766 Sij = Eij
1767 else:
1768 # Beware, orthogonal but not normalized!
1769 Sij = Eij + Eij.transpose()
1770 S.append(Sij)
1771 return tuple(S)
1772
1773
1774 def _complex_hermitian_basis(n, field=QQ):
1775 """
1776 Returns a basis for the space of complex Hermitian n-by-n matrices.
1777
1778 SETUP::
1779
1780 sage: from mjo.eja.eja_algebra import _complex_hermitian_basis
1781
1782 TESTS::
1783
1784 sage: set_random_seed()
1785 sage: n = ZZ.random_element(1,5)
1786 sage: all( M.is_symmetric() for M in _complex_hermitian_basis(n) )
1787 True
1788
1789 """
1790 F = QuadraticField(-1, 'I')
1791 I = F.gen()
1792
1793 # This is like the symmetric case, but we need to be careful:
1794 #
1795 # * We want conjugate-symmetry, not just symmetry.
1796 # * The diagonal will (as a result) be real.
1797 #
1798 S = []
1799 for i in xrange(n):
1800 for j in xrange(i+1):
1801 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1802 if i == j:
1803 Sij = _embed_complex_matrix(Eij)
1804 S.append(Sij)
1805 else:
1806 # Beware, orthogonal but not normalized! The second one
1807 # has a minus because it's conjugated.
1808 Sij_real = _embed_complex_matrix(Eij + Eij.transpose())
1809 S.append(Sij_real)
1810 Sij_imag = _embed_complex_matrix(I*Eij - I*Eij.transpose())
1811 S.append(Sij_imag)
1812 return tuple(S)
1813
1814
1815 def _quaternion_hermitian_basis(n, field=QQ):
1816 """
1817 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
1818
1819 SETUP::
1820
1821 sage: from mjo.eja.eja_algebra import _quaternion_hermitian_basis
1822
1823 TESTS::
1824
1825 sage: set_random_seed()
1826 sage: n = ZZ.random_element(1,5)
1827 sage: all( M.is_symmetric() for M in _quaternion_hermitian_basis(n) )
1828 True
1829
1830 """
1831 Q = QuaternionAlgebra(QQ,-1,-1)
1832 I,J,K = Q.gens()
1833
1834 # This is like the symmetric case, but we need to be careful:
1835 #
1836 # * We want conjugate-symmetry, not just symmetry.
1837 # * The diagonal will (as a result) be real.
1838 #
1839 S = []
1840 for i in xrange(n):
1841 for j in xrange(i+1):
1842 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
1843 if i == j:
1844 Sij = _embed_quaternion_matrix(Eij)
1845 S.append(Sij)
1846 else:
1847 # Beware, orthogonal but not normalized! The second,
1848 # third, and fourth ones have a minus because they're
1849 # conjugated.
1850 Sij_real = _embed_quaternion_matrix(Eij + Eij.transpose())
1851 S.append(Sij_real)
1852 Sij_I = _embed_quaternion_matrix(I*Eij - I*Eij.transpose())
1853 S.append(Sij_I)
1854 Sij_J = _embed_quaternion_matrix(J*Eij - J*Eij.transpose())
1855 S.append(Sij_J)
1856 Sij_K = _embed_quaternion_matrix(K*Eij - K*Eij.transpose())
1857 S.append(Sij_K)
1858 return tuple(S)
1859
1860
1861 def _mat2vec(m):
1862 return vector(m.base_ring(), m.list())
1863
1864 def _vec2mat(v):
1865 return matrix(v.base_ring(), sqrt(v.degree()), v.list())
1866
1867 def _multiplication_table_from_matrix_basis(basis):
1868 """
1869 At least three of the five simple Euclidean Jordan algebras have the
1870 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
1871 multiplication on the right is matrix multiplication. Given a basis
1872 for the underlying matrix space, this function returns a
1873 multiplication table (obtained by looping through the basis
1874 elements) for an algebra of those matrices. A reordered copy
1875 of the basis is also returned to work around the fact that
1876 the ``span()`` in this function will change the order of the basis
1877 from what we think it is, to... something else.
1878 """
1879 # In S^2, for example, we nominally have four coordinates even
1880 # though the space is of dimension three only. The vector space V
1881 # is supposed to hold the entire long vector, and the subspace W
1882 # of V will be spanned by the vectors that arise from symmetric
1883 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
1884 field = basis[0].base_ring()
1885 dimension = basis[0].nrows()
1886
1887 V = VectorSpace(field, dimension**2)
1888 W = V.span( _mat2vec(s) for s in basis )
1889
1890 # Taking the span above reorders our basis (thanks, jerk!) so we
1891 # need to put our "matrix basis" in the same order as the
1892 # (reordered) vector basis.
1893 S = tuple( _vec2mat(b) for b in W.basis() )
1894
1895 Qs = []
1896 for s in S:
1897 # Brute force the multiplication-by-s matrix by looping
1898 # through all elements of the basis and doing the computation
1899 # to find out what the corresponding row should be. BEWARE:
1900 # these multiplication tables won't be symmetric! It therefore
1901 # becomes REALLY IMPORTANT that the underlying algebra
1902 # constructor uses ROW vectors and not COLUMN vectors. That's
1903 # why we're computing rows here and not columns.
1904 Q_rows = []
1905 for t in S:
1906 this_row = _mat2vec((s*t + t*s)/2)
1907 Q_rows.append(W.coordinates(this_row))
1908 Q = matrix(field, W.dimension(), Q_rows)
1909 Qs.append(Q)
1910
1911 return (Qs, S)
1912
1913
1914 def _embed_complex_matrix(M):
1915 """
1916 Embed the n-by-n complex matrix ``M`` into the space of real
1917 matrices of size 2n-by-2n via the map the sends each entry `z = a +
1918 bi` to the block matrix ``[[a,b],[-b,a]]``.
1919
1920 SETUP::
1921
1922 sage: from mjo.eja.eja_algebra import _embed_complex_matrix
1923
1924 EXAMPLES::
1925
1926 sage: F = QuadraticField(-1,'i')
1927 sage: x1 = F(4 - 2*i)
1928 sage: x2 = F(1 + 2*i)
1929 sage: x3 = F(-i)
1930 sage: x4 = F(6)
1931 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
1932 sage: _embed_complex_matrix(M)
1933 [ 4 -2| 1 2]
1934 [ 2 4|-2 1]
1935 [-----+-----]
1936 [ 0 -1| 6 0]
1937 [ 1 0| 0 6]
1938
1939 TESTS:
1940
1941 Embedding is a homomorphism (isomorphism, in fact)::
1942
1943 sage: set_random_seed()
1944 sage: n = ZZ.random_element(5)
1945 sage: F = QuadraticField(-1, 'i')
1946 sage: X = random_matrix(F, n)
1947 sage: Y = random_matrix(F, n)
1948 sage: actual = _embed_complex_matrix(X) * _embed_complex_matrix(Y)
1949 sage: expected = _embed_complex_matrix(X*Y)
1950 sage: actual == expected
1951 True
1952
1953 """
1954 n = M.nrows()
1955 if M.ncols() != n:
1956 raise ValueError("the matrix 'M' must be square")
1957 field = M.base_ring()
1958 blocks = []
1959 for z in M.list():
1960 a = z.real()
1961 b = z.imag()
1962 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1963
1964 # We can drop the imaginaries here.
1965 return matrix.block(field.base_ring(), n, blocks)
1966
1967
1968 def _unembed_complex_matrix(M):
1969 """
1970 The inverse of _embed_complex_matrix().
1971
1972 SETUP::
1973
1974 sage: from mjo.eja.eja_algebra import (_embed_complex_matrix,
1975 ....: _unembed_complex_matrix)
1976
1977 EXAMPLES::
1978
1979 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1980 ....: [-2, 1, -4, 3],
1981 ....: [ 9, 10, 11, 12],
1982 ....: [-10, 9, -12, 11] ])
1983 sage: _unembed_complex_matrix(A)
1984 [ 2*i + 1 4*i + 3]
1985 [ 10*i + 9 12*i + 11]
1986
1987 TESTS:
1988
1989 Unembedding is the inverse of embedding::
1990
1991 sage: set_random_seed()
1992 sage: F = QuadraticField(-1, 'i')
1993 sage: M = random_matrix(F, 3)
1994 sage: _unembed_complex_matrix(_embed_complex_matrix(M)) == M
1995 True
1996
1997 """
1998 n = ZZ(M.nrows())
1999 if M.ncols() != n:
2000 raise ValueError("the matrix 'M' must be square")
2001 if not n.mod(2).is_zero():
2002 raise ValueError("the matrix 'M' must be a complex embedding")
2003
2004 F = QuadraticField(-1, 'i')
2005 i = F.gen()
2006
2007 # Go top-left to bottom-right (reading order), converting every
2008 # 2-by-2 block we see to a single complex element.
2009 elements = []
2010 for k in xrange(n/2):
2011 for j in xrange(n/2):
2012 submat = M[2*k:2*k+2,2*j:2*j+2]
2013 if submat[0,0] != submat[1,1]:
2014 raise ValueError('bad on-diagonal submatrix')
2015 if submat[0,1] != -submat[1,0]:
2016 raise ValueError('bad off-diagonal submatrix')
2017 z = submat[0,0] + submat[0,1]*i
2018 elements.append(z)
2019
2020 return matrix(F, n/2, elements)
2021
2022
2023 def _embed_quaternion_matrix(M):
2024 """
2025 Embed the n-by-n quaternion matrix ``M`` into the space of real
2026 matrices of size 4n-by-4n by first sending each quaternion entry
2027 `z = a + bi + cj + dk` to the block-complex matrix
2028 ``[[a + bi, c+di],[-c + di, a-bi]]`, and then embedding those into
2029 a real matrix.
2030
2031 SETUP::
2032
2033 sage: from mjo.eja.eja_algebra import _embed_quaternion_matrix
2034
2035 EXAMPLES::
2036
2037 sage: Q = QuaternionAlgebra(QQ,-1,-1)
2038 sage: i,j,k = Q.gens()
2039 sage: x = 1 + 2*i + 3*j + 4*k
2040 sage: M = matrix(Q, 1, [[x]])
2041 sage: _embed_quaternion_matrix(M)
2042 [ 1 2 3 4]
2043 [-2 1 -4 3]
2044 [-3 4 1 -2]
2045 [-4 -3 2 1]
2046
2047 Embedding is a homomorphism (isomorphism, in fact)::
2048
2049 sage: set_random_seed()
2050 sage: n = ZZ.random_element(5)
2051 sage: Q = QuaternionAlgebra(QQ,-1,-1)
2052 sage: X = random_matrix(Q, n)
2053 sage: Y = random_matrix(Q, n)
2054 sage: actual = _embed_quaternion_matrix(X)*_embed_quaternion_matrix(Y)
2055 sage: expected = _embed_quaternion_matrix(X*Y)
2056 sage: actual == expected
2057 True
2058
2059 """
2060 quaternions = M.base_ring()
2061 n = M.nrows()
2062 if M.ncols() != n:
2063 raise ValueError("the matrix 'M' must be square")
2064
2065 F = QuadraticField(-1, 'i')
2066 i = F.gen()
2067
2068 blocks = []
2069 for z in M.list():
2070 t = z.coefficient_tuple()
2071 a = t[0]
2072 b = t[1]
2073 c = t[2]
2074 d = t[3]
2075 cplx_matrix = matrix(F, 2, [[ a + b*i, c + d*i],
2076 [-c + d*i, a - b*i]])
2077 blocks.append(_embed_complex_matrix(cplx_matrix))
2078
2079 # We should have real entries by now, so use the realest field
2080 # we've got for the return value.
2081 return matrix.block(quaternions.base_ring(), n, blocks)
2082
2083
2084 def _unembed_quaternion_matrix(M):
2085 """
2086 The inverse of _embed_quaternion_matrix().
2087
2088 SETUP::
2089
2090 sage: from mjo.eja.eja_algebra import (_embed_quaternion_matrix,
2091 ....: _unembed_quaternion_matrix)
2092
2093 EXAMPLES::
2094
2095 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
2096 ....: [-2, 1, -4, 3],
2097 ....: [-3, 4, 1, -2],
2098 ....: [-4, -3, 2, 1]])
2099 sage: _unembed_quaternion_matrix(M)
2100 [1 + 2*i + 3*j + 4*k]
2101
2102 TESTS:
2103
2104 Unembedding is the inverse of embedding::
2105
2106 sage: set_random_seed()
2107 sage: Q = QuaternionAlgebra(QQ, -1, -1)
2108 sage: M = random_matrix(Q, 3)
2109 sage: _unembed_quaternion_matrix(_embed_quaternion_matrix(M)) == M
2110 True
2111
2112 """
2113 n = ZZ(M.nrows())
2114 if M.ncols() != n:
2115 raise ValueError("the matrix 'M' must be square")
2116 if not n.mod(4).is_zero():
2117 raise ValueError("the matrix 'M' must be a complex embedding")
2118
2119 Q = QuaternionAlgebra(QQ,-1,-1)
2120 i,j,k = Q.gens()
2121
2122 # Go top-left to bottom-right (reading order), converting every
2123 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
2124 # quaternion block.
2125 elements = []
2126 for l in xrange(n/4):
2127 for m in xrange(n/4):
2128 submat = _unembed_complex_matrix(M[4*l:4*l+4,4*m:4*m+4])
2129 if submat[0,0] != submat[1,1].conjugate():
2130 raise ValueError('bad on-diagonal submatrix')
2131 if submat[0,1] != -submat[1,0].conjugate():
2132 raise ValueError('bad off-diagonal submatrix')
2133 z = submat[0,0].real() + submat[0,0].imag()*i
2134 z += submat[0,1].real()*j + submat[0,1].imag()*k
2135 elements.append(z)
2136
2137 return matrix(Q, n/4, elements)
2138
2139
2140 # The usual inner product on R^n.
2141 def _usual_ip(x,y):
2142 return x.vector().inner_product(y.vector())
2143
2144 # The inner product used for the real symmetric simple EJA.
2145 # We keep it as a separate function because e.g. the complex
2146 # algebra uses the same inner product, except divided by 2.
2147 def _matrix_ip(X,Y):
2148 X_mat = X.natural_representation()
2149 Y_mat = Y.natural_representation()
2150 return (X_mat*Y_mat).trace()
2151
2152
2153 class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra):
2154 """
2155 The rank-n simple EJA consisting of real symmetric n-by-n
2156 matrices, the usual symmetric Jordan product, and the trace inner
2157 product. It has dimension `(n^2 + n)/2` over the reals.
2158
2159 SETUP::
2160
2161 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
2162
2163 EXAMPLES::
2164
2165 sage: J = RealSymmetricEJA(2)
2166 sage: e0, e1, e2 = J.gens()
2167 sage: e0*e0
2168 e0
2169 sage: e1*e1
2170 e0 + e2
2171 sage: e2*e2
2172 e2
2173
2174 TESTS:
2175
2176 The degree of this algebra is `(n^2 + n) / 2`::
2177
2178 sage: set_random_seed()
2179 sage: n = ZZ.random_element(1,5)
2180 sage: J = RealSymmetricEJA(n)
2181 sage: J.degree() == (n^2 + n)/2
2182 True
2183
2184 The Jordan multiplication is what we think it is::
2185
2186 sage: set_random_seed()
2187 sage: n = ZZ.random_element(1,5)
2188 sage: J = RealSymmetricEJA(n)
2189 sage: x = J.random_element()
2190 sage: y = J.random_element()
2191 sage: actual = (x*y).natural_representation()
2192 sage: X = x.natural_representation()
2193 sage: Y = y.natural_representation()
2194 sage: expected = (X*Y + Y*X)/2
2195 sage: actual == expected
2196 True
2197 sage: J(expected) == x*y
2198 True
2199
2200 """
2201 @staticmethod
2202 def __classcall_private__(cls, n, field=QQ):
2203 S = _real_symmetric_basis(n, field=field)
2204 (Qs, T) = _multiplication_table_from_matrix_basis(S)
2205
2206 fdeja = super(RealSymmetricEJA, cls)
2207 return fdeja.__classcall_private__(cls,
2208 field,
2209 Qs,
2210 rank=n,
2211 natural_basis=T)
2212
2213 def inner_product(self, x, y):
2214 return _matrix_ip(x,y)
2215
2216
2217 class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
2218 """
2219 The rank-n simple EJA consisting of complex Hermitian n-by-n
2220 matrices over the real numbers, the usual symmetric Jordan product,
2221 and the real-part-of-trace inner product. It has dimension `n^2` over
2222 the reals.
2223
2224 SETUP::
2225
2226 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
2227
2228 TESTS:
2229
2230 The degree of this algebra is `n^2`::
2231
2232 sage: set_random_seed()
2233 sage: n = ZZ.random_element(1,5)
2234 sage: J = ComplexHermitianEJA(n)
2235 sage: J.degree() == n^2
2236 True
2237
2238 The Jordan multiplication is what we think it is::
2239
2240 sage: set_random_seed()
2241 sage: n = ZZ.random_element(1,5)
2242 sage: J = ComplexHermitianEJA(n)
2243 sage: x = J.random_element()
2244 sage: y = J.random_element()
2245 sage: actual = (x*y).natural_representation()
2246 sage: X = x.natural_representation()
2247 sage: Y = y.natural_representation()
2248 sage: expected = (X*Y + Y*X)/2
2249 sage: actual == expected
2250 True
2251 sage: J(expected) == x*y
2252 True
2253
2254 """
2255 @staticmethod
2256 def __classcall_private__(cls, n, field=QQ):
2257 S = _complex_hermitian_basis(n)
2258 (Qs, T) = _multiplication_table_from_matrix_basis(S)
2259
2260 fdeja = super(ComplexHermitianEJA, cls)
2261 return fdeja.__classcall_private__(cls,
2262 field,
2263 Qs,
2264 rank=n,
2265 natural_basis=T)
2266
2267 def inner_product(self, x, y):
2268 # Since a+bi on the diagonal is represented as
2269 #
2270 # a + bi = [ a b ]
2271 # [ -b a ],
2272 #
2273 # we'll double-count the "a" entries if we take the trace of
2274 # the embedding.
2275 return _matrix_ip(x,y)/2
2276
2277
2278 class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
2279 """
2280 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
2281 matrices, the usual symmetric Jordan product, and the
2282 real-part-of-trace inner product. It has dimension `2n^2 - n` over
2283 the reals.
2284
2285 SETUP::
2286
2287 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
2288
2289 TESTS:
2290
2291 The degree of this algebra is `n^2`::
2292
2293 sage: set_random_seed()
2294 sage: n = ZZ.random_element(1,5)
2295 sage: J = QuaternionHermitianEJA(n)
2296 sage: J.degree() == 2*(n^2) - n
2297 True
2298
2299 The Jordan multiplication is what we think it is::
2300
2301 sage: set_random_seed()
2302 sage: n = ZZ.random_element(1,5)
2303 sage: J = QuaternionHermitianEJA(n)
2304 sage: x = J.random_element()
2305 sage: y = J.random_element()
2306 sage: actual = (x*y).natural_representation()
2307 sage: X = x.natural_representation()
2308 sage: Y = y.natural_representation()
2309 sage: expected = (X*Y + Y*X)/2
2310 sage: actual == expected
2311 True
2312 sage: J(expected) == x*y
2313 True
2314
2315 """
2316 @staticmethod
2317 def __classcall_private__(cls, n, field=QQ):
2318 S = _quaternion_hermitian_basis(n)
2319 (Qs, T) = _multiplication_table_from_matrix_basis(S)
2320
2321 fdeja = super(QuaternionHermitianEJA, cls)
2322 return fdeja.__classcall_private__(cls,
2323 field,
2324 Qs,
2325 rank=n,
2326 natural_basis=T)
2327
2328 def inner_product(self, x, y):
2329 # Since a+bi+cj+dk on the diagonal is represented as
2330 #
2331 # a + bi +cj + dk = [ a b c d]
2332 # [ -b a -d c]
2333 # [ -c d a -b]
2334 # [ -d -c b a],
2335 #
2336 # we'll quadruple-count the "a" entries if we take the trace of
2337 # the embedding.
2338 return _matrix_ip(x,y)/4
2339
2340
2341 class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra):
2342 """
2343 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2344 with the usual inner product and jordan product ``x*y =
2345 (<x_bar,y_bar>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
2346 the reals.
2347
2348 SETUP::
2349
2350 sage: from mjo.eja.eja_algebra import JordanSpinEJA
2351
2352 EXAMPLES:
2353
2354 This multiplication table can be verified by hand::
2355
2356 sage: J = JordanSpinEJA(4)
2357 sage: e0,e1,e2,e3 = J.gens()
2358 sage: e0*e0
2359 e0
2360 sage: e0*e1
2361 e1
2362 sage: e0*e2
2363 e2
2364 sage: e0*e3
2365 e3
2366 sage: e1*e2
2367 0
2368 sage: e1*e3
2369 0
2370 sage: e2*e3
2371 0
2372
2373 """
2374 @staticmethod
2375 def __classcall_private__(cls, n, field=QQ):
2376 Qs = []
2377 id_matrix = matrix.identity(field, n)
2378 for i in xrange(n):
2379 ei = id_matrix.column(i)
2380 Qi = matrix.zero(field, n)
2381 Qi.set_row(0, ei)
2382 Qi.set_column(0, ei)
2383 Qi += matrix.diagonal(n, [ei[0]]*n)
2384 # The addition of the diagonal matrix adds an extra ei[0] in the
2385 # upper-left corner of the matrix.
2386 Qi[0,0] = Qi[0,0] * ~field(2)
2387 Qs.append(Qi)
2388
2389 # The rank of the spin algebra is two, unless we're in a
2390 # one-dimensional ambient space (because the rank is bounded by
2391 # the ambient dimension).
2392 fdeja = super(JordanSpinEJA, cls)
2393 return fdeja.__classcall_private__(cls, field, Qs, rank=min(n,2))
2394
2395 def inner_product(self, x, y):
2396 return _usual_ip(x,y)