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