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