]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: re-disable base ring to algebra coercion.
[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, 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=QQ)
211 Euclidean Jordan algebra of dimension 2 over Rational 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 1/2*sqrt2] [0 0]
555 [0 0], [1/2*sqrt2 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().left_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(self):
798 """
799 Return the rank of this EJA.
800
801 ALGORITHM:
802
803 The author knows of no algorithm to compute the rank of an EJA
804 where only the multiplication table is known. In lieu of one, we
805 require the rank to be specified when the algebra is created,
806 and simply pass along that number here.
807
808 SETUP::
809
810 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
811 ....: RealSymmetricEJA,
812 ....: ComplexHermitianEJA,
813 ....: QuaternionHermitianEJA,
814 ....: random_eja)
815
816 EXAMPLES:
817
818 The rank of the Jordan spin algebra is always two::
819
820 sage: JordanSpinEJA(2).rank()
821 2
822 sage: JordanSpinEJA(3).rank()
823 2
824 sage: JordanSpinEJA(4).rank()
825 2
826
827 The rank of the `n`-by-`n` Hermitian real, complex, or
828 quaternion matrices is `n`::
829
830 sage: RealSymmetricEJA(4).rank()
831 4
832 sage: ComplexHermitianEJA(3).rank()
833 3
834 sage: QuaternionHermitianEJA(2).rank()
835 2
836
837 TESTS:
838
839 Ensure that every EJA that we know how to construct has a
840 positive integer rank, unless the algebra is trivial in
841 which case its rank will be zero::
842
843 sage: set_random_seed()
844 sage: J = random_eja()
845 sage: r = J.rank()
846 sage: r in ZZ
847 True
848 sage: r > 0 or (r == 0 and J.is_trivial())
849 True
850
851 """
852 return self._rank
853
854
855 def vector_space(self):
856 """
857 Return the vector space that underlies this algebra.
858
859 SETUP::
860
861 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
862
863 EXAMPLES::
864
865 sage: J = RealSymmetricEJA(2)
866 sage: J.vector_space()
867 Vector space of dimension 3 over...
868
869 """
870 return self.zero().to_vector().parent().ambient_vector_space()
871
872
873 Element = FiniteDimensionalEuclideanJordanAlgebraElement
874
875
876 class KnownRankEJA(object):
877 """
878 A class for algebras that we actually know we can construct. The
879 main issue is that, for most of our methods to make sense, we need
880 to know the rank of our algebra. Thus we can't simply generate a
881 "random" algebra, or even check that a given basis and product
882 satisfy the axioms; because even if everything looks OK, we wouldn't
883 know the rank we need to actuallty build the thing.
884
885 Not really a subclass of FDEJA because doing that causes method
886 resolution errors, e.g.
887
888 TypeError: Error when calling the metaclass bases
889 Cannot create a consistent method resolution
890 order (MRO) for bases FiniteDimensionalEuclideanJordanAlgebra,
891 KnownRankEJA
892
893 """
894 @staticmethod
895 def _max_test_case_size():
896 """
897 Return an integer "size" that is an upper bound on the size of
898 this algebra when it is used in a random test
899 case. Unfortunately, the term "size" is quite vague -- when
900 dealing with `R^n` under either the Hadamard or Jordan spin
901 product, the "size" refers to the dimension `n`. When dealing
902 with a matrix algebra (real symmetric or complex/quaternion
903 Hermitian), it refers to the size of the matrix, which is
904 far less than the dimension of the underlying vector space.
905
906 We default to five in this class, which is safe in `R^n`. The
907 matrix algebra subclasses (or any class where the "size" is
908 interpreted to be far less than the dimension) should override
909 with a smaller number.
910 """
911 return 5
912
913 @classmethod
914 def random_instance(cls, field=QQ, **kwargs):
915 """
916 Return a random instance of this type of algebra.
917
918 Beware, this will crash for "most instances" because the
919 constructor below looks wrong.
920 """
921 if cls is TrivialEJA:
922 # The TrivialEJA class doesn't take an "n" argument because
923 # there's only one.
924 return cls(field)
925
926 n = ZZ.random_element(cls._max_test_case_size()) + 1
927 return cls(n, field, **kwargs)
928
929
930 class HadamardEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
931 """
932 Return the Euclidean Jordan Algebra corresponding to the set
933 `R^n` under the Hadamard product.
934
935 Note: this is nothing more than the Cartesian product of ``n``
936 copies of the spin algebra. Once Cartesian product algebras
937 are implemented, this can go.
938
939 SETUP::
940
941 sage: from mjo.eja.eja_algebra import HadamardEJA
942
943 EXAMPLES:
944
945 This multiplication table can be verified by hand::
946
947 sage: J = HadamardEJA(3)
948 sage: e0,e1,e2 = J.gens()
949 sage: e0*e0
950 e0
951 sage: e0*e1
952 0
953 sage: e0*e2
954 0
955 sage: e1*e1
956 e1
957 sage: e1*e2
958 0
959 sage: e2*e2
960 e2
961
962 TESTS:
963
964 We can change the generator prefix::
965
966 sage: HadamardEJA(3, prefix='r').gens()
967 (r0, r1, r2)
968
969 """
970 def __init__(self, n, field=QQ, **kwargs):
971 V = VectorSpace(field, n)
972 mult_table = [ [ V.gen(i)*(i == j) for j in range(n) ]
973 for i in range(n) ]
974
975 fdeja = super(HadamardEJA, self)
976 return fdeja.__init__(field, mult_table, rank=n, **kwargs)
977
978 def inner_product(self, x, y):
979 """
980 Faster to reimplement than to use natural representations.
981
982 SETUP::
983
984 sage: from mjo.eja.eja_algebra import HadamardEJA
985
986 TESTS:
987
988 Ensure that this is the usual inner product for the algebras
989 over `R^n`::
990
991 sage: set_random_seed()
992 sage: J = HadamardEJA.random_instance()
993 sage: x,y = J.random_elements(2)
994 sage: X = x.natural_representation()
995 sage: Y = y.natural_representation()
996 sage: x.inner_product(y) == J.natural_inner_product(X,Y)
997 True
998
999 """
1000 return x.to_vector().inner_product(y.to_vector())
1001
1002
1003 def random_eja(field=QQ, nontrivial=False):
1004 """
1005 Return a "random" finite-dimensional Euclidean Jordan Algebra.
1006
1007 SETUP::
1008
1009 sage: from mjo.eja.eja_algebra import random_eja
1010
1011 TESTS::
1012
1013 sage: random_eja()
1014 Euclidean Jordan algebra of dimension...
1015
1016 """
1017 eja_classes = KnownRankEJA.__subclasses__()
1018 if nontrivial:
1019 eja_classes.remove(TrivialEJA)
1020 classname = choice(eja_classes)
1021 return classname.random_instance(field=field)
1022
1023
1024
1025
1026
1027
1028 class MatrixEuclideanJordanAlgebra(FiniteDimensionalEuclideanJordanAlgebra):
1029 @staticmethod
1030 def _max_test_case_size():
1031 # Play it safe, since this will be squared and the underlying
1032 # field can have dimension 4 (quaternions) too.
1033 return 2
1034
1035 def __init__(self, field, basis, rank, normalize_basis=True, **kwargs):
1036 """
1037 Compared to the superclass constructor, we take a basis instead of
1038 a multiplication table because the latter can be computed in terms
1039 of the former when the product is known (like it is here).
1040 """
1041 # Used in this class's fast _charpoly_coeff() override.
1042 self._basis_normalizers = None
1043
1044 # We're going to loop through this a few times, so now's a good
1045 # time to ensure that it isn't a generator expression.
1046 basis = tuple(basis)
1047
1048 if rank > 1 and normalize_basis:
1049 # We'll need sqrt(2) to normalize the basis, and this
1050 # winds up in the multiplication table, so the whole
1051 # algebra needs to be over the field extension.
1052 R = PolynomialRing(field, 'z')
1053 z = R.gen()
1054 p = z**2 - 2
1055 if p.is_irreducible():
1056 field = field.extension(p, 'sqrt2', embedding=RLF(2).sqrt())
1057 basis = tuple( s.change_ring(field) for s in basis )
1058 self._basis_normalizers = tuple(
1059 ~(self.natural_inner_product(s,s).sqrt()) for s in basis )
1060 basis = tuple(s*c for (s,c) in zip(basis,self._basis_normalizers))
1061
1062 Qs = self.multiplication_table_from_matrix_basis(basis)
1063
1064 fdeja = super(MatrixEuclideanJordanAlgebra, self)
1065 return fdeja.__init__(field,
1066 Qs,
1067 rank=rank,
1068 natural_basis=basis,
1069 **kwargs)
1070
1071
1072 @cached_method
1073 def _charpoly_coeff(self, i):
1074 """
1075 Override the parent method with something that tries to compute
1076 over a faster (non-extension) field.
1077 """
1078 if self._basis_normalizers is None:
1079 # We didn't normalize, so assume that the basis we started
1080 # with had entries in a nice field.
1081 return super(MatrixEuclideanJordanAlgebra, self)._charpoly_coeff(i)
1082 else:
1083 basis = ( (b/n) for (b,n) in zip(self.natural_basis(),
1084 self._basis_normalizers) )
1085
1086 # Do this over the rationals and convert back at the end.
1087 J = MatrixEuclideanJordanAlgebra(QQ,
1088 basis,
1089 self.rank(),
1090 normalize_basis=False)
1091 (_,x,_,_) = J._charpoly_matrix_system()
1092 p = J._charpoly_coeff(i)
1093 # p might be missing some vars, have to substitute "optionally"
1094 pairs = zip(x.base_ring().gens(), self._basis_normalizers)
1095 substitutions = { v: v*c for (v,c) in pairs }
1096 result = p.subs(substitutions)
1097
1098 # The result of "subs" can be either a coefficient-ring
1099 # element or a polynomial. Gotta handle both cases.
1100 if result in QQ:
1101 return self.base_ring()(result)
1102 else:
1103 return result.change_ring(self.base_ring())
1104
1105
1106 @staticmethod
1107 def multiplication_table_from_matrix_basis(basis):
1108 """
1109 At least three of the five simple Euclidean Jordan algebras have the
1110 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
1111 multiplication on the right is matrix multiplication. Given a basis
1112 for the underlying matrix space, this function returns a
1113 multiplication table (obtained by looping through the basis
1114 elements) for an algebra of those matrices.
1115 """
1116 # In S^2, for example, we nominally have four coordinates even
1117 # though the space is of dimension three only. The vector space V
1118 # is supposed to hold the entire long vector, and the subspace W
1119 # of V will be spanned by the vectors that arise from symmetric
1120 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
1121 field = basis[0].base_ring()
1122 dimension = basis[0].nrows()
1123
1124 V = VectorSpace(field, dimension**2)
1125 W = V.span_of_basis( _mat2vec(s) for s in basis )
1126 n = len(basis)
1127 mult_table = [[W.zero() for j in range(n)] for i in range(n)]
1128 for i in range(n):
1129 for j in range(n):
1130 mat_entry = (basis[i]*basis[j] + basis[j]*basis[i])/2
1131 mult_table[i][j] = W.coordinate_vector(_mat2vec(mat_entry))
1132
1133 return mult_table
1134
1135
1136 @staticmethod
1137 def real_embed(M):
1138 """
1139 Embed the matrix ``M`` into a space of real matrices.
1140
1141 The matrix ``M`` can have entries in any field at the moment:
1142 the real numbers, complex numbers, or quaternions. And although
1143 they are not a field, we can probably support octonions at some
1144 point, too. This function returns a real matrix that "acts like"
1145 the original with respect to matrix multiplication; i.e.
1146
1147 real_embed(M*N) = real_embed(M)*real_embed(N)
1148
1149 """
1150 raise NotImplementedError
1151
1152
1153 @staticmethod
1154 def real_unembed(M):
1155 """
1156 The inverse of :meth:`real_embed`.
1157 """
1158 raise NotImplementedError
1159
1160
1161 @classmethod
1162 def natural_inner_product(cls,X,Y):
1163 Xu = cls.real_unembed(X)
1164 Yu = cls.real_unembed(Y)
1165 tr = (Xu*Yu).trace()
1166
1167 if tr in RLF:
1168 # It's real already.
1169 return tr
1170
1171 # Otherwise, try the thing that works for complex numbers; and
1172 # if that doesn't work, the thing that works for quaternions.
1173 try:
1174 return tr.vector()[0] # real part, imag part is index 1
1175 except AttributeError:
1176 # A quaternions doesn't have a vector() method, but does
1177 # have coefficient_tuple() method that returns the
1178 # coefficients of 1, i, j, and k -- in that order.
1179 return tr.coefficient_tuple()[0]
1180
1181
1182 class RealMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1183 @staticmethod
1184 def real_embed(M):
1185 """
1186 The identity function, for embedding real matrices into real
1187 matrices.
1188 """
1189 return M
1190
1191 @staticmethod
1192 def real_unembed(M):
1193 """
1194 The identity function, for unembedding real matrices from real
1195 matrices.
1196 """
1197 return M
1198
1199
1200 class RealSymmetricEJA(RealMatrixEuclideanJordanAlgebra, KnownRankEJA):
1201 """
1202 The rank-n simple EJA consisting of real symmetric n-by-n
1203 matrices, the usual symmetric Jordan product, and the trace inner
1204 product. It has dimension `(n^2 + n)/2` over the reals.
1205
1206 SETUP::
1207
1208 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1209
1210 EXAMPLES::
1211
1212 sage: J = RealSymmetricEJA(2)
1213 sage: e0, e1, e2 = J.gens()
1214 sage: e0*e0
1215 e0
1216 sage: e1*e1
1217 1/2*e0 + 1/2*e2
1218 sage: e2*e2
1219 e2
1220
1221 In theory, our "field" can be any subfield of the reals::
1222
1223 sage: RealSymmetricEJA(2, AA)
1224 Euclidean Jordan algebra of dimension 3 over Algebraic Real Field
1225 sage: RealSymmetricEJA(2, RR)
1226 Euclidean Jordan algebra of dimension 3 over Real Field with
1227 53 bits of precision
1228
1229 TESTS:
1230
1231 The dimension of this algebra is `(n^2 + n) / 2`::
1232
1233 sage: set_random_seed()
1234 sage: n_max = RealSymmetricEJA._max_test_case_size()
1235 sage: n = ZZ.random_element(1, n_max)
1236 sage: J = RealSymmetricEJA(n)
1237 sage: J.dimension() == (n^2 + n)/2
1238 True
1239
1240 The Jordan multiplication is what we think it is::
1241
1242 sage: set_random_seed()
1243 sage: J = RealSymmetricEJA.random_instance()
1244 sage: x,y = J.random_elements(2)
1245 sage: actual = (x*y).natural_representation()
1246 sage: X = x.natural_representation()
1247 sage: Y = y.natural_representation()
1248 sage: expected = (X*Y + Y*X)/2
1249 sage: actual == expected
1250 True
1251 sage: J(expected) == x*y
1252 True
1253
1254 We can change the generator prefix::
1255
1256 sage: RealSymmetricEJA(3, prefix='q').gens()
1257 (q0, q1, q2, q3, q4, q5)
1258
1259 Our natural basis is normalized with respect to the natural inner
1260 product unless we specify otherwise::
1261
1262 sage: set_random_seed()
1263 sage: J = RealSymmetricEJA.random_instance()
1264 sage: all( b.norm() == 1 for b in J.gens() )
1265 True
1266
1267 Since our natural basis is normalized with respect to the natural
1268 inner product, and since we know that this algebra is an EJA, any
1269 left-multiplication operator's matrix will be symmetric because
1270 natural->EJA basis representation is an isometry and within the EJA
1271 the operator is self-adjoint by the Jordan axiom::
1272
1273 sage: set_random_seed()
1274 sage: x = RealSymmetricEJA.random_instance().random_element()
1275 sage: x.operator().matrix().is_symmetric()
1276 True
1277
1278 """
1279 @classmethod
1280 def _denormalized_basis(cls, n, field):
1281 """
1282 Return a basis for the space of real symmetric n-by-n matrices.
1283
1284 SETUP::
1285
1286 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1287
1288 TESTS::
1289
1290 sage: set_random_seed()
1291 sage: n = ZZ.random_element(1,5)
1292 sage: B = RealSymmetricEJA._denormalized_basis(n,QQ)
1293 sage: all( M.is_symmetric() for M in B)
1294 True
1295
1296 """
1297 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1298 # coordinates.
1299 S = []
1300 for i in range(n):
1301 for j in range(i+1):
1302 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1303 if i == j:
1304 Sij = Eij
1305 else:
1306 Sij = Eij + Eij.transpose()
1307 S.append(Sij)
1308 return S
1309
1310
1311 @staticmethod
1312 def _max_test_case_size():
1313 return 4 # Dimension 10
1314
1315
1316 def __init__(self, n, field=QQ, **kwargs):
1317 basis = self._denormalized_basis(n, field)
1318 super(RealSymmetricEJA, self).__init__(field, basis, n, **kwargs)
1319
1320
1321 class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1322 @staticmethod
1323 def real_embed(M):
1324 """
1325 Embed the n-by-n complex matrix ``M`` into the space of real
1326 matrices of size 2n-by-2n via the map the sends each entry `z = a +
1327 bi` to the block matrix ``[[a,b],[-b,a]]``.
1328
1329 SETUP::
1330
1331 sage: from mjo.eja.eja_algebra import \
1332 ....: ComplexMatrixEuclideanJordanAlgebra
1333
1334 EXAMPLES::
1335
1336 sage: F = QuadraticField(-1, 'i')
1337 sage: x1 = F(4 - 2*i)
1338 sage: x2 = F(1 + 2*i)
1339 sage: x3 = F(-i)
1340 sage: x4 = F(6)
1341 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
1342 sage: ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1343 [ 4 -2| 1 2]
1344 [ 2 4|-2 1]
1345 [-----+-----]
1346 [ 0 -1| 6 0]
1347 [ 1 0| 0 6]
1348
1349 TESTS:
1350
1351 Embedding is a homomorphism (isomorphism, in fact)::
1352
1353 sage: set_random_seed()
1354 sage: n_max = ComplexMatrixEuclideanJordanAlgebra._max_test_case_size()
1355 sage: n = ZZ.random_element(n_max)
1356 sage: F = QuadraticField(-1, 'i')
1357 sage: X = random_matrix(F, n)
1358 sage: Y = random_matrix(F, n)
1359 sage: Xe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X)
1360 sage: Ye = ComplexMatrixEuclideanJordanAlgebra.real_embed(Y)
1361 sage: XYe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1362 sage: Xe*Ye == XYe
1363 True
1364
1365 """
1366 n = M.nrows()
1367 if M.ncols() != n:
1368 raise ValueError("the matrix 'M' must be square")
1369
1370 # We don't need any adjoined elements...
1371 field = M.base_ring().base_ring()
1372
1373 blocks = []
1374 for z in M.list():
1375 a = z.list()[0] # real part, I guess
1376 b = z.list()[1] # imag part, I guess
1377 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1378
1379 return matrix.block(field, n, blocks)
1380
1381
1382 @staticmethod
1383 def real_unembed(M):
1384 """
1385 The inverse of _embed_complex_matrix().
1386
1387 SETUP::
1388
1389 sage: from mjo.eja.eja_algebra import \
1390 ....: ComplexMatrixEuclideanJordanAlgebra
1391
1392 EXAMPLES::
1393
1394 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1395 ....: [-2, 1, -4, 3],
1396 ....: [ 9, 10, 11, 12],
1397 ....: [-10, 9, -12, 11] ])
1398 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(A)
1399 [ 2*i + 1 4*i + 3]
1400 [ 10*i + 9 12*i + 11]
1401
1402 TESTS:
1403
1404 Unembedding is the inverse of embedding::
1405
1406 sage: set_random_seed()
1407 sage: F = QuadraticField(-1, 'i')
1408 sage: M = random_matrix(F, 3)
1409 sage: Me = ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1410 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1411 True
1412
1413 """
1414 n = ZZ(M.nrows())
1415 if M.ncols() != n:
1416 raise ValueError("the matrix 'M' must be square")
1417 if not n.mod(2).is_zero():
1418 raise ValueError("the matrix 'M' must be a complex embedding")
1419
1420 # If "M" was normalized, its base ring might have roots
1421 # adjoined and they can stick around after unembedding.
1422 field = M.base_ring()
1423 R = PolynomialRing(field, 'z')
1424 z = R.gen()
1425 F = field.extension(z**2 + 1, 'i', embedding=CLF(-1).sqrt())
1426 i = F.gen()
1427
1428 # Go top-left to bottom-right (reading order), converting every
1429 # 2-by-2 block we see to a single complex element.
1430 elements = []
1431 for k in range(n/2):
1432 for j in range(n/2):
1433 submat = M[2*k:2*k+2,2*j:2*j+2]
1434 if submat[0,0] != submat[1,1]:
1435 raise ValueError('bad on-diagonal submatrix')
1436 if submat[0,1] != -submat[1,0]:
1437 raise ValueError('bad off-diagonal submatrix')
1438 z = submat[0,0] + submat[0,1]*i
1439 elements.append(z)
1440
1441 return matrix(F, n/2, elements)
1442
1443
1444 @classmethod
1445 def natural_inner_product(cls,X,Y):
1446 """
1447 Compute a natural inner product in this algebra directly from
1448 its real embedding.
1449
1450 SETUP::
1451
1452 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1453
1454 TESTS:
1455
1456 This gives the same answer as the slow, default method implemented
1457 in :class:`MatrixEuclideanJordanAlgebra`::
1458
1459 sage: set_random_seed()
1460 sage: J = ComplexHermitianEJA.random_instance()
1461 sage: x,y = J.random_elements(2)
1462 sage: Xe = x.natural_representation()
1463 sage: Ye = y.natural_representation()
1464 sage: X = ComplexHermitianEJA.real_unembed(Xe)
1465 sage: Y = ComplexHermitianEJA.real_unembed(Ye)
1466 sage: expected = (X*Y).trace().vector()[0]
1467 sage: actual = ComplexHermitianEJA.natural_inner_product(Xe,Ye)
1468 sage: actual == expected
1469 True
1470
1471 """
1472 return RealMatrixEuclideanJordanAlgebra.natural_inner_product(X,Y)/2
1473
1474
1475 class ComplexHermitianEJA(ComplexMatrixEuclideanJordanAlgebra, KnownRankEJA):
1476 """
1477 The rank-n simple EJA consisting of complex Hermitian n-by-n
1478 matrices over the real numbers, the usual symmetric Jordan product,
1479 and the real-part-of-trace inner product. It has dimension `n^2` over
1480 the reals.
1481
1482 SETUP::
1483
1484 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1485
1486 EXAMPLES:
1487
1488 In theory, our "field" can be any subfield of the reals::
1489
1490 sage: ComplexHermitianEJA(2, AA)
1491 Euclidean Jordan algebra of dimension 4 over Algebraic Real Field
1492 sage: ComplexHermitianEJA(2, RR)
1493 Euclidean Jordan algebra of dimension 4 over Real Field with
1494 53 bits of precision
1495
1496 TESTS:
1497
1498 The dimension of this algebra is `n^2`::
1499
1500 sage: set_random_seed()
1501 sage: n_max = ComplexHermitianEJA._max_test_case_size()
1502 sage: n = ZZ.random_element(1, n_max)
1503 sage: J = ComplexHermitianEJA(n)
1504 sage: J.dimension() == n^2
1505 True
1506
1507 The Jordan multiplication is what we think it is::
1508
1509 sage: set_random_seed()
1510 sage: J = ComplexHermitianEJA.random_instance()
1511 sage: x,y = J.random_elements(2)
1512 sage: actual = (x*y).natural_representation()
1513 sage: X = x.natural_representation()
1514 sage: Y = y.natural_representation()
1515 sage: expected = (X*Y + Y*X)/2
1516 sage: actual == expected
1517 True
1518 sage: J(expected) == x*y
1519 True
1520
1521 We can change the generator prefix::
1522
1523 sage: ComplexHermitianEJA(2, prefix='z').gens()
1524 (z0, z1, z2, z3)
1525
1526 Our natural basis is normalized with respect to the natural inner
1527 product unless we specify otherwise::
1528
1529 sage: set_random_seed()
1530 sage: J = ComplexHermitianEJA.random_instance()
1531 sage: all( b.norm() == 1 for b in J.gens() )
1532 True
1533
1534 Since our natural basis is normalized with respect to the natural
1535 inner product, and since we know that this algebra is an EJA, any
1536 left-multiplication operator's matrix will be symmetric because
1537 natural->EJA basis representation is an isometry and within the EJA
1538 the operator is self-adjoint by the Jordan axiom::
1539
1540 sage: set_random_seed()
1541 sage: x = ComplexHermitianEJA.random_instance().random_element()
1542 sage: x.operator().matrix().is_symmetric()
1543 True
1544
1545 """
1546
1547 @classmethod
1548 def _denormalized_basis(cls, n, field):
1549 """
1550 Returns a basis for the space of complex Hermitian n-by-n matrices.
1551
1552 Why do we embed these? Basically, because all of numerical linear
1553 algebra assumes that you're working with vectors consisting of `n`
1554 entries from a field and scalars from the same field. There's no way
1555 to tell SageMath that (for example) the vectors contain complex
1556 numbers, while the scalar field is real.
1557
1558 SETUP::
1559
1560 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1561
1562 TESTS::
1563
1564 sage: set_random_seed()
1565 sage: n = ZZ.random_element(1,5)
1566 sage: field = QuadraticField(2, 'sqrt2')
1567 sage: B = ComplexHermitianEJA._denormalized_basis(n, field)
1568 sage: all( M.is_symmetric() for M in B)
1569 True
1570
1571 """
1572 R = PolynomialRing(field, 'z')
1573 z = R.gen()
1574 F = field.extension(z**2 + 1, 'I')
1575 I = F.gen()
1576
1577 # This is like the symmetric case, but we need to be careful:
1578 #
1579 # * We want conjugate-symmetry, not just symmetry.
1580 # * The diagonal will (as a result) be real.
1581 #
1582 S = []
1583 for i in range(n):
1584 for j in range(i+1):
1585 Eij = matrix(F, n, lambda k,l: k==i and l==j)
1586 if i == j:
1587 Sij = cls.real_embed(Eij)
1588 S.append(Sij)
1589 else:
1590 # The second one has a minus because it's conjugated.
1591 Sij_real = cls.real_embed(Eij + Eij.transpose())
1592 S.append(Sij_real)
1593 Sij_imag = cls.real_embed(I*Eij - I*Eij.transpose())
1594 S.append(Sij_imag)
1595
1596 # Since we embedded these, we can drop back to the "field" that we
1597 # started with instead of the complex extension "F".
1598 return ( s.change_ring(field) for s in S )
1599
1600
1601 def __init__(self, n, field=QQ, **kwargs):
1602 basis = self._denormalized_basis(n,field)
1603 super(ComplexHermitianEJA,self).__init__(field, basis, n, **kwargs)
1604
1605
1606 class QuaternionMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1607 @staticmethod
1608 def real_embed(M):
1609 """
1610 Embed the n-by-n quaternion matrix ``M`` into the space of real
1611 matrices of size 4n-by-4n by first sending each quaternion entry `z
1612 = a + bi + cj + dk` to the block-complex matrix ``[[a + bi,
1613 c+di],[-c + di, a-bi]]`, and then embedding those into a real
1614 matrix.
1615
1616 SETUP::
1617
1618 sage: from mjo.eja.eja_algebra import \
1619 ....: QuaternionMatrixEuclideanJordanAlgebra
1620
1621 EXAMPLES::
1622
1623 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1624 sage: i,j,k = Q.gens()
1625 sage: x = 1 + 2*i + 3*j + 4*k
1626 sage: M = matrix(Q, 1, [[x]])
1627 sage: QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1628 [ 1 2 3 4]
1629 [-2 1 -4 3]
1630 [-3 4 1 -2]
1631 [-4 -3 2 1]
1632
1633 Embedding is a homomorphism (isomorphism, in fact)::
1634
1635 sage: set_random_seed()
1636 sage: n_max = QuaternionMatrixEuclideanJordanAlgebra._max_test_case_size()
1637 sage: n = ZZ.random_element(n_max)
1638 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1639 sage: X = random_matrix(Q, n)
1640 sage: Y = random_matrix(Q, n)
1641 sage: Xe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X)
1642 sage: Ye = QuaternionMatrixEuclideanJordanAlgebra.real_embed(Y)
1643 sage: XYe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1644 sage: Xe*Ye == XYe
1645 True
1646
1647 """
1648 quaternions = M.base_ring()
1649 n = M.nrows()
1650 if M.ncols() != n:
1651 raise ValueError("the matrix 'M' must be square")
1652
1653 F = QuadraticField(-1, 'i')
1654 i = F.gen()
1655
1656 blocks = []
1657 for z in M.list():
1658 t = z.coefficient_tuple()
1659 a = t[0]
1660 b = t[1]
1661 c = t[2]
1662 d = t[3]
1663 cplxM = matrix(F, 2, [[ a + b*i, c + d*i],
1664 [-c + d*i, a - b*i]])
1665 realM = ComplexMatrixEuclideanJordanAlgebra.real_embed(cplxM)
1666 blocks.append(realM)
1667
1668 # We should have real entries by now, so use the realest field
1669 # we've got for the return value.
1670 return matrix.block(quaternions.base_ring(), n, blocks)
1671
1672
1673
1674 @staticmethod
1675 def real_unembed(M):
1676 """
1677 The inverse of _embed_quaternion_matrix().
1678
1679 SETUP::
1680
1681 sage: from mjo.eja.eja_algebra import \
1682 ....: QuaternionMatrixEuclideanJordanAlgebra
1683
1684 EXAMPLES::
1685
1686 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
1687 ....: [-2, 1, -4, 3],
1688 ....: [-3, 4, 1, -2],
1689 ....: [-4, -3, 2, 1]])
1690 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(M)
1691 [1 + 2*i + 3*j + 4*k]
1692
1693 TESTS:
1694
1695 Unembedding is the inverse of embedding::
1696
1697 sage: set_random_seed()
1698 sage: Q = QuaternionAlgebra(QQ, -1, -1)
1699 sage: M = random_matrix(Q, 3)
1700 sage: Me = QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1701 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1702 True
1703
1704 """
1705 n = ZZ(M.nrows())
1706 if M.ncols() != n:
1707 raise ValueError("the matrix 'M' must be square")
1708 if not n.mod(4).is_zero():
1709 raise ValueError("the matrix 'M' must be a quaternion embedding")
1710
1711 # Use the base ring of the matrix to ensure that its entries can be
1712 # multiplied by elements of the quaternion algebra.
1713 field = M.base_ring()
1714 Q = QuaternionAlgebra(field,-1,-1)
1715 i,j,k = Q.gens()
1716
1717 # Go top-left to bottom-right (reading order), converting every
1718 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
1719 # quaternion block.
1720 elements = []
1721 for l in range(n/4):
1722 for m in range(n/4):
1723 submat = ComplexMatrixEuclideanJordanAlgebra.real_unembed(
1724 M[4*l:4*l+4,4*m:4*m+4] )
1725 if submat[0,0] != submat[1,1].conjugate():
1726 raise ValueError('bad on-diagonal submatrix')
1727 if submat[0,1] != -submat[1,0].conjugate():
1728 raise ValueError('bad off-diagonal submatrix')
1729 z = submat[0,0].vector()[0] # real part
1730 z += submat[0,0].vector()[1]*i # imag part
1731 z += submat[0,1].vector()[0]*j # real part
1732 z += submat[0,1].vector()[1]*k # imag part
1733 elements.append(z)
1734
1735 return matrix(Q, n/4, elements)
1736
1737
1738 @classmethod
1739 def natural_inner_product(cls,X,Y):
1740 """
1741 Compute a natural inner product in this algebra directly from
1742 its real embedding.
1743
1744 SETUP::
1745
1746 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1747
1748 TESTS:
1749
1750 This gives the same answer as the slow, default method implemented
1751 in :class:`MatrixEuclideanJordanAlgebra`::
1752
1753 sage: set_random_seed()
1754 sage: J = QuaternionHermitianEJA.random_instance()
1755 sage: x,y = J.random_elements(2)
1756 sage: Xe = x.natural_representation()
1757 sage: Ye = y.natural_representation()
1758 sage: X = QuaternionHermitianEJA.real_unembed(Xe)
1759 sage: Y = QuaternionHermitianEJA.real_unembed(Ye)
1760 sage: expected = (X*Y).trace().coefficient_tuple()[0]
1761 sage: actual = QuaternionHermitianEJA.natural_inner_product(Xe,Ye)
1762 sage: actual == expected
1763 True
1764
1765 """
1766 return RealMatrixEuclideanJordanAlgebra.natural_inner_product(X,Y)/4
1767
1768
1769 class QuaternionHermitianEJA(QuaternionMatrixEuclideanJordanAlgebra,
1770 KnownRankEJA):
1771 """
1772 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
1773 matrices, the usual symmetric Jordan product, and the
1774 real-part-of-trace inner product. It has dimension `2n^2 - n` over
1775 the reals.
1776
1777 SETUP::
1778
1779 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1780
1781 EXAMPLES:
1782
1783 In theory, our "field" can be any subfield of the reals::
1784
1785 sage: QuaternionHermitianEJA(2, AA)
1786 Euclidean Jordan algebra of dimension 6 over Algebraic Real Field
1787 sage: QuaternionHermitianEJA(2, RR)
1788 Euclidean Jordan algebra of dimension 6 over Real Field with
1789 53 bits of precision
1790
1791 TESTS:
1792
1793 The dimension of this algebra is `2*n^2 - n`::
1794
1795 sage: set_random_seed()
1796 sage: n_max = QuaternionHermitianEJA._max_test_case_size()
1797 sage: n = ZZ.random_element(1, n_max)
1798 sage: J = QuaternionHermitianEJA(n)
1799 sage: J.dimension() == 2*(n^2) - n
1800 True
1801
1802 The Jordan multiplication is what we think it is::
1803
1804 sage: set_random_seed()
1805 sage: J = QuaternionHermitianEJA.random_instance()
1806 sage: x,y = J.random_elements(2)
1807 sage: actual = (x*y).natural_representation()
1808 sage: X = x.natural_representation()
1809 sage: Y = y.natural_representation()
1810 sage: expected = (X*Y + Y*X)/2
1811 sage: actual == expected
1812 True
1813 sage: J(expected) == x*y
1814 True
1815
1816 We can change the generator prefix::
1817
1818 sage: QuaternionHermitianEJA(2, prefix='a').gens()
1819 (a0, a1, a2, a3, a4, a5)
1820
1821 Our natural basis is normalized with respect to the natural inner
1822 product unless we specify otherwise::
1823
1824 sage: set_random_seed()
1825 sage: J = QuaternionHermitianEJA.random_instance()
1826 sage: all( b.norm() == 1 for b in J.gens() )
1827 True
1828
1829 Since our natural basis is normalized with respect to the natural
1830 inner product, and since we know that this algebra is an EJA, any
1831 left-multiplication operator's matrix will be symmetric because
1832 natural->EJA basis representation is an isometry and within the EJA
1833 the operator is self-adjoint by the Jordan axiom::
1834
1835 sage: set_random_seed()
1836 sage: x = QuaternionHermitianEJA.random_instance().random_element()
1837 sage: x.operator().matrix().is_symmetric()
1838 True
1839
1840 """
1841 @classmethod
1842 def _denormalized_basis(cls, n, field):
1843 """
1844 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
1845
1846 Why do we embed these? Basically, because all of numerical
1847 linear algebra assumes that you're working with vectors consisting
1848 of `n` entries from a field and scalars from the same field. There's
1849 no way to tell SageMath that (for example) the vectors contain
1850 complex numbers, while the scalar field is real.
1851
1852 SETUP::
1853
1854 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1855
1856 TESTS::
1857
1858 sage: set_random_seed()
1859 sage: n = ZZ.random_element(1,5)
1860 sage: B = QuaternionHermitianEJA._denormalized_basis(n,QQ)
1861 sage: all( M.is_symmetric() for M in B )
1862 True
1863
1864 """
1865 Q = QuaternionAlgebra(QQ,-1,-1)
1866 I,J,K = Q.gens()
1867
1868 # This is like the symmetric case, but we need to be careful:
1869 #
1870 # * We want conjugate-symmetry, not just symmetry.
1871 # * The diagonal will (as a result) be real.
1872 #
1873 S = []
1874 for i in range(n):
1875 for j in range(i+1):
1876 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
1877 if i == j:
1878 Sij = cls.real_embed(Eij)
1879 S.append(Sij)
1880 else:
1881 # The second, third, and fourth ones have a minus
1882 # because they're conjugated.
1883 Sij_real = cls.real_embed(Eij + Eij.transpose())
1884 S.append(Sij_real)
1885 Sij_I = cls.real_embed(I*Eij - I*Eij.transpose())
1886 S.append(Sij_I)
1887 Sij_J = cls.real_embed(J*Eij - J*Eij.transpose())
1888 S.append(Sij_J)
1889 Sij_K = cls.real_embed(K*Eij - K*Eij.transpose())
1890 S.append(Sij_K)
1891
1892 # Since we embedded these, we can drop back to the "field" that we
1893 # started with instead of the quaternion algebra "Q".
1894 return ( s.change_ring(field) for s in S )
1895
1896
1897 def __init__(self, n, field=QQ, **kwargs):
1898 basis = self._denormalized_basis(n,field)
1899 super(QuaternionHermitianEJA,self).__init__(field, basis, n, **kwargs)
1900
1901
1902 class BilinearFormEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
1903 r"""
1904 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
1905 with the half-trace inner product and jordan product ``x*y =
1906 (x0*y0 + <B*x_bar,y_bar>, x0*y_bar + y0*x_bar)`` where ``B`` is a
1907 symmetric positive-definite "bilinear form" matrix. It has
1908 dimension `n` over the reals, and reduces to the ``JordanSpinEJA``
1909 when ``B`` is the identity matrix of order ``n-1``.
1910
1911 SETUP::
1912
1913 sage: from mjo.eja.eja_algebra import (BilinearFormEJA,
1914 ....: JordanSpinEJA)
1915
1916 EXAMPLES:
1917
1918 When no bilinear form is specified, the identity matrix is used,
1919 and the resulting algebra is the Jordan spin algebra::
1920
1921 sage: J0 = BilinearFormEJA(3)
1922 sage: J1 = JordanSpinEJA(3)
1923 sage: J0.multiplication_table() == J0.multiplication_table()
1924 True
1925
1926 TESTS:
1927
1928 We can create a zero-dimensional algebra::
1929
1930 sage: J = BilinearFormEJA(0)
1931 sage: J.basis()
1932 Finite family {}
1933
1934 We can check the multiplication condition given in the Jordan, von
1935 Neumann, and Wigner paper (and also discussed on my "On the
1936 symmetry..." paper). Note that this relies heavily on the standard
1937 choice of basis, as does anything utilizing the bilinear form matrix::
1938
1939 sage: set_random_seed()
1940 sage: n = ZZ.random_element(5)
1941 sage: M = matrix.random(QQ, max(0,n-1), algorithm='unimodular')
1942 sage: B = M.transpose()*M
1943 sage: J = BilinearFormEJA(n, B=B)
1944 sage: eis = VectorSpace(M.base_ring(), M.ncols()).basis()
1945 sage: V = J.vector_space()
1946 sage: sis = [ J.from_vector(V([0] + (M.inverse()*ei).list()))
1947 ....: for ei in eis ]
1948 sage: actual = [ sis[i]*sis[j]
1949 ....: for i in range(n-1)
1950 ....: for j in range(n-1) ]
1951 sage: expected = [ J.one() if i == j else J.zero()
1952 ....: for i in range(n-1)
1953 ....: for j in range(n-1) ]
1954 sage: actual == expected
1955 True
1956 """
1957 def __init__(self, n, field=QQ, B=None, **kwargs):
1958 if B is None:
1959 self._B = matrix.identity(field, max(0,n-1))
1960 else:
1961 self._B = B
1962
1963 V = VectorSpace(field, n)
1964 mult_table = [[V.zero() for j in range(n)] for i in range(n)]
1965 for i in range(n):
1966 for j in range(n):
1967 x = V.gen(i)
1968 y = V.gen(j)
1969 x0 = x[0]
1970 xbar = x[1:]
1971 y0 = y[0]
1972 ybar = y[1:]
1973 z0 = x0*y0 + (self._B*xbar).inner_product(ybar)
1974 zbar = y0*xbar + x0*ybar
1975 z = V([z0] + zbar.list())
1976 mult_table[i][j] = z
1977
1978 # The rank of this algebra is two, unless we're in a
1979 # one-dimensional ambient space (because the rank is bounded
1980 # by the ambient dimension).
1981 fdeja = super(BilinearFormEJA, self)
1982 return fdeja.__init__(field, mult_table, rank=min(n,2), **kwargs)
1983
1984 def inner_product(self, x, y):
1985 r"""
1986 Half of the trace inner product.
1987
1988 This is defined so that the special case of the Jordan spin
1989 algebra gets the usual inner product.
1990
1991 SETUP::
1992
1993 sage: from mjo.eja.eja_algebra import BilinearFormEJA
1994
1995 TESTS:
1996
1997 Ensure that this is one-half of the trace inner-product::
1998
1999 sage: set_random_seed()
2000 sage: n = ZZ.random_element(5)
2001 sage: M = matrix.random(QQ, n-1, algorithm='unimodular')
2002 sage: B = M.transpose()*M
2003 sage: J = BilinearFormEJA(n, B=B)
2004 sage: eis = VectorSpace(M.base_ring(), M.ncols()).basis()
2005 sage: V = J.vector_space()
2006 sage: sis = [ J.from_vector(V([0] + (M.inverse()*ei).list()))
2007 ....: for ei in eis ]
2008 sage: actual = [ sis[i]*sis[j]
2009 ....: for i in range(n-1)
2010 ....: for j in range(n-1) ]
2011 sage: expected = [ J.one() if i == j else J.zero()
2012 ....: for i in range(n-1)
2013 ....: for j in range(n-1) ]
2014
2015 """
2016 xvec = x.to_vector()
2017 xbar = xvec[1:]
2018 yvec = y.to_vector()
2019 ybar = yvec[1:]
2020 return x[0]*y[0] + (self._B*xbar).inner_product(ybar)
2021
2022 class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
2023 """
2024 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2025 with the usual inner product and jordan product ``x*y =
2026 (<x,y>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
2027 the reals.
2028
2029 SETUP::
2030
2031 sage: from mjo.eja.eja_algebra import JordanSpinEJA
2032
2033 EXAMPLES:
2034
2035 This multiplication table can be verified by hand::
2036
2037 sage: J = JordanSpinEJA(4)
2038 sage: e0,e1,e2,e3 = J.gens()
2039 sage: e0*e0
2040 e0
2041 sage: e0*e1
2042 e1
2043 sage: e0*e2
2044 e2
2045 sage: e0*e3
2046 e3
2047 sage: e1*e2
2048 0
2049 sage: e1*e3
2050 0
2051 sage: e2*e3
2052 0
2053
2054 We can change the generator prefix::
2055
2056 sage: JordanSpinEJA(2, prefix='B').gens()
2057 (B0, B1)
2058
2059 """
2060 def __init__(self, n, field=QQ, **kwargs):
2061 V = VectorSpace(field, n)
2062 mult_table = [[V.zero() for j in range(n)] for i in range(n)]
2063 for i in range(n):
2064 for j in range(n):
2065 x = V.gen(i)
2066 y = V.gen(j)
2067 x0 = x[0]
2068 xbar = x[1:]
2069 y0 = y[0]
2070 ybar = y[1:]
2071 # z = x*y
2072 z0 = x.inner_product(y)
2073 zbar = y0*xbar + x0*ybar
2074 z = V([z0] + zbar.list())
2075 mult_table[i][j] = z
2076
2077 # The rank of the spin algebra is two, unless we're in a
2078 # one-dimensional ambient space (because the rank is bounded by
2079 # the ambient dimension).
2080 fdeja = super(JordanSpinEJA, self)
2081 return fdeja.__init__(field, mult_table, rank=min(n,2), **kwargs)
2082
2083 def inner_product(self, x, y):
2084 """
2085 Faster to reimplement than to use natural representations.
2086
2087 SETUP::
2088
2089 sage: from mjo.eja.eja_algebra import JordanSpinEJA
2090
2091 TESTS:
2092
2093 Ensure that this is the usual inner product for the algebras
2094 over `R^n`::
2095
2096 sage: set_random_seed()
2097 sage: J = JordanSpinEJA.random_instance()
2098 sage: x,y = J.random_elements(2)
2099 sage: X = x.natural_representation()
2100 sage: Y = y.natural_representation()
2101 sage: x.inner_product(y) == J.natural_inner_product(X,Y)
2102 True
2103
2104 """
2105 return x.to_vector().inner_product(y.to_vector())
2106
2107
2108 class TrivialEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
2109 """
2110 The trivial Euclidean Jordan algebra consisting of only a zero element.
2111
2112 SETUP::
2113
2114 sage: from mjo.eja.eja_algebra import TrivialEJA
2115
2116 EXAMPLES::
2117
2118 sage: J = TrivialEJA()
2119 sage: J.dimension()
2120 0
2121 sage: J.zero()
2122 0
2123 sage: J.one()
2124 0
2125 sage: 7*J.one()*12*J.one()
2126 0
2127 sage: J.one().inner_product(J.one())
2128 0
2129 sage: J.one().norm()
2130 0
2131 sage: J.one().subalgebra_generated_by()
2132 Euclidean Jordan algebra of dimension 0 over Rational Field
2133 sage: J.rank()
2134 0
2135
2136 """
2137 def __init__(self, field=QQ, **kwargs):
2138 mult_table = []
2139 fdeja = super(TrivialEJA, self)
2140 # The rank is zero using my definition, namely the dimension of the
2141 # largest subalgebra generated by any element.
2142 return fdeja.__init__(field, mult_table, rank=0, **kwargs)