00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033 #ifndef EIGEN_SAEIGENSOLVER_MKL_H
00034 #define EIGEN_SAEIGENSOLVER_MKL_H
00035
00036 #include "Eigen/src/Core/util/MKL_support.h"
00037
00038 namespace Eigen {
00039
00042 #define EIGEN_MKL_EIG_SELFADJ(EIGTYPE, MKLTYPE, MKLRTYPE, MKLNAME, EIGCOLROW, MKLCOLROW ) \
00043 template<> inline \
00044 SelfAdjointEigenSolver<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW> >& \
00045 SelfAdjointEigenSolver<Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW> >::compute(const Matrix<EIGTYPE, Dynamic, Dynamic, EIGCOLROW>& matrix, int options) \
00046 { \
00047 eigen_assert(matrix.cols() == matrix.rows()); \
00048 eigen_assert((options&~(EigVecMask|GenEigMask))==0 \
00049 && (options&EigVecMask)!=EigVecMask \
00050 && "invalid option parameter"); \
00051 bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors; \
00052 lapack_int n = matrix.cols(), lda, matrix_order, info; \
00053 m_eivalues.resize(n,1); \
00054 m_subdiag.resize(n-1); \
00055 m_eivec = matrix; \
00056 \
00057 if(n==1) \
00058 { \
00059 m_eivalues.coeffRef(0,0) = internal::real(matrix.coeff(0,0)); \
00060 if(computeEigenvectors) m_eivec.setOnes(n,n); \
00061 m_info = Success; \
00062 m_isInitialized = true; \
00063 m_eigenvectorsOk = computeEigenvectors; \
00064 return *this; \
00065 } \
00066 \
00067 lda = matrix.outerStride(); \
00068 matrix_order=MKLCOLROW; \
00069 char jobz, uplo='L'; \
00070 jobz = computeEigenvectors ? 'V' : 'N'; \
00071 \
00072 info = LAPACKE_##MKLNAME( matrix_order, jobz, uplo, n, (MKLTYPE*)m_eivec.data(), lda, (MKLRTYPE*)m_eivalues.data() ); \
00073 m_info = (info==0) ? Success : NoConvergence; \
00074 m_isInitialized = true; \
00075 m_eigenvectorsOk = computeEigenvectors; \
00076 return *this; \
00077 }
00078
00079
00080 EIGEN_MKL_EIG_SELFADJ(double, double, double, dsyev, ColMajor, LAPACK_COL_MAJOR)
00081 EIGEN_MKL_EIG_SELFADJ(float, float, float, ssyev, ColMajor, LAPACK_COL_MAJOR)
00082 EIGEN_MKL_EIG_SELFADJ(dcomplex, MKL_Complex16, double, zheev, ColMajor, LAPACK_COL_MAJOR)
00083 EIGEN_MKL_EIG_SELFADJ(scomplex, MKL_Complex8, float, cheev, ColMajor, LAPACK_COL_MAJOR)
00084
00085 EIGEN_MKL_EIG_SELFADJ(double, double, double, dsyev, RowMajor, LAPACK_ROW_MAJOR)
00086 EIGEN_MKL_EIG_SELFADJ(float, float, float, ssyev, RowMajor, LAPACK_ROW_MAJOR)
00087 EIGEN_MKL_EIG_SELFADJ(dcomplex, MKL_Complex16, double, zheev, RowMajor, LAPACK_ROW_MAJOR)
00088 EIGEN_MKL_EIG_SELFADJ(scomplex, MKL_Complex8, float, cheev, RowMajor, LAPACK_ROW_MAJOR)
00089
00090 }
00091
00092 #endif // EIGEN_SAEIGENSOLVER_H