MatrixDiag.hh
Go to the documentation of this file.
00001 #ifndef RIVET_MATH_MATRIXDIAG 00002 #define RIVET_MATH_MATRIXDIAG 00003 00004 #include "Rivet/Math/MathHeader.hh" 00005 #include "Rivet/Math/MatrixN.hh" 00006 00007 #include "gsl/gsl_vector.h" 00008 #include "gsl/gsl_matrix.h" 00009 #include "gsl/gsl_eigen.h" 00010 00011 namespace Rivet { 00012 00013 00014 // // GSL forward declarations (avoids need for GSL header files) 00015 // extern "C" { 00016 // struct gsl_vector; 00017 // gsl_vector* gsl_vector_alloc(size_t); 00018 // double gsl_vector_get(gsl_vector*, size_t); 00019 // void gsl_vector_set(gsl_vector*, size_t, double); 00020 // void gsl_vector_free(gsl_vector*); 00021 // struct gsl_matrix; 00022 // gsl_matrix* gsl_matrix_alloc(size_t, size_t); 00023 // double gsl_matrix_get(gsl_matrix*, size_t, size_t); 00024 // void gsl_matrix_set(gsl_matrix*, size_t, size_t, double); 00025 // void gsl_matrix_free(gsl_matrix*); 00026 // struct gsl_eigen_symmv_workspace; 00027 // gsl_eigen_symmv_workspace* gsl_eigen_symmv_alloc(size_t); 00028 // void gsl_eigen_symmv(gsl_matrix*, gsl_vector*, gsl_matrix*, gsl_eigen_symmv_workspace*); 00029 // void gsl_eigen_symmv_free(gsl_eigen_symmv_workspace*); 00030 // typedef enum { 00031 // GSL_EIGEN_SORT_VAL_ASC, 00032 // GSL_EIGEN_SORT_VAL_DESC, 00033 // GSL_EIGEN_SORT_ABS_ASC, 00034 // GSL_EIGEN_SORT_ABS_DESC 00035 // } 00036 // gsl_eigen_sort_t; 00037 // int gsl_eigen_symmv_sort(gsl_vector * eval, gsl_matrix * evec, gsl_eigen_sort_t sort_type); 00038 // } 00039 00040 00041 template <size_t N> 00042 class EigenSystem; 00043 template <size_t N> 00044 EigenSystem<N> diagonalize(const Matrix<N>& m); 00045 00046 00047 /// Handy object containing results of a diagonalization. 00048 template <size_t N> 00049 class EigenSystem { 00050 template <size_t M> 00051 friend EigenSystem<M> diagonalize(const Matrix<M>&); 00052 00053 public: 00054 typedef pair<double, Vector<N> > EigenPair; 00055 typedef vector<EigenPair> EigenPairs; 00056 00057 Vector<N> getDiagVector() const { 00058 assert(_eigenPairs.size() == N); 00059 Vector<N> ret; 00060 for (size_t i = 0; i < N; ++i) { 00061 ret.set(i, _eigenPairs[i].first); 00062 } 00063 return ret; 00064 } 00065 00066 Matrix<N> getDiagMatrix() const { 00067 return Matrix<N>::mkDiag(getDiagVector()); 00068 } 00069 00070 EigenPairs getEigenPairs() const { 00071 return _eigenPairs; 00072 } 00073 00074 vector<double> getEigenValues() const { 00075 assert(_eigenPairs.size() == N); 00076 vector<double> ret; 00077 for (size_t i = 0; i < N; ++i) { 00078 ret.push_back(_eigenPairs[i].first); 00079 } 00080 return ret; 00081 } 00082 00083 vector<Vector<N> > getEigenVectors() const { 00084 assert(_eigenPairs.size() == N); 00085 vector<Vector<N> > ret; 00086 for (size_t i = 0; i < N; ++i) { 00087 ret.push_back(_eigenPairs[i].second); 00088 } 00089 return ret; 00090 } 00091 00092 private: 00093 EigenPairs _eigenPairs; 00094 }; 00095 00096 00097 /// Comparison functor for "eigen-pairs". 00098 template <size_t N> 00099 struct EigenPairCmp : 00100 public std::binary_function<const typename EigenSystem<N>::EigenPair&, 00101 const typename EigenSystem<N>::EigenPair&, bool> { 00102 bool operator()(const typename EigenSystem<N>::EigenPair& a, 00103 const typename EigenSystem<N>::EigenPair& b) { 00104 return a.first < b.first; 00105 } 00106 }; 00107 00108 00109 /// Diagonalize an NxN matrix, returning a collection of pairs of 00110 /// eigenvalues and eigenvectors, ordered decreasing in eigenvalue. 00111 template <size_t N> 00112 EigenSystem<N> diagonalize(const Matrix<N>& m) { 00113 EigenSystem<N> esys; 00114 00115 // Make a GSL matrix. 00116 gsl_matrix* A = gsl_matrix_alloc(N, N); 00117 for (size_t i = 0; i < N; ++i) { 00118 for (size_t j = 0; j < N; ++j) { 00119 gsl_matrix_set(A, i, j, m.get(i, j)); 00120 } 00121 } 00122 00123 // Use GSL diagonalization. 00124 gsl_matrix* vecs = gsl_matrix_alloc(N, N); 00125 gsl_vector* vals = gsl_vector_alloc(N); 00126 gsl_eigen_symmv_workspace* workspace = gsl_eigen_symmv_alloc(N); 00127 gsl_eigen_symmv(A, vals, vecs, workspace); 00128 gsl_eigen_symmv_sort(vals, vecs, GSL_EIGEN_SORT_VAL_DESC); 00129 00130 // Build the vector of "eigen-pairs". 00131 typename EigenSystem<N>::EigenPairs eigensolns; 00132 for (size_t i = 0; i < N; ++i) { 00133 typename EigenSystem<N>::EigenPair ep; 00134 ep.first = gsl_vector_get(vals, i); 00135 Vector<N> ev; 00136 for (size_t j = 0; j < N; ++j) { 00137 ev.set(j, gsl_matrix_get(vecs, j, i)); 00138 } 00139 ep.second = ev; 00140 eigensolns.push_back(ep); 00141 } 00142 00143 // Free GSL memory. 00144 gsl_eigen_symmv_free(workspace); 00145 gsl_matrix_free(A); 00146 gsl_matrix_free(vecs); 00147 gsl_vector_free(vals); 00148 00149 // Populate the returned object. 00150 esys._eigenPairs = eigensolns; 00151 return esys; 00152 } 00153 00154 00155 template <size_t N> 00156 inline const string toString(const typename EigenSystem<N>::EigenPair& e) { 00157 ostringstream ss; 00158 //for (typename EigenSystem<N>::EigenPairs::const_iterator i = e.begin(); i != e.end(); ++i) { 00159 ss << e->first << " -> " << e->second; 00160 // if (i+1 != e.end()) ss << endl; 00161 //} 00162 return ss.str(); 00163 } 00164 00165 template <size_t N> 00166 inline ostream& operator<<(std::ostream& out, const typename EigenSystem<N>::EigenPair& e) { 00167 out << toString(e); 00168 return out; 00169 } 00170 00171 00172 } 00173 00174 #endif Generated on Fri Dec 21 2012 15:03:41 for The Rivet MC analysis system by ![]() |