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