]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: add voodoo _rank_computation() method for EJAs.
[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 itertools import repeat
9
10 from sage.algebras.quatalg.quaternion_algebra import QuaternionAlgebra
11 from sage.categories.magmatic_algebras import MagmaticAlgebras
12 from sage.combinat.free_module import CombinatorialFreeModule
13 from sage.matrix.constructor import matrix
14 from sage.matrix.matrix_space import MatrixSpace
15 from sage.misc.cachefunc import cached_method
16 from sage.misc.lazy_import import lazy_import
17 from sage.misc.prandom import choice
18 from sage.misc.table import table
19 from sage.modules.free_module import FreeModule, VectorSpace
20 from sage.rings.all import (ZZ, QQ, AA, QQbar, RR, RLF, CLF,
21 PolynomialRing,
22 QuadraticField)
23 from mjo.eja.eja_element import FiniteDimensionalEuclideanJordanAlgebraElement
24 lazy_import('mjo.eja.eja_subalgebra',
25 'FiniteDimensionalEuclideanJordanSubalgebra')
26 from mjo.eja.eja_utils import _mat2vec
27
28 class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule):
29
30 def _coerce_map_from_base_ring(self):
31 """
32 Disable the map from the base ring into the algebra.
33
34 Performing a nonsense conversion like this automatically
35 is counterpedagogical. The fallback is to try the usual
36 element constructor, which should also fail.
37
38 SETUP::
39
40 sage: from mjo.eja.eja_algebra import random_eja
41
42 TESTS::
43
44 sage: set_random_seed()
45 sage: J = random_eja()
46 sage: J(1)
47 Traceback (most recent call last):
48 ...
49 ValueError: not a naturally-represented algebra element
50
51 """
52 return None
53
54 def __init__(self,
55 field,
56 mult_table,
57 rank,
58 prefix='e',
59 category=None,
60 natural_basis=None,
61 check=True):
62 """
63 SETUP::
64
65 sage: from mjo.eja.eja_algebra import (JordanSpinEJA, random_eja)
66
67 EXAMPLES:
68
69 By definition, Jordan multiplication commutes::
70
71 sage: set_random_seed()
72 sage: J = random_eja()
73 sage: x,y = J.random_elements(2)
74 sage: x*y == y*x
75 True
76
77 TESTS:
78
79 The ``field`` we're given must be real::
80
81 sage: JordanSpinEJA(2,QQbar)
82 Traceback (most recent call last):
83 ...
84 ValueError: field is not real
85
86 """
87 if check:
88 if not field.is_subring(RR):
89 # Note: this does return true for the real algebraic
90 # field, and any quadratic field where we've specified
91 # a real embedding.
92 raise ValueError('field is not real')
93
94 self._rank = rank
95 self._natural_basis = natural_basis
96
97 if category is None:
98 category = MagmaticAlgebras(field).FiniteDimensional()
99 category = category.WithBasis().Unital()
100
101 fda = super(FiniteDimensionalEuclideanJordanAlgebra, self)
102 fda.__init__(field,
103 range(len(mult_table)),
104 prefix=prefix,
105 category=category)
106 self.print_options(bracket='')
107
108 # The multiplication table we're given is necessarily in terms
109 # of vectors, because we don't have an algebra yet for
110 # anything to be an element of. However, it's faster in the
111 # long run to have the multiplication table be in terms of
112 # algebra elements. We do this after calling the superclass
113 # constructor so that from_vector() knows what to do.
114 self._multiplication_table = [
115 list(map(lambda x: self.from_vector(x), ls))
116 for ls in mult_table
117 ]
118
119
120 def _element_constructor_(self, elt):
121 """
122 Construct an element of this algebra from its natural
123 representation.
124
125 This gets called only after the parent element _call_ method
126 fails to find a coercion for the argument.
127
128 SETUP::
129
130 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
131 ....: HadamardEJA,
132 ....: RealSymmetricEJA)
133
134 EXAMPLES:
135
136 The identity in `S^n` is converted to the identity in the EJA::
137
138 sage: J = RealSymmetricEJA(3)
139 sage: I = matrix.identity(QQ,3)
140 sage: J(I) == J.one()
141 True
142
143 This skew-symmetric matrix can't be represented in the EJA::
144
145 sage: J = RealSymmetricEJA(3)
146 sage: A = matrix(QQ,3, lambda i,j: i-j)
147 sage: J(A)
148 Traceback (most recent call last):
149 ...
150 ArithmeticError: vector is not in free module
151
152 TESTS:
153
154 Ensure that we can convert any element of the two non-matrix
155 simple algebras (whose natural representations are their usual
156 vector representations) back and forth faithfully::
157
158 sage: set_random_seed()
159 sage: J = HadamardEJA.random_instance()
160 sage: x = J.random_element()
161 sage: J(x.to_vector().column()) == x
162 True
163 sage: J = JordanSpinEJA.random_instance()
164 sage: x = J.random_element()
165 sage: J(x.to_vector().column()) == x
166 True
167
168 """
169 msg = "not a naturally-represented algebra element"
170 if elt == 0:
171 # The superclass implementation of random_element()
172 # needs to be able to coerce "0" into the algebra.
173 return self.zero()
174 elif elt in self.base_ring():
175 # Ensure that no base ring -> algebra coercion is performed
176 # by this method. There's some stupidity in sage that would
177 # otherwise propagate to this method; for example, sage thinks
178 # that the integer 3 belongs to the space of 2-by-2 matrices.
179 raise ValueError(msg)
180
181 natural_basis = self.natural_basis()
182 basis_space = natural_basis[0].matrix_space()
183 if elt not in basis_space:
184 raise ValueError(msg)
185
186 # Thanks for nothing! Matrix spaces aren't vector spaces in
187 # Sage, so we have to figure out its natural-basis coordinates
188 # ourselves. We use the basis space's ring instead of the
189 # element's ring because the basis space might be an algebraic
190 # closure whereas the base ring of the 3-by-3 identity matrix
191 # could be QQ instead of QQbar.
192 V = VectorSpace(basis_space.base_ring(), elt.nrows()*elt.ncols())
193 W = V.span_of_basis( _mat2vec(s) for s in natural_basis )
194 coords = W.coordinate_vector(_mat2vec(elt))
195 return self.from_vector(coords)
196
197
198 def _repr_(self):
199 """
200 Return a string representation of ``self``.
201
202 SETUP::
203
204 sage: from mjo.eja.eja_algebra import JordanSpinEJA
205
206 TESTS:
207
208 Ensure that it says what we think it says::
209
210 sage: JordanSpinEJA(2, field=AA)
211 Euclidean Jordan algebra of dimension 2 over Algebraic Real Field
212 sage: JordanSpinEJA(3, field=RDF)
213 Euclidean Jordan algebra of dimension 3 over Real Double Field
214
215 """
216 fmt = "Euclidean Jordan algebra of dimension {} over {}"
217 return fmt.format(self.dimension(), self.base_ring())
218
219 def product_on_basis(self, i, j):
220 return self._multiplication_table[i][j]
221
222 def _a_regular_element(self):
223 """
224 Guess a regular element. Needed to compute the basis for our
225 characteristic polynomial coefficients.
226
227 SETUP::
228
229 sage: from mjo.eja.eja_algebra import random_eja
230
231 TESTS:
232
233 Ensure that this hacky method succeeds for every algebra that we
234 know how to construct::
235
236 sage: set_random_seed()
237 sage: J = random_eja()
238 sage: J._a_regular_element().is_regular()
239 True
240
241 """
242 gs = self.gens()
243 z = self.sum( (i+1)*gs[i] for i in range(len(gs)) )
244 if not z.is_regular():
245 raise ValueError("don't know a regular element")
246 return z
247
248
249 @cached_method
250 def _charpoly_basis_space(self):
251 """
252 Return the vector space spanned by the basis used in our
253 characteristic polynomial coefficients. This is used not only to
254 compute those coefficients, but also any time we need to
255 evaluate the coefficients (like when we compute the trace or
256 determinant).
257 """
258 z = self._a_regular_element()
259 # Don't use the parent vector space directly here in case this
260 # happens to be a subalgebra. In that case, we would be e.g.
261 # two-dimensional but span_of_basis() would expect three
262 # coordinates.
263 V = VectorSpace(self.base_ring(), self.vector_space().dimension())
264 basis = [ (z**k).to_vector() for k in range(self.rank()) ]
265 V1 = V.span_of_basis( basis )
266 b = (V1.basis() + V1.complement().basis())
267 return V.span_of_basis(b)
268
269
270
271 @cached_method
272 def _charpoly_coeff(self, i):
273 """
274 Return the coefficient polynomial "a_{i}" of this algebra's
275 general characteristic polynomial.
276
277 Having this be a separate cached method lets us compute and
278 store the trace/determinant (a_{r-1} and a_{0} respectively)
279 separate from the entire characteristic polynomial.
280 """
281 (A_of_x, x, xr, detA) = self._charpoly_matrix_system()
282 R = A_of_x.base_ring()
283
284 if i == self.rank():
285 return R.one()
286 if i > self.rank():
287 # Guaranteed by theory
288 return R.zero()
289
290 # Danger: the in-place modification is done for performance
291 # reasons (reconstructing a matrix with huge polynomial
292 # entries is slow), but I don't know how cached_method works,
293 # so it's highly possible that we're modifying some global
294 # list variable by reference, here. In other words, you
295 # probably shouldn't call this method twice on the same
296 # algebra, at the same time, in two threads
297 Ai_orig = A_of_x.column(i)
298 A_of_x.set_column(i,xr)
299 numerator = A_of_x.det()
300 A_of_x.set_column(i,Ai_orig)
301
302 # We're relying on the theory here to ensure that each a_i is
303 # indeed back in R, and the added negative signs are to make
304 # the whole charpoly expression sum to zero.
305 return R(-numerator/detA)
306
307
308 @cached_method
309 def _charpoly_matrix_system(self):
310 """
311 Compute the matrix whose entries A_ij are polynomials in
312 X1,...,XN, the vector ``x`` of variables X1,...,XN, the vector
313 corresponding to `x^r` and the determinent of the matrix A =
314 [A_ij]. In other words, all of the fixed (cachable) data needed
315 to compute the coefficients of the characteristic polynomial.
316 """
317 r = self.rank()
318 n = self.dimension()
319
320 # Turn my vector space into a module so that "vectors" can
321 # have multivatiate polynomial entries.
322 names = tuple('X' + str(i) for i in range(1,n+1))
323 R = PolynomialRing(self.base_ring(), names)
324
325 # Using change_ring() on the parent's vector space doesn't work
326 # here because, in a subalgebra, that vector space has a basis
327 # and change_ring() tries to bring the basis along with it. And
328 # that doesn't work unless the new ring is a PID, which it usually
329 # won't be.
330 V = FreeModule(R,n)
331
332 # Now let x = (X1,X2,...,Xn) be the vector whose entries are
333 # indeterminates...
334 x = V(names)
335
336 # And figure out the "left multiplication by x" matrix in
337 # that setting.
338 lmbx_cols = []
339 monomial_matrices = [ self.monomial(i).operator().matrix()
340 for i in range(n) ] # don't recompute these!
341 for k in range(n):
342 ek = self.monomial(k).to_vector()
343 lmbx_cols.append(
344 sum( x[i]*(monomial_matrices[i]*ek)
345 for i in range(n) ) )
346 Lx = matrix.column(R, lmbx_cols)
347
348 # Now we can compute powers of x "symbolically"
349 x_powers = [self.one().to_vector(), x]
350 for d in range(2, r+1):
351 x_powers.append( Lx*(x_powers[-1]) )
352
353 idmat = matrix.identity(R, n)
354
355 W = self._charpoly_basis_space()
356 W = W.change_ring(R.fraction_field())
357
358 # Starting with the standard coordinates x = (X1,X2,...,Xn)
359 # and then converting the entries to W-coordinates allows us
360 # to pass in the standard coordinates to the charpoly and get
361 # back the right answer. Specifically, with x = (X1,X2,...,Xn),
362 # we have
363 #
364 # W.coordinates(x^2) eval'd at (standard z-coords)
365 # =
366 # W-coords of (z^2)
367 # =
368 # W-coords of (standard coords of x^2 eval'd at std-coords of z)
369 #
370 # We want the middle equivalent thing in our matrix, but use
371 # the first equivalent thing instead so that we can pass in
372 # standard coordinates.
373 x_powers = [ W.coordinate_vector(xp) for xp in x_powers ]
374 l2 = [idmat.column(k-1) for k in range(r+1, n+1)]
375 A_of_x = matrix.column(R, n, (x_powers[:r] + l2))
376 return (A_of_x, x, x_powers[r], A_of_x.det())
377
378
379 @cached_method
380 def characteristic_polynomial(self):
381 """
382 Return a characteristic polynomial that works for all elements
383 of this algebra.
384
385 The resulting polynomial has `n+1` variables, where `n` is the
386 dimension of this algebra. The first `n` variables correspond to
387 the coordinates of an algebra element: when evaluated at the
388 coordinates of an algebra element with respect to a certain
389 basis, the result is a univariate polynomial (in the one
390 remaining variable ``t``), namely the characteristic polynomial
391 of that element.
392
393 SETUP::
394
395 sage: from mjo.eja.eja_algebra import JordanSpinEJA, TrivialEJA
396
397 EXAMPLES:
398
399 The characteristic polynomial in the spin algebra is given in
400 Alizadeh, Example 11.11::
401
402 sage: J = JordanSpinEJA(3)
403 sage: p = J.characteristic_polynomial(); p
404 X1^2 - X2^2 - X3^2 + (-2*t)*X1 + t^2
405 sage: xvec = J.one().to_vector()
406 sage: p(*xvec)
407 t^2 - 2*t + 1
408
409 By definition, the characteristic polynomial is a monic
410 degree-zero polynomial in a rank-zero algebra. Note that
411 Cayley-Hamilton is indeed satisfied since the polynomial
412 ``1`` evaluates to the identity element of the algebra on
413 any argument::
414
415 sage: J = TrivialEJA()
416 sage: J.characteristic_polynomial()
417 1
418
419 """
420 r = self.rank()
421 n = self.dimension()
422
423 # The list of coefficient polynomials a_0, a_1, a_2, ..., a_n.
424 a = [ self._charpoly_coeff(i) for i in range(r+1) ]
425
426 # We go to a bit of trouble here to reorder the
427 # indeterminates, so that it's easier to evaluate the
428 # characteristic polynomial at x's coordinates and get back
429 # something in terms of t, which is what we want.
430 R = a[0].parent()
431 S = PolynomialRing(self.base_ring(),'t')
432 t = S.gen(0)
433 S = PolynomialRing(S, R.variable_names())
434 t = S(t)
435
436 return sum( a[k]*(t**k) for k in range(len(a)) )
437
438
439 def inner_product(self, x, y):
440 """
441 The inner product associated with this Euclidean Jordan algebra.
442
443 Defaults to the trace inner product, but can be overridden by
444 subclasses if they are sure that the necessary properties are
445 satisfied.
446
447 SETUP::
448
449 sage: from mjo.eja.eja_algebra import random_eja
450
451 EXAMPLES:
452
453 Our inner product is "associative," which means the following for
454 a symmetric bilinear form::
455
456 sage: set_random_seed()
457 sage: J = random_eja()
458 sage: x,y,z = J.random_elements(3)
459 sage: (x*y).inner_product(z) == y.inner_product(x*z)
460 True
461
462 """
463 X = x.natural_representation()
464 Y = y.natural_representation()
465 return self.natural_inner_product(X,Y)
466
467
468 def is_trivial(self):
469 """
470 Return whether or not this algebra is trivial.
471
472 A trivial algebra contains only the zero element.
473
474 SETUP::
475
476 sage: from mjo.eja.eja_algebra import (ComplexHermitianEJA,
477 ....: TrivialEJA)
478
479 EXAMPLES::
480
481 sage: J = ComplexHermitianEJA(3)
482 sage: J.is_trivial()
483 False
484
485 ::
486
487 sage: J = TrivialEJA()
488 sage: J.is_trivial()
489 True
490
491 """
492 return self.dimension() == 0
493
494
495 def multiplication_table(self):
496 """
497 Return a visual representation of this algebra's multiplication
498 table (on basis elements).
499
500 SETUP::
501
502 sage: from mjo.eja.eja_algebra import JordanSpinEJA
503
504 EXAMPLES::
505
506 sage: J = JordanSpinEJA(4)
507 sage: J.multiplication_table()
508 +----++----+----+----+----+
509 | * || e0 | e1 | e2 | e3 |
510 +====++====+====+====+====+
511 | e0 || e0 | e1 | e2 | e3 |
512 +----++----+----+----+----+
513 | e1 || e1 | e0 | 0 | 0 |
514 +----++----+----+----+----+
515 | e2 || e2 | 0 | e0 | 0 |
516 +----++----+----+----+----+
517 | e3 || e3 | 0 | 0 | e0 |
518 +----++----+----+----+----+
519
520 """
521 M = list(self._multiplication_table) # copy
522 for i in range(len(M)):
523 # M had better be "square"
524 M[i] = [self.monomial(i)] + M[i]
525 M = [["*"] + list(self.gens())] + M
526 return table(M, header_row=True, header_column=True, frame=True)
527
528
529 def natural_basis(self):
530 """
531 Return a more-natural representation of this algebra's basis.
532
533 Every finite-dimensional Euclidean Jordan Algebra is a direct
534 sum of five simple algebras, four of which comprise Hermitian
535 matrices. This method returns the original "natural" basis
536 for our underlying vector space. (Typically, the natural basis
537 is used to construct the multiplication table in the first place.)
538
539 Note that this will always return a matrix. The standard basis
540 in `R^n` will be returned as `n`-by-`1` column matrices.
541
542 SETUP::
543
544 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
545 ....: RealSymmetricEJA)
546
547 EXAMPLES::
548
549 sage: J = RealSymmetricEJA(2)
550 sage: J.basis()
551 Finite family {0: e0, 1: e1, 2: e2}
552 sage: J.natural_basis()
553 (
554 [1 0] [ 0 0.7071067811865475?] [0 0]
555 [0 0], [0.7071067811865475? 0], [0 1]
556 )
557
558 ::
559
560 sage: J = JordanSpinEJA(2)
561 sage: J.basis()
562 Finite family {0: e0, 1: e1}
563 sage: J.natural_basis()
564 (
565 [1] [0]
566 [0], [1]
567 )
568
569 """
570 if self._natural_basis is None:
571 M = self.natural_basis_space()
572 return tuple( M(b.to_vector()) for b in self.basis() )
573 else:
574 return self._natural_basis
575
576
577 def natural_basis_space(self):
578 """
579 Return the matrix space in which this algebra's natural basis
580 elements live.
581 """
582 if self._natural_basis is None or len(self._natural_basis) == 0:
583 return MatrixSpace(self.base_ring(), self.dimension(), 1)
584 else:
585 return self._natural_basis[0].matrix_space()
586
587
588 @staticmethod
589 def natural_inner_product(X,Y):
590 """
591 Compute the inner product of two naturally-represented elements.
592
593 For example in the real symmetric matrix EJA, this will compute
594 the trace inner-product of two n-by-n symmetric matrices. The
595 default should work for the real cartesian product EJA, the
596 Jordan spin EJA, and the real symmetric matrices. The others
597 will have to be overridden.
598 """
599 return (X.conjugate_transpose()*Y).trace()
600
601
602 @cached_method
603 def one(self):
604 """
605 Return the unit element of this algebra.
606
607 SETUP::
608
609 sage: from mjo.eja.eja_algebra import (HadamardEJA,
610 ....: random_eja)
611
612 EXAMPLES::
613
614 sage: J = HadamardEJA(5)
615 sage: J.one()
616 e0 + e1 + e2 + e3 + e4
617
618 TESTS:
619
620 The identity element acts like the identity::
621
622 sage: set_random_seed()
623 sage: J = random_eja()
624 sage: x = J.random_element()
625 sage: J.one()*x == x and x*J.one() == x
626 True
627
628 The matrix of the unit element's operator is the identity::
629
630 sage: set_random_seed()
631 sage: J = random_eja()
632 sage: actual = J.one().operator().matrix()
633 sage: expected = matrix.identity(J.base_ring(), J.dimension())
634 sage: actual == expected
635 True
636
637 """
638 # We can brute-force compute the matrices of the operators
639 # that correspond to the basis elements of this algebra.
640 # If some linear combination of those basis elements is the
641 # algebra identity, then the same linear combination of
642 # their matrices has to be the identity matrix.
643 #
644 # Of course, matrices aren't vectors in sage, so we have to
645 # appeal to the "long vectors" isometry.
646 oper_vecs = [ _mat2vec(g.operator().matrix()) for g in self.gens() ]
647
648 # Now we use basis linear algebra to find the coefficients,
649 # of the matrices-as-vectors-linear-combination, which should
650 # work for the original algebra basis too.
651 A = matrix.column(self.base_ring(), oper_vecs)
652
653 # We used the isometry on the left-hand side already, but we
654 # still need to do it for the right-hand side. Recall that we
655 # wanted something that summed to the identity matrix.
656 b = _mat2vec( matrix.identity(self.base_ring(), self.dimension()) )
657
658 # Now if there's an identity element in the algebra, this should work.
659 coeffs = A.solve_right(b)
660 return self.linear_combination(zip(self.gens(), coeffs))
661
662
663 def peirce_decomposition(self, c):
664 """
665 The Peirce decomposition of this algebra relative to the
666 idempotent ``c``.
667
668 In the future, this can be extended to a complete system of
669 orthogonal idempotents.
670
671 INPUT:
672
673 - ``c`` -- an idempotent of this algebra.
674
675 OUTPUT:
676
677 A triple (J0, J5, J1) containing two subalgebras and one subspace
678 of this algebra,
679
680 - ``J0`` -- the algebra on the eigenspace of ``c.operator()``
681 corresponding to the eigenvalue zero.
682
683 - ``J5`` -- the eigenspace (NOT a subalgebra) of ``c.operator()``
684 corresponding to the eigenvalue one-half.
685
686 - ``J1`` -- the algebra on the eigenspace of ``c.operator()``
687 corresponding to the eigenvalue one.
688
689 These are the only possible eigenspaces for that operator, and this
690 algebra is a direct sum of them. The spaces ``J0`` and ``J1`` are
691 orthogonal, and are subalgebras of this algebra with the appropriate
692 restrictions.
693
694 SETUP::
695
696 sage: from mjo.eja.eja_algebra import random_eja, RealSymmetricEJA
697
698 EXAMPLES:
699
700 The canonical example comes from the symmetric matrices, which
701 decompose into diagonal and off-diagonal parts::
702
703 sage: J = RealSymmetricEJA(3)
704 sage: C = matrix(QQ, [ [1,0,0],
705 ....: [0,1,0],
706 ....: [0,0,0] ])
707 sage: c = J(C)
708 sage: J0,J5,J1 = J.peirce_decomposition(c)
709 sage: J0
710 Euclidean Jordan algebra of dimension 1...
711 sage: J5
712 Vector space of degree 6 and dimension 2...
713 sage: J1
714 Euclidean Jordan algebra of dimension 3...
715
716 TESTS:
717
718 Every algebra decomposes trivially with respect to its identity
719 element::
720
721 sage: set_random_seed()
722 sage: J = random_eja()
723 sage: J0,J5,J1 = J.peirce_decomposition(J.one())
724 sage: J0.dimension() == 0 and J5.dimension() == 0
725 True
726 sage: J1.superalgebra() == J and J1.dimension() == J.dimension()
727 True
728
729 The identity elements in the two subalgebras are the
730 projections onto their respective subspaces of the
731 superalgebra's identity element::
732
733 sage: set_random_seed()
734 sage: J = random_eja()
735 sage: x = J.random_element()
736 sage: if not J.is_trivial():
737 ....: while x.is_nilpotent():
738 ....: x = J.random_element()
739 sage: c = x.subalgebra_idempotent()
740 sage: J0,J5,J1 = J.peirce_decomposition(c)
741 sage: J1(c) == J1.one()
742 True
743 sage: J0(J.one() - c) == J0.one()
744 True
745
746 """
747 if not c.is_idempotent():
748 raise ValueError("element is not idempotent: %s" % c)
749
750 # Default these to what they should be if they turn out to be
751 # trivial, because eigenspaces_left() won't return eigenvalues
752 # corresponding to trivial spaces (e.g. it returns only the
753 # eigenspace corresponding to lambda=1 if you take the
754 # decomposition relative to the identity element).
755 trivial = FiniteDimensionalEuclideanJordanSubalgebra(self, ())
756 J0 = trivial # eigenvalue zero
757 J5 = VectorSpace(self.base_ring(), 0) # eigenvalue one-half
758 J1 = trivial # eigenvalue one
759
760 for (eigval, eigspace) in c.operator().matrix().right_eigenspaces():
761 if eigval == ~(self.base_ring()(2)):
762 J5 = eigspace
763 else:
764 gens = tuple( self.from_vector(b) for b in eigspace.basis() )
765 subalg = FiniteDimensionalEuclideanJordanSubalgebra(self, gens)
766 if eigval == 0:
767 J0 = subalg
768 elif eigval == 1:
769 J1 = subalg
770 else:
771 raise ValueError("unexpected eigenvalue: %s" % eigval)
772
773 return (J0, J5, J1)
774
775
776 def random_elements(self, count):
777 """
778 Return ``count`` random elements as a tuple.
779
780 SETUP::
781
782 sage: from mjo.eja.eja_algebra import JordanSpinEJA
783
784 EXAMPLES::
785
786 sage: J = JordanSpinEJA(3)
787 sage: x,y,z = J.random_elements(3)
788 sage: all( [ x in J, y in J, z in J ])
789 True
790 sage: len( J.random_elements(10) ) == 10
791 True
792
793 """
794 return tuple( self.random_element() for idx in range(count) )
795
796
797 def _rank_computation(self):
798 r"""
799 Compute the rank of this algebra using highly suspicious voodoo.
800
801 ALGORITHM:
802
803 We first compute the basis representation of the operator L_x
804 using polynomial indeterminates are placeholders for the
805 coordinates of "x", which is arbitrary. We then use that
806 matrix to compute the (polynomial) entries of x^0, x^1, ...,
807 x^d,... for increasing values of "d", starting at zero. The
808 idea is that. If we also add "coefficient variables" a_0,
809 a_1,... to the ring, we can form the linear combination
810 a_0*x^0 + ... + a_d*x^d = 0, and ask what dimension the
811 solution space has as an affine variety. When "d" is smaller
812 than the rank, we expect that dimension to be the number of
813 coordinates of "x", since we can set *those* to whatever we
814 want, but linear independence forces the coefficients a_i to
815 be zero. Eventually, when "d" passes the rank, the dimension
816 of the solution space begins to grow, because we can *still*
817 set the coordinates of "x" arbitrarily, but now there are some
818 coefficients that make the sum zero as well. So, when the
819 dimension of the variety jumps, we return the corresponding
820 "d" as the rank of the algebra. This appears to work.
821
822 SETUP::
823
824 sage: from mjo.eja.eja_algebra import (HadamardEJA,
825 ....: JordanSpinEJA,
826 ....: RealSymmetricEJA,
827 ....: ComplexHermitianEJA,
828 ....: QuaternionHermitianEJA)
829
830 EXAMPLES::
831
832 sage: J = HadamardEJA(5)
833 sage: J._rank_computation() == J.rank()
834 True
835 sage: J = JordanSpinEJA(5)
836 sage: J._rank_computation() == J.rank()
837 True
838 sage: J = RealSymmetricEJA(4)
839 sage: J._rank_computation() == J.rank()
840 True
841 sage: J = ComplexHermitianEJA(3)
842 sage: J._rank_computation() == J.rank()
843 True
844 sage: J = QuaternionHermitianEJA(2)
845 sage: J._rank_computation() == J.rank()
846 True
847
848 """
849 n = self.dimension()
850 var_names = [ "X" + str(z) for z in range(1,n+1) ]
851 d = 0
852 ideal_dim = len(var_names)
853 def L_x_i_j(i,j):
854 # From a result in my book, these are the entries of the
855 # basis representation of L_x.
856 return sum( vars[d+k]*self.monomial(k).operator().matrix()[i,j]
857 for k in range(n) )
858
859 while ideal_dim == len(var_names):
860 coeff_names = [ "a" + str(z) for z in range(d) ]
861 R = PolynomialRing(self.base_ring(), coeff_names + var_names)
862 vars = R.gens()
863 L_x = matrix(R, n, n, L_x_i_j)
864 x_powers = [ vars[k]*(L_x**k)*self.one().to_vector()
865 for k in range(d) ]
866 eqs = [ sum(x_powers[k][j] for k in range(d)) for j in range(n) ]
867 ideal_dim = R.ideal(eqs).dimension()
868 d += 1
869
870 # Subtract one because we increment one too many times, and
871 # subtract another one because "d" is one greater than the
872 # answer anyway; when d=3, we go up to x^2.
873 return d-2
874
875 def rank(self):
876 """
877 Return the rank of this EJA.
878
879 ALGORITHM:
880
881 The author knows of no algorithm to compute the rank of an EJA
882 where only the multiplication table is known. In lieu of one, we
883 require the rank to be specified when the algebra is created,
884 and simply pass along that number here.
885
886 SETUP::
887
888 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
889 ....: RealSymmetricEJA,
890 ....: ComplexHermitianEJA,
891 ....: QuaternionHermitianEJA,
892 ....: random_eja)
893
894 EXAMPLES:
895
896 The rank of the Jordan spin algebra is always two::
897
898 sage: JordanSpinEJA(2).rank()
899 2
900 sage: JordanSpinEJA(3).rank()
901 2
902 sage: JordanSpinEJA(4).rank()
903 2
904
905 The rank of the `n`-by-`n` Hermitian real, complex, or
906 quaternion matrices is `n`::
907
908 sage: RealSymmetricEJA(4).rank()
909 4
910 sage: ComplexHermitianEJA(3).rank()
911 3
912 sage: QuaternionHermitianEJA(2).rank()
913 2
914
915 TESTS:
916
917 Ensure that every EJA that we know how to construct has a
918 positive integer rank, unless the algebra is trivial in
919 which case its rank will be zero::
920
921 sage: set_random_seed()
922 sage: J = random_eja()
923 sage: r = J.rank()
924 sage: r in ZZ
925 True
926 sage: r > 0 or (r == 0 and J.is_trivial())
927 True
928
929 """
930 return self._rank
931
932
933 def vector_space(self):
934 """
935 Return the vector space that underlies this algebra.
936
937 SETUP::
938
939 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
940
941 EXAMPLES::
942
943 sage: J = RealSymmetricEJA(2)
944 sage: J.vector_space()
945 Vector space of dimension 3 over...
946
947 """
948 return self.zero().to_vector().parent().ambient_vector_space()
949
950
951 Element = FiniteDimensionalEuclideanJordanAlgebraElement
952
953
954 class KnownRankEJA(object):
955 """
956 A class for algebras that we actually know we can construct. The
957 main issue is that, for most of our methods to make sense, we need
958 to know the rank of our algebra. Thus we can't simply generate a
959 "random" algebra, or even check that a given basis and product
960 satisfy the axioms; because even if everything looks OK, we wouldn't
961 know the rank we need to actuallty build the thing.
962
963 Not really a subclass of FDEJA because doing that causes method
964 resolution errors, e.g.
965
966 TypeError: Error when calling the metaclass bases
967 Cannot create a consistent method resolution
968 order (MRO) for bases FiniteDimensionalEuclideanJordanAlgebra,
969 KnownRankEJA
970
971 """
972 @staticmethod
973 def _max_test_case_size():
974 """
975 Return an integer "size" that is an upper bound on the size of
976 this algebra when it is used in a random test
977 case. Unfortunately, the term "size" is quite vague -- when
978 dealing with `R^n` under either the Hadamard or Jordan spin
979 product, the "size" refers to the dimension `n`. When dealing
980 with a matrix algebra (real symmetric or complex/quaternion
981 Hermitian), it refers to the size of the matrix, which is
982 far less than the dimension of the underlying vector space.
983
984 We default to five in this class, which is safe in `R^n`. The
985 matrix algebra subclasses (or any class where the "size" is
986 interpreted to be far less than the dimension) should override
987 with a smaller number.
988 """
989 return 5
990
991 @classmethod
992 def random_instance(cls, field=AA, **kwargs):
993 """
994 Return a random instance of this type of algebra.
995
996 Beware, this will crash for "most instances" because the
997 constructor below looks wrong.
998 """
999 if cls is TrivialEJA:
1000 # The TrivialEJA class doesn't take an "n" argument because
1001 # there's only one.
1002 return cls(field)
1003
1004 n = ZZ.random_element(cls._max_test_case_size()) + 1
1005 return cls(n, field, **kwargs)
1006
1007
1008 class HadamardEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
1009 """
1010 Return the Euclidean Jordan Algebra corresponding to the set
1011 `R^n` under the Hadamard product.
1012
1013 Note: this is nothing more than the Cartesian product of ``n``
1014 copies of the spin algebra. Once Cartesian product algebras
1015 are implemented, this can go.
1016
1017 SETUP::
1018
1019 sage: from mjo.eja.eja_algebra import HadamardEJA
1020
1021 EXAMPLES:
1022
1023 This multiplication table can be verified by hand::
1024
1025 sage: J = HadamardEJA(3)
1026 sage: e0,e1,e2 = J.gens()
1027 sage: e0*e0
1028 e0
1029 sage: e0*e1
1030 0
1031 sage: e0*e2
1032 0
1033 sage: e1*e1
1034 e1
1035 sage: e1*e2
1036 0
1037 sage: e2*e2
1038 e2
1039
1040 TESTS:
1041
1042 We can change the generator prefix::
1043
1044 sage: HadamardEJA(3, prefix='r').gens()
1045 (r0, r1, r2)
1046
1047 """
1048 def __init__(self, n, field=AA, **kwargs):
1049 V = VectorSpace(field, n)
1050 mult_table = [ [ V.gen(i)*(i == j) for j in range(n) ]
1051 for i in range(n) ]
1052
1053 fdeja = super(HadamardEJA, self)
1054 return fdeja.__init__(field, mult_table, rank=n, **kwargs)
1055
1056 def inner_product(self, x, y):
1057 """
1058 Faster to reimplement than to use natural representations.
1059
1060 SETUP::
1061
1062 sage: from mjo.eja.eja_algebra import HadamardEJA
1063
1064 TESTS:
1065
1066 Ensure that this is the usual inner product for the algebras
1067 over `R^n`::
1068
1069 sage: set_random_seed()
1070 sage: J = HadamardEJA.random_instance()
1071 sage: x,y = J.random_elements(2)
1072 sage: X = x.natural_representation()
1073 sage: Y = y.natural_representation()
1074 sage: x.inner_product(y) == J.natural_inner_product(X,Y)
1075 True
1076
1077 """
1078 return x.to_vector().inner_product(y.to_vector())
1079
1080
1081 def random_eja(field=AA, nontrivial=False):
1082 """
1083 Return a "random" finite-dimensional Euclidean Jordan Algebra.
1084
1085 SETUP::
1086
1087 sage: from mjo.eja.eja_algebra import random_eja
1088
1089 TESTS::
1090
1091 sage: random_eja()
1092 Euclidean Jordan algebra of dimension...
1093
1094 """
1095 eja_classes = KnownRankEJA.__subclasses__()
1096 if nontrivial:
1097 eja_classes.remove(TrivialEJA)
1098 classname = choice(eja_classes)
1099 return classname.random_instance(field=field)
1100
1101
1102
1103
1104
1105
1106 class MatrixEuclideanJordanAlgebra(FiniteDimensionalEuclideanJordanAlgebra):
1107 @staticmethod
1108 def _max_test_case_size():
1109 # Play it safe, since this will be squared and the underlying
1110 # field can have dimension 4 (quaternions) too.
1111 return 2
1112
1113 def __init__(self, field, basis, rank, normalize_basis=True, **kwargs):
1114 """
1115 Compared to the superclass constructor, we take a basis instead of
1116 a multiplication table because the latter can be computed in terms
1117 of the former when the product is known (like it is here).
1118 """
1119 # Used in this class's fast _charpoly_coeff() override.
1120 self._basis_normalizers = None
1121
1122 # We're going to loop through this a few times, so now's a good
1123 # time to ensure that it isn't a generator expression.
1124 basis = tuple(basis)
1125
1126 if rank > 1 and normalize_basis:
1127 # We'll need sqrt(2) to normalize the basis, and this
1128 # winds up in the multiplication table, so the whole
1129 # algebra needs to be over the field extension.
1130 R = PolynomialRing(field, 'z')
1131 z = R.gen()
1132 p = z**2 - 2
1133 if p.is_irreducible():
1134 field = field.extension(p, 'sqrt2', embedding=RLF(2).sqrt())
1135 basis = tuple( s.change_ring(field) for s in basis )
1136 self._basis_normalizers = tuple(
1137 ~(self.natural_inner_product(s,s).sqrt()) for s in basis )
1138 basis = tuple(s*c for (s,c) in zip(basis,self._basis_normalizers))
1139
1140 Qs = self.multiplication_table_from_matrix_basis(basis)
1141
1142 fdeja = super(MatrixEuclideanJordanAlgebra, self)
1143 return fdeja.__init__(field,
1144 Qs,
1145 rank=rank,
1146 natural_basis=basis,
1147 **kwargs)
1148
1149
1150 @cached_method
1151 def _charpoly_coeff(self, i):
1152 """
1153 Override the parent method with something that tries to compute
1154 over a faster (non-extension) field.
1155 """
1156 if self._basis_normalizers is None:
1157 # We didn't normalize, so assume that the basis we started
1158 # with had entries in a nice field.
1159 return super(MatrixEuclideanJordanAlgebra, self)._charpoly_coeff(i)
1160 else:
1161 basis = ( (b/n) for (b,n) in zip(self.natural_basis(),
1162 self._basis_normalizers) )
1163
1164 # Do this over the rationals and convert back at the end.
1165 J = MatrixEuclideanJordanAlgebra(QQ,
1166 basis,
1167 self.rank(),
1168 normalize_basis=False)
1169 (_,x,_,_) = J._charpoly_matrix_system()
1170 p = J._charpoly_coeff(i)
1171 # p might be missing some vars, have to substitute "optionally"
1172 pairs = zip(x.base_ring().gens(), self._basis_normalizers)
1173 substitutions = { v: v*c for (v,c) in pairs }
1174 result = p.subs(substitutions)
1175
1176 # The result of "subs" can be either a coefficient-ring
1177 # element or a polynomial. Gotta handle both cases.
1178 if result in QQ:
1179 return self.base_ring()(result)
1180 else:
1181 return result.change_ring(self.base_ring())
1182
1183
1184 @staticmethod
1185 def multiplication_table_from_matrix_basis(basis):
1186 """
1187 At least three of the five simple Euclidean Jordan algebras have the
1188 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
1189 multiplication on the right is matrix multiplication. Given a basis
1190 for the underlying matrix space, this function returns a
1191 multiplication table (obtained by looping through the basis
1192 elements) for an algebra of those matrices.
1193 """
1194 # In S^2, for example, we nominally have four coordinates even
1195 # though the space is of dimension three only. The vector space V
1196 # is supposed to hold the entire long vector, and the subspace W
1197 # of V will be spanned by the vectors that arise from symmetric
1198 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
1199 field = basis[0].base_ring()
1200 dimension = basis[0].nrows()
1201
1202 V = VectorSpace(field, dimension**2)
1203 W = V.span_of_basis( _mat2vec(s) for s in basis )
1204 n = len(basis)
1205 mult_table = [[W.zero() for j in range(n)] for i in range(n)]
1206 for i in range(n):
1207 for j in range(n):
1208 mat_entry = (basis[i]*basis[j] + basis[j]*basis[i])/2
1209 mult_table[i][j] = W.coordinate_vector(_mat2vec(mat_entry))
1210
1211 return mult_table
1212
1213
1214 @staticmethod
1215 def real_embed(M):
1216 """
1217 Embed the matrix ``M`` into a space of real matrices.
1218
1219 The matrix ``M`` can have entries in any field at the moment:
1220 the real numbers, complex numbers, or quaternions. And although
1221 they are not a field, we can probably support octonions at some
1222 point, too. This function returns a real matrix that "acts like"
1223 the original with respect to matrix multiplication; i.e.
1224
1225 real_embed(M*N) = real_embed(M)*real_embed(N)
1226
1227 """
1228 raise NotImplementedError
1229
1230
1231 @staticmethod
1232 def real_unembed(M):
1233 """
1234 The inverse of :meth:`real_embed`.
1235 """
1236 raise NotImplementedError
1237
1238
1239 @classmethod
1240 def natural_inner_product(cls,X,Y):
1241 Xu = cls.real_unembed(X)
1242 Yu = cls.real_unembed(Y)
1243 tr = (Xu*Yu).trace()
1244
1245 if tr in RLF:
1246 # It's real already.
1247 return tr
1248
1249 # Otherwise, try the thing that works for complex numbers; and
1250 # if that doesn't work, the thing that works for quaternions.
1251 try:
1252 return tr.vector()[0] # real part, imag part is index 1
1253 except AttributeError:
1254 # A quaternions doesn't have a vector() method, but does
1255 # have coefficient_tuple() method that returns the
1256 # coefficients of 1, i, j, and k -- in that order.
1257 return tr.coefficient_tuple()[0]
1258
1259
1260 class RealMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1261 @staticmethod
1262 def real_embed(M):
1263 """
1264 The identity function, for embedding real matrices into real
1265 matrices.
1266 """
1267 return M
1268
1269 @staticmethod
1270 def real_unembed(M):
1271 """
1272 The identity function, for unembedding real matrices from real
1273 matrices.
1274 """
1275 return M
1276
1277
1278 class RealSymmetricEJA(RealMatrixEuclideanJordanAlgebra, KnownRankEJA):
1279 """
1280 The rank-n simple EJA consisting of real symmetric n-by-n
1281 matrices, the usual symmetric Jordan product, and the trace inner
1282 product. It has dimension `(n^2 + n)/2` over the reals.
1283
1284 SETUP::
1285
1286 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1287
1288 EXAMPLES::
1289
1290 sage: J = RealSymmetricEJA(2)
1291 sage: e0, e1, e2 = J.gens()
1292 sage: e0*e0
1293 e0
1294 sage: e1*e1
1295 1/2*e0 + 1/2*e2
1296 sage: e2*e2
1297 e2
1298
1299 In theory, our "field" can be any subfield of the reals::
1300
1301 sage: RealSymmetricEJA(2, RDF)
1302 Euclidean Jordan algebra of dimension 3 over Real Double Field
1303 sage: RealSymmetricEJA(2, RR)
1304 Euclidean Jordan algebra of dimension 3 over Real Field with
1305 53 bits of precision
1306
1307 TESTS:
1308
1309 The dimension of this algebra is `(n^2 + n) / 2`::
1310
1311 sage: set_random_seed()
1312 sage: n_max = RealSymmetricEJA._max_test_case_size()
1313 sage: n = ZZ.random_element(1, n_max)
1314 sage: J = RealSymmetricEJA(n)
1315 sage: J.dimension() == (n^2 + n)/2
1316 True
1317
1318 The Jordan multiplication is what we think it is::
1319
1320 sage: set_random_seed()
1321 sage: J = RealSymmetricEJA.random_instance()
1322 sage: x,y = J.random_elements(2)
1323 sage: actual = (x*y).natural_representation()
1324 sage: X = x.natural_representation()
1325 sage: Y = y.natural_representation()
1326 sage: expected = (X*Y + Y*X)/2
1327 sage: actual == expected
1328 True
1329 sage: J(expected) == x*y
1330 True
1331
1332 We can change the generator prefix::
1333
1334 sage: RealSymmetricEJA(3, prefix='q').gens()
1335 (q0, q1, q2, q3, q4, q5)
1336
1337 Our natural basis is normalized with respect to the natural inner
1338 product unless we specify otherwise::
1339
1340 sage: set_random_seed()
1341 sage: J = RealSymmetricEJA.random_instance()
1342 sage: all( b.norm() == 1 for b in J.gens() )
1343 True
1344
1345 Since our natural basis is normalized with respect to the natural
1346 inner product, and since we know that this algebra is an EJA, any
1347 left-multiplication operator's matrix will be symmetric because
1348 natural->EJA basis representation is an isometry and within the EJA
1349 the operator is self-adjoint by the Jordan axiom::
1350
1351 sage: set_random_seed()
1352 sage: x = RealSymmetricEJA.random_instance().random_element()
1353 sage: x.operator().matrix().is_symmetric()
1354 True
1355
1356 """
1357 @classmethod
1358 def _denormalized_basis(cls, n, field):
1359 """
1360 Return a basis for the space of real symmetric n-by-n matrices.
1361
1362 SETUP::
1363
1364 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1365
1366 TESTS::
1367
1368 sage: set_random_seed()
1369 sage: n = ZZ.random_element(1,5)
1370 sage: B = RealSymmetricEJA._denormalized_basis(n,QQ)
1371 sage: all( M.is_symmetric() for M in B)
1372 True
1373
1374 """
1375 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1376 # coordinates.
1377 S = []
1378 for i in range(n):
1379 for j in range(i+1):
1380 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1381 if i == j:
1382 Sij = Eij
1383 else:
1384 Sij = Eij + Eij.transpose()
1385 S.append(Sij)
1386 return S
1387
1388
1389 @staticmethod
1390 def _max_test_case_size():
1391 return 4 # Dimension 10
1392
1393
1394 def __init__(self, n, field=AA, **kwargs):
1395 basis = self._denormalized_basis(n, field)
1396 super(RealSymmetricEJA, self).__init__(field, basis, n, **kwargs)
1397
1398
1399 class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1400 @staticmethod
1401 def real_embed(M):
1402 """
1403 Embed the n-by-n complex matrix ``M`` into the space of real
1404 matrices of size 2n-by-2n via the map the sends each entry `z = a +
1405 bi` to the block matrix ``[[a,b],[-b,a]]``.
1406
1407 SETUP::
1408
1409 sage: from mjo.eja.eja_algebra import \
1410 ....: ComplexMatrixEuclideanJordanAlgebra
1411
1412 EXAMPLES::
1413
1414 sage: F = QuadraticField(-1, 'I')
1415 sage: x1 = F(4 - 2*i)
1416 sage: x2 = F(1 + 2*i)
1417 sage: x3 = F(-i)
1418 sage: x4 = F(6)
1419 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
1420 sage: ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1421 [ 4 -2| 1 2]
1422 [ 2 4|-2 1]
1423 [-----+-----]
1424 [ 0 -1| 6 0]
1425 [ 1 0| 0 6]
1426
1427 TESTS:
1428
1429 Embedding is a homomorphism (isomorphism, in fact)::
1430
1431 sage: set_random_seed()
1432 sage: n_max = ComplexMatrixEuclideanJordanAlgebra._max_test_case_size()
1433 sage: n = ZZ.random_element(n_max)
1434 sage: F = QuadraticField(-1, 'I')
1435 sage: X = random_matrix(F, n)
1436 sage: Y = random_matrix(F, n)
1437 sage: Xe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X)
1438 sage: Ye = ComplexMatrixEuclideanJordanAlgebra.real_embed(Y)
1439 sage: XYe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1440 sage: Xe*Ye == XYe
1441 True
1442
1443 """
1444 n = M.nrows()
1445 if M.ncols() != n:
1446 raise ValueError("the matrix 'M' must be square")
1447
1448 # We don't need any adjoined elements...
1449 field = M.base_ring().base_ring()
1450
1451 blocks = []
1452 for z in M.list():
1453 a = z.list()[0] # real part, I guess
1454 b = z.list()[1] # imag part, I guess
1455 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1456
1457 return matrix.block(field, n, blocks)
1458
1459
1460 @staticmethod
1461 def real_unembed(M):
1462 """
1463 The inverse of _embed_complex_matrix().
1464
1465 SETUP::
1466
1467 sage: from mjo.eja.eja_algebra import \
1468 ....: ComplexMatrixEuclideanJordanAlgebra
1469
1470 EXAMPLES::
1471
1472 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1473 ....: [-2, 1, -4, 3],
1474 ....: [ 9, 10, 11, 12],
1475 ....: [-10, 9, -12, 11] ])
1476 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(A)
1477 [ 2*I + 1 4*I + 3]
1478 [ 10*I + 9 12*I + 11]
1479
1480 TESTS:
1481
1482 Unembedding is the inverse of embedding::
1483
1484 sage: set_random_seed()
1485 sage: F = QuadraticField(-1, 'I')
1486 sage: M = random_matrix(F, 3)
1487 sage: Me = ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1488 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1489 True
1490
1491 """
1492 n = ZZ(M.nrows())
1493 if M.ncols() != n:
1494 raise ValueError("the matrix 'M' must be square")
1495 if not n.mod(2).is_zero():
1496 raise ValueError("the matrix 'M' must be a complex embedding")
1497
1498 # If "M" was normalized, its base ring might have roots
1499 # adjoined and they can stick around after unembedding.
1500 field = M.base_ring()
1501 R = PolynomialRing(field, 'z')
1502 z = R.gen()
1503 if field is AA:
1504 # Sage doesn't know how to embed AA into QQbar, i.e. how
1505 # to adjoin sqrt(-1) to AA.
1506 F = QQbar
1507 else:
1508 F = field.extension(z**2 + 1, 'I', embedding=CLF(-1).sqrt())
1509 i = F.gen()
1510
1511 # Go top-left to bottom-right (reading order), converting every
1512 # 2-by-2 block we see to a single complex element.
1513 elements = []
1514 for k in range(n/2):
1515 for j in range(n/2):
1516 submat = M[2*k:2*k+2,2*j:2*j+2]
1517 if submat[0,0] != submat[1,1]:
1518 raise ValueError('bad on-diagonal submatrix')
1519 if submat[0,1] != -submat[1,0]:
1520 raise ValueError('bad off-diagonal submatrix')
1521 z = submat[0,0] + submat[0,1]*i
1522 elements.append(z)
1523
1524 return matrix(F, n/2, elements)
1525
1526
1527 @classmethod
1528 def natural_inner_product(cls,X,Y):
1529 """
1530 Compute a natural inner product in this algebra directly from
1531 its real embedding.
1532
1533 SETUP::
1534
1535 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1536
1537 TESTS:
1538
1539 This gives the same answer as the slow, default method implemented
1540 in :class:`MatrixEuclideanJordanAlgebra`::
1541
1542 sage: set_random_seed()
1543 sage: J = ComplexHermitianEJA.random_instance()
1544 sage: x,y = J.random_elements(2)
1545 sage: Xe = x.natural_representation()
1546 sage: Ye = y.natural_representation()
1547 sage: X = ComplexHermitianEJA.real_unembed(Xe)
1548 sage: Y = ComplexHermitianEJA.real_unembed(Ye)
1549 sage: expected = (X*Y).trace().real()
1550 sage: actual = ComplexHermitianEJA.natural_inner_product(Xe,Ye)
1551 sage: actual == expected
1552 True
1553
1554 """
1555 return RealMatrixEuclideanJordanAlgebra.natural_inner_product(X,Y)/2
1556
1557
1558 class ComplexHermitianEJA(ComplexMatrixEuclideanJordanAlgebra, KnownRankEJA):
1559 """
1560 The rank-n simple EJA consisting of complex Hermitian n-by-n
1561 matrices over the real numbers, the usual symmetric Jordan product,
1562 and the real-part-of-trace inner product. It has dimension `n^2` over
1563 the reals.
1564
1565 SETUP::
1566
1567 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1568
1569 EXAMPLES:
1570
1571 In theory, our "field" can be any subfield of the reals::
1572
1573 sage: ComplexHermitianEJA(2, RDF)
1574 Euclidean Jordan algebra of dimension 4 over Real Double Field
1575 sage: ComplexHermitianEJA(2, RR)
1576 Euclidean Jordan algebra of dimension 4 over Real Field with
1577 53 bits of precision
1578
1579 TESTS:
1580
1581 The dimension of this algebra is `n^2`::
1582
1583 sage: set_random_seed()
1584 sage: n_max = ComplexHermitianEJA._max_test_case_size()
1585 sage: n = ZZ.random_element(1, n_max)
1586 sage: J = ComplexHermitianEJA(n)
1587 sage: J.dimension() == n^2
1588 True
1589
1590 The Jordan multiplication is what we think it is::
1591
1592 sage: set_random_seed()
1593 sage: J = ComplexHermitianEJA.random_instance()
1594 sage: x,y = J.random_elements(2)
1595 sage: actual = (x*y).natural_representation()
1596 sage: X = x.natural_representation()
1597 sage: Y = y.natural_representation()
1598 sage: expected = (X*Y + Y*X)/2
1599 sage: actual == expected
1600 True
1601 sage: J(expected) == x*y
1602 True
1603
1604 We can change the generator prefix::
1605
1606 sage: ComplexHermitianEJA(2, prefix='z').gens()
1607 (z0, z1, z2, z3)
1608
1609 Our natural basis is normalized with respect to the natural inner
1610 product unless we specify otherwise::
1611
1612 sage: set_random_seed()
1613 sage: J = ComplexHermitianEJA.random_instance()
1614 sage: all( b.norm() == 1 for b in J.gens() )
1615 True
1616
1617 Since our natural basis is normalized with respect to the natural
1618 inner product, and since we know that this algebra is an EJA, any
1619 left-multiplication operator's matrix will be symmetric because
1620 natural->EJA basis representation is an isometry and within the EJA
1621 the operator is self-adjoint by the Jordan axiom::
1622
1623 sage: set_random_seed()
1624 sage: x = ComplexHermitianEJA.random_instance().random_element()
1625 sage: x.operator().matrix().is_symmetric()
1626 True
1627
1628 """
1629
1630 @classmethod
1631 def _denormalized_basis(cls, n, field):
1632 """
1633 Returns a basis for the space of complex Hermitian n-by-n matrices.
1634
1635 Why do we embed these? Basically, because all of numerical linear
1636 algebra assumes that you're working with vectors consisting of `n`
1637 entries from a field and scalars from the same field. There's no way
1638 to tell SageMath that (for example) the vectors contain complex
1639 numbers, while the scalar field is real.
1640
1641 SETUP::
1642
1643 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1644
1645 TESTS::
1646
1647 sage: set_random_seed()
1648 sage: n = ZZ.random_element(1,5)
1649 sage: field = QuadraticField(2, 'sqrt2')
1650 sage: B = ComplexHermitianEJA._denormalized_basis(n, field)
1651 sage: all( M.is_symmetric() for M in B)
1652 True
1653
1654 """
1655 R = PolynomialRing(field, 'z')
1656 z = R.gen()
1657 F = field.extension(z**2 + 1, 'I')
1658 I = F.gen()
1659
1660 # This is like the symmetric case, but we need to be careful:
1661 #
1662 # * We want conjugate-symmetry, not just symmetry.
1663 # * The diagonal will (as a result) be real.
1664 #
1665 S = []
1666 for i in range(n):
1667 for j in range(i+1):
1668 Eij = matrix(F, n, lambda k,l: k==i and l==j)
1669 if i == j:
1670 Sij = cls.real_embed(Eij)
1671 S.append(Sij)
1672 else:
1673 # The second one has a minus because it's conjugated.
1674 Sij_real = cls.real_embed(Eij + Eij.transpose())
1675 S.append(Sij_real)
1676 Sij_imag = cls.real_embed(I*Eij - I*Eij.transpose())
1677 S.append(Sij_imag)
1678
1679 # Since we embedded these, we can drop back to the "field" that we
1680 # started with instead of the complex extension "F".
1681 return ( s.change_ring(field) for s in S )
1682
1683
1684 def __init__(self, n, field=AA, **kwargs):
1685 basis = self._denormalized_basis(n,field)
1686 super(ComplexHermitianEJA,self).__init__(field, basis, n, **kwargs)
1687
1688
1689 class QuaternionMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1690 @staticmethod
1691 def real_embed(M):
1692 """
1693 Embed the n-by-n quaternion matrix ``M`` into the space of real
1694 matrices of size 4n-by-4n by first sending each quaternion entry `z
1695 = a + bi + cj + dk` to the block-complex matrix ``[[a + bi,
1696 c+di],[-c + di, a-bi]]`, and then embedding those into a real
1697 matrix.
1698
1699 SETUP::
1700
1701 sage: from mjo.eja.eja_algebra import \
1702 ....: QuaternionMatrixEuclideanJordanAlgebra
1703
1704 EXAMPLES::
1705
1706 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1707 sage: i,j,k = Q.gens()
1708 sage: x = 1 + 2*i + 3*j + 4*k
1709 sage: M = matrix(Q, 1, [[x]])
1710 sage: QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1711 [ 1 2 3 4]
1712 [-2 1 -4 3]
1713 [-3 4 1 -2]
1714 [-4 -3 2 1]
1715
1716 Embedding is a homomorphism (isomorphism, in fact)::
1717
1718 sage: set_random_seed()
1719 sage: n_max = QuaternionMatrixEuclideanJordanAlgebra._max_test_case_size()
1720 sage: n = ZZ.random_element(n_max)
1721 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1722 sage: X = random_matrix(Q, n)
1723 sage: Y = random_matrix(Q, n)
1724 sage: Xe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X)
1725 sage: Ye = QuaternionMatrixEuclideanJordanAlgebra.real_embed(Y)
1726 sage: XYe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1727 sage: Xe*Ye == XYe
1728 True
1729
1730 """
1731 quaternions = M.base_ring()
1732 n = M.nrows()
1733 if M.ncols() != n:
1734 raise ValueError("the matrix 'M' must be square")
1735
1736 F = QuadraticField(-1, 'I')
1737 i = F.gen()
1738
1739 blocks = []
1740 for z in M.list():
1741 t = z.coefficient_tuple()
1742 a = t[0]
1743 b = t[1]
1744 c = t[2]
1745 d = t[3]
1746 cplxM = matrix(F, 2, [[ a + b*i, c + d*i],
1747 [-c + d*i, a - b*i]])
1748 realM = ComplexMatrixEuclideanJordanAlgebra.real_embed(cplxM)
1749 blocks.append(realM)
1750
1751 # We should have real entries by now, so use the realest field
1752 # we've got for the return value.
1753 return matrix.block(quaternions.base_ring(), n, blocks)
1754
1755
1756
1757 @staticmethod
1758 def real_unembed(M):
1759 """
1760 The inverse of _embed_quaternion_matrix().
1761
1762 SETUP::
1763
1764 sage: from mjo.eja.eja_algebra import \
1765 ....: QuaternionMatrixEuclideanJordanAlgebra
1766
1767 EXAMPLES::
1768
1769 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
1770 ....: [-2, 1, -4, 3],
1771 ....: [-3, 4, 1, -2],
1772 ....: [-4, -3, 2, 1]])
1773 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(M)
1774 [1 + 2*i + 3*j + 4*k]
1775
1776 TESTS:
1777
1778 Unembedding is the inverse of embedding::
1779
1780 sage: set_random_seed()
1781 sage: Q = QuaternionAlgebra(QQ, -1, -1)
1782 sage: M = random_matrix(Q, 3)
1783 sage: Me = QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1784 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1785 True
1786
1787 """
1788 n = ZZ(M.nrows())
1789 if M.ncols() != n:
1790 raise ValueError("the matrix 'M' must be square")
1791 if not n.mod(4).is_zero():
1792 raise ValueError("the matrix 'M' must be a quaternion embedding")
1793
1794 # Use the base ring of the matrix to ensure that its entries can be
1795 # multiplied by elements of the quaternion algebra.
1796 field = M.base_ring()
1797 Q = QuaternionAlgebra(field,-1,-1)
1798 i,j,k = Q.gens()
1799
1800 # Go top-left to bottom-right (reading order), converting every
1801 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
1802 # quaternion block.
1803 elements = []
1804 for l in range(n/4):
1805 for m in range(n/4):
1806 submat = ComplexMatrixEuclideanJordanAlgebra.real_unembed(
1807 M[4*l:4*l+4,4*m:4*m+4] )
1808 if submat[0,0] != submat[1,1].conjugate():
1809 raise ValueError('bad on-diagonal submatrix')
1810 if submat[0,1] != -submat[1,0].conjugate():
1811 raise ValueError('bad off-diagonal submatrix')
1812 z = submat[0,0].real()
1813 z += submat[0,0].imag()*i
1814 z += submat[0,1].real()*j
1815 z += submat[0,1].imag()*k
1816 elements.append(z)
1817
1818 return matrix(Q, n/4, elements)
1819
1820
1821 @classmethod
1822 def natural_inner_product(cls,X,Y):
1823 """
1824 Compute a natural inner product in this algebra directly from
1825 its real embedding.
1826
1827 SETUP::
1828
1829 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1830
1831 TESTS:
1832
1833 This gives the same answer as the slow, default method implemented
1834 in :class:`MatrixEuclideanJordanAlgebra`::
1835
1836 sage: set_random_seed()
1837 sage: J = QuaternionHermitianEJA.random_instance()
1838 sage: x,y = J.random_elements(2)
1839 sage: Xe = x.natural_representation()
1840 sage: Ye = y.natural_representation()
1841 sage: X = QuaternionHermitianEJA.real_unembed(Xe)
1842 sage: Y = QuaternionHermitianEJA.real_unembed(Ye)
1843 sage: expected = (X*Y).trace().coefficient_tuple()[0]
1844 sage: actual = QuaternionHermitianEJA.natural_inner_product(Xe,Ye)
1845 sage: actual == expected
1846 True
1847
1848 """
1849 return RealMatrixEuclideanJordanAlgebra.natural_inner_product(X,Y)/4
1850
1851
1852 class QuaternionHermitianEJA(QuaternionMatrixEuclideanJordanAlgebra,
1853 KnownRankEJA):
1854 """
1855 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
1856 matrices, the usual symmetric Jordan product, and the
1857 real-part-of-trace inner product. It has dimension `2n^2 - n` over
1858 the reals.
1859
1860 SETUP::
1861
1862 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1863
1864 EXAMPLES:
1865
1866 In theory, our "field" can be any subfield of the reals::
1867
1868 sage: QuaternionHermitianEJA(2, RDF)
1869 Euclidean Jordan algebra of dimension 6 over Real Double Field
1870 sage: QuaternionHermitianEJA(2, RR)
1871 Euclidean Jordan algebra of dimension 6 over Real Field with
1872 53 bits of precision
1873
1874 TESTS:
1875
1876 The dimension of this algebra is `2*n^2 - n`::
1877
1878 sage: set_random_seed()
1879 sage: n_max = QuaternionHermitianEJA._max_test_case_size()
1880 sage: n = ZZ.random_element(1, n_max)
1881 sage: J = QuaternionHermitianEJA(n)
1882 sage: J.dimension() == 2*(n^2) - n
1883 True
1884
1885 The Jordan multiplication is what we think it is::
1886
1887 sage: set_random_seed()
1888 sage: J = QuaternionHermitianEJA.random_instance()
1889 sage: x,y = J.random_elements(2)
1890 sage: actual = (x*y).natural_representation()
1891 sage: X = x.natural_representation()
1892 sage: Y = y.natural_representation()
1893 sage: expected = (X*Y + Y*X)/2
1894 sage: actual == expected
1895 True
1896 sage: J(expected) == x*y
1897 True
1898
1899 We can change the generator prefix::
1900
1901 sage: QuaternionHermitianEJA(2, prefix='a').gens()
1902 (a0, a1, a2, a3, a4, a5)
1903
1904 Our natural basis is normalized with respect to the natural inner
1905 product unless we specify otherwise::
1906
1907 sage: set_random_seed()
1908 sage: J = QuaternionHermitianEJA.random_instance()
1909 sage: all( b.norm() == 1 for b in J.gens() )
1910 True
1911
1912 Since our natural basis is normalized with respect to the natural
1913 inner product, and since we know that this algebra is an EJA, any
1914 left-multiplication operator's matrix will be symmetric because
1915 natural->EJA basis representation is an isometry and within the EJA
1916 the operator is self-adjoint by the Jordan axiom::
1917
1918 sage: set_random_seed()
1919 sage: x = QuaternionHermitianEJA.random_instance().random_element()
1920 sage: x.operator().matrix().is_symmetric()
1921 True
1922
1923 """
1924 @classmethod
1925 def _denormalized_basis(cls, n, field):
1926 """
1927 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
1928
1929 Why do we embed these? Basically, because all of numerical
1930 linear algebra assumes that you're working with vectors consisting
1931 of `n` entries from a field and scalars from the same field. There's
1932 no way to tell SageMath that (for example) the vectors contain
1933 complex numbers, while the scalar field is real.
1934
1935 SETUP::
1936
1937 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1938
1939 TESTS::
1940
1941 sage: set_random_seed()
1942 sage: n = ZZ.random_element(1,5)
1943 sage: B = QuaternionHermitianEJA._denormalized_basis(n,QQ)
1944 sage: all( M.is_symmetric() for M in B )
1945 True
1946
1947 """
1948 Q = QuaternionAlgebra(QQ,-1,-1)
1949 I,J,K = Q.gens()
1950
1951 # This is like the symmetric case, but we need to be careful:
1952 #
1953 # * We want conjugate-symmetry, not just symmetry.
1954 # * The diagonal will (as a result) be real.
1955 #
1956 S = []
1957 for i in range(n):
1958 for j in range(i+1):
1959 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
1960 if i == j:
1961 Sij = cls.real_embed(Eij)
1962 S.append(Sij)
1963 else:
1964 # The second, third, and fourth ones have a minus
1965 # because they're conjugated.
1966 Sij_real = cls.real_embed(Eij + Eij.transpose())
1967 S.append(Sij_real)
1968 Sij_I = cls.real_embed(I*Eij - I*Eij.transpose())
1969 S.append(Sij_I)
1970 Sij_J = cls.real_embed(J*Eij - J*Eij.transpose())
1971 S.append(Sij_J)
1972 Sij_K = cls.real_embed(K*Eij - K*Eij.transpose())
1973 S.append(Sij_K)
1974
1975 # Since we embedded these, we can drop back to the "field" that we
1976 # started with instead of the quaternion algebra "Q".
1977 return ( s.change_ring(field) for s in S )
1978
1979
1980 def __init__(self, n, field=AA, **kwargs):
1981 basis = self._denormalized_basis(n,field)
1982 super(QuaternionHermitianEJA,self).__init__(field, basis, n, **kwargs)
1983
1984
1985 class BilinearFormEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
1986 r"""
1987 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
1988 with the half-trace inner product and jordan product ``x*y =
1989 (x0*y0 + <B*x_bar,y_bar>, x0*y_bar + y0*x_bar)`` where ``B`` is a
1990 symmetric positive-definite "bilinear form" matrix. It has
1991 dimension `n` over the reals, and reduces to the ``JordanSpinEJA``
1992 when ``B`` is the identity matrix of order ``n-1``.
1993
1994 SETUP::
1995
1996 sage: from mjo.eja.eja_algebra import (BilinearFormEJA,
1997 ....: JordanSpinEJA)
1998
1999 EXAMPLES:
2000
2001 When no bilinear form is specified, the identity matrix is used,
2002 and the resulting algebra is the Jordan spin algebra::
2003
2004 sage: J0 = BilinearFormEJA(3)
2005 sage: J1 = JordanSpinEJA(3)
2006 sage: J0.multiplication_table() == J0.multiplication_table()
2007 True
2008
2009 TESTS:
2010
2011 We can create a zero-dimensional algebra::
2012
2013 sage: J = BilinearFormEJA(0)
2014 sage: J.basis()
2015 Finite family {}
2016
2017 We can check the multiplication condition given in the Jordan, von
2018 Neumann, and Wigner paper (and also discussed on my "On the
2019 symmetry..." paper). Note that this relies heavily on the standard
2020 choice of basis, as does anything utilizing the bilinear form matrix::
2021
2022 sage: set_random_seed()
2023 sage: n = ZZ.random_element(5)
2024 sage: M = matrix.random(QQ, max(0,n-1), algorithm='unimodular')
2025 sage: B = M.transpose()*M
2026 sage: J = BilinearFormEJA(n, B=B)
2027 sage: eis = VectorSpace(M.base_ring(), M.ncols()).basis()
2028 sage: V = J.vector_space()
2029 sage: sis = [ J.from_vector(V([0] + (M.inverse()*ei).list()))
2030 ....: for ei in eis ]
2031 sage: actual = [ sis[i]*sis[j]
2032 ....: for i in range(n-1)
2033 ....: for j in range(n-1) ]
2034 sage: expected = [ J.one() if i == j else J.zero()
2035 ....: for i in range(n-1)
2036 ....: for j in range(n-1) ]
2037 sage: actual == expected
2038 True
2039 """
2040 def __init__(self, n, field=AA, B=None, **kwargs):
2041 if B is None:
2042 self._B = matrix.identity(field, max(0,n-1))
2043 else:
2044 self._B = B
2045
2046 V = VectorSpace(field, n)
2047 mult_table = [[V.zero() for j in range(n)] for i in range(n)]
2048 for i in range(n):
2049 for j in range(n):
2050 x = V.gen(i)
2051 y = V.gen(j)
2052 x0 = x[0]
2053 xbar = x[1:]
2054 y0 = y[0]
2055 ybar = y[1:]
2056 z0 = x0*y0 + (self._B*xbar).inner_product(ybar)
2057 zbar = y0*xbar + x0*ybar
2058 z = V([z0] + zbar.list())
2059 mult_table[i][j] = z
2060
2061 # The rank of this algebra is two, unless we're in a
2062 # one-dimensional ambient space (because the rank is bounded
2063 # by the ambient dimension).
2064 fdeja = super(BilinearFormEJA, self)
2065 return fdeja.__init__(field, mult_table, rank=min(n,2), **kwargs)
2066
2067 def inner_product(self, x, y):
2068 r"""
2069 Half of the trace inner product.
2070
2071 This is defined so that the special case of the Jordan spin
2072 algebra gets the usual inner product.
2073
2074 SETUP::
2075
2076 sage: from mjo.eja.eja_algebra import BilinearFormEJA
2077
2078 TESTS:
2079
2080 Ensure that this is one-half of the trace inner-product when
2081 the algebra isn't just the reals (when ``n`` isn't one). This
2082 is in Faraut and Koranyi, and also my "On the symmetry..."
2083 paper::
2084
2085 sage: set_random_seed()
2086 sage: n = ZZ.random_element(2,5)
2087 sage: M = matrix.random(QQ, max(0,n-1), algorithm='unimodular')
2088 sage: B = M.transpose()*M
2089 sage: J = BilinearFormEJA(n, B=B)
2090 sage: x = J.random_element()
2091 sage: y = J.random_element()
2092 sage: x.inner_product(y) == (x*y).trace()/2
2093 True
2094
2095 """
2096 xvec = x.to_vector()
2097 xbar = xvec[1:]
2098 yvec = y.to_vector()
2099 ybar = yvec[1:]
2100 return x[0]*y[0] + (self._B*xbar).inner_product(ybar)
2101
2102
2103 class JordanSpinEJA(BilinearFormEJA):
2104 """
2105 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2106 with the usual inner product and jordan product ``x*y =
2107 (<x,y>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
2108 the reals.
2109
2110 SETUP::
2111
2112 sage: from mjo.eja.eja_algebra import JordanSpinEJA
2113
2114 EXAMPLES:
2115
2116 This multiplication table can be verified by hand::
2117
2118 sage: J = JordanSpinEJA(4)
2119 sage: e0,e1,e2,e3 = J.gens()
2120 sage: e0*e0
2121 e0
2122 sage: e0*e1
2123 e1
2124 sage: e0*e2
2125 e2
2126 sage: e0*e3
2127 e3
2128 sage: e1*e2
2129 0
2130 sage: e1*e3
2131 0
2132 sage: e2*e3
2133 0
2134
2135 We can change the generator prefix::
2136
2137 sage: JordanSpinEJA(2, prefix='B').gens()
2138 (B0, B1)
2139
2140 TESTS:
2141
2142 Ensure that we have the usual inner product on `R^n`::
2143
2144 sage: set_random_seed()
2145 sage: J = JordanSpinEJA.random_instance()
2146 sage: x,y = J.random_elements(2)
2147 sage: X = x.natural_representation()
2148 sage: Y = y.natural_representation()
2149 sage: x.inner_product(y) == J.natural_inner_product(X,Y)
2150 True
2151
2152 """
2153 def __init__(self, n, field=AA, **kwargs):
2154 # This is a special case of the BilinearFormEJA with the identity
2155 # matrix as its bilinear form.
2156 return super(JordanSpinEJA, self).__init__(n, field, **kwargs)
2157
2158
2159 class TrivialEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
2160 """
2161 The trivial Euclidean Jordan algebra consisting of only a zero element.
2162
2163 SETUP::
2164
2165 sage: from mjo.eja.eja_algebra import TrivialEJA
2166
2167 EXAMPLES::
2168
2169 sage: J = TrivialEJA()
2170 sage: J.dimension()
2171 0
2172 sage: J.zero()
2173 0
2174 sage: J.one()
2175 0
2176 sage: 7*J.one()*12*J.one()
2177 0
2178 sage: J.one().inner_product(J.one())
2179 0
2180 sage: J.one().norm()
2181 0
2182 sage: J.one().subalgebra_generated_by()
2183 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
2184 sage: J.rank()
2185 0
2186
2187 """
2188 def __init__(self, field=AA, **kwargs):
2189 mult_table = []
2190 fdeja = super(TrivialEJA, self)
2191 # The rank is zero using my definition, namely the dimension of the
2192 # largest subalgebra generated by any element.
2193 return fdeja.__init__(field, mult_table, rank=0, **kwargs)