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