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