]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: ensure that we can construct RealSymmetricEJA over AA and RR.
[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
943 # Do this over the rationals and convert back at the end.
944 J = MatrixEuclideanJordanAlgebra(QQ,
945 basis,
946 self.rank(),
947 normalize_basis=False)
948 (_,x,_,_) = J._charpoly_matrix_system()
949 p = J._charpoly_coeff(i)
950 # p might be missing some vars, have to substitute "optionally"
951 pairs = izip(x.base_ring().gens(), self._basis_normalizers)
952 substitutions = { v: v*c for (v,c) in pairs }
953 result = p.subs(substitutions)
954
955 # The result of "subs" can be either a coefficient-ring
956 # element or a polynomial. Gotta handle both cases.
957 if result in QQ:
958 return self.base_ring()(result)
959 else:
960 return result.change_ring(self.base_ring())
961
962
963 @staticmethod
964 def multiplication_table_from_matrix_basis(basis):
965 """
966 At least three of the five simple Euclidean Jordan algebras have the
967 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
968 multiplication on the right is matrix multiplication. Given a basis
969 for the underlying matrix space, this function returns a
970 multiplication table (obtained by looping through the basis
971 elements) for an algebra of those matrices.
972 """
973 # In S^2, for example, we nominally have four coordinates even
974 # though the space is of dimension three only. The vector space V
975 # is supposed to hold the entire long vector, and the subspace W
976 # of V will be spanned by the vectors that arise from symmetric
977 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
978 field = basis[0].base_ring()
979 dimension = basis[0].nrows()
980
981 V = VectorSpace(field, dimension**2)
982 W = V.span_of_basis( _mat2vec(s) for s in basis )
983 n = len(basis)
984 mult_table = [[W.zero() for j in xrange(n)] for i in xrange(n)]
985 for i in xrange(n):
986 for j in xrange(n):
987 mat_entry = (basis[i]*basis[j] + basis[j]*basis[i])/2
988 mult_table[i][j] = W.coordinate_vector(_mat2vec(mat_entry))
989
990 return mult_table
991
992
993 @staticmethod
994 def real_embed(M):
995 """
996 Embed the matrix ``M`` into a space of real matrices.
997
998 The matrix ``M`` can have entries in any field at the moment:
999 the real numbers, complex numbers, or quaternions. And although
1000 they are not a field, we can probably support octonions at some
1001 point, too. This function returns a real matrix that "acts like"
1002 the original with respect to matrix multiplication; i.e.
1003
1004 real_embed(M*N) = real_embed(M)*real_embed(N)
1005
1006 """
1007 raise NotImplementedError
1008
1009
1010 @staticmethod
1011 def real_unembed(M):
1012 """
1013 The inverse of :meth:`real_embed`.
1014 """
1015 raise NotImplementedError
1016
1017
1018 @classmethod
1019 def natural_inner_product(cls,X,Y):
1020 Xu = cls.real_unembed(X)
1021 Yu = cls.real_unembed(Y)
1022 tr = (Xu*Yu).trace()
1023
1024 if tr in RLF:
1025 # It's real already.
1026 return tr
1027
1028 # Otherwise, try the thing that works for complex numbers; and
1029 # if that doesn't work, the thing that works for quaternions.
1030 try:
1031 return tr.vector()[0] # real part, imag part is index 1
1032 except AttributeError:
1033 # A quaternions doesn't have a vector() method, but does
1034 # have coefficient_tuple() method that returns the
1035 # coefficients of 1, i, j, and k -- in that order.
1036 return tr.coefficient_tuple()[0]
1037
1038
1039 class RealMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1040 @staticmethod
1041 def real_embed(M):
1042 """
1043 The identity function, for embedding real matrices into real
1044 matrices.
1045 """
1046 return M
1047
1048 @staticmethod
1049 def real_unembed(M):
1050 """
1051 The identity function, for unembedding real matrices from real
1052 matrices.
1053 """
1054 return M
1055
1056
1057 class RealSymmetricEJA(RealMatrixEuclideanJordanAlgebra, KnownRankEJA):
1058 """
1059 The rank-n simple EJA consisting of real symmetric n-by-n
1060 matrices, the usual symmetric Jordan product, and the trace inner
1061 product. It has dimension `(n^2 + n)/2` over the reals.
1062
1063 SETUP::
1064
1065 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1066
1067 EXAMPLES::
1068
1069 sage: J = RealSymmetricEJA(2)
1070 sage: e0, e1, e2 = J.gens()
1071 sage: e0*e0
1072 e0
1073 sage: e1*e1
1074 1/2*e0 + 1/2*e2
1075 sage: e2*e2
1076 e2
1077
1078 In theory, our "field" can be any subfield of the reals::
1079
1080 sage: RealSymmetricEJA(2, AA)
1081 Euclidean Jordan algebra of dimension 3 over Algebraic Real Field
1082 sage: RealSymmetricEJA(2, RR)
1083 Euclidean Jordan algebra of dimension 3 over Real Field with
1084 53 bits of precision
1085
1086 TESTS:
1087
1088 The dimension of this algebra is `(n^2 + n) / 2`::
1089
1090 sage: set_random_seed()
1091 sage: n_max = RealSymmetricEJA._max_test_case_size()
1092 sage: n = ZZ.random_element(1, n_max)
1093 sage: J = RealSymmetricEJA(n)
1094 sage: J.dimension() == (n^2 + n)/2
1095 True
1096
1097 The Jordan multiplication is what we think it is::
1098
1099 sage: set_random_seed()
1100 sage: J = RealSymmetricEJA.random_instance()
1101 sage: x,y = J.random_elements(2)
1102 sage: actual = (x*y).natural_representation()
1103 sage: X = x.natural_representation()
1104 sage: Y = y.natural_representation()
1105 sage: expected = (X*Y + Y*X)/2
1106 sage: actual == expected
1107 True
1108 sage: J(expected) == x*y
1109 True
1110
1111 We can change the generator prefix::
1112
1113 sage: RealSymmetricEJA(3, prefix='q').gens()
1114 (q0, q1, q2, q3, q4, q5)
1115
1116 Our natural basis is normalized with respect to the natural inner
1117 product unless we specify otherwise::
1118
1119 sage: set_random_seed()
1120 sage: J = RealSymmetricEJA.random_instance()
1121 sage: all( b.norm() == 1 for b in J.gens() )
1122 True
1123
1124 Since our natural basis is normalized with respect to the natural
1125 inner product, and since we know that this algebra is an EJA, any
1126 left-multiplication operator's matrix will be symmetric because
1127 natural->EJA basis representation is an isometry and within the EJA
1128 the operator is self-adjoint by the Jordan axiom::
1129
1130 sage: set_random_seed()
1131 sage: x = RealSymmetricEJA.random_instance().random_element()
1132 sage: x.operator().matrix().is_symmetric()
1133 True
1134
1135 """
1136 @classmethod
1137 def _denormalized_basis(cls, n, field):
1138 """
1139 Return a basis for the space of real symmetric n-by-n matrices.
1140
1141 SETUP::
1142
1143 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1144
1145 TESTS::
1146
1147 sage: set_random_seed()
1148 sage: n = ZZ.random_element(1,5)
1149 sage: B = RealSymmetricEJA._denormalized_basis(n,QQ)
1150 sage: all( M.is_symmetric() for M in B)
1151 True
1152
1153 """
1154 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1155 # coordinates.
1156 S = []
1157 for i in xrange(n):
1158 for j in xrange(i+1):
1159 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1160 if i == j:
1161 Sij = Eij
1162 else:
1163 Sij = Eij + Eij.transpose()
1164 S.append(Sij)
1165 return S
1166
1167
1168 @staticmethod
1169 def _max_test_case_size():
1170 return 4 # Dimension 10
1171
1172
1173 def __init__(self, n, field=QQ, **kwargs):
1174 basis = self._denormalized_basis(n, field)
1175 super(RealSymmetricEJA, self).__init__(field, basis, n, **kwargs)
1176
1177
1178 class ComplexMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1179 @staticmethod
1180 def real_embed(M):
1181 """
1182 Embed the n-by-n complex matrix ``M`` into the space of real
1183 matrices of size 2n-by-2n via the map the sends each entry `z = a +
1184 bi` to the block matrix ``[[a,b],[-b,a]]``.
1185
1186 SETUP::
1187
1188 sage: from mjo.eja.eja_algebra import \
1189 ....: ComplexMatrixEuclideanJordanAlgebra
1190
1191 EXAMPLES::
1192
1193 sage: F = QuadraticField(-1, 'i')
1194 sage: x1 = F(4 - 2*i)
1195 sage: x2 = F(1 + 2*i)
1196 sage: x3 = F(-i)
1197 sage: x4 = F(6)
1198 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
1199 sage: ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1200 [ 4 -2| 1 2]
1201 [ 2 4|-2 1]
1202 [-----+-----]
1203 [ 0 -1| 6 0]
1204 [ 1 0| 0 6]
1205
1206 TESTS:
1207
1208 Embedding is a homomorphism (isomorphism, in fact)::
1209
1210 sage: set_random_seed()
1211 sage: n_max = ComplexMatrixEuclideanJordanAlgebra._max_test_case_size()
1212 sage: n = ZZ.random_element(n_max)
1213 sage: F = QuadraticField(-1, 'i')
1214 sage: X = random_matrix(F, n)
1215 sage: Y = random_matrix(F, n)
1216 sage: Xe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X)
1217 sage: Ye = ComplexMatrixEuclideanJordanAlgebra.real_embed(Y)
1218 sage: XYe = ComplexMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1219 sage: Xe*Ye == XYe
1220 True
1221
1222 """
1223 n = M.nrows()
1224 if M.ncols() != n:
1225 raise ValueError("the matrix 'M' must be square")
1226 field = M.base_ring()
1227 blocks = []
1228 for z in M.list():
1229 a = z.vector()[0] # real part, I guess
1230 b = z.vector()[1] # imag part, I guess
1231 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1232
1233 # We can drop the imaginaries here.
1234 return matrix.block(field.base_ring(), n, blocks)
1235
1236
1237 @staticmethod
1238 def real_unembed(M):
1239 """
1240 The inverse of _embed_complex_matrix().
1241
1242 SETUP::
1243
1244 sage: from mjo.eja.eja_algebra import \
1245 ....: ComplexMatrixEuclideanJordanAlgebra
1246
1247 EXAMPLES::
1248
1249 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1250 ....: [-2, 1, -4, 3],
1251 ....: [ 9, 10, 11, 12],
1252 ....: [-10, 9, -12, 11] ])
1253 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(A)
1254 [ 2*i + 1 4*i + 3]
1255 [ 10*i + 9 12*i + 11]
1256
1257 TESTS:
1258
1259 Unembedding is the inverse of embedding::
1260
1261 sage: set_random_seed()
1262 sage: F = QuadraticField(-1, 'i')
1263 sage: M = random_matrix(F, 3)
1264 sage: Me = ComplexMatrixEuclideanJordanAlgebra.real_embed(M)
1265 sage: ComplexMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1266 True
1267
1268 """
1269 n = ZZ(M.nrows())
1270 if M.ncols() != n:
1271 raise ValueError("the matrix 'M' must be square")
1272 if not n.mod(2).is_zero():
1273 raise ValueError("the matrix 'M' must be a complex embedding")
1274
1275 # If "M" was normalized, its base ring might have roots
1276 # adjoined and they can stick around after unembedding.
1277 field = M.base_ring()
1278 R = PolynomialRing(field, 'z')
1279 z = R.gen()
1280 F = NumberField(z**2 + 1,'i', embedding=CLF(-1).sqrt())
1281 i = F.gen()
1282
1283 # Go top-left to bottom-right (reading order), converting every
1284 # 2-by-2 block we see to a single complex element.
1285 elements = []
1286 for k in xrange(n/2):
1287 for j in xrange(n/2):
1288 submat = M[2*k:2*k+2,2*j:2*j+2]
1289 if submat[0,0] != submat[1,1]:
1290 raise ValueError('bad on-diagonal submatrix')
1291 if submat[0,1] != -submat[1,0]:
1292 raise ValueError('bad off-diagonal submatrix')
1293 z = submat[0,0] + submat[0,1]*i
1294 elements.append(z)
1295
1296 return matrix(F, n/2, elements)
1297
1298
1299 @classmethod
1300 def natural_inner_product(cls,X,Y):
1301 """
1302 Compute a natural inner product in this algebra directly from
1303 its real embedding.
1304
1305 SETUP::
1306
1307 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1308
1309 TESTS:
1310
1311 This gives the same answer as the slow, default method implemented
1312 in :class:`MatrixEuclideanJordanAlgebra`::
1313
1314 sage: set_random_seed()
1315 sage: J = ComplexHermitianEJA.random_instance()
1316 sage: x,y = J.random_elements(2)
1317 sage: Xe = x.natural_representation()
1318 sage: Ye = y.natural_representation()
1319 sage: X = ComplexHermitianEJA.real_unembed(Xe)
1320 sage: Y = ComplexHermitianEJA.real_unembed(Ye)
1321 sage: expected = (X*Y).trace().vector()[0]
1322 sage: actual = ComplexHermitianEJA.natural_inner_product(Xe,Ye)
1323 sage: actual == expected
1324 True
1325
1326 """
1327 return RealMatrixEuclideanJordanAlgebra.natural_inner_product(X,Y)/2
1328
1329
1330 class ComplexHermitianEJA(ComplexMatrixEuclideanJordanAlgebra, KnownRankEJA):
1331 """
1332 The rank-n simple EJA consisting of complex Hermitian n-by-n
1333 matrices over the real numbers, the usual symmetric Jordan product,
1334 and the real-part-of-trace inner product. It has dimension `n^2` over
1335 the reals.
1336
1337 SETUP::
1338
1339 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1340
1341 TESTS:
1342
1343 The dimension of this algebra is `n^2`::
1344
1345 sage: set_random_seed()
1346 sage: n_max = ComplexHermitianEJA._max_test_case_size()
1347 sage: n = ZZ.random_element(1, n_max)
1348 sage: J = ComplexHermitianEJA(n)
1349 sage: J.dimension() == n^2
1350 True
1351
1352 The Jordan multiplication is what we think it is::
1353
1354 sage: set_random_seed()
1355 sage: J = ComplexHermitianEJA.random_instance()
1356 sage: x,y = J.random_elements(2)
1357 sage: actual = (x*y).natural_representation()
1358 sage: X = x.natural_representation()
1359 sage: Y = y.natural_representation()
1360 sage: expected = (X*Y + Y*X)/2
1361 sage: actual == expected
1362 True
1363 sage: J(expected) == x*y
1364 True
1365
1366 We can change the generator prefix::
1367
1368 sage: ComplexHermitianEJA(2, prefix='z').gens()
1369 (z0, z1, z2, z3)
1370
1371 Our natural basis is normalized with respect to the natural inner
1372 product unless we specify otherwise::
1373
1374 sage: set_random_seed()
1375 sage: J = ComplexHermitianEJA.random_instance()
1376 sage: all( b.norm() == 1 for b in J.gens() )
1377 True
1378
1379 Since our natural basis is normalized with respect to the natural
1380 inner product, and since we know that this algebra is an EJA, any
1381 left-multiplication operator's matrix will be symmetric because
1382 natural->EJA basis representation is an isometry and within the EJA
1383 the operator is self-adjoint by the Jordan axiom::
1384
1385 sage: set_random_seed()
1386 sage: x = ComplexHermitianEJA.random_instance().random_element()
1387 sage: x.operator().matrix().is_symmetric()
1388 True
1389
1390 """
1391
1392 @classmethod
1393 def _denormalized_basis(cls, n, field):
1394 """
1395 Returns a basis for the space of complex Hermitian n-by-n matrices.
1396
1397 Why do we embed these? Basically, because all of numerical linear
1398 algebra assumes that you're working with vectors consisting of `n`
1399 entries from a field and scalars from the same field. There's no way
1400 to tell SageMath that (for example) the vectors contain complex
1401 numbers, while the scalar field is real.
1402
1403 SETUP::
1404
1405 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1406
1407 TESTS::
1408
1409 sage: set_random_seed()
1410 sage: n = ZZ.random_element(1,5)
1411 sage: field = QuadraticField(2, 'sqrt2')
1412 sage: B = ComplexHermitianEJA._denormalized_basis(n, field)
1413 sage: all( M.is_symmetric() for M in B)
1414 True
1415
1416 """
1417 R = PolynomialRing(field, 'z')
1418 z = R.gen()
1419 F = NumberField(z**2 + 1, 'I', embedding=CLF(-1).sqrt())
1420 I = F.gen()
1421
1422 # This is like the symmetric case, but we need to be careful:
1423 #
1424 # * We want conjugate-symmetry, not just symmetry.
1425 # * The diagonal will (as a result) be real.
1426 #
1427 S = []
1428 for i in xrange(n):
1429 for j in xrange(i+1):
1430 Eij = matrix(F, n, lambda k,l: k==i and l==j)
1431 if i == j:
1432 Sij = cls.real_embed(Eij)
1433 S.append(Sij)
1434 else:
1435 # The second one has a minus because it's conjugated.
1436 Sij_real = cls.real_embed(Eij + Eij.transpose())
1437 S.append(Sij_real)
1438 Sij_imag = cls.real_embed(I*Eij - I*Eij.transpose())
1439 S.append(Sij_imag)
1440
1441 # Since we embedded these, we can drop back to the "field" that we
1442 # started with instead of the complex extension "F".
1443 return ( s.change_ring(field) for s in S )
1444
1445
1446 def __init__(self, n, field=QQ, **kwargs):
1447 basis = self._denormalized_basis(n,field)
1448 super(ComplexHermitianEJA,self).__init__(field, basis, n, **kwargs)
1449
1450
1451 class QuaternionMatrixEuclideanJordanAlgebra(MatrixEuclideanJordanAlgebra):
1452 @staticmethod
1453 def real_embed(M):
1454 """
1455 Embed the n-by-n quaternion matrix ``M`` into the space of real
1456 matrices of size 4n-by-4n by first sending each quaternion entry `z
1457 = a + bi + cj + dk` to the block-complex matrix ``[[a + bi,
1458 c+di],[-c + di, a-bi]]`, and then embedding those into a real
1459 matrix.
1460
1461 SETUP::
1462
1463 sage: from mjo.eja.eja_algebra import \
1464 ....: QuaternionMatrixEuclideanJordanAlgebra
1465
1466 EXAMPLES::
1467
1468 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1469 sage: i,j,k = Q.gens()
1470 sage: x = 1 + 2*i + 3*j + 4*k
1471 sage: M = matrix(Q, 1, [[x]])
1472 sage: QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1473 [ 1 2 3 4]
1474 [-2 1 -4 3]
1475 [-3 4 1 -2]
1476 [-4 -3 2 1]
1477
1478 Embedding is a homomorphism (isomorphism, in fact)::
1479
1480 sage: set_random_seed()
1481 sage: n_max = QuaternionMatrixEuclideanJordanAlgebra._max_test_case_size()
1482 sage: n = ZZ.random_element(n_max)
1483 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1484 sage: X = random_matrix(Q, n)
1485 sage: Y = random_matrix(Q, n)
1486 sage: Xe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X)
1487 sage: Ye = QuaternionMatrixEuclideanJordanAlgebra.real_embed(Y)
1488 sage: XYe = QuaternionMatrixEuclideanJordanAlgebra.real_embed(X*Y)
1489 sage: Xe*Ye == XYe
1490 True
1491
1492 """
1493 quaternions = M.base_ring()
1494 n = M.nrows()
1495 if M.ncols() != n:
1496 raise ValueError("the matrix 'M' must be square")
1497
1498 F = QuadraticField(-1, 'i')
1499 i = F.gen()
1500
1501 blocks = []
1502 for z in M.list():
1503 t = z.coefficient_tuple()
1504 a = t[0]
1505 b = t[1]
1506 c = t[2]
1507 d = t[3]
1508 cplxM = matrix(F, 2, [[ a + b*i, c + d*i],
1509 [-c + d*i, a - b*i]])
1510 realM = ComplexMatrixEuclideanJordanAlgebra.real_embed(cplxM)
1511 blocks.append(realM)
1512
1513 # We should have real entries by now, so use the realest field
1514 # we've got for the return value.
1515 return matrix.block(quaternions.base_ring(), n, blocks)
1516
1517
1518
1519 @staticmethod
1520 def real_unembed(M):
1521 """
1522 The inverse of _embed_quaternion_matrix().
1523
1524 SETUP::
1525
1526 sage: from mjo.eja.eja_algebra import \
1527 ....: QuaternionMatrixEuclideanJordanAlgebra
1528
1529 EXAMPLES::
1530
1531 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
1532 ....: [-2, 1, -4, 3],
1533 ....: [-3, 4, 1, -2],
1534 ....: [-4, -3, 2, 1]])
1535 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(M)
1536 [1 + 2*i + 3*j + 4*k]
1537
1538 TESTS:
1539
1540 Unembedding is the inverse of embedding::
1541
1542 sage: set_random_seed()
1543 sage: Q = QuaternionAlgebra(QQ, -1, -1)
1544 sage: M = random_matrix(Q, 3)
1545 sage: Me = QuaternionMatrixEuclideanJordanAlgebra.real_embed(M)
1546 sage: QuaternionMatrixEuclideanJordanAlgebra.real_unembed(Me) == M
1547 True
1548
1549 """
1550 n = ZZ(M.nrows())
1551 if M.ncols() != n:
1552 raise ValueError("the matrix 'M' must be square")
1553 if not n.mod(4).is_zero():
1554 raise ValueError("the matrix 'M' must be a quaternion embedding")
1555
1556 # Use the base ring of the matrix to ensure that its entries can be
1557 # multiplied by elements of the quaternion algebra.
1558 field = M.base_ring()
1559 Q = QuaternionAlgebra(field,-1,-1)
1560 i,j,k = Q.gens()
1561
1562 # Go top-left to bottom-right (reading order), converting every
1563 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
1564 # quaternion block.
1565 elements = []
1566 for l in xrange(n/4):
1567 for m in xrange(n/4):
1568 submat = ComplexMatrixEuclideanJordanAlgebra.real_unembed(
1569 M[4*l:4*l+4,4*m:4*m+4] )
1570 if submat[0,0] != submat[1,1].conjugate():
1571 raise ValueError('bad on-diagonal submatrix')
1572 if submat[0,1] != -submat[1,0].conjugate():
1573 raise ValueError('bad off-diagonal submatrix')
1574 z = submat[0,0].vector()[0] # real part
1575 z += submat[0,0].vector()[1]*i # imag part
1576 z += submat[0,1].vector()[0]*j # real part
1577 z += submat[0,1].vector()[1]*k # imag part
1578 elements.append(z)
1579
1580 return matrix(Q, n/4, elements)
1581
1582
1583 @classmethod
1584 def natural_inner_product(cls,X,Y):
1585 """
1586 Compute a natural inner product in this algebra directly from
1587 its real embedding.
1588
1589 SETUP::
1590
1591 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1592
1593 TESTS:
1594
1595 This gives the same answer as the slow, default method implemented
1596 in :class:`MatrixEuclideanJordanAlgebra`::
1597
1598 sage: set_random_seed()
1599 sage: J = QuaternionHermitianEJA.random_instance()
1600 sage: x,y = J.random_elements(2)
1601 sage: Xe = x.natural_representation()
1602 sage: Ye = y.natural_representation()
1603 sage: X = QuaternionHermitianEJA.real_unembed(Xe)
1604 sage: Y = QuaternionHermitianEJA.real_unembed(Ye)
1605 sage: expected = (X*Y).trace().coefficient_tuple()[0]
1606 sage: actual = QuaternionHermitianEJA.natural_inner_product(Xe,Ye)
1607 sage: actual == expected
1608 True
1609
1610 """
1611 return RealMatrixEuclideanJordanAlgebra.natural_inner_product(X,Y)/4
1612
1613
1614 class QuaternionHermitianEJA(QuaternionMatrixEuclideanJordanAlgebra,
1615 KnownRankEJA):
1616 """
1617 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
1618 matrices, the usual symmetric Jordan product, and the
1619 real-part-of-trace inner product. It has dimension `2n^2 - n` over
1620 the reals.
1621
1622 SETUP::
1623
1624 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1625
1626 TESTS:
1627
1628 The dimension of this algebra is `2*n^2 - n`::
1629
1630 sage: set_random_seed()
1631 sage: n_max = QuaternionHermitianEJA._max_test_case_size()
1632 sage: n = ZZ.random_element(1, n_max)
1633 sage: J = QuaternionHermitianEJA(n)
1634 sage: J.dimension() == 2*(n^2) - n
1635 True
1636
1637 The Jordan multiplication is what we think it is::
1638
1639 sage: set_random_seed()
1640 sage: J = QuaternionHermitianEJA.random_instance()
1641 sage: x,y = J.random_elements(2)
1642 sage: actual = (x*y).natural_representation()
1643 sage: X = x.natural_representation()
1644 sage: Y = y.natural_representation()
1645 sage: expected = (X*Y + Y*X)/2
1646 sage: actual == expected
1647 True
1648 sage: J(expected) == x*y
1649 True
1650
1651 We can change the generator prefix::
1652
1653 sage: QuaternionHermitianEJA(2, prefix='a').gens()
1654 (a0, a1, a2, a3, a4, a5)
1655
1656 Our natural basis is normalized with respect to the natural inner
1657 product unless we specify otherwise::
1658
1659 sage: set_random_seed()
1660 sage: J = QuaternionHermitianEJA.random_instance()
1661 sage: all( b.norm() == 1 for b in J.gens() )
1662 True
1663
1664 Since our natural basis is normalized with respect to the natural
1665 inner product, and since we know that this algebra is an EJA, any
1666 left-multiplication operator's matrix will be symmetric because
1667 natural->EJA basis representation is an isometry and within the EJA
1668 the operator is self-adjoint by the Jordan axiom::
1669
1670 sage: set_random_seed()
1671 sage: x = QuaternionHermitianEJA.random_instance().random_element()
1672 sage: x.operator().matrix().is_symmetric()
1673 True
1674
1675 """
1676 @classmethod
1677 def _denormalized_basis(cls, n, field):
1678 """
1679 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
1680
1681 Why do we embed these? Basically, because all of numerical
1682 linear algebra assumes that you're working with vectors consisting
1683 of `n` entries from a field and scalars from the same field. There's
1684 no way to tell SageMath that (for example) the vectors contain
1685 complex numbers, while the scalar field is real.
1686
1687 SETUP::
1688
1689 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1690
1691 TESTS::
1692
1693 sage: set_random_seed()
1694 sage: n = ZZ.random_element(1,5)
1695 sage: B = QuaternionHermitianEJA._denormalized_basis(n,QQ)
1696 sage: all( M.is_symmetric() for M in B )
1697 True
1698
1699 """
1700 Q = QuaternionAlgebra(QQ,-1,-1)
1701 I,J,K = Q.gens()
1702
1703 # This is like the symmetric case, but we need to be careful:
1704 #
1705 # * We want conjugate-symmetry, not just symmetry.
1706 # * The diagonal will (as a result) be real.
1707 #
1708 S = []
1709 for i in xrange(n):
1710 for j in xrange(i+1):
1711 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
1712 if i == j:
1713 Sij = cls.real_embed(Eij)
1714 S.append(Sij)
1715 else:
1716 # The second, third, and fourth ones have a minus
1717 # because they're conjugated.
1718 Sij_real = cls.real_embed(Eij + Eij.transpose())
1719 S.append(Sij_real)
1720 Sij_I = cls.real_embed(I*Eij - I*Eij.transpose())
1721 S.append(Sij_I)
1722 Sij_J = cls.real_embed(J*Eij - J*Eij.transpose())
1723 S.append(Sij_J)
1724 Sij_K = cls.real_embed(K*Eij - K*Eij.transpose())
1725 S.append(Sij_K)
1726
1727 # Since we embedded these, we can drop back to the "field" that we
1728 # started with instead of the quaternion algebra "Q".
1729 return ( s.change_ring(field) for s in S )
1730
1731
1732 def __init__(self, n, field=QQ, **kwargs):
1733 basis = self._denormalized_basis(n,field)
1734 super(QuaternionHermitianEJA,self).__init__(field, basis, n, **kwargs)
1735
1736
1737 class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra, KnownRankEJA):
1738 """
1739 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
1740 with the usual inner product and jordan product ``x*y =
1741 (<x_bar,y_bar>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
1742 the reals.
1743
1744 SETUP::
1745
1746 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1747
1748 EXAMPLES:
1749
1750 This multiplication table can be verified by hand::
1751
1752 sage: J = JordanSpinEJA(4)
1753 sage: e0,e1,e2,e3 = J.gens()
1754 sage: e0*e0
1755 e0
1756 sage: e0*e1
1757 e1
1758 sage: e0*e2
1759 e2
1760 sage: e0*e3
1761 e3
1762 sage: e1*e2
1763 0
1764 sage: e1*e3
1765 0
1766 sage: e2*e3
1767 0
1768
1769 We can change the generator prefix::
1770
1771 sage: JordanSpinEJA(2, prefix='B').gens()
1772 (B0, B1)
1773
1774 """
1775 def __init__(self, n, field=QQ, **kwargs):
1776 V = VectorSpace(field, n)
1777 mult_table = [[V.zero() for j in xrange(n)] for i in xrange(n)]
1778 for i in xrange(n):
1779 for j in xrange(n):
1780 x = V.gen(i)
1781 y = V.gen(j)
1782 x0 = x[0]
1783 xbar = x[1:]
1784 y0 = y[0]
1785 ybar = y[1:]
1786 # z = x*y
1787 z0 = x.inner_product(y)
1788 zbar = y0*xbar + x0*ybar
1789 z = V([z0] + zbar.list())
1790 mult_table[i][j] = z
1791
1792 # The rank of the spin algebra is two, unless we're in a
1793 # one-dimensional ambient space (because the rank is bounded by
1794 # the ambient dimension).
1795 fdeja = super(JordanSpinEJA, self)
1796 return fdeja.__init__(field, mult_table, rank=min(n,2), **kwargs)
1797
1798 def inner_product(self, x, y):
1799 """
1800 Faster to reimplement than to use natural representations.
1801
1802 SETUP::
1803
1804 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1805
1806 TESTS:
1807
1808 Ensure that this is the usual inner product for the algebras
1809 over `R^n`::
1810
1811 sage: set_random_seed()
1812 sage: J = JordanSpinEJA.random_instance()
1813 sage: x,y = J.random_elements(2)
1814 sage: X = x.natural_representation()
1815 sage: Y = y.natural_representation()
1816 sage: x.inner_product(y) == J.natural_inner_product(X,Y)
1817 True
1818
1819 """
1820 return x.to_vector().inner_product(y.to_vector())