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