]> gitweb.michael.orlitzky.com - sage.d.git/blob - mjo/eja/euclidean_jordan_algebra.py
eja: fix an erroneous test case.
[sage.d.git] / mjo / eja / euclidean_jordan_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.categories.finite_dimensional_algebras_with_basis import FiniteDimensionalAlgebrasWithBasis
9 from sage.categories.morphism import SetMorphism
10 from sage.structure.element import is_Matrix
11 from sage.structure.category_object import normalize_names
12
13 from sage.algebras.finite_dimensional_algebras.finite_dimensional_algebra import FiniteDimensionalAlgebra
14 from sage.algebras.finite_dimensional_algebras.finite_dimensional_algebra_element import FiniteDimensionalAlgebraElement
15 from sage.algebras.finite_dimensional_algebras.finite_dimensional_algebra_morphism import FiniteDimensionalAlgebraMorphism, FiniteDimensionalAlgebraHomset
16
17
18 class FiniteDimensionalEuclideanJordanAlgebraHomset(FiniteDimensionalAlgebraHomset):
19
20 def has_coerce_map_from(self, S):
21 """
22 EXAMPLES::
23
24 sage: J = RealSymmetricEJA(2)
25 sage: H = J.Hom(J)
26 sage: H.has_coerce_map_from(QQ)
27 True
28
29 """
30 try:
31 # The Homset classes override has_coerce_map_from() with
32 # something that crashes when it's given e.g. QQ.
33 if S.is_subring(self.codomain().base_ring()):
34 return True
35 except:
36 pclass = super(FiniteDimensionalEuclideanJordanAlgebraHomset, self)
37 return pclass.has_coerce_map_from(S)
38
39
40 def _coerce_map_from_(self, S):
41 """
42 EXAMPLES::
43
44 sage: J = RealSymmetricEJA(2)
45 sage: H = J.Hom(J)
46 sage: H.coerce(2)
47 Morphism from Euclidean Jordan algebra of degree 3 over Rational
48 Field to Euclidean Jordan algebra of degree 3 over Rational Field
49 given by matrix
50 [2 0 0]
51 [0 2 0]
52 [0 0 2]
53
54 """
55 C = self.codomain()
56 R = C.base_ring()
57 if S.is_subring(R):
58 h = S.hom(self.codomain())
59 return SetMorphism(Hom(S,C), lambda x: h(x).operator())
60
61
62 def __call__(self, x):
63 """
64 EXAMPLES::
65
66 sage: J = RealSymmetricEJA(2)
67 sage: H = J.Hom(J)
68 sage: H(2)
69 Morphism from Euclidean Jordan algebra of degree 3 over Rational
70 Field to Euclidean Jordan algebra of degree 3 over Rational Field
71 given by matrix
72 [2 0 0]
73 [0 2 0]
74 [0 0 2]
75
76 """
77 if x in self.base_ring():
78 cols = self.domain().dimension()
79 rows = self.codomain().dimension()
80 x = x*identity_matrix(self.codomain().base_ring(), rows, cols)
81 return FiniteDimensionalEuclideanJordanAlgebraMorphism(self, x)
82
83
84 def one(self):
85 """
86 Return the identity morphism, but as a member of the right
87 space (so that we can add it, multiply it, etc.)
88 """
89 cols = self.domain().dimension()
90 rows = self.codomain().dimension()
91 mat = identity_matrix(self.base_ring(), rows, cols)
92 return FiniteDimensionalEuclideanJordanAlgebraMorphism(self, mat)
93
94
95
96 class FiniteDimensionalEuclideanJordanAlgebraMorphism(FiniteDimensionalAlgebraMorphism):
97 """
98 A linear map between two finite-dimensional EJAs.
99
100 This is a very thin wrapper around FiniteDimensionalAlgebraMorphism
101 that does only a few things:
102
103 1. Avoids the ``unitary`` and ``check`` arguments to the constructor
104 that will always be ``False``. This is necessary because these
105 are homomorphisms with respect to ADDITION, but the SageMath
106 machinery wants to check that they're homomorphisms with respect
107 to (Jordan) MULTIPLICATION. That obviously doesn't work.
108
109 2. Inputs and outputs the underlying matrix with respect to COLUMN
110 vectors, unlike the parent class.
111
112 3. Allows us to add, subtract, negate, multiply (compose), and
113 invert morphisms in the obvious way.
114
115 If this seems a bit heavyweight, it is. I would have been happy to
116 use a the ring morphism that underlies the finite-dimensional
117 algebra morphism, but they don't seem to be callable on elements of
118 our EJA, and you can't add/multiply/etc. them.
119 """
120 def _add_(self, other):
121 """
122 Add two EJA morphisms in the obvious way.
123
124 EXAMPLES::
125
126 sage: J = RealSymmetricEJA(3)
127 sage: x = J.zero()
128 sage: y = J.one()
129 sage: x.operator() + y.operator()
130 Morphism from Euclidean Jordan algebra of degree 6 over Rational
131 Field to Euclidean Jordan algebra of degree 6 over Rational Field
132 given by matrix
133 [1 0 0 0 0 0]
134 [0 1 0 0 0 0]
135 [0 0 1 0 0 0]
136 [0 0 0 1 0 0]
137 [0 0 0 0 1 0]
138 [0 0 0 0 0 1]
139
140 TESTS::
141
142 sage: set_random_seed()
143 sage: J = random_eja()
144 sage: x = J.random_element()
145 sage: y = J.random_element()
146 sage: (x.operator() + y.operator()) in J.Hom(J)
147 True
148
149 """
150 P = self.parent()
151 if not other in P:
152 raise ValueError("summands must live in the same space")
153
154 return FiniteDimensionalEuclideanJordanAlgebraMorphism(
155 P,
156 self.matrix() + other.matrix() )
157
158
159 def __init__(self, parent, f):
160 FiniteDimensionalAlgebraMorphism.__init__(self,
161 parent,
162 f.transpose(),
163 unitary=False,
164 check=False)
165
166
167 def __invert__(self):
168 """
169 EXAMPLES::
170
171 sage: J = RealSymmetricEJA(2)
172 sage: x = J.linear_combination(zip(range(len(J.gens())), J.gens()))
173 sage: x.is_invertible()
174 True
175 sage: ~x.operator()
176 Morphism from Euclidean Jordan algebra of degree 3 over Rational
177 Field to Euclidean Jordan algebra of degree 3 over Rational Field
178 given by matrix
179 [-3/2 2 -1/2]
180 [ 1 0 0]
181 [-1/2 0 1/2]
182 sage: x.operator_matrix().inverse()
183 [-3/2 2 -1/2]
184 [ 1 0 0]
185 [-1/2 0 1/2]
186
187 TESTS::
188
189 sage: set_random_seed()
190 sage: J = random_eja()
191 sage: x = J.random_element()
192 sage: not x.is_invertible() or (
193 ....: (~x.operator()).matrix() == x.operator_matrix().inverse() )
194 True
195
196 """
197 A = self.matrix()
198 if not A.is_invertible():
199 raise ValueError("morphism is not invertible")
200
201 P = self.parent()
202 return FiniteDimensionalEuclideanJordanAlgebraMorphism(self.parent(),
203 A.inverse())
204
205 def _lmul_(self, right):
206 """
207 Compose two EJA morphisms using multiplicative notation.
208
209 EXAMPLES::
210
211 sage: J = RealSymmetricEJA(2)
212 sage: x = J.zero()
213 sage: y = J.one()
214 sage: x.operator() * y.operator()
215 Morphism from Euclidean Jordan algebra of degree 3 over Rational
216 Field to Euclidean Jordan algebra of degree 3 over Rational Field
217 given by matrix
218 [0 0 0]
219 [0 0 0]
220 [0 0 0]
221
222 ::
223
224 sage: J = RealSymmetricEJA(2)
225 sage: x = J.linear_combination(zip(range(len(J.gens())), J.gens()))
226 sage: x.operator()
227 Morphism from Euclidean Jordan algebra of degree 3 over Rational
228 Field to Euclidean Jordan algebra of degree 3 over Rational Field
229 given by matrix
230 [ 0 1 0]
231 [1/2 1 1/2]
232 [ 0 1 2]
233 sage: 2*x.operator()
234 Morphism from Euclidean Jordan algebra of degree 3 over Rational
235 Field to Euclidean Jordan algebra of degree 3 over Rational Field
236 given by matrix
237 [0 2 0]
238 [1 2 1]
239 [0 2 4]
240 sage: x.operator()*2
241 Morphism from Euclidean Jordan algebra of degree 3 over Rational
242 Field to Euclidean Jordan algebra of degree 3 over Rational Field
243 given by matrix
244 [0 2 0]
245 [1 2 1]
246 [0 2 4]
247
248 TESTS::
249
250 sage: set_random_seed()
251 sage: J = random_eja()
252 sage: x = J.random_element()
253 sage: y = J.random_element()
254 sage: (x.operator() * y.operator()) in J.Hom(J)
255 True
256
257 """
258 try:
259 # I think the morphism classes break the coercion framework
260 # somewhere along the way, so we have to do this ourselves.
261 right = self.parent().coerce(right)
262 except:
263 pass
264
265 if not right.codomain() is self.domain():
266 raise ValueError("(co)domains must agree for composition")
267
268 return FiniteDimensionalEuclideanJordanAlgebraMorphism(
269 self.parent(),
270 self.matrix()*right.matrix() )
271
272 __mul__ = _lmul_
273
274
275 def __pow__(self, n):
276 """
277
278 TESTS::
279
280 sage: J = JordanSpinEJA(4)
281 sage: e0,e1,e2,e3 = J.gens()
282 sage: x = -5/2*e0 + 1/2*e2 + 20*e3
283 sage: Qx = x.quadratic_representation()
284 sage: Qx^0
285 Morphism from Euclidean Jordan algebra of degree 4 over Rational
286 Field to Euclidean Jordan algebra of degree 4 over Rational Field
287 given by matrix
288 [1 0 0 0]
289 [0 1 0 0]
290 [0 0 1 0]
291 [0 0 0 1]
292 sage: (x^0).quadratic_representation() == Qx^0
293 True
294
295 """
296 if n == 0:
297 # We get back the stupid identity morphism which doesn't
298 # live in the right space.
299 return self.parent().one()
300 elif n == 1:
301 return self
302 else:
303 return FiniteDimensionalAlgebraMorphism.__pow__(self,n)
304
305
306 def _neg_(self):
307 """
308 Negate this morphism.
309
310 EXAMPLES::
311
312 sage: J = RealSymmetricEJA(2)
313 sage: x = J.one()
314 sage: -x.operator()
315 Morphism from Euclidean Jordan algebra of degree 3 over Rational
316 Field to Euclidean Jordan algebra of degree 3 over Rational Field
317 given by matrix
318 [-1 0 0]
319 [ 0 -1 0]
320 [ 0 0 -1]
321
322 TESTS::
323
324 sage: set_random_seed()
325 sage: J = random_eja()
326 sage: x = J.random_element()
327 sage: -x.operator() in J.Hom(J)
328 True
329
330 """
331 return FiniteDimensionalEuclideanJordanAlgebraMorphism(
332 self.parent(),
333 -self.matrix() )
334
335
336 def _repr_(self):
337 """
338 We override only the representation that is shown to the user,
339 because we want the matrix to be with respect to COLUMN vectors.
340
341 TESTS:
342
343 Ensure that we see the transpose of the underlying matrix object:
344
345 sage: J = RealSymmetricEJA(3)
346 sage: x = J.linear_combination(zip(range(len(J.gens())), J.gens()))
347 sage: L = x.operator()
348 sage: L
349 Morphism from Euclidean Jordan algebra of degree 6 over Rational
350 Field to Euclidean Jordan algebra of degree 6 over Rational Field
351 given by matrix
352 [ 0 1 2 0 0 0]
353 [1/2 3/2 2 1/2 1 0]
354 [ 1 2 5/2 0 1/2 1]
355 [ 0 1 0 3 4 0]
356 [ 0 1 1/2 2 4 2]
357 [ 0 0 2 0 4 5]
358 sage: L._matrix
359 [ 0 1/2 1 0 0 0]
360 [ 1 3/2 2 1 1 0]
361 [ 2 2 5/2 0 1/2 2]
362 [ 0 1/2 0 3 2 0]
363 [ 0 1 1/2 4 4 4]
364 [ 0 0 1 0 2 5]
365
366 """
367 return "Morphism from {} to {} given by matrix\n{}".format(
368 self.domain(), self.codomain(), self.matrix())
369
370
371 def __sub__(self, other):
372 """
373 Subtract one morphism from another using addition and negation.
374
375 EXAMPLES::
376
377 sage: J = RealSymmetricEJA(2)
378 sage: L1 = J.one().operator()
379 sage: L1 - L1
380 Morphism from Euclidean Jordan algebra of degree 3 over Rational
381 Field to Euclidean Jordan algebra of degree 3 over Rational
382 Field given by matrix
383 [0 0 0]
384 [0 0 0]
385 [0 0 0]
386
387 TESTS::
388
389 sage: set_random_seed()
390 sage: J = random_eja()
391 sage: x = J.random_element()
392 sage: y = J.random_element()
393 sage: x.operator() - y.operator() in J.Hom(J)
394 True
395
396 """
397 return self + (-other)
398
399
400 def matrix(self):
401 """
402 Return the matrix of this morphism with respect to a left-action
403 on column vectors.
404 """
405 return FiniteDimensionalAlgebraMorphism.matrix(self).transpose()
406
407
408 class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra):
409 @staticmethod
410 def __classcall_private__(cls,
411 field,
412 mult_table,
413 names='e',
414 assume_associative=False,
415 category=None,
416 rank=None,
417 natural_basis=None):
418 n = len(mult_table)
419 mult_table = [b.base_extend(field) for b in mult_table]
420 for b in mult_table:
421 b.set_immutable()
422 if not (is_Matrix(b) and b.dimensions() == (n, n)):
423 raise ValueError("input is not a multiplication table")
424 mult_table = tuple(mult_table)
425
426 cat = FiniteDimensionalAlgebrasWithBasis(field)
427 cat.or_subcategory(category)
428 if assume_associative:
429 cat = cat.Associative()
430
431 names = normalize_names(n, names)
432
433 fda = super(FiniteDimensionalEuclideanJordanAlgebra, cls)
434 return fda.__classcall__(cls,
435 field,
436 mult_table,
437 assume_associative=assume_associative,
438 names=names,
439 category=cat,
440 rank=rank,
441 natural_basis=natural_basis)
442
443
444 def _Hom_(self, B, cat):
445 """
446 Construct a homset of ``self`` and ``B``.
447 """
448 return FiniteDimensionalEuclideanJordanAlgebraHomset(self,
449 B,
450 category=cat)
451
452
453 def __init__(self,
454 field,
455 mult_table,
456 names='e',
457 assume_associative=False,
458 category=None,
459 rank=None,
460 natural_basis=None):
461 """
462 EXAMPLES:
463
464 By definition, Jordan multiplication commutes::
465
466 sage: set_random_seed()
467 sage: J = random_eja()
468 sage: x = J.random_element()
469 sage: y = J.random_element()
470 sage: x*y == y*x
471 True
472
473 """
474 self._rank = rank
475 self._natural_basis = natural_basis
476 self._multiplication_table = mult_table
477 fda = super(FiniteDimensionalEuclideanJordanAlgebra, self)
478 fda.__init__(field,
479 mult_table,
480 names=names,
481 category=category)
482
483
484 def _repr_(self):
485 """
486 Return a string representation of ``self``.
487 """
488 fmt = "Euclidean Jordan algebra of degree {} over {}"
489 return fmt.format(self.degree(), self.base_ring())
490
491
492 def _a_regular_element(self):
493 """
494 Guess a regular element. Needed to compute the basis for our
495 characteristic polynomial coefficients.
496 """
497 gs = self.gens()
498 z = self.sum( (i+1)*gs[i] for i in range(len(gs)) )
499 if not z.is_regular():
500 raise ValueError("don't know a regular element")
501 return z
502
503
504 @cached_method
505 def _charpoly_basis_space(self):
506 """
507 Return the vector space spanned by the basis used in our
508 characteristic polynomial coefficients. This is used not only to
509 compute those coefficients, but also any time we need to
510 evaluate the coefficients (like when we compute the trace or
511 determinant).
512 """
513 z = self._a_regular_element()
514 V = z.vector().parent().ambient_vector_space()
515 V1 = V.span_of_basis( (z**k).vector() for k in range(self.rank()) )
516 b = (V1.basis() + V1.complement().basis())
517 return V.span_of_basis(b)
518
519
520 @cached_method
521 def _charpoly_coeff(self, i):
522 """
523 Return the coefficient polynomial "a_{i}" of this algebra's
524 general characteristic polynomial.
525
526 Having this be a separate cached method lets us compute and
527 store the trace/determinant (a_{r-1} and a_{0} respectively)
528 separate from the entire characteristic polynomial.
529 """
530 (A_of_x, x, xr, detA) = self._charpoly_matrix_system()
531 R = A_of_x.base_ring()
532 if i >= self.rank():
533 # Guaranteed by theory
534 return R.zero()
535
536 # Danger: the in-place modification is done for performance
537 # reasons (reconstructing a matrix with huge polynomial
538 # entries is slow), but I don't know how cached_method works,
539 # so it's highly possible that we're modifying some global
540 # list variable by reference, here. In other words, you
541 # probably shouldn't call this method twice on the same
542 # algebra, at the same time, in two threads
543 Ai_orig = A_of_x.column(i)
544 A_of_x.set_column(i,xr)
545 numerator = A_of_x.det()
546 A_of_x.set_column(i,Ai_orig)
547
548 # We're relying on the theory here to ensure that each a_i is
549 # indeed back in R, and the added negative signs are to make
550 # the whole charpoly expression sum to zero.
551 return R(-numerator/detA)
552
553
554 @cached_method
555 def _charpoly_matrix_system(self):
556 """
557 Compute the matrix whose entries A_ij are polynomials in
558 X1,...,XN, the vector ``x`` of variables X1,...,XN, the vector
559 corresponding to `x^r` and the determinent of the matrix A =
560 [A_ij]. In other words, all of the fixed (cachable) data needed
561 to compute the coefficients of the characteristic polynomial.
562 """
563 r = self.rank()
564 n = self.dimension()
565
566 # Construct a new algebra over a multivariate polynomial ring...
567 names = ['X' + str(i) for i in range(1,n+1)]
568 R = PolynomialRing(self.base_ring(), names)
569 J = FiniteDimensionalEuclideanJordanAlgebra(R,
570 self._multiplication_table,
571 rank=r)
572
573 idmat = identity_matrix(J.base_ring(), n)
574
575 W = self._charpoly_basis_space()
576 W = W.change_ring(R.fraction_field())
577
578 # Starting with the standard coordinates x = (X1,X2,...,Xn)
579 # and then converting the entries to W-coordinates allows us
580 # to pass in the standard coordinates to the charpoly and get
581 # back the right answer. Specifically, with x = (X1,X2,...,Xn),
582 # we have
583 #
584 # W.coordinates(x^2) eval'd at (standard z-coords)
585 # =
586 # W-coords of (z^2)
587 # =
588 # W-coords of (standard coords of x^2 eval'd at std-coords of z)
589 #
590 # We want the middle equivalent thing in our matrix, but use
591 # the first equivalent thing instead so that we can pass in
592 # standard coordinates.
593 x = J(vector(R, R.gens()))
594 l1 = [column_matrix(W.coordinates((x**k).vector())) for k in range(r)]
595 l2 = [idmat.column(k-1).column() for k in range(r+1, n+1)]
596 A_of_x = block_matrix(R, 1, n, (l1 + l2))
597 xr = W.coordinates((x**r).vector())
598 return (A_of_x, x, xr, A_of_x.det())
599
600
601 @cached_method
602 def characteristic_polynomial(self):
603 """
604
605 .. WARNING::
606
607 This implementation doesn't guarantee that the polynomial
608 denominator in the coefficients is not identically zero, so
609 theoretically it could crash. The way that this is handled
610 in e.g. Faraut and Koranyi is to use a basis that guarantees
611 the denominator is non-zero. But, doing so requires knowledge
612 of at least one regular element, and we don't even know how
613 to do that. The trade-off is that, if we use the standard basis,
614 the resulting polynomial will accept the "usual" coordinates. In
615 other words, we don't have to do a change of basis before e.g.
616 computing the trace or determinant.
617
618 EXAMPLES:
619
620 The characteristic polynomial in the spin algebra is given in
621 Alizadeh, Example 11.11::
622
623 sage: J = JordanSpinEJA(3)
624 sage: p = J.characteristic_polynomial(); p
625 X1^2 - X2^2 - X3^2 + (-2*t)*X1 + t^2
626 sage: xvec = J.one().vector()
627 sage: p(*xvec)
628 t^2 - 2*t + 1
629
630 """
631 r = self.rank()
632 n = self.dimension()
633
634 # The list of coefficient polynomials a_1, a_2, ..., a_n.
635 a = [ self._charpoly_coeff(i) for i in range(n) ]
636
637 # We go to a bit of trouble here to reorder the
638 # indeterminates, so that it's easier to evaluate the
639 # characteristic polynomial at x's coordinates and get back
640 # something in terms of t, which is what we want.
641 R = a[0].parent()
642 S = PolynomialRing(self.base_ring(),'t')
643 t = S.gen(0)
644 S = PolynomialRing(S, R.variable_names())
645 t = S(t)
646
647 # Note: all entries past the rth should be zero. The
648 # coefficient of the highest power (x^r) is 1, but it doesn't
649 # appear in the solution vector which contains coefficients
650 # for the other powers (to make them sum to x^r).
651 if (r < n):
652 a[r] = 1 # corresponds to x^r
653 else:
654 # When the rank is equal to the dimension, trying to
655 # assign a[r] goes out-of-bounds.
656 a.append(1) # corresponds to x^r
657
658 return sum( a[k]*(t**k) for k in range(len(a)) )
659
660
661 def inner_product(self, x, y):
662 """
663 The inner product associated with this Euclidean Jordan algebra.
664
665 Defaults to the trace inner product, but can be overridden by
666 subclasses if they are sure that the necessary properties are
667 satisfied.
668
669 EXAMPLES:
670
671 The inner product must satisfy its axiom for this algebra to truly
672 be a Euclidean Jordan Algebra::
673
674 sage: set_random_seed()
675 sage: J = random_eja()
676 sage: x = J.random_element()
677 sage: y = J.random_element()
678 sage: z = J.random_element()
679 sage: (x*y).inner_product(z) == y.inner_product(x*z)
680 True
681
682 """
683 if (not x in self) or (not y in self):
684 raise TypeError("arguments must live in this algebra")
685 return x.trace_inner_product(y)
686
687
688 def natural_basis(self):
689 """
690 Return a more-natural representation of this algebra's basis.
691
692 Every finite-dimensional Euclidean Jordan Algebra is a direct
693 sum of five simple algebras, four of which comprise Hermitian
694 matrices. This method returns the original "natural" basis
695 for our underlying vector space. (Typically, the natural basis
696 is used to construct the multiplication table in the first place.)
697
698 Note that this will always return a matrix. The standard basis
699 in `R^n` will be returned as `n`-by-`1` column matrices.
700
701 EXAMPLES::
702
703 sage: J = RealSymmetricEJA(2)
704 sage: J.basis()
705 Family (e0, e1, e2)
706 sage: J.natural_basis()
707 (
708 [1 0] [0 1] [0 0]
709 [0 0], [1 0], [0 1]
710 )
711
712 ::
713
714 sage: J = JordanSpinEJA(2)
715 sage: J.basis()
716 Family (e0, e1)
717 sage: J.natural_basis()
718 (
719 [1] [0]
720 [0], [1]
721 )
722
723 """
724 if self._natural_basis is None:
725 return tuple( b.vector().column() for b in self.basis() )
726 else:
727 return self._natural_basis
728
729
730 def rank(self):
731 """
732 Return the rank of this EJA.
733 """
734 if self._rank is None:
735 raise ValueError("no rank specified at genesis")
736 else:
737 return self._rank
738
739
740 class Element(FiniteDimensionalAlgebraElement):
741 """
742 An element of a Euclidean Jordan algebra.
743 """
744
745 def __dir__(self):
746 """
747 Oh man, I should not be doing this. This hides the "disabled"
748 methods ``left_matrix`` and ``matrix`` from introspection;
749 in particular it removes them from tab-completion.
750 """
751 return filter(lambda s: s not in ['left_matrix', 'matrix'],
752 dir(self.__class__) )
753
754
755 def __init__(self, A, elt=None):
756 """
757 EXAMPLES:
758
759 The identity in `S^n` is converted to the identity in the EJA::
760
761 sage: J = RealSymmetricEJA(3)
762 sage: I = identity_matrix(QQ,3)
763 sage: J(I) == J.one()
764 True
765
766 This skew-symmetric matrix can't be represented in the EJA::
767
768 sage: J = RealSymmetricEJA(3)
769 sage: A = matrix(QQ,3, lambda i,j: i-j)
770 sage: J(A)
771 Traceback (most recent call last):
772 ...
773 ArithmeticError: vector is not in free module
774
775 """
776 # Goal: if we're given a matrix, and if it lives in our
777 # parent algebra's "natural ambient space," convert it
778 # into an algebra element.
779 #
780 # The catch is, we make a recursive call after converting
781 # the given matrix into a vector that lives in the algebra.
782 # This we need to try the parent class initializer first,
783 # to avoid recursing forever if we're given something that
784 # already fits into the algebra, but also happens to live
785 # in the parent's "natural ambient space" (this happens with
786 # vectors in R^n).
787 try:
788 FiniteDimensionalAlgebraElement.__init__(self, A, elt)
789 except ValueError:
790 natural_basis = A.natural_basis()
791 if elt in natural_basis[0].matrix_space():
792 # Thanks for nothing! Matrix spaces aren't vector
793 # spaces in Sage, so we have to figure out its
794 # natural-basis coordinates ourselves.
795 V = VectorSpace(elt.base_ring(), elt.nrows()**2)
796 W = V.span( _mat2vec(s) for s in natural_basis )
797 coords = W.coordinates(_mat2vec(elt))
798 FiniteDimensionalAlgebraElement.__init__(self, A, coords)
799
800 def __pow__(self, n):
801 """
802 Return ``self`` raised to the power ``n``.
803
804 Jordan algebras are always power-associative; see for
805 example Faraut and Koranyi, Proposition II.1.2 (ii).
806
807 .. WARNING:
808
809 We have to override this because our superclass uses row vectors
810 instead of column vectors! We, on the other hand, assume column
811 vectors everywhere.
812
813 EXAMPLES::
814
815 sage: set_random_seed()
816 sage: x = random_eja().random_element()
817 sage: x.operator_matrix()*x.vector() == (x^2).vector()
818 True
819
820 A few examples of power-associativity::
821
822 sage: set_random_seed()
823 sage: x = random_eja().random_element()
824 sage: x*(x*x)*(x*x) == x^5
825 True
826 sage: (x*x)*(x*x*x) == x^5
827 True
828
829 We also know that powers operator-commute (Koecher, Chapter
830 III, Corollary 1)::
831
832 sage: set_random_seed()
833 sage: x = random_eja().random_element()
834 sage: m = ZZ.random_element(0,10)
835 sage: n = ZZ.random_element(0,10)
836 sage: Lxm = (x^m).operator_matrix()
837 sage: Lxn = (x^n).operator_matrix()
838 sage: Lxm*Lxn == Lxn*Lxm
839 True
840
841 """
842 A = self.parent()
843 if n == 0:
844 return A.one()
845 elif n == 1:
846 return self
847 else:
848 return A( (self.operator_matrix()**(n-1))*self.vector() )
849
850
851 def apply_univariate_polynomial(self, p):
852 """
853 Apply the univariate polynomial ``p`` to this element.
854
855 A priori, SageMath won't allow us to apply a univariate
856 polynomial to an element of an EJA, because we don't know
857 that EJAs are rings (they are usually not associative). Of
858 course, we know that EJAs are power-associative, so the
859 operation is ultimately kosher. This function sidesteps
860 the CAS to get the answer we want and expect.
861
862 EXAMPLES::
863
864 sage: R = PolynomialRing(QQ, 't')
865 sage: t = R.gen(0)
866 sage: p = t^4 - t^3 + 5*t - 2
867 sage: J = RealCartesianProductEJA(5)
868 sage: J.one().apply_univariate_polynomial(p) == 3*J.one()
869 True
870
871 TESTS:
872
873 We should always get back an element of the algebra::
874
875 sage: set_random_seed()
876 sage: p = PolynomialRing(QQ, 't').random_element()
877 sage: J = random_eja()
878 sage: x = J.random_element()
879 sage: x.apply_univariate_polynomial(p) in J
880 True
881
882 """
883 if len(p.variables()) > 1:
884 raise ValueError("not a univariate polynomial")
885 P = self.parent()
886 R = P.base_ring()
887 # Convert the coeficcients to the parent's base ring,
888 # because a priori they might live in an (unnecessarily)
889 # larger ring for which P.sum() would fail below.
890 cs = [ R(c) for c in p.coefficients(sparse=False) ]
891 return P.sum( cs[k]*(self**k) for k in range(len(cs)) )
892
893
894 def characteristic_polynomial(self):
895 """
896 Return the characteristic polynomial of this element.
897
898 EXAMPLES:
899
900 The rank of `R^3` is three, and the minimal polynomial of
901 the identity element is `(t-1)` from which it follows that
902 the characteristic polynomial should be `(t-1)^3`::
903
904 sage: J = RealCartesianProductEJA(3)
905 sage: J.one().characteristic_polynomial()
906 t^3 - 3*t^2 + 3*t - 1
907
908 Likewise, the characteristic of the zero element in the
909 rank-three algebra `R^{n}` should be `t^{3}`::
910
911 sage: J = RealCartesianProductEJA(3)
912 sage: J.zero().characteristic_polynomial()
913 t^3
914
915 The characteristic polynomial of an element should evaluate
916 to zero on that element::
917
918 sage: set_random_seed()
919 sage: x = RealCartesianProductEJA(3).random_element()
920 sage: p = x.characteristic_polynomial()
921 sage: x.apply_univariate_polynomial(p)
922 0
923
924 """
925 p = self.parent().characteristic_polynomial()
926 return p(*self.vector())
927
928
929 def inner_product(self, other):
930 """
931 Return the parent algebra's inner product of myself and ``other``.
932
933 EXAMPLES:
934
935 The inner product in the Jordan spin algebra is the usual
936 inner product on `R^n` (this example only works because the
937 basis for the Jordan algebra is the standard basis in `R^n`)::
938
939 sage: J = JordanSpinEJA(3)
940 sage: x = vector(QQ,[1,2,3])
941 sage: y = vector(QQ,[4,5,6])
942 sage: x.inner_product(y)
943 32
944 sage: J(x).inner_product(J(y))
945 32
946
947 The inner product on `S^n` is `<X,Y> = trace(X*Y)`, where
948 multiplication is the usual matrix multiplication in `S^n`,
949 so the inner product of the identity matrix with itself
950 should be the `n`::
951
952 sage: J = RealSymmetricEJA(3)
953 sage: J.one().inner_product(J.one())
954 3
955
956 Likewise, the inner product on `C^n` is `<X,Y> =
957 Re(trace(X*Y))`, where we must necessarily take the real
958 part because the product of Hermitian matrices may not be
959 Hermitian::
960
961 sage: J = ComplexHermitianEJA(3)
962 sage: J.one().inner_product(J.one())
963 3
964
965 Ditto for the quaternions::
966
967 sage: J = QuaternionHermitianEJA(3)
968 sage: J.one().inner_product(J.one())
969 3
970
971 TESTS:
972
973 Ensure that we can always compute an inner product, and that
974 it gives us back a real number::
975
976 sage: set_random_seed()
977 sage: J = random_eja()
978 sage: x = J.random_element()
979 sage: y = J.random_element()
980 sage: x.inner_product(y) in RR
981 True
982
983 """
984 P = self.parent()
985 if not other in P:
986 raise TypeError("'other' must live in the same algebra")
987
988 return P.inner_product(self, other)
989
990
991 def operator_commutes_with(self, other):
992 """
993 Return whether or not this element operator-commutes
994 with ``other``.
995
996 EXAMPLES:
997
998 The definition of a Jordan algebra says that any element
999 operator-commutes with its square::
1000
1001 sage: set_random_seed()
1002 sage: x = random_eja().random_element()
1003 sage: x.operator_commutes_with(x^2)
1004 True
1005
1006 TESTS:
1007
1008 Test Lemma 1 from Chapter III of Koecher::
1009
1010 sage: set_random_seed()
1011 sage: J = random_eja()
1012 sage: u = J.random_element()
1013 sage: v = J.random_element()
1014 sage: lhs = u.operator_commutes_with(u*v)
1015 sage: rhs = v.operator_commutes_with(u^2)
1016 sage: lhs == rhs
1017 True
1018
1019 """
1020 if not other in self.parent():
1021 raise TypeError("'other' must live in the same algebra")
1022
1023 A = self.operator_matrix()
1024 B = other.operator_matrix()
1025 return (A*B == B*A)
1026
1027
1028 def det(self):
1029 """
1030 Return my determinant, the product of my eigenvalues.
1031
1032 EXAMPLES::
1033
1034 sage: J = JordanSpinEJA(2)
1035 sage: e0,e1 = J.gens()
1036 sage: x = sum( J.gens() )
1037 sage: x.det()
1038 0
1039
1040 ::
1041
1042 sage: J = JordanSpinEJA(3)
1043 sage: e0,e1,e2 = J.gens()
1044 sage: x = sum( J.gens() )
1045 sage: x.det()
1046 -1
1047
1048 TESTS:
1049
1050 An element is invertible if and only if its determinant is
1051 non-zero::
1052
1053 sage: set_random_seed()
1054 sage: x = random_eja().random_element()
1055 sage: x.is_invertible() == (x.det() != 0)
1056 True
1057
1058 """
1059 P = self.parent()
1060 r = P.rank()
1061 p = P._charpoly_coeff(0)
1062 # The _charpoly_coeff function already adds the factor of
1063 # -1 to ensure that _charpoly_coeff(0) is really what
1064 # appears in front of t^{0} in the charpoly. However,
1065 # we want (-1)^r times THAT for the determinant.
1066 return ((-1)**r)*p(*self.vector())
1067
1068
1069 def inverse(self):
1070 """
1071 Return the Jordan-multiplicative inverse of this element.
1072
1073 ALGORITHM:
1074
1075 We appeal to the quadratic representation as in Koecher's
1076 Theorem 12 in Chapter III, Section 5.
1077
1078 EXAMPLES:
1079
1080 The inverse in the spin factor algebra is given in Alizadeh's
1081 Example 11.11::
1082
1083 sage: set_random_seed()
1084 sage: n = ZZ.random_element(1,10)
1085 sage: J = JordanSpinEJA(n)
1086 sage: x = J.random_element()
1087 sage: while not x.is_invertible():
1088 ....: x = J.random_element()
1089 sage: x_vec = x.vector()
1090 sage: x0 = x_vec[0]
1091 sage: x_bar = x_vec[1:]
1092 sage: coeff = ~(x0^2 - x_bar.inner_product(x_bar))
1093 sage: inv_vec = x_vec.parent()([x0] + (-x_bar).list())
1094 sage: x_inverse = coeff*inv_vec
1095 sage: x.inverse() == J(x_inverse)
1096 True
1097
1098 TESTS:
1099
1100 The identity element is its own inverse::
1101
1102 sage: set_random_seed()
1103 sage: J = random_eja()
1104 sage: J.one().inverse() == J.one()
1105 True
1106
1107 If an element has an inverse, it acts like one::
1108
1109 sage: set_random_seed()
1110 sage: J = random_eja()
1111 sage: x = J.random_element()
1112 sage: (not x.is_invertible()) or (x.inverse()*x == J.one())
1113 True
1114
1115 The inverse of the inverse is what we started with::
1116
1117 sage: set_random_seed()
1118 sage: J = random_eja()
1119 sage: x = J.random_element()
1120 sage: (not x.is_invertible()) or (x.inverse().inverse() == x)
1121 True
1122
1123 The zero element is never invertible::
1124
1125 sage: set_random_seed()
1126 sage: J = random_eja().zero().inverse()
1127 Traceback (most recent call last):
1128 ...
1129 ValueError: element is not invertible
1130
1131 """
1132 if not self.is_invertible():
1133 raise ValueError("element is not invertible")
1134
1135 return (~self.quadratic_representation())(self)
1136
1137
1138 def is_invertible(self):
1139 """
1140 Return whether or not this element is invertible.
1141
1142 We can't use the superclass method because it relies on
1143 the algebra being associative.
1144
1145 ALGORITHM:
1146
1147 The usual way to do this is to check if the determinant is
1148 zero, but we need the characteristic polynomial for the
1149 determinant. The minimal polynomial is a lot easier to get,
1150 so we use Corollary 2 in Chapter V of Koecher to check
1151 whether or not the paren't algebra's zero element is a root
1152 of this element's minimal polynomial.
1153
1154 TESTS:
1155
1156 The identity element is always invertible::
1157
1158 sage: set_random_seed()
1159 sage: J = random_eja()
1160 sage: J.one().is_invertible()
1161 True
1162
1163 The zero element is never invertible::
1164
1165 sage: set_random_seed()
1166 sage: J = random_eja()
1167 sage: J.zero().is_invertible()
1168 False
1169
1170 """
1171 zero = self.parent().zero()
1172 p = self.minimal_polynomial()
1173 return not (p(zero) == zero)
1174
1175
1176 def is_nilpotent(self):
1177 """
1178 Return whether or not some power of this element is zero.
1179
1180 The superclass method won't work unless we're in an
1181 associative algebra, and we aren't. However, we generate
1182 an assocoative subalgebra and we're nilpotent there if and
1183 only if we're nilpotent here (probably).
1184
1185 TESTS:
1186
1187 The identity element is never nilpotent::
1188
1189 sage: set_random_seed()
1190 sage: random_eja().one().is_nilpotent()
1191 False
1192
1193 The additive identity is always nilpotent::
1194
1195 sage: set_random_seed()
1196 sage: random_eja().zero().is_nilpotent()
1197 True
1198
1199 """
1200 # The element we're going to call "is_nilpotent()" on.
1201 # Either myself, interpreted as an element of a finite-
1202 # dimensional algebra, or an element of an associative
1203 # subalgebra.
1204 elt = None
1205
1206 if self.parent().is_associative():
1207 elt = FiniteDimensionalAlgebraElement(self.parent(), self)
1208 else:
1209 V = self.span_of_powers()
1210 assoc_subalg = self.subalgebra_generated_by()
1211 # Mis-design warning: the basis used for span_of_powers()
1212 # and subalgebra_generated_by() must be the same, and in
1213 # the same order!
1214 elt = assoc_subalg(V.coordinates(self.vector()))
1215
1216 # Recursive call, but should work since elt lives in an
1217 # associative algebra.
1218 return elt.is_nilpotent()
1219
1220
1221 def is_regular(self):
1222 """
1223 Return whether or not this is a regular element.
1224
1225 EXAMPLES:
1226
1227 The identity element always has degree one, but any element
1228 linearly-independent from it is regular::
1229
1230 sage: J = JordanSpinEJA(5)
1231 sage: J.one().is_regular()
1232 False
1233 sage: e0, e1, e2, e3, e4 = J.gens() # e0 is the identity
1234 sage: for x in J.gens():
1235 ....: (J.one() + x).is_regular()
1236 False
1237 True
1238 True
1239 True
1240 True
1241
1242 """
1243 return self.degree() == self.parent().rank()
1244
1245
1246 def degree(self):
1247 """
1248 Compute the degree of this element the straightforward way
1249 according to the definition; by appending powers to a list
1250 and figuring out its dimension (that is, whether or not
1251 they're linearly dependent).
1252
1253 EXAMPLES::
1254
1255 sage: J = JordanSpinEJA(4)
1256 sage: J.one().degree()
1257 1
1258 sage: e0,e1,e2,e3 = J.gens()
1259 sage: (e0 - e1).degree()
1260 2
1261
1262 In the spin factor algebra (of rank two), all elements that
1263 aren't multiples of the identity are regular::
1264
1265 sage: set_random_seed()
1266 sage: n = ZZ.random_element(1,10)
1267 sage: J = JordanSpinEJA(n)
1268 sage: x = J.random_element()
1269 sage: x == x.coefficient(0)*J.one() or x.degree() == 2
1270 True
1271
1272 """
1273 return self.span_of_powers().dimension()
1274
1275
1276 def left_matrix(self):
1277 """
1278 Our parent class defines ``left_matrix`` and ``matrix``
1279 methods whose names are misleading. We don't want them.
1280 """
1281 raise NotImplementedError("use operator_matrix() instead")
1282
1283 matrix = left_matrix
1284
1285
1286 def minimal_polynomial(self):
1287 """
1288 Return the minimal polynomial of this element,
1289 as a function of the variable `t`.
1290
1291 ALGORITHM:
1292
1293 We restrict ourselves to the associative subalgebra
1294 generated by this element, and then return the minimal
1295 polynomial of this element's operator matrix (in that
1296 subalgebra). This works by Baes Proposition 2.3.16.
1297
1298 TESTS:
1299
1300 The minimal polynomial of the identity and zero elements are
1301 always the same::
1302
1303 sage: set_random_seed()
1304 sage: J = random_eja()
1305 sage: J.one().minimal_polynomial()
1306 t - 1
1307 sage: J.zero().minimal_polynomial()
1308 t
1309
1310 The degree of an element is (by one definition) the degree
1311 of its minimal polynomial::
1312
1313 sage: set_random_seed()
1314 sage: x = random_eja().random_element()
1315 sage: x.degree() == x.minimal_polynomial().degree()
1316 True
1317
1318 The minimal polynomial and the characteristic polynomial coincide
1319 and are known (see Alizadeh, Example 11.11) for all elements of
1320 the spin factor algebra that aren't scalar multiples of the
1321 identity::
1322
1323 sage: set_random_seed()
1324 sage: n = ZZ.random_element(2,10)
1325 sage: J = JordanSpinEJA(n)
1326 sage: y = J.random_element()
1327 sage: while y == y.coefficient(0)*J.one():
1328 ....: y = J.random_element()
1329 sage: y0 = y.vector()[0]
1330 sage: y_bar = y.vector()[1:]
1331 sage: actual = y.minimal_polynomial()
1332 sage: t = PolynomialRing(J.base_ring(),'t').gen(0)
1333 sage: expected = t^2 - 2*y0*t + (y0^2 - norm(y_bar)^2)
1334 sage: bool(actual == expected)
1335 True
1336
1337 The minimal polynomial should always kill its element::
1338
1339 sage: set_random_seed()
1340 sage: x = random_eja().random_element()
1341 sage: p = x.minimal_polynomial()
1342 sage: x.apply_univariate_polynomial(p)
1343 0
1344
1345 """
1346 V = self.span_of_powers()
1347 assoc_subalg = self.subalgebra_generated_by()
1348 # Mis-design warning: the basis used for span_of_powers()
1349 # and subalgebra_generated_by() must be the same, and in
1350 # the same order!
1351 elt = assoc_subalg(V.coordinates(self.vector()))
1352
1353 # We get back a symbolic polynomial in 'x' but want a real
1354 # polynomial in 't'.
1355 p_of_x = elt.operator_matrix().minimal_polynomial()
1356 return p_of_x.change_variable_name('t')
1357
1358
1359 def natural_representation(self):
1360 """
1361 Return a more-natural representation of this element.
1362
1363 Every finite-dimensional Euclidean Jordan Algebra is a
1364 direct sum of five simple algebras, four of which comprise
1365 Hermitian matrices. This method returns the original
1366 "natural" representation of this element as a Hermitian
1367 matrix, if it has one. If not, you get the usual representation.
1368
1369 EXAMPLES::
1370
1371 sage: J = ComplexHermitianEJA(3)
1372 sage: J.one()
1373 e0 + e5 + e8
1374 sage: J.one().natural_representation()
1375 [1 0 0 0 0 0]
1376 [0 1 0 0 0 0]
1377 [0 0 1 0 0 0]
1378 [0 0 0 1 0 0]
1379 [0 0 0 0 1 0]
1380 [0 0 0 0 0 1]
1381
1382 ::
1383
1384 sage: J = QuaternionHermitianEJA(3)
1385 sage: J.one()
1386 e0 + e9 + e14
1387 sage: J.one().natural_representation()
1388 [1 0 0 0 0 0 0 0 0 0 0 0]
1389 [0 1 0 0 0 0 0 0 0 0 0 0]
1390 [0 0 1 0 0 0 0 0 0 0 0 0]
1391 [0 0 0 1 0 0 0 0 0 0 0 0]
1392 [0 0 0 0 1 0 0 0 0 0 0 0]
1393 [0 0 0 0 0 1 0 0 0 0 0 0]
1394 [0 0 0 0 0 0 1 0 0 0 0 0]
1395 [0 0 0 0 0 0 0 1 0 0 0 0]
1396 [0 0 0 0 0 0 0 0 1 0 0 0]
1397 [0 0 0 0 0 0 0 0 0 1 0 0]
1398 [0 0 0 0 0 0 0 0 0 0 1 0]
1399 [0 0 0 0 0 0 0 0 0 0 0 1]
1400
1401 """
1402 B = self.parent().natural_basis()
1403 W = B[0].matrix_space()
1404 return W.linear_combination(zip(self.vector(), B))
1405
1406
1407 def operator(self):
1408 """
1409 Return the left-multiplication-by-this-element
1410 operator on the ambient algebra.
1411
1412 TESTS::
1413
1414 sage: set_random_seed()
1415 sage: J = random_eja()
1416 sage: x = J.random_element()
1417 sage: y = J.random_element()
1418 sage: x.operator()(y) == x*y
1419 True
1420 sage: y.operator()(x) == x*y
1421 True
1422
1423 """
1424 P = self.parent()
1425 return FiniteDimensionalEuclideanJordanAlgebraMorphism(
1426 Hom(P,P),
1427 self.operator_matrix() )
1428
1429
1430
1431 def operator_matrix(self):
1432 """
1433 Return the matrix that represents left- (or right-)
1434 multiplication by this element in the parent algebra.
1435
1436 We implement this ourselves to work around the fact that
1437 our parent class represents everything with row vectors.
1438
1439 EXAMPLES:
1440
1441 Test the first polarization identity from my notes, Koecher Chapter
1442 III, or from Baes (2.3)::
1443
1444 sage: set_random_seed()
1445 sage: J = random_eja()
1446 sage: x = J.random_element()
1447 sage: y = J.random_element()
1448 sage: Lx = x.operator_matrix()
1449 sage: Ly = y.operator_matrix()
1450 sage: Lxx = (x*x).operator_matrix()
1451 sage: Lxy = (x*y).operator_matrix()
1452 sage: bool(2*Lx*Lxy + Ly*Lxx == 2*Lxy*Lx + Lxx*Ly)
1453 True
1454
1455 Test the second polarization identity from my notes or from
1456 Baes (2.4)::
1457
1458 sage: set_random_seed()
1459 sage: J = random_eja()
1460 sage: x = J.random_element()
1461 sage: y = J.random_element()
1462 sage: z = J.random_element()
1463 sage: Lx = x.operator_matrix()
1464 sage: Ly = y.operator_matrix()
1465 sage: Lz = z.operator_matrix()
1466 sage: Lzy = (z*y).operator_matrix()
1467 sage: Lxy = (x*y).operator_matrix()
1468 sage: Lxz = (x*z).operator_matrix()
1469 sage: bool(Lx*Lzy + Lz*Lxy + Ly*Lxz == Lzy*Lx + Lxy*Lz + Lxz*Ly)
1470 True
1471
1472 Test the third polarization identity from my notes or from
1473 Baes (2.5)::
1474
1475 sage: set_random_seed()
1476 sage: J = random_eja()
1477 sage: u = J.random_element()
1478 sage: y = J.random_element()
1479 sage: z = J.random_element()
1480 sage: Lu = u.operator_matrix()
1481 sage: Ly = y.operator_matrix()
1482 sage: Lz = z.operator_matrix()
1483 sage: Lzy = (z*y).operator_matrix()
1484 sage: Luy = (u*y).operator_matrix()
1485 sage: Luz = (u*z).operator_matrix()
1486 sage: Luyz = (u*(y*z)).operator_matrix()
1487 sage: lhs = Lu*Lzy + Lz*Luy + Ly*Luz
1488 sage: rhs = Luyz + Ly*Lu*Lz + Lz*Lu*Ly
1489 sage: bool(lhs == rhs)
1490 True
1491
1492 """
1493 fda_elt = FiniteDimensionalAlgebraElement(self.parent(), self)
1494 return fda_elt.matrix().transpose()
1495
1496
1497 def quadratic_representation(self, other=None):
1498 """
1499 Return the quadratic representation of this element.
1500
1501 EXAMPLES:
1502
1503 The explicit form in the spin factor algebra is given by
1504 Alizadeh's Example 11.12::
1505
1506 sage: set_random_seed()
1507 sage: n = ZZ.random_element(1,10)
1508 sage: J = JordanSpinEJA(n)
1509 sage: x = J.random_element()
1510 sage: x_vec = x.vector()
1511 sage: x0 = x_vec[0]
1512 sage: x_bar = x_vec[1:]
1513 sage: A = matrix(QQ, 1, [x_vec.inner_product(x_vec)])
1514 sage: B = 2*x0*x_bar.row()
1515 sage: C = 2*x0*x_bar.column()
1516 sage: D = identity_matrix(QQ, n-1)
1517 sage: D = (x0^2 - x_bar.inner_product(x_bar))*D
1518 sage: D = D + 2*x_bar.tensor_product(x_bar)
1519 sage: Q = block_matrix(2,2,[A,B,C,D])
1520 sage: Q == x.quadratic_representation().matrix()
1521 True
1522
1523 Test all of the properties from Theorem 11.2 in Alizadeh::
1524
1525 sage: set_random_seed()
1526 sage: J = random_eja()
1527 sage: x = J.random_element()
1528 sage: y = J.random_element()
1529 sage: Lx = x.operator()
1530 sage: Lxx = (x*x).operator()
1531 sage: Qx = x.quadratic_representation()
1532 sage: Qy = y.quadratic_representation()
1533 sage: Qxy = x.quadratic_representation(y)
1534 sage: Qex = J.one().quadratic_representation(x)
1535 sage: n = ZZ.random_element(10)
1536 sage: Qxn = (x^n).quadratic_representation()
1537
1538 Property 1:
1539
1540 sage: 2*Qxy == (x+y).quadratic_representation() - Qx - Qy
1541 True
1542
1543 Property 2:
1544
1545 sage: alpha = QQ.random_element()
1546 sage: (alpha*x).quadratic_representation() == (alpha^2)*Qx
1547 True
1548
1549 Property 3:
1550
1551 sage: not x.is_invertible() or ( Qx(x.inverse()) == x )
1552 True
1553
1554 sage: not x.is_invertible() or (
1555 ....: ~Qx
1556 ....: ==
1557 ....: x.inverse().quadratic_representation() )
1558 True
1559
1560 sage: Qxy(J.one()) == x*y
1561 True
1562
1563 Property 4:
1564
1565 sage: not x.is_invertible() or (
1566 ....: x.quadratic_representation(x.inverse())*Qx
1567 ....: == Qx*x.quadratic_representation(x.inverse()) )
1568 True
1569
1570 sage: not x.is_invertible() or (
1571 ....: x.quadratic_representation(x.inverse())*Qx
1572 ....: ==
1573 ....: 2*x.operator()*Qex - Qx )
1574 True
1575
1576 sage: 2*x.operator()*Qex - Qx == Lxx
1577 True
1578
1579 Property 5:
1580
1581 sage: Qy(x).quadratic_representation() == Qy*Qx*Qy
1582 True
1583
1584 Property 6:
1585
1586 sage: Qxn == (Qx)^n
1587 True
1588
1589 Property 7:
1590
1591 sage: not x.is_invertible() or (
1592 ....: Qx*x.inverse().operator() == Lx )
1593 True
1594
1595 Property 8:
1596
1597 sage: not x.operator_commutes_with(y) or (
1598 ....: Qx(y)^n == Qxn(y^n) )
1599 True
1600
1601 """
1602 if other is None:
1603 other=self
1604 elif not other in self.parent():
1605 raise TypeError("'other' must live in the same algebra")
1606
1607 L = self.operator()
1608 M = other.operator()
1609 return ( L*M + M*L - (self*other).operator() )
1610
1611
1612 def span_of_powers(self):
1613 """
1614 Return the vector space spanned by successive powers of
1615 this element.
1616 """
1617 # The dimension of the subalgebra can't be greater than
1618 # the big algebra, so just put everything into a list
1619 # and let span() get rid of the excess.
1620 #
1621 # We do the extra ambient_vector_space() in case we're messing
1622 # with polynomials and the direct parent is a module.
1623 V = self.vector().parent().ambient_vector_space()
1624 return V.span( (self**d).vector() for d in xrange(V.dimension()) )
1625
1626
1627 def subalgebra_generated_by(self):
1628 """
1629 Return the associative subalgebra of the parent EJA generated
1630 by this element.
1631
1632 TESTS::
1633
1634 sage: set_random_seed()
1635 sage: x = random_eja().random_element()
1636 sage: x.subalgebra_generated_by().is_associative()
1637 True
1638
1639 Squaring in the subalgebra should be the same thing as
1640 squaring in the superalgebra::
1641
1642 sage: set_random_seed()
1643 sage: x = random_eja().random_element()
1644 sage: u = x.subalgebra_generated_by().random_element()
1645 sage: u.operator_matrix()*u.vector() == (u**2).vector()
1646 True
1647
1648 """
1649 # First get the subspace spanned by the powers of myself...
1650 V = self.span_of_powers()
1651 F = self.base_ring()
1652
1653 # Now figure out the entries of the right-multiplication
1654 # matrix for the successive basis elements b0, b1,... of
1655 # that subspace.
1656 mats = []
1657 for b_right in V.basis():
1658 eja_b_right = self.parent()(b_right)
1659 b_right_rows = []
1660 # The first row of the right-multiplication matrix by
1661 # b1 is what we get if we apply that matrix to b1. The
1662 # second row of the right multiplication matrix by b1
1663 # is what we get when we apply that matrix to b2...
1664 #
1665 # IMPORTANT: this assumes that all vectors are COLUMN
1666 # vectors, unlike our superclass (which uses row vectors).
1667 for b_left in V.basis():
1668 eja_b_left = self.parent()(b_left)
1669 # Multiply in the original EJA, but then get the
1670 # coordinates from the subalgebra in terms of its
1671 # basis.
1672 this_row = V.coordinates((eja_b_left*eja_b_right).vector())
1673 b_right_rows.append(this_row)
1674 b_right_matrix = matrix(F, b_right_rows)
1675 mats.append(b_right_matrix)
1676
1677 # It's an algebra of polynomials in one element, and EJAs
1678 # are power-associative.
1679 #
1680 # TODO: choose generator names intelligently.
1681 return FiniteDimensionalEuclideanJordanAlgebra(F, mats, assume_associative=True, names='f')
1682
1683
1684 def subalgebra_idempotent(self):
1685 """
1686 Find an idempotent in the associative subalgebra I generate
1687 using Proposition 2.3.5 in Baes.
1688
1689 TESTS::
1690
1691 sage: set_random_seed()
1692 sage: J = random_eja()
1693 sage: x = J.random_element()
1694 sage: while x.is_nilpotent():
1695 ....: x = J.random_element()
1696 sage: c = x.subalgebra_idempotent()
1697 sage: c^2 == c
1698 True
1699
1700 """
1701 if self.is_nilpotent():
1702 raise ValueError("this only works with non-nilpotent elements!")
1703
1704 V = self.span_of_powers()
1705 J = self.subalgebra_generated_by()
1706 # Mis-design warning: the basis used for span_of_powers()
1707 # and subalgebra_generated_by() must be the same, and in
1708 # the same order!
1709 u = J(V.coordinates(self.vector()))
1710
1711 # The image of the matrix of left-u^m-multiplication
1712 # will be minimal for some natural number s...
1713 s = 0
1714 minimal_dim = V.dimension()
1715 for i in xrange(1, V.dimension()):
1716 this_dim = (u**i).operator_matrix().image().dimension()
1717 if this_dim < minimal_dim:
1718 minimal_dim = this_dim
1719 s = i
1720
1721 # Now minimal_matrix should correspond to the smallest
1722 # non-zero subspace in Baes's (or really, Koecher's)
1723 # proposition.
1724 #
1725 # However, we need to restrict the matrix to work on the
1726 # subspace... or do we? Can't we just solve, knowing that
1727 # A(c) = u^(s+1) should have a solution in the big space,
1728 # too?
1729 #
1730 # Beware, solve_right() means that we're using COLUMN vectors.
1731 # Our FiniteDimensionalAlgebraElement superclass uses rows.
1732 u_next = u**(s+1)
1733 A = u_next.operator_matrix()
1734 c_coordinates = A.solve_right(u_next.vector())
1735
1736 # Now c_coordinates is the idempotent we want, but it's in
1737 # the coordinate system of the subalgebra.
1738 #
1739 # We need the basis for J, but as elements of the parent algebra.
1740 #
1741 basis = [self.parent(v) for v in V.basis()]
1742 return self.parent().linear_combination(zip(c_coordinates, basis))
1743
1744
1745 def trace(self):
1746 """
1747 Return my trace, the sum of my eigenvalues.
1748
1749 EXAMPLES::
1750
1751 sage: J = JordanSpinEJA(3)
1752 sage: x = sum(J.gens())
1753 sage: x.trace()
1754 2
1755
1756 ::
1757
1758 sage: J = RealCartesianProductEJA(5)
1759 sage: J.one().trace()
1760 5
1761
1762 TESTS:
1763
1764 The trace of an element is a real number::
1765
1766 sage: set_random_seed()
1767 sage: J = random_eja()
1768 sage: J.random_element().trace() in J.base_ring()
1769 True
1770
1771 """
1772 P = self.parent()
1773 r = P.rank()
1774 p = P._charpoly_coeff(r-1)
1775 # The _charpoly_coeff function already adds the factor of
1776 # -1 to ensure that _charpoly_coeff(r-1) is really what
1777 # appears in front of t^{r-1} in the charpoly. However,
1778 # we want the negative of THAT for the trace.
1779 return -p(*self.vector())
1780
1781
1782 def trace_inner_product(self, other):
1783 """
1784 Return the trace inner product of myself and ``other``.
1785
1786 TESTS:
1787
1788 The trace inner product is commutative::
1789
1790 sage: set_random_seed()
1791 sage: J = random_eja()
1792 sage: x = J.random_element(); y = J.random_element()
1793 sage: x.trace_inner_product(y) == y.trace_inner_product(x)
1794 True
1795
1796 The trace inner product is bilinear::
1797
1798 sage: set_random_seed()
1799 sage: J = random_eja()
1800 sage: x = J.random_element()
1801 sage: y = J.random_element()
1802 sage: z = J.random_element()
1803 sage: a = QQ.random_element();
1804 sage: actual = (a*(x+z)).trace_inner_product(y)
1805 sage: expected = ( a*x.trace_inner_product(y) +
1806 ....: a*z.trace_inner_product(y) )
1807 sage: actual == expected
1808 True
1809 sage: actual = x.trace_inner_product(a*(y+z))
1810 sage: expected = ( a*x.trace_inner_product(y) +
1811 ....: a*x.trace_inner_product(z) )
1812 sage: actual == expected
1813 True
1814
1815 The trace inner product satisfies the compatibility
1816 condition in the definition of a Euclidean Jordan algebra::
1817
1818 sage: set_random_seed()
1819 sage: J = random_eja()
1820 sage: x = J.random_element()
1821 sage: y = J.random_element()
1822 sage: z = J.random_element()
1823 sage: (x*y).trace_inner_product(z) == y.trace_inner_product(x*z)
1824 True
1825
1826 """
1827 if not other in self.parent():
1828 raise TypeError("'other' must live in the same algebra")
1829
1830 return (self*other).trace()
1831
1832
1833 class RealCartesianProductEJA(FiniteDimensionalEuclideanJordanAlgebra):
1834 """
1835 Return the Euclidean Jordan Algebra corresponding to the set
1836 `R^n` under the Hadamard product.
1837
1838 Note: this is nothing more than the Cartesian product of ``n``
1839 copies of the spin algebra. Once Cartesian product algebras
1840 are implemented, this can go.
1841
1842 EXAMPLES:
1843
1844 This multiplication table can be verified by hand::
1845
1846 sage: J = RealCartesianProductEJA(3)
1847 sage: e0,e1,e2 = J.gens()
1848 sage: e0*e0
1849 e0
1850 sage: e0*e1
1851 0
1852 sage: e0*e2
1853 0
1854 sage: e1*e1
1855 e1
1856 sage: e1*e2
1857 0
1858 sage: e2*e2
1859 e2
1860
1861 """
1862 @staticmethod
1863 def __classcall_private__(cls, n, field=QQ):
1864 # The FiniteDimensionalAlgebra constructor takes a list of
1865 # matrices, the ith representing right multiplication by the ith
1866 # basis element in the vector space. So if e_1 = (1,0,0), then
1867 # right (Hadamard) multiplication of x by e_1 picks out the first
1868 # component of x; and likewise for the ith basis element e_i.
1869 Qs = [ matrix(field, n, n, lambda k,j: 1*(k == j == i))
1870 for i in xrange(n) ]
1871
1872 fdeja = super(RealCartesianProductEJA, cls)
1873 return fdeja.__classcall_private__(cls, field, Qs, rank=n)
1874
1875 def inner_product(self, x, y):
1876 return _usual_ip(x,y)
1877
1878
1879 def random_eja():
1880 """
1881 Return a "random" finite-dimensional Euclidean Jordan Algebra.
1882
1883 ALGORITHM:
1884
1885 For now, we choose a random natural number ``n`` (greater than zero)
1886 and then give you back one of the following:
1887
1888 * The cartesian product of the rational numbers ``n`` times; this is
1889 ``QQ^n`` with the Hadamard product.
1890
1891 * The Jordan spin algebra on ``QQ^n``.
1892
1893 * The ``n``-by-``n`` rational symmetric matrices with the symmetric
1894 product.
1895
1896 * The ``n``-by-``n`` complex-rational Hermitian matrices embedded
1897 in the space of ``2n``-by-``2n`` real symmetric matrices.
1898
1899 * The ``n``-by-``n`` quaternion-rational Hermitian matrices embedded
1900 in the space of ``4n``-by-``4n`` real symmetric matrices.
1901
1902 Later this might be extended to return Cartesian products of the
1903 EJAs above.
1904
1905 TESTS::
1906
1907 sage: random_eja()
1908 Euclidean Jordan algebra of degree...
1909
1910 """
1911
1912 # The max_n component lets us choose different upper bounds on the
1913 # value "n" that gets passed to the constructor. This is needed
1914 # because e.g. R^{10} is reasonable to test, while the Hermitian
1915 # 10-by-10 quaternion matrices are not.
1916 (constructor, max_n) = choice([(RealCartesianProductEJA, 6),
1917 (JordanSpinEJA, 6),
1918 (RealSymmetricEJA, 5),
1919 (ComplexHermitianEJA, 4),
1920 (QuaternionHermitianEJA, 3)])
1921 n = ZZ.random_element(1, max_n)
1922 return constructor(n, field=QQ)
1923
1924
1925
1926 def _real_symmetric_basis(n, field=QQ):
1927 """
1928 Return a basis for the space of real symmetric n-by-n matrices.
1929 """
1930 # The basis of symmetric matrices, as matrices, in their R^(n-by-n)
1931 # coordinates.
1932 S = []
1933 for i in xrange(n):
1934 for j in xrange(i+1):
1935 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1936 if i == j:
1937 Sij = Eij
1938 else:
1939 # Beware, orthogonal but not normalized!
1940 Sij = Eij + Eij.transpose()
1941 S.append(Sij)
1942 return tuple(S)
1943
1944
1945 def _complex_hermitian_basis(n, field=QQ):
1946 """
1947 Returns a basis for the space of complex Hermitian n-by-n matrices.
1948
1949 TESTS::
1950
1951 sage: set_random_seed()
1952 sage: n = ZZ.random_element(1,5)
1953 sage: all( M.is_symmetric() for M in _complex_hermitian_basis(n) )
1954 True
1955
1956 """
1957 F = QuadraticField(-1, 'I')
1958 I = F.gen()
1959
1960 # This is like the symmetric case, but we need to be careful:
1961 #
1962 # * We want conjugate-symmetry, not just symmetry.
1963 # * The diagonal will (as a result) be real.
1964 #
1965 S = []
1966 for i in xrange(n):
1967 for j in xrange(i+1):
1968 Eij = matrix(field, n, lambda k,l: k==i and l==j)
1969 if i == j:
1970 Sij = _embed_complex_matrix(Eij)
1971 S.append(Sij)
1972 else:
1973 # Beware, orthogonal but not normalized! The second one
1974 # has a minus because it's conjugated.
1975 Sij_real = _embed_complex_matrix(Eij + Eij.transpose())
1976 S.append(Sij_real)
1977 Sij_imag = _embed_complex_matrix(I*Eij - I*Eij.transpose())
1978 S.append(Sij_imag)
1979 return tuple(S)
1980
1981
1982 def _quaternion_hermitian_basis(n, field=QQ):
1983 """
1984 Returns a basis for the space of quaternion Hermitian n-by-n matrices.
1985
1986 TESTS::
1987
1988 sage: set_random_seed()
1989 sage: n = ZZ.random_element(1,5)
1990 sage: all( M.is_symmetric() for M in _quaternion_hermitian_basis(n) )
1991 True
1992
1993 """
1994 Q = QuaternionAlgebra(QQ,-1,-1)
1995 I,J,K = Q.gens()
1996
1997 # This is like the symmetric case, but we need to be careful:
1998 #
1999 # * We want conjugate-symmetry, not just symmetry.
2000 # * The diagonal will (as a result) be real.
2001 #
2002 S = []
2003 for i in xrange(n):
2004 for j in xrange(i+1):
2005 Eij = matrix(Q, n, lambda k,l: k==i and l==j)
2006 if i == j:
2007 Sij = _embed_quaternion_matrix(Eij)
2008 S.append(Sij)
2009 else:
2010 # Beware, orthogonal but not normalized! The second,
2011 # third, and fourth ones have a minus because they're
2012 # conjugated.
2013 Sij_real = _embed_quaternion_matrix(Eij + Eij.transpose())
2014 S.append(Sij_real)
2015 Sij_I = _embed_quaternion_matrix(I*Eij - I*Eij.transpose())
2016 S.append(Sij_I)
2017 Sij_J = _embed_quaternion_matrix(J*Eij - J*Eij.transpose())
2018 S.append(Sij_J)
2019 Sij_K = _embed_quaternion_matrix(K*Eij - K*Eij.transpose())
2020 S.append(Sij_K)
2021 return tuple(S)
2022
2023
2024 def _mat2vec(m):
2025 return vector(m.base_ring(), m.list())
2026
2027 def _vec2mat(v):
2028 return matrix(v.base_ring(), sqrt(v.degree()), v.list())
2029
2030 def _multiplication_table_from_matrix_basis(basis):
2031 """
2032 At least three of the five simple Euclidean Jordan algebras have the
2033 symmetric multiplication (A,B) |-> (AB + BA)/2, where the
2034 multiplication on the right is matrix multiplication. Given a basis
2035 for the underlying matrix space, this function returns a
2036 multiplication table (obtained by looping through the basis
2037 elements) for an algebra of those matrices. A reordered copy
2038 of the basis is also returned to work around the fact that
2039 the ``span()`` in this function will change the order of the basis
2040 from what we think it is, to... something else.
2041 """
2042 # In S^2, for example, we nominally have four coordinates even
2043 # though the space is of dimension three only. The vector space V
2044 # is supposed to hold the entire long vector, and the subspace W
2045 # of V will be spanned by the vectors that arise from symmetric
2046 # matrices. Thus for S^2, dim(V) == 4 and dim(W) == 3.
2047 field = basis[0].base_ring()
2048 dimension = basis[0].nrows()
2049
2050 V = VectorSpace(field, dimension**2)
2051 W = V.span( _mat2vec(s) for s in basis )
2052
2053 # Taking the span above reorders our basis (thanks, jerk!) so we
2054 # need to put our "matrix basis" in the same order as the
2055 # (reordered) vector basis.
2056 S = tuple( _vec2mat(b) for b in W.basis() )
2057
2058 Qs = []
2059 for s in S:
2060 # Brute force the multiplication-by-s matrix by looping
2061 # through all elements of the basis and doing the computation
2062 # to find out what the corresponding row should be. BEWARE:
2063 # these multiplication tables won't be symmetric! It therefore
2064 # becomes REALLY IMPORTANT that the underlying algebra
2065 # constructor uses ROW vectors and not COLUMN vectors. That's
2066 # why we're computing rows here and not columns.
2067 Q_rows = []
2068 for t in S:
2069 this_row = _mat2vec((s*t + t*s)/2)
2070 Q_rows.append(W.coordinates(this_row))
2071 Q = matrix(field, W.dimension(), Q_rows)
2072 Qs.append(Q)
2073
2074 return (Qs, S)
2075
2076
2077 def _embed_complex_matrix(M):
2078 """
2079 Embed the n-by-n complex matrix ``M`` into the space of real
2080 matrices of size 2n-by-2n via the map the sends each entry `z = a +
2081 bi` to the block matrix ``[[a,b],[-b,a]]``.
2082
2083 EXAMPLES::
2084
2085 sage: F = QuadraticField(-1,'i')
2086 sage: x1 = F(4 - 2*i)
2087 sage: x2 = F(1 + 2*i)
2088 sage: x3 = F(-i)
2089 sage: x4 = F(6)
2090 sage: M = matrix(F,2,[[x1,x2],[x3,x4]])
2091 sage: _embed_complex_matrix(M)
2092 [ 4 -2| 1 2]
2093 [ 2 4|-2 1]
2094 [-----+-----]
2095 [ 0 -1| 6 0]
2096 [ 1 0| 0 6]
2097
2098 TESTS:
2099
2100 Embedding is a homomorphism (isomorphism, in fact)::
2101
2102 sage: set_random_seed()
2103 sage: n = ZZ.random_element(5)
2104 sage: F = QuadraticField(-1, 'i')
2105 sage: X = random_matrix(F, n)
2106 sage: Y = random_matrix(F, n)
2107 sage: actual = _embed_complex_matrix(X) * _embed_complex_matrix(Y)
2108 sage: expected = _embed_complex_matrix(X*Y)
2109 sage: actual == expected
2110 True
2111
2112 """
2113 n = M.nrows()
2114 if M.ncols() != n:
2115 raise ValueError("the matrix 'M' must be square")
2116 field = M.base_ring()
2117 blocks = []
2118 for z in M.list():
2119 a = z.real()
2120 b = z.imag()
2121 blocks.append(matrix(field, 2, [[a,b],[-b,a]]))
2122
2123 # We can drop the imaginaries here.
2124 return block_matrix(field.base_ring(), n, blocks)
2125
2126
2127 def _unembed_complex_matrix(M):
2128 """
2129 The inverse of _embed_complex_matrix().
2130
2131 EXAMPLES::
2132
2133 sage: A = matrix(QQ,[ [ 1, 2, 3, 4],
2134 ....: [-2, 1, -4, 3],
2135 ....: [ 9, 10, 11, 12],
2136 ....: [-10, 9, -12, 11] ])
2137 sage: _unembed_complex_matrix(A)
2138 [ 2*i + 1 4*i + 3]
2139 [ 10*i + 9 12*i + 11]
2140
2141 TESTS:
2142
2143 Unembedding is the inverse of embedding::
2144
2145 sage: set_random_seed()
2146 sage: F = QuadraticField(-1, 'i')
2147 sage: M = random_matrix(F, 3)
2148 sage: _unembed_complex_matrix(_embed_complex_matrix(M)) == M
2149 True
2150
2151 """
2152 n = ZZ(M.nrows())
2153 if M.ncols() != n:
2154 raise ValueError("the matrix 'M' must be square")
2155 if not n.mod(2).is_zero():
2156 raise ValueError("the matrix 'M' must be a complex embedding")
2157
2158 F = QuadraticField(-1, 'i')
2159 i = F.gen()
2160
2161 # Go top-left to bottom-right (reading order), converting every
2162 # 2-by-2 block we see to a single complex element.
2163 elements = []
2164 for k in xrange(n/2):
2165 for j in xrange(n/2):
2166 submat = M[2*k:2*k+2,2*j:2*j+2]
2167 if submat[0,0] != submat[1,1]:
2168 raise ValueError('bad on-diagonal submatrix')
2169 if submat[0,1] != -submat[1,0]:
2170 raise ValueError('bad off-diagonal submatrix')
2171 z = submat[0,0] + submat[0,1]*i
2172 elements.append(z)
2173
2174 return matrix(F, n/2, elements)
2175
2176
2177 def _embed_quaternion_matrix(M):
2178 """
2179 Embed the n-by-n quaternion matrix ``M`` into the space of real
2180 matrices of size 4n-by-4n by first sending each quaternion entry
2181 `z = a + bi + cj + dk` to the block-complex matrix
2182 ``[[a + bi, c+di],[-c + di, a-bi]]`, and then embedding those into
2183 a real matrix.
2184
2185 EXAMPLES::
2186
2187 sage: Q = QuaternionAlgebra(QQ,-1,-1)
2188 sage: i,j,k = Q.gens()
2189 sage: x = 1 + 2*i + 3*j + 4*k
2190 sage: M = matrix(Q, 1, [[x]])
2191 sage: _embed_quaternion_matrix(M)
2192 [ 1 2 3 4]
2193 [-2 1 -4 3]
2194 [-3 4 1 -2]
2195 [-4 -3 2 1]
2196
2197 Embedding is a homomorphism (isomorphism, in fact)::
2198
2199 sage: set_random_seed()
2200 sage: n = ZZ.random_element(5)
2201 sage: Q = QuaternionAlgebra(QQ,-1,-1)
2202 sage: X = random_matrix(Q, n)
2203 sage: Y = random_matrix(Q, n)
2204 sage: actual = _embed_quaternion_matrix(X)*_embed_quaternion_matrix(Y)
2205 sage: expected = _embed_quaternion_matrix(X*Y)
2206 sage: actual == expected
2207 True
2208
2209 """
2210 quaternions = M.base_ring()
2211 n = M.nrows()
2212 if M.ncols() != n:
2213 raise ValueError("the matrix 'M' must be square")
2214
2215 F = QuadraticField(-1, 'i')
2216 i = F.gen()
2217
2218 blocks = []
2219 for z in M.list():
2220 t = z.coefficient_tuple()
2221 a = t[0]
2222 b = t[1]
2223 c = t[2]
2224 d = t[3]
2225 cplx_matrix = matrix(F, 2, [[ a + b*i, c + d*i],
2226 [-c + d*i, a - b*i]])
2227 blocks.append(_embed_complex_matrix(cplx_matrix))
2228
2229 # We should have real entries by now, so use the realest field
2230 # we've got for the return value.
2231 return block_matrix(quaternions.base_ring(), n, blocks)
2232
2233
2234 def _unembed_quaternion_matrix(M):
2235 """
2236 The inverse of _embed_quaternion_matrix().
2237
2238 EXAMPLES::
2239
2240 sage: M = matrix(QQ, [[ 1, 2, 3, 4],
2241 ....: [-2, 1, -4, 3],
2242 ....: [-3, 4, 1, -2],
2243 ....: [-4, -3, 2, 1]])
2244 sage: _unembed_quaternion_matrix(M)
2245 [1 + 2*i + 3*j + 4*k]
2246
2247 TESTS:
2248
2249 Unembedding is the inverse of embedding::
2250
2251 sage: set_random_seed()
2252 sage: Q = QuaternionAlgebra(QQ, -1, -1)
2253 sage: M = random_matrix(Q, 3)
2254 sage: _unembed_quaternion_matrix(_embed_quaternion_matrix(M)) == M
2255 True
2256
2257 """
2258 n = ZZ(M.nrows())
2259 if M.ncols() != n:
2260 raise ValueError("the matrix 'M' must be square")
2261 if not n.mod(4).is_zero():
2262 raise ValueError("the matrix 'M' must be a complex embedding")
2263
2264 Q = QuaternionAlgebra(QQ,-1,-1)
2265 i,j,k = Q.gens()
2266
2267 # Go top-left to bottom-right (reading order), converting every
2268 # 4-by-4 block we see to a 2-by-2 complex block, to a 1-by-1
2269 # quaternion block.
2270 elements = []
2271 for l in xrange(n/4):
2272 for m in xrange(n/4):
2273 submat = _unembed_complex_matrix(M[4*l:4*l+4,4*m:4*m+4])
2274 if submat[0,0] != submat[1,1].conjugate():
2275 raise ValueError('bad on-diagonal submatrix')
2276 if submat[0,1] != -submat[1,0].conjugate():
2277 raise ValueError('bad off-diagonal submatrix')
2278 z = submat[0,0].real() + submat[0,0].imag()*i
2279 z += submat[0,1].real()*j + submat[0,1].imag()*k
2280 elements.append(z)
2281
2282 return matrix(Q, n/4, elements)
2283
2284
2285 # The usual inner product on R^n.
2286 def _usual_ip(x,y):
2287 return x.vector().inner_product(y.vector())
2288
2289 # The inner product used for the real symmetric simple EJA.
2290 # We keep it as a separate function because e.g. the complex
2291 # algebra uses the same inner product, except divided by 2.
2292 def _matrix_ip(X,Y):
2293 X_mat = X.natural_representation()
2294 Y_mat = Y.natural_representation()
2295 return (X_mat*Y_mat).trace()
2296
2297
2298 class RealSymmetricEJA(FiniteDimensionalEuclideanJordanAlgebra):
2299 """
2300 The rank-n simple EJA consisting of real symmetric n-by-n
2301 matrices, the usual symmetric Jordan product, and the trace inner
2302 product. It has dimension `(n^2 + n)/2` over the reals.
2303
2304 EXAMPLES::
2305
2306 sage: J = RealSymmetricEJA(2)
2307 sage: e0, e1, e2 = J.gens()
2308 sage: e0*e0
2309 e0
2310 sage: e1*e1
2311 e0 + e2
2312 sage: e2*e2
2313 e2
2314
2315 TESTS:
2316
2317 The degree of this algebra is `(n^2 + n) / 2`::
2318
2319 sage: set_random_seed()
2320 sage: n = ZZ.random_element(1,5)
2321 sage: J = RealSymmetricEJA(n)
2322 sage: J.degree() == (n^2 + n)/2
2323 True
2324
2325 The Jordan multiplication is what we think it is::
2326
2327 sage: set_random_seed()
2328 sage: n = ZZ.random_element(1,5)
2329 sage: J = RealSymmetricEJA(n)
2330 sage: x = J.random_element()
2331 sage: y = J.random_element()
2332 sage: actual = (x*y).natural_representation()
2333 sage: X = x.natural_representation()
2334 sage: Y = y.natural_representation()
2335 sage: expected = (X*Y + Y*X)/2
2336 sage: actual == expected
2337 True
2338 sage: J(expected) == x*y
2339 True
2340
2341 """
2342 @staticmethod
2343 def __classcall_private__(cls, n, field=QQ):
2344 S = _real_symmetric_basis(n, field=field)
2345 (Qs, T) = _multiplication_table_from_matrix_basis(S)
2346
2347 fdeja = super(RealSymmetricEJA, cls)
2348 return fdeja.__classcall_private__(cls,
2349 field,
2350 Qs,
2351 rank=n,
2352 natural_basis=T)
2353
2354 def inner_product(self, x, y):
2355 return _matrix_ip(x,y)
2356
2357
2358 class ComplexHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
2359 """
2360 The rank-n simple EJA consisting of complex Hermitian n-by-n
2361 matrices over the real numbers, the usual symmetric Jordan product,
2362 and the real-part-of-trace inner product. It has dimension `n^2` over
2363 the reals.
2364
2365 TESTS:
2366
2367 The degree of this algebra is `n^2`::
2368
2369 sage: set_random_seed()
2370 sage: n = ZZ.random_element(1,5)
2371 sage: J = ComplexHermitianEJA(n)
2372 sage: J.degree() == n^2
2373 True
2374
2375 The Jordan multiplication is what we think it is::
2376
2377 sage: set_random_seed()
2378 sage: n = ZZ.random_element(1,5)
2379 sage: J = ComplexHermitianEJA(n)
2380 sage: x = J.random_element()
2381 sage: y = J.random_element()
2382 sage: actual = (x*y).natural_representation()
2383 sage: X = x.natural_representation()
2384 sage: Y = y.natural_representation()
2385 sage: expected = (X*Y + Y*X)/2
2386 sage: actual == expected
2387 True
2388 sage: J(expected) == x*y
2389 True
2390
2391 """
2392 @staticmethod
2393 def __classcall_private__(cls, n, field=QQ):
2394 S = _complex_hermitian_basis(n)
2395 (Qs, T) = _multiplication_table_from_matrix_basis(S)
2396
2397 fdeja = super(ComplexHermitianEJA, cls)
2398 return fdeja.__classcall_private__(cls,
2399 field,
2400 Qs,
2401 rank=n,
2402 natural_basis=T)
2403
2404 def inner_product(self, x, y):
2405 # Since a+bi on the diagonal is represented as
2406 #
2407 # a + bi = [ a b ]
2408 # [ -b a ],
2409 #
2410 # we'll double-count the "a" entries if we take the trace of
2411 # the embedding.
2412 return _matrix_ip(x,y)/2
2413
2414
2415 class QuaternionHermitianEJA(FiniteDimensionalEuclideanJordanAlgebra):
2416 """
2417 The rank-n simple EJA consisting of self-adjoint n-by-n quaternion
2418 matrices, the usual symmetric Jordan product, and the
2419 real-part-of-trace inner product. It has dimension `2n^2 - n` over
2420 the reals.
2421
2422 TESTS:
2423
2424 The degree of this algebra is `n^2`::
2425
2426 sage: set_random_seed()
2427 sage: n = ZZ.random_element(1,5)
2428 sage: J = QuaternionHermitianEJA(n)
2429 sage: J.degree() == 2*(n^2) - n
2430 True
2431
2432 The Jordan multiplication is what we think it is::
2433
2434 sage: set_random_seed()
2435 sage: n = ZZ.random_element(1,5)
2436 sage: J = QuaternionHermitianEJA(n)
2437 sage: x = J.random_element()
2438 sage: y = J.random_element()
2439 sage: actual = (x*y).natural_representation()
2440 sage: X = x.natural_representation()
2441 sage: Y = y.natural_representation()
2442 sage: expected = (X*Y + Y*X)/2
2443 sage: actual == expected
2444 True
2445 sage: J(expected) == x*y
2446 True
2447
2448 """
2449 @staticmethod
2450 def __classcall_private__(cls, n, field=QQ):
2451 S = _quaternion_hermitian_basis(n)
2452 (Qs, T) = _multiplication_table_from_matrix_basis(S)
2453
2454 fdeja = super(QuaternionHermitianEJA, cls)
2455 return fdeja.__classcall_private__(cls,
2456 field,
2457 Qs,
2458 rank=n,
2459 natural_basis=T)
2460
2461 def inner_product(self, x, y):
2462 # Since a+bi+cj+dk on the diagonal is represented as
2463 #
2464 # a + bi +cj + dk = [ a b c d]
2465 # [ -b a -d c]
2466 # [ -c d a -b]
2467 # [ -d -c b a],
2468 #
2469 # we'll quadruple-count the "a" entries if we take the trace of
2470 # the embedding.
2471 return _matrix_ip(x,y)/4
2472
2473
2474 class JordanSpinEJA(FiniteDimensionalEuclideanJordanAlgebra):
2475 """
2476 The rank-2 simple EJA consisting of real vectors ``x=(x0, x_bar)``
2477 with the usual inner product and jordan product ``x*y =
2478 (<x_bar,y_bar>, x0*y_bar + y0*x_bar)``. It has dimension `n` over
2479 the reals.
2480
2481 EXAMPLES:
2482
2483 This multiplication table can be verified by hand::
2484
2485 sage: J = JordanSpinEJA(4)
2486 sage: e0,e1,e2,e3 = J.gens()
2487 sage: e0*e0
2488 e0
2489 sage: e0*e1
2490 e1
2491 sage: e0*e2
2492 e2
2493 sage: e0*e3
2494 e3
2495 sage: e1*e2
2496 0
2497 sage: e1*e3
2498 0
2499 sage: e2*e3
2500 0
2501
2502 """
2503 @staticmethod
2504 def __classcall_private__(cls, n, field=QQ):
2505 Qs = []
2506 id_matrix = identity_matrix(field, n)
2507 for i in xrange(n):
2508 ei = id_matrix.column(i)
2509 Qi = zero_matrix(field, n)
2510 Qi.set_row(0, ei)
2511 Qi.set_column(0, ei)
2512 Qi += diagonal_matrix(n, [ei[0]]*n)
2513 # The addition of the diagonal matrix adds an extra ei[0] in the
2514 # upper-left corner of the matrix.
2515 Qi[0,0] = Qi[0,0] * ~field(2)
2516 Qs.append(Qi)
2517
2518 # The rank of the spin algebra is two, unless we're in a
2519 # one-dimensional ambient space (because the rank is bounded by
2520 # the ambient dimension).
2521 fdeja = super(JordanSpinEJA, cls)
2522 return fdeja.__classcall_private__(cls, field, Qs, rank=min(n,2))
2523
2524 def inner_product(self, x, y):
2525 return _usual_ip(x,y)