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