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