]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: use izip() instead of zip() in a few places.
[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 itertools import izip, repeat
9
10 from sage.algebras.quatalg.quaternion_algebra import QuaternionAlgebra
11 from sage.categories.magmatic_algebras import MagmaticAlgebras
12 from sage.combinat.free_module import CombinatorialFreeModule
13 from sage.matrix.constructor import matrix
14 from sage.matrix.matrix_space import MatrixSpace
15 from sage.misc.cachefunc import cached_method
16 from sage.misc.prandom import choice
17 from sage.misc.table import table
18 from sage.modules.free_module import FreeModule, VectorSpace
19 from sage.rings.integer_ring import ZZ
20 from sage.rings.number_field.number_field import NumberField, QuadraticField
21 from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing
22 from sage.rings.rational_field import QQ
23 from sage.rings.real_lazy import CLF, RLF
24 from sage.structure.element import is_Matrix
25
26 from mjo.eja.eja_element import FiniteDimensionalEuclideanJordanAlgebraElement
27 from mjo.eja.eja_utils import _mat2vec
28
29 class FiniteDimensionalEuclideanJordanAlgebra(CombinatorialFreeModule):
30 # This is an ugly hack needed to prevent the category framework
31 # from implementing a coercion from our base ring (e.g. the
32 # rationals) into the algebra. First of all -- such a coercion is
33 # nonsense to begin with. But more importantly, it tries to do so
34 # in the category of rings, and since our algebras aren't
35 # associative they generally won't be rings.
36 _no_generic_basering_coercion = True
37
38 def __init__(self,
39 field,
40 mult_table,
41 rank,
42 prefix='e',
43 category=None,
44 natural_basis=None):
45 """
46 SETUP::
47
48 sage: from mjo.eja.eja_algebra import random_eja
49
50 EXAMPLES:
51
52 By definition, Jordan multiplication commutes::
53
54 sage: set_random_seed()
55 sage: J = random_eja()
56 sage: x,y = J.random_elements(2)
57 sage: x*y == y*x
58 True
59
60 """
61 self._rank = rank
62 self._natural_basis = natural_basis
63
64 if category is None:
65 category = MagmaticAlgebras(field).FiniteDimensional()
66 category = category.WithBasis().Unital()
67
68 fda = super(FiniteDimensionalEuclideanJordanAlgebra, self)
69 fda.__init__(field,
70 range(len(mult_table)),
71 prefix=prefix,
72 category=category)
73 self.print_options(bracket='')
74
75 # The multiplication table we're given is necessarily in terms
76 # of vectors, because we don't have an algebra yet for
77 # anything to be an element of. However, it's faster in the
78 # long run to have the multiplication table be in terms of
79 # algebra elements. We do this after calling the superclass
80 # constructor so that from_vector() knows what to do.
81 self._multiplication_table = [ map(lambda x: self.from_vector(x), ls)
82 for ls in mult_table ]
83
84
85 def _element_constructor_(self, elt):
86 """
87 Construct an element of this algebra from its natural
88 representation.
89
90 This gets called only after the parent element _call_ method
91 fails to find a coercion for the argument.
92
93 SETUP::
94
95 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
96 ....: RealCartesianProductEJA,
97 ....: RealSymmetricEJA)
98
99 EXAMPLES:
100
101 The identity in `S^n` is converted to the identity in the EJA::
102
103 sage: J = RealSymmetricEJA(3)
104 sage: I = matrix.identity(QQ,3)
105 sage: J(I) == J.one()
106 True
107
108 This skew-symmetric matrix can't be represented in the EJA::
109
110 sage: J = RealSymmetricEJA(3)
111 sage: A = matrix(QQ,3, lambda i,j: i-j)
112 sage: J(A)
113 Traceback (most recent call last):
114 ...
115 ArithmeticError: vector is not in free module
116
117 TESTS:
118
119 Ensure that we can convert any element of the two non-matrix
120 simple algebras (whose natural representations are their usual
121 vector representations) back and forth faithfully::
122
123 sage: set_random_seed()
124 sage: J = RealCartesianProductEJA.random_instance()
125 sage: x = J.random_element()
126 sage: J(x.to_vector().column()) == x
127 True
128 sage: J = JordanSpinEJA.random_instance()
129 sage: x = J.random_element()
130 sage: J(x.to_vector().column()) == x
131 True
132
133 """
134 if elt == 0:
135 # The superclass implementation of random_element()
136 # needs to be able to coerce "0" into the algebra.
137 return self.zero()
138
139 natural_basis = self.natural_basis()
140 basis_space = natural_basis[0].matrix_space()
141 if elt not in basis_space:
142 raise ValueError("not a naturally-represented algebra element")
143
144 # Thanks for nothing! Matrix spaces aren't vector spaces in
145 # Sage, so we have to figure out its natural-basis coordinates
146 # ourselves. We use the basis space's ring instead of the
147 # element's ring because the basis space might be an algebraic
148 # closure whereas the base ring of the 3-by-3 identity matrix
149 # could be QQ instead of QQbar.
150 V = VectorSpace(basis_space.base_ring(), elt.nrows()*elt.ncols())
151 W = V.span_of_basis( _mat2vec(s) for s in natural_basis )
152 coords = W.coordinate_vector(_mat2vec(elt))
153 return self.from_vector(coords)
154
155
156 @staticmethod
157 def _max_test_case_size():
158 """
159 Return an integer "size" that is an upper bound on the size of
160 this algebra when it is used in a random test
161 case. Unfortunately, the term "size" is quite vague -- when
162 dealing with `R^n` under either the Hadamard or Jordan spin
163 product, the "size" refers to the dimension `n`. When dealing
164 with a matrix algebra (real symmetric or complex/quaternion
165 Hermitian), it refers to the size of the matrix, which is
166 far less than the dimension of the underlying vector space.
167
168 We default to five in this class, which is safe in `R^n`. The
169 matrix algebra subclasses (or any class where the "size" is
170 interpreted to be far less than the dimension) should override
171 with a smaller number.
172 """
173 return 5
174
175
176 def _repr_(self):
177 """
178 Return a string representation of ``self``.
179
180 SETUP::
181
182 sage: from mjo.eja.eja_algebra import JordanSpinEJA
183
184 TESTS:
185
186 Ensure that it says what we think it says::
187
188 sage: JordanSpinEJA(2, field=QQ)
189 Euclidean Jordan algebra of dimension 2 over Rational Field
190 sage: JordanSpinEJA(3, field=RDF)
191 Euclidean Jordan algebra of dimension 3 over Real Double Field
192
193 """
194 fmt = "Euclidean Jordan algebra of dimension {} over {}"
195 return fmt.format(self.dimension(), self.base_ring())
196
197 def product_on_basis(self, i, j):
198 return self._multiplication_table[i][j]
199
200 def _a_regular_element(self):
201 """
202 Guess a regular element. Needed to compute the basis for our
203 characteristic polynomial coefficients.
204
205 SETUP::
206
207 sage: from mjo.eja.eja_algebra import random_eja
208
209 TESTS:
210
211 Ensure that this hacky method succeeds for every algebra that we
212 know how to construct::
213
214 sage: set_random_seed()
215 sage: J = random_eja()
216 sage: J._a_regular_element().is_regular()
217 True
218
219 """
220 gs = self.gens()
221 z = self.sum( (i+1)*gs[i] for i in range(len(gs)) )
222 if not z.is_regular():
223 raise ValueError("don't know a regular element")
224 return z
225
226
227 @cached_method
228 def _charpoly_basis_space(self):
229 """
230 Return the vector space spanned by the basis used in our
231 characteristic polynomial coefficients. This is used not only to
232 compute those coefficients, but also any time we need to
233 evaluate the coefficients (like when we compute the trace or
234 determinant).
235 """
236 z = self._a_regular_element()
237 # Don't use the parent vector space directly here in case this
238 # happens to be a subalgebra. In that case, we would be e.g.
239 # two-dimensional but span_of_basis() would expect three
240 # coordinates.
241 V = VectorSpace(self.base_ring(), self.vector_space().dimension())
242 basis = [ (z**k).to_vector() for k in range(self.rank()) ]
243 V1 = V.span_of_basis( basis )
244 b = (V1.basis() + V1.complement().basis())
245 return V.span_of_basis(b)
246
247
248
249 @cached_method
250 def _charpoly_coeff(self, i):
251 """
252 Return the coefficient polynomial "a_{i}" of this algebra's
253 general characteristic polynomial.
254
255 Having this be a separate cached method lets us compute and
256 store the trace/determinant (a_{r-1} and a_{0} respectively)
257 separate from the entire characteristic polynomial.
258 """
259 (A_of_x, x, xr, detA) = self._charpoly_matrix_system()
260 R = A_of_x.base_ring()
261 if i >= self.rank():
262 # Guaranteed by theory
263 return R.zero()
264
265 # Danger: the in-place modification is done for performance
266 # reasons (reconstructing a matrix with huge polynomial
267 # entries is slow), but I don't know how cached_method works,
268 # so it's highly possible that we're modifying some global
269 # list variable by reference, here. In other words, you
270 # probably shouldn't call this method twice on the same
271 # algebra, at the same time, in two threads
272 Ai_orig = A_of_x.column(i)
273 A_of_x.set_column(i,xr)
274 numerator = A_of_x.det()
275 A_of_x.set_column(i,Ai_orig)
276
277 # We're relying on the theory here to ensure that each a_i is
278 # indeed back in R, and the added negative signs are to make
279 # the whole charpoly expression sum to zero.
280 return R(-numerator/detA)
281
282
283 @cached_method
284 def _charpoly_matrix_system(self):
285 """
286 Compute the matrix whose entries A_ij are polynomials in
287 X1,...,XN, the vector ``x`` of variables X1,...,XN, the vector
288 corresponding to `x^r` and the determinent of the matrix A =
289 [A_ij]. In other words, all of the fixed (cachable) data needed
290 to compute the coefficients of the characteristic polynomial.
291 """
292 r = self.rank()
293 n = self.dimension()
294
295 # Turn my vector space into a module so that "vectors" can
296 # have multivatiate polynomial entries.
297 names = tuple('X' + str(i) for i in range(1,n+1))
298 R = PolynomialRing(self.base_ring(), names)
299
300 # Using change_ring() on the parent's vector space doesn't work
301 # here because, in a subalgebra, that vector space has a basis
302 # and change_ring() tries to bring the basis along with it. And
303 # that doesn't work unless the new ring is a PID, which it usually
304 # won't be.
305 V = FreeModule(R,n)
306
307 # Now let x = (X1,X2,...,Xn) be the vector whose entries are
308 # indeterminates...
309 x = V(names)
310
311 # And figure out the "left multiplication by x" matrix in
312 # that setting.
313 lmbx_cols = []
314 monomial_matrices = [ self.monomial(i).operator().matrix()
315 for i in range(n) ] # don't recompute these!
316 for k in range(n):
317 ek = self.monomial(k).to_vector()
318 lmbx_cols.append(
319 sum( x[i]*(monomial_matrices[i]*ek)
320 for i in range(n) ) )
321 Lx = matrix.column(R, lmbx_cols)
322
323 # Now we can compute powers of x "symbolically"
324 x_powers = [self.one().to_vector(), x]
325 for d in range(2, r+1):
326 x_powers.append( Lx*(x_powers[-1]) )
327
328 idmat = matrix.identity(R, n)
329
330 W = self._charpoly_basis_space()
331 W = W.change_ring(R.fraction_field())
332
333 # Starting with the standard coordinates x = (X1,X2,...,Xn)
334 # and then converting the entries to W-coordinates allows us
335 # to pass in the standard coordinates to the charpoly and get
336 # back the right answer. Specifically, with x = (X1,X2,...,Xn),
337 # we have
338 #
339 # W.coordinates(x^2) eval'd at (standard z-coords)
340 # =
341 # W-coords of (z^2)
342 # =
343 # W-coords of (standard coords of x^2 eval'd at std-coords of z)
344 #
345 # We want the middle equivalent thing in our matrix, but use
346 # the first equivalent thing instead so that we can pass in
347 # standard coordinates.
348 x_powers = [ W.coordinate_vector(xp) for xp in x_powers ]
349 l2 = [idmat.column(k-1) for k in range(r+1, n+1)]
350 A_of_x = matrix.column(R, n, (x_powers[:r] + l2))
351 return (A_of_x, x, x_powers[r], A_of_x.det())
352
353
354 @cached_method
355 def characteristic_polynomial(self):
356 """
357 Return a characteristic polynomial that works for all elements
358 of this algebra.
359
360 The resulting polynomial has `n+1` variables, where `n` is the
361 dimension of this algebra. The first `n` variables correspond to
362 the coordinates of an algebra element: when evaluated at the
363 coordinates of an algebra element with respect to a certain
364 basis, the result is a univariate polynomial (in the one
365 remaining variable ``t``), namely the characteristic polynomial
366 of that element.
367
368 SETUP::
369
370 sage: from mjo.eja.eja_algebra import JordanSpinEJA
371
372 EXAMPLES:
373
374 The characteristic polynomial in the spin algebra is given in
375 Alizadeh, Example 11.11::
376
377 sage: J = JordanSpinEJA(3)
378 sage: p = J.characteristic_polynomial(); p
379 X1^2 - X2^2 - X3^2 + (-2*t)*X1 + t^2
380 sage: xvec = J.one().to_vector()
381 sage: p(*xvec)
382 t^2 - 2*t + 1
383
384 """
385 r = self.rank()
386 n = self.dimension()
387
388 # The list of coefficient polynomials a_1, a_2, ..., a_n.
389 a = [ self._charpoly_coeff(i) for i in range(n) ]
390
391 # We go to a bit of trouble here to reorder the
392 # indeterminates, so that it's easier to evaluate the
393 # characteristic polynomial at x's coordinates and get back
394 # something in terms of t, which is what we want.
395 R = a[0].parent()
396 S = PolynomialRing(self.base_ring(),'t')
397 t = S.gen(0)
398 S = PolynomialRing(S, R.variable_names())
399 t = S(t)
400
401 # Note: all entries past the rth should be zero. The
402 # coefficient of the highest power (x^r) is 1, but it doesn't
403 # appear in the solution vector which contains coefficients
404 # for the other powers (to make them sum to x^r).
405 if (r < n):
406 a[r] = 1 # corresponds to x^r
407 else:
408 # When the rank is equal to the dimension, trying to
409 # assign a[r] goes out-of-bounds.
410 a.append(1) # corresponds to x^r
411
412 return sum( a[k]*(t**k) for k in xrange(len(a)) )
413
414
415 def inner_product(self, x, y):
416 """
417 The inner product associated with this Euclidean Jordan algebra.
418
419 Defaults to the trace inner product, but can be overridden by
420 subclasses if they are sure that the necessary properties are
421 satisfied.
422
423 SETUP::
424
425 sage: from mjo.eja.eja_algebra import random_eja
426
427 EXAMPLES:
428
429 Our inner product satisfies the Jordan axiom, which is also
430 referred to as "associativity" for a symmetric bilinear form::
431
432 sage: set_random_seed()
433 sage: J = random_eja()
434 sage: x,y,z = J.random_elements(3)
435 sage: (x*y).inner_product(z) == y.inner_product(x*z)
436 True
437
438 """
439 X = x.natural_representation()
440 Y = y.natural_representation()
441 return self.natural_inner_product(X,Y)
442
443
444 def is_trivial(self):
445 """
446 Return whether or not this algebra is trivial.
447
448 A trivial algebra contains only the zero element.
449
450 SETUP::
451
452 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
453
454 EXAMPLES::
455
456 sage: J = ComplexHermitianEJA(3)
457 sage: J.is_trivial()
458 False
459 sage: A = J.zero().subalgebra_generated_by()
460 sage: A.is_trivial()
461 True
462
463 """
464 return self.dimension() == 0
465
466
467 def multiplication_table(self):
468 """
469 Return a visual representation of this algebra's multiplication
470 table (on basis elements).
471
472 SETUP::
473
474 sage: from mjo.eja.eja_algebra import JordanSpinEJA
475
476 EXAMPLES::
477
478 sage: J = JordanSpinEJA(4)
479 sage: J.multiplication_table()
480 +----++----+----+----+----+
481 | * || e0 | e1 | e2 | e3 |
482 +====++====+====+====+====+
483 | e0 || e0 | e1 | e2 | e3 |
484 +----++----+----+----+----+
485 | e1 || e1 | e0 | 0 | 0 |
486 +----++----+----+----+----+
487 | e2 || e2 | 0 | e0 | 0 |
488 +----++----+----+----+----+
489 | e3 || e3 | 0 | 0 | e0 |
490 +----++----+----+----+----+
491
492 """
493 M = list(self._multiplication_table) # copy
494 for i in xrange(len(M)):
495 # M had better be "square"
496 M[i] = [self.monomial(i)] + M[i]
497 M = [["*"] + list(self.gens())] + M
498 return table(M, header_row=True, header_column=True, frame=True)
499
500
501 def natural_basis(self):
502 """
503 Return a more-natural representation of this algebra's basis.
504
505 Every finite-dimensional Euclidean Jordan Algebra is a direct
506 sum of five simple algebras, four of which comprise Hermitian
507 matrices. This method returns the original "natural" basis
508 for our underlying vector space. (Typically, the natural basis
509 is used to construct the multiplication table in the first place.)
510
511 Note that this will always return a matrix. The standard basis
512 in `R^n` will be returned as `n`-by-`1` column matrices.
513
514 SETUP::
515
516 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
517 ....: RealSymmetricEJA)
518
519 EXAMPLES::
520
521 sage: J = RealSymmetricEJA(2)
522 sage: J.basis()
523 Finite family {0: e0, 1: e1, 2: e2}
524 sage: J.natural_basis()
525 (
526 [1 0] [ 0 1/2*sqrt2] [0 0]
527 [0 0], [1/2*sqrt2 0], [0 1]
528 )
529
530 ::
531
532 sage: J = JordanSpinEJA(2)
533 sage: J.basis()
534 Finite family {0: e0, 1: e1}
535 sage: J.natural_basis()
536 (
537 [1] [0]
538 [0], [1]
539 )
540
541 """
542 if self._natural_basis is None:
543 M = self.natural_basis_space()
544 return tuple( M(b.to_vector()) for b in self.basis() )
545 else:
546 return self._natural_basis
547
548
549 def natural_basis_space(self):
550 """
551 Return the matrix space in which this algebra's natural basis
552 elements live.
553 """
554 if self._natural_basis is None or len(self._natural_basis) == 0:
555 return MatrixSpace(self.base_ring(), self.dimension(), 1)
556 else:
557 return self._natural_basis[0].matrix_space()
558
559
560 @staticmethod
561 def natural_inner_product(X,Y):
562 """
563 Compute the inner product of two naturally-represented elements.
564
565 For example in the real symmetric matrix EJA, this will compute
566 the trace inner-product of two n-by-n symmetric matrices. The
567 default should work for the real cartesian product EJA, the
568 Jordan spin EJA, and the real symmetric matrices. The others
569 will have to be overridden.
570 """
571 return (X.conjugate_transpose()*Y).trace()
572
573
574 @cached_method
575 def one(self):
576 """
577 Return the unit element of this algebra.
578
579 SETUP::
580
581 sage: from mjo.eja.eja_algebra import (RealCartesianProductEJA,
582 ....: random_eja)
583
584 EXAMPLES::
585
586 sage: J = RealCartesianProductEJA(5)
587 sage: J.one()
588 e0 + e1 + e2 + e3 + e4
589
590 TESTS:
591
592 The identity element acts like the identity::
593
594 sage: set_random_seed()
595 sage: J = random_eja()
596 sage: x = J.random_element()
597 sage: J.one()*x == x and x*J.one() == x
598 True
599
600 The matrix of the unit element's operator is the identity::
601
602 sage: set_random_seed()
603 sage: J = random_eja()
604 sage: actual = J.one().operator().matrix()
605 sage: expected = matrix.identity(J.base_ring(), J.dimension())
606 sage: actual == expected
607 True
608
609 """
610 # We can brute-force compute the matrices of the operators
611 # that correspond to the basis elements of this algebra.
612 # If some linear combination of those basis elements is the
613 # algebra identity, then the same linear combination of
614 # their matrices has to be the identity matrix.
615 #
616 # Of course, matrices aren't vectors in sage, so we have to
617 # appeal to the "long vectors" isometry.
618 oper_vecs = [ _mat2vec(g.operator().matrix()) for g in self.gens() ]
619
620 # Now we use basis linear algebra to find the coefficients,
621 # of the matrices-as-vectors-linear-combination, which should
622 # work for the original algebra basis too.
623 A = matrix.column(self.base_ring(), oper_vecs)
624
625 # We used the isometry on the left-hand side already, but we
626 # still need to do it for the right-hand side. Recall that we
627 # wanted something that summed to the identity matrix.
628 b = _mat2vec( matrix.identity(self.base_ring(), self.dimension()) )
629
630 # Now if there's an identity element in the algebra, this should work.
631 coeffs = A.solve_right(b)
632 return self.linear_combination(zip(self.gens(), coeffs))
633
634
635 def random_element(self):
636 # Temporary workaround for https://trac.sagemath.org/ticket/28327
637 if self.is_trivial():
638 return self.zero()
639 else:
640 s = super(FiniteDimensionalEuclideanJordanAlgebra, self)
641 return s.random_element()
642
643 def random_elements(self, count):
644 """
645 Return ``count`` random elements as a tuple.
646
647 SETUP::
648
649 sage: from mjo.eja.eja_algebra import JordanSpinEJA
650
651 EXAMPLES::
652
653 sage: J = JordanSpinEJA(3)
654 sage: x,y,z = J.random_elements(3)
655 sage: all( [ x in J, y in J, z in J ])
656 True
657 sage: len( J.random_elements(10) ) == 10
658 True
659
660 """
661 return tuple( self.random_element() for idx in xrange(count) )
662
663 @classmethod
664 def random_instance(cls, field=QQ, **kwargs):
665 """
666 Return a random instance of this type of algebra.
667
668 In subclasses for algebras that we know how to construct, this
669 is a shortcut for constructing test cases and examples.
670 """
671 if cls is FiniteDimensionalEuclideanJordanAlgebra:
672 # Red flag! But in theory we could do this I guess. The
673 # only finite-dimensional exceptional EJA is the
674 # octononions. So, we could just create an EJA from an
675 # associative matrix algebra (generated by a subset of
676 # elements) with the symmetric product. Or, we could punt
677 # to random_eja() here, override it in our subclasses, and
678 # not worry about it.
679 raise NotImplementedError
680
681 n = ZZ.random_element(cls._max_test_case_size()) + 1
682 return cls(n, field, **kwargs)
683
684
685 def rank(self):
686 """
687 Return the rank of this EJA.
688
689 ALGORITHM:
690
691 The author knows of no algorithm to compute the rank of an EJA
692 where only the multiplication table is known. In lieu of one, we
693 require the rank to be specified when the algebra is created,
694 and simply pass along that number here.
695
696 SETUP::
697
698 sage: from mjo.eja.eja_algebra import (JordanSpinEJA,
699 ....: RealSymmetricEJA,
700 ....: ComplexHermitianEJA,
701 ....: QuaternionHermitianEJA,
702 ....: random_eja)
703
704 EXAMPLES:
705
706 The rank of the Jordan spin algebra is always two::
707
708 sage: JordanSpinEJA(2).rank()
709 2
710 sage: JordanSpinEJA(3).rank()
711 2
712 sage: JordanSpinEJA(4).rank()
713 2
714
715 The rank of the `n`-by-`n` Hermitian real, complex, or
716 quaternion matrices is `n`::
717
718 sage: RealSymmetricEJA(4).rank()
719 4
720 sage: ComplexHermitianEJA(3).rank()
721 3
722 sage: QuaternionHermitianEJA(2).rank()
723 2
724
725 TESTS:
726
727 Ensure that every EJA that we know how to construct has a
728 positive integer rank::
729
730 sage: set_random_seed()
731 sage: r = random_eja().rank()
732 sage: r in ZZ and r > 0
733 True
734
735 """
736 return self._rank
737
738
739 def vector_space(self):
740 """
741 Return the vector space that underlies this algebra.
742
743 SETUP::
744
745 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
746
747 EXAMPLES::
748
749 sage: J = RealSymmetricEJA(2)
750 sage: J.vector_space()
751 Vector space of dimension 3 over...
752
753 """
754 return self.zero().to_vector().parent().ambient_vector_space()
755
756
757 Element = FiniteDimensionalEuclideanJordanAlgebraElement
758
759
760 class RealCartesianProductEJA(FiniteDimensionalEuclideanJordanAlgebra):
761 """
762 Return the Euclidean Jordan Algebra corresponding to the set
763 `R^n` under the Hadamard product.
764
765 Note: this is nothing more than the Cartesian product of ``n``
766 copies of the spin algebra. Once Cartesian product algebras
767 are implemented, this can go.
768
769 SETUP::
770
771 sage: from mjo.eja.eja_algebra import RealCartesianProductEJA
772
773 EXAMPLES:
774
775 This multiplication table can be verified by hand::
776
777 sage: J = RealCartesianProductEJA(3)
778 sage: e0,e1,e2 = J.gens()
779 sage: e0*e0
780 e0
781 sage: e0*e1
782 0
783 sage: e0*e2
784 0
785 sage: e1*e1
786 e1
787 sage: e1*e2
788 0
789 sage: e2*e2
790 e2
791
792 TESTS:
793
794 We can change the generator prefix::
795
796 sage: RealCartesianProductEJA(3, prefix='r').gens()
797 (r0, r1, r2)
798
799 """
800 def __init__(self, n, field=QQ, **kwargs):
801 V = VectorSpace(field, n)
802 mult_table = [ [ V.gen(i)*(i == j) for j in xrange(n) ]
803 for i in xrange(n) ]
804
805 fdeja = super(RealCartesianProductEJA, self)
806 return fdeja.__init__(field, mult_table, rank=n, **kwargs)
807
808 def inner_product(self, x, y):
809 """
810 Faster to reimplement than to use natural representations.
811
812 SETUP::
813
814 sage: from mjo.eja.eja_algebra import RealCartesianProductEJA
815
816 TESTS:
817
818 Ensure that this is the usual inner product for the algebras
819 over `R^n`::
820
821 sage: set_random_seed()
822 sage: J = RealCartesianProductEJA.random_instance()
823 sage: x,y = J.random_elements(2)
824 sage: X = x.natural_representation()
825 sage: Y = y.natural_representation()
826 sage: x.inner_product(y) == J.natural_inner_product(X,Y)
827 True
828
829 """
830 return x.to_vector().inner_product(y.to_vector())
831
832
833 def random_eja():
834 """
835 Return a "random" finite-dimensional Euclidean Jordan Algebra.
836
837 ALGORITHM:
838
839 For now, we choose a random natural number ``n`` (greater than zero)
840 and then give you back one of the following:
841
842 * The cartesian product of the rational numbers ``n`` times; this is
843 ``QQ^n`` with the Hadamard product.
844
845 * The Jordan spin algebra on ``QQ^n``.
846
847 * The ``n``-by-``n`` rational symmetric matrices with the symmetric
848 product.
849
850 * The ``n``-by-``n`` complex-rational Hermitian matrices embedded
851 in the space of ``2n``-by-``2n`` real symmetric matrices.
852
853 * The ``n``-by-``n`` quaternion-rational Hermitian matrices embedded
854 in the space of ``4n``-by-``4n`` real symmetric matrices.
855
856 Later this might be extended to return Cartesian products of the
857 EJAs above.
858
859 SETUP::
860
861 sage: from mjo.eja.eja_algebra import random_eja
862
863 TESTS::
864
865 sage: random_eja()
866 Euclidean Jordan algebra of dimension...
867
868 """
869 classname = choice([RealCartesianProductEJA,
870 JordanSpinEJA,
871 RealSymmetricEJA,
872 ComplexHermitianEJA,
873 QuaternionHermitianEJA])
874 return classname.random_instance()
875
876
877
878
879
880
881 class MatrixEuclideanJordanAlgebra(FiniteDimensionalEuclideanJordanAlgebra):
882 @staticmethod
883 def _max_test_case_size():
884 # Play it safe, since this will be squared and the underlying
885 # field can have dimension 4 (quaternions) too.
886 return 2
887
888 @classmethod
889 def _denormalized_basis(cls, n, field):
890 raise NotImplementedError
891
892 def __init__(self, n, field=QQ, normalize_basis=True, **kwargs):
893 S = self._denormalized_basis(n, field)
894
895 # Used in this class's fast _charpoly_coeff() override.
896 self._basis_normalizers = None
897
898 if n > 1 and normalize_basis:
899 # We'll need sqrt(2) to normalize the basis, and this
900 # winds up in the multiplication table, so the whole
901 # algebra needs to be over the field extension.
902 R = PolynomialRing(field, 'z')
903 z = R.gen()
904 p = z**2 - 2
905 if p.is_irreducible():
906 field = NumberField(p, 'sqrt2', embedding=RLF(2).sqrt())
907 S = [ s.change_ring(field) for s in S ]
908 self._basis_normalizers = tuple(
909 ~(self.natural_inner_product(s,s).sqrt()) for s in S )
910 S = tuple( s*c for (s,c) in zip(S,self._basis_normalizers) )
911
912 Qs = self.multiplication_table_from_matrix_basis(S)
913
914 fdeja = super(MatrixEuclideanJordanAlgebra, self)
915 return fdeja.__init__(field,
916 Qs,
917 rank=n,
918 natural_basis=S,
919 **kwargs)
920
921
922 @cached_method
923 def _charpoly_coeff(self, i):
924 """
925 Override the parent method with something that tries to compute
926 over a faster (non-extension) field.
927 """
928 if self._basis_normalizers is None:
929 # We didn't normalize, so assume that the basis we started
930 # with had entries in a nice field.
931 return super(MatrixEuclideanJordanAlgebra, self)._charpoly_coeff(i)
932 else:
933 # If we didn't unembed first, this number would be wrong
934 # by a power-of-two factor for complex/quaternion matrices.
935 n = self.real_unembed(self.natural_basis_space().zero()).nrows()
936 field = self.base_ring().base_ring() # yeeeeaaaahhh
937 J = self.__class__(n, field, False)
938 (_,x,_,_) = J._charpoly_matrix_system()
939 p = J._charpoly_coeff(i)
940 # p might be missing some vars, have to substitute "optionally"
941 pairs = izip(x.base_ring().gens(), self._basis_normalizers)
942 substitutions = { v: v*c for (v,c) in pairs }
943 return p.subs(substitutions)
944
945
946 @staticmethod
947 def multiplication_table_from_matrix_basis(basis):
948 """
949 At least three of the five simple Euclidean Jordan algebras have the
950 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
951 multiplication on the right is matrix multiplication. Given a basis
952 for the underlying matrix space, this function returns a
953 multiplication table (obtained by looping through the basis
954 elements) for an algebra of those matrices.
955 """
956 # In S^2, for example, we nominally have four coordinates even
957 # though the space is of dimension three only. The vector space V
958 # is supposed to hold the entire long vector, and the subspace W
959 # of V will be spanned by the vectors that arise from symmetric
960 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
961 field = basis[0].base_ring()
962 dimension = basis[0].nrows()
963
964 V = VectorSpace(field, dimension**2)
965 W = V.span_of_basis( _mat2vec(s) for s in basis )
966 n = len(basis)
967 mult_table = [[W.zero() for j in xrange(n)] for i in xrange(n)]
968 for i in xrange(n):
969 for j in xrange(n):
970 mat_entry = (basis[i]*basis[j] + basis[j]*basis[i])/2
971 mult_table[i][j] = W.coordinate_vector(_mat2vec(mat_entry))
972
973 return mult_table
974
975
976 @staticmethod
977 def real_embed(M):
978 """
979 Embed the matrix ``M`` into a space of real matrices.
980
981 The matrix ``M`` can have entries in any field at the moment:
982 the real numbers, complex numbers, or quaternions. And although
983 they are not a field, we can probably support octonions at some
984 point, too. This function returns a real matrix that "acts like"
985 the original with respect to matrix multiplication; i.e.
986
987 real_embed(M*N) = real_embed(M)*real_embed(N)
988
989 """
990 raise NotImplementedError
991
992
993 @staticmethod
994 def real_unembed(M):
995 """
996 The inverse of :meth:`real_embed`.
997 """
998 raise NotImplementedError
999
1000
1001 @classmethod
1002 def natural_inner_product(cls,X,Y):
1003 Xu = cls.real_unembed(X)
1004 Yu = cls.real_unembed(Y)
1005 tr = (Xu*Yu).trace()
1006 if tr in RLF:
1007 # It's real already.
1008 return tr
1009
1010 # Otherwise, try the thing that works for complex numbers; and
1011 # if that doesn't work, the thing that works for quaternions.
1012 try:
1013 return tr.vector()[0] # real part, imag part is index 1
1014 except AttributeError:
1015 # A quaternions doesn't have a vector() method, but does
1016 # have coefficient_tuple() method that returns the
1017 # coefficients of 1, i, j, and k -- in that order.
1018 return tr.coefficient_tuple()[0]
1019
1020
1021 class RealMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1022 @staticmethod
1023 def real_embed(M):
1024 """
1025 Embed the matrix ``M`` into a space of real matrices.
1026
1027 The matrix ``M`` can have entries in any field at the moment:
1028 the real numbers, complex numbers, or quaternions. And although
1029 they are not a field, we can probably support octonions at some
1030 point, too. This function returns a real matrix that "acts like"
1031 the original with respect to matrix multiplication; i.e.
1032
1033 real_embed(M*N) = real_embed(M)*real_embed(N)
1034
1035 """
1036 return M
1037
1038
1039 @staticmethod
1040 def real_unembed(M):
1041 """
1042 The inverse of :meth:`real_embed`.
1043 """
1044 return M
1045
1046
1047 class RealSymmetricEJA(RealMatrixEuclideanJordanAlgebra):
1048 """
1049 The rank-n simple EJA consisting of real symmetric n-by-n
1050 matrices, the usual symmetric Jordan product, and the trace inner
1051 product. It has dimension `(n^2 + n)/2` over the reals.
1052
1053 SETUP::
1054
1055 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1056
1057 EXAMPLES::
1058
1059 sage: J = RealSymmetricEJA(2)
1060 sage: e0, e1, e2 = J.gens()
1061 sage: e0*e0
1062 e0
1063 sage: e1*e1
1064 1/2*e0 + 1/2*e2
1065 sage: e2*e2
1066 e2
1067
1068 TESTS:
1069
1070 The dimension of this algebra is `(n^2 + n) / 2`::
1071
1072 sage: set_random_seed()
1073 sage: n_max = RealSymmetricEJA._max_test_case_size()
1074 sage: n = ZZ.random_element(1, n_max)
1075 sage: J = RealSymmetricEJA(n)
1076 sage: J.dimension() == (n^2 + n)/2
1077 True
1078
1079 The Jordan multiplication is what we think it is::
1080
1081 sage: set_random_seed()
1082 sage: J = RealSymmetricEJA.random_instance()
1083 sage: x,y = J.random_elements(2)
1084 sage: actual = (x*y).natural_representation()
1085 sage: X = x.natural_representation()
1086 sage: Y = y.natural_representation()
1087 sage: expected = (X*Y + Y*X)/2
1088 sage: actual == expected
1089 True
1090 sage: J(expected) == x*y
1091 True
1092
1093 We can change the generator prefix::
1094
1095 sage: RealSymmetricEJA(3, prefix='q').gens()
1096 (q0, q1, q2, q3, q4, q5)
1097
1098 Our natural basis is normalized with respect to the natural inner
1099 product unless we specify otherwise::
1100
1101 sage: set_random_seed()
1102 sage: J = RealSymmetricEJA.random_instance()
1103 sage: all( b.norm() == 1 for b in J.gens() )
1104 True
1105
1106 Since our natural basis is normalized with respect to the natural
1107 inner product, and since we know that this algebra is an EJA, any
1108 left-multiplication operator's matrix will be symmetric because
1109 natural->EJA basis representation is an isometry and within the EJA
1110 the operator is self-adjoint by the Jordan axiom::
1111
1112 sage: set_random_seed()
1113 sage: x = RealSymmetricEJA.random_instance().random_element()
1114 sage: x.operator().matrix().is_symmetric()
1115 True
1116
1117 """
1118 @classmethod
1119 def _denormalized_basis(cls, n, field):
1120 """
1121 Return a basis for the space of real symmetric n-by-n matrices.
1122
1123 SETUP::
1124
1125 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1126
1127 TESTS::
1128
1129 sage: set_random_seed()
1130 sage: n = ZZ.random_element(1,5)
1131 sage: B = RealSymmetricEJA._denormalized_basis(n,QQ)
1132 sage: all( M.is_symmetric() for M in B)
1133 True
1134
1135 """
1136 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1137 # coordinates.
1138 S = []
1139 for i in xrange(n):
1140 for j in xrange(i+1):
1141 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1142 if i == j:
1143 Sij = Eij
1144 else:
1145 Sij = Eij + Eij.transpose()
1146 S.append(Sij)
1147 return tuple(S)
1148
1149
1150 @staticmethod
1151 def _max_test_case_size():
1152 return 4 # Dimension 10
1153
1154
1155
1156 class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1157 @staticmethod
1158 def real_embed(M):
1159 """
1160 Embed the n-by-n complex matrix ``M`` into the space of real
1161 matrices of size 2n-by-2n via the map the sends each entry `z = a +
1162 bi` to the block matrix ``[[a,b],[-b,a]]``.
1163
1164 SETUP::
1165
1166 sage: from mjo.eja.eja_algebra import \
1167 ....: ComplexMatrixEuclideanJordanAlgebra
1168
1169 EXAMPLES::
1170
1171 sage: F = QuadraticField(-1, 'i')
1172 sage: x1 = F(4 - 2*i)
1173 sage: x2 = F(1 + 2*i)
1174 sage: x3 = F(-i)
1175 sage: x4 = F(6)
1176 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
1177 sage: ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1178 [ 4 -2| 1 2]
1179 [ 2 4|-2 1]
1180 [-----+-----]
1181 [ 0 -1| 6 0]
1182 [ 1 0| 0 6]
1183
1184 TESTS:
1185
1186 Embedding is a homomorphism (isomorphism, in fact)::
1187
1188 sage: set_random_seed()
1189 sage: n_max = ComplexMatrixEuclideanJordanAlgebra._max_test_case_size()
1190 sage: n = ZZ.random_element(n_max)
1191 sage: F = QuadraticField(-1, 'i')
1192 sage: X = random_matrix(F, n)
1193 sage: Y = random_matrix(F, n)
1194 sage: Xe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X)
1195 sage: Ye = ComplexMatrixEuclideanJordanAlgebra.real_embed(Y)
1196 sage: XYe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1197 sage: Xe*Ye == XYe
1198 True
1199
1200 """
1201 n = M.nrows()
1202 if M.ncols() != n:
1203 raise ValueError("the matrix 'M' must be square")
1204 field = M.base_ring()
1205 blocks = []
1206 for z in M.list():
1207 a = z.vector()[0] # real part, I guess
1208 b = z.vector()[1] # imag part, I guess
1209 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1210
1211 # We can drop the imaginaries here.
1212 return matrix.block(field.base_ring(), n, blocks)
1213
1214
1215 @staticmethod
1216 def real_unembed(M):
1217 """
1218 The inverse of _embed_complex_matrix().
1219
1220 SETUP::
1221
1222 sage: from mjo.eja.eja_algebra import \
1223 ....: ComplexMatrixEuclideanJordanAlgebra
1224
1225 EXAMPLES::
1226
1227 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1228 ....: [-2, 1, -4, 3],
1229 ....: [ 9, 10, 11, 12],
1230 ....: [-10, 9, -12, 11] ])
1231 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(A)
1232 [ 2*i + 1 4*i + 3]
1233 [ 10*i + 9 12*i + 11]
1234
1235 TESTS:
1236
1237 Unembedding is the inverse of embedding::
1238
1239 sage: set_random_seed()
1240 sage: F = QuadraticField(-1, 'i')
1241 sage: M = random_matrix(F, 3)
1242 sage: Me = ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1243 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1244 True
1245
1246 """
1247 n = ZZ(M.nrows())
1248 if M.ncols() != n:
1249 raise ValueError("the matrix 'M' must be square")
1250 if not n.mod(2).is_zero():
1251 raise ValueError("the matrix 'M' must be a complex embedding")
1252
1253 field = M.base_ring() # This should already have sqrt2
1254 R = PolynomialRing(field, 'z')
1255 z = R.gen()
1256 F = NumberField(z**2 + 1,'i', embedding=CLF(-1).sqrt())
1257 i = F.gen()
1258
1259 # Go top-left to bottom-right (reading order), converting every
1260 # 2-by-2 block we see to a single complex element.
1261 elements = []
1262 for k in xrange(n/2):
1263 for j in xrange(n/2):
1264 submat = M[2*k:2*k+2,2*j:2*j+2]
1265 if submat[0,0] != submat[1,1]:
1266 raise ValueError('bad on-diagonal submatrix')
1267 if submat[0,1] != -submat[1,0]:
1268 raise ValueError('bad off-diagonal submatrix')
1269 z = submat[0,0] + submat[0,1]*i
1270 elements.append(z)
1271
1272 return matrix(F, n/2, elements)
1273
1274
1275 @classmethod
1276 def natural_inner_product(cls,X,Y):
1277 """
1278 Compute a natural inner product in this algebra directly from
1279 its real embedding.
1280
1281 SETUP::
1282
1283 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1284
1285 TESTS:
1286
1287 This gives the same answer as the slow, default method implemented
1288 in :class:`MatrixEuclideanJordanAlgebra`::
1289
1290 sage: set_random_seed()
1291 sage: J = ComplexHermitianEJA.random_instance()
1292 sage: x,y = J.random_elements(2)
1293 sage: Xe = x.natural_representation()
1294 sage: Ye = y.natural_representation()
1295 sage: X = ComplexHermitianEJA.real_unembed(Xe)
1296 sage: Y = ComplexHermitianEJA.real_unembed(Ye)
1297 sage: expected = (X*Y).trace().vector()[0]
1298 sage: actual = ComplexHermitianEJA.natural_inner_product(Xe,Ye)
1299 sage: actual == expected
1300 True
1301
1302 """
1303 return RealMatrixEuclideanJordanAlgebra.natural_inner_product(X,Y)/2
1304
1305
1306 class ComplexHermitianEJA(ComplexMatrixEuclideanJordanAlgebra):
1307 """
1308 The rank-n simple EJA consisting of complex Hermitian n-by-n
1309 matrices over the real numbers, the usual symmetric Jordan product,
1310 and the real-part-of-trace inner product. It has dimension `n^2` over
1311 the reals.
1312
1313 SETUP::
1314
1315 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1316
1317 TESTS:
1318
1319 The dimension of this algebra is `n^2`::
1320
1321 sage: set_random_seed()
1322 sage: n_max = ComplexHermitianEJA._max_test_case_size()
1323 sage: n = ZZ.random_element(1, n_max)
1324 sage: J = ComplexHermitianEJA(n)
1325 sage: J.dimension() == n^2
1326 True
1327
1328 The Jordan multiplication is what we think it is::
1329
1330 sage: set_random_seed()
1331 sage: J = ComplexHermitianEJA.random_instance()
1332 sage: x,y = J.random_elements(2)
1333 sage: actual = (x*y).natural_representation()
1334 sage: X = x.natural_representation()
1335 sage: Y = y.natural_representation()
1336 sage: expected = (X*Y + Y*X)/2
1337 sage: actual == expected
1338 True
1339 sage: J(expected) == x*y
1340 True
1341
1342 We can change the generator prefix::
1343
1344 sage: ComplexHermitianEJA(2, prefix='z').gens()
1345 (z0, z1, z2, z3)
1346
1347 Our natural basis is normalized with respect to the natural inner
1348 product unless we specify otherwise::
1349
1350 sage: set_random_seed()
1351 sage: J = ComplexHermitianEJA.random_instance()
1352 sage: all( b.norm() == 1 for b in J.gens() )
1353 True
1354
1355 Since our natural basis is normalized with respect to the natural
1356 inner product, and since we know that this algebra is an EJA, any
1357 left-multiplication operator's matrix will be symmetric because
1358 natural->EJA basis representation is an isometry and within the EJA
1359 the operator is self-adjoint by the Jordan axiom::
1360
1361 sage: set_random_seed()
1362 sage: x = ComplexHermitianEJA.random_instance().random_element()
1363 sage: x.operator().matrix().is_symmetric()
1364 True
1365
1366 """
1367 @classmethod
1368 def _denormalized_basis(cls, n, field):
1369 """
1370 Returns a basis for the space of complex Hermitian n-by-n matrices.
1371
1372 Why do we embed these? Basically, because all of numerical linear
1373 algebra assumes that you're working with vectors consisting of `n`
1374 entries from a field and scalars from the same field. There's no way
1375 to tell SageMath that (for example) the vectors contain complex
1376 numbers, while the scalar field is real.
1377
1378 SETUP::
1379
1380 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1381
1382 TESTS::
1383
1384 sage: set_random_seed()
1385 sage: n = ZZ.random_element(1,5)
1386 sage: field = QuadraticField(2, 'sqrt2')
1387 sage: B = ComplexHermitianEJA._denormalized_basis(n, field)
1388 sage: all( M.is_symmetric() for M in B)
1389 True
1390
1391 """
1392 R = PolynomialRing(field, 'z')
1393 z = R.gen()
1394 F = NumberField(z**2 + 1, 'I', embedding=CLF(-1).sqrt())
1395 I = F.gen()
1396
1397 # This is like the symmetric case, but we need to be careful:
1398 #
1399 # * We want conjugate-symmetry, not just symmetry.
1400 # * The diagonal will (as a result) be real.
1401 #
1402 S = []
1403 for i in xrange(n):
1404 for j in xrange(i+1):
1405 Eij = matrix(F, n, lambda k,l: k==i and l==j)
1406 if i == j:
1407 Sij = cls.real_embed(Eij)
1408 S.append(Sij)
1409 else:
1410 # The second one has a minus because it's conjugated.
1411 Sij_real = cls.real_embed(Eij + Eij.transpose())
1412 S.append(Sij_real)
1413 Sij_imag = cls.real_embed(I*Eij - I*Eij.transpose())
1414 S.append(Sij_imag)
1415
1416 # Since we embedded these, we can drop back to the "field" that we
1417 # started with instead of the complex extension "F".
1418 return tuple( s.change_ring(field) for s in S )
1419
1420
1421
1422 class QuaternionMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1423 @staticmethod
1424 def real_embed(M):
1425 """
1426 Embed the n-by-n quaternion matrix ``M`` into the space of real
1427 matrices of size 4n-by-4n by first sending each quaternion entry `z
1428 = a + bi + cj + dk` to the block-complex matrix ``[[a + bi,
1429 c+di],[-c + di, a-bi]]`, and then embedding those into a real
1430 matrix.
1431
1432 SETUP::
1433
1434 sage: from mjo.eja.eja_algebra import \
1435 ....: QuaternionMatrixEuclideanJordanAlgebra
1436
1437 EXAMPLES::
1438
1439 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1440 sage: i,j,k = Q.gens()
1441 sage: x = 1 + 2*i + 3*j + 4*k
1442 sage: M = matrix(Q, 1, [[x]])
1443 sage: QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1444 [ 1 2 3 4]
1445 [-2 1 -4 3]
1446 [-3 4 1 -2]
1447 [-4 -3 2 1]
1448
1449 Embedding is a homomorphism (isomorphism, in fact)::
1450
1451 sage: set_random_seed()
1452 sage: n_max = QuaternionMatrixEuclideanJordanAlgebra._max_test_case_size()
1453 sage: n = ZZ.random_element(n_max)
1454 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1455 sage: X = random_matrix(Q, n)
1456 sage: Y = random_matrix(Q, n)
1457 sage: Xe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X)
1458 sage: Ye = QuaternionMatrixEuclideanJordanAlgebra.real_embed(Y)
1459 sage: XYe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1460 sage: Xe*Ye == XYe
1461 True
1462
1463 """
1464 quaternions = M.base_ring()
1465 n = M.nrows()
1466 if M.ncols() != n:
1467 raise ValueError("the matrix 'M' must be square")
1468
1469 F = QuadraticField(-1, 'i')
1470 i = F.gen()
1471
1472 blocks = []
1473 for z in M.list():
1474 t = z.coefficient_tuple()
1475 a = t[0]
1476 b = t[1]
1477 c = t[2]
1478 d = t[3]
1479 cplxM = matrix(F, 2, [[ a + b*i, c + d*i],
1480 [-c + d*i, a - b*i]])
1481 realM = ComplexMatrixEuclideanJordanAlgebra.real_embed(cplxM)
1482 blocks.append(realM)
1483
1484 # We should have real entries by now, so use the realest field
1485 # we've got for the return value.
1486 return matrix.block(quaternions.base_ring(), n, blocks)
1487
1488
1489
1490 @staticmethod
1491 def real_unembed(M):
1492 """
1493 The inverse of _embed_quaternion_matrix().
1494
1495 SETUP::
1496
1497 sage: from mjo.eja.eja_algebra import \
1498 ....: QuaternionMatrixEuclideanJordanAlgebra
1499
1500 EXAMPLES::
1501
1502 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
1503 ....: [-2, 1, -4, 3],
1504 ....: [-3, 4, 1, -2],
1505 ....: [-4, -3, 2, 1]])
1506 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(M)
1507 [1 + 2*i + 3*j + 4*k]
1508
1509 TESTS:
1510
1511 Unembedding is the inverse of embedding::
1512
1513 sage: set_random_seed()
1514 sage: Q = QuaternionAlgebra(QQ, -1, -1)
1515 sage: M = random_matrix(Q, 3)
1516 sage: Me = QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1517 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1518 True
1519
1520 """
1521 n = ZZ(M.nrows())
1522 if M.ncols() != n:
1523 raise ValueError("the matrix 'M' must be square")
1524 if not n.mod(4).is_zero():
1525 raise ValueError("the matrix 'M' must be a complex embedding")
1526
1527 # Use the base ring of the matrix to ensure that its entries can be
1528 # multiplied by elements of the quaternion algebra.
1529 field = M.base_ring()
1530 Q = QuaternionAlgebra(field,-1,-1)
1531 i,j,k = Q.gens()
1532
1533 # Go top-left to bottom-right (reading order), converting every
1534 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
1535 # quaternion block.
1536 elements = []
1537 for l in xrange(n/4):
1538 for m in xrange(n/4):
1539 submat = ComplexMatrixEuclideanJordanAlgebra.real_unembed(
1540 M[4*l:4*l+4,4*m:4*m+4] )
1541 if submat[0,0] != submat[1,1].conjugate():
1542 raise ValueError('bad on-diagonal submatrix')
1543 if submat[0,1] != -submat[1,0].conjugate():
1544 raise ValueError('bad off-diagonal submatrix')
1545 z = submat[0,0].vector()[0] # real part
1546 z += submat[0,0].vector()[1]*i # imag part
1547 z += submat[0,1].vector()[0]*j # real part
1548 z += submat[0,1].vector()[1]*k # imag part
1549 elements.append(z)
1550
1551 return matrix(Q, n/4, elements)
1552
1553
1554 @classmethod
1555 def natural_inner_product(cls,X,Y):
1556 """
1557 Compute a natural inner product in this algebra directly from
1558 its real embedding.
1559
1560 SETUP::
1561
1562 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1563
1564 TESTS:
1565
1566 This gives the same answer as the slow, default method implemented
1567 in :class:`MatrixEuclideanJordanAlgebra`::
1568
1569 sage: set_random_seed()
1570 sage: J = QuaternionHermitianEJA.random_instance()
1571 sage: x,y = J.random_elements(2)
1572 sage: Xe = x.natural_representation()
1573 sage: Ye = y.natural_representation()
1574 sage: X = QuaternionHermitianEJA.real_unembed(Xe)
1575 sage: Y = QuaternionHermitianEJA.real_unembed(Ye)
1576 sage: expected = (X*Y).trace().coefficient_tuple()[0]
1577 sage: actual = QuaternionHermitianEJA.natural_inner_product(Xe,Ye)
1578 sage: actual == expected
1579 True
1580
1581 """
1582 return RealMatrixEuclideanJordanAlgebra.natural_inner_product(X,Y)/4
1583
1584
1585 class QuaternionHermitianEJA(QuaternionMatrixEuclideanJordanAlgebra):
1586 """
1587 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
1588 matrices, the usual symmetric Jordan product, and the
1589 real-part-of-trace inner product. It has dimension `2n^2 - n` over
1590 the reals.
1591
1592 SETUP::
1593
1594 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1595
1596 TESTS:
1597
1598 The dimension of this algebra is `2*n^2 - n`::
1599
1600 sage: set_random_seed()
1601 sage: n_max = QuaternionHermitianEJA._max_test_case_size()
1602 sage: n = ZZ.random_element(1, n_max)
1603 sage: J = QuaternionHermitianEJA(n)
1604 sage: J.dimension() == 2*(n^2) - n
1605 True
1606
1607 The Jordan multiplication is what we think it is::
1608
1609 sage: set_random_seed()
1610 sage: J = QuaternionHermitianEJA.random_instance()
1611 sage: x,y = J.random_elements(2)
1612 sage: actual = (x*y).natural_representation()
1613 sage: X = x.natural_representation()
1614 sage: Y = y.natural_representation()
1615 sage: expected = (X*Y + Y*X)/2
1616 sage: actual == expected
1617 True
1618 sage: J(expected) == x*y
1619 True
1620
1621 We can change the generator prefix::
1622
1623 sage: QuaternionHermitianEJA(2, prefix='a').gens()
1624 (a0, a1, a2, a3, a4, a5)
1625
1626 Our natural basis is normalized with respect to the natural inner
1627 product unless we specify otherwise::
1628
1629 sage: set_random_seed()
1630 sage: J = QuaternionHermitianEJA.random_instance()
1631 sage: all( b.norm() == 1 for b in J.gens() )
1632 True
1633
1634 Since our natural basis is normalized with respect to the natural
1635 inner product, and since we know that this algebra is an EJA, any
1636 left-multiplication operator's matrix will be symmetric because
1637 natural->EJA basis representation is an isometry and within the EJA
1638 the operator is self-adjoint by the Jordan axiom::
1639
1640 sage: set_random_seed()
1641 sage: x = QuaternionHermitianEJA.random_instance().random_element()
1642 sage: x.operator().matrix().is_symmetric()
1643 True
1644
1645 """
1646 @classmethod
1647 def _denormalized_basis(cls, n, field):
1648 """
1649 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
1650
1651 Why do we embed these? Basically, because all of numerical
1652 linear algebra assumes that you're working with vectors consisting
1653 of `n` entries from a field and scalars from the same field. There's
1654 no way to tell SageMath that (for example) the vectors contain
1655 complex numbers, while the scalar field is real.
1656
1657 SETUP::
1658
1659 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1660
1661 TESTS::
1662
1663 sage: set_random_seed()
1664 sage: n = ZZ.random_element(1,5)
1665 sage: B = QuaternionHermitianEJA._denormalized_basis(n,QQ)
1666 sage: all( M.is_symmetric() for M in B )
1667 True
1668
1669 """
1670 Q = QuaternionAlgebra(QQ,-1,-1)
1671 I,J,K = Q.gens()
1672
1673 # This is like the symmetric case, but we need to be careful:
1674 #
1675 # * We want conjugate-symmetry, not just symmetry.
1676 # * The diagonal will (as a result) be real.
1677 #
1678 S = []
1679 for i in xrange(n):
1680 for j in xrange(i+1):
1681 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
1682 if i == j:
1683 Sij = cls.real_embed(Eij)
1684 S.append(Sij)
1685 else:
1686 # The second, third, and fourth ones have a minus
1687 # because they're conjugated.
1688 Sij_real = cls.real_embed(Eij + Eij.transpose())
1689 S.append(Sij_real)
1690 Sij_I = cls.real_embed(I*Eij - I*Eij.transpose())
1691 S.append(Sij_I)
1692 Sij_J = cls.real_embed(J*Eij - J*Eij.transpose())
1693 S.append(Sij_J)
1694 Sij_K = cls.real_embed(K*Eij - K*Eij.transpose())
1695 S.append(Sij_K)
1696 return tuple(S)
1697
1698
1699
1700 class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra):
1701 """
1702 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
1703 with the usual inner product and jordan product ``x*y =
1704 (<x_bar,y_bar>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
1705 the reals.
1706
1707 SETUP::
1708
1709 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1710
1711 EXAMPLES:
1712
1713 This multiplication table can be verified by hand::
1714
1715 sage: J = JordanSpinEJA(4)
1716 sage: e0,e1,e2,e3 = J.gens()
1717 sage: e0*e0
1718 e0
1719 sage: e0*e1
1720 e1
1721 sage: e0*e2
1722 e2
1723 sage: e0*e3
1724 e3
1725 sage: e1*e2
1726 0
1727 sage: e1*e3
1728 0
1729 sage: e2*e3
1730 0
1731
1732 We can change the generator prefix::
1733
1734 sage: JordanSpinEJA(2, prefix='B').gens()
1735 (B0, B1)
1736
1737 """
1738 def __init__(self, n, field=QQ, **kwargs):
1739 V = VectorSpace(field, n)
1740 mult_table = [[V.zero() for j in xrange(n)] for i in xrange(n)]
1741 for i in xrange(n):
1742 for j in xrange(n):
1743 x = V.gen(i)
1744 y = V.gen(j)
1745 x0 = x[0]
1746 xbar = x[1:]
1747 y0 = y[0]
1748 ybar = y[1:]
1749 # z = x*y
1750 z0 = x.inner_product(y)
1751 zbar = y0*xbar + x0*ybar
1752 z = V([z0] + zbar.list())
1753 mult_table[i][j] = z
1754
1755 # The rank of the spin algebra is two, unless we're in a
1756 # one-dimensional ambient space (because the rank is bounded by
1757 # the ambient dimension).
1758 fdeja = super(JordanSpinEJA, self)
1759 return fdeja.__init__(field, mult_table, rank=min(n,2), **kwargs)
1760
1761 def inner_product(self, x, y):
1762 """
1763 Faster to reimplement than to use natural representations.
1764
1765 SETUP::
1766
1767 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1768
1769 TESTS:
1770
1771 Ensure that this is the usual inner product for the algebras
1772 over `R^n`::
1773
1774 sage: set_random_seed()
1775 sage: J = JordanSpinEJA.random_instance()
1776 sage: x,y = J.random_elements(2)
1777 sage: X = x.natural_representation()
1778 sage: Y = y.natural_representation()
1779 sage: x.inner_product(y) == J.natural_inner_product(X,Y)
1780 True
1781
1782 """
1783 return x.to_vector().inner_product(y.to_vector())