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