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