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