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