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