]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/euclidean_jordan_algebra.py
eja: add vector_space() method for EJAs to get the ambient space.
[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 = self.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 def vector_space(self):
744 """
745 Return the vector space that underlies this algebra.
746
747 EXAMPLES::
748
749 sage: J = RealSymmetricEJA(2)
750 sage: J.vector_space()
751 Vector space of dimension 3 over Rational Field
752
753 """
754 return self.zero().vector().parent().ambient_vector_space()
755
756
757 class Element(FiniteDimensionalAlgebraElement):
758 """
759 An element of a Euclidean Jordan algebra.
760 """
761
762 def __dir__(self):
763 """
764 Oh man, I should not be doing this. This hides the "disabled"
765 methods ``left_matrix`` and ``matrix`` from introspection;
766 in particular it removes them from tab-completion.
767 """
768 return filter(lambda s: s not in ['left_matrix', 'matrix'],
769 dir(self.__class__) )
770
771
772 def __init__(self, A, elt=None):
773 """
774 EXAMPLES:
775
776 The identity in `S^n` is converted to the identity in the EJA::
777
778 sage: J = RealSymmetricEJA(3)
779 sage: I = identity_matrix(QQ,3)
780 sage: J(I) == J.one()
781 True
782
783 This skew-symmetric matrix can't be represented in the EJA::
784
785 sage: J = RealSymmetricEJA(3)
786 sage: A = matrix(QQ,3, lambda i,j: i-j)
787 sage: J(A)
788 Traceback (most recent call last):
789 ...
790 ArithmeticError: vector is not in free module
791
792 """
793 # Goal: if we're given a matrix, and if it lives in our
794 # parent algebra's "natural ambient space," convert it
795 # into an algebra element.
796 #
797 # The catch is, we make a recursive call after converting
798 # the given matrix into a vector that lives in the algebra.
799 # This we need to try the parent class initializer first,
800 # to avoid recursing forever if we're given something that
801 # already fits into the algebra, but also happens to live
802 # in the parent's "natural ambient space" (this happens with
803 # vectors in R^n).
804 try:
805 FiniteDimensionalAlgebraElement.__init__(self, A, elt)
806 except ValueError:
807 natural_basis = A.natural_basis()
808 if elt in natural_basis[0].matrix_space():
809 # Thanks for nothing! Matrix spaces aren't vector
810 # spaces in Sage, so we have to figure out its
811 # natural-basis coordinates ourselves.
812 V = VectorSpace(elt.base_ring(), elt.nrows()**2)
813 W = V.span( _mat2vec(s) for s in natural_basis )
814 coords = W.coordinates(_mat2vec(elt))
815 FiniteDimensionalAlgebraElement.__init__(self, A, coords)
816
817 def __pow__(self, n):
818 """
819 Return ``self`` raised to the power ``n``.
820
821 Jordan algebras are always power-associative; see for
822 example Faraut and Koranyi, Proposition II.1.2 (ii).
823
824 .. WARNING:
825
826 We have to override this because our superclass uses row vectors
827 instead of column vectors! We, on the other hand, assume column
828 vectors everywhere.
829
830 EXAMPLES::
831
832 sage: set_random_seed()
833 sage: x = random_eja().random_element()
834 sage: x.operator_matrix()*x.vector() == (x^2).vector()
835 True
836
837 A few examples of power-associativity::
838
839 sage: set_random_seed()
840 sage: x = random_eja().random_element()
841 sage: x*(x*x)*(x*x) == x^5
842 True
843 sage: (x*x)*(x*x*x) == x^5
844 True
845
846 We also know that powers operator-commute (Koecher, Chapter
847 III, Corollary 1)::
848
849 sage: set_random_seed()
850 sage: x = random_eja().random_element()
851 sage: m = ZZ.random_element(0,10)
852 sage: n = ZZ.random_element(0,10)
853 sage: Lxm = (x^m).operator_matrix()
854 sage: Lxn = (x^n).operator_matrix()
855 sage: Lxm*Lxn == Lxn*Lxm
856 True
857
858 """
859 A = self.parent()
860 if n == 0:
861 return A.one()
862 elif n == 1:
863 return self
864 else:
865 return A( (self.operator_matrix()**(n-1))*self.vector() )
866
867
868 def apply_univariate_polynomial(self, p):
869 """
870 Apply the univariate polynomial ``p`` to this element.
871
872 A priori, SageMath won't allow us to apply a univariate
873 polynomial to an element of an EJA, because we don't know
874 that EJAs are rings (they are usually not associative). Of
875 course, we know that EJAs are power-associative, so the
876 operation is ultimately kosher. This function sidesteps
877 the CAS to get the answer we want and expect.
878
879 EXAMPLES::
880
881 sage: R = PolynomialRing(QQ, 't')
882 sage: t = R.gen(0)
883 sage: p = t^4 - t^3 + 5*t - 2
884 sage: J = RealCartesianProductEJA(5)
885 sage: J.one().apply_univariate_polynomial(p) == 3*J.one()
886 True
887
888 TESTS:
889
890 We should always get back an element of the algebra::
891
892 sage: set_random_seed()
893 sage: p = PolynomialRing(QQ, 't').random_element()
894 sage: J = random_eja()
895 sage: x = J.random_element()
896 sage: x.apply_univariate_polynomial(p) in J
897 True
898
899 """
900 if len(p.variables()) > 1:
901 raise ValueError("not a univariate polynomial")
902 P = self.parent()
903 R = P.base_ring()
904 # Convert the coeficcients to the parent's base ring,
905 # because a priori they might live in an (unnecessarily)
906 # larger ring for which P.sum() would fail below.
907 cs = [ R(c) for c in p.coefficients(sparse=False) ]
908 return P.sum( cs[k]*(self**k) for k in range(len(cs)) )
909
910
911 def characteristic_polynomial(self):
912 """
913 Return the characteristic polynomial of this element.
914
915 EXAMPLES:
916
917 The rank of `R^3` is three, and the minimal polynomial of
918 the identity element is `(t-1)` from which it follows that
919 the characteristic polynomial should be `(t-1)^3`::
920
921 sage: J = RealCartesianProductEJA(3)
922 sage: J.one().characteristic_polynomial()
923 t^3 - 3*t^2 + 3*t - 1
924
925 Likewise, the characteristic of the zero element in the
926 rank-three algebra `R^{n}` should be `t^{3}`::
927
928 sage: J = RealCartesianProductEJA(3)
929 sage: J.zero().characteristic_polynomial()
930 t^3
931
932 The characteristic polynomial of an element should evaluate
933 to zero on that element::
934
935 sage: set_random_seed()
936 sage: x = RealCartesianProductEJA(3).random_element()
937 sage: p = x.characteristic_polynomial()
938 sage: x.apply_univariate_polynomial(p)
939 0
940
941 """
942 p = self.parent().characteristic_polynomial()
943 return p(*self.vector())
944
945
946 def inner_product(self, other):
947 """
948 Return the parent algebra's inner product of myself and ``other``.
949
950 EXAMPLES:
951
952 The inner product in the Jordan spin algebra is the usual
953 inner product on `R^n` (this example only works because the
954 basis for the Jordan algebra is the standard basis in `R^n`)::
955
956 sage: J = JordanSpinEJA(3)
957 sage: x = vector(QQ,[1,2,3])
958 sage: y = vector(QQ,[4,5,6])
959 sage: x.inner_product(y)
960 32
961 sage: J(x).inner_product(J(y))
962 32
963
964 The inner product on `S^n` is `<X,Y> = trace(X*Y)`, where
965 multiplication is the usual matrix multiplication in `S^n`,
966 so the inner product of the identity matrix with itself
967 should be the `n`::
968
969 sage: J = RealSymmetricEJA(3)
970 sage: J.one().inner_product(J.one())
971 3
972
973 Likewise, the inner product on `C^n` is `<X,Y> =
974 Re(trace(X*Y))`, where we must necessarily take the real
975 part because the product of Hermitian matrices may not be
976 Hermitian::
977
978 sage: J = ComplexHermitianEJA(3)
979 sage: J.one().inner_product(J.one())
980 3
981
982 Ditto for the quaternions::
983
984 sage: J = QuaternionHermitianEJA(3)
985 sage: J.one().inner_product(J.one())
986 3
987
988 TESTS:
989
990 Ensure that we can always compute an inner product, and that
991 it gives us back a real number::
992
993 sage: set_random_seed()
994 sage: J = random_eja()
995 sage: x = J.random_element()
996 sage: y = J.random_element()
997 sage: x.inner_product(y) in RR
998 True
999
1000 """
1001 P = self.parent()
1002 if not other in P:
1003 raise TypeError("'other' must live in the same algebra")
1004
1005 return P.inner_product(self, other)
1006
1007
1008 def operator_commutes_with(self, other):
1009 """
1010 Return whether or not this element operator-commutes
1011 with ``other``.
1012
1013 EXAMPLES:
1014
1015 The definition of a Jordan algebra says that any element
1016 operator-commutes with its square::
1017
1018 sage: set_random_seed()
1019 sage: x = random_eja().random_element()
1020 sage: x.operator_commutes_with(x^2)
1021 True
1022
1023 TESTS:
1024
1025 Test Lemma 1 from Chapter III of Koecher::
1026
1027 sage: set_random_seed()
1028 sage: J = random_eja()
1029 sage: u = J.random_element()
1030 sage: v = J.random_element()
1031 sage: lhs = u.operator_commutes_with(u*v)
1032 sage: rhs = v.operator_commutes_with(u^2)
1033 sage: lhs == rhs
1034 True
1035
1036 """
1037 if not other in self.parent():
1038 raise TypeError("'other' must live in the same algebra")
1039
1040 A = self.operator_matrix()
1041 B = other.operator_matrix()
1042 return (A*B == B*A)
1043
1044
1045 def det(self):
1046 """
1047 Return my determinant, the product of my eigenvalues.
1048
1049 EXAMPLES::
1050
1051 sage: J = JordanSpinEJA(2)
1052 sage: e0,e1 = J.gens()
1053 sage: x = sum( J.gens() )
1054 sage: x.det()
1055 0
1056
1057 ::
1058
1059 sage: J = JordanSpinEJA(3)
1060 sage: e0,e1,e2 = J.gens()
1061 sage: x = sum( J.gens() )
1062 sage: x.det()
1063 -1
1064
1065 TESTS:
1066
1067 An element is invertible if and only if its determinant is
1068 non-zero::
1069
1070 sage: set_random_seed()
1071 sage: x = random_eja().random_element()
1072 sage: x.is_invertible() == (x.det() != 0)
1073 True
1074
1075 """
1076 P = self.parent()
1077 r = P.rank()
1078 p = P._charpoly_coeff(0)
1079 # The _charpoly_coeff function already adds the factor of
1080 # -1 to ensure that _charpoly_coeff(0) is really what
1081 # appears in front of t^{0} in the charpoly. However,
1082 # we want (-1)^r times THAT for the determinant.
1083 return ((-1)**r)*p(*self.vector())
1084
1085
1086 def inverse(self):
1087 """
1088 Return the Jordan-multiplicative inverse of this element.
1089
1090 ALGORITHM:
1091
1092 We appeal to the quadratic representation as in Koecher's
1093 Theorem 12 in Chapter III, Section 5.
1094
1095 EXAMPLES:
1096
1097 The inverse in the spin factor algebra is given in Alizadeh's
1098 Example 11.11::
1099
1100 sage: set_random_seed()
1101 sage: n = ZZ.random_element(1,10)
1102 sage: J = JordanSpinEJA(n)
1103 sage: x = J.random_element()
1104 sage: while not x.is_invertible():
1105 ....: x = J.random_element()
1106 sage: x_vec = x.vector()
1107 sage: x0 = x_vec[0]
1108 sage: x_bar = x_vec[1:]
1109 sage: coeff = ~(x0^2 - x_bar.inner_product(x_bar))
1110 sage: inv_vec = x_vec.parent()([x0] + (-x_bar).list())
1111 sage: x_inverse = coeff*inv_vec
1112 sage: x.inverse() == J(x_inverse)
1113 True
1114
1115 TESTS:
1116
1117 The identity element is its own inverse::
1118
1119 sage: set_random_seed()
1120 sage: J = random_eja()
1121 sage: J.one().inverse() == J.one()
1122 True
1123
1124 If an element has an inverse, it acts like one::
1125
1126 sage: set_random_seed()
1127 sage: J = random_eja()
1128 sage: x = J.random_element()
1129 sage: (not x.is_invertible()) or (x.inverse()*x == J.one())
1130 True
1131
1132 The inverse of the inverse is what we started with::
1133
1134 sage: set_random_seed()
1135 sage: J = random_eja()
1136 sage: x = J.random_element()
1137 sage: (not x.is_invertible()) or (x.inverse().inverse() == x)
1138 True
1139
1140 The zero element is never invertible::
1141
1142 sage: set_random_seed()
1143 sage: J = random_eja().zero().inverse()
1144 Traceback (most recent call last):
1145 ...
1146 ValueError: element is not invertible
1147
1148 """
1149 if not self.is_invertible():
1150 raise ValueError("element is not invertible")
1151
1152 return (~self.quadratic_representation())(self)
1153
1154
1155 def is_invertible(self):
1156 """
1157 Return whether or not this element is invertible.
1158
1159 We can't use the superclass method because it relies on
1160 the algebra being associative.
1161
1162 ALGORITHM:
1163
1164 The usual way to do this is to check if the determinant is
1165 zero, but we need the characteristic polynomial for the
1166 determinant. The minimal polynomial is a lot easier to get,
1167 so we use Corollary 2 in Chapter V of Koecher to check
1168 whether or not the paren't algebra's zero element is a root
1169 of this element's minimal polynomial.
1170
1171 TESTS:
1172
1173 The identity element is always invertible::
1174
1175 sage: set_random_seed()
1176 sage: J = random_eja()
1177 sage: J.one().is_invertible()
1178 True
1179
1180 The zero element is never invertible::
1181
1182 sage: set_random_seed()
1183 sage: J = random_eja()
1184 sage: J.zero().is_invertible()
1185 False
1186
1187 """
1188 zero = self.parent().zero()
1189 p = self.minimal_polynomial()
1190 return not (p(zero) == zero)
1191
1192
1193 def is_nilpotent(self):
1194 """
1195 Return whether or not some power of this element is zero.
1196
1197 The superclass method won't work unless we're in an
1198 associative algebra, and we aren't. However, we generate
1199 an assocoative subalgebra and we're nilpotent there if and
1200 only if we're nilpotent here (probably).
1201
1202 TESTS:
1203
1204 The identity element is never nilpotent::
1205
1206 sage: set_random_seed()
1207 sage: random_eja().one().is_nilpotent()
1208 False
1209
1210 The additive identity is always nilpotent::
1211
1212 sage: set_random_seed()
1213 sage: random_eja().zero().is_nilpotent()
1214 True
1215
1216 """
1217 # The element we're going to call "is_nilpotent()" on.
1218 # Either myself, interpreted as an element of a finite-
1219 # dimensional algebra, or an element of an associative
1220 # subalgebra.
1221 elt = None
1222
1223 if self.parent().is_associative():
1224 elt = FiniteDimensionalAlgebraElement(self.parent(), self)
1225 else:
1226 V = self.span_of_powers()
1227 assoc_subalg = self.subalgebra_generated_by()
1228 # Mis-design warning: the basis used for span_of_powers()
1229 # and subalgebra_generated_by() must be the same, and in
1230 # the same order!
1231 elt = assoc_subalg(V.coordinates(self.vector()))
1232
1233 # Recursive call, but should work since elt lives in an
1234 # associative algebra.
1235 return elt.is_nilpotent()
1236
1237
1238 def is_regular(self):
1239 """
1240 Return whether or not this is a regular element.
1241
1242 EXAMPLES:
1243
1244 The identity element always has degree one, but any element
1245 linearly-independent from it is regular::
1246
1247 sage: J = JordanSpinEJA(5)
1248 sage: J.one().is_regular()
1249 False
1250 sage: e0, e1, e2, e3, e4 = J.gens() # e0 is the identity
1251 sage: for x in J.gens():
1252 ....: (J.one() + x).is_regular()
1253 False
1254 True
1255 True
1256 True
1257 True
1258
1259 """
1260 return self.degree() == self.parent().rank()
1261
1262
1263 def degree(self):
1264 """
1265 Compute the degree of this element the straightforward way
1266 according to the definition; by appending powers to a list
1267 and figuring out its dimension (that is, whether or not
1268 they're linearly dependent).
1269
1270 EXAMPLES::
1271
1272 sage: J = JordanSpinEJA(4)
1273 sage: J.one().degree()
1274 1
1275 sage: e0,e1,e2,e3 = J.gens()
1276 sage: (e0 - e1).degree()
1277 2
1278
1279 In the spin factor algebra (of rank two), all elements that
1280 aren't multiples of the identity are regular::
1281
1282 sage: set_random_seed()
1283 sage: n = ZZ.random_element(1,10)
1284 sage: J = JordanSpinEJA(n)
1285 sage: x = J.random_element()
1286 sage: x == x.coefficient(0)*J.one() or x.degree() == 2
1287 True
1288
1289 """
1290 return self.span_of_powers().dimension()
1291
1292
1293 def left_matrix(self):
1294 """
1295 Our parent class defines ``left_matrix`` and ``matrix``
1296 methods whose names are misleading. We don't want them.
1297 """
1298 raise NotImplementedError("use operator_matrix() instead")
1299
1300 matrix = left_matrix
1301
1302
1303 def minimal_polynomial(self):
1304 """
1305 Return the minimal polynomial of this element,
1306 as a function of the variable `t`.
1307
1308 ALGORITHM:
1309
1310 We restrict ourselves to the associative subalgebra
1311 generated by this element, and then return the minimal
1312 polynomial of this element's operator matrix (in that
1313 subalgebra). This works by Baes Proposition 2.3.16.
1314
1315 TESTS:
1316
1317 The minimal polynomial of the identity and zero elements are
1318 always the same::
1319
1320 sage: set_random_seed()
1321 sage: J = random_eja()
1322 sage: J.one().minimal_polynomial()
1323 t - 1
1324 sage: J.zero().minimal_polynomial()
1325 t
1326
1327 The degree of an element is (by one definition) the degree
1328 of its minimal polynomial::
1329
1330 sage: set_random_seed()
1331 sage: x = random_eja().random_element()
1332 sage: x.degree() == x.minimal_polynomial().degree()
1333 True
1334
1335 The minimal polynomial and the characteristic polynomial coincide
1336 and are known (see Alizadeh, Example 11.11) for all elements of
1337 the spin factor algebra that aren't scalar multiples of the
1338 identity::
1339
1340 sage: set_random_seed()
1341 sage: n = ZZ.random_element(2,10)
1342 sage: J = JordanSpinEJA(n)
1343 sage: y = J.random_element()
1344 sage: while y == y.coefficient(0)*J.one():
1345 ....: y = J.random_element()
1346 sage: y0 = y.vector()[0]
1347 sage: y_bar = y.vector()[1:]
1348 sage: actual = y.minimal_polynomial()
1349 sage: t = PolynomialRing(J.base_ring(),'t').gen(0)
1350 sage: expected = t^2 - 2*y0*t + (y0^2 - norm(y_bar)^2)
1351 sage: bool(actual == expected)
1352 True
1353
1354 The minimal polynomial should always kill its element::
1355
1356 sage: set_random_seed()
1357 sage: x = random_eja().random_element()
1358 sage: p = x.minimal_polynomial()
1359 sage: x.apply_univariate_polynomial(p)
1360 0
1361
1362 """
1363 V = self.span_of_powers()
1364 assoc_subalg = self.subalgebra_generated_by()
1365 # Mis-design warning: the basis used for span_of_powers()
1366 # and subalgebra_generated_by() must be the same, and in
1367 # the same order!
1368 elt = assoc_subalg(V.coordinates(self.vector()))
1369
1370 # We get back a symbolic polynomial in 'x' but want a real
1371 # polynomial in 't'.
1372 p_of_x = elt.operator_matrix().minimal_polynomial()
1373 return p_of_x.change_variable_name('t')
1374
1375
1376 def natural_representation(self):
1377 """
1378 Return a more-natural representation of this element.
1379
1380 Every finite-dimensional Euclidean Jordan Algebra is a
1381 direct sum of five simple algebras, four of which comprise
1382 Hermitian matrices. This method returns the original
1383 "natural" representation of this element as a Hermitian
1384 matrix, if it has one. If not, you get the usual representation.
1385
1386 EXAMPLES::
1387
1388 sage: J = ComplexHermitianEJA(3)
1389 sage: J.one()
1390 e0 + e5 + e8
1391 sage: J.one().natural_representation()
1392 [1 0 0 0 0 0]
1393 [0 1 0 0 0 0]
1394 [0 0 1 0 0 0]
1395 [0 0 0 1 0 0]
1396 [0 0 0 0 1 0]
1397 [0 0 0 0 0 1]
1398
1399 ::
1400
1401 sage: J = QuaternionHermitianEJA(3)
1402 sage: J.one()
1403 e0 + e9 + e14
1404 sage: J.one().natural_representation()
1405 [1 0 0 0 0 0 0 0 0 0 0 0]
1406 [0 1 0 0 0 0 0 0 0 0 0 0]
1407 [0 0 1 0 0 0 0 0 0 0 0 0]
1408 [0 0 0 1 0 0 0 0 0 0 0 0]
1409 [0 0 0 0 1 0 0 0 0 0 0 0]
1410 [0 0 0 0 0 1 0 0 0 0 0 0]
1411 [0 0 0 0 0 0 1 0 0 0 0 0]
1412 [0 0 0 0 0 0 0 1 0 0 0 0]
1413 [0 0 0 0 0 0 0 0 1 0 0 0]
1414 [0 0 0 0 0 0 0 0 0 1 0 0]
1415 [0 0 0 0 0 0 0 0 0 0 1 0]
1416 [0 0 0 0 0 0 0 0 0 0 0 1]
1417
1418 """
1419 B = self.parent().natural_basis()
1420 W = B[0].matrix_space()
1421 return W.linear_combination(zip(self.vector(), B))
1422
1423
1424 def operator(self):
1425 """
1426 Return the left-multiplication-by-this-element
1427 operator on the ambient algebra.
1428
1429 TESTS::
1430
1431 sage: set_random_seed()
1432 sage: J = random_eja()
1433 sage: x = J.random_element()
1434 sage: y = J.random_element()
1435 sage: x.operator()(y) == x*y
1436 True
1437 sage: y.operator()(x) == x*y
1438 True
1439
1440 """
1441 P = self.parent()
1442 return FiniteDimensionalEuclideanJordanAlgebraMorphism(
1443 Hom(P,P),
1444 self.operator_matrix() )
1445
1446
1447
1448 def operator_matrix(self):
1449 """
1450 Return the matrix that represents left- (or right-)
1451 multiplication by this element in the parent algebra.
1452
1453 We implement this ourselves to work around the fact that
1454 our parent class represents everything with row vectors.
1455
1456 EXAMPLES:
1457
1458 Test the first polarization identity from my notes, Koecher Chapter
1459 III, or from Baes (2.3)::
1460
1461 sage: set_random_seed()
1462 sage: J = random_eja()
1463 sage: x = J.random_element()
1464 sage: y = J.random_element()
1465 sage: Lx = x.operator_matrix()
1466 sage: Ly = y.operator_matrix()
1467 sage: Lxx = (x*x).operator_matrix()
1468 sage: Lxy = (x*y).operator_matrix()
1469 sage: bool(2*Lx*Lxy + Ly*Lxx == 2*Lxy*Lx + Lxx*Ly)
1470 True
1471
1472 Test the second polarization identity from my notes or from
1473 Baes (2.4)::
1474
1475 sage: set_random_seed()
1476 sage: J = random_eja()
1477 sage: x = J.random_element()
1478 sage: y = J.random_element()
1479 sage: z = J.random_element()
1480 sage: Lx = x.operator_matrix()
1481 sage: Ly = y.operator_matrix()
1482 sage: Lz = z.operator_matrix()
1483 sage: Lzy = (z*y).operator_matrix()
1484 sage: Lxy = (x*y).operator_matrix()
1485 sage: Lxz = (x*z).operator_matrix()
1486 sage: bool(Lx*Lzy + Lz*Lxy + Ly*Lxz == Lzy*Lx + Lxy*Lz + Lxz*Ly)
1487 True
1488
1489 Test the third polarization identity from my notes or from
1490 Baes (2.5)::
1491
1492 sage: set_random_seed()
1493 sage: J = random_eja()
1494 sage: u = J.random_element()
1495 sage: y = J.random_element()
1496 sage: z = J.random_element()
1497 sage: Lu = u.operator_matrix()
1498 sage: Ly = y.operator_matrix()
1499 sage: Lz = z.operator_matrix()
1500 sage: Lzy = (z*y).operator_matrix()
1501 sage: Luy = (u*y).operator_matrix()
1502 sage: Luz = (u*z).operator_matrix()
1503 sage: Luyz = (u*(y*z)).operator_matrix()
1504 sage: lhs = Lu*Lzy + Lz*Luy + Ly*Luz
1505 sage: rhs = Luyz + Ly*Lu*Lz + Lz*Lu*Ly
1506 sage: bool(lhs == rhs)
1507 True
1508
1509 """
1510 fda_elt = FiniteDimensionalAlgebraElement(self.parent(), self)
1511 return fda_elt.matrix().transpose()
1512
1513
1514 def quadratic_representation(self, other=None):
1515 """
1516 Return the quadratic representation of this element.
1517
1518 EXAMPLES:
1519
1520 The explicit form in the spin factor algebra is given by
1521 Alizadeh's Example 11.12::
1522
1523 sage: set_random_seed()
1524 sage: n = ZZ.random_element(1,10)
1525 sage: J = JordanSpinEJA(n)
1526 sage: x = J.random_element()
1527 sage: x_vec = x.vector()
1528 sage: x0 = x_vec[0]
1529 sage: x_bar = x_vec[1:]
1530 sage: A = matrix(QQ, 1, [x_vec.inner_product(x_vec)])
1531 sage: B = 2*x0*x_bar.row()
1532 sage: C = 2*x0*x_bar.column()
1533 sage: D = identity_matrix(QQ, n-1)
1534 sage: D = (x0^2 - x_bar.inner_product(x_bar))*D
1535 sage: D = D + 2*x_bar.tensor_product(x_bar)
1536 sage: Q = block_matrix(2,2,[A,B,C,D])
1537 sage: Q == x.quadratic_representation().matrix()
1538 True
1539
1540 Test all of the properties from Theorem 11.2 in Alizadeh::
1541
1542 sage: set_random_seed()
1543 sage: J = random_eja()
1544 sage: x = J.random_element()
1545 sage: y = J.random_element()
1546 sage: Lx = x.operator()
1547 sage: Lxx = (x*x).operator()
1548 sage: Qx = x.quadratic_representation()
1549 sage: Qy = y.quadratic_representation()
1550 sage: Qxy = x.quadratic_representation(y)
1551 sage: Qex = J.one().quadratic_representation(x)
1552 sage: n = ZZ.random_element(10)
1553 sage: Qxn = (x^n).quadratic_representation()
1554
1555 Property 1:
1556
1557 sage: 2*Qxy == (x+y).quadratic_representation() - Qx - Qy
1558 True
1559
1560 Property 2:
1561
1562 sage: alpha = QQ.random_element()
1563 sage: (alpha*x).quadratic_representation() == (alpha^2)*Qx
1564 True
1565
1566 Property 3:
1567
1568 sage: not x.is_invertible() or ( Qx(x.inverse()) == x )
1569 True
1570
1571 sage: not x.is_invertible() or (
1572 ....: ~Qx
1573 ....: ==
1574 ....: x.inverse().quadratic_representation() )
1575 True
1576
1577 sage: Qxy(J.one()) == x*y
1578 True
1579
1580 Property 4:
1581
1582 sage: not x.is_invertible() or (
1583 ....: x.quadratic_representation(x.inverse())*Qx
1584 ....: == Qx*x.quadratic_representation(x.inverse()) )
1585 True
1586
1587 sage: not x.is_invertible() or (
1588 ....: x.quadratic_representation(x.inverse())*Qx
1589 ....: ==
1590 ....: 2*x.operator()*Qex - Qx )
1591 True
1592
1593 sage: 2*x.operator()*Qex - Qx == Lxx
1594 True
1595
1596 Property 5:
1597
1598 sage: Qy(x).quadratic_representation() == Qy*Qx*Qy
1599 True
1600
1601 Property 6:
1602
1603 sage: Qxn == (Qx)^n
1604 True
1605
1606 Property 7:
1607
1608 sage: not x.is_invertible() or (
1609 ....: Qx*x.inverse().operator() == Lx )
1610 True
1611
1612 Property 8:
1613
1614 sage: not x.operator_commutes_with(y) or (
1615 ....: Qx(y)^n == Qxn(y^n) )
1616 True
1617
1618 """
1619 if other is None:
1620 other=self
1621 elif not other in self.parent():
1622 raise TypeError("'other' must live in the same algebra")
1623
1624 L = self.operator()
1625 M = other.operator()
1626 return ( L*M + M*L - (self*other).operator() )
1627
1628
1629 def span_of_powers(self):
1630 """
1631 Return the vector space spanned by successive powers of
1632 this element.
1633 """
1634 # The dimension of the subalgebra can't be greater than
1635 # the big algebra, so just put everything into a list
1636 # and let span() get rid of the excess.
1637 #
1638 # We do the extra ambient_vector_space() in case we're messing
1639 # with polynomials and the direct parent is a module.
1640 V = self.parent().vector_space()
1641 return V.span( (self**d).vector() for d in xrange(V.dimension()) )
1642
1643
1644 def subalgebra_generated_by(self):
1645 """
1646 Return the associative subalgebra of the parent EJA generated
1647 by this element.
1648
1649 TESTS::
1650
1651 sage: set_random_seed()
1652 sage: x = random_eja().random_element()
1653 sage: x.subalgebra_generated_by().is_associative()
1654 True
1655
1656 Squaring in the subalgebra should be the same thing as
1657 squaring in the superalgebra::
1658
1659 sage: set_random_seed()
1660 sage: x = random_eja().random_element()
1661 sage: u = x.subalgebra_generated_by().random_element()
1662 sage: u.operator_matrix()*u.vector() == (u**2).vector()
1663 True
1664
1665 """
1666 # First get the subspace spanned by the powers of myself...
1667 V = self.span_of_powers()
1668 F = self.base_ring()
1669
1670 # Now figure out the entries of the right-multiplication
1671 # matrix for the successive basis elements b0, b1,... of
1672 # that subspace.
1673 mats = []
1674 for b_right in V.basis():
1675 eja_b_right = self.parent()(b_right)
1676 b_right_rows = []
1677 # The first row of the right-multiplication matrix by
1678 # b1 is what we get if we apply that matrix to b1. The
1679 # second row of the right multiplication matrix by b1
1680 # is what we get when we apply that matrix to b2...
1681 #
1682 # IMPORTANT: this assumes that all vectors are COLUMN
1683 # vectors, unlike our superclass (which uses row vectors).
1684 for b_left in V.basis():
1685 eja_b_left = self.parent()(b_left)
1686 # Multiply in the original EJA, but then get the
1687 # coordinates from the subalgebra in terms of its
1688 # basis.
1689 this_row = V.coordinates((eja_b_left*eja_b_right).vector())
1690 b_right_rows.append(this_row)
1691 b_right_matrix = matrix(F, b_right_rows)
1692 mats.append(b_right_matrix)
1693
1694 # It's an algebra of polynomials in one element, and EJAs
1695 # are power-associative.
1696 #
1697 # TODO: choose generator names intelligently.
1698 return FiniteDimensionalEuclideanJordanAlgebra(F, mats, assume_associative=True, names='f')
1699
1700
1701 def subalgebra_idempotent(self):
1702 """
1703 Find an idempotent in the associative subalgebra I generate
1704 using Proposition 2.3.5 in Baes.
1705
1706 TESTS::
1707
1708 sage: set_random_seed()
1709 sage: J = random_eja()
1710 sage: x = J.random_element()
1711 sage: while x.is_nilpotent():
1712 ....: x = J.random_element()
1713 sage: c = x.subalgebra_idempotent()
1714 sage: c^2 == c
1715 True
1716
1717 """
1718 if self.is_nilpotent():
1719 raise ValueError("this only works with non-nilpotent elements!")
1720
1721 V = self.span_of_powers()
1722 J = self.subalgebra_generated_by()
1723 # Mis-design warning: the basis used for span_of_powers()
1724 # and subalgebra_generated_by() must be the same, and in
1725 # the same order!
1726 u = J(V.coordinates(self.vector()))
1727
1728 # The image of the matrix of left-u^m-multiplication
1729 # will be minimal for some natural number s...
1730 s = 0
1731 minimal_dim = V.dimension()
1732 for i in xrange(1, V.dimension()):
1733 this_dim = (u**i).operator_matrix().image().dimension()
1734 if this_dim < minimal_dim:
1735 minimal_dim = this_dim
1736 s = i
1737
1738 # Now minimal_matrix should correspond to the smallest
1739 # non-zero subspace in Baes's (or really, Koecher's)
1740 # proposition.
1741 #
1742 # However, we need to restrict the matrix to work on the
1743 # subspace... or do we? Can't we just solve, knowing that
1744 # A(c) = u^(s+1) should have a solution in the big space,
1745 # too?
1746 #
1747 # Beware, solve_right() means that we're using COLUMN vectors.
1748 # Our FiniteDimensionalAlgebraElement superclass uses rows.
1749 u_next = u**(s+1)
1750 A = u_next.operator_matrix()
1751 c_coordinates = A.solve_right(u_next.vector())
1752
1753 # Now c_coordinates is the idempotent we want, but it's in
1754 # the coordinate system of the subalgebra.
1755 #
1756 # We need the basis for J, but as elements of the parent algebra.
1757 #
1758 basis = [self.parent(v) for v in V.basis()]
1759 return self.parent().linear_combination(zip(c_coordinates, basis))
1760
1761
1762 def trace(self):
1763 """
1764 Return my trace, the sum of my eigenvalues.
1765
1766 EXAMPLES::
1767
1768 sage: J = JordanSpinEJA(3)
1769 sage: x = sum(J.gens())
1770 sage: x.trace()
1771 2
1772
1773 ::
1774
1775 sage: J = RealCartesianProductEJA(5)
1776 sage: J.one().trace()
1777 5
1778
1779 TESTS:
1780
1781 The trace of an element is a real number::
1782
1783 sage: set_random_seed()
1784 sage: J = random_eja()
1785 sage: J.random_element().trace() in J.base_ring()
1786 True
1787
1788 """
1789 P = self.parent()
1790 r = P.rank()
1791 p = P._charpoly_coeff(r-1)
1792 # The _charpoly_coeff function already adds the factor of
1793 # -1 to ensure that _charpoly_coeff(r-1) is really what
1794 # appears in front of t^{r-1} in the charpoly. However,
1795 # we want the negative of THAT for the trace.
1796 return -p(*self.vector())
1797
1798
1799 def trace_inner_product(self, other):
1800 """
1801 Return the trace inner product of myself and ``other``.
1802
1803 TESTS:
1804
1805 The trace inner product is commutative::
1806
1807 sage: set_random_seed()
1808 sage: J = random_eja()
1809 sage: x = J.random_element(); y = J.random_element()
1810 sage: x.trace_inner_product(y) == y.trace_inner_product(x)
1811 True
1812
1813 The trace inner product is bilinear::
1814
1815 sage: set_random_seed()
1816 sage: J = random_eja()
1817 sage: x = J.random_element()
1818 sage: y = J.random_element()
1819 sage: z = J.random_element()
1820 sage: a = QQ.random_element();
1821 sage: actual = (a*(x+z)).trace_inner_product(y)
1822 sage: expected = ( a*x.trace_inner_product(y) +
1823 ....: a*z.trace_inner_product(y) )
1824 sage: actual == expected
1825 True
1826 sage: actual = x.trace_inner_product(a*(y+z))
1827 sage: expected = ( a*x.trace_inner_product(y) +
1828 ....: a*x.trace_inner_product(z) )
1829 sage: actual == expected
1830 True
1831
1832 The trace inner product satisfies the compatibility
1833 condition in the definition of a Euclidean Jordan algebra::
1834
1835 sage: set_random_seed()
1836 sage: J = random_eja()
1837 sage: x = J.random_element()
1838 sage: y = J.random_element()
1839 sage: z = J.random_element()
1840 sage: (x*y).trace_inner_product(z) == y.trace_inner_product(x*z)
1841 True
1842
1843 """
1844 if not other in self.parent():
1845 raise TypeError("'other' must live in the same algebra")
1846
1847 return (self*other).trace()
1848
1849
1850 class RealCartesianProductEJA(FiniteDimensionalEuclideanJordanAlgebra):
1851 """
1852 Return the Euclidean Jordan Algebra corresponding to the set
1853 `R^n` under the Hadamard product.
1854
1855 Note: this is nothing more than the Cartesian product of ``n``
1856 copies of the spin algebra. Once Cartesian product algebras
1857 are implemented, this can go.
1858
1859 EXAMPLES:
1860
1861 This multiplication table can be verified by hand::
1862
1863 sage: J = RealCartesianProductEJA(3)
1864 sage: e0,e1,e2 = J.gens()
1865 sage: e0*e0
1866 e0
1867 sage: e0*e1
1868 0
1869 sage: e0*e2
1870 0
1871 sage: e1*e1
1872 e1
1873 sage: e1*e2
1874 0
1875 sage: e2*e2
1876 e2
1877
1878 """
1879 @staticmethod
1880 def __classcall_private__(cls, n, field=QQ):
1881 # The FiniteDimensionalAlgebra constructor takes a list of
1882 # matrices, the ith representing right multiplication by the ith
1883 # basis element in the vector space. So if e_1 = (1,0,0), then
1884 # right (Hadamard) multiplication of x by e_1 picks out the first
1885 # component of x; and likewise for the ith basis element e_i.
1886 Qs = [ matrix(field, n, n, lambda k,j: 1*(k == j == i))
1887 for i in xrange(n) ]
1888
1889 fdeja = super(RealCartesianProductEJA, cls)
1890 return fdeja.__classcall_private__(cls, field, Qs, rank=n)
1891
1892 def inner_product(self, x, y):
1893 return _usual_ip(x,y)
1894
1895
1896 def random_eja():
1897 """
1898 Return a "random" finite-dimensional Euclidean Jordan Algebra.
1899
1900 ALGORITHM:
1901
1902 For now, we choose a random natural number ``n`` (greater than zero)
1903 and then give you back one of the following:
1904
1905 * The cartesian product of the rational numbers ``n`` times; this is
1906 ``QQ^n`` with the Hadamard product.
1907
1908 * The Jordan spin algebra on ``QQ^n``.
1909
1910 * The ``n``-by-``n`` rational symmetric matrices with the symmetric
1911 product.
1912
1913 * The ``n``-by-``n`` complex-rational Hermitian matrices embedded
1914 in the space of ``2n``-by-``2n`` real symmetric matrices.
1915
1916 * The ``n``-by-``n`` quaternion-rational Hermitian matrices embedded
1917 in the space of ``4n``-by-``4n`` real symmetric matrices.
1918
1919 Later this might be extended to return Cartesian products of the
1920 EJAs above.
1921
1922 TESTS::
1923
1924 sage: random_eja()
1925 Euclidean Jordan algebra of degree...
1926
1927 """
1928
1929 # The max_n component lets us choose different upper bounds on the
1930 # value "n" that gets passed to the constructor. This is needed
1931 # because e.g. R^{10} is reasonable to test, while the Hermitian
1932 # 10-by-10 quaternion matrices are not.
1933 (constructor, max_n) = choice([(RealCartesianProductEJA, 6),
1934 (JordanSpinEJA, 6),
1935 (RealSymmetricEJA, 5),
1936 (ComplexHermitianEJA, 4),
1937 (QuaternionHermitianEJA, 3)])
1938 n = ZZ.random_element(1, max_n)
1939 return constructor(n, field=QQ)
1940
1941
1942
1943 def _real_symmetric_basis(n, field=QQ):
1944 """
1945 Return a basis for the space of real symmetric n-by-n matrices.
1946 """
1947 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1948 # coordinates.
1949 S = []
1950 for i in xrange(n):
1951 for j in xrange(i+1):
1952 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1953 if i == j:
1954 Sij = Eij
1955 else:
1956 # Beware, orthogonal but not normalized!
1957 Sij = Eij + Eij.transpose()
1958 S.append(Sij)
1959 return tuple(S)
1960
1961
1962 def _complex_hermitian_basis(n, field=QQ):
1963 """
1964 Returns a basis for the space of complex Hermitian n-by-n matrices.
1965
1966 TESTS::
1967
1968 sage: set_random_seed()
1969 sage: n = ZZ.random_element(1,5)
1970 sage: all( M.is_symmetric() for M in _complex_hermitian_basis(n) )
1971 True
1972
1973 """
1974 F = QuadraticField(-1, 'I')
1975 I = F.gen()
1976
1977 # This is like the symmetric case, but we need to be careful:
1978 #
1979 # * We want conjugate-symmetry, not just symmetry.
1980 # * The diagonal will (as a result) be real.
1981 #
1982 S = []
1983 for i in xrange(n):
1984 for j in xrange(i+1):
1985 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1986 if i == j:
1987 Sij = _embed_complex_matrix(Eij)
1988 S.append(Sij)
1989 else:
1990 # Beware, orthogonal but not normalized! The second one
1991 # has a minus because it's conjugated.
1992 Sij_real = _embed_complex_matrix(Eij + Eij.transpose())
1993 S.append(Sij_real)
1994 Sij_imag = _embed_complex_matrix(I*Eij - I*Eij.transpose())
1995 S.append(Sij_imag)
1996 return tuple(S)
1997
1998
1999 def _quaternion_hermitian_basis(n, field=QQ):
2000 """
2001 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
2002
2003 TESTS::
2004
2005 sage: set_random_seed()
2006 sage: n = ZZ.random_element(1,5)
2007 sage: all( M.is_symmetric() for M in _quaternion_hermitian_basis(n) )
2008 True
2009
2010 """
2011 Q = QuaternionAlgebra(QQ,-1,-1)
2012 I,J,K = Q.gens()
2013
2014 # This is like the symmetric case, but we need to be careful:
2015 #
2016 # * We want conjugate-symmetry, not just symmetry.
2017 # * The diagonal will (as a result) be real.
2018 #
2019 S = []
2020 for i in xrange(n):
2021 for j in xrange(i+1):
2022 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
2023 if i == j:
2024 Sij = _embed_quaternion_matrix(Eij)
2025 S.append(Sij)
2026 else:
2027 # Beware, orthogonal but not normalized! The second,
2028 # third, and fourth ones have a minus because they're
2029 # conjugated.
2030 Sij_real = _embed_quaternion_matrix(Eij + Eij.transpose())
2031 S.append(Sij_real)
2032 Sij_I = _embed_quaternion_matrix(I*Eij - I*Eij.transpose())
2033 S.append(Sij_I)
2034 Sij_J = _embed_quaternion_matrix(J*Eij - J*Eij.transpose())
2035 S.append(Sij_J)
2036 Sij_K = _embed_quaternion_matrix(K*Eij - K*Eij.transpose())
2037 S.append(Sij_K)
2038 return tuple(S)
2039
2040
2041 def _mat2vec(m):
2042 return vector(m.base_ring(), m.list())
2043
2044 def _vec2mat(v):
2045 return matrix(v.base_ring(), sqrt(v.degree()), v.list())
2046
2047 def _multiplication_table_from_matrix_basis(basis):
2048 """
2049 At least three of the five simple Euclidean Jordan algebras have the
2050 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
2051 multiplication on the right is matrix multiplication. Given a basis
2052 for the underlying matrix space, this function returns a
2053 multiplication table (obtained by looping through the basis
2054 elements) for an algebra of those matrices. A reordered copy
2055 of the basis is also returned to work around the fact that
2056 the ``span()`` in this function will change the order of the basis
2057 from what we think it is, to... something else.
2058 """
2059 # In S^2, for example, we nominally have four coordinates even
2060 # though the space is of dimension three only. The vector space V
2061 # is supposed to hold the entire long vector, and the subspace W
2062 # of V will be spanned by the vectors that arise from symmetric
2063 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
2064 field = basis[0].base_ring()
2065 dimension = basis[0].nrows()
2066
2067 V = VectorSpace(field, dimension**2)
2068 W = V.span( _mat2vec(s) for s in basis )
2069
2070 # Taking the span above reorders our basis (thanks, jerk!) so we
2071 # need to put our "matrix basis" in the same order as the
2072 # (reordered) vector basis.
2073 S = tuple( _vec2mat(b) for b in W.basis() )
2074
2075 Qs = []
2076 for s in S:
2077 # Brute force the multiplication-by-s matrix by looping
2078 # through all elements of the basis and doing the computation
2079 # to find out what the corresponding row should be. BEWARE:
2080 # these multiplication tables won't be symmetric! It therefore
2081 # becomes REALLY IMPORTANT that the underlying algebra
2082 # constructor uses ROW vectors and not COLUMN vectors. That's
2083 # why we're computing rows here and not columns.
2084 Q_rows = []
2085 for t in S:
2086 this_row = _mat2vec((s*t + t*s)/2)
2087 Q_rows.append(W.coordinates(this_row))
2088 Q = matrix(field, W.dimension(), Q_rows)
2089 Qs.append(Q)
2090
2091 return (Qs, S)
2092
2093
2094 def _embed_complex_matrix(M):
2095 """
2096 Embed the n-by-n complex matrix ``M`` into the space of real
2097 matrices of size 2n-by-2n via the map the sends each entry `z = a +
2098 bi` to the block matrix ``[[a,b],[-b,a]]``.
2099
2100 EXAMPLES::
2101
2102 sage: F = QuadraticField(-1,'i')
2103 sage: x1 = F(4 - 2*i)
2104 sage: x2 = F(1 + 2*i)
2105 sage: x3 = F(-i)
2106 sage: x4 = F(6)
2107 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
2108 sage: _embed_complex_matrix(M)
2109 [ 4 -2| 1 2]
2110 [ 2 4|-2 1]
2111 [-----+-----]
2112 [ 0 -1| 6 0]
2113 [ 1 0| 0 6]
2114
2115 TESTS:
2116
2117 Embedding is a homomorphism (isomorphism, in fact)::
2118
2119 sage: set_random_seed()
2120 sage: n = ZZ.random_element(5)
2121 sage: F = QuadraticField(-1, 'i')
2122 sage: X = random_matrix(F, n)
2123 sage: Y = random_matrix(F, n)
2124 sage: actual = _embed_complex_matrix(X) * _embed_complex_matrix(Y)
2125 sage: expected = _embed_complex_matrix(X*Y)
2126 sage: actual == expected
2127 True
2128
2129 """
2130 n = M.nrows()
2131 if M.ncols() != n:
2132 raise ValueError("the matrix 'M' must be square")
2133 field = M.base_ring()
2134 blocks = []
2135 for z in M.list():
2136 a = z.real()
2137 b = z.imag()
2138 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
2139
2140 # We can drop the imaginaries here.
2141 return block_matrix(field.base_ring(), n, blocks)
2142
2143
2144 def _unembed_complex_matrix(M):
2145 """
2146 The inverse of _embed_complex_matrix().
2147
2148 EXAMPLES::
2149
2150 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
2151 ....: [-2, 1, -4, 3],
2152 ....: [ 9, 10, 11, 12],
2153 ....: [-10, 9, -12, 11] ])
2154 sage: _unembed_complex_matrix(A)
2155 [ 2*i + 1 4*i + 3]
2156 [ 10*i + 9 12*i + 11]
2157
2158 TESTS:
2159
2160 Unembedding is the inverse of embedding::
2161
2162 sage: set_random_seed()
2163 sage: F = QuadraticField(-1, 'i')
2164 sage: M = random_matrix(F, 3)
2165 sage: _unembed_complex_matrix(_embed_complex_matrix(M)) == M
2166 True
2167
2168 """
2169 n = ZZ(M.nrows())
2170 if M.ncols() != n:
2171 raise ValueError("the matrix 'M' must be square")
2172 if not n.mod(2).is_zero():
2173 raise ValueError("the matrix 'M' must be a complex embedding")
2174
2175 F = QuadraticField(-1, 'i')
2176 i = F.gen()
2177
2178 # Go top-left to bottom-right (reading order), converting every
2179 # 2-by-2 block we see to a single complex element.
2180 elements = []
2181 for k in xrange(n/2):
2182 for j in xrange(n/2):
2183 submat = M[2*k:2*k+2,2*j:2*j+2]
2184 if submat[0,0] != submat[1,1]:
2185 raise ValueError('bad on-diagonal submatrix')
2186 if submat[0,1] != -submat[1,0]:
2187 raise ValueError('bad off-diagonal submatrix')
2188 z = submat[0,0] + submat[0,1]*i
2189 elements.append(z)
2190
2191 return matrix(F, n/2, elements)
2192
2193
2194 def _embed_quaternion_matrix(M):
2195 """
2196 Embed the n-by-n quaternion matrix ``M`` into the space of real
2197 matrices of size 4n-by-4n by first sending each quaternion entry
2198 `z = a + bi + cj + dk` to the block-complex matrix
2199 ``[[a + bi, c+di],[-c + di, a-bi]]`, and then embedding those into
2200 a real matrix.
2201
2202 EXAMPLES::
2203
2204 sage: Q = QuaternionAlgebra(QQ,-1,-1)
2205 sage: i,j,k = Q.gens()
2206 sage: x = 1 + 2*i + 3*j + 4*k
2207 sage: M = matrix(Q, 1, [[x]])
2208 sage: _embed_quaternion_matrix(M)
2209 [ 1 2 3 4]
2210 [-2 1 -4 3]
2211 [-3 4 1 -2]
2212 [-4 -3 2 1]
2213
2214 Embedding is a homomorphism (isomorphism, in fact)::
2215
2216 sage: set_random_seed()
2217 sage: n = ZZ.random_element(5)
2218 sage: Q = QuaternionAlgebra(QQ,-1,-1)
2219 sage: X = random_matrix(Q, n)
2220 sage: Y = random_matrix(Q, n)
2221 sage: actual = _embed_quaternion_matrix(X)*_embed_quaternion_matrix(Y)
2222 sage: expected = _embed_quaternion_matrix(X*Y)
2223 sage: actual == expected
2224 True
2225
2226 """
2227 quaternions = M.base_ring()
2228 n = M.nrows()
2229 if M.ncols() != n:
2230 raise ValueError("the matrix 'M' must be square")
2231
2232 F = QuadraticField(-1, 'i')
2233 i = F.gen()
2234
2235 blocks = []
2236 for z in M.list():
2237 t = z.coefficient_tuple()
2238 a = t[0]
2239 b = t[1]
2240 c = t[2]
2241 d = t[3]
2242 cplx_matrix = matrix(F, 2, [[ a + b*i, c + d*i],
2243 [-c + d*i, a - b*i]])
2244 blocks.append(_embed_complex_matrix(cplx_matrix))
2245
2246 # We should have real entries by now, so use the realest field
2247 # we've got for the return value.
2248 return block_matrix(quaternions.base_ring(), n, blocks)
2249
2250
2251 def _unembed_quaternion_matrix(M):
2252 """
2253 The inverse of _embed_quaternion_matrix().
2254
2255 EXAMPLES::
2256
2257 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
2258 ....: [-2, 1, -4, 3],
2259 ....: [-3, 4, 1, -2],
2260 ....: [-4, -3, 2, 1]])
2261 sage: _unembed_quaternion_matrix(M)
2262 [1 + 2*i + 3*j + 4*k]
2263
2264 TESTS:
2265
2266 Unembedding is the inverse of embedding::
2267
2268 sage: set_random_seed()
2269 sage: Q = QuaternionAlgebra(QQ, -1, -1)
2270 sage: M = random_matrix(Q, 3)
2271 sage: _unembed_quaternion_matrix(_embed_quaternion_matrix(M)) == M
2272 True
2273
2274 """
2275 n = ZZ(M.nrows())
2276 if M.ncols() != n:
2277 raise ValueError("the matrix 'M' must be square")
2278 if not n.mod(4).is_zero():
2279 raise ValueError("the matrix 'M' must be a complex embedding")
2280
2281 Q = QuaternionAlgebra(QQ,-1,-1)
2282 i,j,k = Q.gens()
2283
2284 # Go top-left to bottom-right (reading order), converting every
2285 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
2286 # quaternion block.
2287 elements = []
2288 for l in xrange(n/4):
2289 for m in xrange(n/4):
2290 submat = _unembed_complex_matrix(M[4*l:4*l+4,4*m:4*m+4])
2291 if submat[0,0] != submat[1,1].conjugate():
2292 raise ValueError('bad on-diagonal submatrix')
2293 if submat[0,1] != -submat[1,0].conjugate():
2294 raise ValueError('bad off-diagonal submatrix')
2295 z = submat[0,0].real() + submat[0,0].imag()*i
2296 z += submat[0,1].real()*j + submat[0,1].imag()*k
2297 elements.append(z)
2298
2299 return matrix(Q, n/4, elements)
2300
2301
2302 # The usual inner product on R^n.
2303 def _usual_ip(x,y):
2304 return x.vector().inner_product(y.vector())
2305
2306 # The inner product used for the real symmetric simple EJA.
2307 # We keep it as a separate function because e.g. the complex
2308 # algebra uses the same inner product, except divided by 2.
2309 def _matrix_ip(X,Y):
2310 X_mat = X.natural_representation()
2311 Y_mat = Y.natural_representation()
2312 return (X_mat*Y_mat).trace()
2313
2314
2315 class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra):
2316 """
2317 The rank-n simple EJA consisting of real symmetric n-by-n
2318 matrices, the usual symmetric Jordan product, and the trace inner
2319 product. It has dimension `(n^2 + n)/2` over the reals.
2320
2321 EXAMPLES::
2322
2323 sage: J = RealSymmetricEJA(2)
2324 sage: e0, e1, e2 = J.gens()
2325 sage: e0*e0
2326 e0
2327 sage: e1*e1
2328 e0 + e2
2329 sage: e2*e2
2330 e2
2331
2332 TESTS:
2333
2334 The degree of this algebra is `(n^2 + n) / 2`::
2335
2336 sage: set_random_seed()
2337 sage: n = ZZ.random_element(1,5)
2338 sage: J = RealSymmetricEJA(n)
2339 sage: J.degree() == (n^2 + n)/2
2340 True
2341
2342 The Jordan multiplication is what we think it is::
2343
2344 sage: set_random_seed()
2345 sage: n = ZZ.random_element(1,5)
2346 sage: J = RealSymmetricEJA(n)
2347 sage: x = J.random_element()
2348 sage: y = J.random_element()
2349 sage: actual = (x*y).natural_representation()
2350 sage: X = x.natural_representation()
2351 sage: Y = y.natural_representation()
2352 sage: expected = (X*Y + Y*X)/2
2353 sage: actual == expected
2354 True
2355 sage: J(expected) == x*y
2356 True
2357
2358 """
2359 @staticmethod
2360 def __classcall_private__(cls, n, field=QQ):
2361 S = _real_symmetric_basis(n, field=field)
2362 (Qs, T) = _multiplication_table_from_matrix_basis(S)
2363
2364 fdeja = super(RealSymmetricEJA, cls)
2365 return fdeja.__classcall_private__(cls,
2366 field,
2367 Qs,
2368 rank=n,
2369 natural_basis=T)
2370
2371 def inner_product(self, x, y):
2372 return _matrix_ip(x,y)
2373
2374
2375 class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
2376 """
2377 The rank-n simple EJA consisting of complex Hermitian n-by-n
2378 matrices over the real numbers, the usual symmetric Jordan product,
2379 and the real-part-of-trace inner product. It has dimension `n^2` over
2380 the reals.
2381
2382 TESTS:
2383
2384 The degree of this algebra is `n^2`::
2385
2386 sage: set_random_seed()
2387 sage: n = ZZ.random_element(1,5)
2388 sage: J = ComplexHermitianEJA(n)
2389 sage: J.degree() == n^2
2390 True
2391
2392 The Jordan multiplication is what we think it is::
2393
2394 sage: set_random_seed()
2395 sage: n = ZZ.random_element(1,5)
2396 sage: J = ComplexHermitianEJA(n)
2397 sage: x = J.random_element()
2398 sage: y = J.random_element()
2399 sage: actual = (x*y).natural_representation()
2400 sage: X = x.natural_representation()
2401 sage: Y = y.natural_representation()
2402 sage: expected = (X*Y + Y*X)/2
2403 sage: actual == expected
2404 True
2405 sage: J(expected) == x*y
2406 True
2407
2408 """
2409 @staticmethod
2410 def __classcall_private__(cls, n, field=QQ):
2411 S = _complex_hermitian_basis(n)
2412 (Qs, T) = _multiplication_table_from_matrix_basis(S)
2413
2414 fdeja = super(ComplexHermitianEJA, cls)
2415 return fdeja.__classcall_private__(cls,
2416 field,
2417 Qs,
2418 rank=n,
2419 natural_basis=T)
2420
2421 def inner_product(self, x, y):
2422 # Since a+bi on the diagonal is represented as
2423 #
2424 # a + bi = [ a b ]
2425 # [ -b a ],
2426 #
2427 # we'll double-count the "a" entries if we take the trace of
2428 # the embedding.
2429 return _matrix_ip(x,y)/2
2430
2431
2432 class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
2433 """
2434 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
2435 matrices, the usual symmetric Jordan product, and the
2436 real-part-of-trace inner product. It has dimension `2n^2 - n` over
2437 the reals.
2438
2439 TESTS:
2440
2441 The degree of this algebra is `n^2`::
2442
2443 sage: set_random_seed()
2444 sage: n = ZZ.random_element(1,5)
2445 sage: J = QuaternionHermitianEJA(n)
2446 sage: J.degree() == 2*(n^2) - n
2447 True
2448
2449 The Jordan multiplication is what we think it is::
2450
2451 sage: set_random_seed()
2452 sage: n = ZZ.random_element(1,5)
2453 sage: J = QuaternionHermitianEJA(n)
2454 sage: x = J.random_element()
2455 sage: y = J.random_element()
2456 sage: actual = (x*y).natural_representation()
2457 sage: X = x.natural_representation()
2458 sage: Y = y.natural_representation()
2459 sage: expected = (X*Y + Y*X)/2
2460 sage: actual == expected
2461 True
2462 sage: J(expected) == x*y
2463 True
2464
2465 """
2466 @staticmethod
2467 def __classcall_private__(cls, n, field=QQ):
2468 S = _quaternion_hermitian_basis(n)
2469 (Qs, T) = _multiplication_table_from_matrix_basis(S)
2470
2471 fdeja = super(QuaternionHermitianEJA, cls)
2472 return fdeja.__classcall_private__(cls,
2473 field,
2474 Qs,
2475 rank=n,
2476 natural_basis=T)
2477
2478 def inner_product(self, x, y):
2479 # Since a+bi+cj+dk on the diagonal is represented as
2480 #
2481 # a + bi +cj + dk = [ a b c d]
2482 # [ -b a -d c]
2483 # [ -c d a -b]
2484 # [ -d -c b a],
2485 #
2486 # we'll quadruple-count the "a" entries if we take the trace of
2487 # the embedding.
2488 return _matrix_ip(x,y)/4
2489
2490
2491 class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra):
2492 """
2493 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2494 with the usual inner product and jordan product ``x*y =
2495 (<x_bar,y_bar>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
2496 the reals.
2497
2498 EXAMPLES:
2499
2500 This multiplication table can be verified by hand::
2501
2502 sage: J = JordanSpinEJA(4)
2503 sage: e0,e1,e2,e3 = J.gens()
2504 sage: e0*e0
2505 e0
2506 sage: e0*e1
2507 e1
2508 sage: e0*e2
2509 e2
2510 sage: e0*e3
2511 e3
2512 sage: e1*e2
2513 0
2514 sage: e1*e3
2515 0
2516 sage: e2*e3
2517 0
2518
2519 """
2520 @staticmethod
2521 def __classcall_private__(cls, n, field=QQ):
2522 Qs = []
2523 id_matrix = identity_matrix(field, n)
2524 for i in xrange(n):
2525 ei = id_matrix.column(i)
2526 Qi = zero_matrix(field, n)
2527 Qi.set_row(0, ei)
2528 Qi.set_column(0, ei)
2529 Qi += diagonal_matrix(n, [ei[0]]*n)
2530 # The addition of the diagonal matrix adds an extra ei[0] in the
2531 # upper-left corner of the matrix.
2532 Qi[0,0] = Qi[0,0] * ~field(2)
2533 Qs.append(Qi)
2534
2535 # The rank of the spin algebra is two, unless we're in a
2536 # one-dimensional ambient space (because the rank is bounded by
2537 # the ambient dimension).
2538 fdeja = super(JordanSpinEJA, cls)
2539 return fdeja.__classcall_private__(cls, field, Qs, rank=min(n,2))
2540
2541 def inner_product(self, x, y):
2542 return _usual_ip(x,y)