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