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