]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: move away from using matrices as our "multiplication table."
[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 from sage.algebras.quatalg.quaternion_algebra import QuaternionAlgebra
9 from sage.categories.finite_dimensional_algebras_with_basis import FiniteDimensionalAlgebrasWithBasis
10 from sage.combinat.free_module import CombinatorialFreeModule
11 from sage.matrix.constructor import matrix
12 from sage.misc.cachefunc import cached_method
13 from sage.misc.prandom import choice
14 from sage.modules.free_module import VectorSpace
15 from sage.rings.integer_ring import ZZ
16 from sage.rings.number_field.number_field import QuadraticField
17 from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing
18 from sage.rings.rational_field import QQ
19 from sage.structure.element import is_Matrix
20
21 from mjo.eja.eja_element import FiniteDimensionalEuclideanJordanAlgebraElement
22 from mjo.eja.eja_utils import _mat2vec
23
24 class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule):
25 def __init__(self,
26 field,
27 mult_table,
28 rank,
29 prefix='e',
30 category=None,
31 natural_basis=None):
32 """
33 SETUP::
34
35 sage: from mjo.eja.eja_algebra import random_eja
36
37 EXAMPLES:
38
39 By definition, Jordan multiplication commutes::
40
41 sage: set_random_seed()
42 sage: J = random_eja()
43 sage: x = J.random_element()
44 sage: y = J.random_element()
45 sage: x*y == y*x
46 True
47
48 """
49 self._rank = rank
50 self._natural_basis = natural_basis
51
52 if category is None:
53 category = FiniteDimensionalAlgebrasWithBasis(field).Unital()
54 fda = super(FiniteDimensionalEuclideanJordanAlgebra, self)
55 fda.__init__(field,
56 range(len(mult_table)),
57 prefix=prefix,
58 category=category)
59 self.print_options(bracket='')
60
61 # The multiplication table we're given is necessarily in terms
62 # of vectors, because we don't have an algebra yet for
63 # anything to be an element of. However, it's faster in the
64 # long run to have the multiplication table be in terms of
65 # algebra elements. We do this after calling the superclass
66 # constructor so that from_vector() knows what to do.
67 self._multiplication_table = [ map(lambda x: self.from_vector(x), ls)
68 for ls in mult_table ]
69
70
71 def _element_constructor_(self, elt):
72 """
73 Construct an element of this algebra from its natural
74 representation.
75
76 This gets called only after the parent element _call_ method
77 fails to find a coercion for the argument.
78
79 SETUP::
80
81 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
82 ....: RealCartesianProductEJA,
83 ....: RealSymmetricEJA)
84
85 EXAMPLES:
86
87 The identity in `S^n` is converted to the identity in the EJA::
88
89 sage: J = RealSymmetricEJA(3)
90 sage: I = matrix.identity(QQ,3)
91 sage: J(I) == J.one()
92 True
93
94 This skew-symmetric matrix can't be represented in the EJA::
95
96 sage: J = RealSymmetricEJA(3)
97 sage: A = matrix(QQ,3, lambda i,j: i-j)
98 sage: J(A)
99 Traceback (most recent call last):
100 ...
101 ArithmeticError: vector is not in free module
102
103 TESTS:
104
105 Ensure that we can convert any element of the two non-matrix
106 simple algebras (whose natural representations are their usual
107 vector representations) back and forth faithfully::
108
109 sage: set_random_seed()
110 sage: J = RealCartesianProductEJA(5)
111 sage: x = J.random_element()
112 sage: J(x.to_vector().column()) == x
113 True
114 sage: J = JordanSpinEJA(5)
115 sage: x = J.random_element()
116 sage: J(x.to_vector().column()) == x
117 True
118
119 """
120 natural_basis = self.natural_basis()
121 if elt not in natural_basis[0].matrix_space():
122 raise ValueError("not a naturally-represented algebra element")
123
124 # Thanks for nothing! Matrix spaces aren't vector
125 # spaces in Sage, so we have to figure out its
126 # natural-basis coordinates ourselves.
127 V = VectorSpace(elt.base_ring(), elt.nrows()*elt.ncols())
128 W = V.span_of_basis( _mat2vec(s) for s in natural_basis )
129 coords = W.coordinate_vector(_mat2vec(elt))
130 return self.from_vector(coords)
131
132
133 def _repr_(self):
134 """
135 Return a string representation of ``self``.
136
137 SETUP::
138
139 sage: from mjo.eja.eja_algebra import JordanSpinEJA
140
141 TESTS:
142
143 Ensure that it says what we think it says::
144
145 sage: JordanSpinEJA(2, field=QQ)
146 Euclidean Jordan algebra of degree 2 over Rational Field
147 sage: JordanSpinEJA(3, field=RDF)
148 Euclidean Jordan algebra of degree 3 over Real Double Field
149
150 """
151 # TODO: change this to say "dimension" and fix all the tests.
152 fmt = "Euclidean Jordan algebra of degree {} over {}"
153 return fmt.format(self.dimension(), self.base_ring())
154
155 def product_on_basis(self, i, j):
156 return self._multiplication_table[i][j]
157
158 def _a_regular_element(self):
159 """
160 Guess a regular element. Needed to compute the basis for our
161 characteristic polynomial coefficients.
162
163 SETUP::
164
165 sage: from mjo.eja.eja_algebra import random_eja
166
167 TESTS:
168
169 Ensure that this hacky method succeeds for every algebra that we
170 know how to construct::
171
172 sage: set_random_seed()
173 sage: J = random_eja()
174 sage: J._a_regular_element().is_regular()
175 True
176
177 """
178 gs = self.gens()
179 z = self.sum( (i+1)*gs[i] for i in range(len(gs)) )
180 if not z.is_regular():
181 raise ValueError("don't know a regular element")
182 return z
183
184
185 @cached_method
186 def _charpoly_basis_space(self):
187 """
188 Return the vector space spanned by the basis used in our
189 characteristic polynomial coefficients. This is used not only to
190 compute those coefficients, but also any time we need to
191 evaluate the coefficients (like when we compute the trace or
192 determinant).
193 """
194 z = self._a_regular_element()
195 V = self.vector_space()
196 V1 = V.span_of_basis( (z**k).to_vector() for k in range(self.rank()) )
197 b = (V1.basis() + V1.complement().basis())
198 return V.span_of_basis(b)
199
200
201 @cached_method
202 def _charpoly_coeff(self, i):
203 """
204 Return the coefficient polynomial "a_{i}" of this algebra's
205 general characteristic polynomial.
206
207 Having this be a separate cached method lets us compute and
208 store the trace/determinant (a_{r-1} and a_{0} respectively)
209 separate from the entire characteristic polynomial.
210 """
211 (A_of_x, x, xr, detA) = self._charpoly_matrix_system()
212 R = A_of_x.base_ring()
213 if i >= self.rank():
214 # Guaranteed by theory
215 return R.zero()
216
217 # Danger: the in-place modification is done for performance
218 # reasons (reconstructing a matrix with huge polynomial
219 # entries is slow), but I don't know how cached_method works,
220 # so it's highly possible that we're modifying some global
221 # list variable by reference, here. In other words, you
222 # probably shouldn't call this method twice on the same
223 # algebra, at the same time, in two threads
224 Ai_orig = A_of_x.column(i)
225 A_of_x.set_column(i,xr)
226 numerator = A_of_x.det()
227 A_of_x.set_column(i,Ai_orig)
228
229 # We're relying on the theory here to ensure that each a_i is
230 # indeed back in R, and the added negative signs are to make
231 # the whole charpoly expression sum to zero.
232 return R(-numerator/detA)
233
234
235 @cached_method
236 def _charpoly_matrix_system(self):
237 """
238 Compute the matrix whose entries A_ij are polynomials in
239 X1,...,XN, the vector ``x`` of variables X1,...,XN, the vector
240 corresponding to `x^r` and the determinent of the matrix A =
241 [A_ij]. In other words, all of the fixed (cachable) data needed
242 to compute the coefficients of the characteristic polynomial.
243 """
244 r = self.rank()
245 n = self.dimension()
246
247 # Construct a new algebra over a multivariate polynomial ring...
248 names = tuple('X' + str(i) for i in range(1,n+1))
249 R = PolynomialRing(self.base_ring(), names)
250 # Hack around the fact that our multiplication table is in terms of
251 # algebra elements but the constructor wants it in terms of vectors.
252 vmt = [ tuple(map(lambda x: x.to_vector(), ls))
253 for ls in self._multiplication_table ]
254 J = FiniteDimensionalEuclideanJordanAlgebra(R, tuple(vmt), r)
255
256 idmat = matrix.identity(J.base_ring(), n)
257
258 W = self._charpoly_basis_space()
259 W = W.change_ring(R.fraction_field())
260
261 # Starting with the standard coordinates x = (X1,X2,...,Xn)
262 # and then converting the entries to W-coordinates allows us
263 # to pass in the standard coordinates to the charpoly and get
264 # back the right answer. Specifically, with x = (X1,X2,...,Xn),
265 # we have
266 #
267 # W.coordinates(x^2) eval'd at (standard z-coords)
268 # =
269 # W-coords of (z^2)
270 # =
271 # W-coords of (standard coords of x^2 eval'd at std-coords of z)
272 #
273 # We want the middle equivalent thing in our matrix, but use
274 # the first equivalent thing instead so that we can pass in
275 # standard coordinates.
276 x = J.from_vector(W(R.gens()))
277
278 # Handle the zeroth power separately, because computing
279 # the unit element in J is mathematically suspect.
280 x0 = W.coordinate_vector(self.one().to_vector())
281 l1 = [ x0.column() ]
282 l1 += [ W.coordinate_vector((x**k).to_vector()).column()
283 for k in range(1,r) ]
284 l2 = [idmat.column(k-1).column() for k in range(r+1, n+1)]
285 A_of_x = matrix.block(R, 1, n, (l1 + l2))
286 xr = W.coordinate_vector((x**r).to_vector())
287 return (A_of_x, x, xr, A_of_x.det())
288
289
290 @cached_method
291 def characteristic_polynomial(self):
292 """
293 Return a characteristic polynomial that works for all elements
294 of this algebra.
295
296 The resulting polynomial has `n+1` variables, where `n` is the
297 dimension of this algebra. The first `n` variables correspond to
298 the coordinates of an algebra element: when evaluated at the
299 coordinates of an algebra element with respect to a certain
300 basis, the result is a univariate polynomial (in the one
301 remaining variable ``t``), namely the characteristic polynomial
302 of that element.
303
304 SETUP::
305
306 sage: from mjo.eja.eja_algebra import JordanSpinEJA
307
308 EXAMPLES:
309
310 The characteristic polynomial in the spin algebra is given in
311 Alizadeh, Example 11.11::
312
313 sage: J = JordanSpinEJA(3)
314 sage: p = J.characteristic_polynomial(); p
315 X1^2 - X2^2 - X3^2 + (-2*t)*X1 + t^2
316 sage: xvec = J.one().to_vector()
317 sage: p(*xvec)
318 t^2 - 2*t + 1
319
320 """
321 r = self.rank()
322 n = self.dimension()
323
324 # The list of coefficient polynomials a_1, a_2, ..., a_n.
325 a = [ self._charpoly_coeff(i) for i in range(n) ]
326
327 # We go to a bit of trouble here to reorder the
328 # indeterminates, so that it's easier to evaluate the
329 # characteristic polynomial at x's coordinates and get back
330 # something in terms of t, which is what we want.
331 R = a[0].parent()
332 S = PolynomialRing(self.base_ring(),'t')
333 t = S.gen(0)
334 S = PolynomialRing(S, R.variable_names())
335 t = S(t)
336
337 # Note: all entries past the rth should be zero. The
338 # coefficient of the highest power (x^r) is 1, but it doesn't
339 # appear in the solution vector which contains coefficients
340 # for the other powers (to make them sum to x^r).
341 if (r < n):
342 a[r] = 1 # corresponds to x^r
343 else:
344 # When the rank is equal to the dimension, trying to
345 # assign a[r] goes out-of-bounds.
346 a.append(1) # corresponds to x^r
347
348 return sum( a[k]*(t**k) for k in range(len(a)) )
349
350
351 def inner_product(self, x, y):
352 """
353 The inner product associated with this Euclidean Jordan algebra.
354
355 Defaults to the trace inner product, but can be overridden by
356 subclasses if they are sure that the necessary properties are
357 satisfied.
358
359 SETUP::
360
361 sage: from mjo.eja.eja_algebra import random_eja
362
363 EXAMPLES:
364
365 The inner product must satisfy its axiom for this algebra to truly
366 be a Euclidean Jordan Algebra::
367
368 sage: set_random_seed()
369 sage: J = random_eja()
370 sage: x = J.random_element()
371 sage: y = J.random_element()
372 sage: z = J.random_element()
373 sage: (x*y).inner_product(z) == y.inner_product(x*z)
374 True
375
376 """
377 if (not x in self) or (not y in self):
378 raise TypeError("arguments must live in this algebra")
379 return x.trace_inner_product(y)
380
381
382 def natural_basis(self):
383 """
384 Return a more-natural representation of this algebra's basis.
385
386 Every finite-dimensional Euclidean Jordan Algebra is a direct
387 sum of five simple algebras, four of which comprise Hermitian
388 matrices. This method returns the original "natural" basis
389 for our underlying vector space. (Typically, the natural basis
390 is used to construct the multiplication table in the first place.)
391
392 Note that this will always return a matrix. The standard basis
393 in `R^n` will be returned as `n`-by-`1` column matrices.
394
395 SETUP::
396
397 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
398 ....: RealSymmetricEJA)
399
400 EXAMPLES::
401
402 sage: J = RealSymmetricEJA(2)
403 sage: J.basis()
404 Finite family {0: e0, 1: e1, 2: e2}
405 sage: J.natural_basis()
406 (
407 [1 0] [0 1] [0 0]
408 [0 0], [1 0], [0 1]
409 )
410
411 ::
412
413 sage: J = JordanSpinEJA(2)
414 sage: J.basis()
415 Finite family {0: e0, 1: e1}
416 sage: J.natural_basis()
417 (
418 [1] [0]
419 [0], [1]
420 )
421
422 """
423 if self._natural_basis is None:
424 return tuple( b.to_vector().column() for b in self.basis() )
425 else:
426 return self._natural_basis
427
428
429 @cached_method
430 def one(self):
431 """
432 Return the unit element of this algebra.
433
434 SETUP::
435
436 sage: from mjo.eja.eja_algebra import (RealCartesianProductEJA,
437 ....: random_eja)
438
439 EXAMPLES::
440
441 sage: J = RealCartesianProductEJA(5)
442 sage: J.one()
443 e0 + e1 + e2 + e3 + e4
444
445 TESTS::
446
447 The identity element acts like the identity::
448
449 sage: set_random_seed()
450 sage: J = random_eja()
451 sage: x = J.random_element()
452 sage: J.one()*x == x and x*J.one() == x
453 True
454
455 The matrix of the unit element's operator is the identity::
456
457 sage: set_random_seed()
458 sage: J = random_eja()
459 sage: actual = J.one().operator().matrix()
460 sage: expected = matrix.identity(J.base_ring(), J.dimension())
461 sage: actual == expected
462 True
463
464 """
465 # We can brute-force compute the matrices of the operators
466 # that correspond to the basis elements of this algebra.
467 # If some linear combination of those basis elements is the
468 # algebra identity, then the same linear combination of
469 # their matrices has to be the identity matrix.
470 #
471 # Of course, matrices aren't vectors in sage, so we have to
472 # appeal to the "long vectors" isometry.
473 oper_vecs = [ _mat2vec(g.operator().matrix()) for g in self.gens() ]
474
475 # Now we use basis linear algebra to find the coefficients,
476 # of the matrices-as-vectors-linear-combination, which should
477 # work for the original algebra basis too.
478 A = matrix.column(self.base_ring(), oper_vecs)
479
480 # We used the isometry on the left-hand side already, but we
481 # still need to do it for the right-hand side. Recall that we
482 # wanted something that summed to the identity matrix.
483 b = _mat2vec( matrix.identity(self.base_ring(), self.dimension()) )
484
485 # Now if there's an identity element in the algebra, this should work.
486 coeffs = A.solve_right(b)
487 return self.linear_combination(zip(self.gens(), coeffs))
488
489
490 def rank(self):
491 """
492 Return the rank of this EJA.
493
494 ALGORITHM:
495
496 The author knows of no algorithm to compute the rank of an EJA
497 where only the multiplication table is known. In lieu of one, we
498 require the rank to be specified when the algebra is created,
499 and simply pass along that number here.
500
501 SETUP::
502
503 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
504 ....: RealSymmetricEJA,
505 ....: ComplexHermitianEJA,
506 ....: QuaternionHermitianEJA,
507 ....: random_eja)
508
509 EXAMPLES:
510
511 The rank of the Jordan spin algebra is always two::
512
513 sage: JordanSpinEJA(2).rank()
514 2
515 sage: JordanSpinEJA(3).rank()
516 2
517 sage: JordanSpinEJA(4).rank()
518 2
519
520 The rank of the `n`-by-`n` Hermitian real, complex, or
521 quaternion matrices is `n`::
522
523 sage: RealSymmetricEJA(2).rank()
524 2
525 sage: ComplexHermitianEJA(2).rank()
526 2
527 sage: QuaternionHermitianEJA(2).rank()
528 2
529 sage: RealSymmetricEJA(5).rank()
530 5
531 sage: ComplexHermitianEJA(5).rank()
532 5
533 sage: QuaternionHermitianEJA(5).rank()
534 5
535
536 TESTS:
537
538 Ensure that every EJA that we know how to construct has a
539 positive integer rank::
540
541 sage: set_random_seed()
542 sage: r = random_eja().rank()
543 sage: r in ZZ and r > 0
544 True
545
546 """
547 return self._rank
548
549
550 def vector_space(self):
551 """
552 Return the vector space that underlies this algebra.
553
554 SETUP::
555
556 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
557
558 EXAMPLES::
559
560 sage: J = RealSymmetricEJA(2)
561 sage: J.vector_space()
562 Vector space of dimension 3 over Rational Field
563
564 """
565 return self.zero().to_vector().parent().ambient_vector_space()
566
567
568 Element = FiniteDimensionalEuclideanJordanAlgebraElement
569
570
571 class RealCartesianProductEJA(FiniteDimensionalEuclideanJordanAlgebra):
572 """
573 Return the Euclidean Jordan Algebra corresponding to the set
574 `R^n` under the Hadamard product.
575
576 Note: this is nothing more than the Cartesian product of ``n``
577 copies of the spin algebra. Once Cartesian product algebras
578 are implemented, this can go.
579
580 SETUP::
581
582 sage: from mjo.eja.eja_algebra import RealCartesianProductEJA
583
584 EXAMPLES:
585
586 This multiplication table can be verified by hand::
587
588 sage: J = RealCartesianProductEJA(3)
589 sage: e0,e1,e2 = J.gens()
590 sage: e0*e0
591 e0
592 sage: e0*e1
593 0
594 sage: e0*e2
595 0
596 sage: e1*e1
597 e1
598 sage: e1*e2
599 0
600 sage: e2*e2
601 e2
602
603 """
604 def __init__(self, n, field=QQ):
605 V = VectorSpace(field, n)
606 mult_table = [ [ V.basis()[i]*(i == j) for i in range(n) ]
607 for j in range(n) ]
608
609 fdeja = super(RealCartesianProductEJA, self)
610 return fdeja.__init__(field, mult_table, rank=n)
611
612 def inner_product(self, x, y):
613 return _usual_ip(x,y)
614
615
616 def random_eja():
617 """
618 Return a "random" finite-dimensional Euclidean Jordan Algebra.
619
620 ALGORITHM:
621
622 For now, we choose a random natural number ``n`` (greater than zero)
623 and then give you back one of the following:
624
625 * The cartesian product of the rational numbers ``n`` times; this is
626 ``QQ^n`` with the Hadamard product.
627
628 * The Jordan spin algebra on ``QQ^n``.
629
630 * The ``n``-by-``n`` rational symmetric matrices with the symmetric
631 product.
632
633 * The ``n``-by-``n`` complex-rational Hermitian matrices embedded
634 in the space of ``2n``-by-``2n`` real symmetric matrices.
635
636 * The ``n``-by-``n`` quaternion-rational Hermitian matrices embedded
637 in the space of ``4n``-by-``4n`` real symmetric matrices.
638
639 Later this might be extended to return Cartesian products of the
640 EJAs above.
641
642 SETUP::
643
644 sage: from mjo.eja.eja_algebra import random_eja
645
646 TESTS::
647
648 sage: random_eja()
649 Euclidean Jordan algebra of degree...
650
651 """
652
653 # The max_n component lets us choose different upper bounds on the
654 # value "n" that gets passed to the constructor. This is needed
655 # because e.g. R^{10} is reasonable to test, while the Hermitian
656 # 10-by-10 quaternion matrices are not.
657 (constructor, max_n) = choice([(RealCartesianProductEJA, 6),
658 (JordanSpinEJA, 6),
659 (RealSymmetricEJA, 5),
660 (ComplexHermitianEJA, 4),
661 (QuaternionHermitianEJA, 3)])
662 n = ZZ.random_element(1, max_n)
663 return constructor(n, field=QQ)
664
665
666
667 def _real_symmetric_basis(n, field=QQ):
668 """
669 Return a basis for the space of real symmetric n-by-n matrices.
670 """
671 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
672 # coordinates.
673 S = []
674 for i in xrange(n):
675 for j in xrange(i+1):
676 Eij = matrix(field, n, lambda k,l: k==i and l==j)
677 if i == j:
678 Sij = Eij
679 else:
680 # Beware, orthogonal but not normalized!
681 Sij = Eij + Eij.transpose()
682 S.append(Sij)
683 return tuple(S)
684
685
686 def _complex_hermitian_basis(n, field=QQ):
687 """
688 Returns a basis for the space of complex Hermitian n-by-n matrices.
689
690 SETUP::
691
692 sage: from mjo.eja.eja_algebra import _complex_hermitian_basis
693
694 TESTS::
695
696 sage: set_random_seed()
697 sage: n = ZZ.random_element(1,5)
698 sage: all( M.is_symmetric() for M in _complex_hermitian_basis(n) )
699 True
700
701 """
702 F = QuadraticField(-1, 'I')
703 I = F.gen()
704
705 # This is like the symmetric case, but we need to be careful:
706 #
707 # * We want conjugate-symmetry, not just symmetry.
708 # * The diagonal will (as a result) be real.
709 #
710 S = []
711 for i in xrange(n):
712 for j in xrange(i+1):
713 Eij = matrix(field, n, lambda k,l: k==i and l==j)
714 if i == j:
715 Sij = _embed_complex_matrix(Eij)
716 S.append(Sij)
717 else:
718 # Beware, orthogonal but not normalized! The second one
719 # has a minus because it's conjugated.
720 Sij_real = _embed_complex_matrix(Eij + Eij.transpose())
721 S.append(Sij_real)
722 Sij_imag = _embed_complex_matrix(I*Eij - I*Eij.transpose())
723 S.append(Sij_imag)
724 return tuple(S)
725
726
727 def _quaternion_hermitian_basis(n, field=QQ):
728 """
729 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
730
731 SETUP::
732
733 sage: from mjo.eja.eja_algebra import _quaternion_hermitian_basis
734
735 TESTS::
736
737 sage: set_random_seed()
738 sage: n = ZZ.random_element(1,5)
739 sage: all( M.is_symmetric() for M in _quaternion_hermitian_basis(n) )
740 True
741
742 """
743 Q = QuaternionAlgebra(QQ,-1,-1)
744 I,J,K = Q.gens()
745
746 # This is like the symmetric case, but we need to be careful:
747 #
748 # * We want conjugate-symmetry, not just symmetry.
749 # * The diagonal will (as a result) be real.
750 #
751 S = []
752 for i in xrange(n):
753 for j in xrange(i+1):
754 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
755 if i == j:
756 Sij = _embed_quaternion_matrix(Eij)
757 S.append(Sij)
758 else:
759 # Beware, orthogonal but not normalized! The second,
760 # third, and fourth ones have a minus because they're
761 # conjugated.
762 Sij_real = _embed_quaternion_matrix(Eij + Eij.transpose())
763 S.append(Sij_real)
764 Sij_I = _embed_quaternion_matrix(I*Eij - I*Eij.transpose())
765 S.append(Sij_I)
766 Sij_J = _embed_quaternion_matrix(J*Eij - J*Eij.transpose())
767 S.append(Sij_J)
768 Sij_K = _embed_quaternion_matrix(K*Eij - K*Eij.transpose())
769 S.append(Sij_K)
770 return tuple(S)
771
772
773
774 def _multiplication_table_from_matrix_basis(basis):
775 """
776 At least three of the five simple Euclidean Jordan algebras have the
777 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
778 multiplication on the right is matrix multiplication. Given a basis
779 for the underlying matrix space, this function returns a
780 multiplication table (obtained by looping through the basis
781 elements) for an algebra of those matrices.
782 """
783 # In S^2, for example, we nominally have four coordinates even
784 # though the space is of dimension three only. The vector space V
785 # is supposed to hold the entire long vector, and the subspace W
786 # of V will be spanned by the vectors that arise from symmetric
787 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
788 field = basis[0].base_ring()
789 dimension = basis[0].nrows()
790
791 V = VectorSpace(field, dimension**2)
792 W = V.span_of_basis( _mat2vec(s) for s in basis )
793 n = len(basis)
794 mult_table = [[W.zero() for i in range(n)] for j in range(n)]
795 for i in range(n):
796 for j in range(n):
797 mat_entry = (basis[i]*basis[j] + basis[j]*basis[i])/2
798 mult_table[i][j] = W.coordinate_vector(_mat2vec(mat_entry))
799
800 return mult_table
801
802
803 def _embed_complex_matrix(M):
804 """
805 Embed the n-by-n complex matrix ``M`` into the space of real
806 matrices of size 2n-by-2n via the map the sends each entry `z = a +
807 bi` to the block matrix ``[[a,b],[-b,a]]``.
808
809 SETUP::
810
811 sage: from mjo.eja.eja_algebra import _embed_complex_matrix
812
813 EXAMPLES::
814
815 sage: F = QuadraticField(-1,'i')
816 sage: x1 = F(4 - 2*i)
817 sage: x2 = F(1 + 2*i)
818 sage: x3 = F(-i)
819 sage: x4 = F(6)
820 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
821 sage: _embed_complex_matrix(M)
822 [ 4 -2| 1 2]
823 [ 2 4|-2 1]
824 [-----+-----]
825 [ 0 -1| 6 0]
826 [ 1 0| 0 6]
827
828 TESTS:
829
830 Embedding is a homomorphism (isomorphism, in fact)::
831
832 sage: set_random_seed()
833 sage: n = ZZ.random_element(5)
834 sage: F = QuadraticField(-1, 'i')
835 sage: X = random_matrix(F, n)
836 sage: Y = random_matrix(F, n)
837 sage: actual = _embed_complex_matrix(X) * _embed_complex_matrix(Y)
838 sage: expected = _embed_complex_matrix(X*Y)
839 sage: actual == expected
840 True
841
842 """
843 n = M.nrows()
844 if M.ncols() != n:
845 raise ValueError("the matrix 'M' must be square")
846 field = M.base_ring()
847 blocks = []
848 for z in M.list():
849 a = z.real()
850 b = z.imag()
851 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
852
853 # We can drop the imaginaries here.
854 return matrix.block(field.base_ring(), n, blocks)
855
856
857 def _unembed_complex_matrix(M):
858 """
859 The inverse of _embed_complex_matrix().
860
861 SETUP::
862
863 sage: from mjo.eja.eja_algebra import (_embed_complex_matrix,
864 ....: _unembed_complex_matrix)
865
866 EXAMPLES::
867
868 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
869 ....: [-2, 1, -4, 3],
870 ....: [ 9, 10, 11, 12],
871 ....: [-10, 9, -12, 11] ])
872 sage: _unembed_complex_matrix(A)
873 [ 2*i + 1 4*i + 3]
874 [ 10*i + 9 12*i + 11]
875
876 TESTS:
877
878 Unembedding is the inverse of embedding::
879
880 sage: set_random_seed()
881 sage: F = QuadraticField(-1, 'i')
882 sage: M = random_matrix(F, 3)
883 sage: _unembed_complex_matrix(_embed_complex_matrix(M)) == M
884 True
885
886 """
887 n = ZZ(M.nrows())
888 if M.ncols() != n:
889 raise ValueError("the matrix 'M' must be square")
890 if not n.mod(2).is_zero():
891 raise ValueError("the matrix 'M' must be a complex embedding")
892
893 F = QuadraticField(-1, 'i')
894 i = F.gen()
895
896 # Go top-left to bottom-right (reading order), converting every
897 # 2-by-2 block we see to a single complex element.
898 elements = []
899 for k in xrange(n/2):
900 for j in xrange(n/2):
901 submat = M[2*k:2*k+2,2*j:2*j+2]
902 if submat[0,0] != submat[1,1]:
903 raise ValueError('bad on-diagonal submatrix')
904 if submat[0,1] != -submat[1,0]:
905 raise ValueError('bad off-diagonal submatrix')
906 z = submat[0,0] + submat[0,1]*i
907 elements.append(z)
908
909 return matrix(F, n/2, elements)
910
911
912 def _embed_quaternion_matrix(M):
913 """
914 Embed the n-by-n quaternion matrix ``M`` into the space of real
915 matrices of size 4n-by-4n by first sending each quaternion entry
916 `z = a + bi + cj + dk` to the block-complex matrix
917 ``[[a + bi, c+di],[-c + di, a-bi]]`, and then embedding those into
918 a real matrix.
919
920 SETUP::
921
922 sage: from mjo.eja.eja_algebra import _embed_quaternion_matrix
923
924 EXAMPLES::
925
926 sage: Q = QuaternionAlgebra(QQ,-1,-1)
927 sage: i,j,k = Q.gens()
928 sage: x = 1 + 2*i + 3*j + 4*k
929 sage: M = matrix(Q, 1, [[x]])
930 sage: _embed_quaternion_matrix(M)
931 [ 1 2 3 4]
932 [-2 1 -4 3]
933 [-3 4 1 -2]
934 [-4 -3 2 1]
935
936 Embedding is a homomorphism (isomorphism, in fact)::
937
938 sage: set_random_seed()
939 sage: n = ZZ.random_element(5)
940 sage: Q = QuaternionAlgebra(QQ,-1,-1)
941 sage: X = random_matrix(Q, n)
942 sage: Y = random_matrix(Q, n)
943 sage: actual = _embed_quaternion_matrix(X)*_embed_quaternion_matrix(Y)
944 sage: expected = _embed_quaternion_matrix(X*Y)
945 sage: actual == expected
946 True
947
948 """
949 quaternions = M.base_ring()
950 n = M.nrows()
951 if M.ncols() != n:
952 raise ValueError("the matrix 'M' must be square")
953
954 F = QuadraticField(-1, 'i')
955 i = F.gen()
956
957 blocks = []
958 for z in M.list():
959 t = z.coefficient_tuple()
960 a = t[0]
961 b = t[1]
962 c = t[2]
963 d = t[3]
964 cplx_matrix = matrix(F, 2, [[ a + b*i, c + d*i],
965 [-c + d*i, a - b*i]])
966 blocks.append(_embed_complex_matrix(cplx_matrix))
967
968 # We should have real entries by now, so use the realest field
969 # we've got for the return value.
970 return matrix.block(quaternions.base_ring(), n, blocks)
971
972
973 def _unembed_quaternion_matrix(M):
974 """
975 The inverse of _embed_quaternion_matrix().
976
977 SETUP::
978
979 sage: from mjo.eja.eja_algebra import (_embed_quaternion_matrix,
980 ....: _unembed_quaternion_matrix)
981
982 EXAMPLES::
983
984 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
985 ....: [-2, 1, -4, 3],
986 ....: [-3, 4, 1, -2],
987 ....: [-4, -3, 2, 1]])
988 sage: _unembed_quaternion_matrix(M)
989 [1 + 2*i + 3*j + 4*k]
990
991 TESTS:
992
993 Unembedding is the inverse of embedding::
994
995 sage: set_random_seed()
996 sage: Q = QuaternionAlgebra(QQ, -1, -1)
997 sage: M = random_matrix(Q, 3)
998 sage: _unembed_quaternion_matrix(_embed_quaternion_matrix(M)) == M
999 True
1000
1001 """
1002 n = ZZ(M.nrows())
1003 if M.ncols() != n:
1004 raise ValueError("the matrix 'M' must be square")
1005 if not n.mod(4).is_zero():
1006 raise ValueError("the matrix 'M' must be a complex embedding")
1007
1008 Q = QuaternionAlgebra(QQ,-1,-1)
1009 i,j,k = Q.gens()
1010
1011 # Go top-left to bottom-right (reading order), converting every
1012 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
1013 # quaternion block.
1014 elements = []
1015 for l in xrange(n/4):
1016 for m in xrange(n/4):
1017 submat = _unembed_complex_matrix(M[4*l:4*l+4,4*m:4*m+4])
1018 if submat[0,0] != submat[1,1].conjugate():
1019 raise ValueError('bad on-diagonal submatrix')
1020 if submat[0,1] != -submat[1,0].conjugate():
1021 raise ValueError('bad off-diagonal submatrix')
1022 z = submat[0,0].real() + submat[0,0].imag()*i
1023 z += submat[0,1].real()*j + submat[0,1].imag()*k
1024 elements.append(z)
1025
1026 return matrix(Q, n/4, elements)
1027
1028
1029 # The usual inner product on R^n.
1030 def _usual_ip(x,y):
1031 return x.to_vector().inner_product(y.to_vector())
1032
1033 # The inner product used for the real symmetric simple EJA.
1034 # We keep it as a separate function because e.g. the complex
1035 # algebra uses the same inner product, except divided by 2.
1036 def _matrix_ip(X,Y):
1037 X_mat = X.natural_representation()
1038 Y_mat = Y.natural_representation()
1039 return (X_mat*Y_mat).trace()
1040
1041
1042 class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra):
1043 """
1044 The rank-n simple EJA consisting of real symmetric n-by-n
1045 matrices, the usual symmetric Jordan product, and the trace inner
1046 product. It has dimension `(n^2 + n)/2` over the reals.
1047
1048 SETUP::
1049
1050 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1051
1052 EXAMPLES::
1053
1054 sage: J = RealSymmetricEJA(2)
1055 sage: e0, e1, e2 = J.gens()
1056 sage: e0*e0
1057 e0
1058 sage: e1*e1
1059 e0 + e2
1060 sage: e2*e2
1061 e2
1062
1063 TESTS:
1064
1065 The dimension of this algebra is `(n^2 + n) / 2`::
1066
1067 sage: set_random_seed()
1068 sage: n = ZZ.random_element(1,5)
1069 sage: J = RealSymmetricEJA(n)
1070 sage: J.dimension() == (n^2 + n)/2
1071 True
1072
1073 The Jordan multiplication is what we think it is::
1074
1075 sage: set_random_seed()
1076 sage: n = ZZ.random_element(1,5)
1077 sage: J = RealSymmetricEJA(n)
1078 sage: x = J.random_element()
1079 sage: y = J.random_element()
1080 sage: actual = (x*y).natural_representation()
1081 sage: X = x.natural_representation()
1082 sage: Y = y.natural_representation()
1083 sage: expected = (X*Y + Y*X)/2
1084 sage: actual == expected
1085 True
1086 sage: J(expected) == x*y
1087 True
1088
1089 """
1090 def __init__(self, n, field=QQ):
1091 S = _real_symmetric_basis(n, field=field)
1092 Qs = _multiplication_table_from_matrix_basis(S)
1093
1094 fdeja = super(RealSymmetricEJA, self)
1095 return fdeja.__init__(field,
1096 Qs,
1097 rank=n,
1098 natural_basis=S)
1099
1100 def inner_product(self, x, y):
1101 return _matrix_ip(x,y)
1102
1103
1104 class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
1105 """
1106 The rank-n simple EJA consisting of complex Hermitian n-by-n
1107 matrices over the real numbers, the usual symmetric Jordan product,
1108 and the real-part-of-trace inner product. It has dimension `n^2` over
1109 the reals.
1110
1111 SETUP::
1112
1113 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1114
1115 TESTS:
1116
1117 The dimension of this algebra is `n^2`::
1118
1119 sage: set_random_seed()
1120 sage: n = ZZ.random_element(1,5)
1121 sage: J = ComplexHermitianEJA(n)
1122 sage: J.dimension() == n^2
1123 True
1124
1125 The Jordan multiplication is what we think it is::
1126
1127 sage: set_random_seed()
1128 sage: n = ZZ.random_element(1,5)
1129 sage: J = ComplexHermitianEJA(n)
1130 sage: x = J.random_element()
1131 sage: y = J.random_element()
1132 sage: actual = (x*y).natural_representation()
1133 sage: X = x.natural_representation()
1134 sage: Y = y.natural_representation()
1135 sage: expected = (X*Y + Y*X)/2
1136 sage: actual == expected
1137 True
1138 sage: J(expected) == x*y
1139 True
1140
1141 """
1142 def __init__(self, n, field=QQ):
1143 S = _complex_hermitian_basis(n)
1144 Qs = _multiplication_table_from_matrix_basis(S)
1145
1146 fdeja = super(ComplexHermitianEJA, self)
1147 return fdeja.__init__(field,
1148 Qs,
1149 rank=n,
1150 natural_basis=S)
1151
1152
1153 def inner_product(self, x, y):
1154 # Since a+bi on the diagonal is represented as
1155 #
1156 # a + bi = [ a b ]
1157 # [ -b a ],
1158 #
1159 # we'll double-count the "a" entries if we take the trace of
1160 # the embedding.
1161 return _matrix_ip(x,y)/2
1162
1163
1164 class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
1165 """
1166 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
1167 matrices, the usual symmetric Jordan product, and the
1168 real-part-of-trace inner product. It has dimension `2n^2 - n` over
1169 the reals.
1170
1171 SETUP::
1172
1173 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1174
1175 TESTS:
1176
1177 The dimension of this algebra is `n^2`::
1178
1179 sage: set_random_seed()
1180 sage: n = ZZ.random_element(1,5)
1181 sage: J = QuaternionHermitianEJA(n)
1182 sage: J.dimension() == 2*(n^2) - n
1183 True
1184
1185 The Jordan multiplication is what we think it is::
1186
1187 sage: set_random_seed()
1188 sage: n = ZZ.random_element(1,5)
1189 sage: J = QuaternionHermitianEJA(n)
1190 sage: x = J.random_element()
1191 sage: y = J.random_element()
1192 sage: actual = (x*y).natural_representation()
1193 sage: X = x.natural_representation()
1194 sage: Y = y.natural_representation()
1195 sage: expected = (X*Y + Y*X)/2
1196 sage: actual == expected
1197 True
1198 sage: J(expected) == x*y
1199 True
1200
1201 """
1202 def __init__(self, n, field=QQ):
1203 S = _quaternion_hermitian_basis(n)
1204 Qs = _multiplication_table_from_matrix_basis(S)
1205
1206 fdeja = super(QuaternionHermitianEJA, self)
1207 return fdeja.__init__(field,
1208 Qs,
1209 rank=n,
1210 natural_basis=S)
1211
1212 def inner_product(self, x, y):
1213 # Since a+bi+cj+dk on the diagonal is represented as
1214 #
1215 # a + bi +cj + dk = [ a b c d]
1216 # [ -b a -d c]
1217 # [ -c d a -b]
1218 # [ -d -c b a],
1219 #
1220 # we'll quadruple-count the "a" entries if we take the trace of
1221 # the embedding.
1222 return _matrix_ip(x,y)/4
1223
1224
1225 class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra):
1226 """
1227 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
1228 with the usual inner product and jordan product ``x*y =
1229 (<x_bar,y_bar>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
1230 the reals.
1231
1232 SETUP::
1233
1234 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1235
1236 EXAMPLES:
1237
1238 This multiplication table can be verified by hand::
1239
1240 sage: J = JordanSpinEJA(4)
1241 sage: e0,e1,e2,e3 = J.gens()
1242 sage: e0*e0
1243 e0
1244 sage: e0*e1
1245 e1
1246 sage: e0*e2
1247 e2
1248 sage: e0*e3
1249 e3
1250 sage: e1*e2
1251 0
1252 sage: e1*e3
1253 0
1254 sage: e2*e3
1255 0
1256
1257 """
1258 def __init__(self, n, field=QQ):
1259 V = VectorSpace(field, n)
1260 mult_table = [[V.zero() for i in range(n)] for j in range(n)]
1261 for i in range(n):
1262 for j in range(n):
1263 x = V.basis()[i]
1264 y = V.basis()[j]
1265 x0 = x[0]
1266 xbar = x[1:]
1267 y0 = y[0]
1268 ybar = y[1:]
1269 # z = x*y
1270 z0 = x.inner_product(y)
1271 zbar = y0*xbar + x0*ybar
1272 z = V([z0] + zbar.list())
1273 mult_table[i][j] = z
1274
1275 # The rank of the spin algebra is two, unless we're in a
1276 # one-dimensional ambient space (because the rank is bounded by
1277 # the ambient dimension).
1278 fdeja = super(JordanSpinEJA, self)
1279 return fdeja.__init__(field, mult_table, rank=min(n,2))
1280
1281 def inner_product(self, x, y):
1282 return _usual_ip(x,y)