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