]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: normalize the real symmetric matrix basis.
[sage.d.git] / mjo / eja / eja_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.algebras.quatalg.quaternion_algebra import QuaternionAlgebra
9 from sage.categories.magmatic_algebras import MagmaticAlgebras
10 from sage.combinat.free_module import CombinatorialFreeModule
11 from sage.matrix.constructor import matrix
12 from sage.matrix.matrix_space import MatrixSpace
13 from sage.misc.cachefunc import cached_method
14 from sage.misc.prandom import choice
15 from sage.misc.table import table
16 from sage.modules.free_module import FreeModule, VectorSpace
17 from sage.rings.integer_ring import ZZ
18 from sage.rings.number_field.number_field import QuadraticField
19 from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing
20 from sage.rings.rational_field import QQ
21 from sage.structure.element import is_Matrix
22
23 from mjo.eja.eja_element import FiniteDimensionalEuclideanJordanAlgebraElement
24 from mjo.eja.eja_utils import _mat2vec
25
26 class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule):
27 # This is an ugly hack needed to prevent the category framework
28 # from implementing a coercion from our base ring (e.g. the
29 # rationals) into the algebra. First of all -- such a coercion is
30 # nonsense to begin with. But more importantly, it tries to do so
31 # in the category of rings, and since our algebras aren't
32 # associative they generally won't be rings.
33 _no_generic_basering_coercion = True
34
35 def __init__(self,
36 field,
37 mult_table,
38 rank,
39 prefix='e',
40 category=None,
41 natural_basis=None):
42 """
43 SETUP::
44
45 sage: from mjo.eja.eja_algebra import random_eja
46
47 EXAMPLES:
48
49 By definition, Jordan multiplication commutes::
50
51 sage: set_random_seed()
52 sage: J = random_eja()
53 sage: x = J.random_element()
54 sage: y = J.random_element()
55 sage: x*y == y*x
56 True
57
58 """
59 self._rank = rank
60 self._natural_basis = natural_basis
61
62 if category is None:
63 category = MagmaticAlgebras(field).FiniteDimensional()
64 category = category.WithBasis().Unital()
65
66 fda = super(FiniteDimensionalEuclideanJordanAlgebra, self)
67 fda.__init__(field,
68 range(len(mult_table)),
69 prefix=prefix,
70 category=category)
71 self.print_options(bracket='')
72
73 # The multiplication table we're given is necessarily in terms
74 # of vectors, because we don't have an algebra yet for
75 # anything to be an element of. However, it's faster in the
76 # long run to have the multiplication table be in terms of
77 # algebra elements. We do this after calling the superclass
78 # constructor so that from_vector() knows what to do.
79 self._multiplication_table = [ map(lambda x: self.from_vector(x), ls)
80 for ls in mult_table ]
81
82
83 def _element_constructor_(self, elt):
84 """
85 Construct an element of this algebra from its natural
86 representation.
87
88 This gets called only after the parent element _call_ method
89 fails to find a coercion for the argument.
90
91 SETUP::
92
93 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
94 ....: RealCartesianProductEJA,
95 ....: RealSymmetricEJA)
96
97 EXAMPLES:
98
99 The identity in `S^n` is converted to the identity in the EJA::
100
101 sage: J = RealSymmetricEJA(3)
102 sage: I = matrix.identity(QQ,3)
103 sage: J(I) == J.one()
104 True
105
106 This skew-symmetric matrix can't be represented in the EJA::
107
108 sage: J = RealSymmetricEJA(3)
109 sage: A = matrix(QQ,3, lambda i,j: i-j)
110 sage: J(A)
111 Traceback (most recent call last):
112 ...
113 ArithmeticError: vector is not in free module
114
115 TESTS:
116
117 Ensure that we can convert any element of the two non-matrix
118 simple algebras (whose natural representations are their usual
119 vector representations) back and forth faithfully::
120
121 sage: set_random_seed()
122 sage: J = RealCartesianProductEJA(5)
123 sage: x = J.random_element()
124 sage: J(x.to_vector().column()) == x
125 True
126 sage: J = JordanSpinEJA(5)
127 sage: x = J.random_element()
128 sage: J(x.to_vector().column()) == x
129 True
130
131 """
132 if elt == 0:
133 # The superclass implementation of random_element()
134 # needs to be able to coerce "0" into the algebra.
135 return self.zero()
136
137 natural_basis = self.natural_basis()
138 basis_space = natural_basis[0].matrix_space()
139 if elt not in basis_space:
140 raise ValueError("not a naturally-represented algebra element")
141
142 # Thanks for nothing! Matrix spaces aren't vector spaces in
143 # Sage, so we have to figure out its natural-basis coordinates
144 # ourselves. We use the basis space's ring instead of the
145 # element's ring because the basis space might be an algebraic
146 # closure whereas the base ring of the 3-by-3 identity matrix
147 # could be QQ instead of QQbar.
148 V = VectorSpace(basis_space.base_ring(), elt.nrows()*elt.ncols())
149 W = V.span_of_basis( _mat2vec(s) for s in natural_basis )
150 coords = W.coordinate_vector(_mat2vec(elt))
151 return self.from_vector(coords)
152
153
154 def _repr_(self):
155 """
156 Return a string representation of ``self``.
157
158 SETUP::
159
160 sage: from mjo.eja.eja_algebra import JordanSpinEJA
161
162 TESTS:
163
164 Ensure that it says what we think it says::
165
166 sage: JordanSpinEJA(2, field=QQ)
167 Euclidean Jordan algebra of dimension 2 over Rational Field
168 sage: JordanSpinEJA(3, field=RDF)
169 Euclidean Jordan algebra of dimension 3 over Real Double Field
170
171 """
172 fmt = "Euclidean Jordan algebra of dimension {} over {}"
173 return fmt.format(self.dimension(), self.base_ring())
174
175 def product_on_basis(self, i, j):
176 return self._multiplication_table[i][j]
177
178 def _a_regular_element(self):
179 """
180 Guess a regular element. Needed to compute the basis for our
181 characteristic polynomial coefficients.
182
183 SETUP::
184
185 sage: from mjo.eja.eja_algebra import random_eja
186
187 TESTS:
188
189 Ensure that this hacky method succeeds for every algebra that we
190 know how to construct::
191
192 sage: set_random_seed()
193 sage: J = random_eja()
194 sage: J._a_regular_element().is_regular()
195 True
196
197 """
198 gs = self.gens()
199 z = self.sum( (i+1)*gs[i] for i in range(len(gs)) )
200 if not z.is_regular():
201 raise ValueError("don't know a regular element")
202 return z
203
204
205 @cached_method
206 def _charpoly_basis_space(self):
207 """
208 Return the vector space spanned by the basis used in our
209 characteristic polynomial coefficients. This is used not only to
210 compute those coefficients, but also any time we need to
211 evaluate the coefficients (like when we compute the trace or
212 determinant).
213 """
214 z = self._a_regular_element()
215 # Don't use the parent vector space directly here in case this
216 # happens to be a subalgebra. In that case, we would be e.g.
217 # two-dimensional but span_of_basis() would expect three
218 # coordinates.
219 V = VectorSpace(self.base_ring(), self.vector_space().dimension())
220 basis = [ (z**k).to_vector() for k in range(self.rank()) ]
221 V1 = V.span_of_basis( basis )
222 b = (V1.basis() + V1.complement().basis())
223 return V.span_of_basis(b)
224
225
226 @cached_method
227 def _charpoly_coeff(self, i):
228 """
229 Return the coefficient polynomial "a_{i}" of this algebra's
230 general characteristic polynomial.
231
232 Having this be a separate cached method lets us compute and
233 store the trace/determinant (a_{r-1} and a_{0} respectively)
234 separate from the entire characteristic polynomial.
235 """
236 (A_of_x, x, xr, detA) = self._charpoly_matrix_system()
237 R = A_of_x.base_ring()
238 if i >= self.rank():
239 # Guaranteed by theory
240 return R.zero()
241
242 # Danger: the in-place modification is done for performance
243 # reasons (reconstructing a matrix with huge polynomial
244 # entries is slow), but I don't know how cached_method works,
245 # so it's highly possible that we're modifying some global
246 # list variable by reference, here. In other words, you
247 # probably shouldn't call this method twice on the same
248 # algebra, at the same time, in two threads
249 Ai_orig = A_of_x.column(i)
250 A_of_x.set_column(i,xr)
251 numerator = A_of_x.det()
252 A_of_x.set_column(i,Ai_orig)
253
254 # We're relying on the theory here to ensure that each a_i is
255 # indeed back in R, and the added negative signs are to make
256 # the whole charpoly expression sum to zero.
257 return R(-numerator/detA)
258
259
260 @cached_method
261 def _charpoly_matrix_system(self):
262 """
263 Compute the matrix whose entries A_ij are polynomials in
264 X1,...,XN, the vector ``x`` of variables X1,...,XN, the vector
265 corresponding to `x^r` and the determinent of the matrix A =
266 [A_ij]. In other words, all of the fixed (cachable) data needed
267 to compute the coefficients of the characteristic polynomial.
268 """
269 r = self.rank()
270 n = self.dimension()
271
272 # Turn my vector space into a module so that "vectors" can
273 # have multivatiate polynomial entries.
274 names = tuple('X' + str(i) for i in range(1,n+1))
275 R = PolynomialRing(self.base_ring(), names)
276
277 # Using change_ring() on the parent's vector space doesn't work
278 # here because, in a subalgebra, that vector space has a basis
279 # and change_ring() tries to bring the basis along with it. And
280 # that doesn't work unless the new ring is a PID, which it usually
281 # won't be.
282 V = FreeModule(R,n)
283
284 # Now let x = (X1,X2,...,Xn) be the vector whose entries are
285 # indeterminates...
286 x = V(names)
287
288 # And figure out the "left multiplication by x" matrix in
289 # that setting.
290 lmbx_cols = []
291 monomial_matrices = [ self.monomial(i).operator().matrix()
292 for i in range(n) ] # don't recompute these!
293 for k in range(n):
294 ek = self.monomial(k).to_vector()
295 lmbx_cols.append(
296 sum( x[i]*(monomial_matrices[i]*ek)
297 for i in range(n) ) )
298 Lx = matrix.column(R, lmbx_cols)
299
300 # Now we can compute powers of x "symbolically"
301 x_powers = [self.one().to_vector(), x]
302 for d in range(2, r+1):
303 x_powers.append( Lx*(x_powers[-1]) )
304
305 idmat = matrix.identity(R, n)
306
307 W = self._charpoly_basis_space()
308 W = W.change_ring(R.fraction_field())
309
310 # Starting with the standard coordinates x = (X1,X2,...,Xn)
311 # and then converting the entries to W-coordinates allows us
312 # to pass in the standard coordinates to the charpoly and get
313 # back the right answer. Specifically, with x = (X1,X2,...,Xn),
314 # we have
315 #
316 # W.coordinates(x^2) eval'd at (standard z-coords)
317 # =
318 # W-coords of (z^2)
319 # =
320 # W-coords of (standard coords of x^2 eval'd at std-coords of z)
321 #
322 # We want the middle equivalent thing in our matrix, but use
323 # the first equivalent thing instead so that we can pass in
324 # standard coordinates.
325 x_powers = [ W.coordinate_vector(xp) for xp in x_powers ]
326 l2 = [idmat.column(k-1) for k in range(r+1, n+1)]
327 A_of_x = matrix.column(R, n, (x_powers[:r] + l2))
328 return (A_of_x, x, x_powers[r], A_of_x.det())
329
330
331 @cached_method
332 def characteristic_polynomial(self):
333 """
334 Return a characteristic polynomial that works for all elements
335 of this algebra.
336
337 The resulting polynomial has `n+1` variables, where `n` is the
338 dimension of this algebra. The first `n` variables correspond to
339 the coordinates of an algebra element: when evaluated at the
340 coordinates of an algebra element with respect to a certain
341 basis, the result is a univariate polynomial (in the one
342 remaining variable ``t``), namely the characteristic polynomial
343 of that element.
344
345 SETUP::
346
347 sage: from mjo.eja.eja_algebra import JordanSpinEJA
348
349 EXAMPLES:
350
351 The characteristic polynomial in the spin algebra is given in
352 Alizadeh, Example 11.11::
353
354 sage: J = JordanSpinEJA(3)
355 sage: p = J.characteristic_polynomial(); p
356 X1^2 - X2^2 - X3^2 + (-2*t)*X1 + t^2
357 sage: xvec = J.one().to_vector()
358 sage: p(*xvec)
359 t^2 - 2*t + 1
360
361 """
362 r = self.rank()
363 n = self.dimension()
364
365 # The list of coefficient polynomials a_1, a_2, ..., a_n.
366 a = [ self._charpoly_coeff(i) for i in range(n) ]
367
368 # We go to a bit of trouble here to reorder the
369 # indeterminates, so that it's easier to evaluate the
370 # characteristic polynomial at x's coordinates and get back
371 # something in terms of t, which is what we want.
372 R = a[0].parent()
373 S = PolynomialRing(self.base_ring(),'t')
374 t = S.gen(0)
375 S = PolynomialRing(S, R.variable_names())
376 t = S(t)
377
378 # Note: all entries past the rth should be zero. The
379 # coefficient of the highest power (x^r) is 1, but it doesn't
380 # appear in the solution vector which contains coefficients
381 # for the other powers (to make them sum to x^r).
382 if (r < n):
383 a[r] = 1 # corresponds to x^r
384 else:
385 # When the rank is equal to the dimension, trying to
386 # assign a[r] goes out-of-bounds.
387 a.append(1) # corresponds to x^r
388
389 return sum( a[k]*(t**k) for k in range(len(a)) )
390
391
392 def inner_product(self, x, y):
393 """
394 The inner product associated with this Euclidean Jordan algebra.
395
396 Defaults to the trace inner product, but can be overridden by
397 subclasses if they are sure that the necessary properties are
398 satisfied.
399
400 SETUP::
401
402 sage: from mjo.eja.eja_algebra import random_eja
403
404 EXAMPLES:
405
406 The inner product must satisfy its axiom for this algebra to truly
407 be a Euclidean Jordan Algebra::
408
409 sage: set_random_seed()
410 sage: J = random_eja()
411 sage: x = J.random_element()
412 sage: y = J.random_element()
413 sage: z = J.random_element()
414 sage: (x*y).inner_product(z) == y.inner_product(x*z)
415 True
416
417 """
418 if (not x in self) or (not y in self):
419 raise TypeError("arguments must live in this algebra")
420 return x.trace_inner_product(y)
421
422
423 def is_trivial(self):
424 """
425 Return whether or not this algebra is trivial.
426
427 A trivial algebra contains only the zero element.
428
429 SETUP::
430
431 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
432
433 EXAMPLES::
434
435 sage: J = ComplexHermitianEJA(3)
436 sage: J.is_trivial()
437 False
438 sage: A = J.zero().subalgebra_generated_by()
439 sage: A.is_trivial()
440 True
441
442 """
443 return self.dimension() == 0
444
445
446 def multiplication_table(self):
447 """
448 Return a visual representation of this algebra's multiplication
449 table (on basis elements).
450
451 SETUP::
452
453 sage: from mjo.eja.eja_algebra import JordanSpinEJA
454
455 EXAMPLES::
456
457 sage: J = JordanSpinEJA(4)
458 sage: J.multiplication_table()
459 +----++----+----+----+----+
460 | * || e0 | e1 | e2 | e3 |
461 +====++====+====+====+====+
462 | e0 || e0 | e1 | e2 | e3 |
463 +----++----+----+----+----+
464 | e1 || e1 | e0 | 0 | 0 |
465 +----++----+----+----+----+
466 | e2 || e2 | 0 | e0 | 0 |
467 +----++----+----+----+----+
468 | e3 || e3 | 0 | 0 | e0 |
469 +----++----+----+----+----+
470
471 """
472 M = list(self._multiplication_table) # copy
473 for i in range(len(M)):
474 # M had better be "square"
475 M[i] = [self.monomial(i)] + M[i]
476 M = [["*"] + list(self.gens())] + M
477 return table(M, header_row=True, header_column=True, frame=True)
478
479
480 def natural_basis(self):
481 """
482 Return a more-natural representation of this algebra's basis.
483
484 Every finite-dimensional Euclidean Jordan Algebra is a direct
485 sum of five simple algebras, four of which comprise Hermitian
486 matrices. This method returns the original "natural" basis
487 for our underlying vector space. (Typically, the natural basis
488 is used to construct the multiplication table in the first place.)
489
490 Note that this will always return a matrix. The standard basis
491 in `R^n` will be returned as `n`-by-`1` column matrices.
492
493 SETUP::
494
495 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
496 ....: RealSymmetricEJA)
497
498 EXAMPLES::
499
500 sage: J = RealSymmetricEJA(2)
501 sage: J.basis()
502 Finite family {0: e0, 1: e1, 2: e2}
503 sage: J.natural_basis()
504 (
505 [1 0] [ 0 1/2*sqrt2] [0 0]
506 [0 0], [1/2*sqrt2 0], [0 1]
507 )
508
509 ::
510
511 sage: J = JordanSpinEJA(2)
512 sage: J.basis()
513 Finite family {0: e0, 1: e1}
514 sage: J.natural_basis()
515 (
516 [1] [0]
517 [0], [1]
518 )
519
520 """
521 if self._natural_basis is None:
522 M = self.natural_basis_space()
523 return tuple( M(b.to_vector()) for b in self.basis() )
524 else:
525 return self._natural_basis
526
527
528 def natural_basis_space(self):
529 """
530 Return the matrix space in which this algebra's natural basis
531 elements live.
532 """
533 if self._natural_basis is None or len(self._natural_basis) == 0:
534 return MatrixSpace(self.base_ring(), self.dimension(), 1)
535 else:
536 return self._natural_basis[0].matrix_space()
537
538
539 @cached_method
540 def one(self):
541 """
542 Return the unit element of this algebra.
543
544 SETUP::
545
546 sage: from mjo.eja.eja_algebra import (RealCartesianProductEJA,
547 ....: random_eja)
548
549 EXAMPLES::
550
551 sage: J = RealCartesianProductEJA(5)
552 sage: J.one()
553 e0 + e1 + e2 + e3 + e4
554
555 TESTS:
556
557 The identity element acts like the identity::
558
559 sage: set_random_seed()
560 sage: J = random_eja()
561 sage: x = J.random_element()
562 sage: J.one()*x == x and x*J.one() == x
563 True
564
565 The matrix of the unit element's operator is the identity::
566
567 sage: set_random_seed()
568 sage: J = random_eja()
569 sage: actual = J.one().operator().matrix()
570 sage: expected = matrix.identity(J.base_ring(), J.dimension())
571 sage: actual == expected
572 True
573
574 """
575 # We can brute-force compute the matrices of the operators
576 # that correspond to the basis elements of this algebra.
577 # If some linear combination of those basis elements is the
578 # algebra identity, then the same linear combination of
579 # their matrices has to be the identity matrix.
580 #
581 # Of course, matrices aren't vectors in sage, so we have to
582 # appeal to the "long vectors" isometry.
583 oper_vecs = [ _mat2vec(g.operator().matrix()) for g in self.gens() ]
584
585 # Now we use basis linear algebra to find the coefficients,
586 # of the matrices-as-vectors-linear-combination, which should
587 # work for the original algebra basis too.
588 A = matrix.column(self.base_ring(), oper_vecs)
589
590 # We used the isometry on the left-hand side already, but we
591 # still need to do it for the right-hand side. Recall that we
592 # wanted something that summed to the identity matrix.
593 b = _mat2vec( matrix.identity(self.base_ring(), self.dimension()) )
594
595 # Now if there's an identity element in the algebra, this should work.
596 coeffs = A.solve_right(b)
597 return self.linear_combination(zip(self.gens(), coeffs))
598
599
600 def random_element(self):
601 # Temporary workaround for https://trac.sagemath.org/ticket/28327
602 if self.is_trivial():
603 return self.zero()
604 else:
605 s = super(FiniteDimensionalEuclideanJordanAlgebra, self)
606 return s.random_element()
607
608
609 def rank(self):
610 """
611 Return the rank of this EJA.
612
613 ALGORITHM:
614
615 The author knows of no algorithm to compute the rank of an EJA
616 where only the multiplication table is known. In lieu of one, we
617 require the rank to be specified when the algebra is created,
618 and simply pass along that number here.
619
620 SETUP::
621
622 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
623 ....: RealSymmetricEJA,
624 ....: ComplexHermitianEJA,
625 ....: QuaternionHermitianEJA,
626 ....: random_eja)
627
628 EXAMPLES:
629
630 The rank of the Jordan spin algebra is always two::
631
632 sage: JordanSpinEJA(2).rank()
633 2
634 sage: JordanSpinEJA(3).rank()
635 2
636 sage: JordanSpinEJA(4).rank()
637 2
638
639 The rank of the `n`-by-`n` Hermitian real, complex, or
640 quaternion matrices is `n`::
641
642 sage: RealSymmetricEJA(2).rank()
643 2
644 sage: ComplexHermitianEJA(2).rank()
645 2
646 sage: QuaternionHermitianEJA(2).rank()
647 2
648 sage: RealSymmetricEJA(5).rank()
649 5
650 sage: ComplexHermitianEJA(5).rank()
651 5
652 sage: QuaternionHermitianEJA(5).rank()
653 5
654
655 TESTS:
656
657 Ensure that every EJA that we know how to construct has a
658 positive integer rank::
659
660 sage: set_random_seed()
661 sage: r = random_eja().rank()
662 sage: r in ZZ and r > 0
663 True
664
665 """
666 return self._rank
667
668
669 def vector_space(self):
670 """
671 Return the vector space that underlies this algebra.
672
673 SETUP::
674
675 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
676
677 EXAMPLES::
678
679 sage: J = RealSymmetricEJA(2)
680 sage: J.vector_space()
681 Vector space of dimension 3 over...
682
683 """
684 return self.zero().to_vector().parent().ambient_vector_space()
685
686
687 Element = FiniteDimensionalEuclideanJordanAlgebraElement
688
689
690 class RealCartesianProductEJA(FiniteDimensionalEuclideanJordanAlgebra):
691 """
692 Return the Euclidean Jordan Algebra corresponding to the set
693 `R^n` under the Hadamard product.
694
695 Note: this is nothing more than the Cartesian product of ``n``
696 copies of the spin algebra. Once Cartesian product algebras
697 are implemented, this can go.
698
699 SETUP::
700
701 sage: from mjo.eja.eja_algebra import RealCartesianProductEJA
702
703 EXAMPLES:
704
705 This multiplication table can be verified by hand::
706
707 sage: J = RealCartesianProductEJA(3)
708 sage: e0,e1,e2 = J.gens()
709 sage: e0*e0
710 e0
711 sage: e0*e1
712 0
713 sage: e0*e2
714 0
715 sage: e1*e1
716 e1
717 sage: e1*e2
718 0
719 sage: e2*e2
720 e2
721
722 TESTS:
723
724 We can change the generator prefix::
725
726 sage: RealCartesianProductEJA(3, prefix='r').gens()
727 (r0, r1, r2)
728
729 Our inner product satisfies the Jordan axiom::
730
731 sage: set_random_seed()
732 sage: n = ZZ.random_element(1,5)
733 sage: J = RealCartesianProductEJA(n)
734 sage: x = J.random_element()
735 sage: y = J.random_element()
736 sage: z = J.random_element()
737 sage: (x*y).inner_product(z) == y.inner_product(x*z)
738 True
739
740 """
741 def __init__(self, n, field=QQ, **kwargs):
742 V = VectorSpace(field, n)
743 mult_table = [ [ V.gen(i)*(i == j) for j in range(n) ]
744 for i in range(n) ]
745
746 fdeja = super(RealCartesianProductEJA, self)
747 return fdeja.__init__(field, mult_table, rank=n, **kwargs)
748
749 def inner_product(self, x, y):
750 return _usual_ip(x,y)
751
752
753 def random_eja():
754 """
755 Return a "random" finite-dimensional Euclidean Jordan Algebra.
756
757 ALGORITHM:
758
759 For now, we choose a random natural number ``n`` (greater than zero)
760 and then give you back one of the following:
761
762 * The cartesian product of the rational numbers ``n`` times; this is
763 ``QQ^n`` with the Hadamard product.
764
765 * The Jordan spin algebra on ``QQ^n``.
766
767 * The ``n``-by-``n`` rational symmetric matrices with the symmetric
768 product.
769
770 * The ``n``-by-``n`` complex-rational Hermitian matrices embedded
771 in the space of ``2n``-by-``2n`` real symmetric matrices.
772
773 * The ``n``-by-``n`` quaternion-rational Hermitian matrices embedded
774 in the space of ``4n``-by-``4n`` real symmetric matrices.
775
776 Later this might be extended to return Cartesian products of the
777 EJAs above.
778
779 SETUP::
780
781 sage: from mjo.eja.eja_algebra import random_eja
782
783 TESTS::
784
785 sage: random_eja()
786 Euclidean Jordan algebra of dimension...
787
788 """
789
790 # The max_n component lets us choose different upper bounds on the
791 # value "n" that gets passed to the constructor. This is needed
792 # because e.g. R^{10} is reasonable to test, while the Hermitian
793 # 10-by-10 quaternion matrices are not.
794 (constructor, max_n) = choice([(RealCartesianProductEJA, 6),
795 (JordanSpinEJA, 6),
796 (RealSymmetricEJA, 5),
797 (ComplexHermitianEJA, 4),
798 (QuaternionHermitianEJA, 3)])
799 n = ZZ.random_element(1, max_n)
800 return constructor(n, field=QQ)
801
802
803
804 def _real_symmetric_basis(n, field):
805 """
806 Return a basis for the space of real symmetric n-by-n matrices.
807
808 SETUP::
809
810 sage: from mjo.eja.eja_algebra import _real_symmetric_basis
811
812 TESTS::
813
814 sage: set_random_seed()
815 sage: n = ZZ.random_element(1,5)
816 sage: B = _real_symmetric_basis(n, QQbar)
817 sage: all( M.is_symmetric() for M in B)
818 True
819
820 """
821 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
822 # coordinates.
823 S = []
824 for i in xrange(n):
825 for j in xrange(i+1):
826 Eij = matrix(field, n, lambda k,l: k==i and l==j)
827 if i == j:
828 Sij = Eij
829 else:
830 Sij = Eij + Eij.transpose()
831 # Now normalize it.
832 Sij = Sij / _real_symmetric_matrix_ip(Sij,Sij).sqrt()
833 S.append(Sij)
834 return tuple(S)
835
836
837 def _complex_hermitian_basis(n, field):
838 """
839 Returns a basis for the space of complex Hermitian n-by-n matrices.
840
841 SETUP::
842
843 sage: from mjo.eja.eja_algebra import _complex_hermitian_basis
844
845 TESTS::
846
847 sage: set_random_seed()
848 sage: n = ZZ.random_element(1,5)
849 sage: B = _complex_hermitian_basis(n, QQ)
850 sage: all( M.is_symmetric() for M in B)
851 True
852
853 """
854 F = QuadraticField(-1, 'I')
855 I = F.gen()
856
857 # This is like the symmetric case, but we need to be careful:
858 #
859 # * We want conjugate-symmetry, not just symmetry.
860 # * The diagonal will (as a result) be real.
861 #
862 S = []
863 for i in xrange(n):
864 for j in xrange(i+1):
865 Eij = matrix(field, n, lambda k,l: k==i and l==j)
866 if i == j:
867 Sij = _embed_complex_matrix(Eij)
868 S.append(Sij)
869 else:
870 # Beware, orthogonal but not normalized! The second one
871 # has a minus because it's conjugated.
872 Sij_real = _embed_complex_matrix(Eij + Eij.transpose())
873 S.append(Sij_real)
874 Sij_imag = _embed_complex_matrix(I*Eij - I*Eij.transpose())
875 S.append(Sij_imag)
876 return tuple(S)
877
878
879 def _quaternion_hermitian_basis(n, field):
880 """
881 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
882
883 SETUP::
884
885 sage: from mjo.eja.eja_algebra import _quaternion_hermitian_basis
886
887 TESTS::
888
889 sage: set_random_seed()
890 sage: n = ZZ.random_element(1,5)
891 sage: B = _quaternion_hermitian_basis(n, QQbar)
892 sage: all( M.is_symmetric() for M in B )
893 True
894
895 """
896 Q = QuaternionAlgebra(QQ,-1,-1)
897 I,J,K = Q.gens()
898
899 # This is like the symmetric case, but we need to be careful:
900 #
901 # * We want conjugate-symmetry, not just symmetry.
902 # * The diagonal will (as a result) be real.
903 #
904 S = []
905 for i in xrange(n):
906 for j in xrange(i+1):
907 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
908 if i == j:
909 Sij = _embed_quaternion_matrix(Eij)
910 S.append(Sij)
911 else:
912 # Beware, orthogonal but not normalized! The second,
913 # third, and fourth ones have a minus because they're
914 # conjugated.
915 Sij_real = _embed_quaternion_matrix(Eij + Eij.transpose())
916 S.append(Sij_real)
917 Sij_I = _embed_quaternion_matrix(I*Eij - I*Eij.transpose())
918 S.append(Sij_I)
919 Sij_J = _embed_quaternion_matrix(J*Eij - J*Eij.transpose())
920 S.append(Sij_J)
921 Sij_K = _embed_quaternion_matrix(K*Eij - K*Eij.transpose())
922 S.append(Sij_K)
923 return tuple(S)
924
925
926
927 def _multiplication_table_from_matrix_basis(basis):
928 """
929 At least three of the five simple Euclidean Jordan algebras have the
930 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
931 multiplication on the right is matrix multiplication. Given a basis
932 for the underlying matrix space, this function returns a
933 multiplication table (obtained by looping through the basis
934 elements) for an algebra of those matrices.
935 """
936 # In S^2, for example, we nominally have four coordinates even
937 # though the space is of dimension three only. The vector space V
938 # is supposed to hold the entire long vector, and the subspace W
939 # of V will be spanned by the vectors that arise from symmetric
940 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
941 field = basis[0].base_ring()
942 dimension = basis[0].nrows()
943
944 V = VectorSpace(field, dimension**2)
945 W = V.span_of_basis( _mat2vec(s) for s in basis )
946 n = len(basis)
947 mult_table = [[W.zero() for j in range(n)] for i in range(n)]
948 for i in range(n):
949 for j in range(n):
950 mat_entry = (basis[i]*basis[j] + basis[j]*basis[i])/2
951 mult_table[i][j] = W.coordinate_vector(_mat2vec(mat_entry))
952
953 return mult_table
954
955
956 def _embed_complex_matrix(M):
957 """
958 Embed the n-by-n complex matrix ``M`` into the space of real
959 matrices of size 2n-by-2n via the map the sends each entry `z = a +
960 bi` to the block matrix ``[[a,b],[-b,a]]``.
961
962 SETUP::
963
964 sage: from mjo.eja.eja_algebra import _embed_complex_matrix
965
966 EXAMPLES::
967
968 sage: F = QuadraticField(-1,'i')
969 sage: x1 = F(4 - 2*i)
970 sage: x2 = F(1 + 2*i)
971 sage: x3 = F(-i)
972 sage: x4 = F(6)
973 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
974 sage: _embed_complex_matrix(M)
975 [ 4 -2| 1 2]
976 [ 2 4|-2 1]
977 [-----+-----]
978 [ 0 -1| 6 0]
979 [ 1 0| 0 6]
980
981 TESTS:
982
983 Embedding is a homomorphism (isomorphism, in fact)::
984
985 sage: set_random_seed()
986 sage: n = ZZ.random_element(5)
987 sage: F = QuadraticField(-1, 'i')
988 sage: X = random_matrix(F, n)
989 sage: Y = random_matrix(F, n)
990 sage: actual = _embed_complex_matrix(X) * _embed_complex_matrix(Y)
991 sage: expected = _embed_complex_matrix(X*Y)
992 sage: actual == expected
993 True
994
995 """
996 n = M.nrows()
997 if M.ncols() != n:
998 raise ValueError("the matrix 'M' must be square")
999 field = M.base_ring()
1000 blocks = []
1001 for z in M.list():
1002 a = z.real()
1003 b = z.imag()
1004 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1005
1006 # We can drop the imaginaries here.
1007 return matrix.block(field.base_ring(), n, blocks)
1008
1009
1010 def _unembed_complex_matrix(M):
1011 """
1012 The inverse of _embed_complex_matrix().
1013
1014 SETUP::
1015
1016 sage: from mjo.eja.eja_algebra import (_embed_complex_matrix,
1017 ....: _unembed_complex_matrix)
1018
1019 EXAMPLES::
1020
1021 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1022 ....: [-2, 1, -4, 3],
1023 ....: [ 9, 10, 11, 12],
1024 ....: [-10, 9, -12, 11] ])
1025 sage: _unembed_complex_matrix(A)
1026 [ 2*i + 1 4*i + 3]
1027 [ 10*i + 9 12*i + 11]
1028
1029 TESTS:
1030
1031 Unembedding is the inverse of embedding::
1032
1033 sage: set_random_seed()
1034 sage: F = QuadraticField(-1, 'i')
1035 sage: M = random_matrix(F, 3)
1036 sage: _unembed_complex_matrix(_embed_complex_matrix(M)) == M
1037 True
1038
1039 """
1040 n = ZZ(M.nrows())
1041 if M.ncols() != n:
1042 raise ValueError("the matrix 'M' must be square")
1043 if not n.mod(2).is_zero():
1044 raise ValueError("the matrix 'M' must be a complex embedding")
1045
1046 F = QuadraticField(-1, 'i')
1047 i = F.gen()
1048
1049 # Go top-left to bottom-right (reading order), converting every
1050 # 2-by-2 block we see to a single complex element.
1051 elements = []
1052 for k in xrange(n/2):
1053 for j in xrange(n/2):
1054 submat = M[2*k:2*k+2,2*j:2*j+2]
1055 if submat[0,0] != submat[1,1]:
1056 raise ValueError('bad on-diagonal submatrix')
1057 if submat[0,1] != -submat[1,0]:
1058 raise ValueError('bad off-diagonal submatrix')
1059 z = submat[0,0] + submat[0,1]*i
1060 elements.append(z)
1061
1062 return matrix(F, n/2, elements)
1063
1064
1065 def _embed_quaternion_matrix(M):
1066 """
1067 Embed the n-by-n quaternion matrix ``M`` into the space of real
1068 matrices of size 4n-by-4n by first sending each quaternion entry
1069 `z = a + bi + cj + dk` to the block-complex matrix
1070 ``[[a + bi, c+di],[-c + di, a-bi]]`, and then embedding those into
1071 a real matrix.
1072
1073 SETUP::
1074
1075 sage: from mjo.eja.eja_algebra import _embed_quaternion_matrix
1076
1077 EXAMPLES::
1078
1079 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1080 sage: i,j,k = Q.gens()
1081 sage: x = 1 + 2*i + 3*j + 4*k
1082 sage: M = matrix(Q, 1, [[x]])
1083 sage: _embed_quaternion_matrix(M)
1084 [ 1 2 3 4]
1085 [-2 1 -4 3]
1086 [-3 4 1 -2]
1087 [-4 -3 2 1]
1088
1089 Embedding is a homomorphism (isomorphism, in fact)::
1090
1091 sage: set_random_seed()
1092 sage: n = ZZ.random_element(5)
1093 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1094 sage: X = random_matrix(Q, n)
1095 sage: Y = random_matrix(Q, n)
1096 sage: actual = _embed_quaternion_matrix(X)*_embed_quaternion_matrix(Y)
1097 sage: expected = _embed_quaternion_matrix(X*Y)
1098 sage: actual == expected
1099 True
1100
1101 """
1102 quaternions = M.base_ring()
1103 n = M.nrows()
1104 if M.ncols() != n:
1105 raise ValueError("the matrix 'M' must be square")
1106
1107 F = QuadraticField(-1, 'i')
1108 i = F.gen()
1109
1110 blocks = []
1111 for z in M.list():
1112 t = z.coefficient_tuple()
1113 a = t[0]
1114 b = t[1]
1115 c = t[2]
1116 d = t[3]
1117 cplx_matrix = matrix(F, 2, [[ a + b*i, c + d*i],
1118 [-c + d*i, a - b*i]])
1119 blocks.append(_embed_complex_matrix(cplx_matrix))
1120
1121 # We should have real entries by now, so use the realest field
1122 # we've got for the return value.
1123 return matrix.block(quaternions.base_ring(), n, blocks)
1124
1125
1126 def _unembed_quaternion_matrix(M):
1127 """
1128 The inverse of _embed_quaternion_matrix().
1129
1130 SETUP::
1131
1132 sage: from mjo.eja.eja_algebra import (_embed_quaternion_matrix,
1133 ....: _unembed_quaternion_matrix)
1134
1135 EXAMPLES::
1136
1137 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
1138 ....: [-2, 1, -4, 3],
1139 ....: [-3, 4, 1, -2],
1140 ....: [-4, -3, 2, 1]])
1141 sage: _unembed_quaternion_matrix(M)
1142 [1 + 2*i + 3*j + 4*k]
1143
1144 TESTS:
1145
1146 Unembedding is the inverse of embedding::
1147
1148 sage: set_random_seed()
1149 sage: Q = QuaternionAlgebra(QQ, -1, -1)
1150 sage: M = random_matrix(Q, 3)
1151 sage: _unembed_quaternion_matrix(_embed_quaternion_matrix(M)) == M
1152 True
1153
1154 """
1155 n = ZZ(M.nrows())
1156 if M.ncols() != n:
1157 raise ValueError("the matrix 'M' must be square")
1158 if not n.mod(4).is_zero():
1159 raise ValueError("the matrix 'M' must be a complex embedding")
1160
1161 Q = QuaternionAlgebra(QQ,-1,-1)
1162 i,j,k = Q.gens()
1163
1164 # Go top-left to bottom-right (reading order), converting every
1165 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
1166 # quaternion block.
1167 elements = []
1168 for l in xrange(n/4):
1169 for m in xrange(n/4):
1170 submat = _unembed_complex_matrix(M[4*l:4*l+4,4*m:4*m+4])
1171 if submat[0,0] != submat[1,1].conjugate():
1172 raise ValueError('bad on-diagonal submatrix')
1173 if submat[0,1] != -submat[1,0].conjugate():
1174 raise ValueError('bad off-diagonal submatrix')
1175 z = submat[0,0].real() + submat[0,0].imag()*i
1176 z += submat[0,1].real()*j + submat[0,1].imag()*k
1177 elements.append(z)
1178
1179 return matrix(Q, n/4, elements)
1180
1181
1182 # The usual inner product on R^n.
1183 def _usual_ip(x,y):
1184 return x.to_vector().inner_product(y.to_vector())
1185
1186 # The inner product used for the real symmetric simple EJA.
1187 # We keep it as a separate function because e.g. the complex
1188 # algebra uses the same inner product, except divided by 2.
1189 def _matrix_ip(X,Y):
1190 X_mat = X.natural_representation()
1191 Y_mat = Y.natural_representation()
1192 return (X_mat*Y_mat).trace()
1193
1194 def _real_symmetric_matrix_ip(X,Y):
1195 return (X*Y).trace()
1196
1197
1198 class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra):
1199 """
1200 The rank-n simple EJA consisting of real symmetric n-by-n
1201 matrices, the usual symmetric Jordan product, and the trace inner
1202 product. It has dimension `(n^2 + n)/2` over the reals.
1203
1204 SETUP::
1205
1206 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1207
1208 EXAMPLES::
1209
1210 sage: J = RealSymmetricEJA(2)
1211 sage: e0, e1, e2 = J.gens()
1212 sage: e0*e0
1213 e0
1214 sage: e1*e1
1215 1/2*e0 + 1/2*e2
1216 sage: e2*e2
1217 e2
1218
1219 TESTS:
1220
1221 The dimension of this algebra is `(n^2 + n) / 2`::
1222
1223 sage: set_random_seed()
1224 sage: n = ZZ.random_element(1,5)
1225 sage: J = RealSymmetricEJA(n)
1226 sage: J.dimension() == (n^2 + n)/2
1227 True
1228
1229 The Jordan multiplication is what we think it is::
1230
1231 sage: set_random_seed()
1232 sage: n = ZZ.random_element(1,5)
1233 sage: J = RealSymmetricEJA(n)
1234 sage: x = J.random_element()
1235 sage: y = J.random_element()
1236 sage: actual = (x*y).natural_representation()
1237 sage: X = x.natural_representation()
1238 sage: Y = y.natural_representation()
1239 sage: expected = (X*Y + Y*X)/2
1240 sage: actual == expected
1241 True
1242 sage: J(expected) == x*y
1243 True
1244
1245 We can change the generator prefix::
1246
1247 sage: RealSymmetricEJA(3, prefix='q').gens()
1248 (q0, q1, q2, q3, q4, q5)
1249
1250 Our inner product satisfies the Jordan axiom::
1251
1252 sage: set_random_seed()
1253 sage: n = ZZ.random_element(1,5)
1254 sage: J = RealSymmetricEJA(n)
1255 sage: x = J.random_element()
1256 sage: y = J.random_element()
1257 sage: z = J.random_element()
1258 sage: (x*y).inner_product(z) == y.inner_product(x*z)
1259 True
1260
1261 Our basis is normalized with respect to the natural inner product::
1262
1263 sage: set_random_seed()
1264 sage: n = ZZ.random_element(1,5)
1265 sage: J = RealSymmetricEJA(n)
1266 sage: all( b.norm() == 1 for b in J.gens() )
1267 True
1268
1269 Left-multiplication operators are symmetric because they satisfy
1270 the Jordan axiom::
1271
1272 sage: set_random_seed()
1273 sage: n = ZZ.random_element(1,5)
1274 sage: x = RealSymmetricEJA(n).random_element()
1275 sage: x.operator().matrix().is_symmetric()
1276 True
1277
1278 """
1279 def __init__(self, n, field=QQ, **kwargs):
1280 if n > 1 and field is QQ:
1281 # We'll need sqrt(2) to normalize the basis, and this
1282 # winds up in the multiplication table, so the whole
1283 # algebra needs to be over the field extension.
1284 field = QuadraticField(2, 'sqrt2')
1285
1286 S = _real_symmetric_basis(n, field)
1287 Qs = _multiplication_table_from_matrix_basis(S)
1288
1289 fdeja = super(RealSymmetricEJA, self)
1290 return fdeja.__init__(field,
1291 Qs,
1292 rank=n,
1293 natural_basis=S,
1294 **kwargs)
1295
1296 def inner_product(self, x, y):
1297 X = x.natural_representation()
1298 Y = y.natural_representation()
1299 return _real_symmetric_matrix_ip(X,Y)
1300
1301
1302 class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
1303 """
1304 The rank-n simple EJA consisting of complex Hermitian n-by-n
1305 matrices over the real numbers, the usual symmetric Jordan product,
1306 and the real-part-of-trace inner product. It has dimension `n^2` over
1307 the reals.
1308
1309 SETUP::
1310
1311 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1312
1313 TESTS:
1314
1315 The dimension of this algebra is `n^2`::
1316
1317 sage: set_random_seed()
1318 sage: n = ZZ.random_element(1,5)
1319 sage: J = ComplexHermitianEJA(n)
1320 sage: J.dimension() == n^2
1321 True
1322
1323 The Jordan multiplication is what we think it is::
1324
1325 sage: set_random_seed()
1326 sage: n = ZZ.random_element(1,5)
1327 sage: J = ComplexHermitianEJA(n)
1328 sage: x = J.random_element()
1329 sage: y = J.random_element()
1330 sage: actual = (x*y).natural_representation()
1331 sage: X = x.natural_representation()
1332 sage: Y = y.natural_representation()
1333 sage: expected = (X*Y + Y*X)/2
1334 sage: actual == expected
1335 True
1336 sage: J(expected) == x*y
1337 True
1338
1339 We can change the generator prefix::
1340
1341 sage: ComplexHermitianEJA(2, prefix='z').gens()
1342 (z0, z1, z2, z3)
1343
1344 Our inner product satisfies the Jordan axiom::
1345
1346 sage: set_random_seed()
1347 sage: n = ZZ.random_element(1,5)
1348 sage: J = ComplexHermitianEJA(n)
1349 sage: x = J.random_element()
1350 sage: y = J.random_element()
1351 sage: z = J.random_element()
1352 sage: (x*y).inner_product(z) == y.inner_product(x*z)
1353 True
1354
1355 """
1356 def __init__(self, n, field=QQ, **kwargs):
1357 S = _complex_hermitian_basis(n, field)
1358 Qs = _multiplication_table_from_matrix_basis(S)
1359
1360 fdeja = super(ComplexHermitianEJA, self)
1361 return fdeja.__init__(field,
1362 Qs,
1363 rank=n,
1364 natural_basis=S,
1365 **kwargs)
1366
1367
1368 def inner_product(self, x, y):
1369 # Since a+bi on the diagonal is represented as
1370 #
1371 # a + bi = [ a b ]
1372 # [ -b a ],
1373 #
1374 # we'll double-count the "a" entries if we take the trace of
1375 # the embedding.
1376 return _matrix_ip(x,y)/2
1377
1378
1379 class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
1380 """
1381 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
1382 matrices, the usual symmetric Jordan product, and the
1383 real-part-of-trace inner product. It has dimension `2n^2 - n` over
1384 the reals.
1385
1386 SETUP::
1387
1388 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1389
1390 TESTS:
1391
1392 The dimension of this algebra is `n^2`::
1393
1394 sage: set_random_seed()
1395 sage: n = ZZ.random_element(1,5)
1396 sage: J = QuaternionHermitianEJA(n)
1397 sage: J.dimension() == 2*(n^2) - n
1398 True
1399
1400 The Jordan multiplication is what we think it is::
1401
1402 sage: set_random_seed()
1403 sage: n = ZZ.random_element(1,5)
1404 sage: J = QuaternionHermitianEJA(n)
1405 sage: x = J.random_element()
1406 sage: y = J.random_element()
1407 sage: actual = (x*y).natural_representation()
1408 sage: X = x.natural_representation()
1409 sage: Y = y.natural_representation()
1410 sage: expected = (X*Y + Y*X)/2
1411 sage: actual == expected
1412 True
1413 sage: J(expected) == x*y
1414 True
1415
1416 We can change the generator prefix::
1417
1418 sage: QuaternionHermitianEJA(2, prefix='a').gens()
1419 (a0, a1, a2, a3, a4, a5)
1420
1421 Our inner product satisfies the Jordan axiom::
1422
1423 sage: set_random_seed()
1424 sage: n = ZZ.random_element(1,5)
1425 sage: J = QuaternionHermitianEJA(n)
1426 sage: x = J.random_element()
1427 sage: y = J.random_element()
1428 sage: z = J.random_element()
1429 sage: (x*y).inner_product(z) == y.inner_product(x*z)
1430 True
1431
1432 """
1433 def __init__(self, n, field=QQ, **kwargs):
1434 S = _quaternion_hermitian_basis(n, field)
1435 Qs = _multiplication_table_from_matrix_basis(S)
1436
1437 fdeja = super(QuaternionHermitianEJA, self)
1438 return fdeja.__init__(field,
1439 Qs,
1440 rank=n,
1441 natural_basis=S,
1442 **kwargs)
1443
1444 def inner_product(self, x, y):
1445 # Since a+bi+cj+dk on the diagonal is represented as
1446 #
1447 # a + bi +cj + dk = [ a b c d]
1448 # [ -b a -d c]
1449 # [ -c d a -b]
1450 # [ -d -c b a],
1451 #
1452 # we'll quadruple-count the "a" entries if we take the trace of
1453 # the embedding.
1454 return _matrix_ip(x,y)/4
1455
1456
1457 class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra):
1458 """
1459 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
1460 with the usual inner product and jordan product ``x*y =
1461 (<x_bar,y_bar>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
1462 the reals.
1463
1464 SETUP::
1465
1466 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1467
1468 EXAMPLES:
1469
1470 This multiplication table can be verified by hand::
1471
1472 sage: J = JordanSpinEJA(4)
1473 sage: e0,e1,e2,e3 = J.gens()
1474 sage: e0*e0
1475 e0
1476 sage: e0*e1
1477 e1
1478 sage: e0*e2
1479 e2
1480 sage: e0*e3
1481 e3
1482 sage: e1*e2
1483 0
1484 sage: e1*e3
1485 0
1486 sage: e2*e3
1487 0
1488
1489 We can change the generator prefix::
1490
1491 sage: JordanSpinEJA(2, prefix='B').gens()
1492 (B0, B1)
1493
1494 Our inner product satisfies the Jordan axiom::
1495
1496 sage: set_random_seed()
1497 sage: n = ZZ.random_element(1,5)
1498 sage: J = JordanSpinEJA(n)
1499 sage: x = J.random_element()
1500 sage: y = J.random_element()
1501 sage: z = J.random_element()
1502 sage: (x*y).inner_product(z) == y.inner_product(x*z)
1503 True
1504
1505 """
1506 def __init__(self, n, field=QQ, **kwargs):
1507 V = VectorSpace(field, n)
1508 mult_table = [[V.zero() for j in range(n)] for i in range(n)]
1509 for i in range(n):
1510 for j in range(n):
1511 x = V.gen(i)
1512 y = V.gen(j)
1513 x0 = x[0]
1514 xbar = x[1:]
1515 y0 = y[0]
1516 ybar = y[1:]
1517 # z = x*y
1518 z0 = x.inner_product(y)
1519 zbar = y0*xbar + x0*ybar
1520 z = V([z0] + zbar.list())
1521 mult_table[i][j] = z
1522
1523 # The rank of the spin algebra is two, unless we're in a
1524 # one-dimensional ambient space (because the rank is bounded by
1525 # the ambient dimension).
1526 fdeja = super(JordanSpinEJA, self)
1527 return fdeja.__init__(field, mult_table, rank=min(n,2), **kwargs)
1528
1529 def inner_product(self, x, y):
1530 return _usual_ip(x,y)