]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: unbreak DirectSumEJA, poorly.
[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 SETUP::
9
10 sage: from mjo.eja.eja_algebra import random_eja
11
12 EXAMPLES::
13
14 sage: random_eja()
15 Euclidean Jordan algebra of dimension...
16
17 """
18
19 from itertools import repeat
20
21 from sage.algebras.quatalg.quaternion_algebra import QuaternionAlgebra
22 from sage.categories.magmatic_algebras import MagmaticAlgebras
23 from sage.combinat.free_module import CombinatorialFreeModule
24 from sage.matrix.constructor import matrix
25 from sage.matrix.matrix_space import MatrixSpace
26 from sage.misc.cachefunc import cached_method
27 from sage.misc.table import table
28 from sage.modules.free_module import FreeModule, VectorSpace
29 from sage.rings.all import (ZZ, QQ, AA, QQbar, RR, RLF, CLF,
30 PolynomialRing,
31 QuadraticField)
32 from mjo.eja.eja_element import FiniteDimensionalEJAElement
33 from mjo.eja.eja_operator import FiniteDimensionalEJAOperator
34 from mjo.eja.eja_utils import _mat2vec
35
36 class FiniteDimensionalEJA(CombinatorialFreeModule):
37 r"""
38 A finite-dimensional Euclidean Jordan algebra.
39
40 INPUT:
41
42 - basis -- a tuple of basis elements in their matrix form.
43
44 - jordan_product -- function of two elements (in matrix form)
45 that returns their jordan product in this algebra; this will
46 be applied to ``basis`` to compute a multiplication table for
47 the algebra.
48
49 - inner_product -- function of two elements (in matrix form) that
50 returns their inner product. This will be applied to ``basis`` to
51 compute an inner-product table (basically a matrix) for this algebra.
52
53 """
54 Element = FiniteDimensionalEJAElement
55
56 def __init__(self,
57 basis,
58 jordan_product,
59 inner_product,
60 field=AA,
61 orthonormalize=True,
62 associative=False,
63 check_field=True,
64 check_axioms=True,
65 prefix='e'):
66
67 if check_field:
68 if not field.is_subring(RR):
69 # Note: this does return true for the real algebraic
70 # field, the rationals, and any quadratic field where
71 # we've specified a real embedding.
72 raise ValueError("scalar field is not real")
73
74 # If the basis given to us wasn't over the field that it's
75 # supposed to be over, fix that. Or, you know, crash.
76 basis = tuple( b.change_ring(field) for b in basis )
77
78 if check_axioms:
79 # Check commutativity of the Jordan and inner-products.
80 # This has to be done before we build the multiplication
81 # and inner-product tables/matrices, because we take
82 # advantage of symmetry in the process.
83 if not all( jordan_product(bi,bj) == jordan_product(bj,bi)
84 for bi in basis
85 for bj in basis ):
86 raise ValueError("Jordan product is not commutative")
87
88 if not all( inner_product(bi,bj) == inner_product(bj,bi)
89 for bi in basis
90 for bj in basis ):
91 raise ValueError("inner-product is not commutative")
92
93
94 category = MagmaticAlgebras(field).FiniteDimensional()
95 category = category.WithBasis().Unital()
96 if associative:
97 # Element subalgebras can take advantage of this.
98 category = category.Associative()
99
100 # Call the superclass constructor so that we can use its from_vector()
101 # method to build our multiplication table.
102 n = len(basis)
103 super().__init__(field,
104 range(n),
105 prefix=prefix,
106 category=category,
107 bracket=False)
108
109 # Now comes all of the hard work. We'll be constructing an
110 # ambient vector space V that our (vectorized) basis lives in,
111 # as well as a subspace W of V spanned by those (vectorized)
112 # basis elements. The W-coordinates are the coefficients that
113 # we see in things like x = 1*e1 + 2*e2.
114 vector_basis = basis
115
116 degree = 0
117 if n > 0:
118 # Works on both column and square matrices...
119 degree = len(basis[0].list())
120
121 # Build an ambient space that fits our matrix basis when
122 # written out as "long vectors."
123 V = VectorSpace(field, degree)
124
125 # The matrix that will hole the orthonormal -> unorthonormal
126 # coordinate transformation.
127 self._deortho_matrix = None
128
129 if orthonormalize:
130 # Save a copy of the un-orthonormalized basis for later.
131 # Convert it to ambient V (vector) coordinates while we're
132 # at it, because we'd have to do it later anyway.
133 deortho_vector_basis = tuple( V(b.list()) for b in basis )
134
135 from mjo.eja.eja_utils import gram_schmidt
136 basis = tuple(gram_schmidt(basis, inner_product))
137
138 # Save the (possibly orthonormalized) matrix basis for
139 # later...
140 self._matrix_basis = basis
141
142 # Now create the vector space for the algebra, which will have
143 # its own set of non-ambient coordinates (in terms of the
144 # supplied basis).
145 vector_basis = tuple( V(b.list()) for b in basis )
146 W = V.span_of_basis( vector_basis, check=check_axioms)
147
148 if orthonormalize:
149 # Now "W" is the vector space of our algebra coordinates. The
150 # variables "X1", "X2",... refer to the entries of vectors in
151 # W. Thus to convert back and forth between the orthonormal
152 # coordinates and the given ones, we need to stick the original
153 # basis in W.
154 U = V.span_of_basis( deortho_vector_basis, check=check_axioms)
155 self._deortho_matrix = matrix( U.coordinate_vector(q)
156 for q in vector_basis )
157
158
159 # Now we actually compute the multiplication and inner-product
160 # tables/matrices using the possibly-orthonormalized basis.
161 self._inner_product_matrix = matrix.identity(field, n)
162 self._multiplication_table = [ [0 for j in range(i+1)]
163 for i in range(n) ]
164
165 # Note: the Jordan and inner-products are defined in terms
166 # of the ambient basis. It's important that their arguments
167 # are in ambient coordinates as well.
168 for i in range(n):
169 for j in range(i+1):
170 # ortho basis w.r.t. ambient coords
171 q_i = basis[i]
172 q_j = basis[j]
173
174 # The jordan product returns a matrixy answer, so we
175 # have to convert it to the algebra coordinates.
176 elt = jordan_product(q_i, q_j)
177 elt = W.coordinate_vector(V(elt.list()))
178 self._multiplication_table[i][j] = self.from_vector(elt)
179
180 if not orthonormalize:
181 # If we're orthonormalizing the basis with respect
182 # to an inner-product, then the inner-product
183 # matrix with respect to the resulting basis is
184 # just going to be the identity.
185 ip = inner_product(q_i, q_j)
186 self._inner_product_matrix[i,j] = ip
187 self._inner_product_matrix[j,i] = ip
188
189 self._inner_product_matrix._cache = {'hermitian': True}
190 self._inner_product_matrix.set_immutable()
191
192 if check_axioms:
193 if not self._is_jordanian():
194 raise ValueError("Jordan identity does not hold")
195 if not self._inner_product_is_associative():
196 raise ValueError("inner product is not associative")
197
198
199 def _coerce_map_from_base_ring(self):
200 """
201 Disable the map from the base ring into the algebra.
202
203 Performing a nonsense conversion like this automatically
204 is counterpedagogical. The fallback is to try the usual
205 element constructor, which should also fail.
206
207 SETUP::
208
209 sage: from mjo.eja.eja_algebra import random_eja
210
211 TESTS::
212
213 sage: set_random_seed()
214 sage: J = random_eja()
215 sage: J(1)
216 Traceback (most recent call last):
217 ...
218 ValueError: not an element of this algebra
219
220 """
221 return None
222
223
224 def product_on_basis(self, i, j):
225 # We only stored the lower-triangular portion of the
226 # multiplication table.
227 if j <= i:
228 return self._multiplication_table[i][j]
229 else:
230 return self._multiplication_table[j][i]
231
232 def inner_product(self, x, y):
233 """
234 The inner product associated with this Euclidean Jordan algebra.
235
236 Defaults to the trace inner product, but can be overridden by
237 subclasses if they are sure that the necessary properties are
238 satisfied.
239
240 SETUP::
241
242 sage: from mjo.eja.eja_algebra import (random_eja,
243 ....: HadamardEJA,
244 ....: BilinearFormEJA)
245
246 EXAMPLES:
247
248 Our inner product is "associative," which means the following for
249 a symmetric bilinear form::
250
251 sage: set_random_seed()
252 sage: J = random_eja()
253 sage: x,y,z = J.random_elements(3)
254 sage: (x*y).inner_product(z) == y.inner_product(x*z)
255 True
256
257 TESTS:
258
259 Ensure that this is the usual inner product for the algebras
260 over `R^n`::
261
262 sage: set_random_seed()
263 sage: J = HadamardEJA.random_instance()
264 sage: x,y = J.random_elements(2)
265 sage: actual = x.inner_product(y)
266 sage: expected = x.to_vector().inner_product(y.to_vector())
267 sage: actual == expected
268 True
269
270 Ensure that this is one-half of the trace inner-product in a
271 BilinearFormEJA that isn't just the reals (when ``n`` isn't
272 one). This is in Faraut and Koranyi, and also my "On the
273 symmetry..." paper::
274
275 sage: set_random_seed()
276 sage: J = BilinearFormEJA.random_instance()
277 sage: n = J.dimension()
278 sage: x = J.random_element()
279 sage: y = J.random_element()
280 sage: (n == 1) or (x.inner_product(y) == (x*y).trace()/2)
281 True
282 """
283 B = self._inner_product_matrix
284 return (B*x.to_vector()).inner_product(y.to_vector())
285
286
287 def _is_commutative(self):
288 r"""
289 Whether or not this algebra's multiplication table is commutative.
290
291 This method should of course always return ``True``, unless
292 this algebra was constructed with ``check_axioms=False`` and
293 passed an invalid multiplication table.
294 """
295 return all( self.product_on_basis(i,j) == self.product_on_basis(i,j)
296 for i in range(self.dimension())
297 for j in range(self.dimension()) )
298
299 def _is_jordanian(self):
300 r"""
301 Whether or not this algebra's multiplication table respects the
302 Jordan identity `(x^{2})(xy) = x(x^{2}y)`.
303
304 We only check one arrangement of `x` and `y`, so for a
305 ``True`` result to be truly true, you should also check
306 :meth:`_is_commutative`. This method should of course always
307 return ``True``, unless this algebra was constructed with
308 ``check_axioms=False`` and passed an invalid multiplication table.
309 """
310 return all( (self.monomial(i)**2)*(self.monomial(i)*self.monomial(j))
311 ==
312 (self.monomial(i))*((self.monomial(i)**2)*self.monomial(j))
313 for i in range(self.dimension())
314 for j in range(self.dimension()) )
315
316 def _inner_product_is_associative(self):
317 r"""
318 Return whether or not this algebra's inner product `B` is
319 associative; that is, whether or not `B(xy,z) = B(x,yz)`.
320
321 This method should of course always return ``True``, unless
322 this algebra was constructed with ``check_axioms=False`` and
323 passed an invalid Jordan or inner-product.
324 """
325
326 # Used to check whether or not something is zero in an inexact
327 # ring. This number is sufficient to allow the construction of
328 # QuaternionHermitianEJA(2, field=RDF) with check_axioms=True.
329 epsilon = 1e-16
330
331 for i in range(self.dimension()):
332 for j in range(self.dimension()):
333 for k in range(self.dimension()):
334 x = self.monomial(i)
335 y = self.monomial(j)
336 z = self.monomial(k)
337 diff = (x*y).inner_product(z) - x.inner_product(y*z)
338
339 if self.base_ring().is_exact():
340 if diff != 0:
341 return False
342 else:
343 if diff.abs() > epsilon:
344 return False
345
346 return True
347
348 def _element_constructor_(self, elt):
349 """
350 Construct an element of this algebra from its vector or matrix
351 representation.
352
353 This gets called only after the parent element _call_ method
354 fails to find a coercion for the argument.
355
356 SETUP::
357
358 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
359 ....: HadamardEJA,
360 ....: RealSymmetricEJA)
361
362 EXAMPLES:
363
364 The identity in `S^n` is converted to the identity in the EJA::
365
366 sage: J = RealSymmetricEJA(3)
367 sage: I = matrix.identity(QQ,3)
368 sage: J(I) == J.one()
369 True
370
371 This skew-symmetric matrix can't be represented in the EJA::
372
373 sage: J = RealSymmetricEJA(3)
374 sage: A = matrix(QQ,3, lambda i,j: i-j)
375 sage: J(A)
376 Traceback (most recent call last):
377 ...
378 ValueError: not an element of this algebra
379
380 TESTS:
381
382 Ensure that we can convert any element of the two non-matrix
383 simple algebras (whose matrix representations are columns)
384 back and forth faithfully::
385
386 sage: set_random_seed()
387 sage: J = HadamardEJA.random_instance()
388 sage: x = J.random_element()
389 sage: J(x.to_vector().column()) == x
390 True
391 sage: J = JordanSpinEJA.random_instance()
392 sage: x = J.random_element()
393 sage: J(x.to_vector().column()) == x
394 True
395
396 """
397 msg = "not an element of this algebra"
398 if elt == 0:
399 # The superclass implementation of random_element()
400 # needs to be able to coerce "0" into the algebra.
401 return self.zero()
402 elif elt in self.base_ring():
403 # Ensure that no base ring -> algebra coercion is performed
404 # by this method. There's some stupidity in sage that would
405 # otherwise propagate to this method; for example, sage thinks
406 # that the integer 3 belongs to the space of 2-by-2 matrices.
407 raise ValueError(msg)
408
409 try:
410 elt = elt.column()
411 except (AttributeError, TypeError):
412 # Try to convert a vector into a column-matrix
413 pass
414
415 if elt not in self.matrix_space():
416 raise ValueError(msg)
417
418 # Thanks for nothing! Matrix spaces aren't vector spaces in
419 # Sage, so we have to figure out its matrix-basis coordinates
420 # ourselves. We use the basis space's ring instead of the
421 # element's ring because the basis space might be an algebraic
422 # closure whereas the base ring of the 3-by-3 identity matrix
423 # could be QQ instead of QQbar.
424 #
425 # We pass check=False because the matrix basis is "guaranteed"
426 # to be linearly independent... right? Ha ha.
427 V = VectorSpace(self.base_ring(), elt.nrows()*elt.ncols())
428 W = V.span_of_basis( (_mat2vec(s) for s in self.matrix_basis()),
429 check=False)
430
431 try:
432 coords = W.coordinate_vector(_mat2vec(elt))
433 except ArithmeticError: # vector is not in free module
434 raise ValueError(msg)
435
436 return self.from_vector(coords)
437
438 def _repr_(self):
439 """
440 Return a string representation of ``self``.
441
442 SETUP::
443
444 sage: from mjo.eja.eja_algebra import JordanSpinEJA
445
446 TESTS:
447
448 Ensure that it says what we think it says::
449
450 sage: JordanSpinEJA(2, field=AA)
451 Euclidean Jordan algebra of dimension 2 over Algebraic Real Field
452 sage: JordanSpinEJA(3, field=RDF)
453 Euclidean Jordan algebra of dimension 3 over Real Double Field
454
455 """
456 fmt = "Euclidean Jordan algebra of dimension {} over {}"
457 return fmt.format(self.dimension(), self.base_ring())
458
459
460 @cached_method
461 def characteristic_polynomial_of(self):
462 """
463 Return the algebra's "characteristic polynomial of" function,
464 which is itself a multivariate polynomial that, when evaluated
465 at the coordinates of some algebra element, returns that
466 element's characteristic polynomial.
467
468 The resulting polynomial has `n+1` variables, where `n` is the
469 dimension of this algebra. The first `n` variables correspond to
470 the coordinates of an algebra element: when evaluated at the
471 coordinates of an algebra element with respect to a certain
472 basis, the result is a univariate polynomial (in the one
473 remaining variable ``t``), namely the characteristic polynomial
474 of that element.
475
476 SETUP::
477
478 sage: from mjo.eja.eja_algebra import JordanSpinEJA, TrivialEJA
479
480 EXAMPLES:
481
482 The characteristic polynomial in the spin algebra is given in
483 Alizadeh, Example 11.11::
484
485 sage: J = JordanSpinEJA(3)
486 sage: p = J.characteristic_polynomial_of(); p
487 X1^2 - X2^2 - X3^2 + (-2*t)*X1 + t^2
488 sage: xvec = J.one().to_vector()
489 sage: p(*xvec)
490 t^2 - 2*t + 1
491
492 By definition, the characteristic polynomial is a monic
493 degree-zero polynomial in a rank-zero algebra. Note that
494 Cayley-Hamilton is indeed satisfied since the polynomial
495 ``1`` evaluates to the identity element of the algebra on
496 any argument::
497
498 sage: J = TrivialEJA()
499 sage: J.characteristic_polynomial_of()
500 1
501
502 """
503 r = self.rank()
504 n = self.dimension()
505
506 # The list of coefficient polynomials a_0, a_1, a_2, ..., a_(r-1).
507 a = self._charpoly_coefficients()
508
509 # We go to a bit of trouble here to reorder the
510 # indeterminates, so that it's easier to evaluate the
511 # characteristic polynomial at x's coordinates and get back
512 # something in terms of t, which is what we want.
513 S = PolynomialRing(self.base_ring(),'t')
514 t = S.gen(0)
515 if r > 0:
516 R = a[0].parent()
517 S = PolynomialRing(S, R.variable_names())
518 t = S(t)
519
520 return (t**r + sum( a[k]*(t**k) for k in range(r) ))
521
522 def coordinate_polynomial_ring(self):
523 r"""
524 The multivariate polynomial ring in which this algebra's
525 :meth:`characteristic_polynomial_of` lives.
526
527 SETUP::
528
529 sage: from mjo.eja.eja_algebra import (HadamardEJA,
530 ....: RealSymmetricEJA)
531
532 EXAMPLES::
533
534 sage: J = HadamardEJA(2)
535 sage: J.coordinate_polynomial_ring()
536 Multivariate Polynomial Ring in X1, X2...
537 sage: J = RealSymmetricEJA(3,field=QQ,orthonormalize=False)
538 sage: J.coordinate_polynomial_ring()
539 Multivariate Polynomial Ring in X1, X2, X3, X4, X5, X6...
540
541 """
542 var_names = tuple( "X%d" % z for z in range(1, self.dimension()+1) )
543 return PolynomialRing(self.base_ring(), var_names)
544
545 def inner_product(self, x, y):
546 """
547 The inner product associated with this Euclidean Jordan algebra.
548
549 Defaults to the trace inner product, but can be overridden by
550 subclasses if they are sure that the necessary properties are
551 satisfied.
552
553 SETUP::
554
555 sage: from mjo.eja.eja_algebra import (random_eja,
556 ....: HadamardEJA,
557 ....: BilinearFormEJA)
558
559 EXAMPLES:
560
561 Our inner product is "associative," which means the following for
562 a symmetric bilinear form::
563
564 sage: set_random_seed()
565 sage: J = random_eja()
566 sage: x,y,z = J.random_elements(3)
567 sage: (x*y).inner_product(z) == y.inner_product(x*z)
568 True
569
570 TESTS:
571
572 Ensure that this is the usual inner product for the algebras
573 over `R^n`::
574
575 sage: set_random_seed()
576 sage: J = HadamardEJA.random_instance()
577 sage: x,y = J.random_elements(2)
578 sage: actual = x.inner_product(y)
579 sage: expected = x.to_vector().inner_product(y.to_vector())
580 sage: actual == expected
581 True
582
583 Ensure that this is one-half of the trace inner-product in a
584 BilinearFormEJA that isn't just the reals (when ``n`` isn't
585 one). This is in Faraut and Koranyi, and also my "On the
586 symmetry..." paper::
587
588 sage: set_random_seed()
589 sage: J = BilinearFormEJA.random_instance()
590 sage: n = J.dimension()
591 sage: x = J.random_element()
592 sage: y = J.random_element()
593 sage: (n == 1) or (x.inner_product(y) == (x*y).trace()/2)
594 True
595 """
596 B = self._inner_product_matrix
597 return (B*x.to_vector()).inner_product(y.to_vector())
598
599
600 def is_trivial(self):
601 """
602 Return whether or not this algebra is trivial.
603
604 A trivial algebra contains only the zero element.
605
606 SETUP::
607
608 sage: from mjo.eja.eja_algebra import (ComplexHermitianEJA,
609 ....: TrivialEJA)
610
611 EXAMPLES::
612
613 sage: J = ComplexHermitianEJA(3)
614 sage: J.is_trivial()
615 False
616
617 ::
618
619 sage: J = TrivialEJA()
620 sage: J.is_trivial()
621 True
622
623 """
624 return self.dimension() == 0
625
626
627 def multiplication_table(self):
628 """
629 Return a visual representation of this algebra's multiplication
630 table (on basis elements).
631
632 SETUP::
633
634 sage: from mjo.eja.eja_algebra import JordanSpinEJA
635
636 EXAMPLES::
637
638 sage: J = JordanSpinEJA(4)
639 sage: J.multiplication_table()
640 +----++----+----+----+----+
641 | * || e0 | e1 | e2 | e3 |
642 +====++====+====+====+====+
643 | e0 || e0 | e1 | e2 | e3 |
644 +----++----+----+----+----+
645 | e1 || e1 | e0 | 0 | 0 |
646 +----++----+----+----+----+
647 | e2 || e2 | 0 | e0 | 0 |
648 +----++----+----+----+----+
649 | e3 || e3 | 0 | 0 | e0 |
650 +----++----+----+----+----+
651
652 """
653 n = self.dimension()
654 # Prepend the header row.
655 M = [["*"] + list(self.gens())]
656
657 # And to each subsequent row, prepend an entry that belongs to
658 # the left-side "header column."
659 M += [ [self.monomial(i)] + [ self.product_on_basis(i,j)
660 for j in range(n) ]
661 for i in range(n) ]
662
663 return table(M, header_row=True, header_column=True, frame=True)
664
665
666 def matrix_basis(self):
667 """
668 Return an (often more natural) representation of this algebras
669 basis as an ordered tuple of matrices.
670
671 Every finite-dimensional Euclidean Jordan Algebra is a, up to
672 Jordan isomorphism, a direct sum of five simple
673 algebras---four of which comprise Hermitian matrices. And the
674 last type of algebra can of course be thought of as `n`-by-`1`
675 column matrices (ambiguusly called column vectors) to avoid
676 special cases. As a result, matrices (and column vectors) are
677 a natural representation format for Euclidean Jordan algebra
678 elements.
679
680 But, when we construct an algebra from a basis of matrices,
681 those matrix representations are lost in favor of coordinate
682 vectors *with respect to* that basis. We could eventually
683 convert back if we tried hard enough, but having the original
684 representations handy is valuable enough that we simply store
685 them and return them from this method.
686
687 Why implement this for non-matrix algebras? Avoiding special
688 cases for the :class:`BilinearFormEJA` pays with simplicity in
689 its own right. But mainly, we would like to be able to assume
690 that elements of a :class:`DirectSumEJA` can be displayed
691 nicely, without having to have special classes for direct sums
692 one of whose components was a matrix algebra.
693
694 SETUP::
695
696 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
697 ....: RealSymmetricEJA)
698
699 EXAMPLES::
700
701 sage: J = RealSymmetricEJA(2)
702 sage: J.basis()
703 Finite family {0: e0, 1: e1, 2: e2}
704 sage: J.matrix_basis()
705 (
706 [1 0] [ 0 0.7071067811865475?] [0 0]
707 [0 0], [0.7071067811865475? 0], [0 1]
708 )
709
710 ::
711
712 sage: J = JordanSpinEJA(2)
713 sage: J.basis()
714 Finite family {0: e0, 1: e1}
715 sage: J.matrix_basis()
716 (
717 [1] [0]
718 [0], [1]
719 )
720 """
721 return self._matrix_basis
722
723
724 def matrix_space(self):
725 """
726 Return the matrix space in which this algebra's elements live, if
727 we think of them as matrices (including column vectors of the
728 appropriate size).
729
730 Generally this will be an `n`-by-`1` column-vector space,
731 except when the algebra is trivial. There it's `n`-by-`n`
732 (where `n` is zero), to ensure that two elements of the matrix
733 space (empty matrices) can be multiplied.
734
735 Matrix algebras override this with something more useful.
736 """
737 if self.is_trivial():
738 return MatrixSpace(self.base_ring(), 0)
739 else:
740 return self.matrix_basis()[0].parent()
741
742
743 @cached_method
744 def one(self):
745 """
746 Return the unit element of this algebra.
747
748 SETUP::
749
750 sage: from mjo.eja.eja_algebra import (HadamardEJA,
751 ....: random_eja)
752
753 EXAMPLES:
754
755 We can compute unit element in the Hadamard EJA::
756
757 sage: J = HadamardEJA(5)
758 sage: J.one()
759 e0 + e1 + e2 + e3 + e4
760
761 The unit element in the Hadamard EJA is inherited in the
762 subalgebras generated by its elements::
763
764 sage: J = HadamardEJA(5)
765 sage: J.one()
766 e0 + e1 + e2 + e3 + e4
767 sage: x = sum(J.gens())
768 sage: A = x.subalgebra_generated_by(orthonormalize=False)
769 sage: A.one()
770 f0
771 sage: A.one().superalgebra_element()
772 e0 + e1 + e2 + e3 + e4
773
774 TESTS:
775
776 The identity element acts like the identity, regardless of
777 whether or not we orthonormalize::
778
779 sage: set_random_seed()
780 sage: J = random_eja()
781 sage: x = J.random_element()
782 sage: J.one()*x == x and x*J.one() == x
783 True
784 sage: A = x.subalgebra_generated_by()
785 sage: y = A.random_element()
786 sage: A.one()*y == y and y*A.one() == y
787 True
788
789 ::
790
791 sage: set_random_seed()
792 sage: J = random_eja(field=QQ, orthonormalize=False)
793 sage: x = J.random_element()
794 sage: J.one()*x == x and x*J.one() == x
795 True
796 sage: A = x.subalgebra_generated_by(orthonormalize=False)
797 sage: y = A.random_element()
798 sage: A.one()*y == y and y*A.one() == y
799 True
800
801 The matrix of the unit element's operator is the identity,
802 regardless of the base field and whether or not we
803 orthonormalize::
804
805 sage: set_random_seed()
806 sage: J = random_eja()
807 sage: actual = J.one().operator().matrix()
808 sage: expected = matrix.identity(J.base_ring(), J.dimension())
809 sage: actual == expected
810 True
811 sage: x = J.random_element()
812 sage: A = x.subalgebra_generated_by()
813 sage: actual = A.one().operator().matrix()
814 sage: expected = matrix.identity(A.base_ring(), A.dimension())
815 sage: actual == expected
816 True
817
818 ::
819
820 sage: set_random_seed()
821 sage: J = random_eja(field=QQ, orthonormalize=False)
822 sage: actual = J.one().operator().matrix()
823 sage: expected = matrix.identity(J.base_ring(), J.dimension())
824 sage: actual == expected
825 True
826 sage: x = J.random_element()
827 sage: A = x.subalgebra_generated_by(orthonormalize=False)
828 sage: actual = A.one().operator().matrix()
829 sage: expected = matrix.identity(A.base_ring(), A.dimension())
830 sage: actual == expected
831 True
832
833 Ensure that the cached unit element (often precomputed by
834 hand) agrees with the computed one::
835
836 sage: set_random_seed()
837 sage: J = random_eja()
838 sage: cached = J.one()
839 sage: J.one.clear_cache()
840 sage: J.one() == cached
841 True
842
843 ::
844
845 sage: set_random_seed()
846 sage: J = random_eja(field=QQ, orthonormalize=False)
847 sage: cached = J.one()
848 sage: J.one.clear_cache()
849 sage: J.one() == cached
850 True
851
852 """
853 # We can brute-force compute the matrices of the operators
854 # that correspond to the basis elements of this algebra.
855 # If some linear combination of those basis elements is the
856 # algebra identity, then the same linear combination of
857 # their matrices has to be the identity matrix.
858 #
859 # Of course, matrices aren't vectors in sage, so we have to
860 # appeal to the "long vectors" isometry.
861 oper_vecs = [ _mat2vec(g.operator().matrix()) for g in self.gens() ]
862
863 # Now we use basic linear algebra to find the coefficients,
864 # of the matrices-as-vectors-linear-combination, which should
865 # work for the original algebra basis too.
866 A = matrix(self.base_ring(), oper_vecs)
867
868 # We used the isometry on the left-hand side already, but we
869 # still need to do it for the right-hand side. Recall that we
870 # wanted something that summed to the identity matrix.
871 b = _mat2vec( matrix.identity(self.base_ring(), self.dimension()) )
872
873 # Now if there's an identity element in the algebra, this
874 # should work. We solve on the left to avoid having to
875 # transpose the matrix "A".
876 return self.from_vector(A.solve_left(b))
877
878
879 def peirce_decomposition(self, c):
880 """
881 The Peirce decomposition of this algebra relative to the
882 idempotent ``c``.
883
884 In the future, this can be extended to a complete system of
885 orthogonal idempotents.
886
887 INPUT:
888
889 - ``c`` -- an idempotent of this algebra.
890
891 OUTPUT:
892
893 A triple (J0, J5, J1) containing two subalgebras and one subspace
894 of this algebra,
895
896 - ``J0`` -- the algebra on the eigenspace of ``c.operator()``
897 corresponding to the eigenvalue zero.
898
899 - ``J5`` -- the eigenspace (NOT a subalgebra) of ``c.operator()``
900 corresponding to the eigenvalue one-half.
901
902 - ``J1`` -- the algebra on the eigenspace of ``c.operator()``
903 corresponding to the eigenvalue one.
904
905 These are the only possible eigenspaces for that operator, and this
906 algebra is a direct sum of them. The spaces ``J0`` and ``J1`` are
907 orthogonal, and are subalgebras of this algebra with the appropriate
908 restrictions.
909
910 SETUP::
911
912 sage: from mjo.eja.eja_algebra import random_eja, RealSymmetricEJA
913
914 EXAMPLES:
915
916 The canonical example comes from the symmetric matrices, which
917 decompose into diagonal and off-diagonal parts::
918
919 sage: J = RealSymmetricEJA(3)
920 sage: C = matrix(QQ, [ [1,0,0],
921 ....: [0,1,0],
922 ....: [0,0,0] ])
923 sage: c = J(C)
924 sage: J0,J5,J1 = J.peirce_decomposition(c)
925 sage: J0
926 Euclidean Jordan algebra of dimension 1...
927 sage: J5
928 Vector space of degree 6 and dimension 2...
929 sage: J1
930 Euclidean Jordan algebra of dimension 3...
931 sage: J0.one().to_matrix()
932 [0 0 0]
933 [0 0 0]
934 [0 0 1]
935 sage: orig_df = AA.options.display_format
936 sage: AA.options.display_format = 'radical'
937 sage: J.from_vector(J5.basis()[0]).to_matrix()
938 [ 0 0 1/2*sqrt(2)]
939 [ 0 0 0]
940 [1/2*sqrt(2) 0 0]
941 sage: J.from_vector(J5.basis()[1]).to_matrix()
942 [ 0 0 0]
943 [ 0 0 1/2*sqrt(2)]
944 [ 0 1/2*sqrt(2) 0]
945 sage: AA.options.display_format = orig_df
946 sage: J1.one().to_matrix()
947 [1 0 0]
948 [0 1 0]
949 [0 0 0]
950
951 TESTS:
952
953 Every algebra decomposes trivially with respect to its identity
954 element::
955
956 sage: set_random_seed()
957 sage: J = random_eja()
958 sage: J0,J5,J1 = J.peirce_decomposition(J.one())
959 sage: J0.dimension() == 0 and J5.dimension() == 0
960 True
961 sage: J1.superalgebra() == J and J1.dimension() == J.dimension()
962 True
963
964 The decomposition is into eigenspaces, and its components are
965 therefore necessarily orthogonal. Moreover, the identity
966 elements in the two subalgebras are the projections onto their
967 respective subspaces of the superalgebra's identity element::
968
969 sage: set_random_seed()
970 sage: J = random_eja()
971 sage: x = J.random_element()
972 sage: if not J.is_trivial():
973 ....: while x.is_nilpotent():
974 ....: x = J.random_element()
975 sage: c = x.subalgebra_idempotent()
976 sage: J0,J5,J1 = J.peirce_decomposition(c)
977 sage: ipsum = 0
978 sage: for (w,y,z) in zip(J0.basis(), J5.basis(), J1.basis()):
979 ....: w = w.superalgebra_element()
980 ....: y = J.from_vector(y)
981 ....: z = z.superalgebra_element()
982 ....: ipsum += w.inner_product(y).abs()
983 ....: ipsum += w.inner_product(z).abs()
984 ....: ipsum += y.inner_product(z).abs()
985 sage: ipsum
986 0
987 sage: J1(c) == J1.one()
988 True
989 sage: J0(J.one() - c) == J0.one()
990 True
991
992 """
993 if not c.is_idempotent():
994 raise ValueError("element is not idempotent: %s" % c)
995
996 from mjo.eja.eja_subalgebra import FiniteDimensionalEJASubalgebra
997
998 # Default these to what they should be if they turn out to be
999 # trivial, because eigenspaces_left() won't return eigenvalues
1000 # corresponding to trivial spaces (e.g. it returns only the
1001 # eigenspace corresponding to lambda=1 if you take the
1002 # decomposition relative to the identity element).
1003 trivial = FiniteDimensionalEJASubalgebra(self, ())
1004 J0 = trivial # eigenvalue zero
1005 J5 = VectorSpace(self.base_ring(), 0) # eigenvalue one-half
1006 J1 = trivial # eigenvalue one
1007
1008 for (eigval, eigspace) in c.operator().matrix().right_eigenspaces():
1009 if eigval == ~(self.base_ring()(2)):
1010 J5 = eigspace
1011 else:
1012 gens = tuple( self.from_vector(b) for b in eigspace.basis() )
1013 subalg = FiniteDimensionalEJASubalgebra(self,
1014 gens,
1015 check_axioms=False)
1016 if eigval == 0:
1017 J0 = subalg
1018 elif eigval == 1:
1019 J1 = subalg
1020 else:
1021 raise ValueError("unexpected eigenvalue: %s" % eigval)
1022
1023 return (J0, J5, J1)
1024
1025
1026 def random_element(self, thorough=False):
1027 r"""
1028 Return a random element of this algebra.
1029
1030 Our algebra superclass method only returns a linear
1031 combination of at most two basis elements. We instead
1032 want the vector space "random element" method that
1033 returns a more diverse selection.
1034
1035 INPUT:
1036
1037 - ``thorough`` -- (boolean; default False) whether or not we
1038 should generate irrational coefficients for the random
1039 element when our base ring is irrational; this slows the
1040 algebra operations to a crawl, but any truly random method
1041 should include them
1042
1043 """
1044 # For a general base ring... maybe we can trust this to do the
1045 # right thing? Unlikely, but.
1046 V = self.vector_space()
1047 v = V.random_element()
1048
1049 if self.base_ring() is AA:
1050 # The "random element" method of the algebraic reals is
1051 # stupid at the moment, and only returns integers between
1052 # -2 and 2, inclusive:
1053 #
1054 # https://trac.sagemath.org/ticket/30875
1055 #
1056 # Instead, we implement our own "random vector" method,
1057 # and then coerce that into the algebra. We use the vector
1058 # space degree here instead of the dimension because a
1059 # subalgebra could (for example) be spanned by only two
1060 # vectors, each with five coordinates. We need to
1061 # generate all five coordinates.
1062 if thorough:
1063 v *= QQbar.random_element().real()
1064 else:
1065 v *= QQ.random_element()
1066
1067 return self.from_vector(V.coordinate_vector(v))
1068
1069 def random_elements(self, count, thorough=False):
1070 """
1071 Return ``count`` random elements as a tuple.
1072
1073 INPUT:
1074
1075 - ``thorough`` -- (boolean; default False) whether or not we
1076 should generate irrational coefficients for the random
1077 elements when our base ring is irrational; this slows the
1078 algebra operations to a crawl, but any truly random method
1079 should include them
1080
1081 SETUP::
1082
1083 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1084
1085 EXAMPLES::
1086
1087 sage: J = JordanSpinEJA(3)
1088 sage: x,y,z = J.random_elements(3)
1089 sage: all( [ x in J, y in J, z in J ])
1090 True
1091 sage: len( J.random_elements(10) ) == 10
1092 True
1093
1094 """
1095 return tuple( self.random_element(thorough)
1096 for idx in range(count) )
1097
1098
1099 @cached_method
1100 def _charpoly_coefficients(self):
1101 r"""
1102 The `r` polynomial coefficients of the "characteristic polynomial
1103 of" function.
1104 """
1105 n = self.dimension()
1106 R = self.coordinate_polynomial_ring()
1107 vars = R.gens()
1108 F = R.fraction_field()
1109
1110 def L_x_i_j(i,j):
1111 # From a result in my book, these are the entries of the
1112 # basis representation of L_x.
1113 return sum( vars[k]*self.monomial(k).operator().matrix()[i,j]
1114 for k in range(n) )
1115
1116 L_x = matrix(F, n, n, L_x_i_j)
1117
1118 r = None
1119 if self.rank.is_in_cache():
1120 r = self.rank()
1121 # There's no need to pad the system with redundant
1122 # columns if we *know* they'll be redundant.
1123 n = r
1124
1125 # Compute an extra power in case the rank is equal to
1126 # the dimension (otherwise, we would stop at x^(r-1)).
1127 x_powers = [ (L_x**k)*self.one().to_vector()
1128 for k in range(n+1) ]
1129 A = matrix.column(F, x_powers[:n])
1130 AE = A.extended_echelon_form()
1131 E = AE[:,n:]
1132 A_rref = AE[:,:n]
1133 if r is None:
1134 r = A_rref.rank()
1135 b = x_powers[r]
1136
1137 # The theory says that only the first "r" coefficients are
1138 # nonzero, and they actually live in the original polynomial
1139 # ring and not the fraction field. We negate them because
1140 # in the actual characteristic polynomial, they get moved
1141 # to the other side where x^r lives.
1142 return -A_rref.solve_right(E*b).change_ring(R)[:r]
1143
1144 @cached_method
1145 def rank(self):
1146 r"""
1147 Return the rank of this EJA.
1148
1149 This is a cached method because we know the rank a priori for
1150 all of the algebras we can construct. Thus we can avoid the
1151 expensive ``_charpoly_coefficients()`` call unless we truly
1152 need to compute the whole characteristic polynomial.
1153
1154 SETUP::
1155
1156 sage: from mjo.eja.eja_algebra import (HadamardEJA,
1157 ....: JordanSpinEJA,
1158 ....: RealSymmetricEJA,
1159 ....: ComplexHermitianEJA,
1160 ....: QuaternionHermitianEJA,
1161 ....: random_eja)
1162
1163 EXAMPLES:
1164
1165 The rank of the Jordan spin algebra is always two::
1166
1167 sage: JordanSpinEJA(2).rank()
1168 2
1169 sage: JordanSpinEJA(3).rank()
1170 2
1171 sage: JordanSpinEJA(4).rank()
1172 2
1173
1174 The rank of the `n`-by-`n` Hermitian real, complex, or
1175 quaternion matrices is `n`::
1176
1177 sage: RealSymmetricEJA(4).rank()
1178 4
1179 sage: ComplexHermitianEJA(3).rank()
1180 3
1181 sage: QuaternionHermitianEJA(2).rank()
1182 2
1183
1184 TESTS:
1185
1186 Ensure that every EJA that we know how to construct has a
1187 positive integer rank, unless the algebra is trivial in
1188 which case its rank will be zero::
1189
1190 sage: set_random_seed()
1191 sage: J = random_eja()
1192 sage: r = J.rank()
1193 sage: r in ZZ
1194 True
1195 sage: r > 0 or (r == 0 and J.is_trivial())
1196 True
1197
1198 Ensure that computing the rank actually works, since the ranks
1199 of all simple algebras are known and will be cached by default::
1200
1201 sage: set_random_seed() # long time
1202 sage: J = random_eja() # long time
1203 sage: caches = J.rank() # long time
1204 sage: J.rank.clear_cache() # long time
1205 sage: J.rank() == cached # long time
1206 True
1207
1208 """
1209 return len(self._charpoly_coefficients())
1210
1211
1212 def vector_space(self):
1213 """
1214 Return the vector space that underlies this algebra.
1215
1216 SETUP::
1217
1218 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1219
1220 EXAMPLES::
1221
1222 sage: J = RealSymmetricEJA(2)
1223 sage: J.vector_space()
1224 Vector space of dimension 3 over...
1225
1226 """
1227 return self.zero().to_vector().parent().ambient_vector_space()
1228
1229
1230 Element = FiniteDimensionalEJAElement
1231
1232 class RationalBasisEJA(FiniteDimensionalEJA):
1233 r"""
1234 New class for algebras whose supplied basis elements have all rational entries.
1235
1236 SETUP::
1237
1238 sage: from mjo.eja.eja_algebra import BilinearFormEJA
1239
1240 EXAMPLES:
1241
1242 The supplied basis is orthonormalized by default::
1243
1244 sage: B = matrix(QQ, [[1, 0, 0], [0, 25, -32], [0, -32, 41]])
1245 sage: J = BilinearFormEJA(B)
1246 sage: J.matrix_basis()
1247 (
1248 [1] [ 0] [ 0]
1249 [0] [1/5] [32/5]
1250 [0], [ 0], [ 5]
1251 )
1252
1253 """
1254 def __init__(self,
1255 basis,
1256 jordan_product,
1257 inner_product,
1258 field=AA,
1259 check_field=True,
1260 **kwargs):
1261
1262 if check_field:
1263 # Abuse the check_field parameter to check that the entries of
1264 # out basis (in ambient coordinates) are in the field QQ.
1265 if not all( all(b_i in QQ for b_i in b.list()) for b in basis ):
1266 raise TypeError("basis not rational")
1267
1268 self._rational_algebra = None
1269 if field is not QQ:
1270 # There's no point in constructing the extra algebra if this
1271 # one is already rational.
1272 #
1273 # Note: the same Jordan and inner-products work here,
1274 # because they are necessarily defined with respect to
1275 # ambient coordinates and not any particular basis.
1276 self._rational_algebra = FiniteDimensionalEJA(
1277 basis,
1278 jordan_product,
1279 inner_product,
1280 field=QQ,
1281 orthonormalize=False,
1282 check_field=False,
1283 check_axioms=False)
1284
1285 super().__init__(basis,
1286 jordan_product,
1287 inner_product,
1288 field=field,
1289 check_field=check_field,
1290 **kwargs)
1291
1292 @cached_method
1293 def _charpoly_coefficients(self):
1294 r"""
1295 SETUP::
1296
1297 sage: from mjo.eja.eja_algebra import (BilinearFormEJA,
1298 ....: JordanSpinEJA)
1299
1300 EXAMPLES:
1301
1302 The base ring of the resulting polynomial coefficients is what
1303 it should be, and not the rationals (unless the algebra was
1304 already over the rationals)::
1305
1306 sage: J = JordanSpinEJA(3)
1307 sage: J._charpoly_coefficients()
1308 (X1^2 - X2^2 - X3^2, -2*X1)
1309 sage: a0 = J._charpoly_coefficients()[0]
1310 sage: J.base_ring()
1311 Algebraic Real Field
1312 sage: a0.base_ring()
1313 Algebraic Real Field
1314
1315 """
1316 if self._rational_algebra is None:
1317 # There's no need to construct *another* algebra over the
1318 # rationals if this one is already over the
1319 # rationals. Likewise, if we never orthonormalized our
1320 # basis, we might as well just use the given one.
1321 return super()._charpoly_coefficients()
1322
1323 # Do the computation over the rationals. The answer will be
1324 # the same, because all we've done is a change of basis.
1325 # Then, change back from QQ to our real base ring
1326 a = ( a_i.change_ring(self.base_ring())
1327 for a_i in self._rational_algebra._charpoly_coefficients() )
1328
1329 if self._deortho_matrix is None:
1330 # This can happen if our base ring was, say, AA and we
1331 # chose not to (or didn't need to) orthonormalize. It's
1332 # still faster to do the computations over QQ even if
1333 # the numbers in the boxes stay the same.
1334 return tuple(a)
1335
1336 # Otherwise, convert the coordinate variables back to the
1337 # deorthonormalized ones.
1338 R = self.coordinate_polynomial_ring()
1339 from sage.modules.free_module_element import vector
1340 X = vector(R, R.gens())
1341 BX = self._deortho_matrix*X
1342
1343 subs_dict = { X[i]: BX[i] for i in range(len(X)) }
1344 return tuple( a_i.subs(subs_dict) for a_i in a )
1345
1346 class ConcreteEJA(RationalBasisEJA):
1347 r"""
1348 A class for the Euclidean Jordan algebras that we know by name.
1349
1350 These are the Jordan algebras whose basis, multiplication table,
1351 rank, and so on are known a priori. More to the point, they are
1352 the Euclidean Jordan algebras for which we are able to conjure up
1353 a "random instance."
1354
1355 SETUP::
1356
1357 sage: from mjo.eja.eja_algebra import ConcreteEJA
1358
1359 TESTS:
1360
1361 Our basis is normalized with respect to the algebra's inner
1362 product, unless we specify otherwise::
1363
1364 sage: set_random_seed()
1365 sage: J = ConcreteEJA.random_instance()
1366 sage: all( b.norm() == 1 for b in J.gens() )
1367 True
1368
1369 Since our basis is orthonormal with respect to the algebra's inner
1370 product, and since we know that this algebra is an EJA, any
1371 left-multiplication operator's matrix will be symmetric because
1372 natural->EJA basis representation is an isometry and within the
1373 EJA the operator is self-adjoint by the Jordan axiom::
1374
1375 sage: set_random_seed()
1376 sage: J = ConcreteEJA.random_instance()
1377 sage: x = J.random_element()
1378 sage: x.operator().is_self_adjoint()
1379 True
1380 """
1381
1382 @staticmethod
1383 def _max_random_instance_size():
1384 """
1385 Return an integer "size" that is an upper bound on the size of
1386 this algebra when it is used in a random test
1387 case. Unfortunately, the term "size" is ambiguous -- when
1388 dealing with `R^n` under either the Hadamard or Jordan spin
1389 product, the "size" refers to the dimension `n`. When dealing
1390 with a matrix algebra (real symmetric or complex/quaternion
1391 Hermitian), it refers to the size of the matrix, which is far
1392 less than the dimension of the underlying vector space.
1393
1394 This method must be implemented in each subclass.
1395 """
1396 raise NotImplementedError
1397
1398 @classmethod
1399 def random_instance(cls, *args, **kwargs):
1400 """
1401 Return a random instance of this type of algebra.
1402
1403 This method should be implemented in each subclass.
1404 """
1405 from sage.misc.prandom import choice
1406 eja_class = choice(cls.__subclasses__())
1407
1408 # These all bubble up to the RationalBasisEJA superclass
1409 # constructor, so any (kw)args valid there are also valid
1410 # here.
1411 return eja_class.random_instance(*args, **kwargs)
1412
1413
1414 class MatrixEJA:
1415 @staticmethod
1416 def dimension_over_reals():
1417 r"""
1418 The dimension of this matrix's base ring over the reals.
1419
1420 The reals are dimension one over themselves, obviously; that's
1421 just `\mathbb{R}^{1}`. Likewise, the complex numbers `a + bi`
1422 have dimension two. Finally, the quaternions have dimension
1423 four over the reals.
1424
1425 This is used to determine the size of the matrix returned from
1426 :meth:`real_embed`, among other things.
1427 """
1428 raise NotImplementedError
1429
1430 @classmethod
1431 def real_embed(cls,M):
1432 """
1433 Embed the matrix ``M`` into a space of real matrices.
1434
1435 The matrix ``M`` can have entries in any field at the moment:
1436 the real numbers, complex numbers, or quaternions. And although
1437 they are not a field, we can probably support octonions at some
1438 point, too. This function returns a real matrix that "acts like"
1439 the original with respect to matrix multiplication; i.e.
1440
1441 real_embed(M*N) = real_embed(M)*real_embed(N)
1442
1443 """
1444 if M.ncols() != M.nrows():
1445 raise ValueError("the matrix 'M' must be square")
1446 return M
1447
1448
1449 @classmethod
1450 def real_unembed(cls,M):
1451 """
1452 The inverse of :meth:`real_embed`.
1453 """
1454 if M.ncols() != M.nrows():
1455 raise ValueError("the matrix 'M' must be square")
1456 if not ZZ(M.nrows()).mod(cls.dimension_over_reals()).is_zero():
1457 raise ValueError("the matrix 'M' must be a real embedding")
1458 return M
1459
1460 @staticmethod
1461 def jordan_product(X,Y):
1462 return (X*Y + Y*X)/2
1463
1464 @classmethod
1465 def trace_inner_product(cls,X,Y):
1466 r"""
1467 Compute the trace inner-product of two real-embeddings.
1468
1469 SETUP::
1470
1471 sage: from mjo.eja.eja_algebra import (RealSymmetricEJA,
1472 ....: ComplexHermitianEJA,
1473 ....: QuaternionHermitianEJA)
1474
1475 EXAMPLES::
1476
1477 This gives the same answer as it would if we computed the trace
1478 from the unembedded (original) matrices::
1479
1480 sage: set_random_seed()
1481 sage: J = RealSymmetricEJA.random_instance()
1482 sage: x,y = J.random_elements(2)
1483 sage: Xe = x.to_matrix()
1484 sage: Ye = y.to_matrix()
1485 sage: X = J.real_unembed(Xe)
1486 sage: Y = J.real_unembed(Ye)
1487 sage: expected = (X*Y).trace()
1488 sage: actual = J.trace_inner_product(Xe,Ye)
1489 sage: actual == expected
1490 True
1491
1492 ::
1493
1494 sage: set_random_seed()
1495 sage: J = ComplexHermitianEJA.random_instance()
1496 sage: x,y = J.random_elements(2)
1497 sage: Xe = x.to_matrix()
1498 sage: Ye = y.to_matrix()
1499 sage: X = J.real_unembed(Xe)
1500 sage: Y = J.real_unembed(Ye)
1501 sage: expected = (X*Y).trace().real()
1502 sage: actual = J.trace_inner_product(Xe,Ye)
1503 sage: actual == expected
1504 True
1505
1506 ::
1507
1508 sage: set_random_seed()
1509 sage: J = QuaternionHermitianEJA.random_instance()
1510 sage: x,y = J.random_elements(2)
1511 sage: Xe = x.to_matrix()
1512 sage: Ye = y.to_matrix()
1513 sage: X = J.real_unembed(Xe)
1514 sage: Y = J.real_unembed(Ye)
1515 sage: expected = (X*Y).trace().coefficient_tuple()[0]
1516 sage: actual = J.trace_inner_product(Xe,Ye)
1517 sage: actual == expected
1518 True
1519
1520 """
1521 Xu = cls.real_unembed(X)
1522 Yu = cls.real_unembed(Y)
1523 tr = (Xu*Yu).trace()
1524
1525 try:
1526 # Works in QQ, AA, RDF, et cetera.
1527 return tr.real()
1528 except AttributeError:
1529 # A quaternion doesn't have a real() method, but does
1530 # have coefficient_tuple() method that returns the
1531 # coefficients of 1, i, j, and k -- in that order.
1532 return tr.coefficient_tuple()[0]
1533
1534
1535 class RealMatrixEJA(MatrixEJA):
1536 @staticmethod
1537 def dimension_over_reals():
1538 return 1
1539
1540
1541 class RealSymmetricEJA(ConcreteEJA, RealMatrixEJA):
1542 """
1543 The rank-n simple EJA consisting of real symmetric n-by-n
1544 matrices, the usual symmetric Jordan product, and the trace inner
1545 product. It has dimension `(n^2 + n)/2` over the reals.
1546
1547 SETUP::
1548
1549 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1550
1551 EXAMPLES::
1552
1553 sage: J = RealSymmetricEJA(2)
1554 sage: e0, e1, e2 = J.gens()
1555 sage: e0*e0
1556 e0
1557 sage: e1*e1
1558 1/2*e0 + 1/2*e2
1559 sage: e2*e2
1560 e2
1561
1562 In theory, our "field" can be any subfield of the reals::
1563
1564 sage: RealSymmetricEJA(2, field=RDF)
1565 Euclidean Jordan algebra of dimension 3 over Real Double Field
1566 sage: RealSymmetricEJA(2, field=RR)
1567 Euclidean Jordan algebra of dimension 3 over Real Field with
1568 53 bits of precision
1569
1570 TESTS:
1571
1572 The dimension of this algebra is `(n^2 + n) / 2`::
1573
1574 sage: set_random_seed()
1575 sage: n_max = RealSymmetricEJA._max_random_instance_size()
1576 sage: n = ZZ.random_element(1, n_max)
1577 sage: J = RealSymmetricEJA(n)
1578 sage: J.dimension() == (n^2 + n)/2
1579 True
1580
1581 The Jordan multiplication is what we think it is::
1582
1583 sage: set_random_seed()
1584 sage: J = RealSymmetricEJA.random_instance()
1585 sage: x,y = J.random_elements(2)
1586 sage: actual = (x*y).to_matrix()
1587 sage: X = x.to_matrix()
1588 sage: Y = y.to_matrix()
1589 sage: expected = (X*Y + Y*X)/2
1590 sage: actual == expected
1591 True
1592 sage: J(expected) == x*y
1593 True
1594
1595 We can change the generator prefix::
1596
1597 sage: RealSymmetricEJA(3, prefix='q').gens()
1598 (q0, q1, q2, q3, q4, q5)
1599
1600 We can construct the (trivial) algebra of rank zero::
1601
1602 sage: RealSymmetricEJA(0)
1603 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
1604
1605 """
1606 @classmethod
1607 def _denormalized_basis(cls, n):
1608 """
1609 Return a basis for the space of real symmetric n-by-n matrices.
1610
1611 SETUP::
1612
1613 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1614
1615 TESTS::
1616
1617 sage: set_random_seed()
1618 sage: n = ZZ.random_element(1,5)
1619 sage: B = RealSymmetricEJA._denormalized_basis(n)
1620 sage: all( M.is_symmetric() for M in B)
1621 True
1622
1623 """
1624 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1625 # coordinates.
1626 S = []
1627 for i in range(n):
1628 for j in range(i+1):
1629 Eij = matrix(ZZ, n, lambda k,l: k==i and l==j)
1630 if i == j:
1631 Sij = Eij
1632 else:
1633 Sij = Eij + Eij.transpose()
1634 S.append(Sij)
1635 return tuple(S)
1636
1637
1638 @staticmethod
1639 def _max_random_instance_size():
1640 return 4 # Dimension 10
1641
1642 @classmethod
1643 def random_instance(cls, **kwargs):
1644 """
1645 Return a random instance of this type of algebra.
1646 """
1647 n = ZZ.random_element(cls._max_random_instance_size() + 1)
1648 return cls(n, **kwargs)
1649
1650 def __init__(self, n, **kwargs):
1651 # We know this is a valid EJA, but will double-check
1652 # if the user passes check_axioms=True.
1653 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
1654
1655 super(RealSymmetricEJA, self).__init__(self._denormalized_basis(n),
1656 self.jordan_product,
1657 self.trace_inner_product,
1658 **kwargs)
1659
1660 # TODO: this could be factored out somehow, but is left here
1661 # because the MatrixEJA is not presently a subclass of the
1662 # FDEJA class that defines rank() and one().
1663 self.rank.set_cache(n)
1664 idV = matrix.identity(ZZ, self.dimension_over_reals()*n)
1665 self.one.set_cache(self(idV))
1666
1667
1668
1669 class ComplexMatrixEJA(MatrixEJA):
1670 # A manual dictionary-cache for the complex_extension() method,
1671 # since apparently @classmethods can't also be @cached_methods.
1672 _complex_extension = {}
1673
1674 @classmethod
1675 def complex_extension(cls,field):
1676 r"""
1677 The complex field that we embed/unembed, as an extension
1678 of the given ``field``.
1679 """
1680 if field in cls._complex_extension:
1681 return cls._complex_extension[field]
1682
1683 # Sage doesn't know how to adjoin the complex "i" (the root of
1684 # x^2 + 1) to a field in a general way. Here, we just enumerate
1685 # all of the cases that I have cared to support so far.
1686 if field is AA:
1687 # Sage doesn't know how to embed AA into QQbar, i.e. how
1688 # to adjoin sqrt(-1) to AA.
1689 F = QQbar
1690 elif not field.is_exact():
1691 # RDF or RR
1692 F = field.complex_field()
1693 else:
1694 # Works for QQ and... maybe some other fields.
1695 R = PolynomialRing(field, 'z')
1696 z = R.gen()
1697 F = field.extension(z**2 + 1, 'I', embedding=CLF(-1).sqrt())
1698
1699 cls._complex_extension[field] = F
1700 return F
1701
1702 @staticmethod
1703 def dimension_over_reals():
1704 return 2
1705
1706 @classmethod
1707 def real_embed(cls,M):
1708 """
1709 Embed the n-by-n complex matrix ``M`` into the space of real
1710 matrices of size 2n-by-2n via the map the sends each entry `z = a +
1711 bi` to the block matrix ``[[a,b],[-b,a]]``.
1712
1713 SETUP::
1714
1715 sage: from mjo.eja.eja_algebra import ComplexMatrixEJA
1716
1717 EXAMPLES::
1718
1719 sage: F = QuadraticField(-1, 'I')
1720 sage: x1 = F(4 - 2*i)
1721 sage: x2 = F(1 + 2*i)
1722 sage: x3 = F(-i)
1723 sage: x4 = F(6)
1724 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
1725 sage: ComplexMatrixEJA.real_embed(M)
1726 [ 4 -2| 1 2]
1727 [ 2 4|-2 1]
1728 [-----+-----]
1729 [ 0 -1| 6 0]
1730 [ 1 0| 0 6]
1731
1732 TESTS:
1733
1734 Embedding is a homomorphism (isomorphism, in fact)::
1735
1736 sage: set_random_seed()
1737 sage: n = ZZ.random_element(3)
1738 sage: F = QuadraticField(-1, 'I')
1739 sage: X = random_matrix(F, n)
1740 sage: Y = random_matrix(F, n)
1741 sage: Xe = ComplexMatrixEJA.real_embed(X)
1742 sage: Ye = ComplexMatrixEJA.real_embed(Y)
1743 sage: XYe = ComplexMatrixEJA.real_embed(X*Y)
1744 sage: Xe*Ye == XYe
1745 True
1746
1747 """
1748 super(ComplexMatrixEJA,cls).real_embed(M)
1749 n = M.nrows()
1750
1751 # We don't need any adjoined elements...
1752 field = M.base_ring().base_ring()
1753
1754 blocks = []
1755 for z in M.list():
1756 a = z.real()
1757 b = z.imag()
1758 blocks.append(matrix(field, 2, [ [ a, b],
1759 [-b, a] ]))
1760
1761 return matrix.block(field, n, blocks)
1762
1763
1764 @classmethod
1765 def real_unembed(cls,M):
1766 """
1767 The inverse of _embed_complex_matrix().
1768
1769 SETUP::
1770
1771 sage: from mjo.eja.eja_algebra import ComplexMatrixEJA
1772
1773 EXAMPLES::
1774
1775 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1776 ....: [-2, 1, -4, 3],
1777 ....: [ 9, 10, 11, 12],
1778 ....: [-10, 9, -12, 11] ])
1779 sage: ComplexMatrixEJA.real_unembed(A)
1780 [ 2*I + 1 4*I + 3]
1781 [ 10*I + 9 12*I + 11]
1782
1783 TESTS:
1784
1785 Unembedding is the inverse of embedding::
1786
1787 sage: set_random_seed()
1788 sage: F = QuadraticField(-1, 'I')
1789 sage: M = random_matrix(F, 3)
1790 sage: Me = ComplexMatrixEJA.real_embed(M)
1791 sage: ComplexMatrixEJA.real_unembed(Me) == M
1792 True
1793
1794 """
1795 super(ComplexMatrixEJA,cls).real_unembed(M)
1796 n = ZZ(M.nrows())
1797 d = cls.dimension_over_reals()
1798 F = cls.complex_extension(M.base_ring())
1799 i = F.gen()
1800
1801 # Go top-left to bottom-right (reading order), converting every
1802 # 2-by-2 block we see to a single complex element.
1803 elements = []
1804 for k in range(n/d):
1805 for j in range(n/d):
1806 submat = M[d*k:d*k+d,d*j:d*j+d]
1807 if submat[0,0] != submat[1,1]:
1808 raise ValueError('bad on-diagonal submatrix')
1809 if submat[0,1] != -submat[1,0]:
1810 raise ValueError('bad off-diagonal submatrix')
1811 z = submat[0,0] + submat[0,1]*i
1812 elements.append(z)
1813
1814 return matrix(F, n/d, elements)
1815
1816
1817 class ComplexHermitianEJA(ConcreteEJA, ComplexMatrixEJA):
1818 """
1819 The rank-n simple EJA consisting of complex Hermitian n-by-n
1820 matrices over the real numbers, the usual symmetric Jordan product,
1821 and the real-part-of-trace inner product. It has dimension `n^2` over
1822 the reals.
1823
1824 SETUP::
1825
1826 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1827
1828 EXAMPLES:
1829
1830 In theory, our "field" can be any subfield of the reals::
1831
1832 sage: ComplexHermitianEJA(2, field=RDF)
1833 Euclidean Jordan algebra of dimension 4 over Real Double Field
1834 sage: ComplexHermitianEJA(2, field=RR)
1835 Euclidean Jordan algebra of dimension 4 over Real Field with
1836 53 bits of precision
1837
1838 TESTS:
1839
1840 The dimension of this algebra is `n^2`::
1841
1842 sage: set_random_seed()
1843 sage: n_max = ComplexHermitianEJA._max_random_instance_size()
1844 sage: n = ZZ.random_element(1, n_max)
1845 sage: J = ComplexHermitianEJA(n)
1846 sage: J.dimension() == n^2
1847 True
1848
1849 The Jordan multiplication is what we think it is::
1850
1851 sage: set_random_seed()
1852 sage: J = ComplexHermitianEJA.random_instance()
1853 sage: x,y = J.random_elements(2)
1854 sage: actual = (x*y).to_matrix()
1855 sage: X = x.to_matrix()
1856 sage: Y = y.to_matrix()
1857 sage: expected = (X*Y + Y*X)/2
1858 sage: actual == expected
1859 True
1860 sage: J(expected) == x*y
1861 True
1862
1863 We can change the generator prefix::
1864
1865 sage: ComplexHermitianEJA(2, prefix='z').gens()
1866 (z0, z1, z2, z3)
1867
1868 We can construct the (trivial) algebra of rank zero::
1869
1870 sage: ComplexHermitianEJA(0)
1871 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
1872
1873 """
1874
1875 @classmethod
1876 def _denormalized_basis(cls, n):
1877 """
1878 Returns a basis for the space of complex Hermitian n-by-n matrices.
1879
1880 Why do we embed these? Basically, because all of numerical linear
1881 algebra assumes that you're working with vectors consisting of `n`
1882 entries from a field and scalars from the same field. There's no way
1883 to tell SageMath that (for example) the vectors contain complex
1884 numbers, while the scalar field is real.
1885
1886 SETUP::
1887
1888 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1889
1890 TESTS::
1891
1892 sage: set_random_seed()
1893 sage: n = ZZ.random_element(1,5)
1894 sage: B = ComplexHermitianEJA._denormalized_basis(n)
1895 sage: all( M.is_symmetric() for M in B)
1896 True
1897
1898 """
1899 field = ZZ
1900 R = PolynomialRing(field, 'z')
1901 z = R.gen()
1902 F = field.extension(z**2 + 1, 'I')
1903 I = F.gen(1)
1904
1905 # This is like the symmetric case, but we need to be careful:
1906 #
1907 # * We want conjugate-symmetry, not just symmetry.
1908 # * The diagonal will (as a result) be real.
1909 #
1910 S = []
1911 Eij = matrix.zero(F,n)
1912 for i in range(n):
1913 for j in range(i+1):
1914 # "build" E_ij
1915 Eij[i,j] = 1
1916 if i == j:
1917 Sij = cls.real_embed(Eij)
1918 S.append(Sij)
1919 else:
1920 # The second one has a minus because it's conjugated.
1921 Eij[j,i] = 1 # Eij = Eij + Eij.transpose()
1922 Sij_real = cls.real_embed(Eij)
1923 S.append(Sij_real)
1924 # Eij = I*Eij - I*Eij.transpose()
1925 Eij[i,j] = I
1926 Eij[j,i] = -I
1927 Sij_imag = cls.real_embed(Eij)
1928 S.append(Sij_imag)
1929 Eij[j,i] = 0
1930 # "erase" E_ij
1931 Eij[i,j] = 0
1932
1933 # Since we embedded these, we can drop back to the "field" that we
1934 # started with instead of the complex extension "F".
1935 return tuple( s.change_ring(field) for s in S )
1936
1937
1938 def __init__(self, n, **kwargs):
1939 # We know this is a valid EJA, but will double-check
1940 # if the user passes check_axioms=True.
1941 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
1942
1943 super(ComplexHermitianEJA, self).__init__(self._denormalized_basis(n),
1944 self.jordan_product,
1945 self.trace_inner_product,
1946 **kwargs)
1947 # TODO: this could be factored out somehow, but is left here
1948 # because the MatrixEJA is not presently a subclass of the
1949 # FDEJA class that defines rank() and one().
1950 self.rank.set_cache(n)
1951 idV = matrix.identity(ZZ, self.dimension_over_reals()*n)
1952 self.one.set_cache(self(idV))
1953
1954 @staticmethod
1955 def _max_random_instance_size():
1956 return 3 # Dimension 9
1957
1958 @classmethod
1959 def random_instance(cls, **kwargs):
1960 """
1961 Return a random instance of this type of algebra.
1962 """
1963 n = ZZ.random_element(cls._max_random_instance_size() + 1)
1964 return cls(n, **kwargs)
1965
1966 class QuaternionMatrixEJA(MatrixEJA):
1967
1968 # A manual dictionary-cache for the quaternion_extension() method,
1969 # since apparently @classmethods can't also be @cached_methods.
1970 _quaternion_extension = {}
1971
1972 @classmethod
1973 def quaternion_extension(cls,field):
1974 r"""
1975 The quaternion field that we embed/unembed, as an extension
1976 of the given ``field``.
1977 """
1978 if field in cls._quaternion_extension:
1979 return cls._quaternion_extension[field]
1980
1981 Q = QuaternionAlgebra(field,-1,-1)
1982
1983 cls._quaternion_extension[field] = Q
1984 return Q
1985
1986 @staticmethod
1987 def dimension_over_reals():
1988 return 4
1989
1990 @classmethod
1991 def real_embed(cls,M):
1992 """
1993 Embed the n-by-n quaternion matrix ``M`` into the space of real
1994 matrices of size 4n-by-4n by first sending each quaternion entry `z
1995 = a + bi + cj + dk` to the block-complex matrix ``[[a + bi,
1996 c+di],[-c + di, a-bi]]`, and then embedding those into a real
1997 matrix.
1998
1999 SETUP::
2000
2001 sage: from mjo.eja.eja_algebra import QuaternionMatrixEJA
2002
2003 EXAMPLES::
2004
2005 sage: Q = QuaternionAlgebra(QQ,-1,-1)
2006 sage: i,j,k = Q.gens()
2007 sage: x = 1 + 2*i + 3*j + 4*k
2008 sage: M = matrix(Q, 1, [[x]])
2009 sage: QuaternionMatrixEJA.real_embed(M)
2010 [ 1 2 3 4]
2011 [-2 1 -4 3]
2012 [-3 4 1 -2]
2013 [-4 -3 2 1]
2014
2015 Embedding is a homomorphism (isomorphism, in fact)::
2016
2017 sage: set_random_seed()
2018 sage: n = ZZ.random_element(2)
2019 sage: Q = QuaternionAlgebra(QQ,-1,-1)
2020 sage: X = random_matrix(Q, n)
2021 sage: Y = random_matrix(Q, n)
2022 sage: Xe = QuaternionMatrixEJA.real_embed(X)
2023 sage: Ye = QuaternionMatrixEJA.real_embed(Y)
2024 sage: XYe = QuaternionMatrixEJA.real_embed(X*Y)
2025 sage: Xe*Ye == XYe
2026 True
2027
2028 """
2029 super(QuaternionMatrixEJA,cls).real_embed(M)
2030 quaternions = M.base_ring()
2031 n = M.nrows()
2032
2033 F = QuadraticField(-1, 'I')
2034 i = F.gen()
2035
2036 blocks = []
2037 for z in M.list():
2038 t = z.coefficient_tuple()
2039 a = t[0]
2040 b = t[1]
2041 c = t[2]
2042 d = t[3]
2043 cplxM = matrix(F, 2, [[ a + b*i, c + d*i],
2044 [-c + d*i, a - b*i]])
2045 realM = ComplexMatrixEJA.real_embed(cplxM)
2046 blocks.append(realM)
2047
2048 # We should have real entries by now, so use the realest field
2049 # we've got for the return value.
2050 return matrix.block(quaternions.base_ring(), n, blocks)
2051
2052
2053
2054 @classmethod
2055 def real_unembed(cls,M):
2056 """
2057 The inverse of _embed_quaternion_matrix().
2058
2059 SETUP::
2060
2061 sage: from mjo.eja.eja_algebra import QuaternionMatrixEJA
2062
2063 EXAMPLES::
2064
2065 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
2066 ....: [-2, 1, -4, 3],
2067 ....: [-3, 4, 1, -2],
2068 ....: [-4, -3, 2, 1]])
2069 sage: QuaternionMatrixEJA.real_unembed(M)
2070 [1 + 2*i + 3*j + 4*k]
2071
2072 TESTS:
2073
2074 Unembedding is the inverse of embedding::
2075
2076 sage: set_random_seed()
2077 sage: Q = QuaternionAlgebra(QQ, -1, -1)
2078 sage: M = random_matrix(Q, 3)
2079 sage: Me = QuaternionMatrixEJA.real_embed(M)
2080 sage: QuaternionMatrixEJA.real_unembed(Me) == M
2081 True
2082
2083 """
2084 super(QuaternionMatrixEJA,cls).real_unembed(M)
2085 n = ZZ(M.nrows())
2086 d = cls.dimension_over_reals()
2087
2088 # Use the base ring of the matrix to ensure that its entries can be
2089 # multiplied by elements of the quaternion algebra.
2090 Q = cls.quaternion_extension(M.base_ring())
2091 i,j,k = Q.gens()
2092
2093 # Go top-left to bottom-right (reading order), converting every
2094 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
2095 # quaternion block.
2096 elements = []
2097 for l in range(n/d):
2098 for m in range(n/d):
2099 submat = ComplexMatrixEJA.real_unembed(
2100 M[d*l:d*l+d,d*m:d*m+d] )
2101 if submat[0,0] != submat[1,1].conjugate():
2102 raise ValueError('bad on-diagonal submatrix')
2103 if submat[0,1] != -submat[1,0].conjugate():
2104 raise ValueError('bad off-diagonal submatrix')
2105 z = submat[0,0].real()
2106 z += submat[0,0].imag()*i
2107 z += submat[0,1].real()*j
2108 z += submat[0,1].imag()*k
2109 elements.append(z)
2110
2111 return matrix(Q, n/d, elements)
2112
2113
2114 class QuaternionHermitianEJA(ConcreteEJA, QuaternionMatrixEJA):
2115 r"""
2116 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
2117 matrices, the usual symmetric Jordan product, and the
2118 real-part-of-trace inner product. It has dimension `2n^2 - n` over
2119 the reals.
2120
2121 SETUP::
2122
2123 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
2124
2125 EXAMPLES:
2126
2127 In theory, our "field" can be any subfield of the reals::
2128
2129 sage: QuaternionHermitianEJA(2, field=RDF)
2130 Euclidean Jordan algebra of dimension 6 over Real Double Field
2131 sage: QuaternionHermitianEJA(2, field=RR)
2132 Euclidean Jordan algebra of dimension 6 over Real Field with
2133 53 bits of precision
2134
2135 TESTS:
2136
2137 The dimension of this algebra is `2*n^2 - n`::
2138
2139 sage: set_random_seed()
2140 sage: n_max = QuaternionHermitianEJA._max_random_instance_size()
2141 sage: n = ZZ.random_element(1, n_max)
2142 sage: J = QuaternionHermitianEJA(n)
2143 sage: J.dimension() == 2*(n^2) - n
2144 True
2145
2146 The Jordan multiplication is what we think it is::
2147
2148 sage: set_random_seed()
2149 sage: J = QuaternionHermitianEJA.random_instance()
2150 sage: x,y = J.random_elements(2)
2151 sage: actual = (x*y).to_matrix()
2152 sage: X = x.to_matrix()
2153 sage: Y = y.to_matrix()
2154 sage: expected = (X*Y + Y*X)/2
2155 sage: actual == expected
2156 True
2157 sage: J(expected) == x*y
2158 True
2159
2160 We can change the generator prefix::
2161
2162 sage: QuaternionHermitianEJA(2, prefix='a').gens()
2163 (a0, a1, a2, a3, a4, a5)
2164
2165 We can construct the (trivial) algebra of rank zero::
2166
2167 sage: QuaternionHermitianEJA(0)
2168 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
2169
2170 """
2171 @classmethod
2172 def _denormalized_basis(cls, n):
2173 """
2174 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
2175
2176 Why do we embed these? Basically, because all of numerical
2177 linear algebra assumes that you're working with vectors consisting
2178 of `n` entries from a field and scalars from the same field. There's
2179 no way to tell SageMath that (for example) the vectors contain
2180 complex numbers, while the scalar field is real.
2181
2182 SETUP::
2183
2184 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
2185
2186 TESTS::
2187
2188 sage: set_random_seed()
2189 sage: n = ZZ.random_element(1,5)
2190 sage: B = QuaternionHermitianEJA._denormalized_basis(n)
2191 sage: all( M.is_symmetric() for M in B )
2192 True
2193
2194 """
2195 field = ZZ
2196 Q = QuaternionAlgebra(QQ,-1,-1)
2197 I,J,K = Q.gens()
2198
2199 # This is like the symmetric case, but we need to be careful:
2200 #
2201 # * We want conjugate-symmetry, not just symmetry.
2202 # * The diagonal will (as a result) be real.
2203 #
2204 S = []
2205 Eij = matrix.zero(Q,n)
2206 for i in range(n):
2207 for j in range(i+1):
2208 # "build" E_ij
2209 Eij[i,j] = 1
2210 if i == j:
2211 Sij = cls.real_embed(Eij)
2212 S.append(Sij)
2213 else:
2214 # The second, third, and fourth ones have a minus
2215 # because they're conjugated.
2216 # Eij = Eij + Eij.transpose()
2217 Eij[j,i] = 1
2218 Sij_real = cls.real_embed(Eij)
2219 S.append(Sij_real)
2220 # Eij = I*(Eij - Eij.transpose())
2221 Eij[i,j] = I
2222 Eij[j,i] = -I
2223 Sij_I = cls.real_embed(Eij)
2224 S.append(Sij_I)
2225 # Eij = J*(Eij - Eij.transpose())
2226 Eij[i,j] = J
2227 Eij[j,i] = -J
2228 Sij_J = cls.real_embed(Eij)
2229 S.append(Sij_J)
2230 # Eij = K*(Eij - Eij.transpose())
2231 Eij[i,j] = K
2232 Eij[j,i] = -K
2233 Sij_K = cls.real_embed(Eij)
2234 S.append(Sij_K)
2235 Eij[j,i] = 0
2236 # "erase" E_ij
2237 Eij[i,j] = 0
2238
2239 # Since we embedded these, we can drop back to the "field" that we
2240 # started with instead of the quaternion algebra "Q".
2241 return tuple( s.change_ring(field) for s in S )
2242
2243
2244 def __init__(self, n, **kwargs):
2245 # We know this is a valid EJA, but will double-check
2246 # if the user passes check_axioms=True.
2247 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
2248
2249 super(QuaternionHermitianEJA, self).__init__(self._denormalized_basis(n),
2250 self.jordan_product,
2251 self.trace_inner_product,
2252 **kwargs)
2253 # TODO: this could be factored out somehow, but is left here
2254 # because the MatrixEJA is not presently a subclass of the
2255 # FDEJA class that defines rank() and one().
2256 self.rank.set_cache(n)
2257 idV = matrix.identity(ZZ, self.dimension_over_reals()*n)
2258 self.one.set_cache(self(idV))
2259
2260
2261 @staticmethod
2262 def _max_random_instance_size():
2263 r"""
2264 The maximum rank of a random QuaternionHermitianEJA.
2265 """
2266 return 2 # Dimension 6
2267
2268 @classmethod
2269 def random_instance(cls, **kwargs):
2270 """
2271 Return a random instance of this type of algebra.
2272 """
2273 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2274 return cls(n, **kwargs)
2275
2276
2277 class HadamardEJA(ConcreteEJA):
2278 """
2279 Return the Euclidean Jordan Algebra corresponding to the set
2280 `R^n` under the Hadamard product.
2281
2282 Note: this is nothing more than the Cartesian product of ``n``
2283 copies of the spin algebra. Once Cartesian product algebras
2284 are implemented, this can go.
2285
2286 SETUP::
2287
2288 sage: from mjo.eja.eja_algebra import HadamardEJA
2289
2290 EXAMPLES:
2291
2292 This multiplication table can be verified by hand::
2293
2294 sage: J = HadamardEJA(3)
2295 sage: e0,e1,e2 = J.gens()
2296 sage: e0*e0
2297 e0
2298 sage: e0*e1
2299 0
2300 sage: e0*e2
2301 0
2302 sage: e1*e1
2303 e1
2304 sage: e1*e2
2305 0
2306 sage: e2*e2
2307 e2
2308
2309 TESTS:
2310
2311 We can change the generator prefix::
2312
2313 sage: HadamardEJA(3, prefix='r').gens()
2314 (r0, r1, r2)
2315
2316 """
2317 def __init__(self, n, **kwargs):
2318 if n == 0:
2319 jordan_product = lambda x,y: x
2320 inner_product = lambda x,y: x
2321 else:
2322 def jordan_product(x,y):
2323 P = x.parent()
2324 return P( xi*yi for (xi,yi) in zip(x,y) )
2325
2326 def inner_product(x,y):
2327 return (x.T*y)[0,0]
2328
2329 # New defaults for keyword arguments. Don't orthonormalize
2330 # because our basis is already orthonormal with respect to our
2331 # inner-product. Don't check the axioms, because we know this
2332 # is a valid EJA... but do double-check if the user passes
2333 # check_axioms=True. Note: we DON'T override the "check_field"
2334 # default here, because the user can pass in a field!
2335 if "orthonormalize" not in kwargs: kwargs["orthonormalize"] = False
2336 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
2337
2338 column_basis = tuple( b.column() for b in FreeModule(ZZ, n).basis() )
2339 super().__init__(column_basis, jordan_product, inner_product, **kwargs)
2340 self.rank.set_cache(n)
2341
2342 if n == 0:
2343 self.one.set_cache( self.zero() )
2344 else:
2345 self.one.set_cache( sum(self.gens()) )
2346
2347 @staticmethod
2348 def _max_random_instance_size():
2349 r"""
2350 The maximum dimension of a random HadamardEJA.
2351 """
2352 return 5
2353
2354 @classmethod
2355 def random_instance(cls, **kwargs):
2356 """
2357 Return a random instance of this type of algebra.
2358 """
2359 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2360 return cls(n, **kwargs)
2361
2362
2363 class BilinearFormEJA(ConcreteEJA):
2364 r"""
2365 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2366 with the half-trace inner product and jordan product ``x*y =
2367 (<Bx,y>,y_bar>, x0*y_bar + y0*x_bar)`` where `B = 1 \times B22` is
2368 a symmetric positive-definite "bilinear form" matrix. Its
2369 dimension is the size of `B`, and it has rank two in dimensions
2370 larger than two. It reduces to the ``JordanSpinEJA`` when `B` is
2371 the identity matrix of order ``n``.
2372
2373 We insist that the one-by-one upper-left identity block of `B` be
2374 passed in as well so that we can be passed a matrix of size zero
2375 to construct a trivial algebra.
2376
2377 SETUP::
2378
2379 sage: from mjo.eja.eja_algebra import (BilinearFormEJA,
2380 ....: JordanSpinEJA)
2381
2382 EXAMPLES:
2383
2384 When no bilinear form is specified, the identity matrix is used,
2385 and the resulting algebra is the Jordan spin algebra::
2386
2387 sage: B = matrix.identity(AA,3)
2388 sage: J0 = BilinearFormEJA(B)
2389 sage: J1 = JordanSpinEJA(3)
2390 sage: J0.multiplication_table() == J0.multiplication_table()
2391 True
2392
2393 An error is raised if the matrix `B` does not correspond to a
2394 positive-definite bilinear form::
2395
2396 sage: B = matrix.random(QQ,2,3)
2397 sage: J = BilinearFormEJA(B)
2398 Traceback (most recent call last):
2399 ...
2400 ValueError: bilinear form is not positive-definite
2401 sage: B = matrix.zero(QQ,3)
2402 sage: J = BilinearFormEJA(B)
2403 Traceback (most recent call last):
2404 ...
2405 ValueError: bilinear form is not positive-definite
2406
2407 TESTS:
2408
2409 We can create a zero-dimensional algebra::
2410
2411 sage: B = matrix.identity(AA,0)
2412 sage: J = BilinearFormEJA(B)
2413 sage: J.basis()
2414 Finite family {}
2415
2416 We can check the multiplication condition given in the Jordan, von
2417 Neumann, and Wigner paper (and also discussed on my "On the
2418 symmetry..." paper). Note that this relies heavily on the standard
2419 choice of basis, as does anything utilizing the bilinear form
2420 matrix. We opt not to orthonormalize the basis, because if we
2421 did, we would have to normalize the `s_{i}` in a similar manner::
2422
2423 sage: set_random_seed()
2424 sage: n = ZZ.random_element(5)
2425 sage: M = matrix.random(QQ, max(0,n-1), algorithm='unimodular')
2426 sage: B11 = matrix.identity(QQ,1)
2427 sage: B22 = M.transpose()*M
2428 sage: B = block_matrix(2,2,[ [B11,0 ],
2429 ....: [0, B22 ] ])
2430 sage: J = BilinearFormEJA(B, orthonormalize=False)
2431 sage: eis = VectorSpace(M.base_ring(), M.ncols()).basis()
2432 sage: V = J.vector_space()
2433 sage: sis = [ J( V([0] + (M.inverse()*ei).list()).column() )
2434 ....: for ei in eis ]
2435 sage: actual = [ sis[i]*sis[j]
2436 ....: for i in range(n-1)
2437 ....: for j in range(n-1) ]
2438 sage: expected = [ J.one() if i == j else J.zero()
2439 ....: for i in range(n-1)
2440 ....: for j in range(n-1) ]
2441 sage: actual == expected
2442 True
2443
2444 """
2445 def __init__(self, B, **kwargs):
2446 # The matrix "B" is supplied by the user in most cases,
2447 # so it makes sense to check whether or not its positive-
2448 # definite unless we are specifically asked not to...
2449 if ("check_axioms" not in kwargs) or kwargs["check_axioms"]:
2450 if not B.is_positive_definite():
2451 raise ValueError("bilinear form is not positive-definite")
2452
2453 # However, all of the other data for this EJA is computed
2454 # by us in manner that guarantees the axioms are
2455 # satisfied. So, again, unless we are specifically asked to
2456 # verify things, we'll skip the rest of the checks.
2457 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
2458
2459 def inner_product(x,y):
2460 return (y.T*B*x)[0,0]
2461
2462 def jordan_product(x,y):
2463 P = x.parent()
2464 x0 = x[0,0]
2465 xbar = x[1:,0]
2466 y0 = y[0,0]
2467 ybar = y[1:,0]
2468 z0 = inner_product(y,x)
2469 zbar = y0*xbar + x0*ybar
2470 return P([z0] + zbar.list())
2471
2472 n = B.nrows()
2473 column_basis = tuple( b.column() for b in FreeModule(ZZ, n).basis() )
2474 super(BilinearFormEJA, self).__init__(column_basis,
2475 jordan_product,
2476 inner_product,
2477 **kwargs)
2478
2479 # The rank of this algebra is two, unless we're in a
2480 # one-dimensional ambient space (because the rank is bounded
2481 # by the ambient dimension).
2482 self.rank.set_cache(min(n,2))
2483
2484 if n == 0:
2485 self.one.set_cache( self.zero() )
2486 else:
2487 self.one.set_cache( self.monomial(0) )
2488
2489 @staticmethod
2490 def _max_random_instance_size():
2491 r"""
2492 The maximum dimension of a random BilinearFormEJA.
2493 """
2494 return 5
2495
2496 @classmethod
2497 def random_instance(cls, **kwargs):
2498 """
2499 Return a random instance of this algebra.
2500 """
2501 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2502 if n.is_zero():
2503 B = matrix.identity(ZZ, n)
2504 return cls(B, **kwargs)
2505
2506 B11 = matrix.identity(ZZ, 1)
2507 M = matrix.random(ZZ, n-1)
2508 I = matrix.identity(ZZ, n-1)
2509 alpha = ZZ.zero()
2510 while alpha.is_zero():
2511 alpha = ZZ.random_element().abs()
2512 B22 = M.transpose()*M + alpha*I
2513
2514 from sage.matrix.special import block_matrix
2515 B = block_matrix(2,2, [ [B11, ZZ(0) ],
2516 [ZZ(0), B22 ] ])
2517
2518 return cls(B, **kwargs)
2519
2520
2521 class JordanSpinEJA(BilinearFormEJA):
2522 """
2523 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2524 with the usual inner product and jordan product ``x*y =
2525 (<x,y>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
2526 the reals.
2527
2528 SETUP::
2529
2530 sage: from mjo.eja.eja_algebra import JordanSpinEJA
2531
2532 EXAMPLES:
2533
2534 This multiplication table can be verified by hand::
2535
2536 sage: J = JordanSpinEJA(4)
2537 sage: e0,e1,e2,e3 = J.gens()
2538 sage: e0*e0
2539 e0
2540 sage: e0*e1
2541 e1
2542 sage: e0*e2
2543 e2
2544 sage: e0*e3
2545 e3
2546 sage: e1*e2
2547 0
2548 sage: e1*e3
2549 0
2550 sage: e2*e3
2551 0
2552
2553 We can change the generator prefix::
2554
2555 sage: JordanSpinEJA(2, prefix='B').gens()
2556 (B0, B1)
2557
2558 TESTS:
2559
2560 Ensure that we have the usual inner product on `R^n`::
2561
2562 sage: set_random_seed()
2563 sage: J = JordanSpinEJA.random_instance()
2564 sage: x,y = J.random_elements(2)
2565 sage: actual = x.inner_product(y)
2566 sage: expected = x.to_vector().inner_product(y.to_vector())
2567 sage: actual == expected
2568 True
2569
2570 """
2571 def __init__(self, n, **kwargs):
2572 # This is a special case of the BilinearFormEJA with the
2573 # identity matrix as its bilinear form.
2574 B = matrix.identity(ZZ, n)
2575
2576 # Don't orthonormalize because our basis is already
2577 # orthonormal with respect to our inner-product.
2578 if "orthonormalize" not in kwargs: kwargs["orthonormalize"] = False
2579
2580 # But also don't pass check_field=False here, because the user
2581 # can pass in a field!
2582 super(JordanSpinEJA, self).__init__(B, **kwargs)
2583
2584 @staticmethod
2585 def _max_random_instance_size():
2586 r"""
2587 The maximum dimension of a random JordanSpinEJA.
2588 """
2589 return 5
2590
2591 @classmethod
2592 def random_instance(cls, **kwargs):
2593 """
2594 Return a random instance of this type of algebra.
2595
2596 Needed here to override the implementation for ``BilinearFormEJA``.
2597 """
2598 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2599 return cls(n, **kwargs)
2600
2601
2602 class TrivialEJA(ConcreteEJA):
2603 """
2604 The trivial Euclidean Jordan algebra consisting of only a zero element.
2605
2606 SETUP::
2607
2608 sage: from mjo.eja.eja_algebra import TrivialEJA
2609
2610 EXAMPLES::
2611
2612 sage: J = TrivialEJA()
2613 sage: J.dimension()
2614 0
2615 sage: J.zero()
2616 0
2617 sage: J.one()
2618 0
2619 sage: 7*J.one()*12*J.one()
2620 0
2621 sage: J.one().inner_product(J.one())
2622 0
2623 sage: J.one().norm()
2624 0
2625 sage: J.one().subalgebra_generated_by()
2626 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
2627 sage: J.rank()
2628 0
2629
2630 """
2631 def __init__(self, **kwargs):
2632 jordan_product = lambda x,y: x
2633 inner_product = lambda x,y: 0
2634 basis = ()
2635
2636 # New defaults for keyword arguments
2637 if "orthonormalize" not in kwargs: kwargs["orthonormalize"] = False
2638 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
2639
2640 super(TrivialEJA, self).__init__(basis,
2641 jordan_product,
2642 inner_product,
2643 **kwargs)
2644 # The rank is zero using my definition, namely the dimension of the
2645 # largest subalgebra generated by any element.
2646 self.rank.set_cache(0)
2647 self.one.set_cache( self.zero() )
2648
2649 @classmethod
2650 def random_instance(cls, **kwargs):
2651 # We don't take a "size" argument so the superclass method is
2652 # inappropriate for us.
2653 return cls(**kwargs)
2654
2655
2656 class DirectSumEJA(FiniteDimensionalEJA):
2657 r"""
2658 The external (orthogonal) direct sum of two other Euclidean Jordan
2659 algebras. Essentially the Cartesian product of its two factors.
2660 Every Euclidean Jordan algebra decomposes into an orthogonal
2661 direct sum of simple Euclidean Jordan algebras, so no generality
2662 is lost by providing only this construction.
2663
2664 SETUP::
2665
2666 sage: from mjo.eja.eja_algebra import (random_eja,
2667 ....: HadamardEJA,
2668 ....: RealSymmetricEJA,
2669 ....: DirectSumEJA)
2670
2671 EXAMPLES::
2672
2673 sage: J1 = HadamardEJA(2)
2674 sage: J2 = RealSymmetricEJA(3)
2675 sage: J = DirectSumEJA(J1,J2)
2676 sage: J.dimension()
2677 8
2678 sage: J.rank()
2679 5
2680 sage: J.matrix_space()
2681 The Cartesian product of (Full MatrixSpace of 2 by 1 dense matrices
2682 over Algebraic Real Field, Full MatrixSpace of 3 by 3 dense matrices
2683 over Algebraic Real Field)
2684
2685 TESTS:
2686
2687 The external direct sum construction is only valid when the two factors
2688 have the same base ring; an error is raised otherwise::
2689
2690 sage: set_random_seed()
2691 sage: J1 = random_eja(field=AA)
2692 sage: J2 = random_eja(field=QQ,orthonormalize=False)
2693 sage: J = DirectSumEJA(J1,J2)
2694 Traceback (most recent call last):
2695 ...
2696 ValueError: algebras must share the same base field
2697
2698 """
2699 def __init__(self, J1, J2, **kwargs):
2700 if J1.base_ring() != J2.base_ring():
2701 raise ValueError("algebras must share the same base field")
2702 field = J1.base_ring()
2703
2704 M = J1.matrix_space().cartesian_product(J2.matrix_space())
2705 self._cartprod_algebra = J1.cartesian_product(J2)
2706
2707 self._matrix_basis = tuple( [M((a,0)) for a in J1.matrix_basis()] +
2708 [M((0,b)) for b in J2.matrix_basis()] )
2709
2710 n = len(self._matrix_basis)
2711 self._sets = None
2712 CombinatorialFreeModule.__init__(
2713 self,
2714 field,
2715 range(n),
2716 category=self._cartprod_algebra.category(),
2717 bracket=False,
2718 **kwargs)
2719 self.rank.set_cache(J1.rank() + J2.rank())
2720
2721
2722
2723 def product(self,x,y):
2724 r"""
2725 SETUP::
2726
2727 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
2728 ....: ComplexHermitianEJA,
2729 ....: DirectSumEJA)
2730
2731 TESTS::
2732
2733 sage: set_random_seed()
2734 sage: J1 = JordanSpinEJA(3, field=QQ)
2735 sage: J2 = ComplexHermitianEJA(2, field=QQ, orthonormalize=False)
2736 sage: J = DirectSumEJA(J1,J2)
2737 sage: J.random_element()*J.random_element() in J
2738 True
2739
2740 """
2741 xv = self._cartprod_algebra.from_vector(x.to_vector())
2742 yv = self._cartprod_algebra.from_vector(y.to_vector())
2743 return self.from_vector((xv*yv).to_vector())
2744
2745
2746 def cartesian_factors(self):
2747 r"""
2748 Return the pair of this algebra's factors.
2749
2750 SETUP::
2751
2752 sage: from mjo.eja.eja_algebra import (HadamardEJA,
2753 ....: JordanSpinEJA,
2754 ....: DirectSumEJA)
2755
2756 EXAMPLES::
2757
2758 sage: J1 = HadamardEJA(2, field=QQ)
2759 sage: J2 = JordanSpinEJA(3, field=QQ)
2760 sage: J = DirectSumEJA(J1,J2)
2761 sage: J.cartesian_factors()
2762 (Euclidean Jordan algebra of dimension 2 over Rational Field,
2763 Euclidean Jordan algebra of dimension 3 over Rational Field)
2764
2765 """
2766 return self._cartprod_algebra.cartesian_factors()
2767
2768
2769 # def projections(self):
2770 # r"""
2771 # Return a pair of projections onto this algebra's factors.
2772
2773 # SETUP::
2774
2775 # sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
2776 # ....: ComplexHermitianEJA,
2777 # ....: DirectSumEJA)
2778
2779 # EXAMPLES::
2780
2781 # sage: J1 = JordanSpinEJA(2)
2782 # sage: J2 = ComplexHermitianEJA(2)
2783 # sage: J = DirectSumEJA(J1,J2)
2784 # sage: (pi_left, pi_right) = J.projections()
2785 # sage: J.one().to_vector()
2786 # (1, 0, 1, 0, 0, 1)
2787 # sage: pi_left(J.one()).to_vector()
2788 # (1, 0)
2789 # sage: pi_right(J.one()).to_vector()
2790 # (1, 0, 0, 1)
2791
2792 # """
2793 # (J1,J2) = self.factors()
2794 # m = J1.dimension()
2795 # n = J2.dimension()
2796 # V_basis = self.vector_space().basis()
2797 # # Need to specify the dimensions explicitly so that we don't
2798 # # wind up with a zero-by-zero matrix when we want e.g. a
2799 # # zero-by-two matrix (important for composing things).
2800 # P1 = matrix(self.base_ring(), m, m+n, V_basis[:m])
2801 # P2 = matrix(self.base_ring(), n, m+n, V_basis[m:])
2802 # pi_left = FiniteDimensionalEJAOperator(self,J1,P1)
2803 # pi_right = FiniteDimensionalEJAOperator(self,J2,P2)
2804 # return (pi_left, pi_right)
2805
2806 # def inclusions(self):
2807 # r"""
2808 # Return the pair of inclusion maps from our factors into us.
2809
2810 # SETUP::
2811
2812 # sage: from mjo.eja.eja_algebra import (random_eja,
2813 # ....: JordanSpinEJA,
2814 # ....: RealSymmetricEJA,
2815 # ....: DirectSumEJA)
2816
2817 # EXAMPLES::
2818
2819 # sage: J1 = JordanSpinEJA(3)
2820 # sage: J2 = RealSymmetricEJA(2)
2821 # sage: J = DirectSumEJA(J1,J2)
2822 # sage: (iota_left, iota_right) = J.inclusions()
2823 # sage: iota_left(J1.zero()) == J.zero()
2824 # True
2825 # sage: iota_right(J2.zero()) == J.zero()
2826 # True
2827 # sage: J1.one().to_vector()
2828 # (1, 0, 0)
2829 # sage: iota_left(J1.one()).to_vector()
2830 # (1, 0, 0, 0, 0, 0)
2831 # sage: J2.one().to_vector()
2832 # (1, 0, 1)
2833 # sage: iota_right(J2.one()).to_vector()
2834 # (0, 0, 0, 1, 0, 1)
2835 # sage: J.one().to_vector()
2836 # (1, 0, 0, 1, 0, 1)
2837
2838 # TESTS:
2839
2840 # Composing a projection with the corresponding inclusion should
2841 # produce the identity map, and mismatching them should produce
2842 # the zero map::
2843
2844 # sage: set_random_seed()
2845 # sage: J1 = random_eja()
2846 # sage: J2 = random_eja()
2847 # sage: J = DirectSumEJA(J1,J2)
2848 # sage: (iota_left, iota_right) = J.inclusions()
2849 # sage: (pi_left, pi_right) = J.projections()
2850 # sage: pi_left*iota_left == J1.one().operator()
2851 # True
2852 # sage: pi_right*iota_right == J2.one().operator()
2853 # True
2854 # sage: (pi_left*iota_right).is_zero()
2855 # True
2856 # sage: (pi_right*iota_left).is_zero()
2857 # True
2858
2859 # """
2860 # (J1,J2) = self.factors()
2861 # m = J1.dimension()
2862 # n = J2.dimension()
2863 # V_basis = self.vector_space().basis()
2864 # # Need to specify the dimensions explicitly so that we don't
2865 # # wind up with a zero-by-zero matrix when we want e.g. a
2866 # # two-by-zero matrix (important for composing things).
2867 # I1 = matrix.column(self.base_ring(), m, m+n, V_basis[:m])
2868 # I2 = matrix.column(self.base_ring(), n, m+n, V_basis[m:])
2869 # iota_left = FiniteDimensionalEJAOperator(J1,self,I1)
2870 # iota_right = FiniteDimensionalEJAOperator(J2,self,I2)
2871 # return (iota_left, iota_right)
2872
2873 # def inner_product(self, x, y):
2874 # r"""
2875 # The standard Cartesian inner-product.
2876
2877 # We project ``x`` and ``y`` onto our factors, and add up the
2878 # inner-products from the subalgebras.
2879
2880 # SETUP::
2881
2882
2883 # sage: from mjo.eja.eja_algebra import (HadamardEJA,
2884 # ....: QuaternionHermitianEJA,
2885 # ....: DirectSumEJA)
2886
2887 # EXAMPLE::
2888
2889 # sage: J1 = HadamardEJA(3,field=QQ)
2890 # sage: J2 = QuaternionHermitianEJA(2,field=QQ,orthonormalize=False)
2891 # sage: J = DirectSumEJA(J1,J2)
2892 # sage: x1 = J1.one()
2893 # sage: x2 = x1
2894 # sage: y1 = J2.one()
2895 # sage: y2 = y1
2896 # sage: x1.inner_product(x2)
2897 # 3
2898 # sage: y1.inner_product(y2)
2899 # 2
2900 # sage: J.one().inner_product(J.one())
2901 # 5
2902
2903 # """
2904 # (pi_left, pi_right) = self.projections()
2905 # x1 = pi_left(x)
2906 # x2 = pi_right(x)
2907 # y1 = pi_left(y)
2908 # y2 = pi_right(y)
2909
2910 # return (x1.inner_product(y1) + x2.inner_product(y2))
2911
2912
2913
2914 random_eja = ConcreteEJA.random_instance