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