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