]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/eja_algebra.py
eja: add "normalize" argument to matrix algebra constructors.
[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, normalize):
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, False)
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 if normalize:
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, normalize):
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, False)
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 # Since we embedded these, we can drop back to the "field" that we
888 # started with instead of the complex extension "F".
889 S = [ s.change_ring(field) for s in S ]
890 if normalize:
891 S = [ s / _complex_hermitian_matrix_ip(s,s).sqrt() for s in S ]
892
893 return tuple(S)
894
895
896
897 def _quaternion_hermitian_basis(n, field, normalize):
898 """
899 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
900
901 Why do we embed these? Basically, because all of numerical linear
902 algebra assumes that you're working with vectors consisting of `n`
903 entries from a field and scalars from the same field. There's no way
904 to tell SageMath that (for example) the vectors contain complex
905 numbers, while the scalar field is real.
906
907 SETUP::
908
909 sage: from mjo.eja.eja_algebra import _quaternion_hermitian_basis
910
911 TESTS::
912
913 sage: set_random_seed()
914 sage: n = ZZ.random_element(1,5)
915 sage: B = _quaternion_hermitian_basis(n, QQ, False)
916 sage: all( M.is_symmetric() for M in B )
917 True
918
919 """
920 Q = QuaternionAlgebra(QQ,-1,-1)
921 I,J,K = Q.gens()
922
923 # This is like the symmetric case, but we need to be careful:
924 #
925 # * We want conjugate-symmetry, not just symmetry.
926 # * The diagonal will (as a result) be real.
927 #
928 S = []
929 for i in xrange(n):
930 for j in xrange(i+1):
931 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
932 if i == j:
933 Sij = _embed_quaternion_matrix(Eij)
934 S.append(Sij)
935 else:
936 # Beware, orthogonal but not normalized! The second,
937 # third, and fourth ones have a minus because they're
938 # conjugated.
939 Sij_real = _embed_quaternion_matrix(Eij + Eij.transpose())
940 S.append(Sij_real)
941 Sij_I = _embed_quaternion_matrix(I*Eij - I*Eij.transpose())
942 S.append(Sij_I)
943 Sij_J = _embed_quaternion_matrix(J*Eij - J*Eij.transpose())
944 S.append(Sij_J)
945 Sij_K = _embed_quaternion_matrix(K*Eij - K*Eij.transpose())
946 S.append(Sij_K)
947 return tuple(S)
948
949
950
951 def _multiplication_table_from_matrix_basis(basis):
952 """
953 At least three of the five simple Euclidean Jordan algebras have the
954 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
955 multiplication on the right is matrix multiplication. Given a basis
956 for the underlying matrix space, this function returns a
957 multiplication table (obtained by looping through the basis
958 elements) for an algebra of those matrices.
959 """
960 # In S^2, for example, we nominally have four coordinates even
961 # though the space is of dimension three only. The vector space V
962 # is supposed to hold the entire long vector, and the subspace W
963 # of V will be spanned by the vectors that arise from symmetric
964 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
965 field = basis[0].base_ring()
966 dimension = basis[0].nrows()
967
968 V = VectorSpace(field, dimension**2)
969 W = V.span_of_basis( _mat2vec(s) for s in basis )
970 n = len(basis)
971 mult_table = [[W.zero() for j in range(n)] for i in range(n)]
972 for i in range(n):
973 for j in range(n):
974 mat_entry = (basis[i]*basis[j] + basis[j]*basis[i])/2
975 mult_table[i][j] = W.coordinate_vector(_mat2vec(mat_entry))
976
977 return mult_table
978
979
980 def _embed_complex_matrix(M):
981 """
982 Embed the n-by-n complex matrix ``M`` into the space of real
983 matrices of size 2n-by-2n via the map the sends each entry `z = a +
984 bi` to the block matrix ``[[a,b],[-b,a]]``.
985
986 SETUP::
987
988 sage: from mjo.eja.eja_algebra import _embed_complex_matrix
989
990 EXAMPLES::
991
992 sage: F = QuadraticField(-1, 'i')
993 sage: x1 = F(4 - 2*i)
994 sage: x2 = F(1 + 2*i)
995 sage: x3 = F(-i)
996 sage: x4 = F(6)
997 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
998 sage: _embed_complex_matrix(M)
999 [ 4 -2| 1 2]
1000 [ 2 4|-2 1]
1001 [-----+-----]
1002 [ 0 -1| 6 0]
1003 [ 1 0| 0 6]
1004
1005 TESTS:
1006
1007 Embedding is a homomorphism (isomorphism, in fact)::
1008
1009 sage: set_random_seed()
1010 sage: n = ZZ.random_element(5)
1011 sage: F = QuadraticField(-1, 'i')
1012 sage: X = random_matrix(F, n)
1013 sage: Y = random_matrix(F, n)
1014 sage: actual = _embed_complex_matrix(X) * _embed_complex_matrix(Y)
1015 sage: expected = _embed_complex_matrix(X*Y)
1016 sage: actual == expected
1017 True
1018
1019 """
1020 n = M.nrows()
1021 if M.ncols() != n:
1022 raise ValueError("the matrix 'M' must be square")
1023 field = M.base_ring()
1024 blocks = []
1025 for z in M.list():
1026 a = z.vector()[0] # real part, I guess
1027 b = z.vector()[1] # imag part, I guess
1028 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
1029
1030 # We can drop the imaginaries here.
1031 return matrix.block(field.base_ring(), n, blocks)
1032
1033
1034 def _unembed_complex_matrix(M):
1035 """
1036 The inverse of _embed_complex_matrix().
1037
1038 SETUP::
1039
1040 sage: from mjo.eja.eja_algebra import (_embed_complex_matrix,
1041 ....: _unembed_complex_matrix)
1042
1043 EXAMPLES::
1044
1045 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
1046 ....: [-2, 1, -4, 3],
1047 ....: [ 9, 10, 11, 12],
1048 ....: [-10, 9, -12, 11] ])
1049 sage: _unembed_complex_matrix(A)
1050 [ 2*i + 1 4*i + 3]
1051 [ 10*i + 9 12*i + 11]
1052
1053 TESTS:
1054
1055 Unembedding is the inverse of embedding::
1056
1057 sage: set_random_seed()
1058 sage: F = QuadraticField(-1, 'i')
1059 sage: M = random_matrix(F, 3)
1060 sage: _unembed_complex_matrix(_embed_complex_matrix(M)) == M
1061 True
1062
1063 """
1064 n = ZZ(M.nrows())
1065 if M.ncols() != n:
1066 raise ValueError("the matrix 'M' must be square")
1067 if not n.mod(2).is_zero():
1068 raise ValueError("the matrix 'M' must be a complex embedding")
1069
1070 field = M.base_ring() # This should already have sqrt2
1071 R = PolynomialRing(field, 'z')
1072 z = R.gen()
1073 F = NumberField(z**2 + 1,'i', embedding=CLF(-1).sqrt())
1074 i = F.gen()
1075
1076 # Go top-left to bottom-right (reading order), converting every
1077 # 2-by-2 block we see to a single complex element.
1078 elements = []
1079 for k in xrange(n/2):
1080 for j in xrange(n/2):
1081 submat = M[2*k:2*k+2,2*j:2*j+2]
1082 if submat[0,0] != submat[1,1]:
1083 raise ValueError('bad on-diagonal submatrix')
1084 if submat[0,1] != -submat[1,0]:
1085 raise ValueError('bad off-diagonal submatrix')
1086 z = submat[0,0] + submat[0,1]*i
1087 elements.append(z)
1088
1089 return matrix(F, n/2, elements)
1090
1091
1092 def _embed_quaternion_matrix(M):
1093 """
1094 Embed the n-by-n quaternion matrix ``M`` into the space of real
1095 matrices of size 4n-by-4n by first sending each quaternion entry
1096 `z = a + bi + cj + dk` to the block-complex matrix
1097 ``[[a + bi, c+di],[-c + di, a-bi]]`, and then embedding those into
1098 a real matrix.
1099
1100 SETUP::
1101
1102 sage: from mjo.eja.eja_algebra import _embed_quaternion_matrix
1103
1104 EXAMPLES::
1105
1106 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1107 sage: i,j,k = Q.gens()
1108 sage: x = 1 + 2*i + 3*j + 4*k
1109 sage: M = matrix(Q, 1, [[x]])
1110 sage: _embed_quaternion_matrix(M)
1111 [ 1 2 3 4]
1112 [-2 1 -4 3]
1113 [-3 4 1 -2]
1114 [-4 -3 2 1]
1115
1116 Embedding is a homomorphism (isomorphism, in fact)::
1117
1118 sage: set_random_seed()
1119 sage: n = ZZ.random_element(5)
1120 sage: Q = QuaternionAlgebra(QQ,-1,-1)
1121 sage: X = random_matrix(Q, n)
1122 sage: Y = random_matrix(Q, n)
1123 sage: actual = _embed_quaternion_matrix(X)*_embed_quaternion_matrix(Y)
1124 sage: expected = _embed_quaternion_matrix(X*Y)
1125 sage: actual == expected
1126 True
1127
1128 """
1129 quaternions = M.base_ring()
1130 n = M.nrows()
1131 if M.ncols() != n:
1132 raise ValueError("the matrix 'M' must be square")
1133
1134 F = QuadraticField(-1, 'i')
1135 i = F.gen()
1136
1137 blocks = []
1138 for z in M.list():
1139 t = z.coefficient_tuple()
1140 a = t[0]
1141 b = t[1]
1142 c = t[2]
1143 d = t[3]
1144 cplx_matrix = matrix(F, 2, [[ a + b*i, c + d*i],
1145 [-c + d*i, a - b*i]])
1146 blocks.append(_embed_complex_matrix(cplx_matrix))
1147
1148 # We should have real entries by now, so use the realest field
1149 # we've got for the return value.
1150 return matrix.block(quaternions.base_ring(), n, blocks)
1151
1152
1153 def _unembed_quaternion_matrix(M):
1154 """
1155 The inverse of _embed_quaternion_matrix().
1156
1157 SETUP::
1158
1159 sage: from mjo.eja.eja_algebra import (_embed_quaternion_matrix,
1160 ....: _unembed_quaternion_matrix)
1161
1162 EXAMPLES::
1163
1164 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
1165 ....: [-2, 1, -4, 3],
1166 ....: [-3, 4, 1, -2],
1167 ....: [-4, -3, 2, 1]])
1168 sage: _unembed_quaternion_matrix(M)
1169 [1 + 2*i + 3*j + 4*k]
1170
1171 TESTS:
1172
1173 Unembedding is the inverse of embedding::
1174
1175 sage: set_random_seed()
1176 sage: Q = QuaternionAlgebra(QQ, -1, -1)
1177 sage: M = random_matrix(Q, 3)
1178 sage: _unembed_quaternion_matrix(_embed_quaternion_matrix(M)) == M
1179 True
1180
1181 """
1182 n = ZZ(M.nrows())
1183 if M.ncols() != n:
1184 raise ValueError("the matrix 'M' must be square")
1185 if not n.mod(4).is_zero():
1186 raise ValueError("the matrix 'M' must be a complex embedding")
1187
1188 Q = QuaternionAlgebra(QQ,-1,-1)
1189 i,j,k = Q.gens()
1190
1191 # Go top-left to bottom-right (reading order), converting every
1192 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
1193 # quaternion block.
1194 elements = []
1195 for l in xrange(n/4):
1196 for m in xrange(n/4):
1197 submat = _unembed_complex_matrix(M[4*l:4*l+4,4*m:4*m+4])
1198 if submat[0,0] != submat[1,1].conjugate():
1199 raise ValueError('bad on-diagonal submatrix')
1200 if submat[0,1] != -submat[1,0].conjugate():
1201 raise ValueError('bad off-diagonal submatrix')
1202 z = submat[0,0].real() + submat[0,0].imag()*i
1203 z += submat[0,1].real()*j + submat[0,1].imag()*k
1204 elements.append(z)
1205
1206 return matrix(Q, n/4, elements)
1207
1208
1209 # The usual inner product on R^n.
1210 def _usual_ip(x,y):
1211 return x.to_vector().inner_product(y.to_vector())
1212
1213 # The inner product used for the real symmetric simple EJA.
1214 # We keep it as a separate function because e.g. the complex
1215 # algebra uses the same inner product, except divided by 2.
1216 def _matrix_ip(X,Y):
1217 X_mat = X.natural_representation()
1218 Y_mat = Y.natural_representation()
1219 return (X_mat*Y_mat).trace()
1220
1221 def _real_symmetric_matrix_ip(X,Y):
1222 return (X*Y).trace()
1223
1224 def _complex_hermitian_matrix_ip(X,Y):
1225 # This takes EMBEDDED matrices.
1226 Xu = _unembed_complex_matrix(X)
1227 Yu = _unembed_complex_matrix(Y)
1228 # The trace need not be real; consider Xu = (i*I) and Yu = I.
1229 return ((Xu*Yu).trace()).vector()[0] # real part, I guess
1230
1231 class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra):
1232 """
1233 The rank-n simple EJA consisting of real symmetric n-by-n
1234 matrices, the usual symmetric Jordan product, and the trace inner
1235 product. It has dimension `(n^2 + n)/2` over the reals.
1236
1237 SETUP::
1238
1239 sage: from mjo.eja.eja_algebra import RealSymmetricEJA
1240
1241 EXAMPLES::
1242
1243 sage: J = RealSymmetricEJA(2)
1244 sage: e0, e1, e2 = J.gens()
1245 sage: e0*e0
1246 e0
1247 sage: e1*e1
1248 1/2*e0 + 1/2*e2
1249 sage: e2*e2
1250 e2
1251
1252 TESTS:
1253
1254 The dimension of this algebra is `(n^2 + n) / 2`::
1255
1256 sage: set_random_seed()
1257 sage: n = ZZ.random_element(1,5)
1258 sage: J = RealSymmetricEJA(n)
1259 sage: J.dimension() == (n^2 + n)/2
1260 True
1261
1262 The Jordan multiplication is what we think it is::
1263
1264 sage: set_random_seed()
1265 sage: n = ZZ.random_element(1,5)
1266 sage: J = RealSymmetricEJA(n)
1267 sage: x = J.random_element()
1268 sage: y = J.random_element()
1269 sage: actual = (x*y).natural_representation()
1270 sage: X = x.natural_representation()
1271 sage: Y = y.natural_representation()
1272 sage: expected = (X*Y + Y*X)/2
1273 sage: actual == expected
1274 True
1275 sage: J(expected) == x*y
1276 True
1277
1278 We can change the generator prefix::
1279
1280 sage: RealSymmetricEJA(3, prefix='q').gens()
1281 (q0, q1, q2, q3, q4, q5)
1282
1283 Our inner product satisfies the Jordan axiom::
1284
1285 sage: set_random_seed()
1286 sage: n = ZZ.random_element(1,5)
1287 sage: J = RealSymmetricEJA(n)
1288 sage: x = J.random_element()
1289 sage: y = J.random_element()
1290 sage: z = J.random_element()
1291 sage: (x*y).inner_product(z) == y.inner_product(x*z)
1292 True
1293
1294 Our basis is normalized with respect to the natural inner product::
1295
1296 sage: set_random_seed()
1297 sage: n = ZZ.random_element(1,5)
1298 sage: J = RealSymmetricEJA(n)
1299 sage: all( b.norm() == 1 for b in J.gens() )
1300 True
1301
1302 Left-multiplication operators are symmetric because they satisfy
1303 the Jordan axiom::
1304
1305 sage: set_random_seed()
1306 sage: n = ZZ.random_element(1,5)
1307 sage: x = RealSymmetricEJA(n).random_element()
1308 sage: x.operator().matrix().is_symmetric()
1309 True
1310
1311 """
1312 def __init__(self, n, field=QQ, normalize_basis=True, **kwargs):
1313 if n > 1:
1314 # We'll need sqrt(2) to normalize the basis, and this
1315 # winds up in the multiplication table, so the whole
1316 # algebra needs to be over the field extension.
1317 R = PolynomialRing(field, 'z')
1318 z = R.gen()
1319 p = z**2 - 2
1320 if p.is_irreducible():
1321 field = NumberField(p, 'sqrt2', embedding=RLF(2).sqrt())
1322
1323 S = _real_symmetric_basis(n, field, normalize_basis)
1324 Qs = _multiplication_table_from_matrix_basis(S)
1325
1326 fdeja = super(RealSymmetricEJA, self)
1327 return fdeja.__init__(field,
1328 Qs,
1329 rank=n,
1330 natural_basis=S,
1331 **kwargs)
1332
1333 def inner_product(self, x, y):
1334 X = x.natural_representation()
1335 Y = y.natural_representation()
1336 return _real_symmetric_matrix_ip(X,Y)
1337
1338
1339 class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
1340 """
1341 The rank-n simple EJA consisting of complex Hermitian n-by-n
1342 matrices over the real numbers, the usual symmetric Jordan product,
1343 and the real-part-of-trace inner product. It has dimension `n^2` over
1344 the reals.
1345
1346 SETUP::
1347
1348 sage: from mjo.eja.eja_algebra import ComplexHermitianEJA
1349
1350 TESTS:
1351
1352 The dimension of this algebra is `n^2`::
1353
1354 sage: set_random_seed()
1355 sage: n = ZZ.random_element(1,5)
1356 sage: J = ComplexHermitianEJA(n)
1357 sage: J.dimension() == n^2
1358 True
1359
1360 The Jordan multiplication is what we think it is::
1361
1362 sage: set_random_seed()
1363 sage: n = ZZ.random_element(1,5)
1364 sage: J = ComplexHermitianEJA(n)
1365 sage: x = J.random_element()
1366 sage: y = J.random_element()
1367 sage: actual = (x*y).natural_representation()
1368 sage: X = x.natural_representation()
1369 sage: Y = y.natural_representation()
1370 sage: expected = (X*Y + Y*X)/2
1371 sage: actual == expected
1372 True
1373 sage: J(expected) == x*y
1374 True
1375
1376 We can change the generator prefix::
1377
1378 sage: ComplexHermitianEJA(2, prefix='z').gens()
1379 (z0, z1, z2, z3)
1380
1381 Our inner product satisfies the Jordan axiom::
1382
1383 sage: set_random_seed()
1384 sage: n = ZZ.random_element(1,5)
1385 sage: J = ComplexHermitianEJA(n)
1386 sage: x = J.random_element()
1387 sage: y = J.random_element()
1388 sage: z = J.random_element()
1389 sage: (x*y).inner_product(z) == y.inner_product(x*z)
1390 True
1391
1392 Our basis is normalized with respect to the natural inner product::
1393
1394 sage: set_random_seed()
1395 sage: n = ZZ.random_element(1,4)
1396 sage: J = ComplexHermitianEJA(n)
1397 sage: all( b.norm() == 1 for b in J.gens() )
1398 True
1399
1400 Left-multiplication operators are symmetric because they satisfy
1401 the Jordan axiom::
1402
1403 sage: set_random_seed()
1404 sage: n = ZZ.random_element(1,5)
1405 sage: x = ComplexHermitianEJA(n).random_element()
1406 sage: x.operator().matrix().is_symmetric()
1407 True
1408
1409 """
1410 def __init__(self, n, field=QQ, normalize_basis=True, **kwargs):
1411 if n > 1:
1412 # We'll need sqrt(2) to normalize the basis, and this
1413 # winds up in the multiplication table, so the whole
1414 # algebra needs to be over the field extension.
1415 R = PolynomialRing(field, 'z')
1416 z = R.gen()
1417 p = z**2 - 2
1418 if p.is_irreducible():
1419 field = NumberField(p, 'sqrt2', embedding=RLF(2).sqrt())
1420
1421 S = _complex_hermitian_basis(n, field, normalize_basis)
1422 Qs = _multiplication_table_from_matrix_basis(S)
1423
1424 fdeja = super(ComplexHermitianEJA, self)
1425 return fdeja.__init__(field,
1426 Qs,
1427 rank=n,
1428 natural_basis=S,
1429 **kwargs)
1430
1431
1432 def inner_product(self, x, y):
1433 X = x.natural_representation()
1434 Y = y.natural_representation()
1435 return _complex_hermitian_matrix_ip(X,Y)
1436
1437
1438 class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
1439 """
1440 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
1441 matrices, the usual symmetric Jordan product, and the
1442 real-part-of-trace inner product. It has dimension `2n^2 - n` over
1443 the reals.
1444
1445 SETUP::
1446
1447 sage: from mjo.eja.eja_algebra import QuaternionHermitianEJA
1448
1449 TESTS:
1450
1451 The dimension of this algebra is `n^2`::
1452
1453 sage: set_random_seed()
1454 sage: n = ZZ.random_element(1,5)
1455 sage: J = QuaternionHermitianEJA(n)
1456 sage: J.dimension() == 2*(n^2) - n
1457 True
1458
1459 The Jordan multiplication is what we think it is::
1460
1461 sage: set_random_seed()
1462 sage: n = ZZ.random_element(1,5)
1463 sage: J = QuaternionHermitianEJA(n)
1464 sage: x = J.random_element()
1465 sage: y = J.random_element()
1466 sage: actual = (x*y).natural_representation()
1467 sage: X = x.natural_representation()
1468 sage: Y = y.natural_representation()
1469 sage: expected = (X*Y + Y*X)/2
1470 sage: actual == expected
1471 True
1472 sage: J(expected) == x*y
1473 True
1474
1475 We can change the generator prefix::
1476
1477 sage: QuaternionHermitianEJA(2, prefix='a').gens()
1478 (a0, a1, a2, a3, a4, a5)
1479
1480 Our inner product satisfies the Jordan axiom::
1481
1482 sage: set_random_seed()
1483 sage: n = ZZ.random_element(1,5)
1484 sage: J = QuaternionHermitianEJA(n)
1485 sage: x = J.random_element()
1486 sage: y = J.random_element()
1487 sage: z = J.random_element()
1488 sage: (x*y).inner_product(z) == y.inner_product(x*z)
1489 True
1490
1491 """
1492 def __init__(self, n, field=QQ, normalize_basis=True, **kwargs):
1493 S = _quaternion_hermitian_basis(n, field, normalize_basis)
1494 Qs = _multiplication_table_from_matrix_basis(S)
1495
1496 fdeja = super(QuaternionHermitianEJA, self)
1497 return fdeja.__init__(field,
1498 Qs,
1499 rank=n,
1500 natural_basis=S,
1501 **kwargs)
1502
1503 def inner_product(self, x, y):
1504 # Since a+bi+cj+dk on the diagonal is represented as
1505 #
1506 # a + bi +cj + dk = [ a b c d]
1507 # [ -b a -d c]
1508 # [ -c d a -b]
1509 # [ -d -c b a],
1510 #
1511 # we'll quadruple-count the "a" entries if we take the trace of
1512 # the embedding.
1513 return _matrix_ip(x,y)/4
1514
1515
1516 class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra):
1517 """
1518 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
1519 with the usual inner product and jordan product ``x*y =
1520 (<x_bar,y_bar>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
1521 the reals.
1522
1523 SETUP::
1524
1525 sage: from mjo.eja.eja_algebra import JordanSpinEJA
1526
1527 EXAMPLES:
1528
1529 This multiplication table can be verified by hand::
1530
1531 sage: J = JordanSpinEJA(4)
1532 sage: e0,e1,e2,e3 = J.gens()
1533 sage: e0*e0
1534 e0
1535 sage: e0*e1
1536 e1
1537 sage: e0*e2
1538 e2
1539 sage: e0*e3
1540 e3
1541 sage: e1*e2
1542 0
1543 sage: e1*e3
1544 0
1545 sage: e2*e3
1546 0
1547
1548 We can change the generator prefix::
1549
1550 sage: JordanSpinEJA(2, prefix='B').gens()
1551 (B0, B1)
1552
1553 Our inner product satisfies the Jordan axiom::
1554
1555 sage: set_random_seed()
1556 sage: n = ZZ.random_element(1,5)
1557 sage: J = JordanSpinEJA(n)
1558 sage: x = J.random_element()
1559 sage: y = J.random_element()
1560 sage: z = J.random_element()
1561 sage: (x*y).inner_product(z) == y.inner_product(x*z)
1562 True
1563
1564 """
1565 def __init__(self, n, field=QQ, **kwargs):
1566 V = VectorSpace(field, n)
1567 mult_table = [[V.zero() for j in range(n)] for i in range(n)]
1568 for i in range(n):
1569 for j in range(n):
1570 x = V.gen(i)
1571 y = V.gen(j)
1572 x0 = x[0]
1573 xbar = x[1:]
1574 y0 = y[0]
1575 ybar = y[1:]
1576 # z = x*y
1577 z0 = x.inner_product(y)
1578 zbar = y0*xbar + x0*ybar
1579 z = V([z0] + zbar.list())
1580 mult_table[i][j] = z
1581
1582 # The rank of the spin algebra is two, unless we're in a
1583 # one-dimensional ambient space (because the rank is bounded by
1584 # the ambient dimension).
1585 fdeja = super(JordanSpinEJA, self)
1586 return fdeja.__init__(field, mult_table, rank=min(n,2), **kwargs)
1587
1588 def inner_product(self, x, y):
1589 return _usual_ip(x,y)