]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: fix fast charpoly in one case.
[sage.d.git] / mjo / eja / eja_algebra.py
1 """
2 Euclidean Jordan Algebras. These are formally-real Jordan Algebras;
3 specifically those where u^2 + v^2 = 0 implies that u = v = 0. They
4 are used in optimization, and have some additional nice methods beyond
5 what can be supported in a general Jordan Algebra.
6
7
8 SETUP::
9
10 sage: from mjo.eja.eja_algebra import random_eja
11
12 EXAMPLES::
13
14 sage: random_eja()
15 Euclidean Jordan algebra of dimension...
16
17 """
18
19 from itertools import repeat
20
21 from sage.algebras.quatalg.quaternion_algebra import QuaternionAlgebra
22 from sage.categories.magmatic_algebras import MagmaticAlgebras
23 from sage.combinat.free_module import CombinatorialFreeModule
24 from sage.matrix.constructor import matrix
25 from sage.matrix.matrix_space import MatrixSpace
26 from sage.misc.cachefunc import cached_method
27 from sage.misc.table import table
28 from sage.modules.free_module import FreeModule, VectorSpace
29 from sage.rings.all import (ZZ, QQ, AA, QQbar, RR, RLF, CLF,
30 PolynomialRing,
31 QuadraticField)
32 from mjo.eja.eja_element import FiniteDimensionalEJAElement
33 from mjo.eja.eja_operator import FiniteDimensionalEJAOperator
34 from mjo.eja.eja_utils import _mat2vec
35
36 class FiniteDimensionalEJA(CombinatorialFreeModule):
37 r"""
38 A finite-dimensional Euclidean Jordan algebra.
39
40 INPUT:
41
42 - basis -- a tuple of basis elements in their matrix form.
43
44 - jordan_product -- function of two elements (in matrix form)
45 that returns their jordan product in this algebra; this will
46 be applied to ``basis`` to compute a multiplication table for
47 the algebra.
48
49 - inner_product -- function of two elements (in matrix form) that
50 returns their inner product. This will be applied to ``basis`` to
51 compute an inner-product table (basically a matrix) for this algebra.
52
53 """
54 Element = FiniteDimensionalEJAElement
55
56 def __init__(self,
57 basis,
58 jordan_product,
59 inner_product,
60 field=AA,
61 orthonormalize=True,
62 associative=False,
63 check_field=True,
64 check_axioms=True,
65 prefix='e'):
66
67 if check_field:
68 if not field.is_subring(RR):
69 # Note: this does return true for the real algebraic
70 # field, the rationals, and any quadratic field where
71 # we've specified a real embedding.
72 raise ValueError("scalar field is not real")
73
74 # If the basis given to us wasn't over the field that it's
75 # supposed to be over, fix that. Or, you know, crash.
76 basis = tuple( b.change_ring(field) for b in basis )
77
78 if check_axioms:
79 # Check commutativity of the Jordan and inner-products.
80 # This has to be done before we build the multiplication
81 # and inner-product tables/matrices, because we take
82 # advantage of symmetry in the process.
83 if not all( jordan_product(bi,bj) == jordan_product(bj,bi)
84 for bi in basis
85 for bj in basis ):
86 raise ValueError("Jordan product is not commutative")
87
88 if not all( inner_product(bi,bj) == inner_product(bj,bi)
89 for bi in basis
90 for bj in basis ):
91 raise ValueError("inner-product is not commutative")
92
93
94 category = MagmaticAlgebras(field).FiniteDimensional()
95 category = category.WithBasis().Unital()
96 if associative:
97 # Element subalgebras can take advantage of this.
98 category = category.Associative()
99
100 # Call the superclass constructor so that we can use its from_vector()
101 # method to build our multiplication table.
102 n = len(basis)
103 super().__init__(field,
104 range(n),
105 prefix=prefix,
106 category=category,
107 bracket=False)
108
109 # Now comes all of the hard work. We'll be constructing an
110 # ambient vector space V that our (vectorized) basis lives in,
111 # as well as a subspace W of V spanned by those (vectorized)
112 # basis elements. The W-coordinates are the coefficients that
113 # we see in things like x = 1*e1 + 2*e2.
114 vector_basis = basis
115
116 degree = 0
117 if n > 0:
118 # Works on both column and square matrices...
119 degree = len(basis[0].list())
120
121 # Build an ambient space that fits our matrix basis when
122 # written out as "long vectors."
123 V = VectorSpace(field, degree)
124
125 # The matrix that will hole the orthonormal -> unorthonormal
126 # coordinate transformation.
127 self._deortho_matrix = None
128
129 if orthonormalize:
130 # Save a copy of the un-orthonormalized basis for later.
131 # Convert it to ambient V (vector) coordinates while we're
132 # at it, because we'd have to do it later anyway.
133 deortho_vector_basis = tuple( V(b.list()) for b in basis )
134
135 from mjo.eja.eja_utils import gram_schmidt
136 basis = tuple(gram_schmidt(basis, inner_product))
137
138 # Save the (possibly orthonormalized) matrix basis for
139 # later...
140 self._matrix_basis = basis
141
142 # Now create the vector space for the algebra, which will have
143 # its own set of non-ambient coordinates (in terms of the
144 # supplied basis).
145 vector_basis = tuple( V(b.list()) for b in basis )
146 W = V.span_of_basis( vector_basis, check=check_axioms)
147
148 if orthonormalize:
149 # Now "W" is the vector space of our algebra coordinates. The
150 # variables "X1", "X2",... refer to the entries of vectors in
151 # W. Thus to convert back and forth between the orthonormal
152 # coordinates and the given ones, we need to stick the original
153 # basis in W.
154 U = V.span_of_basis( deortho_vector_basis, check=check_axioms)
155 self._deortho_matrix = matrix( U.coordinate_vector(q)
156 for q in vector_basis )
157
158
159 # Now we actually compute the multiplication and inner-product
160 # tables/matrices using the possibly-orthonormalized basis.
161 self._inner_product_matrix = matrix.identity(field, n)
162 self._multiplication_table = [ [0 for j in range(i+1)]
163 for i in range(n) ]
164
165 # Note: the Jordan and inner-products are defined in terms
166 # of the ambient basis. It's important that their arguments
167 # are in ambient coordinates as well.
168 for i in range(n):
169 for j in range(i+1):
170 # ortho basis w.r.t. ambient coords
171 q_i = basis[i]
172 q_j = basis[j]
173
174 # The jordan product returns a matrixy answer, so we
175 # have to convert it to the algebra coordinates.
176 elt = jordan_product(q_i, q_j)
177 elt = W.coordinate_vector(V(elt.list()))
178 self._multiplication_table[i][j] = self.from_vector(elt)
179
180 if not orthonormalize:
181 # If we're orthonormalizing the basis with respect
182 # to an inner-product, then the inner-product
183 # matrix with respect to the resulting basis is
184 # just going to be the identity.
185 ip = inner_product(q_i, q_j)
186 self._inner_product_matrix[i,j] = ip
187 self._inner_product_matrix[j,i] = ip
188
189 self._inner_product_matrix._cache = {'hermitian': True}
190 self._inner_product_matrix.set_immutable()
191
192 if check_axioms:
193 if not self._is_jordanian():
194 raise ValueError("Jordan identity does not hold")
195 if not self._inner_product_is_associative():
196 raise ValueError("inner product is not associative")
197
198
199 def _coerce_map_from_base_ring(self):
200 """
201 Disable the map from the base ring into the algebra.
202
203 Performing a nonsense conversion like this automatically
204 is counterpedagogical. The fallback is to try the usual
205 element constructor, which should also fail.
206
207 SETUP::
208
209 sage: from mjo.eja.eja_algebra import random_eja
210
211 TESTS::
212
213 sage: set_random_seed()
214 sage: J = random_eja()
215 sage: J(1)
216 Traceback (most recent call last):
217 ...
218 ValueError: not an element of this algebra
219
220 """
221 return None
222
223
224 def product_on_basis(self, i, j):
225 # We only stored the lower-triangular portion of the
226 # multiplication table.
227 if j <= i:
228 return self._multiplication_table[i][j]
229 else:
230 return self._multiplication_table[j][i]
231
232 def inner_product(self, x, y):
233 """
234 The inner product associated with this Euclidean Jordan algebra.
235
236 Defaults to the trace inner product, but can be overridden by
237 subclasses if they are sure that the necessary properties are
238 satisfied.
239
240 SETUP::
241
242 sage: from mjo.eja.eja_algebra import (random_eja,
243 ....: HadamardEJA,
244 ....: BilinearFormEJA)
245
246 EXAMPLES:
247
248 Our inner product is "associative," which means the following for
249 a symmetric bilinear form::
250
251 sage: set_random_seed()
252 sage: J = random_eja()
253 sage: x,y,z = J.random_elements(3)
254 sage: (x*y).inner_product(z) == y.inner_product(x*z)
255 True
256
257 TESTS:
258
259 Ensure that this is the usual inner product for the algebras
260 over `R^n`::
261
262 sage: set_random_seed()
263 sage: J = HadamardEJA.random_instance()
264 sage: x,y = J.random_elements(2)
265 sage: actual = x.inner_product(y)
266 sage: expected = x.to_vector().inner_product(y.to_vector())
267 sage: actual == expected
268 True
269
270 Ensure that this is one-half of the trace inner-product in a
271 BilinearFormEJA that isn't just the reals (when ``n`` isn't
272 one). This is in Faraut and Koranyi, and also my "On the
273 symmetry..." paper::
274
275 sage: set_random_seed()
276 sage: J = BilinearFormEJA.random_instance()
277 sage: n = J.dimension()
278 sage: x = J.random_element()
279 sage: y = J.random_element()
280 sage: (n == 1) or (x.inner_product(y) == (x*y).trace()/2)
281 True
282 """
283 B = self._inner_product_matrix
284 return (B*x.to_vector()).inner_product(y.to_vector())
285
286
287 def _is_commutative(self):
288 r"""
289 Whether or not this algebra's multiplication table is commutative.
290
291 This method should of course always return ``True``, unless
292 this algebra was constructed with ``check_axioms=False`` and
293 passed an invalid multiplication table.
294 """
295 return all( self.product_on_basis(i,j) == self.product_on_basis(i,j)
296 for i in range(self.dimension())
297 for j in range(self.dimension()) )
298
299 def _is_jordanian(self):
300 r"""
301 Whether or not this algebra's multiplication table respects the
302 Jordan identity `(x^{2})(xy) = x(x^{2}y)`.
303
304 We only check one arrangement of `x` and `y`, so for a
305 ``True`` result to be truly true, you should also check
306 :meth:`_is_commutative`. This method should of course always
307 return ``True``, unless this algebra was constructed with
308 ``check_axioms=False`` and passed an invalid multiplication table.
309 """
310 return all( (self.monomial(i)**2)*(self.monomial(i)*self.monomial(j))
311 ==
312 (self.monomial(i))*((self.monomial(i)**2)*self.monomial(j))
313 for i in range(self.dimension())
314 for j in range(self.dimension()) )
315
316 def _inner_product_is_associative(self):
317 r"""
318 Return whether or not this algebra's inner product `B` is
319 associative; that is, whether or not `B(xy,z) = B(x,yz)`.
320
321 This method should of course always return ``True``, unless
322 this algebra was constructed with ``check_axioms=False`` and
323 passed an invalid Jordan or inner-product.
324 """
325
326 # Used to check whether or not something is zero in an inexact
327 # ring. This number is sufficient to allow the construction of
328 # QuaternionHermitianEJA(2, field=RDF) with check_axioms=True.
329 epsilon = 1e-16
330
331 for i in range(self.dimension()):
332 for j in range(self.dimension()):
333 for k in range(self.dimension()):
334 x = self.monomial(i)
335 y = self.monomial(j)
336 z = self.monomial(k)
337 diff = (x*y).inner_product(z) - x.inner_product(y*z)
338
339 if self.base_ring().is_exact():
340 if diff != 0:
341 return False
342 else:
343 if diff.abs() > epsilon:
344 return False
345
346 return True
347
348 def _element_constructor_(self, elt):
349 """
350 Construct an element of this algebra from its vector or matrix
351 representation.
352
353 This gets called only after the parent element _call_ method
354 fails to find a coercion for the argument.
355
356 SETUP::
357
358 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
359 ....: HadamardEJA,
360 ....: RealSymmetricEJA)
361
362 EXAMPLES:
363
364 The identity in `S^n` is converted to the identity in the EJA::
365
366 sage: J = RealSymmetricEJA(3)
367 sage: I = matrix.identity(QQ,3)
368 sage: J(I) == J.one()
369 True
370
371 This skew-symmetric matrix can't be represented in the EJA::
372
373 sage: J = RealSymmetricEJA(3)
374 sage: A = matrix(QQ,3, lambda i,j: i-j)
375 sage: J(A)
376 Traceback (most recent call last):
377 ...
378 ValueError: not an element of this algebra
379
380 TESTS:
381
382 Ensure that we can convert any element of the two non-matrix
383 simple algebras (whose matrix representations are columns)
384 back and forth faithfully::
385
386 sage: set_random_seed()
387 sage: J = HadamardEJA.random_instance()
388 sage: x = J.random_element()
389 sage: J(x.to_vector().column()) == x
390 True
391 sage: J = JordanSpinEJA.random_instance()
392 sage: x = J.random_element()
393 sage: J(x.to_vector().column()) == x
394 True
395
396 """
397 msg = "not an element of this algebra"
398 if elt == 0:
399 # The superclass implementation of random_element()
400 # needs to be able to coerce "0" into the algebra.
401 return self.zero()
402 elif elt in self.base_ring():
403 # Ensure that no base ring -> algebra coercion is performed
404 # by this method. There's some stupidity in sage that would
405 # otherwise propagate to this method; for example, sage thinks
406 # that the integer 3 belongs to the space of 2-by-2 matrices.
407 raise ValueError(msg)
408
409 try:
410 elt = elt.column()
411 except (AttributeError, TypeError):
412 # Try to convert a vector into a column-matrix
413 pass
414
415 if elt not in self.matrix_space():
416 raise ValueError(msg)
417
418 # Thanks for nothing! Matrix spaces aren't vector spaces in
419 # Sage, so we have to figure out its matrix-basis coordinates
420 # ourselves. We use the basis space's ring instead of the
421 # element's ring because the basis space might be an algebraic
422 # closure whereas the base ring of the 3-by-3 identity matrix
423 # could be QQ instead of QQbar.
424 #
425 # We pass check=False because the matrix basis is "guaranteed"
426 # to be linearly independent... right? Ha ha.
427 V = VectorSpace(self.base_ring(), elt.nrows()*elt.ncols())
428 W = V.span_of_basis( (_mat2vec(s) for s in self.matrix_basis()),
429 check=False)
430
431 try:
432 coords = W.coordinate_vector(_mat2vec(elt))
433 except ArithmeticError: # vector is not in free module
434 raise ValueError(msg)
435
436 return self.from_vector(coords)
437
438 def _repr_(self):
439 """
440 Return a string representation of ``self``.
441
442 SETUP::
443
444 sage: from mjo.eja.eja_algebra import JordanSpinEJA
445
446 TESTS:
447
448 Ensure that it says what we think it says::
449
450 sage: JordanSpinEJA(2, field=AA)
451 Euclidean Jordan algebra of dimension 2 over Algebraic Real Field
452 sage: JordanSpinEJA(3, field=RDF)
453 Euclidean Jordan algebra of dimension 3 over Real Double Field
454
455 """
456 fmt = "Euclidean Jordan algebra of dimension {} over {}"
457 return fmt.format(self.dimension(), self.base_ring())
458
459
460 @cached_method
461 def characteristic_polynomial_of(self):
462 """
463 Return the algebra's "characteristic polynomial of" function,
464 which is itself a multivariate polynomial that, when evaluated
465 at the coordinates of some algebra element, returns that
466 element's characteristic polynomial.
467
468 The resulting polynomial has `n+1` variables, where `n` is the
469 dimension of this algebra. The first `n` variables correspond to
470 the coordinates of an algebra element: when evaluated at the
471 coordinates of an algebra element with respect to a certain
472 basis, the result is a univariate polynomial (in the one
473 remaining variable ``t``), namely the characteristic polynomial
474 of that element.
475
476 SETUP::
477
478 sage: from mjo.eja.eja_algebra import JordanSpinEJA, TrivialEJA
479
480 EXAMPLES:
481
482 The characteristic polynomial in the spin algebra is given in
483 Alizadeh, Example 11.11::
484
485 sage: J = JordanSpinEJA(3)
486 sage: p = J.characteristic_polynomial_of(); p
487 X1^2 - X2^2 - X3^2 + (-2*t)*X1 + t^2
488 sage: xvec = J.one().to_vector()
489 sage: p(*xvec)
490 t^2 - 2*t + 1
491
492 By definition, the characteristic polynomial is a monic
493 degree-zero polynomial in a rank-zero algebra. Note that
494 Cayley-Hamilton is indeed satisfied since the polynomial
495 ``1`` evaluates to the identity element of the algebra on
496 any argument::
497
498 sage: J = TrivialEJA()
499 sage: J.characteristic_polynomial_of()
500 1
501
502 """
503 r = self.rank()
504 n = self.dimension()
505
506 # The list of coefficient polynomials a_0, a_1, a_2, ..., a_(r-1).
507 a = self._charpoly_coefficients()
508
509 # We go to a bit of trouble here to reorder the
510 # indeterminates, so that it's easier to evaluate the
511 # characteristic polynomial at x's coordinates and get back
512 # something in terms of t, which is what we want.
513 S = PolynomialRing(self.base_ring(),'t')
514 t = S.gen(0)
515 if r > 0:
516 R = a[0].parent()
517 S = PolynomialRing(S, R.variable_names())
518 t = S(t)
519
520 return (t**r + sum( a[k]*(t**k) for k in range(r) ))
521
522 def coordinate_polynomial_ring(self):
523 r"""
524 The multivariate polynomial ring in which this algebra's
525 :meth:`characteristic_polynomial_of` lives.
526
527 SETUP::
528
529 sage: from mjo.eja.eja_algebra import (HadamardEJA,
530 ....: RealSymmetricEJA)
531
532 EXAMPLES::
533
534 sage: J = HadamardEJA(2)
535 sage: J.coordinate_polynomial_ring()
536 Multivariate Polynomial Ring in X1, X2...
537 sage: J = RealSymmetricEJA(3,field=QQ,orthonormalize=False)
538 sage: J.coordinate_polynomial_ring()
539 Multivariate Polynomial Ring in X1, X2, X3, X4, X5, X6...
540
541 """
542 var_names = tuple( "X%d" % z for z in range(1, self.dimension()+1) )
543 return PolynomialRing(self.base_ring(), var_names)
544
545 def inner_product(self, x, y):
546 """
547 The inner product associated with this Euclidean Jordan algebra.
548
549 Defaults to the trace inner product, but can be overridden by
550 subclasses if they are sure that the necessary properties are
551 satisfied.
552
553 SETUP::
554
555 sage: from mjo.eja.eja_algebra import (random_eja,
556 ....: HadamardEJA,
557 ....: BilinearFormEJA)
558
559 EXAMPLES:
560
561 Our inner product is "associative," which means the following for
562 a symmetric bilinear form::
563
564 sage: set_random_seed()
565 sage: J = random_eja()
566 sage: x,y,z = J.random_elements(3)
567 sage: (x*y).inner_product(z) == y.inner_product(x*z)
568 True
569
570 TESTS:
571
572 Ensure that this is the usual inner product for the algebras
573 over `R^n`::
574
575 sage: set_random_seed()
576 sage: J = HadamardEJA.random_instance()
577 sage: x,y = J.random_elements(2)
578 sage: actual = x.inner_product(y)
579 sage: expected = x.to_vector().inner_product(y.to_vector())
580 sage: actual == expected
581 True
582
583 Ensure that this is one-half of the trace inner-product in a
584 BilinearFormEJA that isn't just the reals (when ``n`` isn't
585 one). This is in Faraut and Koranyi, and also my "On the
586 symmetry..." paper::
587
588 sage: set_random_seed()
589 sage: J = BilinearFormEJA.random_instance()
590 sage: n = J.dimension()
591 sage: x = J.random_element()
592 sage: y = J.random_element()
593 sage: (n == 1) or (x.inner_product(y) == (x*y).trace()/2)
594 True
595 """
596 B = self._inner_product_matrix
597 return (B*x.to_vector()).inner_product(y.to_vector())
598
599
600 def is_trivial(self):
601 """
602 Return whether or not this algebra is trivial.
603
604 A trivial algebra contains only the zero element.
605
606 SETUP::
607
608 sage: from mjo.eja.eja_algebra import (ComplexHermitianEJA,
609 ....: TrivialEJA)
610
611 EXAMPLES::
612
613 sage: J = ComplexHermitianEJA(3)
614 sage: J.is_trivial()
615 False
616
617 ::
618
619 sage: J = TrivialEJA()
620 sage: J.is_trivial()
621 True
622
623 """
624 return self.dimension() == 0
625
626
627 def multiplication_table(self):
628 """
629 Return a visual representation of this algebra's multiplication
630 table (on basis elements).
631
632 SETUP::
633
634 sage: from mjo.eja.eja_algebra import JordanSpinEJA
635
636 EXAMPLES::
637
638 sage: J = JordanSpinEJA(4)
639 sage: J.multiplication_table()
640 +----++----+----+----+----+
641 | * || e0 | e1 | e2 | e3 |
642 +====++====+====+====+====+
643 | e0 || e0 | e1 | e2 | e3 |
644 +----++----+----+----+----+
645 | e1 || e1 | e0 | 0 | 0 |
646 +----++----+----+----+----+
647 | e2 || e2 | 0 | e0 | 0 |
648 +----++----+----+----+----+
649 | e3 || e3 | 0 | 0 | e0 |
650 +----++----+----+----+----+
651
652 """
653 n = self.dimension()
654 # Prepend the header row.
655 M = [["*"] + list(self.gens())]
656
657 # And to each subsequent row, prepend an entry that belongs to
658 # the left-side "header column."
659 M += [ [self.monomial(i)] + [ self.product_on_basis(i,j)
660 for j in range(n) ]
661 for i in range(n) ]
662
663 return table(M, header_row=True, header_column=True, frame=True)
664
665
666 def matrix_basis(self):
667 """
668 Return an (often more natural) representation of this algebras
669 basis as an ordered tuple of matrices.
670
671 Every finite-dimensional Euclidean Jordan Algebra is a, up to
672 Jordan isomorphism, a direct sum of five simple
673 algebras---four of which comprise Hermitian matrices. And the
674 last type of algebra can of course be thought of as `n`-by-`1`
675 column matrices (ambiguusly called column vectors) to avoid
676 special cases. As a result, matrices (and column vectors) are
677 a natural representation format for Euclidean Jordan algebra
678 elements.
679
680 But, when we construct an algebra from a basis of matrices,
681 those matrix representations are lost in favor of coordinate
682 vectors *with respect to* that basis. We could eventually
683 convert back if we tried hard enough, but having the original
684 representations handy is valuable enough that we simply store
685 them and return them from this method.
686
687 Why implement this for non-matrix algebras? Avoiding special
688 cases for the :class:`BilinearFormEJA` pays with simplicity in
689 its own right. But mainly, we would like to be able to assume
690 that elements of a :class:`DirectSumEJA` can be displayed
691 nicely, without having to have special classes for direct sums
692 one of whose components was a matrix algebra.
693
694 SETUP::
695
696 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
697 ....: RealSymmetricEJA)
698
699 EXAMPLES::
700
701 sage: J = RealSymmetricEJA(2)
702 sage: J.basis()
703 Finite family {0: e0, 1: e1, 2: e2}
704 sage: J.matrix_basis()
705 (
706 [1 0] [ 0 0.7071067811865475?] [0 0]
707 [0 0], [0.7071067811865475? 0], [0 1]
708 )
709
710 ::
711
712 sage: J = JordanSpinEJA(2)
713 sage: J.basis()
714 Finite family {0: e0, 1: e1}
715 sage: J.matrix_basis()
716 (
717 [1] [0]
718 [0], [1]
719 )
720 """
721 return self._matrix_basis
722
723
724 def matrix_space(self):
725 """
726 Return the matrix space in which this algebra's elements live, if
727 we think of them as matrices (including column vectors of the
728 appropriate size).
729
730 Generally this will be an `n`-by-`1` column-vector space,
731 except when the algebra is trivial. There it's `n`-by-`n`
732 (where `n` is zero), to ensure that two elements of the matrix
733 space (empty matrices) can be multiplied.
734
735 Matrix algebras override this with something more useful.
736 """
737 if self.is_trivial():
738 return MatrixSpace(self.base_ring(), 0)
739 else:
740 return self._matrix_basis[0].matrix_space()
741
742
743 @cached_method
744 def one(self):
745 """
746 Return the unit element of this algebra.
747
748 SETUP::
749
750 sage: from mjo.eja.eja_algebra import (HadamardEJA,
751 ....: random_eja)
752
753 EXAMPLES::
754
755 sage: J = HadamardEJA(5)
756 sage: J.one()
757 e0 + e1 + e2 + e3 + e4
758
759 TESTS:
760
761 The identity element acts like the identity::
762
763 sage: set_random_seed()
764 sage: J = random_eja()
765 sage: x = J.random_element()
766 sage: J.one()*x == x and x*J.one() == x
767 True
768
769 The matrix of the unit element's operator is the identity::
770
771 sage: set_random_seed()
772 sage: J = random_eja()
773 sage: actual = J.one().operator().matrix()
774 sage: expected = matrix.identity(J.base_ring(), J.dimension())
775 sage: actual == expected
776 True
777
778 Ensure that the cached unit element (often precomputed by
779 hand) agrees with the computed one::
780
781 sage: set_random_seed()
782 sage: J = random_eja()
783 sage: cached = J.one()
784 sage: J.one.clear_cache()
785 sage: J.one() == cached
786 True
787
788 """
789 # We can brute-force compute the matrices of the operators
790 # that correspond to the basis elements of this algebra.
791 # If some linear combination of those basis elements is the
792 # algebra identity, then the same linear combination of
793 # their matrices has to be the identity matrix.
794 #
795 # Of course, matrices aren't vectors in sage, so we have to
796 # appeal to the "long vectors" isometry.
797 oper_vecs = [ _mat2vec(g.operator().matrix()) for g in self.gens() ]
798
799 # Now we use basic linear algebra to find the coefficients,
800 # of the matrices-as-vectors-linear-combination, which should
801 # work for the original algebra basis too.
802 A = matrix(self.base_ring(), oper_vecs)
803
804 # We used the isometry on the left-hand side already, but we
805 # still need to do it for the right-hand side. Recall that we
806 # wanted something that summed to the identity matrix.
807 b = _mat2vec( matrix.identity(self.base_ring(), self.dimension()) )
808
809 # Now if there's an identity element in the algebra, this
810 # should work. We solve on the left to avoid having to
811 # transpose the matrix "A".
812 return self.from_vector(A.solve_left(b))
813
814
815 def peirce_decomposition(self, c):
816 """
817 The Peirce decomposition of this algebra relative to the
818 idempotent ``c``.
819
820 In the future, this can be extended to a complete system of
821 orthogonal idempotents.
822
823 INPUT:
824
825 - ``c`` -- an idempotent of this algebra.
826
827 OUTPUT:
828
829 A triple (J0, J5, J1) containing two subalgebras and one subspace
830 of this algebra,
831
832 - ``J0`` -- the algebra on the eigenspace of ``c.operator()``
833 corresponding to the eigenvalue zero.
834
835 - ``J5`` -- the eigenspace (NOT a subalgebra) of ``c.operator()``
836 corresponding to the eigenvalue one-half.
837
838 - ``J1`` -- the algebra on the eigenspace of ``c.operator()``
839 corresponding to the eigenvalue one.
840
841 These are the only possible eigenspaces for that operator, and this
842 algebra is a direct sum of them. The spaces ``J0`` and ``J1`` are
843 orthogonal, and are subalgebras of this algebra with the appropriate
844 restrictions.
845
846 SETUP::
847
848 sage: from mjo.eja.eja_algebra import random_eja, RealSymmetricEJA
849
850 EXAMPLES:
851
852 The canonical example comes from the symmetric matrices, which
853 decompose into diagonal and off-diagonal parts::
854
855 sage: J = RealSymmetricEJA(3)
856 sage: C = matrix(QQ, [ [1,0,0],
857 ....: [0,1,0],
858 ....: [0,0,0] ])
859 sage: c = J(C)
860 sage: J0,J5,J1 = J.peirce_decomposition(c)
861 sage: J0
862 Euclidean Jordan algebra of dimension 1...
863 sage: J5
864 Vector space of degree 6 and dimension 2...
865 sage: J1
866 Euclidean Jordan algebra of dimension 3...
867 sage: J0.one().to_matrix()
868 [0 0 0]
869 [0 0 0]
870 [0 0 1]
871 sage: orig_df = AA.options.display_format
872 sage: AA.options.display_format = 'radical'
873 sage: J.from_vector(J5.basis()[0]).to_matrix()
874 [ 0 0 1/2*sqrt(2)]
875 [ 0 0 0]
876 [1/2*sqrt(2) 0 0]
877 sage: J.from_vector(J5.basis()[1]).to_matrix()
878 [ 0 0 0]
879 [ 0 0 1/2*sqrt(2)]
880 [ 0 1/2*sqrt(2) 0]
881 sage: AA.options.display_format = orig_df
882 sage: J1.one().to_matrix()
883 [1 0 0]
884 [0 1 0]
885 [0 0 0]
886
887 TESTS:
888
889 Every algebra decomposes trivially with respect to its identity
890 element::
891
892 sage: set_random_seed()
893 sage: J = random_eja()
894 sage: J0,J5,J1 = J.peirce_decomposition(J.one())
895 sage: J0.dimension() == 0 and J5.dimension() == 0
896 True
897 sage: J1.superalgebra() == J and J1.dimension() == J.dimension()
898 True
899
900 The decomposition is into eigenspaces, and its components are
901 therefore necessarily orthogonal. Moreover, the identity
902 elements in the two subalgebras are the projections onto their
903 respective subspaces of the superalgebra's identity element::
904
905 sage: set_random_seed()
906 sage: J = random_eja()
907 sage: x = J.random_element()
908 sage: if not J.is_trivial():
909 ....: while x.is_nilpotent():
910 ....: x = J.random_element()
911 sage: c = x.subalgebra_idempotent()
912 sage: J0,J5,J1 = J.peirce_decomposition(c)
913 sage: ipsum = 0
914 sage: for (w,y,z) in zip(J0.basis(), J5.basis(), J1.basis()):
915 ....: w = w.superalgebra_element()
916 ....: y = J.from_vector(y)
917 ....: z = z.superalgebra_element()
918 ....: ipsum += w.inner_product(y).abs()
919 ....: ipsum += w.inner_product(z).abs()
920 ....: ipsum += y.inner_product(z).abs()
921 sage: ipsum
922 0
923 sage: J1(c) == J1.one()
924 True
925 sage: J0(J.one() - c) == J0.one()
926 True
927
928 """
929 if not c.is_idempotent():
930 raise ValueError("element is not idempotent: %s" % c)
931
932 from mjo.eja.eja_subalgebra import FiniteDimensionalEJASubalgebra
933
934 # Default these to what they should be if they turn out to be
935 # trivial, because eigenspaces_left() won't return eigenvalues
936 # corresponding to trivial spaces (e.g. it returns only the
937 # eigenspace corresponding to lambda=1 if you take the
938 # decomposition relative to the identity element).
939 trivial = FiniteDimensionalEJASubalgebra(self, ())
940 J0 = trivial # eigenvalue zero
941 J5 = VectorSpace(self.base_ring(), 0) # eigenvalue one-half
942 J1 = trivial # eigenvalue one
943
944 for (eigval, eigspace) in c.operator().matrix().right_eigenspaces():
945 if eigval == ~(self.base_ring()(2)):
946 J5 = eigspace
947 else:
948 gens = tuple( self.from_vector(b) for b in eigspace.basis() )
949 subalg = FiniteDimensionalEJASubalgebra(self,
950 gens,
951 check_axioms=False)
952 if eigval == 0:
953 J0 = subalg
954 elif eigval == 1:
955 J1 = subalg
956 else:
957 raise ValueError("unexpected eigenvalue: %s" % eigval)
958
959 return (J0, J5, J1)
960
961
962 def random_element(self, thorough=False):
963 r"""
964 Return a random element of this algebra.
965
966 Our algebra superclass method only returns a linear
967 combination of at most two basis elements. We instead
968 want the vector space "random element" method that
969 returns a more diverse selection.
970
971 INPUT:
972
973 - ``thorough`` -- (boolean; default False) whether or not we
974 should generate irrational coefficients for the random
975 element when our base ring is irrational; this slows the
976 algebra operations to a crawl, but any truly random method
977 should include them
978
979 """
980 # For a general base ring... maybe we can trust this to do the
981 # right thing? Unlikely, but.
982 V = self.vector_space()
983 v = V.random_element()
984
985 if self.base_ring() is AA:
986 # The "random element" method of the algebraic reals is
987 # stupid at the moment, and only returns integers between
988 # -2 and 2, inclusive:
989 #
990 # https://trac.sagemath.org/ticket/30875
991 #
992 # Instead, we implement our own "random vector" method,
993 # and then coerce that into the algebra. We use the vector
994 # space degree here instead of the dimension because a
995 # subalgebra could (for example) be spanned by only two
996 # vectors, each with five coordinates. We need to
997 # generate all five coordinates.
998 if thorough:
999 v *= QQbar.random_element().real()
1000 else:
1001 v *= QQ.random_element()
1002
1003 return self.from_vector(V.coordinate_vector(v))
1004
1005 def random_elements(self, count, thorough=False):
1006 """
1007 Return ``count`` random elements as a tuple.
1008
1009 INPUT:
1010
1011 - ``thorough`` -- (boolean; default False) whether or not we
1012 should generate irrational coefficients for the random
1013 elements when our base ring is irrational; this slows the
1014 algebra operations to a crawl, but any truly random method
1015 should include them
1016
1017 SETUP::
1018
1019 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1020
1021 EXAMPLES::
1022
1023 sage: J = JordanSpinEJA(3)
1024 sage: x,y,z = J.random_elements(3)
1025 sage: all( [ x in J, y in J, z in J ])
1026 True
1027 sage: len( J.random_elements(10) ) == 10
1028 True
1029
1030 """
1031 return tuple( self.random_element(thorough)
1032 for idx in range(count) )
1033
1034
1035 @cached_method
1036 def _charpoly_coefficients(self):
1037 r"""
1038 The `r` polynomial coefficients of the "characteristic polynomial
1039 of" function.
1040 """
1041 n = self.dimension()
1042 R = self.coordinate_polynomial_ring()
1043 vars = R.gens()
1044 F = R.fraction_field()
1045
1046 def L_x_i_j(i,j):
1047 # From a result in my book, these are the entries of the
1048 # basis representation of L_x.
1049 return sum( vars[k]*self.monomial(k).operator().matrix()[i,j]
1050 for k in range(n) )
1051
1052 L_x = matrix(F, n, n, L_x_i_j)
1053
1054 r = None
1055 if self.rank.is_in_cache():
1056 r = self.rank()
1057 # There's no need to pad the system with redundant
1058 # columns if we *know* they'll be redundant.
1059 n = r
1060
1061 # Compute an extra power in case the rank is equal to
1062 # the dimension (otherwise, we would stop at x^(r-1)).
1063 x_powers = [ (L_x**k)*self.one().to_vector()
1064 for k in range(n+1) ]
1065 A = matrix.column(F, x_powers[:n])
1066 AE = A.extended_echelon_form()
1067 E = AE[:,n:]
1068 A_rref = AE[:,:n]
1069 if r is None:
1070 r = A_rref.rank()
1071 b = x_powers[r]
1072
1073 # The theory says that only the first "r" coefficients are
1074 # nonzero, and they actually live in the original polynomial
1075 # ring and not the fraction field. We negate them because
1076 # in the actual characteristic polynomial, they get moved
1077 # to the other side where x^r lives.
1078 return -A_rref.solve_right(E*b).change_ring(R)[:r]
1079
1080 @cached_method
1081 def rank(self):
1082 r"""
1083 Return the rank of this EJA.
1084
1085 This is a cached method because we know the rank a priori for
1086 all of the algebras we can construct. Thus we can avoid the
1087 expensive ``_charpoly_coefficients()`` call unless we truly
1088 need to compute the whole characteristic polynomial.
1089
1090 SETUP::
1091
1092 sage: from mjo.eja.eja_algebra import (HadamardEJA,
1093 ....: JordanSpinEJA,
1094 ....: RealSymmetricEJA,
1095 ....: ComplexHermitianEJA,
1096 ....: QuaternionHermitianEJA,
1097 ....: random_eja)
1098
1099 EXAMPLES:
1100
1101 The rank of the Jordan spin algebra is always two::
1102
1103 sage: JordanSpinEJA(2).rank()
1104 2
1105 sage: JordanSpinEJA(3).rank()
1106 2
1107 sage: JordanSpinEJA(4).rank()
1108 2
1109
1110 The rank of the `n`-by-`n` Hermitian real, complex, or
1111 quaternion matrices is `n`::
1112
1113 sage: RealSymmetricEJA(4).rank()
1114 4
1115 sage: ComplexHermitianEJA(3).rank()
1116 3
1117 sage: QuaternionHermitianEJA(2).rank()
1118 2
1119
1120 TESTS:
1121
1122 Ensure that every EJA that we know how to construct has a
1123 positive integer rank, unless the algebra is trivial in
1124 which case its rank will be zero::
1125
1126 sage: set_random_seed()
1127 sage: J = random_eja()
1128 sage: r = J.rank()
1129 sage: r in ZZ
1130 True
1131 sage: r > 0 or (r == 0 and J.is_trivial())
1132 True
1133
1134 Ensure that computing the rank actually works, since the ranks
1135 of all simple algebras are known and will be cached by default::
1136
1137 sage: set_random_seed() # long time
1138 sage: J = random_eja() # long time
1139 sage: caches = J.rank() # long time
1140 sage: J.rank.clear_cache() # long time
1141 sage: J.rank() == cached # long time
1142 True
1143
1144 """
1145 return len(self._charpoly_coefficients())
1146
1147
1148 def vector_space(self):
1149 """
1150 Return the vector space that underlies this algebra.
1151
1152 SETUP::
1153
1154 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1155
1156 EXAMPLES::
1157
1158 sage: J = RealSymmetricEJA(2)
1159 sage: J.vector_space()
1160 Vector space of dimension 3 over...
1161
1162 """
1163 return self.zero().to_vector().parent().ambient_vector_space()
1164
1165
1166 Element = FiniteDimensionalEJAElement
1167
1168 class RationalBasisEJA(FiniteDimensionalEJA):
1169 r"""
1170 New class for algebras whose supplied basis elements have all rational entries.
1171
1172 SETUP::
1173
1174 sage: from mjo.eja.eja_algebra import BilinearFormEJA
1175
1176 EXAMPLES:
1177
1178 The supplied basis is orthonormalized by default::
1179
1180 sage: B = matrix(QQ, [[1, 0, 0], [0, 25, -32], [0, -32, 41]])
1181 sage: J = BilinearFormEJA(B)
1182 sage: J.matrix_basis()
1183 (
1184 [1] [ 0] [ 0]
1185 [0] [1/5] [32/5]
1186 [0], [ 0], [ 5]
1187 )
1188
1189 """
1190 def __init__(self,
1191 basis,
1192 jordan_product,
1193 inner_product,
1194 field=AA,
1195 check_field=True,
1196 **kwargs):
1197
1198 if check_field:
1199 # Abuse the check_field parameter to check that the entries of
1200 # out basis (in ambient coordinates) are in the field QQ.
1201 if not all( all(b_i in QQ for b_i in b.list()) for b in basis ):
1202 raise TypeError("basis not rational")
1203
1204 if field is not QQ:
1205 # There's no point in constructing the extra algebra if this
1206 # one is already rational.
1207 #
1208 # Note: the same Jordan and inner-products work here,
1209 # because they are necessarily defined with respect to
1210 # ambient coordinates and not any particular basis.
1211 self._rational_algebra = FiniteDimensionalEJA(
1212 basis,
1213 jordan_product,
1214 inner_product,
1215 field=QQ,
1216 orthonormalize=False,
1217 check_field=False,
1218 check_axioms=False)
1219
1220 super().__init__(basis,
1221 jordan_product,
1222 inner_product,
1223 field=field,
1224 check_field=check_field,
1225 **kwargs)
1226
1227 @cached_method
1228 def _charpoly_coefficients(self):
1229 r"""
1230 SETUP::
1231
1232 sage: from mjo.eja.eja_algebra import (BilinearFormEJA,
1233 ....: JordanSpinEJA)
1234
1235 EXAMPLES:
1236
1237 The base ring of the resulting polynomial coefficients is what
1238 it should be, and not the rationals (unless the algebra was
1239 already over the rationals)::
1240
1241 sage: J = JordanSpinEJA(3)
1242 sage: J._charpoly_coefficients()
1243 (X1^2 - X2^2 - X3^2, -2*X1)
1244 sage: a0 = J._charpoly_coefficients()[0]
1245 sage: J.base_ring()
1246 Algebraic Real Field
1247 sage: a0.base_ring()
1248 Algebraic Real Field
1249
1250 """
1251 if self._rational_algebra is None:
1252 # There's no need to construct *another* algebra over the
1253 # rationals if this one is already over the
1254 # rationals. Likewise, if we never orthonormalized our
1255 # basis, we might as well just use the given one.
1256 return super()._charpoly_coefficients()
1257
1258 # Do the computation over the rationals. The answer will be
1259 # the same, because all we've done is a change of basis.
1260 # Then, change back from QQ to our real base ring
1261 a = ( a_i.change_ring(self.base_ring())
1262 for a_i in self._rational_algebra._charpoly_coefficients() )
1263
1264 if self._deortho_matrix is None:
1265 # This can happen if our base ring was, say, AA and we
1266 # chose not to (or didn't need to) orthonormalize. It's
1267 # still faster to do the computations over QQ even if
1268 # the numbers in the boxes stay the same.
1269 return tuple(a)
1270
1271 # Otherwise, convert the coordinate variables back to the
1272 # deorthonormalized ones.
1273 R = self.coordinate_polynomial_ring()
1274 from sage.modules.free_module_element import vector
1275 X = vector(R, R.gens())
1276 BX = self._deortho_matrix*X
1277
1278 subs_dict = { X[i]: BX[i] for i in range(len(X)) }
1279 return tuple( a_i.subs(subs_dict) for a_i in a )
1280
1281 class ConcreteEJA(RationalBasisEJA):
1282 r"""
1283 A class for the Euclidean Jordan algebras that we know by name.
1284
1285 These are the Jordan algebras whose basis, multiplication table,
1286 rank, and so on are known a priori. More to the point, they are
1287 the Euclidean Jordan algebras for which we are able to conjure up
1288 a "random instance."
1289
1290 SETUP::
1291
1292 sage: from mjo.eja.eja_algebra import ConcreteEJA
1293
1294 TESTS:
1295
1296 Our basis is normalized with respect to the algebra's inner
1297 product, unless we specify otherwise::
1298
1299 sage: set_random_seed()
1300 sage: J = ConcreteEJA.random_instance()
1301 sage: all( b.norm() == 1 for b in J.gens() )
1302 True
1303
1304 Since our basis is orthonormal with respect to the algebra's inner
1305 product, and since we know that this algebra is an EJA, any
1306 left-multiplication operator's matrix will be symmetric because
1307 natural->EJA basis representation is an isometry and within the
1308 EJA the operator is self-adjoint by the Jordan axiom::
1309
1310 sage: set_random_seed()
1311 sage: J = ConcreteEJA.random_instance()
1312 sage: x = J.random_element()
1313 sage: x.operator().is_self_adjoint()
1314 True
1315 """
1316
1317 @staticmethod
1318 def _max_random_instance_size():
1319 """
1320 Return an integer "size" that is an upper bound on the size of
1321 this algebra when it is used in a random test
1322 case. Unfortunately, the term "size" is ambiguous -- when
1323 dealing with `R^n` under either the Hadamard or Jordan spin
1324 product, the "size" refers to the dimension `n`. When dealing
1325 with a matrix algebra (real symmetric or complex/quaternion
1326 Hermitian), it refers to the size of the matrix, which is far
1327 less than the dimension of the underlying vector space.
1328
1329 This method must be implemented in each subclass.
1330 """
1331 raise NotImplementedError
1332
1333 @classmethod
1334 def random_instance(cls, *args, **kwargs):
1335 """
1336 Return a random instance of this type of algebra.
1337
1338 This method should be implemented in each subclass.
1339 """
1340 from sage.misc.prandom import choice
1341 eja_class = choice(cls.__subclasses__())
1342
1343 # These all bubble up to the RationalBasisEJA superclass
1344 # constructor, so any (kw)args valid there are also valid
1345 # here.
1346 return eja_class.random_instance(*args, **kwargs)
1347
1348
1349 class MatrixEJA:
1350 @staticmethod
1351 def dimension_over_reals():
1352 r"""
1353 The dimension of this matrix's base ring over the reals.
1354
1355 The reals are dimension one over themselves, obviously; that's
1356 just `\mathbb{R}^{1}`. Likewise, the complex numbers `a + bi`
1357 have dimension two. Finally, the quaternions have dimension
1358 four over the reals.
1359
1360 This is used to determine the size of the matrix returned from
1361 :meth:`real_embed`, among other things.
1362 """
1363 raise NotImplementedError
1364
1365 @classmethod
1366 def real_embed(cls,M):
1367 """
1368 Embed the matrix ``M`` into a space of real matrices.
1369
1370 The matrix ``M`` can have entries in any field at the moment:
1371 the real numbers, complex numbers, or quaternions. And although
1372 they are not a field, we can probably support octonions at some
1373 point, too. This function returns a real matrix that "acts like"
1374 the original with respect to matrix multiplication; i.e.
1375
1376 real_embed(M*N) = real_embed(M)*real_embed(N)
1377
1378 """
1379 if M.ncols() != M.nrows():
1380 raise ValueError("the matrix 'M' must be square")
1381 return M
1382
1383
1384 @classmethod
1385 def real_unembed(cls,M):
1386 """
1387 The inverse of :meth:`real_embed`.
1388 """
1389 if M.ncols() != M.nrows():
1390 raise ValueError("the matrix 'M' must be square")
1391 if not ZZ(M.nrows()).mod(cls.dimension_over_reals()).is_zero():
1392 raise ValueError("the matrix 'M' must be a real embedding")
1393 return M
1394
1395 @staticmethod
1396 def jordan_product(X,Y):
1397 return (X*Y + Y*X)/2
1398
1399 @classmethod
1400 def trace_inner_product(cls,X,Y):
1401 r"""
1402 Compute the trace inner-product of two real-embeddings.
1403
1404 SETUP::
1405
1406 sage: from mjo.eja.eja_algebra import (RealSymmetricEJA,
1407 ....: ComplexHermitianEJA,
1408 ....: QuaternionHermitianEJA)
1409
1410 EXAMPLES::
1411
1412 This gives the same answer as it would if we computed the trace
1413 from the unembedded (original) matrices::
1414
1415 sage: set_random_seed()
1416 sage: J = RealSymmetricEJA.random_instance()
1417 sage: x,y = J.random_elements(2)
1418 sage: Xe = x.to_matrix()
1419 sage: Ye = y.to_matrix()
1420 sage: X = J.real_unembed(Xe)
1421 sage: Y = J.real_unembed(Ye)
1422 sage: expected = (X*Y).trace()
1423 sage: actual = J.trace_inner_product(Xe,Ye)
1424 sage: actual == expected
1425 True
1426
1427 ::
1428
1429 sage: set_random_seed()
1430 sage: J = ComplexHermitianEJA.random_instance()
1431 sage: x,y = J.random_elements(2)
1432 sage: Xe = x.to_matrix()
1433 sage: Ye = y.to_matrix()
1434 sage: X = J.real_unembed(Xe)
1435 sage: Y = J.real_unembed(Ye)
1436 sage: expected = (X*Y).trace().real()
1437 sage: actual = J.trace_inner_product(Xe,Ye)
1438 sage: actual == expected
1439 True
1440
1441 ::
1442
1443 sage: set_random_seed()
1444 sage: J = QuaternionHermitianEJA.random_instance()
1445 sage: x,y = J.random_elements(2)
1446 sage: Xe = x.to_matrix()
1447 sage: Ye = y.to_matrix()
1448 sage: X = J.real_unembed(Xe)
1449 sage: Y = J.real_unembed(Ye)
1450 sage: expected = (X*Y).trace().coefficient_tuple()[0]
1451 sage: actual = J.trace_inner_product(Xe,Ye)
1452 sage: actual == expected
1453 True
1454
1455 """
1456 Xu = cls.real_unembed(X)
1457 Yu = cls.real_unembed(Y)
1458 tr = (Xu*Yu).trace()
1459
1460 try:
1461 # Works in QQ, AA, RDF, et cetera.
1462 return tr.real()
1463 except AttributeError:
1464 # A quaternion doesn't have a real() method, but does
1465 # have coefficient_tuple() method that returns the
1466 # coefficients of 1, i, j, and k -- in that order.
1467 return tr.coefficient_tuple()[0]
1468
1469
1470 class RealMatrixEJA(MatrixEJA):
1471 @staticmethod
1472 def dimension_over_reals():
1473 return 1
1474
1475
1476 class RealSymmetricEJA(ConcreteEJA, RealMatrixEJA):
1477 """
1478 The rank-n simple EJA consisting of real symmetric n-by-n
1479 matrices, the usual symmetric Jordan product, and the trace inner
1480 product. It has dimension `(n^2 + n)/2` over the reals.
1481
1482 SETUP::
1483
1484 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1485
1486 EXAMPLES::
1487
1488 sage: J = RealSymmetricEJA(2)
1489 sage: e0, e1, e2 = J.gens()
1490 sage: e0*e0
1491 e0
1492 sage: e1*e1
1493 1/2*e0 + 1/2*e2
1494 sage: e2*e2
1495 e2
1496
1497 In theory, our "field" can be any subfield of the reals::
1498
1499 sage: RealSymmetricEJA(2, field=RDF)
1500 Euclidean Jordan algebra of dimension 3 over Real Double Field
1501 sage: RealSymmetricEJA(2, field=RR)
1502 Euclidean Jordan algebra of dimension 3 over Real Field with
1503 53 bits of precision
1504
1505 TESTS:
1506
1507 The dimension of this algebra is `(n^2 + n) / 2`::
1508
1509 sage: set_random_seed()
1510 sage: n_max = RealSymmetricEJA._max_random_instance_size()
1511 sage: n = ZZ.random_element(1, n_max)
1512 sage: J = RealSymmetricEJA(n)
1513 sage: J.dimension() == (n^2 + n)/2
1514 True
1515
1516 The Jordan multiplication is what we think it is::
1517
1518 sage: set_random_seed()
1519 sage: J = RealSymmetricEJA.random_instance()
1520 sage: x,y = J.random_elements(2)
1521 sage: actual = (x*y).to_matrix()
1522 sage: X = x.to_matrix()
1523 sage: Y = y.to_matrix()
1524 sage: expected = (X*Y + Y*X)/2
1525 sage: actual == expected
1526 True
1527 sage: J(expected) == x*y
1528 True
1529
1530 We can change the generator prefix::
1531
1532 sage: RealSymmetricEJA(3, prefix='q').gens()
1533 (q0, q1, q2, q3, q4, q5)
1534
1535 We can construct the (trivial) algebra of rank zero::
1536
1537 sage: RealSymmetricEJA(0)
1538 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
1539
1540 """
1541 @classmethod
1542 def _denormalized_basis(cls, n):
1543 """
1544 Return a basis for the space of real symmetric n-by-n matrices.
1545
1546 SETUP::
1547
1548 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1549
1550 TESTS::
1551
1552 sage: set_random_seed()
1553 sage: n = ZZ.random_element(1,5)
1554 sage: B = RealSymmetricEJA._denormalized_basis(n)
1555 sage: all( M.is_symmetric() for M in B)
1556 True
1557
1558 """
1559 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1560 # coordinates.
1561 S = []
1562 for i in range(n):
1563 for j in range(i+1):
1564 Eij = matrix(ZZ, n, lambda k,l: k==i and l==j)
1565 if i == j:
1566 Sij = Eij
1567 else:
1568 Sij = Eij + Eij.transpose()
1569 S.append(Sij)
1570 return tuple(S)
1571
1572
1573 @staticmethod
1574 def _max_random_instance_size():
1575 return 4 # Dimension 10
1576
1577 @classmethod
1578 def random_instance(cls, **kwargs):
1579 """
1580 Return a random instance of this type of algebra.
1581 """
1582 n = ZZ.random_element(cls._max_random_instance_size() + 1)
1583 return cls(n, **kwargs)
1584
1585 def __init__(self, n, **kwargs):
1586 # We know this is a valid EJA, but will double-check
1587 # if the user passes check_axioms=True.
1588 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
1589
1590 super(RealSymmetricEJA, self).__init__(self._denormalized_basis(n),
1591 self.jordan_product,
1592 self.trace_inner_product,
1593 **kwargs)
1594
1595 # TODO: this could be factored out somehow, but is left here
1596 # because the MatrixEJA is not presently a subclass of the
1597 # FDEJA class that defines rank() and one().
1598 self.rank.set_cache(n)
1599 idV = matrix.identity(ZZ, self.dimension_over_reals()*n)
1600 self.one.set_cache(self(idV))
1601
1602
1603
1604 class ComplexMatrixEJA(MatrixEJA):
1605 @staticmethod
1606 def dimension_over_reals():
1607 return 2
1608
1609 @classmethod
1610 def real_embed(cls,M):
1611 """
1612 Embed the n-by-n complex matrix ``M`` into the space of real
1613 matrices of size 2n-by-2n via the map the sends each entry `z = a +
1614 bi` to the block matrix ``[[a,b],[-b,a]]``.
1615
1616 SETUP::
1617
1618 sage: from mjo.eja.eja_algebra import ComplexMatrixEJA
1619
1620 EXAMPLES::
1621
1622 sage: F = QuadraticField(-1, 'I')
1623 sage: x1 = F(4 - 2*i)
1624 sage: x2 = F(1 + 2*i)
1625 sage: x3 = F(-i)
1626 sage: x4 = F(6)
1627 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
1628 sage: ComplexMatrixEJA.real_embed(M)
1629 [ 4 -2| 1 2]
1630 [ 2 4|-2 1]
1631 [-----+-----]
1632 [ 0 -1| 6 0]
1633 [ 1 0| 0 6]
1634
1635 TESTS:
1636
1637 Embedding is a homomorphism (isomorphism, in fact)::
1638
1639 sage: set_random_seed()
1640 sage: n = ZZ.random_element(3)
1641 sage: F = QuadraticField(-1, 'I')
1642 sage: X = random_matrix(F, n)
1643 sage: Y = random_matrix(F, n)
1644 sage: Xe = ComplexMatrixEJA.real_embed(X)
1645 sage: Ye = ComplexMatrixEJA.real_embed(Y)
1646 sage: XYe = ComplexMatrixEJA.real_embed(X*Y)
1647 sage: Xe*Ye == XYe
1648 True
1649
1650 """
1651 super(ComplexMatrixEJA,cls).real_embed(M)
1652 n = M.nrows()
1653
1654 # We don't need any adjoined elements...
1655 field = M.base_ring().base_ring()
1656
1657 blocks = []
1658 for z in M.list():
1659 a = z.list()[0] # real part, I guess
1660 b = z.list()[1] # imag part, I guess
1661 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1662
1663 return matrix.block(field, n, blocks)
1664
1665
1666 @classmethod
1667 def real_unembed(cls,M):
1668 """
1669 The inverse of _embed_complex_matrix().
1670
1671 SETUP::
1672
1673 sage: from mjo.eja.eja_algebra import ComplexMatrixEJA
1674
1675 EXAMPLES::
1676
1677 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1678 ....: [-2, 1, -4, 3],
1679 ....: [ 9, 10, 11, 12],
1680 ....: [-10, 9, -12, 11] ])
1681 sage: ComplexMatrixEJA.real_unembed(A)
1682 [ 2*I + 1 4*I + 3]
1683 [ 10*I + 9 12*I + 11]
1684
1685 TESTS:
1686
1687 Unembedding is the inverse of embedding::
1688
1689 sage: set_random_seed()
1690 sage: F = QuadraticField(-1, 'I')
1691 sage: M = random_matrix(F, 3)
1692 sage: Me = ComplexMatrixEJA.real_embed(M)
1693 sage: ComplexMatrixEJA.real_unembed(Me) == M
1694 True
1695
1696 """
1697 super(ComplexMatrixEJA,cls).real_unembed(M)
1698 n = ZZ(M.nrows())
1699 d = cls.dimension_over_reals()
1700
1701 # If "M" was normalized, its base ring might have roots
1702 # adjoined and they can stick around after unembedding.
1703 field = M.base_ring()
1704 R = PolynomialRing(field, 'z')
1705 z = R.gen()
1706
1707 # Sage doesn't know how to adjoin the complex "i" (the root of
1708 # x^2 + 1) to a field in a general way. Here, we just enumerate
1709 # all of the cases that I have cared to support so far.
1710 if field is AA:
1711 # Sage doesn't know how to embed AA into QQbar, i.e. how
1712 # to adjoin sqrt(-1) to AA.
1713 F = QQbar
1714 elif not field.is_exact():
1715 # RDF or RR
1716 F = field.complex_field()
1717 else:
1718 # Works for QQ and... maybe some other fields.
1719 F = field.extension(z**2 + 1, 'I', embedding=CLF(-1).sqrt())
1720 i = F.gen()
1721
1722 # Go top-left to bottom-right (reading order), converting every
1723 # 2-by-2 block we see to a single complex element.
1724 elements = []
1725 for k in range(n/d):
1726 for j in range(n/d):
1727 submat = M[d*k:d*k+d,d*j:d*j+d]
1728 if submat[0,0] != submat[1,1]:
1729 raise ValueError('bad on-diagonal submatrix')
1730 if submat[0,1] != -submat[1,0]:
1731 raise ValueError('bad off-diagonal submatrix')
1732 z = submat[0,0] + submat[0,1]*i
1733 elements.append(z)
1734
1735 return matrix(F, n/d, elements)
1736
1737
1738 class ComplexHermitianEJA(ConcreteEJA, ComplexMatrixEJA):
1739 """
1740 The rank-n simple EJA consisting of complex Hermitian n-by-n
1741 matrices over the real numbers, the usual symmetric Jordan product,
1742 and the real-part-of-trace inner product. It has dimension `n^2` over
1743 the reals.
1744
1745 SETUP::
1746
1747 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1748
1749 EXAMPLES:
1750
1751 In theory, our "field" can be any subfield of the reals::
1752
1753 sage: ComplexHermitianEJA(2, field=RDF)
1754 Euclidean Jordan algebra of dimension 4 over Real Double Field
1755 sage: ComplexHermitianEJA(2, field=RR)
1756 Euclidean Jordan algebra of dimension 4 over Real Field with
1757 53 bits of precision
1758
1759 TESTS:
1760
1761 The dimension of this algebra is `n^2`::
1762
1763 sage: set_random_seed()
1764 sage: n_max = ComplexHermitianEJA._max_random_instance_size()
1765 sage: n = ZZ.random_element(1, n_max)
1766 sage: J = ComplexHermitianEJA(n)
1767 sage: J.dimension() == n^2
1768 True
1769
1770 The Jordan multiplication is what we think it is::
1771
1772 sage: set_random_seed()
1773 sage: J = ComplexHermitianEJA.random_instance()
1774 sage: x,y = J.random_elements(2)
1775 sage: actual = (x*y).to_matrix()
1776 sage: X = x.to_matrix()
1777 sage: Y = y.to_matrix()
1778 sage: expected = (X*Y + Y*X)/2
1779 sage: actual == expected
1780 True
1781 sage: J(expected) == x*y
1782 True
1783
1784 We can change the generator prefix::
1785
1786 sage: ComplexHermitianEJA(2, prefix='z').gens()
1787 (z0, z1, z2, z3)
1788
1789 We can construct the (trivial) algebra of rank zero::
1790
1791 sage: ComplexHermitianEJA(0)
1792 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
1793
1794 """
1795
1796 @classmethod
1797 def _denormalized_basis(cls, n):
1798 """
1799 Returns a basis for the space of complex Hermitian n-by-n matrices.
1800
1801 Why do we embed these? Basically, because all of numerical linear
1802 algebra assumes that you're working with vectors consisting of `n`
1803 entries from a field and scalars from the same field. There's no way
1804 to tell SageMath that (for example) the vectors contain complex
1805 numbers, while the scalar field is real.
1806
1807 SETUP::
1808
1809 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1810
1811 TESTS::
1812
1813 sage: set_random_seed()
1814 sage: n = ZZ.random_element(1,5)
1815 sage: field = QuadraticField(2, 'sqrt2')
1816 sage: B = ComplexHermitianEJA._denormalized_basis(n)
1817 sage: all( M.is_symmetric() for M in B)
1818 True
1819
1820 """
1821 field = ZZ
1822 R = PolynomialRing(field, 'z')
1823 z = R.gen()
1824 F = field.extension(z**2 + 1, 'I')
1825 I = F.gen(1)
1826
1827 # This is like the symmetric case, but we need to be careful:
1828 #
1829 # * We want conjugate-symmetry, not just symmetry.
1830 # * The diagonal will (as a result) be real.
1831 #
1832 S = []
1833 for i in range(n):
1834 for j in range(i+1):
1835 Eij = matrix(F, n, lambda k,l: k==i and l==j)
1836 if i == j:
1837 Sij = cls.real_embed(Eij)
1838 S.append(Sij)
1839 else:
1840 # The second one has a minus because it's conjugated.
1841 Sij_real = cls.real_embed(Eij + Eij.transpose())
1842 S.append(Sij_real)
1843 Sij_imag = cls.real_embed(I*Eij - I*Eij.transpose())
1844 S.append(Sij_imag)
1845
1846 # Since we embedded these, we can drop back to the "field" that we
1847 # started with instead of the complex extension "F".
1848 return tuple( s.change_ring(field) for s in S )
1849
1850
1851 def __init__(self, n, **kwargs):
1852 # We know this is a valid EJA, but will double-check
1853 # if the user passes check_axioms=True.
1854 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
1855
1856 super(ComplexHermitianEJA, self).__init__(self._denormalized_basis(n),
1857 self.jordan_product,
1858 self.trace_inner_product,
1859 **kwargs)
1860 # TODO: this could be factored out somehow, but is left here
1861 # because the MatrixEJA is not presently a subclass of the
1862 # FDEJA class that defines rank() and one().
1863 self.rank.set_cache(n)
1864 idV = matrix.identity(ZZ, self.dimension_over_reals()*n)
1865 self.one.set_cache(self(idV))
1866
1867 @staticmethod
1868 def _max_random_instance_size():
1869 return 3 # Dimension 9
1870
1871 @classmethod
1872 def random_instance(cls, **kwargs):
1873 """
1874 Return a random instance of this type of algebra.
1875 """
1876 n = ZZ.random_element(cls._max_random_instance_size() + 1)
1877 return cls(n, **kwargs)
1878
1879 class QuaternionMatrixEJA(MatrixEJA):
1880 @staticmethod
1881 def dimension_over_reals():
1882 return 4
1883
1884 @classmethod
1885 def real_embed(cls,M):
1886 """
1887 Embed the n-by-n quaternion matrix ``M`` into the space of real
1888 matrices of size 4n-by-4n by first sending each quaternion entry `z
1889 = a + bi + cj + dk` to the block-complex matrix ``[[a + bi,
1890 c+di],[-c + di, a-bi]]`, and then embedding those into a real
1891 matrix.
1892
1893 SETUP::
1894
1895 sage: from mjo.eja.eja_algebra import QuaternionMatrixEJA
1896
1897 EXAMPLES::
1898
1899 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1900 sage: i,j,k = Q.gens()
1901 sage: x = 1 + 2*i + 3*j + 4*k
1902 sage: M = matrix(Q, 1, [[x]])
1903 sage: QuaternionMatrixEJA.real_embed(M)
1904 [ 1 2 3 4]
1905 [-2 1 -4 3]
1906 [-3 4 1 -2]
1907 [-4 -3 2 1]
1908
1909 Embedding is a homomorphism (isomorphism, in fact)::
1910
1911 sage: set_random_seed()
1912 sage: n = ZZ.random_element(2)
1913 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1914 sage: X = random_matrix(Q, n)
1915 sage: Y = random_matrix(Q, n)
1916 sage: Xe = QuaternionMatrixEJA.real_embed(X)
1917 sage: Ye = QuaternionMatrixEJA.real_embed(Y)
1918 sage: XYe = QuaternionMatrixEJA.real_embed(X*Y)
1919 sage: Xe*Ye == XYe
1920 True
1921
1922 """
1923 super(QuaternionMatrixEJA,cls).real_embed(M)
1924 quaternions = M.base_ring()
1925 n = M.nrows()
1926
1927 F = QuadraticField(-1, 'I')
1928 i = F.gen()
1929
1930 blocks = []
1931 for z in M.list():
1932 t = z.coefficient_tuple()
1933 a = t[0]
1934 b = t[1]
1935 c = t[2]
1936 d = t[3]
1937 cplxM = matrix(F, 2, [[ a + b*i, c + d*i],
1938 [-c + d*i, a - b*i]])
1939 realM = ComplexMatrixEJA.real_embed(cplxM)
1940 blocks.append(realM)
1941
1942 # We should have real entries by now, so use the realest field
1943 # we've got for the return value.
1944 return matrix.block(quaternions.base_ring(), n, blocks)
1945
1946
1947
1948 @classmethod
1949 def real_unembed(cls,M):
1950 """
1951 The inverse of _embed_quaternion_matrix().
1952
1953 SETUP::
1954
1955 sage: from mjo.eja.eja_algebra import QuaternionMatrixEJA
1956
1957 EXAMPLES::
1958
1959 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
1960 ....: [-2, 1, -4, 3],
1961 ....: [-3, 4, 1, -2],
1962 ....: [-4, -3, 2, 1]])
1963 sage: QuaternionMatrixEJA.real_unembed(M)
1964 [1 + 2*i + 3*j + 4*k]
1965
1966 TESTS:
1967
1968 Unembedding is the inverse of embedding::
1969
1970 sage: set_random_seed()
1971 sage: Q = QuaternionAlgebra(QQ, -1, -1)
1972 sage: M = random_matrix(Q, 3)
1973 sage: Me = QuaternionMatrixEJA.real_embed(M)
1974 sage: QuaternionMatrixEJA.real_unembed(Me) == M
1975 True
1976
1977 """
1978 super(QuaternionMatrixEJA,cls).real_unembed(M)
1979 n = ZZ(M.nrows())
1980 d = cls.dimension_over_reals()
1981
1982 # Use the base ring of the matrix to ensure that its entries can be
1983 # multiplied by elements of the quaternion algebra.
1984 field = M.base_ring()
1985 Q = QuaternionAlgebra(field,-1,-1)
1986 i,j,k = Q.gens()
1987
1988 # Go top-left to bottom-right (reading order), converting every
1989 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
1990 # quaternion block.
1991 elements = []
1992 for l in range(n/d):
1993 for m in range(n/d):
1994 submat = ComplexMatrixEJA.real_unembed(
1995 M[d*l:d*l+d,d*m:d*m+d] )
1996 if submat[0,0] != submat[1,1].conjugate():
1997 raise ValueError('bad on-diagonal submatrix')
1998 if submat[0,1] != -submat[1,0].conjugate():
1999 raise ValueError('bad off-diagonal submatrix')
2000 z = submat[0,0].real()
2001 z += submat[0,0].imag()*i
2002 z += submat[0,1].real()*j
2003 z += submat[0,1].imag()*k
2004 elements.append(z)
2005
2006 return matrix(Q, n/d, elements)
2007
2008
2009 class QuaternionHermitianEJA(ConcreteEJA, QuaternionMatrixEJA):
2010 r"""
2011 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
2012 matrices, the usual symmetric Jordan product, and the
2013 real-part-of-trace inner product. It has dimension `2n^2 - n` over
2014 the reals.
2015
2016 SETUP::
2017
2018 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
2019
2020 EXAMPLES:
2021
2022 In theory, our "field" can be any subfield of the reals::
2023
2024 sage: QuaternionHermitianEJA(2, field=RDF)
2025 Euclidean Jordan algebra of dimension 6 over Real Double Field
2026 sage: QuaternionHermitianEJA(2, field=RR)
2027 Euclidean Jordan algebra of dimension 6 over Real Field with
2028 53 bits of precision
2029
2030 TESTS:
2031
2032 The dimension of this algebra is `2*n^2 - n`::
2033
2034 sage: set_random_seed()
2035 sage: n_max = QuaternionHermitianEJA._max_random_instance_size()
2036 sage: n = ZZ.random_element(1, n_max)
2037 sage: J = QuaternionHermitianEJA(n)
2038 sage: J.dimension() == 2*(n^2) - n
2039 True
2040
2041 The Jordan multiplication is what we think it is::
2042
2043 sage: set_random_seed()
2044 sage: J = QuaternionHermitianEJA.random_instance()
2045 sage: x,y = J.random_elements(2)
2046 sage: actual = (x*y).to_matrix()
2047 sage: X = x.to_matrix()
2048 sage: Y = y.to_matrix()
2049 sage: expected = (X*Y + Y*X)/2
2050 sage: actual == expected
2051 True
2052 sage: J(expected) == x*y
2053 True
2054
2055 We can change the generator prefix::
2056
2057 sage: QuaternionHermitianEJA(2, prefix='a').gens()
2058 (a0, a1, a2, a3, a4, a5)
2059
2060 We can construct the (trivial) algebra of rank zero::
2061
2062 sage: QuaternionHermitianEJA(0)
2063 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
2064
2065 """
2066 @classmethod
2067 def _denormalized_basis(cls, n):
2068 """
2069 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
2070
2071 Why do we embed these? Basically, because all of numerical
2072 linear algebra assumes that you're working with vectors consisting
2073 of `n` entries from a field and scalars from the same field. There's
2074 no way to tell SageMath that (for example) the vectors contain
2075 complex numbers, while the scalar field is real.
2076
2077 SETUP::
2078
2079 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
2080
2081 TESTS::
2082
2083 sage: set_random_seed()
2084 sage: n = ZZ.random_element(1,5)
2085 sage: B = QuaternionHermitianEJA._denormalized_basis(n)
2086 sage: all( M.is_symmetric() for M in B )
2087 True
2088
2089 """
2090 field = ZZ
2091 Q = QuaternionAlgebra(QQ,-1,-1)
2092 I,J,K = Q.gens()
2093
2094 # This is like the symmetric case, but we need to be careful:
2095 #
2096 # * We want conjugate-symmetry, not just symmetry.
2097 # * The diagonal will (as a result) be real.
2098 #
2099 S = []
2100 for i in range(n):
2101 for j in range(i+1):
2102 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
2103 if i == j:
2104 Sij = cls.real_embed(Eij)
2105 S.append(Sij)
2106 else:
2107 # The second, third, and fourth ones have a minus
2108 # because they're conjugated.
2109 Sij_real = cls.real_embed(Eij + Eij.transpose())
2110 S.append(Sij_real)
2111 Sij_I = cls.real_embed(I*Eij - I*Eij.transpose())
2112 S.append(Sij_I)
2113 Sij_J = cls.real_embed(J*Eij - J*Eij.transpose())
2114 S.append(Sij_J)
2115 Sij_K = cls.real_embed(K*Eij - K*Eij.transpose())
2116 S.append(Sij_K)
2117
2118 # Since we embedded these, we can drop back to the "field" that we
2119 # started with instead of the quaternion algebra "Q".
2120 return tuple( s.change_ring(field) for s in S )
2121
2122
2123 def __init__(self, n, **kwargs):
2124 # We know this is a valid EJA, but will double-check
2125 # if the user passes check_axioms=True.
2126 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
2127
2128 super(QuaternionHermitianEJA, self).__init__(self._denormalized_basis(n),
2129 self.jordan_product,
2130 self.trace_inner_product,
2131 **kwargs)
2132 # TODO: this could be factored out somehow, but is left here
2133 # because the MatrixEJA is not presently a subclass of the
2134 # FDEJA class that defines rank() and one().
2135 self.rank.set_cache(n)
2136 idV = matrix.identity(ZZ, self.dimension_over_reals()*n)
2137 self.one.set_cache(self(idV))
2138
2139
2140 @staticmethod
2141 def _max_random_instance_size():
2142 r"""
2143 The maximum rank of a random QuaternionHermitianEJA.
2144 """
2145 return 2 # Dimension 6
2146
2147 @classmethod
2148 def random_instance(cls, **kwargs):
2149 """
2150 Return a random instance of this type of algebra.
2151 """
2152 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2153 return cls(n, **kwargs)
2154
2155
2156 class HadamardEJA(ConcreteEJA):
2157 """
2158 Return the Euclidean Jordan Algebra corresponding to the set
2159 `R^n` under the Hadamard product.
2160
2161 Note: this is nothing more than the Cartesian product of ``n``
2162 copies of the spin algebra. Once Cartesian product algebras
2163 are implemented, this can go.
2164
2165 SETUP::
2166
2167 sage: from mjo.eja.eja_algebra import HadamardEJA
2168
2169 EXAMPLES:
2170
2171 This multiplication table can be verified by hand::
2172
2173 sage: J = HadamardEJA(3)
2174 sage: e0,e1,e2 = J.gens()
2175 sage: e0*e0
2176 e0
2177 sage: e0*e1
2178 0
2179 sage: e0*e2
2180 0
2181 sage: e1*e1
2182 e1
2183 sage: e1*e2
2184 0
2185 sage: e2*e2
2186 e2
2187
2188 TESTS:
2189
2190 We can change the generator prefix::
2191
2192 sage: HadamardEJA(3, prefix='r').gens()
2193 (r0, r1, r2)
2194
2195 """
2196 def __init__(self, n, **kwargs):
2197 if n == 0:
2198 jordan_product = lambda x,y: x
2199 inner_product = lambda x,y: x
2200 else:
2201 def jordan_product(x,y):
2202 P = x.parent()
2203 return P( xi*yi for (xi,yi) in zip(x,y) )
2204
2205 def inner_product(x,y):
2206 return (x.T*y)[0,0]
2207
2208 # New defaults for keyword arguments. Don't orthonormalize
2209 # because our basis is already orthonormal with respect to our
2210 # inner-product. Don't check the axioms, because we know this
2211 # is a valid EJA... but do double-check if the user passes
2212 # check_axioms=True. Note: we DON'T override the "check_field"
2213 # default here, because the user can pass in a field!
2214 if "orthonormalize" not in kwargs: kwargs["orthonormalize"] = False
2215 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
2216
2217 column_basis = tuple( b.column() for b in FreeModule(ZZ, n).basis() )
2218 super().__init__(column_basis, jordan_product, inner_product, **kwargs)
2219 self.rank.set_cache(n)
2220
2221 if n == 0:
2222 self.one.set_cache( self.zero() )
2223 else:
2224 self.one.set_cache( sum(self.gens()) )
2225
2226 @staticmethod
2227 def _max_random_instance_size():
2228 r"""
2229 The maximum dimension of a random HadamardEJA.
2230 """
2231 return 5
2232
2233 @classmethod
2234 def random_instance(cls, **kwargs):
2235 """
2236 Return a random instance of this type of algebra.
2237 """
2238 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2239 return cls(n, **kwargs)
2240
2241
2242 class BilinearFormEJA(ConcreteEJA):
2243 r"""
2244 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2245 with the half-trace inner product and jordan product ``x*y =
2246 (<Bx,y>,y_bar>, x0*y_bar + y0*x_bar)`` where `B = 1 \times B22` is
2247 a symmetric positive-definite "bilinear form" matrix. Its
2248 dimension is the size of `B`, and it has rank two in dimensions
2249 larger than two. It reduces to the ``JordanSpinEJA`` when `B` is
2250 the identity matrix of order ``n``.
2251
2252 We insist that the one-by-one upper-left identity block of `B` be
2253 passed in as well so that we can be passed a matrix of size zero
2254 to construct a trivial algebra.
2255
2256 SETUP::
2257
2258 sage: from mjo.eja.eja_algebra import (BilinearFormEJA,
2259 ....: JordanSpinEJA)
2260
2261 EXAMPLES:
2262
2263 When no bilinear form is specified, the identity matrix is used,
2264 and the resulting algebra is the Jordan spin algebra::
2265
2266 sage: B = matrix.identity(AA,3)
2267 sage: J0 = BilinearFormEJA(B)
2268 sage: J1 = JordanSpinEJA(3)
2269 sage: J0.multiplication_table() == J0.multiplication_table()
2270 True
2271
2272 An error is raised if the matrix `B` does not correspond to a
2273 positive-definite bilinear form::
2274
2275 sage: B = matrix.random(QQ,2,3)
2276 sage: J = BilinearFormEJA(B)
2277 Traceback (most recent call last):
2278 ...
2279 ValueError: bilinear form is not positive-definite
2280 sage: B = matrix.zero(QQ,3)
2281 sage: J = BilinearFormEJA(B)
2282 Traceback (most recent call last):
2283 ...
2284 ValueError: bilinear form is not positive-definite
2285
2286 TESTS:
2287
2288 We can create a zero-dimensional algebra::
2289
2290 sage: B = matrix.identity(AA,0)
2291 sage: J = BilinearFormEJA(B)
2292 sage: J.basis()
2293 Finite family {}
2294
2295 We can check the multiplication condition given in the Jordan, von
2296 Neumann, and Wigner paper (and also discussed on my "On the
2297 symmetry..." paper). Note that this relies heavily on the standard
2298 choice of basis, as does anything utilizing the bilinear form
2299 matrix. We opt not to orthonormalize the basis, because if we
2300 did, we would have to normalize the `s_{i}` in a similar manner::
2301
2302 sage: set_random_seed()
2303 sage: n = ZZ.random_element(5)
2304 sage: M = matrix.random(QQ, max(0,n-1), algorithm='unimodular')
2305 sage: B11 = matrix.identity(QQ,1)
2306 sage: B22 = M.transpose()*M
2307 sage: B = block_matrix(2,2,[ [B11,0 ],
2308 ....: [0, B22 ] ])
2309 sage: J = BilinearFormEJA(B, orthonormalize=False)
2310 sage: eis = VectorSpace(M.base_ring(), M.ncols()).basis()
2311 sage: V = J.vector_space()
2312 sage: sis = [ J( V([0] + (M.inverse()*ei).list()).column() )
2313 ....: for ei in eis ]
2314 sage: actual = [ sis[i]*sis[j]
2315 ....: for i in range(n-1)
2316 ....: for j in range(n-1) ]
2317 sage: expected = [ J.one() if i == j else J.zero()
2318 ....: for i in range(n-1)
2319 ....: for j in range(n-1) ]
2320 sage: actual == expected
2321 True
2322
2323 """
2324 def __init__(self, B, **kwargs):
2325 # The matrix "B" is supplied by the user in most cases,
2326 # so it makes sense to check whether or not its positive-
2327 # definite unless we are specifically asked not to...
2328 if ("check_axioms" not in kwargs) or kwargs["check_axioms"]:
2329 if not B.is_positive_definite():
2330 raise ValueError("bilinear form is not positive-definite")
2331
2332 # However, all of the other data for this EJA is computed
2333 # by us in manner that guarantees the axioms are
2334 # satisfied. So, again, unless we are specifically asked to
2335 # verify things, we'll skip the rest of the checks.
2336 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
2337
2338 def inner_product(x,y):
2339 return (y.T*B*x)[0,0]
2340
2341 def jordan_product(x,y):
2342 P = x.parent()
2343 x0 = x[0,0]
2344 xbar = x[1:,0]
2345 y0 = y[0,0]
2346 ybar = y[1:,0]
2347 z0 = inner_product(y,x)
2348 zbar = y0*xbar + x0*ybar
2349 return P([z0] + zbar.list())
2350
2351 n = B.nrows()
2352 column_basis = tuple( b.column() for b in FreeModule(ZZ, n).basis() )
2353 super(BilinearFormEJA, self).__init__(column_basis,
2354 jordan_product,
2355 inner_product,
2356 **kwargs)
2357
2358 # The rank of this algebra is two, unless we're in a
2359 # one-dimensional ambient space (because the rank is bounded
2360 # by the ambient dimension).
2361 self.rank.set_cache(min(n,2))
2362
2363 if n == 0:
2364 self.one.set_cache( self.zero() )
2365 else:
2366 self.one.set_cache( self.monomial(0) )
2367
2368 @staticmethod
2369 def _max_random_instance_size():
2370 r"""
2371 The maximum dimension of a random BilinearFormEJA.
2372 """
2373 return 5
2374
2375 @classmethod
2376 def random_instance(cls, **kwargs):
2377 """
2378 Return a random instance of this algebra.
2379 """
2380 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2381 if n.is_zero():
2382 B = matrix.identity(ZZ, n)
2383 return cls(B, **kwargs)
2384
2385 B11 = matrix.identity(ZZ, 1)
2386 M = matrix.random(ZZ, n-1)
2387 I = matrix.identity(ZZ, n-1)
2388 alpha = ZZ.zero()
2389 while alpha.is_zero():
2390 alpha = ZZ.random_element().abs()
2391 B22 = M.transpose()*M + alpha*I
2392
2393 from sage.matrix.special import block_matrix
2394 B = block_matrix(2,2, [ [B11, ZZ(0) ],
2395 [ZZ(0), B22 ] ])
2396
2397 return cls(B, **kwargs)
2398
2399
2400 class JordanSpinEJA(BilinearFormEJA):
2401 """
2402 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2403 with the usual inner product and jordan product ``x*y =
2404 (<x,y>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
2405 the reals.
2406
2407 SETUP::
2408
2409 sage: from mjo.eja.eja_algebra import JordanSpinEJA
2410
2411 EXAMPLES:
2412
2413 This multiplication table can be verified by hand::
2414
2415 sage: J = JordanSpinEJA(4)
2416 sage: e0,e1,e2,e3 = J.gens()
2417 sage: e0*e0
2418 e0
2419 sage: e0*e1
2420 e1
2421 sage: e0*e2
2422 e2
2423 sage: e0*e3
2424 e3
2425 sage: e1*e2
2426 0
2427 sage: e1*e3
2428 0
2429 sage: e2*e3
2430 0
2431
2432 We can change the generator prefix::
2433
2434 sage: JordanSpinEJA(2, prefix='B').gens()
2435 (B0, B1)
2436
2437 TESTS:
2438
2439 Ensure that we have the usual inner product on `R^n`::
2440
2441 sage: set_random_seed()
2442 sage: J = JordanSpinEJA.random_instance()
2443 sage: x,y = J.random_elements(2)
2444 sage: actual = x.inner_product(y)
2445 sage: expected = x.to_vector().inner_product(y.to_vector())
2446 sage: actual == expected
2447 True
2448
2449 """
2450 def __init__(self, n, **kwargs):
2451 # This is a special case of the BilinearFormEJA with the
2452 # identity matrix as its bilinear form.
2453 B = matrix.identity(ZZ, n)
2454
2455 # Don't orthonormalize because our basis is already
2456 # orthonormal with respect to our inner-product.
2457 if "orthonormalize" not in kwargs: kwargs["orthonormalize"] = False
2458
2459 # But also don't pass check_field=False here, because the user
2460 # can pass in a field!
2461 super(JordanSpinEJA, self).__init__(B, **kwargs)
2462
2463 @staticmethod
2464 def _max_random_instance_size():
2465 r"""
2466 The maximum dimension of a random JordanSpinEJA.
2467 """
2468 return 5
2469
2470 @classmethod
2471 def random_instance(cls, **kwargs):
2472 """
2473 Return a random instance of this type of algebra.
2474
2475 Needed here to override the implementation for ``BilinearFormEJA``.
2476 """
2477 n = ZZ.random_element(cls._max_random_instance_size() + 1)
2478 return cls(n, **kwargs)
2479
2480
2481 class TrivialEJA(ConcreteEJA):
2482 """
2483 The trivial Euclidean Jordan algebra consisting of only a zero element.
2484
2485 SETUP::
2486
2487 sage: from mjo.eja.eja_algebra import TrivialEJA
2488
2489 EXAMPLES::
2490
2491 sage: J = TrivialEJA()
2492 sage: J.dimension()
2493 0
2494 sage: J.zero()
2495 0
2496 sage: J.one()
2497 0
2498 sage: 7*J.one()*12*J.one()
2499 0
2500 sage: J.one().inner_product(J.one())
2501 0
2502 sage: J.one().norm()
2503 0
2504 sage: J.one().subalgebra_generated_by()
2505 Euclidean Jordan algebra of dimension 0 over Algebraic Real Field
2506 sage: J.rank()
2507 0
2508
2509 """
2510 def __init__(self, **kwargs):
2511 jordan_product = lambda x,y: x
2512 inner_product = lambda x,y: 0
2513 basis = ()
2514
2515 # New defaults for keyword arguments
2516 if "orthonormalize" not in kwargs: kwargs["orthonormalize"] = False
2517 if "check_axioms" not in kwargs: kwargs["check_axioms"] = False
2518
2519 super(TrivialEJA, self).__init__(basis,
2520 jordan_product,
2521 inner_product,
2522 **kwargs)
2523 # The rank is zero using my definition, namely the dimension of the
2524 # largest subalgebra generated by any element.
2525 self.rank.set_cache(0)
2526 self.one.set_cache( self.zero() )
2527
2528 @classmethod
2529 def random_instance(cls, **kwargs):
2530 # We don't take a "size" argument so the superclass method is
2531 # inappropriate for us.
2532 return cls(**kwargs)
2533
2534 # class DirectSumEJA(ConcreteEJA):
2535 # r"""
2536 # The external (orthogonal) direct sum of two other Euclidean Jordan
2537 # algebras. Essentially the Cartesian product of its two factors.
2538 # Every Euclidean Jordan algebra decomposes into an orthogonal
2539 # direct sum of simple Euclidean Jordan algebras, so no generality
2540 # is lost by providing only this construction.
2541
2542 # SETUP::
2543
2544 # sage: from mjo.eja.eja_algebra import (random_eja,
2545 # ....: HadamardEJA,
2546 # ....: RealSymmetricEJA,
2547 # ....: DirectSumEJA)
2548
2549 # EXAMPLES::
2550
2551 # sage: J1 = HadamardEJA(2)
2552 # sage: J2 = RealSymmetricEJA(3)
2553 # sage: J = DirectSumEJA(J1,J2)
2554 # sage: J.dimension()
2555 # 8
2556 # sage: J.rank()
2557 # 5
2558
2559 # TESTS:
2560
2561 # The external direct sum construction is only valid when the two factors
2562 # have the same base ring; an error is raised otherwise::
2563
2564 # sage: set_random_seed()
2565 # sage: J1 = random_eja(field=AA)
2566 # sage: J2 = random_eja(field=QQ,orthonormalize=False)
2567 # sage: J = DirectSumEJA(J1,J2)
2568 # Traceback (most recent call last):
2569 # ...
2570 # ValueError: algebras must share the same base field
2571
2572 # """
2573 # def __init__(self, J1, J2, **kwargs):
2574 # if J1.base_ring() != J2.base_ring():
2575 # raise ValueError("algebras must share the same base field")
2576 # field = J1.base_ring()
2577
2578 # self._factors = (J1, J2)
2579 # n1 = J1.dimension()
2580 # n2 = J2.dimension()
2581 # n = n1+n2
2582 # V = VectorSpace(field, n)
2583 # mult_table = [ [ V.zero() for j in range(i+1) ]
2584 # for i in range(n) ]
2585 # for i in range(n1):
2586 # for j in range(i+1):
2587 # p = (J1.monomial(i)*J1.monomial(j)).to_vector()
2588 # mult_table[i][j] = V(p.list() + [field.zero()]*n2)
2589
2590 # for i in range(n2):
2591 # for j in range(i+1):
2592 # p = (J2.monomial(i)*J2.monomial(j)).to_vector()
2593 # mult_table[n1+i][n1+j] = V([field.zero()]*n1 + p.list())
2594
2595 # # TODO: build the IP table here from the two constituent IP
2596 # # matrices (it'll be block diagonal, I think).
2597 # ip_table = [ [ field.zero() for j in range(i+1) ]
2598 # for i in range(n) ]
2599 # super(DirectSumEJA, self).__init__(field,
2600 # mult_table,
2601 # ip_table,
2602 # check_axioms=False,
2603 # **kwargs)
2604 # self.rank.set_cache(J1.rank() + J2.rank())
2605
2606
2607 # def factors(self):
2608 # r"""
2609 # Return the pair of this algebra's factors.
2610
2611 # SETUP::
2612
2613 # sage: from mjo.eja.eja_algebra import (HadamardEJA,
2614 # ....: JordanSpinEJA,
2615 # ....: DirectSumEJA)
2616
2617 # EXAMPLES::
2618
2619 # sage: J1 = HadamardEJA(2, field=QQ)
2620 # sage: J2 = JordanSpinEJA(3, field=QQ)
2621 # sage: J = DirectSumEJA(J1,J2)
2622 # sage: J.factors()
2623 # (Euclidean Jordan algebra of dimension 2 over Rational Field,
2624 # Euclidean Jordan algebra of dimension 3 over Rational Field)
2625
2626 # """
2627 # return self._factors
2628
2629 # def projections(self):
2630 # r"""
2631 # Return a pair of projections onto this algebra's factors.
2632
2633 # SETUP::
2634
2635 # sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
2636 # ....: ComplexHermitianEJA,
2637 # ....: DirectSumEJA)
2638
2639 # EXAMPLES::
2640
2641 # sage: J1 = JordanSpinEJA(2)
2642 # sage: J2 = ComplexHermitianEJA(2)
2643 # sage: J = DirectSumEJA(J1,J2)
2644 # sage: (pi_left, pi_right) = J.projections()
2645 # sage: J.one().to_vector()
2646 # (1, 0, 1, 0, 0, 1)
2647 # sage: pi_left(J.one()).to_vector()
2648 # (1, 0)
2649 # sage: pi_right(J.one()).to_vector()
2650 # (1, 0, 0, 1)
2651
2652 # """
2653 # (J1,J2) = self.factors()
2654 # m = J1.dimension()
2655 # n = J2.dimension()
2656 # V_basis = self.vector_space().basis()
2657 # # Need to specify the dimensions explicitly so that we don't
2658 # # wind up with a zero-by-zero matrix when we want e.g. a
2659 # # zero-by-two matrix (important for composing things).
2660 # P1 = matrix(self.base_ring(), m, m+n, V_basis[:m])
2661 # P2 = matrix(self.base_ring(), n, m+n, V_basis[m:])
2662 # pi_left = FiniteDimensionalEJAOperator(self,J1,P1)
2663 # pi_right = FiniteDimensionalEJAOperator(self,J2,P2)
2664 # return (pi_left, pi_right)
2665
2666 # def inclusions(self):
2667 # r"""
2668 # Return the pair of inclusion maps from our factors into us.
2669
2670 # SETUP::
2671
2672 # sage: from mjo.eja.eja_algebra import (random_eja,
2673 # ....: JordanSpinEJA,
2674 # ....: RealSymmetricEJA,
2675 # ....: DirectSumEJA)
2676
2677 # EXAMPLES::
2678
2679 # sage: J1 = JordanSpinEJA(3)
2680 # sage: J2 = RealSymmetricEJA(2)
2681 # sage: J = DirectSumEJA(J1,J2)
2682 # sage: (iota_left, iota_right) = J.inclusions()
2683 # sage: iota_left(J1.zero()) == J.zero()
2684 # True
2685 # sage: iota_right(J2.zero()) == J.zero()
2686 # True
2687 # sage: J1.one().to_vector()
2688 # (1, 0, 0)
2689 # sage: iota_left(J1.one()).to_vector()
2690 # (1, 0, 0, 0, 0, 0)
2691 # sage: J2.one().to_vector()
2692 # (1, 0, 1)
2693 # sage: iota_right(J2.one()).to_vector()
2694 # (0, 0, 0, 1, 0, 1)
2695 # sage: J.one().to_vector()
2696 # (1, 0, 0, 1, 0, 1)
2697
2698 # TESTS:
2699
2700 # Composing a projection with the corresponding inclusion should
2701 # produce the identity map, and mismatching them should produce
2702 # the zero map::
2703
2704 # sage: set_random_seed()
2705 # sage: J1 = random_eja()
2706 # sage: J2 = random_eja()
2707 # sage: J = DirectSumEJA(J1,J2)
2708 # sage: (iota_left, iota_right) = J.inclusions()
2709 # sage: (pi_left, pi_right) = J.projections()
2710 # sage: pi_left*iota_left == J1.one().operator()
2711 # True
2712 # sage: pi_right*iota_right == J2.one().operator()
2713 # True
2714 # sage: (pi_left*iota_right).is_zero()
2715 # True
2716 # sage: (pi_right*iota_left).is_zero()
2717 # True
2718
2719 # """
2720 # (J1,J2) = self.factors()
2721 # m = J1.dimension()
2722 # n = J2.dimension()
2723 # V_basis = self.vector_space().basis()
2724 # # Need to specify the dimensions explicitly so that we don't
2725 # # wind up with a zero-by-zero matrix when we want e.g. a
2726 # # two-by-zero matrix (important for composing things).
2727 # I1 = matrix.column(self.base_ring(), m, m+n, V_basis[:m])
2728 # I2 = matrix.column(self.base_ring(), n, m+n, V_basis[m:])
2729 # iota_left = FiniteDimensionalEJAOperator(J1,self,I1)
2730 # iota_right = FiniteDimensionalEJAOperator(J2,self,I2)
2731 # return (iota_left, iota_right)
2732
2733 # def inner_product(self, x, y):
2734 # r"""
2735 # The standard Cartesian inner-product.
2736
2737 # We project ``x`` and ``y`` onto our factors, and add up the
2738 # inner-products from the subalgebras.
2739
2740 # SETUP::
2741
2742
2743 # sage: from mjo.eja.eja_algebra import (HadamardEJA,
2744 # ....: QuaternionHermitianEJA,
2745 # ....: DirectSumEJA)
2746
2747 # EXAMPLE::
2748
2749 # sage: J1 = HadamardEJA(3,field=QQ)
2750 # sage: J2 = QuaternionHermitianEJA(2,field=QQ,orthonormalize=False)
2751 # sage: J = DirectSumEJA(J1,J2)
2752 # sage: x1 = J1.one()
2753 # sage: x2 = x1
2754 # sage: y1 = J2.one()
2755 # sage: y2 = y1
2756 # sage: x1.inner_product(x2)
2757 # 3
2758 # sage: y1.inner_product(y2)
2759 # 2
2760 # sage: J.one().inner_product(J.one())
2761 # 5
2762
2763 # """
2764 # (pi_left, pi_right) = self.projections()
2765 # x1 = pi_left(x)
2766 # x2 = pi_right(x)
2767 # y1 = pi_left(y)
2768 # y2 = pi_right(y)
2769
2770 # return (x1.inner_product(y1) + x2.inner_product(y2))
2771
2772
2773
2774 random_eja = ConcreteEJA.random_instance