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