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