18 #ifndef itkSymmetricEigenAnalysis_h
19 #define itkSymmetricEigenAnalysis_h
22 #include "itk_eigen.h"
23 #include ITK_EIGEN(Eigenvalues)
27 #include "vnl/vnl_matrix.h"
28 #include "vnl/vnl_matrix_fixed.h"
36 template <
typename TValueType,
unsigned int NRows,
unsigned int NCols>
40 return inputMatrix.data_block();
42 template <
typename TValueType>
46 return inputMatrix.data_block();
49 template <
typename TValueType,
unsigned int NRows,
unsigned int NCols>
71 template <
typename TArray>
75 std::vector<int> indicesSortPermutations(numberOfElements, 0);
76 std::iota(std::begin(indicesSortPermutations), std::end(indicesSortPermutations), 0);
80 std::begin(indicesSortPermutations),
81 std::end(indicesSortPermutations),
82 [&eigenValues](
unsigned int a,
unsigned int b) {
return std::abs(eigenValues[a]) < std::abs(eigenValues[b]); });
83 auto tmpCopy = eigenValues;
84 for (
unsigned int i = 0; i < numberOfElements; ++i)
86 eigenValues[i] = tmpCopy[indicesSortPermutations[i]];
88 return indicesSortPermutations;
99 template <
typename QMatrix>
103 using EigenLibPermutationMatrix = Eigen::PermutationMatrix<Eigen::Dynamic, Eigen::Dynamic>;
104 auto numberOfElements = indicesSortPermutations.size();
107 EigenLibPermutationMatrix perm(numberOfElements);
109 std::copy(indicesSortPermutations.begin(), indicesSortPermutations.end(), perm.indices().data());
111 eigenVectors = eigenVectors * perm;
138 extern ITKCommon_EXPORT std::ostream &
157 itkGenericExceptionMacro(<<
"Error: Invalid value for conversion.");
160 #if !defined(ITK_LEGACY_REMOVE)
202 template <
typename TMatrix,
typename TVector,
typename TEigenMatrix = TMatrix>
207 #if !defined(ITK_LEGACY_REMOVE)
215 : m_Dimension(dimension)
239 ComputeEigenValues(
const TMatrix & A, TVector & D)
const;
262 ComputeEigenValuesAndVectors(
const TMatrix & A, TVector & EigenValues, TEigenMatrix & EigenVectors)
const;
335 m_Order = m_Dimension;
352 m_UseEigenLibrary = input;
357 m_UseEigenLibrary =
true;
362 m_UseEigenLibrary =
false;
367 return m_UseEigenLibrary;
372 bool m_UseEigenLibrary{
false };
373 unsigned int m_Dimension{ 0 };
374 unsigned int m_Order{ 0 };
397 ReduceToTridiagonalMatrix(
double * a,
double * d,
double *
e,
double * e2)
const;
421 ReduceToTridiagonalMatrixAndGetTransformation(
const double * a,
double * d,
double *
e,
double * z)
const;
453 ComputeEigenValuesUsingQL(
double * d,
double *
e)
const;
493 ComputeEigenValuesAndVectorsUsingQL(
double * d,
double *
e,
double * z)
const;
499 ComputeEigenValuesLegacy(
const TMatrix & A, TVector & D)
const;
502 ComputeEigenValuesAndVectorsLegacy(
const TMatrix & A, TVector & EigenValues, TEigenMatrix & EigenVectors)
const;
518 template <
typename QMatrix = TMatrix>
522 return QMatrix::element_type();
524 template <
typename QMatrix = TMatrix>
528 return QMatrix::ValueType();
534 TVector & EigenValues,
535 TEigenMatrix & EigenVectors)
const
537 return ComputeEigenValuesAndVectorsWithEigenLibraryImpl(A, EigenValues, EigenVectors,
true);
546 template <
typename QMatrix>
549 TVector & EigenValues,
550 TEigenMatrix & EigenVectors,
551 long)
const -> decltype(static_cast<unsigned int>(1))
553 using ValueType = decltype(GetMatrixValueType(
true));
554 using EigenLibMatrixType = Eigen::Matrix<ValueType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
555 EigenLibMatrixType inputMatrix(m_Dimension, m_Dimension);
556 for (
unsigned int row = 0; row < m_Dimension; ++row)
558 for (
unsigned int col = 0; col < m_Dimension; ++col)
560 inputMatrix(row, col) = A(row, col);
563 using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
564 EigenSolverType solver(inputMatrix);
565 const auto & eigenValues = solver.eigenvalues();
569 const auto & eigenVectors = solver.eigenvectors();
573 auto copyEigenValues = eigenValues;
574 auto copyEigenVectors = eigenVectors;
578 for (
unsigned int row = 0; row < m_Dimension; ++row)
580 EigenValues[row] = copyEigenValues[row];
581 for (
unsigned int col = 0; col < m_Dimension; ++col)
583 EigenVectors[row][col] = copyEigenVectors(col, row);
589 for (
unsigned int row = 0; row < m_Dimension; ++row)
591 EigenValues[row] = eigenValues[row];
592 for (
unsigned int col = 0; col < m_Dimension; ++col)
594 EigenVectors[row][col] = eigenVectors(col, row);
612 template <
typename QMatrix>
615 TVector & EigenValues,
616 TEigenMatrix & EigenVectors,
621 using PointerType = decltype(pointerToData);
622 using ValueTypeCV =
typename std::remove_pointer<PointerType>::type;
623 using ValueType =
typename std::remove_cv<ValueTypeCV>::type;
624 using EigenLibMatrixType = Eigen::Matrix<ValueType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
625 using EigenConstMatrixMap = Eigen::Map<const EigenLibMatrixType>;
626 EigenConstMatrixMap inputMatrix(pointerToData, m_Dimension, m_Dimension);
627 using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
628 EigenSolverType solver(inputMatrix);
629 const auto & eigenValues = solver.eigenvalues();
633 const auto & eigenVectors = solver.eigenvectors();
636 auto copyEigenValues = eigenValues;
637 auto copyEigenVectors = eigenVectors;
640 for (
unsigned int row = 0; row < m_Dimension; ++row)
642 EigenValues[row] = copyEigenValues[row];
643 for (
unsigned int col = 0; col < m_Dimension; ++col)
645 EigenVectors[row][col] = copyEigenVectors(col, row);
651 for (
unsigned int row = 0; row < m_Dimension; ++row)
653 EigenValues[row] = eigenValues[row];
654 for (
unsigned int col = 0; col < m_Dimension; ++col)
656 EigenVectors[row][col] = eigenVectors(col, row);
668 return ComputeEigenValuesWithEigenLibraryImpl(A, EigenValues,
true);
677 template <
typename QMatrix>
680 -> decltype(static_cast<unsigned int>(1))
682 using ValueType = decltype(GetMatrixValueType(
true));
683 using EigenLibMatrixType = Eigen::Matrix<ValueType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
684 EigenLibMatrixType inputMatrix(m_Dimension, m_Dimension);
685 for (
unsigned int row = 0; row < m_Dimension; ++row)
687 for (
unsigned int col = 0; col < m_Dimension; ++col)
689 inputMatrix(row, col) = A(row, col);
692 using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
693 EigenSolverType solver(inputMatrix, Eigen::EigenvaluesOnly);
694 auto eigenValues = solver.eigenvalues();
699 for (
unsigned int i = 0; i < m_Dimension; ++i)
701 EigenValues[i] = eigenValues[i];
717 template <
typename QMatrix>
723 using PointerType = decltype(pointerToData);
724 using ValueTypeCV =
typename std::remove_pointer<PointerType>::type;
725 using ValueType =
typename std::remove_cv<ValueTypeCV>::type;
726 using EigenLibMatrixType = Eigen::Matrix<ValueType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
727 using EigenConstMatrixMap = Eigen::Map<const EigenLibMatrixType>;
728 EigenConstMatrixMap inputMatrix(pointerToData, m_Dimension, m_Dimension);
729 using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
730 EigenSolverType solver(inputMatrix, Eigen::EigenvaluesOnly);
731 auto eigenValues = solver.eigenvalues();
736 for (
unsigned int i = 0; i < m_Dimension; ++i)
738 EigenValues[i] = eigenValues[i];
745 template <
typename TMatrix,
typename TVector,
typename TEigenMatrix>
749 os <<
"[ClassType: SymmetricEigenAnalysis]" << std::endl;
751 os <<
" Order : " << s.
GetOrder() << std::endl;
758 template <
unsigned int VDimension,
typename TMatrix,
typename TVector,
typename TEigenMatrix = TMatrix>
762 #if !defined(ITK_LEGACY_REMOVE)
766 #if !defined(ITK_LEGACY_REMOVE)
797 return ComputeEigenValuesWithEigenLibraryImpl(A, EigenValues,
true);
823 return ComputeEigenValuesAndVectorsWithEigenLibraryImpl(A, EigenValues, EigenVectors,
true);
860 constexpr
unsigned int
865 constexpr
unsigned int
887 template <
typename QMatrix = TMatrix>
891 return QMatrix::element_type();
893 template <
typename QMatrix = TMatrix>
897 return QMatrix::ValueType();
909 template <
typename QMatrix>
912 TVector & EigenValues,
913 TEigenMatrix & EigenVectors,
918 using PointerType = decltype(pointerToData);
919 using ValueTypeCV =
typename std::remove_pointer<PointerType>::type;
920 using ValueType =
typename std::remove_cv<ValueTypeCV>::type;
921 using EigenLibMatrixType = Eigen::Matrix<ValueType, VDimension, VDimension, Eigen::RowMajor>;
922 using EigenConstMatrixMap = Eigen::Map<const EigenLibMatrixType>;
923 EigenConstMatrixMap inputMatrix(pointerToData);
924 using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
925 EigenSolverType solver(inputMatrix);
926 const auto & eigenValues = solver.eigenvalues();
930 const auto & eigenVectors = solver.eigenvectors();
933 auto copyEigenValues = eigenValues;
934 auto copyEigenVectors = eigenVectors;
937 for (
unsigned int row = 0; row < VDimension; ++row)
939 EigenValues[row] = copyEigenValues[row];
940 for (
unsigned int col = 0; col < VDimension; ++col)
942 EigenVectors[row][col] = copyEigenVectors(col, row);
948 for (
unsigned int row = 0; row < VDimension; ++row)
950 EigenValues[row] = eigenValues[row];
951 for (
unsigned int col = 0; col < VDimension; ++col)
953 EigenVectors[row][col] = eigenVectors(col, row);
967 template <
typename QMatrix>
970 TVector & EigenValues,
971 TEigenMatrix & EigenVectors,
972 long)
const -> decltype(static_cast<unsigned int>(1))
974 using ValueType = decltype(GetMatrixValueType(
true));
975 using EigenLibMatrixType = Eigen::Matrix<ValueType, VDimension, VDimension, Eigen::RowMajor>;
976 EigenLibMatrixType inputMatrix;
977 for (
unsigned int row = 0; row < VDimension; ++row)
979 for (
unsigned int col = 0; col < VDimension; ++col)
981 inputMatrix(row, col) = A(row, col);
984 using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
985 EigenSolverType solver(inputMatrix);
986 const auto & eigenValues = solver.eigenvalues();
990 const auto & eigenVectors = solver.eigenvectors();
994 auto copyEigenValues = eigenValues;
995 auto copyEigenVectors = eigenVectors;
999 for (
unsigned int row = 0; row < VDimension; ++row)
1001 EigenValues[row] = copyEigenValues[row];
1002 for (
unsigned int col = 0; col < VDimension; ++col)
1004 EigenVectors[row][col] = copyEigenVectors(col, row);
1010 for (
unsigned int row = 0; row < VDimension; ++row)
1012 EigenValues[row] = eigenValues[row];
1013 for (
unsigned int col = 0; col < VDimension; ++col)
1015 EigenVectors[row][col] = eigenVectors(col, row);
1029 template <
typename QMatrix>
1032 -> decltype(static_cast<unsigned int>(1))
1034 using ValueType = decltype(GetMatrixValueType(
true));
1035 using EigenLibMatrixType = Eigen::Matrix<ValueType, VDimension, VDimension, Eigen::RowMajor>;
1036 EigenLibMatrixType inputMatrix;
1037 for (
unsigned int row = 0; row < VDimension; ++row)
1039 for (
unsigned int col = 0; col < VDimension; ++col)
1041 inputMatrix(row, col) = A(row, col);
1044 using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
1045 EigenSolverType solver(inputMatrix, Eigen::EigenvaluesOnly);
1046 auto eigenValues = solver.eigenvalues();
1051 for (
unsigned int i = 0; i < VDimension; ++i)
1053 EigenValues[i] = eigenValues[i];
1069 template <
typename QMatrix>
1075 using PointerType = decltype(pointerToData);
1076 using ValueTypeCV =
typename std::remove_pointer<PointerType>::type;
1077 using ValueType =
typename std::remove_cv<ValueTypeCV>::type;
1078 using EigenLibMatrixType = Eigen::Matrix<ValueType, VDimension, VDimension, Eigen::RowMajor>;
1079 using EigenConstMatrixMap = Eigen::Map<const EigenLibMatrixType>;
1080 EigenConstMatrixMap inputMatrix(pointerToData);
1081 using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
1082 EigenSolverType solver(inputMatrix, Eigen::EigenvaluesOnly);
1083 auto eigenValues = solver.eigenvalues();
1088 for (
unsigned int i = 0; i < VDimension; ++i)
1090 EigenValues[i] = eigenValues[i];
1097 template <
unsigned int VDimension,
typename TMatrix,
typename TVector,
typename TEigenMatrix>
1102 os <<
"[ClassType: SymmetricEigenAnalysisFixedDimension]" << std::endl;
1103 os <<
" Dimension : " << s.
GetDimension() << std::endl;
1104 os <<
" Order : " << s.
GetOrder() << std::endl;
1112 #ifndef ITK_MANUAL_INSTANTIATION
1113 # include "itkSymmetricEigenAnalysis.hxx"