X-Git-Url: http://gitweb.michael.orlitzky.com/?a=blobdiff_plain;f=mjo%2Feja%2Feja_element.py;h=f4d5995ce9302d2fdc518dac2b0ca20f0fadf6d1;hb=HEAD;hp=0953b2f27d67d059e88588fc9adb1f3081fdecc6;hpb=ffabf130fc8fa05e68d9b3d55a1f3ccf89f5e390;p=sage.d.git diff --git a/mjo/eja/eja_element.py b/mjo/eja/eja_element.py index 0953b2f..f4d5995 100644 --- a/mjo/eja/eja_element.py +++ b/mjo/eja/eja_element.py @@ -1,18 +1,13 @@ -# -*- coding: utf-8 -*- - from sage.matrix.constructor import matrix +from sage.misc.cachefunc import cached_method from sage.modules.free_module import VectorSpace from sage.modules.with_basis.indexed_element import IndexedFreeModuleElement -# TODO: make this unnecessary somehow. -from sage.misc.lazy_import import lazy_import -lazy_import('mjo.eja.eja_algebra', 'FiniteDimensionalEuclideanJordanAlgebra') -lazy_import('mjo.eja.eja_element_subalgebra', - 'FiniteDimensionalEuclideanJordanElementSubalgebra') -from mjo.eja.eja_operator import FiniteDimensionalEuclideanJordanAlgebraOperator -from mjo.eja.eja_utils import _mat2vec +from mjo.eja.eja_operator import EJAOperator +from mjo.eja.eja_utils import _scale + -class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): +class EJAElement(IndexedFreeModuleElement): """ An element of a Euclidean Jordan algebra. """ @@ -48,14 +43,12 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The definition of `x^2` is the unambiguous `x*x`:: - sage: set_random_seed() sage: x = random_eja().random_element() sage: x*x == (x^2) True A few examples of power-associativity:: - sage: set_random_seed() sage: x = random_eja().random_element() sage: x*(x*x)*(x*x) == x^5 True @@ -65,7 +58,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): We also know that powers operator-commute (Koecher, Chapter III, Corollary 1):: - sage: set_random_seed() sage: x = random_eja().random_element() sage: m = ZZ.random_element(0,10) sage: n = ZZ.random_element(0,10) @@ -112,7 +104,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): We should always get back an element of the algebra:: - sage: set_random_seed() sage: p = PolynomialRing(AA, 't').random_element() sage: J = random_eja() sage: x = J.random_element() @@ -137,7 +128,8 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): SETUP:: - sage: from mjo.eja.eja_algebra import HadamardEJA + sage: from mjo.eja.eja_algebra import (random_eja, + ....: HadamardEJA) EXAMPLES: @@ -161,11 +153,10 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The characteristic polynomial of an element should evaluate to zero on that element:: - sage: set_random_seed() - sage: x = HadamardEJA(3).random_element() + sage: x = random_eja().random_element() sage: p = x.characteristic_polynomial() - sage: x.apply_univariate_polynomial(p) - 0 + sage: x.apply_univariate_polynomial(p).is_zero() + True The characteristic polynomials of the zero and unit elements should be what we think they are in a subalgebra, too:: @@ -173,8 +164,8 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): sage: J = HadamardEJA(3) sage: p1 = J.one().characteristic_polynomial() sage: q1 = J.zero().characteristic_polynomial() - sage: e0,e1,e2 = J.gens() - sage: A = (e0 + 2*e1 + 3*e2).subalgebra_generated_by() # dim 3 + sage: b0,b1,b2 = J.gens() + sage: A = (b0 + 2*b1 + 3*b2).subalgebra_generated_by() # dim 3 sage: p2 = A.one().characteristic_polynomial() sage: q2 = A.zero().characteristic_polynomial() sage: p1 == p2 @@ -183,7 +174,7 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): True """ - p = self.parent().characteristic_polynomial() + p = self.parent().characteristic_polynomial_of() return p(*self.to_vector()) @@ -234,16 +225,15 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): Ditto for the quaternions:: - sage: J = QuaternionHermitianEJA(3) + sage: J = QuaternionHermitianEJA(2) sage: J.one().inner_product(J.one()) - 3 + 2 TESTS: Ensure that we can always compute an inner product, and that it gives us back a real number:: - sage: set_random_seed() sage: J = random_eja() sage: x,y = J.random_elements(2) sage: x.inner_product(y) in RLF @@ -271,7 +261,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The definition of a Jordan algebra says that any element operator-commutes with its square:: - sage: set_random_seed() sage: x = random_eja().random_element() sage: x.operator_commutes_with(x^2) True @@ -280,7 +269,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): Test Lemma 1 from Chapter III of Koecher:: - sage: set_random_seed() sage: u,v = random_eja().random_elements(2) sage: lhs = u.operator_commutes_with(u*v) sage: rhs = v.operator_commutes_with(u^2) @@ -290,7 +278,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): Test the first polarization identity from my notes, Koecher Chapter III, or from Baes (2.3):: - sage: set_random_seed() sage: x,y = random_eja().random_elements(2) sage: Lx = x.operator() sage: Ly = y.operator() @@ -302,32 +289,32 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): Test the second polarization identity from my notes or from Baes (2.4):: - sage: set_random_seed() - sage: x,y,z = random_eja().random_elements(3) - sage: Lx = x.operator() - sage: Ly = y.operator() - sage: Lz = z.operator() - sage: Lzy = (z*y).operator() - sage: Lxy = (x*y).operator() - sage: Lxz = (x*z).operator() - sage: bool(Lx*Lzy + Lz*Lxy + Ly*Lxz == Lzy*Lx + Lxy*Lz + Lxz*Ly) + sage: x,y,z = random_eja().random_elements(3) # long time + sage: Lx = x.operator() # long time + sage: Ly = y.operator() # long time + sage: Lz = z.operator() # long time + sage: Lzy = (z*y).operator() # long time + sage: Lxy = (x*y).operator() # long time + sage: Lxz = (x*z).operator() # long time + sage: lhs = Lx*Lzy + Lz*Lxy + Ly*Lxz # long time + sage: rhs = Lzy*Lx + Lxy*Lz + Lxz*Ly # long time + sage: bool(lhs == rhs) # long time True Test the third polarization identity from my notes or from Baes (2.5):: - sage: set_random_seed() - sage: u,y,z = random_eja().random_elements(3) - sage: Lu = u.operator() - sage: Ly = y.operator() - sage: Lz = z.operator() - sage: Lzy = (z*y).operator() - sage: Luy = (u*y).operator() - sage: Luz = (u*z).operator() - sage: Luyz = (u*(y*z)).operator() - sage: lhs = Lu*Lzy + Lz*Luy + Ly*Luz - sage: rhs = Luyz + Ly*Lu*Lz + Lz*Lu*Ly - sage: bool(lhs == rhs) + sage: u,y,z = random_eja().random_elements(3) # long time + sage: Lu = u.operator() # long time + sage: Ly = y.operator() # long time + sage: Lz = z.operator() # long time + sage: Lzy = (z*y).operator() # long time + sage: Luy = (u*y).operator() # long time + sage: Luz = (u*z).operator() # long time + sage: Luyz = (u*(y*z)).operator() # long time + sage: lhs = Lu*Lzy + Lz*Luy + Ly*Luz # long time + sage: rhs = Luyz + Ly*Lu*Lz + Lz*Lu*Ly # long time + sage: bool(lhs == rhs) # long time True """ @@ -345,14 +332,16 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): SETUP:: - sage: from mjo.eja.eja_algebra import (JordanSpinEJA, + sage: from mjo.eja.eja_algebra import (AlbertEJA, + ....: JordanSpinEJA, ....: TrivialEJA, + ....: RealSymmetricEJA, + ....: ComplexHermitianEJA, ....: random_eja) EXAMPLES:: sage: J = JordanSpinEJA(2) - sage: e0,e1 = J.gens() sage: x = sum( J.gens() ) sage: x.det() 0 @@ -360,7 +349,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): :: sage: J = JordanSpinEJA(3) - sage: e0,e1,e2 = J.gens() sage: x = sum( J.gens() ) sage: x.det() -1 @@ -381,7 +369,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): An element is invertible if and only if its determinant is non-zero:: - sage: set_random_seed() sage: x = random_eja().random_element() sage: x.is_invertible() == (x.det() != 0) True @@ -389,11 +376,42 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): Ensure that the determinant is multiplicative on an associative subalgebra as in Faraut and Korányi's Proposition II.2.2:: - sage: set_random_seed() - sage: J = random_eja().random_element().subalgebra_generated_by() + sage: x0 = random_eja().random_element() + sage: J = x0.subalgebra_generated_by(orthonormalize=False) sage: x,y = J.random_elements(2) sage: (x*y).det() == x.det()*y.det() True + + The determinant in real matrix algebras is the usual determinant:: + + sage: X = matrix.random(QQ,3) + sage: X = X + X.T + sage: J1 = RealSymmetricEJA(3) + sage: J2 = RealSymmetricEJA(3,field=QQ,orthonormalize=False) + sage: expected = X.det() + sage: actual1 = J1(X).det() + sage: actual2 = J2(X).det() + sage: actual1 == expected + True + sage: actual2 == expected + True + + There's a formula for the determinant of the Albert algebra + (Yokota, Section 2.1):: + + sage: def albert_det(x): + ....: X = x.to_matrix() + ....: res = X[0,0]*X[1,1]*X[2,2] + ....: res += 2*(X[1,2]*X[2,0]*X[0,1]).real() + ....: res -= X[0,0]*X[1,2]*X[2,1] + ....: res -= X[1,1]*X[2,0]*X[0,2] + ....: res -= X[2,2]*X[0,1]*X[1,0] + ....: return res.leading_coefficient() + sage: J = AlbertEJA(field=QQ, orthonormalize=False) + sage: xs = J.random_elements(10) + sage: all( albert_det(x) == x.det() for x in xs ) + True + """ P = self.parent() r = P.rank() @@ -412,14 +430,20 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): return ((-1)**r)*p(*self.to_vector()) + @cached_method def inverse(self): """ Return the Jordan-multiplicative inverse of this element. ALGORITHM: - We appeal to the quadratic representation as in Koecher's - Theorem 12 in Chapter III, Section 5. + In general we appeal to the quadratic representation as in + Koecher's Theorem 12 in Chapter III, Section 5. But if the + parent algebra's "characteristic polynomial of" coefficients + happen to be cached, then we use Proposition II.2.4 in Faraut + and Korányi which gives a formula for the inverse based on the + characteristic polynomial and the Cayley-Hamilton theorem for + Euclidean Jordan algebras:: SETUP:: @@ -432,17 +456,16 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The inverse in the spin factor algebra is given in Alizadeh's Example 11.11:: - sage: set_random_seed() sage: J = JordanSpinEJA.random_instance() sage: x = J.random_element() sage: while not x.is_invertible(): ....: x = J.random_element() sage: x_vec = x.to_vector() - sage: x0 = x_vec[0] + sage: x0 = x_vec[:1] sage: x_bar = x_vec[1:] - sage: coeff = ~(x0^2 - x_bar.inner_product(x_bar)) - sage: inv_vec = x_vec.parent()([x0] + (-x_bar).list()) - sage: x_inverse = coeff*inv_vec + sage: coeff = x0.inner_product(x0) - x_bar.inner_product(x_bar) + sage: x_inverse = x_vec.parent()(x0.list() + (-x_bar).list()) + sage: if not coeff.is_zero(): x_inverse = x_inverse/coeff sage: x.inverse() == J.from_vector(x_inverse) True @@ -451,20 +474,18 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): sage: JordanSpinEJA(3).zero().inverse() Traceback (most recent call last): ... - ValueError: element is not invertible + ZeroDivisionError: element is not invertible TESTS: The identity element is its own inverse:: - sage: set_random_seed() sage: J = random_eja() sage: J.one().inverse() == J.one() True If an element has an inverse, it acts like one:: - sage: set_random_seed() sage: J = random_eja() sage: x = J.random_element() sage: (not x.is_invertible()) or (x.inverse()*x == J.one()) @@ -472,7 +493,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The inverse of the inverse is what we started with:: - sage: set_random_seed() sage: J = random_eja() sage: x = J.random_element() sage: (not x.is_invertible()) or (x.inverse().inverse() == x) @@ -482,51 +502,62 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): of an element is the inverse of its left-multiplication operator applied to the algebra's identity, when that inverse exists:: - sage: set_random_seed() - sage: J = random_eja() - sage: x = J.random_element() - sage: (not x.operator().is_invertible()) or ( - ....: x.operator().inverse()(J.one()) == x.inverse() ) + sage: J = random_eja() # long time + sage: x = J.random_element() # long time + sage: (not x.operator().is_invertible()) or ( # long time + ....: x.operator().inverse()(J.one()) # long time + ....: == # long time + ....: x.inverse() ) # long time True - Proposition II.2.4 in Faraut and Korányi gives a formula for - the inverse based on the characteristic polynomial and the - Cayley-Hamilton theorem for Euclidean Jordan algebras:: + Check that the fast (cached) and slow algorithms give the same + answer:: - sage: set_random_seed() - sage: J = ComplexHermitianEJA(3) - sage: x = J.random_element() - sage: while not x.is_invertible(): - ....: x = J.random_element() - sage: r = J.rank() - sage: a = x.characteristic_polynomial().coefficients(sparse=False) - sage: expected = (-1)^(r+1)/x.det() - sage: expected *= sum( a[i+1]*x^i for i in range(r) ) - sage: x.inverse() == expected + sage: J = random_eja(field=QQ, orthonormalize=False) # long time + sage: x = J.random_element() # long time + sage: while not x.is_invertible(): # long time + ....: x = J.random_element() # long time + sage: slow = x.inverse() # long time + sage: _ = J._charpoly_coefficients() # long time + sage: fast = x.inverse() # long time + sage: slow == fast # long time True - """ - if not self.is_invertible(): - raise ValueError("element is not invertible") - - return (~self.quadratic_representation())(self) - - + not_invertible_msg = "element is not invertible" + + algebra = self.parent() + if algebra._charpoly_coefficients.is_in_cache(): + # We can invert using our charpoly if it will be fast to + # compute. If the coefficients are cached, our rank had + # better be too! + if self.det().is_zero(): + raise ZeroDivisionError(not_invertible_msg) + r = algebra.rank() + a = self.characteristic_polynomial().coefficients(sparse=False) + return (-1)**(r+1)*algebra.sum(a[i+1]*self**i + for i in range(r))/self.det() + + try: + inv = (~self.quadratic_representation())(self) + self.is_invertible.set_cache(True) + return inv + except ZeroDivisionError: + self.is_invertible.set_cache(False) + raise ZeroDivisionError(not_invertible_msg) + + + @cached_method def is_invertible(self): """ Return whether or not this element is invertible. ALGORITHM: - The usual way to do this is to check if the determinant is - zero, but we need the characteristic polynomial for the - determinant. The minimal polynomial is a lot easier to get, - so we use Corollary 2 in Chapter V of Koecher to check - whether or not the paren't algebra's zero element is a root - of this element's minimal polynomial. - - Beware that we can't use the superclass method, because it - relies on the algebra being associative. + If computing my determinant will be fast, we do so and compare + with zero (Proposition II.2.4 in Faraut and + Korányi). Otherwise, Proposition II.3.2 in Faraut and Korányi + reduces the problem to the invertibility of my quadratic + representation. SETUP:: @@ -536,18 +567,26 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The identity element is always invertible:: - sage: set_random_seed() sage: J = random_eja() sage: J.one().is_invertible() True The zero element is never invertible in a non-trivial algebra:: - sage: set_random_seed() sage: J = random_eja() sage: (not J.is_trivial()) and J.zero().is_invertible() False + Test that the fast (cached) and slow algorithms give the same + answer:: + + sage: J = random_eja(field=QQ, orthonormalize=False) # long time + sage: x = J.random_element() # long time + sage: slow = x.is_invertible() # long time + sage: _ = J._charpoly_coefficients() # long time + sage: fast = x.is_invertible() # long time + sage: slow == fast # long time + True """ if self.is_zero(): if self.parent().is_trivial(): @@ -555,11 +594,18 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): else: return False - # In fact, we only need to know if the constant term is non-zero, - # so we can pass in the field's zero element instead. - zero = self.base_ring().zero() - p = self.minimal_polynomial() - return not (p(zero) == zero) + if self.parent()._charpoly_coefficients.is_in_cache(): + # The determinant will be quicker than inverting the + # quadratic representation, most likely. + return (not self.det().is_zero()) + + # The easiest way to determine if I'm invertible is to try. + try: + inv = (~self.quadratic_representation())(self) + self.inverse.set_cache(inv) + return True + except ZeroDivisionError: + return False def is_primitive_idempotent(self): @@ -605,7 +651,7 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): element should always be in terms of minimal idempotents:: sage: J = JordanSpinEJA(4) - sage: x = sum( i*J.gens()[i] for i in range(len(J.gens())) ) + sage: x = sum( i*J.monomial(i) for i in range(len(J.gens())) ) sage: x.is_regular() True sage: [ c.is_primitive_idempotent() @@ -616,14 +662,12 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The identity element is minimal only in an EJA of rank one:: - sage: set_random_seed() sage: J = random_eja() sage: J.rank() == 1 or not J.one().is_primitive_idempotent() True A non-idempotent cannot be a minimal idempotent:: - sage: set_random_seed() sage: J = JordanSpinEJA(4) sage: x = J.random_element() sage: (not x.is_idempotent()) and x.is_primitive_idempotent() @@ -633,7 +677,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): idempotent if and only if it's idempotent with trace equal to unity:: - sage: set_random_seed() sage: J = JordanSpinEJA(4) sage: x = J.random_element() sage: expected = (x.is_idempotent() and x.trace() == 1) @@ -643,7 +686,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): Primitive idempotents must be non-zero:: - sage: set_random_seed() sage: J = random_eja() sage: J.zero().is_idempotent() True @@ -700,14 +742,12 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The identity element is never nilpotent, except in a trivial EJA:: - sage: set_random_seed() sage: J = random_eja() sage: J.one().is_nilpotent() and not J.is_trivial() False The additive identity is always nilpotent:: - sage: set_random_seed() sage: random_eja().zero().is_nilpotent() True @@ -734,7 +774,9 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): sage: J = JordanSpinEJA(5) sage: J.one().is_regular() False - sage: e0, e1, e2, e3, e4 = J.gens() # e0 is the identity + sage: b0, b1, b2, b3, b4 = J.gens() + sage: b0 == J.one() + True sage: for x in J.gens(): ....: (J.one() + x).is_regular() False @@ -748,7 +790,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The zero element should never be regular, unless the parent algebra has dimension less than or equal to one:: - sage: set_random_seed() sage: J = random_eja() sage: J.dimension() <= 1 or not J.zero().is_regular() True @@ -756,7 +797,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The unit element isn't regular unless the algebra happens to consist of only its scalar multiples:: - sage: set_random_seed() sage: J = random_eja() sage: J.dimension() <= 1 or not J.one().is_regular() True @@ -772,10 +812,23 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): ALGORITHM: - For now, we skip the messy minimal polynomial computation - and instead return the dimension of the vector space spanned - by the powers of this element. The latter is a bit more - straightforward to compute. + First we handle the special cases where the algebra is + trivial, this element is zero, or the dimension of the algebra + is one and this element is not zero. With those out of the + way, we may assume that ``self`` is nonzero, the algebra is + nontrivial, and that the dimension of the algebra is at least + two. + + Beginning with the algebra's unit element (power zero), we add + successive (basis representations of) powers of this element + to a matrix, row-reducing at each step. After row-reducing, we + check the rank of the matrix. If adding a row and row-reducing + does not increase the rank of the matrix at any point, the row + we've just added lives in the span of the previous ones; thus + the corresponding power of ``self`` lives in the span of its + lesser powers. When that happens, the degree of the minimal + polynomial is the rank of the matrix; if it never happens, the + degree must be the dimension of the entire space. SETUP:: @@ -787,17 +840,17 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): sage: J = JordanSpinEJA(4) sage: J.one().degree() 1 - sage: e0,e1,e2,e3 = J.gens() - sage: (e0 - e1).degree() + sage: b0,b1,b2,b3 = J.gens() + sage: (b0 - b1).degree() 2 In the spin factor algebra (of rank two), all elements that aren't multiples of the identity are regular:: - sage: set_random_seed() sage: J = JordanSpinEJA.random_instance() + sage: n = J.dimension() sage: x = J.random_element() - sage: x == x.coefficient(0)*J.one() or x.degree() == 2 + sage: x.degree() == min(n,2) or (x == x.coefficient(0)*J.one()) True TESTS: @@ -805,7 +858,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The zero and unit elements are both of degree one in nontrivial algebras:: - sage: set_random_seed() sage: J = random_eja() sage: d = J.zero().degree() sage: (J.is_trivial() and d == 0) or d == 1 @@ -816,18 +868,63 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): Our implementation agrees with the definition:: - sage: set_random_seed() sage: x = random_eja().random_element() sage: x.degree() == x.minimal_polynomial().degree() True - """ - if self.is_zero() and not self.parent().is_trivial(): + n = self.parent().dimension() + + if n == 0: + # The minimal polynomial is an empty product, i.e. the + # constant polynomial "1" having degree zero. + return 0 + elif self.is_zero(): # The minimal polynomial of zero in a nontrivial algebra - # is "t"; in a trivial algebra it's "1" by convention - # (it's an empty product). + # is "t", and is of degree one. return 1 - return self.subalgebra_generated_by().dimension() + elif n == 1: + # If this is a nonzero element of a nontrivial algebra, it + # has degree at least one. It follows that, in an algebra + # of dimension one, the degree must be actually one. + return 1 + + # BEWARE: The subalgebra_generated_by() method uses the result + # of this method to construct a basis for the subalgebra. That + # means, in particular, that we cannot implement this method + # as ``self.subalgebra_generated_by().dimension()``. + + # Algorithm: keep appending (vector representations of) powers + # self as rows to a matrix and echelonizing it. When its rank + # stops increasing, we've reached a redundancy. + + # Given the special cases above, we can assume that "self" is + # nonzero, the algebra is nontrivial, and that its dimension + # is at least two. + M = matrix([(self.parent().one()).to_vector()]) + old_rank = 1 + + # Specifying the row-reduction algorithm can e.g. help over + # AA because it avoids the RecursionError that gets thrown + # when we have to look too hard for a root. + # + # Beware: QQ supports an entirely different set of "algorithm" + # keywords than do AA and RR. + algo = None + from sage.rings.all import QQ + if self.parent().base_ring() is not QQ: + algo = "scaled_partial_pivoting" + + for d in range(1,n): + M = matrix(M.rows() + [(self**d).to_vector()]) + M.echelonize(algo) + new_rank = M.rank() + if new_rank == old_rank: + return new_rank + else: + old_rank = new_rank + + return n + def left_matrix(self): @@ -875,19 +972,22 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): TESTS: The minimal polynomial of the identity and zero elements are - always the same:: + always the same, except in trivial algebras where the minimal + polynomial of the unit/zero element is ``1``:: - sage: set_random_seed() - sage: J = random_eja(nontrivial=True) - sage: J.one().minimal_polynomial() + sage: J = random_eja() + sage: mu = J.one().minimal_polynomial() + sage: t = mu.parent().gen() + sage: mu + int(J.is_trivial())*(t-2) t - 1 - sage: J.zero().minimal_polynomial() + sage: mu = J.zero().minimal_polynomial() + sage: t = mu.parent().gen() + sage: mu + int(J.is_trivial())*(t-1) t The degree of an element is (by one definition) the degree of its minimal polynomial:: - sage: set_random_seed() sage: x = random_eja().random_element() sage: x.degree() == x.minimal_polynomial().degree() True @@ -898,9 +998,8 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): identity. We require the dimension of the algebra to be at least two here so that said elements actually exist:: - sage: set_random_seed() - sage: n_max = max(2, JordanSpinEJA._max_test_case_size()) - sage: n = ZZ.random_element(2, n_max) + sage: d_max = JordanSpinEJA._max_random_instance_dimension() + sage: n = ZZ.random_element(2, max(2,d_max)) sage: J = JordanSpinEJA(n) sage: y = J.random_element() sage: while y == y.coefficient(0)*J.one(): @@ -915,20 +1014,19 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The minimal polynomial should always kill its element:: - sage: set_random_seed() - sage: x = random_eja().random_element() - sage: p = x.minimal_polynomial() - sage: x.apply_univariate_polynomial(p) + sage: x = random_eja().random_element() # long time + sage: p = x.minimal_polynomial() # long time + sage: x.apply_univariate_polynomial(p) # long time 0 The minimal polynomial is invariant under a change of basis, and in particular, a re-scaling of the basis:: - sage: set_random_seed() - sage: n_max = RealSymmetricEJA._max_test_case_size() - sage: n = ZZ.random_element(1, n_max) + sage: d_max = RealSymmetricEJA._max_random_instance_dimension() + sage: d = ZZ.random_element(1, d_max) + sage: n = RealSymmetricEJA._max_random_instance_size(d) sage: J1 = RealSymmetricEJA(n) - sage: J2 = RealSymmetricEJA(n,normalize_basis=False) + sage: J2 = RealSymmetricEJA(n,orthonormalize=False) sage: X = random_matrix(AA,n) sage: X = X*X.transpose() sage: x1 = J1(X) @@ -938,73 +1036,98 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): """ if self.is_zero(): - # We would generate a zero-dimensional subalgebra - # where the minimal polynomial would be constant. - # That might be correct, but only if *this* algebra - # is trivial too. - if not self.parent().is_trivial(): - # Pretty sure we know what the minimal polynomial of - # the zero operator is going to be. This ensures - # consistency of e.g. the polynomial variable returned - # in the "normal" case without us having to think about it. - return self.operator().minimal_polynomial() - - A = self.subalgebra_generated_by() - return A(self).operator().minimal_polynomial() + # Pretty sure we know what the minimal polynomial of + # the zero operator is going to be. This ensures + # consistency of e.g. the polynomial variable returned + # in the "normal" case without us having to think about it. + return self.operator().minimal_polynomial() + + # If we don't orthonormalize the subalgebra's basis, then the + # first two monomials in the subalgebra will be self^0 and + # self^1... assuming that self^1 is not a scalar multiple of + # self^0 (the unit element). We special case these to avoid + # having to solve a system to coerce self into the subalgebra. + A = self.subalgebra_generated_by(orthonormalize=False) + + if A.dimension() == 1: + # Does a solve to find the scalar multiple alpha such that + # alpha*unit = self. We have to do this because the basis + # for the subalgebra will be [ self^0 ], and not [ self^1 ]! + unit = self.parent().one() + alpha = self.to_vector() / unit.to_vector() + return (unit.operator()*alpha).minimal_polynomial() + else: + # If the dimension of the subalgebra is >= 2, then we just + # use the second basis element. + return A.monomial(1).operator().minimal_polynomial() - def natural_representation(self): + def to_matrix(self): """ - Return a more-natural representation of this element. + Return an (often more natural) representation of this element as a + matrix. - Every finite-dimensional Euclidean Jordan Algebra is a - direct sum of five simple algebras, four of which comprise - Hermitian matrices. This method returns the original - "natural" representation of this element as a Hermitian - matrix, if it has one. If not, you get the usual representation. + Every finite-dimensional Euclidean Jordan Algebra is a direct + sum of five simple algebras, four of which comprise Hermitian + matrices. This method returns a "natural" matrix + representation of this element as either a Hermitian matrix or + column vector. SETUP:: sage: from mjo.eja.eja_algebra import (ComplexHermitianEJA, - ....: QuaternionHermitianEJA) + ....: HadamardEJA, + ....: QuaternionHermitianEJA, + ....: RealSymmetricEJA) EXAMPLES:: sage: J = ComplexHermitianEJA(3) sage: J.one() - e0 + e3 + e8 - sage: J.one().natural_representation() - [1 0 0 0 0 0] - [0 1 0 0 0 0] - [0 0 1 0 0 0] - [0 0 0 1 0 0] - [0 0 0 0 1 0] - [0 0 0 0 0 1] + b0 + b3 + b8 + sage: J.one().to_matrix() + +---+---+---+ + | 1 | 0 | 0 | + +---+---+---+ + | 0 | 1 | 0 | + +---+---+---+ + | 0 | 0 | 1 | + +---+---+---+ :: - sage: J = QuaternionHermitianEJA(3) + sage: J = QuaternionHermitianEJA(2) sage: J.one() - e0 + e5 + e14 - sage: J.one().natural_representation() - [1 0 0 0 0 0 0 0 0 0 0 0] - [0 1 0 0 0 0 0 0 0 0 0 0] - [0 0 1 0 0 0 0 0 0 0 0 0] - [0 0 0 1 0 0 0 0 0 0 0 0] - [0 0 0 0 1 0 0 0 0 0 0 0] - [0 0 0 0 0 1 0 0 0 0 0 0] - [0 0 0 0 0 0 1 0 0 0 0 0] - [0 0 0 0 0 0 0 1 0 0 0 0] - [0 0 0 0 0 0 0 0 1 0 0 0] - [0 0 0 0 0 0 0 0 0 1 0 0] - [0 0 0 0 0 0 0 0 0 0 1 0] - [0 0 0 0 0 0 0 0 0 0 0 1] + b0 + b5 + sage: J.one().to_matrix() + +---+---+ + | 1 | 0 | + +---+---+ + | 0 | 1 | + +---+---+ + + This also works in Cartesian product algebras:: + + sage: J1 = HadamardEJA(1) + sage: J2 = RealSymmetricEJA(2) + sage: J = cartesian_product([J1,J2]) + sage: x = sum(J.gens()) + sage: x.to_matrix()[0] + [1] + sage: x.to_matrix()[1] + [ 1 0.7071067811865475?] + [0.7071067811865475? 1] """ - B = self.parent().natural_basis() - W = self.parent().natural_basis_space() - return W.linear_combination(zip(B,self.to_vector())) + B = self.parent().matrix_basis() + W = self.parent().matrix_space() + + # This is just a manual "from_vector()", but of course + # matrix spaces aren't vector spaces in sage, so they + # don't have a from_vector() method. + return W.linear_combination( zip(B, self.to_vector()) ) + def norm(self): @@ -1045,7 +1168,6 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): TESTS:: - sage: set_random_seed() sage: J = random_eja() sage: x,y = J.random_elements(2) sage: x.operator()(y) == x*y @@ -1057,10 +1179,7 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): P = self.parent() left_mult_by_self = lambda y: self*y L = P.module_morphism(function=left_mult_by_self, codomain=P) - return FiniteDimensionalEuclideanJordanAlgebraOperator( - P, - P, - L.matrix() ) + return EJAOperator(P, P, L.matrix() ) def quadratic_representation(self, other=None): @@ -1077,25 +1196,25 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The explicit form in the spin factor algebra is given by Alizadeh's Example 11.12:: - sage: set_random_seed() sage: x = JordanSpinEJA.random_instance().random_element() sage: x_vec = x.to_vector() + sage: Q = matrix.identity(x.base_ring(), 0) sage: n = x_vec.degree() - sage: x0 = x_vec[0] - sage: x_bar = x_vec[1:] - sage: A = matrix(AA, 1, [x_vec.inner_product(x_vec)]) - sage: B = 2*x0*x_bar.row() - sage: C = 2*x0*x_bar.column() - sage: D = matrix.identity(AA, n-1) - sage: D = (x0^2 - x_bar.inner_product(x_bar))*D - sage: D = D + 2*x_bar.tensor_product(x_bar) - sage: Q = matrix.block(2,2,[A,B,C,D]) + sage: if n > 0: + ....: x0 = x_vec[0] + ....: x_bar = x_vec[1:] + ....: A = matrix(x.base_ring(), 1, [x_vec.inner_product(x_vec)]) + ....: B = 2*x0*x_bar.row() + ....: C = 2*x0*x_bar.column() + ....: D = matrix.identity(x.base_ring(), n-1) + ....: D = (x0^2 - x_bar.inner_product(x_bar))*D + ....: D = D + 2*x_bar.tensor_product(x_bar) + ....: Q = matrix.block(2,2,[A,B,C,D]) sage: Q == x.quadratic_representation().matrix() True Test all of the properties from Theorem 11.2 in Alizadeh:: - sage: set_random_seed() sage: J = random_eja() sage: x,y = J.random_elements(2) sage: Lx = x.operator() @@ -1212,11 +1331,11 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): sage: J = RealSymmetricEJA(3) sage: J.one() - e0 + e2 + e5 + b0 + b2 + b5 sage: J.one().spectral_decomposition() - [(1, e0 + e2 + e5)] + [(1, b0 + b2 + b5)] sage: J.zero().spectral_decomposition() - [(0, e0 + e2 + e5)] + [(0, b0 + b2 + b5)] TESTS:: @@ -1241,22 +1360,22 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The spectral decomposition should work in subalgebras, too:: sage: J = RealSymmetricEJA(4) - sage: (e0, e1, e2, e3, e4, e5, e6, e7, e8, e9) = J.gens() - sage: A = 2*e5 - 2*e8 + sage: (b0, b1, b2, b3, b4, b5, b6, b7, b8, b9) = J.gens() + sage: A = 2*b5 - 2*b8 sage: (lambda1, c1) = A.spectral_decomposition()[1] sage: (J0, J5, J1) = J.peirce_decomposition(c1) sage: (f0, f1, f2) = J1.gens() sage: f0.spectral_decomposition() - [(0, 1.000000000000000?*f2), (1, 1.000000000000000?*f0)] + [(0, 1.000000000000000?*c2), (1, 1.000000000000000?*c0)] """ - A = self.subalgebra_generated_by(orthonormalize_basis=True) + A = self.subalgebra_generated_by(orthonormalize=True) result = [] for (evalue, proj) in A(self).operator().spectral_decomposition(): result.append( (evalue, proj(A.one()).superalgebra_element()) ) return result - def subalgebra_generated_by(self, orthonormalize_basis=False): + def subalgebra_generated_by(self, **kwargs): """ Return the associative subalgebra of the parent EJA generated by this element. @@ -1270,15 +1389,27 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): SETUP:: - sage: from mjo.eja.eja_algebra import random_eja + sage: from mjo.eja.eja_algebra import (random_eja, + ....: HadamardEJA, + ....: RealSymmetricEJA) + + EXAMPLES: + + We can create subalgebras of Cartesian product EJAs that are not + themselves Cartesian product EJAs (they're just "regular" EJAs):: + + sage: J1 = HadamardEJA(3) + sage: J2 = RealSymmetricEJA(2) + sage: J = cartesian_product([J1,J2]) + sage: J.one().subalgebra_generated_by() + Euclidean Jordan algebra of dimension 1 over Algebraic Real Field TESTS: This subalgebra, being composed of only powers, is associative:: - sage: set_random_seed() sage: x0 = random_eja().random_element() - sage: A = x0.subalgebra_generated_by() + sage: A = x0.subalgebra_generated_by(orthonormalize=False) sage: x,y,z = A.random_elements(3) sage: (x*y)*z == x*(y*z) True @@ -1286,9 +1417,8 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): Squaring in the subalgebra should work the same as in the superalgebra:: - sage: set_random_seed() sage: x = random_eja().random_element() - sage: A = x.subalgebra_generated_by() + sage: A = x.subalgebra_generated_by(orthonormalize=False) sage: A(x^2) == A(x)*A(x) True @@ -1297,13 +1427,19 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): element... unless the original algebra was trivial, in which case the subalgebra is trivial too:: - sage: set_random_seed() sage: A = random_eja().zero().subalgebra_generated_by() sage: (A.is_trivial() and A.dimension() == 0) or A.dimension() == 1 True """ - return FiniteDimensionalEuclideanJordanElementSubalgebra(self, orthonormalize_basis) + powers = tuple( self**k for k in range(self.degree()) ) + A = self.parent().subalgebra(powers, + associative=True, + check_field=False, + check_axioms=False, + **kwargs) + A.one.set_cache(A(self.parent().one())) + return A def subalgebra_idempotent(self): @@ -1318,12 +1454,12 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): TESTS: Ensure that we can find an idempotent in a non-trivial algebra - where there are non-nilpotent elements:: + where there are non-nilpotent elements, or that we get the dumb + solution in the trivial algebra:: - sage: set_random_seed() - sage: J = random_eja(nontrivial=True) + sage: J = random_eja(field=QQ, orthonormalize=False) sage: x = J.random_element() - sage: while x.is_nilpotent(): + sage: while x.is_nilpotent() and not J.is_trivial(): ....: x = J.random_element() sage: c = x.subalgebra_idempotent() sage: c^2 == c @@ -1336,7 +1472,10 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): if self.is_nilpotent(): raise ValueError("this only works with non-nilpotent elements!") - J = self.subalgebra_generated_by() + # The subalgebra is transient (we return an element of the + # superalgebra, i.e. this algebra) so why bother + # orthonormalizing? + J = self.subalgebra_generated_by(orthonormalize=False) u = J(self) # The image of the matrix of left-u^m-multiplication @@ -1357,14 +1496,12 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): # subspace... or do we? Can't we just solve, knowing that # A(c) = u^(s+1) should have a solution in the big space, # too? - # - # Beware, solve_right() means that we're using COLUMN vectors. - # Our FiniteDimensionalAlgebraElement superclass uses rows. u_next = u**(s+1) A = u_next.operator().matrix() c = J.from_vector(A.solve_right(u_next.to_vector())) - # Now c is the idempotent we want, but it still lives in the subalgebra. + # Now c is the idempotent we want, but it still lives in + # the subalgebra. return c.superalgebra_element() @@ -1404,11 +1541,24 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The trace of an element is a real number:: - sage: set_random_seed() sage: J = random_eja() sage: J.random_element().trace() in RLF True + The trace is linear:: + + sage: J = random_eja() + sage: x,y = J.random_elements(2) + sage: alpha = J.base_ring().random_element() + sage: (alpha*x + y).trace() == alpha*x.trace() + y.trace() + True + + The trace of a square is nonnegative:: + + sage: x = random_eja().random_element() + sage: (x*x).trace() >= 0 + True + """ P = self.parent() r = P.rank() @@ -1438,14 +1588,13 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): The trace inner product is commutative, bilinear, and associative:: - sage: set_random_seed() sage: J = random_eja() sage: x,y,z = J.random_elements(3) sage: # commutative sage: x.trace_inner_product(y) == y.trace_inner_product(x) True sage: # bilinear - sage: a = J.base_ring().random_element(); + sage: a = J.base_ring().random_element() sage: actual = (a*(x+z)).trace_inner_product(y) sage: expected = ( a*x.trace_inner_product(y) + ....: a*z.trace_inner_product(y) ) @@ -1492,3 +1641,187 @@ class FiniteDimensionalEuclideanJordanAlgebraElement(IndexedFreeModuleElement): """ return self.trace_inner_product(self).sqrt() + + + def operator_trace_inner_product(self, other): + r""" + Return the operator inner product of myself and ``other``. + + The "operator inner product," whose name is not standard, is + defined be the usual linear-algebraic trace of the + ``(x*y).operator()``. + + Proposition III.1.5 in Faraut and Korányi shows that on any + Euclidean Jordan algebra, this is another associative inner + product under which the cone of squares is symmetric. + + This works even if the basis hasn't been orthonormalized + because the eigenvalues of the corresponding matrix don't + change when the basis does (they're preserved by any + similarity transformation). + + SETUP:: + + sage: from mjo.eja.eja_algebra import (JordanSpinEJA, + ....: RealSymmetricEJA, + ....: ComplexHermitianEJA, + ....: random_eja) + + EXAMPLES: + + Proposition III.4.2 of Faraut and Korányi shows that on a + simple algebra of rank `r` and dimension `n`, this inner + product is `n/r` times the canonical + :meth:`trace_inner_product`:: + + sage: J = JordanSpinEJA(4, field=QQ) + sage: x,y = J.random_elements(2) + sage: n = J.dimension() + sage: r = J.rank() + sage: actual = x.operator_trace_inner_product(y) + sage: expected = (n/r)*x.trace_inner_product(y) + sage: actual == expected + True + + :: + + sage: J = RealSymmetricEJA(3) + sage: x,y = J.random_elements(2) + sage: n = J.dimension() + sage: r = J.rank() + sage: actual = x.operator_trace_inner_product(y) + sage: expected = (n/r)*x.trace_inner_product(y) + sage: actual == expected + True + + :: + + sage: J = ComplexHermitianEJA(3, field=QQ, orthonormalize=False) + sage: x,y = J.random_elements(2) + sage: n = J.dimension() + sage: r = J.rank() + sage: actual = x.operator_trace_inner_product(y) + sage: expected = (n/r)*x.trace_inner_product(y) + sage: actual == expected + True + + TESTS: + + The operator inner product is commutative, bilinear, and + associative:: + + sage: J = random_eja() + sage: x,y,z = J.random_elements(3) + sage: # commutative + sage: actual = x.operator_trace_inner_product(y) + sage: expected = y.operator_trace_inner_product(x) + sage: actual == expected + True + sage: # bilinear + sage: a = J.base_ring().random_element() + sage: actual = (a*(x+z)).operator_trace_inner_product(y) + sage: expected = ( a*x.operator_trace_inner_product(y) + + ....: a*z.operator_trace_inner_product(y) ) + sage: actual == expected + True + sage: actual = x.operator_trace_inner_product(a*(y+z)) + sage: expected = ( a*x.operator_trace_inner_product(y) + + ....: a*x.operator_trace_inner_product(z) ) + sage: actual == expected + True + sage: # associative + sage: actual = (x*y).operator_trace_inner_product(z) + sage: expected = y.operator_trace_inner_product(x*z) + sage: actual == expected + True + + Despite the fact that the implementation uses a matrix representation, + the answer is independent of the basis used:: + + sage: J = RealSymmetricEJA(3, field=QQ, orthonormalize=False) + sage: V = RealSymmetricEJA(3) + sage: x,y = J.random_elements(2) + sage: w = V(x.to_matrix()) + sage: z = V(y.to_matrix()) + sage: expected = x.operator_trace_inner_product(y) + sage: actual = w.operator_trace_inner_product(z) + sage: actual == expected + True + + """ + if not other in self.parent(): + raise TypeError("'other' must live in the same algebra") + + return (self*other).operator().matrix().trace() + + + def operator_trace_norm(self): + """ + The norm of this element with respect to + :meth:`operator_trace_inner_product`. + + SETUP:: + + sage: from mjo.eja.eja_algebra import (JordanSpinEJA, + ....: HadamardEJA) + + EXAMPLES: + + On a simple algebra, this will differ from :meth:`trace_norm` + by the scalar factor ``(n/r).sqrt()``, where `n` is the + dimension of the algebra and `r` its rank. This follows from + the corresponding result (Proposition III.4.2 of Faraut and + Korányi) for the trace inner product:: + + sage: J = HadamardEJA(2) + sage: x = sum(J.gens()) + sage: x.operator_trace_norm() + 1.414213562373095? + + :: + + sage: J = JordanSpinEJA(4) + sage: x = sum(J.gens()) + sage: x.operator_trace_norm() + 4 + + """ + return self.operator_trace_inner_product(self).sqrt() + + +class CartesianProductParentEJAElement(EJAElement): + r""" + An intermediate class for elements that have a Cartesian + product as their parent algebra. + + This is needed because the ``to_matrix`` method (which gives you a + representation from the superalgebra) needs to do special stuff + for Cartesian products. Specifically, an EJA subalgebra of a + Cartesian product EJA will not itself be a Cartesian product (it + has its own basis) -- but we want ``to_matrix()`` to be able to + give us a Cartesian product representation. + """ + def to_matrix(self): + # An override is necessary to call our custom _scale(). + B = self.parent().matrix_basis() + W = self.parent().matrix_space() + + # Aaaaand linear combinations don't work in Cartesian + # product spaces, even though they provide a method with + # that name. This is hidden in a subclass because the + # _scale() function is slow. + pairs = zip(B, self.to_vector()) + return W.sum( _scale(b, alpha) for (b,alpha) in pairs ) + +class CartesianProductEJAElement(CartesianProductParentEJAElement): + def det(self): + r""" + Compute the determinant of this product-element using the + determianants of its factors. + + This result Follows from the spectral decomposition of (say) + the pair `(x,y)` in terms of the Jordan frame `\left\{ (c_1, + 0),(c_2, 0),...,(0,d_1),(0,d_2),... \right\}. + """ + from sage.misc.misc_c import prod + return prod( f.det() for f in self.cartesian_factors() )