X-Git-Url: http://gitweb.michael.orlitzky.com/?a=blobdiff_plain;f=mjo%2Feja%2Feuclidean_jordan_algebra.py;h=a0ba1c68bb30bf54383807c52961ad82e4a69f03;hb=cc75c4e093a920c44f1eaaef4a1baa525b5c5727;hp=421c70bfa42d1e0ae8f837086437c491a9ab98f9;hpb=0e576704d235f1818b4ecd9f6ac0bdd5823fd745;p=sage.d.git diff --git a/mjo/eja/euclidean_jordan_algebra.py b/mjo/eja/euclidean_jordan_algebra.py index 421c70b..a0ba1c6 100644 --- a/mjo/eja/euclidean_jordan_algebra.py +++ b/mjo/eja/euclidean_jordan_algebra.py @@ -69,7 +69,6 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra): True """ - self._charpoly = None # for caching self._rank = rank self._natural_basis = natural_basis self._multiplication_table = mult_table @@ -89,6 +88,56 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra): + @cached_method + def _charpoly_coeff(self, i): + """ + Return the coefficient polynomial "a_{i}" of this algebra's + general characteristic polynomial. + + Having this be a separate cached method lets us compute and + store the trace/determinant (a_{r-1} and a_{0} respectively) + separate from the entire characteristic polynomial. + """ + (A_of_x, x) = self._charpoly_matrix() + R = A_of_x.base_ring() + A_cols = A_of_x.columns() + A_cols[i] = (x**self.rank()).vector() + numerator = column_matrix(A_of_x.base_ring(), A_cols).det() + denominator = A_of_x.det() + + # We're relying on the theory here to ensure that each a_i is + # indeed back in R, and the added negative signs are to make + # the whole charpoly expression sum to zero. + return R(-numerator/denominator) + + + @cached_method + def _charpoly_matrix(self): + """ + Compute the matrix whose entries A_ij are polynomials in + X1,...,XN. This same matrix is used in more than one method and + it's not so fast to construct. + """ + r = self.rank() + n = self.dimension() + + # Construct a new algebra over a multivariate polynomial ring... + names = ['X' + str(i) for i in range(1,n+1)] + R = PolynomialRing(self.base_ring(), names) + J = FiniteDimensionalEuclideanJordanAlgebra(R, + self._multiplication_table, + rank=r) + + idmat = identity_matrix(J.base_ring(), n) + + x = J(vector(R, R.gens())) + l1 = [column_matrix((x**k).vector()) for k in range(r)] + l2 = [idmat.column(k-1).column() for k in range(r+1, n+1)] + A_of_x = block_matrix(R, 1, n, (l1 + l2)) + return (A_of_x, x) + + + @cached_method def characteristic_polynomial(self): """ EXAMPLES: @@ -104,72 +153,22 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra): t^2 - 2*t + 1 """ - if self._charpoly is not None: - return self._charpoly - r = self.rank() n = self.dimension() - # First, compute the basis B... - x0 = self.zero() - c = 1 - for g in self.gens(): - x0 += c*g - c +=1 - if not x0.is_regular(): - raise ValueError("don't know a regular element") - - V = x0.vector().parent().ambient_vector_space() - V1 = V.span_of_basis( (x0**k).vector() for k in range(self.rank()) ) - B = (V1.basis() + V1.complement().basis()) - - # Now switch to the polynomial rings. - - names = ['X' + str(i) for i in range(1,n+1)] - R = PolynomialRing(self.base_ring(), names) - J = FiniteDimensionalEuclideanJordanAlgebra(R, - self._multiplication_table, - rank=r) - B = [ b.change_ring(R.fraction_field()) for b in B ] - # Get the vector space (as opposed to module) so that - # span_of_basis() works. - V = J.zero().vector().parent().ambient_vector_space() - W = V.span_of_basis(B) - - def e(k): - # The coordinates of e_k with respect to the basis B. - # But, the e_k are elements of B... - return identity_matrix(J.base_ring(), n).column(k-1).column() - - # A matrix implementation 1 - x = J(vector(R, R.gens())) - l1 = [column_matrix(W.coordinates((x**k).vector())) for k in range(r)] - l2 = [e(k) for k in range(r+1, n+1)] - A_of_x = block_matrix(1, n, (l1 + l2)) - xr = W.coordinates((x**r).vector()) - a = [] - denominator = A_of_x.det() # This is constant - for i in range(n): - A_cols = A_of_x.columns() - A_cols[i] = xr - numerator = column_matrix(A_of_x.base_ring(), A_cols).det() - ai = numerator/denominator - a.append(ai) + # The list of coefficient polynomials a_1, a_2, ..., a_n. + a = [ self._charpoly_coeff(i) for i in range(n) ] # We go to a bit of trouble here to reorder the # indeterminates, so that it's easier to evaluate the # characteristic polynomial at x's coordinates and get back # something in terms of t, which is what we want. + R = a[0].parent() S = PolynomialRing(self.base_ring(),'t') t = S.gen(0) S = PolynomialRing(S, R.variable_names()) t = S(t) - # We're relying on the theory here to ensure that each entry - # a[i] is indeed back in R, and the added negative signs are - # to make the whole expression sum to zero. - a = [R(-ai) for ai in a] # corresponds to powerx x^0 through x^(r-1) - # Note: all entries past the rth should be zero. The # coefficient of the highest power (x^r) is 1, but it doesn't # appear in the solution vector which contains coefficients @@ -181,8 +180,7 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra): # assign a[r] goes out-of-bounds. a.append(1) # corresponds to x^r - self._charpoly = sum( a[k]*(t**k) for k in range(len(a)) ) - return self._charpoly + return sum( a[k]*(t**k) for k in range(len(a)) ) def inner_product(self, x, y): @@ -1192,17 +1190,34 @@ class FiniteDimensionalEuclideanJordanAlgebra(FiniteDimensionalAlgebra): EXAMPLES:: sage: J = JordanSpinEJA(3) - sage: e0,e1,e2 = J.gens() - sage: x = e0 + e1 + e2 + sage: x = sum(J.gens()) sage: x.trace() 2 + :: + + sage: J = RealCartesianProductEJA(5) + sage: J.one().trace() + 5 + + TESTS: + + The trace of an element is a real number:: + + sage: set_random_seed() + sage: J = random_eja() + sage: J.random_element().trace() in J.base_ring() + True + """ - cs = self.characteristic_polynomial().coefficients(sparse=False) - if len(cs) >= 2: - return -1*cs[-2] - else: - raise ValueError('charpoly had fewer than 2 coefficients') + P = self.parent() + r = P.rank() + p = P._charpoly_coeff(r-1) + # The _charpoly_coeff function already adds the factor of + # -1 to ensure that _charpoly_coeff(r-1) is really what + # appears in front of t^{r-1} in the charpoly. However, + # we want the negative of THAT for the trace. + return -p(*self.vector()) def trace_inner_product(self, other):