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