SelfAdjointEigenSolver.h 33 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. // Copyright (C) 2010 Jitse Niesen <jitse@maths.leeds.ac.uk>
  6. //
  7. // This Source Code Form is subject to the terms of the Mozilla
  8. // Public License v. 2.0. If a copy of the MPL was not distributed
  9. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  10. #ifndef EIGEN_SELFADJOINTEIGENSOLVER_H
  11. #define EIGEN_SELFADJOINTEIGENSOLVER_H
  12. #include "./Tridiagonalization.h"
  13. namespace Eigen {
  14. template<typename _MatrixType>
  15. class GeneralizedSelfAdjointEigenSolver;
  16. namespace internal {
  17. template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues;
  18. template<typename MatrixType, typename DiagType, typename SubDiagType>
  19. ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType& eivec);
  20. }
  21. /** \eigenvalues_module \ingroup Eigenvalues_Module
  22. *
  23. *
  24. * \class SelfAdjointEigenSolver
  25. *
  26. * \brief Computes eigenvalues and eigenvectors of selfadjoint matrices
  27. *
  28. * \tparam _MatrixType the type of the matrix of which we are computing the
  29. * eigendecomposition; this is expected to be an instantiation of the Matrix
  30. * class template.
  31. *
  32. * A matrix \f$ A \f$ is selfadjoint if it equals its adjoint. For real
  33. * matrices, this means that the matrix is symmetric: it equals its
  34. * transpose. This class computes the eigenvalues and eigenvectors of a
  35. * selfadjoint matrix. These are the scalars \f$ \lambda \f$ and vectors
  36. * \f$ v \f$ such that \f$ Av = \lambda v \f$. The eigenvalues of a
  37. * selfadjoint matrix are always real. If \f$ D \f$ is a diagonal matrix with
  38. * the eigenvalues on the diagonal, and \f$ V \f$ is a matrix with the
  39. * eigenvectors as its columns, then \f$ A = V D V^{-1} \f$ (for selfadjoint
  40. * matrices, the matrix \f$ V \f$ is always invertible). This is called the
  41. * eigendecomposition.
  42. *
  43. * The algorithm exploits the fact that the matrix is selfadjoint, making it
  44. * faster and more accurate than the general purpose eigenvalue algorithms
  45. * implemented in EigenSolver and ComplexEigenSolver.
  46. *
  47. * Only the \b lower \b triangular \b part of the input matrix is referenced.
  48. *
  49. * Call the function compute() to compute the eigenvalues and eigenvectors of
  50. * a given matrix. Alternatively, you can use the
  51. * SelfAdjointEigenSolver(const MatrixType&, int) constructor which computes
  52. * the eigenvalues and eigenvectors at construction time. Once the eigenvalue
  53. * and eigenvectors are computed, they can be retrieved with the eigenvalues()
  54. * and eigenvectors() functions.
  55. *
  56. * The documentation for SelfAdjointEigenSolver(const MatrixType&, int)
  57. * contains an example of the typical use of this class.
  58. *
  59. * To solve the \em generalized eigenvalue problem \f$ Av = \lambda Bv \f$ and
  60. * the likes, see the class GeneralizedSelfAdjointEigenSolver.
  61. *
  62. * \sa MatrixBase::eigenvalues(), class EigenSolver, class ComplexEigenSolver
  63. */
  64. template<typename _MatrixType> class SelfAdjointEigenSolver
  65. {
  66. public:
  67. typedef _MatrixType MatrixType;
  68. enum {
  69. Size = MatrixType::RowsAtCompileTime,
  70. ColsAtCompileTime = MatrixType::ColsAtCompileTime,
  71. Options = MatrixType::Options,
  72. MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
  73. };
  74. /** \brief Scalar type for matrices of type \p _MatrixType. */
  75. typedef typename MatrixType::Scalar Scalar;
  76. typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
  77. typedef Matrix<Scalar,Size,Size,ColMajor,MaxColsAtCompileTime,MaxColsAtCompileTime> EigenvectorsType;
  78. /** \brief Real scalar type for \p _MatrixType.
  79. *
  80. * This is just \c Scalar if #Scalar is real (e.g., \c float or
  81. * \c double), and the type of the real part of \c Scalar if #Scalar is
  82. * complex.
  83. */
  84. typedef typename NumTraits<Scalar>::Real RealScalar;
  85. friend struct internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>;
  86. /** \brief Type for vector of eigenvalues as returned by eigenvalues().
  87. *
  88. * This is a column vector with entries of type #RealScalar.
  89. * The length of the vector is the size of \p _MatrixType.
  90. */
  91. typedef typename internal::plain_col_type<MatrixType, RealScalar>::type RealVectorType;
  92. typedef Tridiagonalization<MatrixType> TridiagonalizationType;
  93. typedef typename TridiagonalizationType::SubDiagonalType SubDiagonalType;
  94. /** \brief Default constructor for fixed-size matrices.
  95. *
  96. * The default constructor is useful in cases in which the user intends to
  97. * perform decompositions via compute(). This constructor
  98. * can only be used if \p _MatrixType is a fixed-size matrix; use
  99. * SelfAdjointEigenSolver(Index) for dynamic-size matrices.
  100. *
  101. * Example: \include SelfAdjointEigenSolver_SelfAdjointEigenSolver.cpp
  102. * Output: \verbinclude SelfAdjointEigenSolver_SelfAdjointEigenSolver.out
  103. */
  104. EIGEN_DEVICE_FUNC
  105. SelfAdjointEigenSolver()
  106. : m_eivec(),
  107. m_eivalues(),
  108. m_subdiag(),
  109. m_isInitialized(false)
  110. { }
  111. /** \brief Constructor, pre-allocates memory for dynamic-size matrices.
  112. *
  113. * \param [in] size Positive integer, size of the matrix whose
  114. * eigenvalues and eigenvectors will be computed.
  115. *
  116. * This constructor is useful for dynamic-size matrices, when the user
  117. * intends to perform decompositions via compute(). The \p size
  118. * parameter is only used as a hint. It is not an error to give a wrong
  119. * \p size, but it may impair performance.
  120. *
  121. * \sa compute() for an example
  122. */
  123. EIGEN_DEVICE_FUNC
  124. explicit SelfAdjointEigenSolver(Index size)
  125. : m_eivec(size, size),
  126. m_eivalues(size),
  127. m_subdiag(size > 1 ? size - 1 : 1),
  128. m_isInitialized(false)
  129. {}
  130. /** \brief Constructor; computes eigendecomposition of given matrix.
  131. *
  132. * \param[in] matrix Selfadjoint matrix whose eigendecomposition is to
  133. * be computed. Only the lower triangular part of the matrix is referenced.
  134. * \param[in] options Can be #ComputeEigenvectors (default) or #EigenvaluesOnly.
  135. *
  136. * This constructor calls compute(const MatrixType&, int) to compute the
  137. * eigenvalues of the matrix \p matrix. The eigenvectors are computed if
  138. * \p options equals #ComputeEigenvectors.
  139. *
  140. * Example: \include SelfAdjointEigenSolver_SelfAdjointEigenSolver_MatrixType.cpp
  141. * Output: \verbinclude SelfAdjointEigenSolver_SelfAdjointEigenSolver_MatrixType.out
  142. *
  143. * \sa compute(const MatrixType&, int)
  144. */
  145. template<typename InputType>
  146. EIGEN_DEVICE_FUNC
  147. explicit SelfAdjointEigenSolver(const EigenBase<InputType>& matrix, int options = ComputeEigenvectors)
  148. : m_eivec(matrix.rows(), matrix.cols()),
  149. m_eivalues(matrix.cols()),
  150. m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
  151. m_isInitialized(false)
  152. {
  153. compute(matrix.derived(), options);
  154. }
  155. /** \brief Computes eigendecomposition of given matrix.
  156. *
  157. * \param[in] matrix Selfadjoint matrix whose eigendecomposition is to
  158. * be computed. Only the lower triangular part of the matrix is referenced.
  159. * \param[in] options Can be #ComputeEigenvectors (default) or #EigenvaluesOnly.
  160. * \returns Reference to \c *this
  161. *
  162. * This function computes the eigenvalues of \p matrix. The eigenvalues()
  163. * function can be used to retrieve them. If \p options equals #ComputeEigenvectors,
  164. * then the eigenvectors are also computed and can be retrieved by
  165. * calling eigenvectors().
  166. *
  167. * This implementation uses a symmetric QR algorithm. The matrix is first
  168. * reduced to tridiagonal form using the Tridiagonalization class. The
  169. * tridiagonal matrix is then brought to diagonal form with implicit
  170. * symmetric QR steps with Wilkinson shift. Details can be found in
  171. * Section 8.3 of Golub \& Van Loan, <i>%Matrix Computations</i>.
  172. *
  173. * The cost of the computation is about \f$ 9n^3 \f$ if the eigenvectors
  174. * are required and \f$ 4n^3/3 \f$ if they are not required.
  175. *
  176. * This method reuses the memory in the SelfAdjointEigenSolver object that
  177. * was allocated when the object was constructed, if the size of the
  178. * matrix does not change.
  179. *
  180. * Example: \include SelfAdjointEigenSolver_compute_MatrixType.cpp
  181. * Output: \verbinclude SelfAdjointEigenSolver_compute_MatrixType.out
  182. *
  183. * \sa SelfAdjointEigenSolver(const MatrixType&, int)
  184. */
  185. template<typename InputType>
  186. EIGEN_DEVICE_FUNC
  187. SelfAdjointEigenSolver& compute(const EigenBase<InputType>& matrix, int options = ComputeEigenvectors);
  188. /** \brief Computes eigendecomposition of given matrix using a closed-form algorithm
  189. *
  190. * This is a variant of compute(const MatrixType&, int options) which
  191. * directly solves the underlying polynomial equation.
  192. *
  193. * Currently only 2x2 and 3x3 matrices for which the sizes are known at compile time are supported (e.g., Matrix3d).
  194. *
  195. * This method is usually significantly faster than the QR iterative algorithm
  196. * but it might also be less accurate. It is also worth noting that
  197. * for 3x3 matrices it involves trigonometric operations which are
  198. * not necessarily available for all scalar types.
  199. *
  200. * For the 3x3 case, we observed the following worst case relative error regarding the eigenvalues:
  201. * - double: 1e-8
  202. * - float: 1e-3
  203. *
  204. * \sa compute(const MatrixType&, int options)
  205. */
  206. EIGEN_DEVICE_FUNC
  207. SelfAdjointEigenSolver& computeDirect(const MatrixType& matrix, int options = ComputeEigenvectors);
  208. /**
  209. *\brief Computes the eigen decomposition from a tridiagonal symmetric matrix
  210. *
  211. * \param[in] diag The vector containing the diagonal of the matrix.
  212. * \param[in] subdiag The subdiagonal of the matrix.
  213. * \param[in] options Can be #ComputeEigenvectors (default) or #EigenvaluesOnly.
  214. * \returns Reference to \c *this
  215. *
  216. * This function assumes that the matrix has been reduced to tridiagonal form.
  217. *
  218. * \sa compute(const MatrixType&, int) for more information
  219. */
  220. SelfAdjointEigenSolver& computeFromTridiagonal(const RealVectorType& diag, const SubDiagonalType& subdiag , int options=ComputeEigenvectors);
  221. /** \brief Returns the eigenvectors of given matrix.
  222. *
  223. * \returns A const reference to the matrix whose columns are the eigenvectors.
  224. *
  225. * \pre The eigenvectors have been computed before.
  226. *
  227. * Column \f$ k \f$ of the returned matrix is an eigenvector corresponding
  228. * to eigenvalue number \f$ k \f$ as returned by eigenvalues(). The
  229. * eigenvectors are normalized to have (Euclidean) norm equal to one. If
  230. * this object was used to solve the eigenproblem for the selfadjoint
  231. * matrix \f$ A \f$, then the matrix returned by this function is the
  232. * matrix \f$ V \f$ in the eigendecomposition \f$ A = V D V^{-1} \f$.
  233. *
  234. * Example: \include SelfAdjointEigenSolver_eigenvectors.cpp
  235. * Output: \verbinclude SelfAdjointEigenSolver_eigenvectors.out
  236. *
  237. * \sa eigenvalues()
  238. */
  239. EIGEN_DEVICE_FUNC
  240. const EigenvectorsType& eigenvectors() const
  241. {
  242. eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
  243. eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
  244. return m_eivec;
  245. }
  246. /** \brief Returns the eigenvalues of given matrix.
  247. *
  248. * \returns A const reference to the column vector containing the eigenvalues.
  249. *
  250. * \pre The eigenvalues have been computed before.
  251. *
  252. * The eigenvalues are repeated according to their algebraic multiplicity,
  253. * so there are as many eigenvalues as rows in the matrix. The eigenvalues
  254. * are sorted in increasing order.
  255. *
  256. * Example: \include SelfAdjointEigenSolver_eigenvalues.cpp
  257. * Output: \verbinclude SelfAdjointEigenSolver_eigenvalues.out
  258. *
  259. * \sa eigenvectors(), MatrixBase::eigenvalues()
  260. */
  261. EIGEN_DEVICE_FUNC
  262. const RealVectorType& eigenvalues() const
  263. {
  264. eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
  265. return m_eivalues;
  266. }
  267. /** \brief Computes the positive-definite square root of the matrix.
  268. *
  269. * \returns the positive-definite square root of the matrix
  270. *
  271. * \pre The eigenvalues and eigenvectors of a positive-definite matrix
  272. * have been computed before.
  273. *
  274. * The square root of a positive-definite matrix \f$ A \f$ is the
  275. * positive-definite matrix whose square equals \f$ A \f$. This function
  276. * uses the eigendecomposition \f$ A = V D V^{-1} \f$ to compute the
  277. * square root as \f$ A^{1/2} = V D^{1/2} V^{-1} \f$.
  278. *
  279. * Example: \include SelfAdjointEigenSolver_operatorSqrt.cpp
  280. * Output: \verbinclude SelfAdjointEigenSolver_operatorSqrt.out
  281. *
  282. * \sa operatorInverseSqrt(), <a href="unsupported/group__MatrixFunctions__Module.html">MatrixFunctions Module</a>
  283. */
  284. EIGEN_DEVICE_FUNC
  285. MatrixType operatorSqrt() const
  286. {
  287. eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
  288. eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
  289. return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
  290. }
  291. /** \brief Computes the inverse square root of the matrix.
  292. *
  293. * \returns the inverse positive-definite square root of the matrix
  294. *
  295. * \pre The eigenvalues and eigenvectors of a positive-definite matrix
  296. * have been computed before.
  297. *
  298. * This function uses the eigendecomposition \f$ A = V D V^{-1} \f$ to
  299. * compute the inverse square root as \f$ V D^{-1/2} V^{-1} \f$. This is
  300. * cheaper than first computing the square root with operatorSqrt() and
  301. * then its inverse with MatrixBase::inverse().
  302. *
  303. * Example: \include SelfAdjointEigenSolver_operatorInverseSqrt.cpp
  304. * Output: \verbinclude SelfAdjointEigenSolver_operatorInverseSqrt.out
  305. *
  306. * \sa operatorSqrt(), MatrixBase::inverse(), <a href="unsupported/group__MatrixFunctions__Module.html">MatrixFunctions Module</a>
  307. */
  308. EIGEN_DEVICE_FUNC
  309. MatrixType operatorInverseSqrt() const
  310. {
  311. eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
  312. eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
  313. return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
  314. }
  315. /** \brief Reports whether previous computation was successful.
  316. *
  317. * \returns \c Success if computation was succesful, \c NoConvergence otherwise.
  318. */
  319. EIGEN_DEVICE_FUNC
  320. ComputationInfo info() const
  321. {
  322. eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
  323. return m_info;
  324. }
  325. /** \brief Maximum number of iterations.
  326. *
  327. * The algorithm terminates if it does not converge within m_maxIterations * n iterations, where n
  328. * denotes the size of the matrix. This value is currently set to 30 (copied from LAPACK).
  329. */
  330. static const int m_maxIterations = 30;
  331. protected:
  332. static void check_template_parameters()
  333. {
  334. EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
  335. }
  336. EigenvectorsType m_eivec;
  337. RealVectorType m_eivalues;
  338. typename TridiagonalizationType::SubDiagonalType m_subdiag;
  339. ComputationInfo m_info;
  340. bool m_isInitialized;
  341. bool m_eigenvectorsOk;
  342. };
  343. namespace internal {
  344. /** \internal
  345. *
  346. * \eigenvalues_module \ingroup Eigenvalues_Module
  347. *
  348. * Performs a QR step on a tridiagonal symmetric matrix represented as a
  349. * pair of two vectors \a diag and \a subdiag.
  350. *
  351. * \param diag the diagonal part of the input selfadjoint tridiagonal matrix
  352. * \param subdiag the sub-diagonal part of the input selfadjoint tridiagonal matrix
  353. * \param start starting index of the submatrix to work on
  354. * \param end last+1 index of the submatrix to work on
  355. * \param matrixQ pointer to the column-major matrix holding the eigenvectors, can be 0
  356. * \param n size of the input matrix
  357. *
  358. * For compilation efficiency reasons, this procedure does not use eigen expression
  359. * for its arguments.
  360. *
  361. * Implemented from Golub's "Matrix Computations", algorithm 8.3.2:
  362. * "implicit symmetric QR step with Wilkinson shift"
  363. */
  364. template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
  365. EIGEN_DEVICE_FUNC
  366. static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n);
  367. }
  368. template<typename MatrixType>
  369. template<typename InputType>
  370. EIGEN_DEVICE_FUNC
  371. SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
  372. ::compute(const EigenBase<InputType>& a_matrix, int options)
  373. {
  374. check_template_parameters();
  375. const InputType &matrix(a_matrix.derived());
  376. using std::abs;
  377. eigen_assert(matrix.cols() == matrix.rows());
  378. eigen_assert((options&~(EigVecMask|GenEigMask))==0
  379. && (options&EigVecMask)!=EigVecMask
  380. && "invalid option parameter");
  381. bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
  382. Index n = matrix.cols();
  383. m_eivalues.resize(n,1);
  384. if(n==1)
  385. {
  386. m_eivec = matrix;
  387. m_eivalues.coeffRef(0,0) = numext::real(m_eivec.coeff(0,0));
  388. if(computeEigenvectors)
  389. m_eivec.setOnes(n,n);
  390. m_info = Success;
  391. m_isInitialized = true;
  392. m_eigenvectorsOk = computeEigenvectors;
  393. return *this;
  394. }
  395. // declare some aliases
  396. RealVectorType& diag = m_eivalues;
  397. EigenvectorsType& mat = m_eivec;
  398. // map the matrix coefficients to [-1:1] to avoid over- and underflow.
  399. mat = matrix.template triangularView<Lower>();
  400. RealScalar scale = mat.cwiseAbs().maxCoeff();
  401. if(scale==RealScalar(0)) scale = RealScalar(1);
  402. mat.template triangularView<Lower>() /= scale;
  403. m_subdiag.resize(n-1);
  404. internal::tridiagonalization_inplace(mat, diag, m_subdiag, computeEigenvectors);
  405. m_info = internal::computeFromTridiagonal_impl(diag, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
  406. // scale back the eigen values
  407. m_eivalues *= scale;
  408. m_isInitialized = true;
  409. m_eigenvectorsOk = computeEigenvectors;
  410. return *this;
  411. }
  412. template<typename MatrixType>
  413. SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
  414. ::computeFromTridiagonal(const RealVectorType& diag, const SubDiagonalType& subdiag , int options)
  415. {
  416. //TODO : Add an option to scale the values beforehand
  417. bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
  418. m_eivalues = diag;
  419. m_subdiag = subdiag;
  420. if (computeEigenvectors)
  421. {
  422. m_eivec.setIdentity(diag.size(), diag.size());
  423. }
  424. m_info = internal::computeFromTridiagonal_impl(m_eivalues, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
  425. m_isInitialized = true;
  426. m_eigenvectorsOk = computeEigenvectors;
  427. return *this;
  428. }
  429. namespace internal {
  430. /**
  431. * \internal
  432. * \brief Compute the eigendecomposition from a tridiagonal matrix
  433. *
  434. * \param[in,out] diag : On input, the diagonal of the matrix, on output the eigenvalues
  435. * \param[in,out] subdiag : The subdiagonal part of the matrix (entries are modified during the decomposition)
  436. * \param[in] maxIterations : the maximum number of iterations
  437. * \param[in] computeEigenvectors : whether the eigenvectors have to be computed or not
  438. * \param[out] eivec : The matrix to store the eigenvectors if computeEigenvectors==true. Must be allocated on input.
  439. * \returns \c Success or \c NoConvergence
  440. */
  441. template<typename MatrixType, typename DiagType, typename SubDiagType>
  442. ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType& eivec)
  443. {
  444. using std::abs;
  445. ComputationInfo info;
  446. typedef typename MatrixType::Scalar Scalar;
  447. Index n = diag.size();
  448. Index end = n-1;
  449. Index start = 0;
  450. Index iter = 0; // total number of iterations
  451. typedef typename DiagType::RealScalar RealScalar;
  452. const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
  453. const RealScalar precision = RealScalar(2)*NumTraits<RealScalar>::epsilon();
  454. while (end>0)
  455. {
  456. for (Index i = start; i<end; ++i)
  457. if (internal::isMuchSmallerThan(abs(subdiag[i]),(abs(diag[i])+abs(diag[i+1])),precision) || abs(subdiag[i]) <= considerAsZero)
  458. subdiag[i] = 0;
  459. // find the largest unreduced block
  460. while (end>0 && subdiag[end-1]==RealScalar(0))
  461. {
  462. end--;
  463. }
  464. if (end<=0)
  465. break;
  466. // if we spent too many iterations, we give up
  467. iter++;
  468. if(iter > maxIterations * n) break;
  469. start = end - 1;
  470. while (start>0 && subdiag[start-1]!=0)
  471. start--;
  472. internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), subdiag.data(), start, end, computeEigenvectors ? eivec.data() : (Scalar*)0, n);
  473. }
  474. if (iter <= maxIterations * n)
  475. info = Success;
  476. else
  477. info = NoConvergence;
  478. // Sort eigenvalues and corresponding vectors.
  479. // TODO make the sort optional ?
  480. // TODO use a better sort algorithm !!
  481. if (info == Success)
  482. {
  483. for (Index i = 0; i < n-1; ++i)
  484. {
  485. Index k;
  486. diag.segment(i,n-i).minCoeff(&k);
  487. if (k > 0)
  488. {
  489. std::swap(diag[i], diag[k+i]);
  490. if(computeEigenvectors)
  491. eivec.col(i).swap(eivec.col(k+i));
  492. }
  493. }
  494. }
  495. return info;
  496. }
  497. template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues
  498. {
  499. EIGEN_DEVICE_FUNC
  500. static inline void run(SolverType& eig, const typename SolverType::MatrixType& A, int options)
  501. { eig.compute(A,options); }
  502. };
  503. template<typename SolverType> struct direct_selfadjoint_eigenvalues<SolverType,3,false>
  504. {
  505. typedef typename SolverType::MatrixType MatrixType;
  506. typedef typename SolverType::RealVectorType VectorType;
  507. typedef typename SolverType::Scalar Scalar;
  508. typedef typename SolverType::EigenvectorsType EigenvectorsType;
  509. /** \internal
  510. * Computes the roots of the characteristic polynomial of \a m.
  511. * For numerical stability m.trace() should be near zero and to avoid over- or underflow m should be normalized.
  512. */
  513. EIGEN_DEVICE_FUNC
  514. static inline void computeRoots(const MatrixType& m, VectorType& roots)
  515. {
  516. EIGEN_USING_STD_MATH(sqrt)
  517. EIGEN_USING_STD_MATH(atan2)
  518. EIGEN_USING_STD_MATH(cos)
  519. EIGEN_USING_STD_MATH(sin)
  520. const Scalar s_inv3 = Scalar(1)/Scalar(3);
  521. const Scalar s_sqrt3 = sqrt(Scalar(3));
  522. // The characteristic equation is x^3 - c2*x^2 + c1*x - c0 = 0. The
  523. // eigenvalues are the roots to this equation, all guaranteed to be
  524. // real-valued, because the matrix is symmetric.
  525. Scalar c0 = m(0,0)*m(1,1)*m(2,2) + Scalar(2)*m(1,0)*m(2,0)*m(2,1) - m(0,0)*m(2,1)*m(2,1) - m(1,1)*m(2,0)*m(2,0) - m(2,2)*m(1,0)*m(1,0);
  526. Scalar c1 = m(0,0)*m(1,1) - m(1,0)*m(1,0) + m(0,0)*m(2,2) - m(2,0)*m(2,0) + m(1,1)*m(2,2) - m(2,1)*m(2,1);
  527. Scalar c2 = m(0,0) + m(1,1) + m(2,2);
  528. // Construct the parameters used in classifying the roots of the equation
  529. // and in solving the equation for the roots in closed form.
  530. Scalar c2_over_3 = c2*s_inv3;
  531. Scalar a_over_3 = (c2*c2_over_3 - c1)*s_inv3;
  532. a_over_3 = numext::maxi(a_over_3, Scalar(0));
  533. Scalar half_b = Scalar(0.5)*(c0 + c2_over_3*(Scalar(2)*c2_over_3*c2_over_3 - c1));
  534. Scalar q = a_over_3*a_over_3*a_over_3 - half_b*half_b;
  535. q = numext::maxi(q, Scalar(0));
  536. // Compute the eigenvalues by solving for the roots of the polynomial.
  537. Scalar rho = sqrt(a_over_3);
  538. Scalar theta = atan2(sqrt(q),half_b)*s_inv3; // since sqrt(q) > 0, atan2 is in [0, pi] and theta is in [0, pi/3]
  539. Scalar cos_theta = cos(theta);
  540. Scalar sin_theta = sin(theta);
  541. // roots are already sorted, since cos is monotonically decreasing on [0, pi]
  542. roots(0) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta); // == 2*rho*cos(theta+2pi/3)
  543. roots(1) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta); // == 2*rho*cos(theta+ pi/3)
  544. roots(2) = c2_over_3 + Scalar(2)*rho*cos_theta;
  545. }
  546. EIGEN_DEVICE_FUNC
  547. static inline bool extract_kernel(MatrixType& mat, Ref<VectorType> res, Ref<VectorType> representative)
  548. {
  549. EIGEN_USING_STD_MATH(sqrt)
  550. EIGEN_USING_STD_MATH(abs)
  551. Index i0;
  552. // Find non-zero column i0 (by construction, there must exist a non zero coefficient on the diagonal):
  553. mat.diagonal().cwiseAbs().maxCoeff(&i0);
  554. // mat.col(i0) is a good candidate for an orthogonal vector to the current eigenvector,
  555. // so let's save it:
  556. representative = mat.col(i0);
  557. Scalar n0, n1;
  558. VectorType c0, c1;
  559. n0 = (c0 = representative.cross(mat.col((i0+1)%3))).squaredNorm();
  560. n1 = (c1 = representative.cross(mat.col((i0+2)%3))).squaredNorm();
  561. if(n0>n1) res = c0/sqrt(n0);
  562. else res = c1/sqrt(n1);
  563. return true;
  564. }
  565. EIGEN_DEVICE_FUNC
  566. static inline void run(SolverType& solver, const MatrixType& mat, int options)
  567. {
  568. eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
  569. eigen_assert((options&~(EigVecMask|GenEigMask))==0
  570. && (options&EigVecMask)!=EigVecMask
  571. && "invalid option parameter");
  572. bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
  573. EigenvectorsType& eivecs = solver.m_eivec;
  574. VectorType& eivals = solver.m_eivalues;
  575. // Shift the matrix to the mean eigenvalue and map the matrix coefficients to [-1:1] to avoid over- and underflow.
  576. Scalar shift = mat.trace() / Scalar(3);
  577. // TODO Avoid this copy. Currently it is necessary to suppress bogus values when determining maxCoeff and for computing the eigenvectors later
  578. MatrixType scaledMat = mat.template selfadjointView<Lower>();
  579. scaledMat.diagonal().array() -= shift;
  580. Scalar scale = scaledMat.cwiseAbs().maxCoeff();
  581. if(scale > 0) scaledMat /= scale; // TODO for scale==0 we could save the remaining operations
  582. // compute the eigenvalues
  583. computeRoots(scaledMat,eivals);
  584. // compute the eigenvectors
  585. if(computeEigenvectors)
  586. {
  587. if((eivals(2)-eivals(0))<=Eigen::NumTraits<Scalar>::epsilon())
  588. {
  589. // All three eigenvalues are numerically the same
  590. eivecs.setIdentity();
  591. }
  592. else
  593. {
  594. MatrixType tmp;
  595. tmp = scaledMat;
  596. // Compute the eigenvector of the most distinct eigenvalue
  597. Scalar d0 = eivals(2) - eivals(1);
  598. Scalar d1 = eivals(1) - eivals(0);
  599. Index k(0), l(2);
  600. if(d0 > d1)
  601. {
  602. numext::swap(k,l);
  603. d0 = d1;
  604. }
  605. // Compute the eigenvector of index k
  606. {
  607. tmp.diagonal().array () -= eivals(k);
  608. // By construction, 'tmp' is of rank 2, and its kernel corresponds to the respective eigenvector.
  609. extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
  610. }
  611. // Compute eigenvector of index l
  612. if(d0<=2*Eigen::NumTraits<Scalar>::epsilon()*d1)
  613. {
  614. // If d0 is too small, then the two other eigenvalues are numerically the same,
  615. // and thus we only have to ortho-normalize the near orthogonal vector we saved above.
  616. eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l))*eivecs.col(l);
  617. eivecs.col(l).normalize();
  618. }
  619. else
  620. {
  621. tmp = scaledMat;
  622. tmp.diagonal().array () -= eivals(l);
  623. VectorType dummy;
  624. extract_kernel(tmp, eivecs.col(l), dummy);
  625. }
  626. // Compute last eigenvector from the other two
  627. eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
  628. }
  629. }
  630. // Rescale back to the original size.
  631. eivals *= scale;
  632. eivals.array() += shift;
  633. solver.m_info = Success;
  634. solver.m_isInitialized = true;
  635. solver.m_eigenvectorsOk = computeEigenvectors;
  636. }
  637. };
  638. // 2x2 direct eigenvalues decomposition, code from Hauke Heibel
  639. template<typename SolverType>
  640. struct direct_selfadjoint_eigenvalues<SolverType,2,false>
  641. {
  642. typedef typename SolverType::MatrixType MatrixType;
  643. typedef typename SolverType::RealVectorType VectorType;
  644. typedef typename SolverType::Scalar Scalar;
  645. typedef typename SolverType::EigenvectorsType EigenvectorsType;
  646. EIGEN_DEVICE_FUNC
  647. static inline void computeRoots(const MatrixType& m, VectorType& roots)
  648. {
  649. using std::sqrt;
  650. const Scalar t0 = Scalar(0.5) * sqrt( numext::abs2(m(0,0)-m(1,1)) + Scalar(4)*numext::abs2(m(1,0)));
  651. const Scalar t1 = Scalar(0.5) * (m(0,0) + m(1,1));
  652. roots(0) = t1 - t0;
  653. roots(1) = t1 + t0;
  654. }
  655. EIGEN_DEVICE_FUNC
  656. static inline void run(SolverType& solver, const MatrixType& mat, int options)
  657. {
  658. EIGEN_USING_STD_MATH(sqrt);
  659. EIGEN_USING_STD_MATH(abs);
  660. eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
  661. eigen_assert((options&~(EigVecMask|GenEigMask))==0
  662. && (options&EigVecMask)!=EigVecMask
  663. && "invalid option parameter");
  664. bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
  665. EigenvectorsType& eivecs = solver.m_eivec;
  666. VectorType& eivals = solver.m_eivalues;
  667. // Shift the matrix to the mean eigenvalue and map the matrix coefficients to [-1:1] to avoid over- and underflow.
  668. Scalar shift = mat.trace() / Scalar(2);
  669. MatrixType scaledMat = mat;
  670. scaledMat.coeffRef(0,1) = mat.coeff(1,0);
  671. scaledMat.diagonal().array() -= shift;
  672. Scalar scale = scaledMat.cwiseAbs().maxCoeff();
  673. if(scale > Scalar(0))
  674. scaledMat /= scale;
  675. // Compute the eigenvalues
  676. computeRoots(scaledMat,eivals);
  677. // compute the eigen vectors
  678. if(computeEigenvectors)
  679. {
  680. if((eivals(1)-eivals(0))<=abs(eivals(1))*Eigen::NumTraits<Scalar>::epsilon())
  681. {
  682. eivecs.setIdentity();
  683. }
  684. else
  685. {
  686. scaledMat.diagonal().array () -= eivals(1);
  687. Scalar a2 = numext::abs2(scaledMat(0,0));
  688. Scalar c2 = numext::abs2(scaledMat(1,1));
  689. Scalar b2 = numext::abs2(scaledMat(1,0));
  690. if(a2>c2)
  691. {
  692. eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
  693. eivecs.col(1) /= sqrt(a2+b2);
  694. }
  695. else
  696. {
  697. eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
  698. eivecs.col(1) /= sqrt(c2+b2);
  699. }
  700. eivecs.col(0) << eivecs.col(1).unitOrthogonal();
  701. }
  702. }
  703. // Rescale back to the original size.
  704. eivals *= scale;
  705. eivals.array() += shift;
  706. solver.m_info = Success;
  707. solver.m_isInitialized = true;
  708. solver.m_eigenvectorsOk = computeEigenvectors;
  709. }
  710. };
  711. }
  712. template<typename MatrixType>
  713. EIGEN_DEVICE_FUNC
  714. SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
  715. ::computeDirect(const MatrixType& matrix, int options)
  716. {
  717. internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>::run(*this,matrix,options);
  718. return *this;
  719. }
  720. namespace internal {
  721. template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
  722. EIGEN_DEVICE_FUNC
  723. static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n)
  724. {
  725. using std::abs;
  726. RealScalar td = (diag[end-1] - diag[end])*RealScalar(0.5);
  727. RealScalar e = subdiag[end-1];
  728. // Note that thanks to scaling, e^2 or td^2 cannot overflow, however they can still
  729. // underflow thus leading to inf/NaN values when using the following commented code:
  730. // RealScalar e2 = numext::abs2(subdiag[end-1]);
  731. // RealScalar mu = diag[end] - e2 / (td + (td>0 ? 1 : -1) * sqrt(td*td + e2));
  732. // This explain the following, somewhat more complicated, version:
  733. RealScalar mu = diag[end];
  734. if(td==RealScalar(0))
  735. mu -= abs(e);
  736. else
  737. {
  738. RealScalar e2 = numext::abs2(subdiag[end-1]);
  739. RealScalar h = numext::hypot(td,e);
  740. if(e2==RealScalar(0)) mu -= (e / (td + (td>RealScalar(0) ? RealScalar(1) : RealScalar(-1)))) * (e / h);
  741. else mu -= e2 / (td + (td>RealScalar(0) ? h : -h));
  742. }
  743. RealScalar x = diag[start] - mu;
  744. RealScalar z = subdiag[start];
  745. for (Index k = start; k < end; ++k)
  746. {
  747. JacobiRotation<RealScalar> rot;
  748. rot.makeGivens(x, z);
  749. // do T = G' T G
  750. RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
  751. RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k+1];
  752. diag[k] = rot.c() * (rot.c() * diag[k] - rot.s() * subdiag[k]) - rot.s() * (rot.c() * subdiag[k] - rot.s() * diag[k+1]);
  753. diag[k+1] = rot.s() * sdk + rot.c() * dkp1;
  754. subdiag[k] = rot.c() * sdk - rot.s() * dkp1;
  755. if (k > start)
  756. subdiag[k - 1] = rot.c() * subdiag[k-1] - rot.s() * z;
  757. x = subdiag[k];
  758. if (k < end - 1)
  759. {
  760. z = -rot.s() * subdiag[k+1];
  761. subdiag[k + 1] = rot.c() * subdiag[k+1];
  762. }
  763. // apply the givens rotation to the unit matrix Q = Q * G
  764. if (matrixQ)
  765. {
  766. // FIXME if StorageOrder == RowMajor this operation is not very efficient
  767. Map<Matrix<Scalar,Dynamic,Dynamic,StorageOrder> > q(matrixQ,n,n);
  768. q.applyOnTheRight(k,k+1,rot);
  769. }
  770. }
  771. }
  772. } // end namespace internal
  773. } // end namespace Eigen
  774. #endif // EIGEN_SELFADJOINTEIGENSOLVER_H