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