]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: add a TODO about factoring out the Jordan axiom / norm tests.
[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] [0 0]
506 [0 0], [1 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 Rational Field
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 # Beware, orthogonal but not normalized!
831 Sij = Eij + Eij.transpose()
832 S.append(Sij)
833 return tuple(S)
834
835
836 def _complex_hermitian_basis(n, field):
837 """
838 Returns a basis for the space of complex Hermitian n-by-n matrices.
839
840 SETUP::
841
842 sage: from mjo.eja.eja_algebra import _complex_hermitian_basis
843
844 TESTS::
845
846 sage: set_random_seed()
847 sage: n = ZZ.random_element(1,5)
848 sage: all( M.is_symmetric() for M in _complex_hermitian_basis(n) )
849 True
850
851 """
852 F = QuadraticField(-1, 'I')
853 I = F.gen()
854
855 # This is like the symmetric case, but we need to be careful:
856 #
857 # * We want conjugate-symmetry, not just symmetry.
858 # * The diagonal will (as a result) be real.
859 #
860 S = []
861 for i in xrange(n):
862 for j in xrange(i+1):
863 Eij = matrix(field, n, lambda k,l: k==i and l==j)
864 if i == j:
865 Sij = _embed_complex_matrix(Eij)
866 S.append(Sij)
867 else:
868 # Beware, orthogonal but not normalized! The second one
869 # has a minus because it's conjugated.
870 Sij_real = _embed_complex_matrix(Eij + Eij.transpose())
871 S.append(Sij_real)
872 Sij_imag = _embed_complex_matrix(I*Eij - I*Eij.transpose())
873 S.append(Sij_imag)
874 return tuple(S)
875
876
877 def _quaternion_hermitian_basis(n, field):
878 """
879 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
880
881 SETUP::
882
883 sage: from mjo.eja.eja_algebra import _quaternion_hermitian_basis
884
885 TESTS::
886
887 sage: set_random_seed()
888 sage: n = ZZ.random_element(1,5)
889 sage: all( M.is_symmetric() for M in _quaternion_hermitian_basis(n) )
890 True
891
892 """
893 Q = QuaternionAlgebra(QQ,-1,-1)
894 I,J,K = Q.gens()
895
896 # This is like the symmetric case, but we need to be careful:
897 #
898 # * We want conjugate-symmetry, not just symmetry.
899 # * The diagonal will (as a result) be real.
900 #
901 S = []
902 for i in xrange(n):
903 for j in xrange(i+1):
904 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
905 if i == j:
906 Sij = _embed_quaternion_matrix(Eij)
907 S.append(Sij)
908 else:
909 # Beware, orthogonal but not normalized! The second,
910 # third, and fourth ones have a minus because they're
911 # conjugated.
912 Sij_real = _embed_quaternion_matrix(Eij + Eij.transpose())
913 S.append(Sij_real)
914 Sij_I = _embed_quaternion_matrix(I*Eij - I*Eij.transpose())
915 S.append(Sij_I)
916 Sij_J = _embed_quaternion_matrix(J*Eij - J*Eij.transpose())
917 S.append(Sij_J)
918 Sij_K = _embed_quaternion_matrix(K*Eij - K*Eij.transpose())
919 S.append(Sij_K)
920 return tuple(S)
921
922
923
924 def _multiplication_table_from_matrix_basis(basis):
925 """
926 At least three of the five simple Euclidean Jordan algebras have the
927 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
928 multiplication on the right is matrix multiplication. Given a basis
929 for the underlying matrix space, this function returns a
930 multiplication table (obtained by looping through the basis
931 elements) for an algebra of those matrices.
932 """
933 # In S^2, for example, we nominally have four coordinates even
934 # though the space is of dimension three only. The vector space V
935 # is supposed to hold the entire long vector, and the subspace W
936 # of V will be spanned by the vectors that arise from symmetric
937 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
938 field = basis[0].base_ring()
939 dimension = basis[0].nrows()
940
941 V = VectorSpace(field, dimension**2)
942 W = V.span_of_basis( _mat2vec(s) for s in basis )
943 n = len(basis)
944 mult_table = [[W.zero() for j in range(n)] for i in range(n)]
945 for i in range(n):
946 for j in range(n):
947 mat_entry = (basis[i]*basis[j] + basis[j]*basis[i])/2
948 mult_table[i][j] = W.coordinate_vector(_mat2vec(mat_entry))
949
950 return mult_table
951
952
953 def _embed_complex_matrix(M):
954 """
955 Embed the n-by-n complex matrix ``M`` into the space of real
956 matrices of size 2n-by-2n via the map the sends each entry `z = a +
957 bi` to the block matrix ``[[a,b],[-b,a]]``.
958
959 SETUP::
960
961 sage: from mjo.eja.eja_algebra import _embed_complex_matrix
962
963 EXAMPLES::
964
965 sage: F = QuadraticField(-1,'i')
966 sage: x1 = F(4 - 2*i)
967 sage: x2 = F(1 + 2*i)
968 sage: x3 = F(-i)
969 sage: x4 = F(6)
970 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
971 sage: _embed_complex_matrix(M)
972 [ 4 -2| 1 2]
973 [ 2 4|-2 1]
974 [-----+-----]
975 [ 0 -1| 6 0]
976 [ 1 0| 0 6]
977
978 TESTS:
979
980 Embedding is a homomorphism (isomorphism, in fact)::
981
982 sage: set_random_seed()
983 sage: n = ZZ.random_element(5)
984 sage: F = QuadraticField(-1, 'i')
985 sage: X = random_matrix(F, n)
986 sage: Y = random_matrix(F, n)
987 sage: actual = _embed_complex_matrix(X) * _embed_complex_matrix(Y)
988 sage: expected = _embed_complex_matrix(X*Y)
989 sage: actual == expected
990 True
991
992 """
993 n = M.nrows()
994 if M.ncols() != n:
995 raise ValueError("the matrix 'M' must be square")
996 field = M.base_ring()
997 blocks = []
998 for z in M.list():
999 a = z.real()
1000 b = z.imag()
1001 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1002
1003 # We can drop the imaginaries here.
1004 return matrix.block(field.base_ring(), n, blocks)
1005
1006
1007 def _unembed_complex_matrix(M):
1008 """
1009 The inverse of _embed_complex_matrix().
1010
1011 SETUP::
1012
1013 sage: from mjo.eja.eja_algebra import (_embed_complex_matrix,
1014 ....: _unembed_complex_matrix)
1015
1016 EXAMPLES::
1017
1018 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1019 ....: [-2, 1, -4, 3],
1020 ....: [ 9, 10, 11, 12],
1021 ....: [-10, 9, -12, 11] ])
1022 sage: _unembed_complex_matrix(A)
1023 [ 2*i + 1 4*i + 3]
1024 [ 10*i + 9 12*i + 11]
1025
1026 TESTS:
1027
1028 Unembedding is the inverse of embedding::
1029
1030 sage: set_random_seed()
1031 sage: F = QuadraticField(-1, 'i')
1032 sage: M = random_matrix(F, 3)
1033 sage: _unembed_complex_matrix(_embed_complex_matrix(M)) == M
1034 True
1035
1036 """
1037 n = ZZ(M.nrows())
1038 if M.ncols() != n:
1039 raise ValueError("the matrix 'M' must be square")
1040 if not n.mod(2).is_zero():
1041 raise ValueError("the matrix 'M' must be a complex embedding")
1042
1043 F = QuadraticField(-1, 'i')
1044 i = F.gen()
1045
1046 # Go top-left to bottom-right (reading order), converting every
1047 # 2-by-2 block we see to a single complex element.
1048 elements = []
1049 for k in xrange(n/2):
1050 for j in xrange(n/2):
1051 submat = M[2*k:2*k+2,2*j:2*j+2]
1052 if submat[0,0] != submat[1,1]:
1053 raise ValueError('bad on-diagonal submatrix')
1054 if submat[0,1] != -submat[1,0]:
1055 raise ValueError('bad off-diagonal submatrix')
1056 z = submat[0,0] + submat[0,1]*i
1057 elements.append(z)
1058
1059 return matrix(F, n/2, elements)
1060
1061
1062 def _embed_quaternion_matrix(M):
1063 """
1064 Embed the n-by-n quaternion matrix ``M`` into the space of real
1065 matrices of size 4n-by-4n by first sending each quaternion entry
1066 `z = a + bi + cj + dk` to the block-complex matrix
1067 ``[[a + bi, c+di],[-c + di, a-bi]]`, and then embedding those into
1068 a real matrix.
1069
1070 SETUP::
1071
1072 sage: from mjo.eja.eja_algebra import _embed_quaternion_matrix
1073
1074 EXAMPLES::
1075
1076 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1077 sage: i,j,k = Q.gens()
1078 sage: x = 1 + 2*i + 3*j + 4*k
1079 sage: M = matrix(Q, 1, [[x]])
1080 sage: _embed_quaternion_matrix(M)
1081 [ 1 2 3 4]
1082 [-2 1 -4 3]
1083 [-3 4 1 -2]
1084 [-4 -3 2 1]
1085
1086 Embedding is a homomorphism (isomorphism, in fact)::
1087
1088 sage: set_random_seed()
1089 sage: n = ZZ.random_element(5)
1090 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1091 sage: X = random_matrix(Q, n)
1092 sage: Y = random_matrix(Q, n)
1093 sage: actual = _embed_quaternion_matrix(X)*_embed_quaternion_matrix(Y)
1094 sage: expected = _embed_quaternion_matrix(X*Y)
1095 sage: actual == expected
1096 True
1097
1098 """
1099 quaternions = M.base_ring()
1100 n = M.nrows()
1101 if M.ncols() != n:
1102 raise ValueError("the matrix 'M' must be square")
1103
1104 F = QuadraticField(-1, 'i')
1105 i = F.gen()
1106
1107 blocks = []
1108 for z in M.list():
1109 t = z.coefficient_tuple()
1110 a = t[0]
1111 b = t[1]
1112 c = t[2]
1113 d = t[3]
1114 cplx_matrix = matrix(F, 2, [[ a + b*i, c + d*i],
1115 [-c + d*i, a - b*i]])
1116 blocks.append(_embed_complex_matrix(cplx_matrix))
1117
1118 # We should have real entries by now, so use the realest field
1119 # we've got for the return value.
1120 return matrix.block(quaternions.base_ring(), n, blocks)
1121
1122
1123 def _unembed_quaternion_matrix(M):
1124 """
1125 The inverse of _embed_quaternion_matrix().
1126
1127 SETUP::
1128
1129 sage: from mjo.eja.eja_algebra import (_embed_quaternion_matrix,
1130 ....: _unembed_quaternion_matrix)
1131
1132 EXAMPLES::
1133
1134 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
1135 ....: [-2, 1, -4, 3],
1136 ....: [-3, 4, 1, -2],
1137 ....: [-4, -3, 2, 1]])
1138 sage: _unembed_quaternion_matrix(M)
1139 [1 + 2*i + 3*j + 4*k]
1140
1141 TESTS:
1142
1143 Unembedding is the inverse of embedding::
1144
1145 sage: set_random_seed()
1146 sage: Q = QuaternionAlgebra(QQ, -1, -1)
1147 sage: M = random_matrix(Q, 3)
1148 sage: _unembed_quaternion_matrix(_embed_quaternion_matrix(M)) == M
1149 True
1150
1151 """
1152 n = ZZ(M.nrows())
1153 if M.ncols() != n:
1154 raise ValueError("the matrix 'M' must be square")
1155 if not n.mod(4).is_zero():
1156 raise ValueError("the matrix 'M' must be a complex embedding")
1157
1158 Q = QuaternionAlgebra(QQ,-1,-1)
1159 i,j,k = Q.gens()
1160
1161 # Go top-left to bottom-right (reading order), converting every
1162 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
1163 # quaternion block.
1164 elements = []
1165 for l in xrange(n/4):
1166 for m in xrange(n/4):
1167 submat = _unembed_complex_matrix(M[4*l:4*l+4,4*m:4*m+4])
1168 if submat[0,0] != submat[1,1].conjugate():
1169 raise ValueError('bad on-diagonal submatrix')
1170 if submat[0,1] != -submat[1,0].conjugate():
1171 raise ValueError('bad off-diagonal submatrix')
1172 z = submat[0,0].real() + submat[0,0].imag()*i
1173 z += submat[0,1].real()*j + submat[0,1].imag()*k
1174 elements.append(z)
1175
1176 return matrix(Q, n/4, elements)
1177
1178
1179 # The usual inner product on R^n.
1180 def _usual_ip(x,y):
1181 return x.to_vector().inner_product(y.to_vector())
1182
1183 # The inner product used for the real symmetric simple EJA.
1184 # We keep it as a separate function because e.g. the complex
1185 # algebra uses the same inner product, except divided by 2.
1186 def _matrix_ip(X,Y):
1187 X_mat = X.natural_representation()
1188 Y_mat = Y.natural_representation()
1189 return (X_mat*Y_mat).trace()
1190
1191
1192 class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra):
1193 """
1194 The rank-n simple EJA consisting of real symmetric n-by-n
1195 matrices, the usual symmetric Jordan product, and the trace inner
1196 product. It has dimension `(n^2 + n)/2` over the reals.
1197
1198 SETUP::
1199
1200 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1201
1202 EXAMPLES::
1203
1204 sage: J = RealSymmetricEJA(2)
1205 sage: e0, e1, e2 = J.gens()
1206 sage: e0*e0
1207 e0
1208 sage: e1*e1
1209 e0 + e2
1210 sage: e2*e2
1211 e2
1212
1213 TESTS:
1214
1215 The dimension of this algebra is `(n^2 + n) / 2`::
1216
1217 sage: set_random_seed()
1218 sage: n = ZZ.random_element(1,5)
1219 sage: J = RealSymmetricEJA(n)
1220 sage: J.dimension() == (n^2 + n)/2
1221 True
1222
1223 The Jordan multiplication is what we think it is::
1224
1225 sage: set_random_seed()
1226 sage: n = ZZ.random_element(1,5)
1227 sage: J = RealSymmetricEJA(n)
1228 sage: x = J.random_element()
1229 sage: y = J.random_element()
1230 sage: actual = (x*y).natural_representation()
1231 sage: X = x.natural_representation()
1232 sage: Y = y.natural_representation()
1233 sage: expected = (X*Y + Y*X)/2
1234 sage: actual == expected
1235 True
1236 sage: J(expected) == x*y
1237 True
1238
1239 We can change the generator prefix::
1240
1241 sage: RealSymmetricEJA(3, prefix='q').gens()
1242 (q0, q1, q2, q3, q4, q5)
1243
1244 Our inner product satisfies the Jordan axiom::
1245
1246 sage: set_random_seed()
1247 sage: n = ZZ.random_element(1,5)
1248 sage: J = RealSymmetricEJA(n)
1249 sage: x = J.random_element()
1250 sage: y = J.random_element()
1251 sage: z = J.random_element()
1252 sage: (x*y).inner_product(z) == y.inner_product(x*z)
1253 True
1254
1255 """
1256 def __init__(self, n, field=QQ, **kwargs):
1257 S = _real_symmetric_basis(n, field)
1258 Qs = _multiplication_table_from_matrix_basis(S)
1259
1260 fdeja = super(RealSymmetricEJA, self)
1261 return fdeja.__init__(field,
1262 Qs,
1263 rank=n,
1264 natural_basis=S,
1265 **kwargs)
1266
1267 def inner_product(self, x, y):
1268 return _matrix_ip(x,y)
1269
1270
1271 class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
1272 """
1273 The rank-n simple EJA consisting of complex Hermitian n-by-n
1274 matrices over the real numbers, the usual symmetric Jordan product,
1275 and the real-part-of-trace inner product. It has dimension `n^2` over
1276 the reals.
1277
1278 SETUP::
1279
1280 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1281
1282 TESTS:
1283
1284 The dimension of this algebra is `n^2`::
1285
1286 sage: set_random_seed()
1287 sage: n = ZZ.random_element(1,5)
1288 sage: J = ComplexHermitianEJA(n)
1289 sage: J.dimension() == n^2
1290 True
1291
1292 The Jordan multiplication is what we think it is::
1293
1294 sage: set_random_seed()
1295 sage: n = ZZ.random_element(1,5)
1296 sage: J = ComplexHermitianEJA(n)
1297 sage: x = J.random_element()
1298 sage: y = J.random_element()
1299 sage: actual = (x*y).natural_representation()
1300 sage: X = x.natural_representation()
1301 sage: Y = y.natural_representation()
1302 sage: expected = (X*Y + Y*X)/2
1303 sage: actual == expected
1304 True
1305 sage: J(expected) == x*y
1306 True
1307
1308 We can change the generator prefix::
1309
1310 sage: ComplexHermitianEJA(2, prefix='z').gens()
1311 (z0, z1, z2, z3)
1312
1313 Our inner product satisfies the Jordan axiom::
1314
1315 sage: set_random_seed()
1316 sage: n = ZZ.random_element(1,5)
1317 sage: J = ComplexHermitianEJA(n)
1318 sage: x = J.random_element()
1319 sage: y = J.random_element()
1320 sage: z = J.random_element()
1321 sage: (x*y).inner_product(z) == y.inner_product(x*z)
1322 True
1323
1324 """
1325 def __init__(self, n, field=QQ, **kwargs):
1326 S = _complex_hermitian_basis(n, field)
1327 Qs = _multiplication_table_from_matrix_basis(S)
1328
1329 fdeja = super(ComplexHermitianEJA, self)
1330 return fdeja.__init__(field,
1331 Qs,
1332 rank=n,
1333 natural_basis=S,
1334 **kwargs)
1335
1336
1337 def inner_product(self, x, y):
1338 # Since a+bi on the diagonal is represented as
1339 #
1340 # a + bi = [ a b ]
1341 # [ -b a ],
1342 #
1343 # we'll double-count the "a" entries if we take the trace of
1344 # the embedding.
1345 return _matrix_ip(x,y)/2
1346
1347
1348 class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
1349 """
1350 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
1351 matrices, the usual symmetric Jordan product, and the
1352 real-part-of-trace inner product. It has dimension `2n^2 - n` over
1353 the reals.
1354
1355 SETUP::
1356
1357 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1358
1359 TESTS:
1360
1361 The dimension of this algebra is `n^2`::
1362
1363 sage: set_random_seed()
1364 sage: n = ZZ.random_element(1,5)
1365 sage: J = QuaternionHermitianEJA(n)
1366 sage: J.dimension() == 2*(n^2) - n
1367 True
1368
1369 The Jordan multiplication is what we think it is::
1370
1371 sage: set_random_seed()
1372 sage: n = ZZ.random_element(1,5)
1373 sage: J = QuaternionHermitianEJA(n)
1374 sage: x = J.random_element()
1375 sage: y = J.random_element()
1376 sage: actual = (x*y).natural_representation()
1377 sage: X = x.natural_representation()
1378 sage: Y = y.natural_representation()
1379 sage: expected = (X*Y + Y*X)/2
1380 sage: actual == expected
1381 True
1382 sage: J(expected) == x*y
1383 True
1384
1385 We can change the generator prefix::
1386
1387 sage: QuaternionHermitianEJA(2, prefix='a').gens()
1388 (a0, a1, a2, a3, a4, a5)
1389
1390 Our inner product satisfies the Jordan axiom::
1391
1392 sage: set_random_seed()
1393 sage: n = ZZ.random_element(1,5)
1394 sage: J = QuaternionHermitianEJA(n)
1395 sage: x = J.random_element()
1396 sage: y = J.random_element()
1397 sage: z = J.random_element()
1398 sage: (x*y).inner_product(z) == y.inner_product(x*z)
1399 True
1400
1401 """
1402 def __init__(self, n, field=QQ, **kwargs):
1403 S = _quaternion_hermitian_basis(n, field)
1404 Qs = _multiplication_table_from_matrix_basis(S)
1405
1406 fdeja = super(QuaternionHermitianEJA, self)
1407 return fdeja.__init__(field,
1408 Qs,
1409 rank=n,
1410 natural_basis=S,
1411 **kwargs)
1412
1413 def inner_product(self, x, y):
1414 # Since a+bi+cj+dk on the diagonal is represented as
1415 #
1416 # a + bi +cj + dk = [ a b c d]
1417 # [ -b a -d c]
1418 # [ -c d a -b]
1419 # [ -d -c b a],
1420 #
1421 # we'll quadruple-count the "a" entries if we take the trace of
1422 # the embedding.
1423 return _matrix_ip(x,y)/4
1424
1425
1426 class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra):
1427 """
1428 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
1429 with the usual inner product and jordan product ``x*y =
1430 (<x_bar,y_bar>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
1431 the reals.
1432
1433 SETUP::
1434
1435 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1436
1437 EXAMPLES:
1438
1439 This multiplication table can be verified by hand::
1440
1441 sage: J = JordanSpinEJA(4)
1442 sage: e0,e1,e2,e3 = J.gens()
1443 sage: e0*e0
1444 e0
1445 sage: e0*e1
1446 e1
1447 sage: e0*e2
1448 e2
1449 sage: e0*e3
1450 e3
1451 sage: e1*e2
1452 0
1453 sage: e1*e3
1454 0
1455 sage: e2*e3
1456 0
1457
1458 We can change the generator prefix::
1459
1460 sage: JordanSpinEJA(2, prefix='B').gens()
1461 (B0, B1)
1462
1463 Our inner product satisfies the Jordan axiom::
1464
1465 sage: set_random_seed()
1466 sage: n = ZZ.random_element(1,5)
1467 sage: J = JordanSpinEJA(n)
1468 sage: x = J.random_element()
1469 sage: y = J.random_element()
1470 sage: z = J.random_element()
1471 sage: (x*y).inner_product(z) == y.inner_product(x*z)
1472 True
1473
1474 """
1475 def __init__(self, n, field=QQ, **kwargs):
1476 V = VectorSpace(field, n)
1477 mult_table = [[V.zero() for j in range(n)] for i in range(n)]
1478 for i in range(n):
1479 for j in range(n):
1480 x = V.gen(i)
1481 y = V.gen(j)
1482 x0 = x[0]
1483 xbar = x[1:]
1484 y0 = y[0]
1485 ybar = y[1:]
1486 # z = x*y
1487 z0 = x.inner_product(y)
1488 zbar = y0*xbar + x0*ybar
1489 z = V([z0] + zbar.list())
1490 mult_table[i][j] = z
1491
1492 # The rank of the spin algebra is two, unless we're in a
1493 # one-dimensional ambient space (because the rank is bounded by
1494 # the ambient dimension).
1495 fdeja = super(JordanSpinEJA, self)
1496 return fdeja.__init__(field, mult_table, rank=min(n,2), **kwargs)
1497
1498 def inner_product(self, x, y):
1499 return _usual_ip(x,y)