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