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