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