EigenSolver.h 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. // Copyright (C) 2010,2012 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_EIGENSOLVER_H
  11. #define EIGEN_EIGENSOLVER_H
  12. #include "./RealSchur.h"
  13. namespace Eigen {
  14. /** \eigenvalues_module \ingroup Eigenvalues_Module
  15. *
  16. *
  17. * \class EigenSolver
  18. *
  19. * \brief Computes eigenvalues and eigenvectors of general matrices
  20. *
  21. * \tparam _MatrixType the type of the matrix of which we are computing the
  22. * eigendecomposition; this is expected to be an instantiation of the Matrix
  23. * class template. Currently, only real matrices are supported.
  24. *
  25. * The eigenvalues and eigenvectors of a matrix \f$ A \f$ are scalars
  26. * \f$ \lambda \f$ and vectors \f$ v \f$ such that \f$ Av = \lambda v \f$. If
  27. * \f$ D \f$ is a diagonal matrix with the eigenvalues on the diagonal, and
  28. * \f$ V \f$ is a matrix with the eigenvectors as its columns, then \f$ A V =
  29. * V D \f$. The matrix \f$ V \f$ is almost always invertible, in which case we
  30. * have \f$ A = V D V^{-1} \f$. This is called the eigendecomposition.
  31. *
  32. * The eigenvalues and eigenvectors of a matrix may be complex, even when the
  33. * matrix is real. However, we can choose real matrices \f$ V \f$ and \f$ D
  34. * \f$ satisfying \f$ A V = V D \f$, just like the eigendecomposition, if the
  35. * matrix \f$ D \f$ is not required to be diagonal, but if it is allowed to
  36. * have blocks of the form
  37. * \f[ \begin{bmatrix} u & v \\ -v & u \end{bmatrix} \f]
  38. * (where \f$ u \f$ and \f$ v \f$ are real numbers) on the diagonal. These
  39. * blocks correspond to complex eigenvalue pairs \f$ u \pm iv \f$. We call
  40. * this variant of the eigendecomposition the pseudo-eigendecomposition.
  41. *
  42. * Call the function compute() to compute the eigenvalues and eigenvectors of
  43. * a given matrix. Alternatively, you can use the
  44. * EigenSolver(const MatrixType&, bool) constructor which computes the
  45. * eigenvalues and eigenvectors at construction time. Once the eigenvalue and
  46. * eigenvectors are computed, they can be retrieved with the eigenvalues() and
  47. * eigenvectors() functions. The pseudoEigenvalueMatrix() and
  48. * pseudoEigenvectors() methods allow the construction of the
  49. * pseudo-eigendecomposition.
  50. *
  51. * The documentation for EigenSolver(const MatrixType&, bool) contains an
  52. * example of the typical use of this class.
  53. *
  54. * \note The implementation is adapted from
  55. * <a href="http://math.nist.gov/javanumerics/jama/">JAMA</a> (public domain).
  56. * Their code is based on EISPACK.
  57. *
  58. * \sa MatrixBase::eigenvalues(), class ComplexEigenSolver, class SelfAdjointEigenSolver
  59. */
  60. template<typename _MatrixType> class EigenSolver
  61. {
  62. public:
  63. /** \brief Synonym for the template parameter \p _MatrixType. */
  64. typedef _MatrixType MatrixType;
  65. enum {
  66. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  67. ColsAtCompileTime = MatrixType::ColsAtCompileTime,
  68. Options = MatrixType::Options,
  69. MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  70. MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
  71. };
  72. /** \brief Scalar type for matrices of type #MatrixType. */
  73. typedef typename MatrixType::Scalar Scalar;
  74. typedef typename NumTraits<Scalar>::Real RealScalar;
  75. typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
  76. /** \brief Complex scalar type for #MatrixType.
  77. *
  78. * This is \c std::complex<Scalar> if #Scalar is real (e.g.,
  79. * \c float or \c double) and just \c Scalar if #Scalar is
  80. * complex.
  81. */
  82. typedef std::complex<RealScalar> ComplexScalar;
  83. /** \brief Type for vector of eigenvalues as returned by eigenvalues().
  84. *
  85. * This is a column vector with entries of type #ComplexScalar.
  86. * The length of the vector is the size of #MatrixType.
  87. */
  88. typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> EigenvalueType;
  89. /** \brief Type for matrix of eigenvectors as returned by eigenvectors().
  90. *
  91. * This is a square matrix with entries of type #ComplexScalar.
  92. * The size is the same as the size of #MatrixType.
  93. */
  94. typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> EigenvectorsType;
  95. /** \brief Default constructor.
  96. *
  97. * The default constructor is useful in cases in which the user intends to
  98. * perform decompositions via EigenSolver::compute(const MatrixType&, bool).
  99. *
  100. * \sa compute() for an example.
  101. */
  102. EigenSolver() : m_eivec(), m_eivalues(), m_isInitialized(false), m_realSchur(), m_matT(), m_tmp() {}
  103. /** \brief Default constructor with memory preallocation
  104. *
  105. * Like the default constructor but with preallocation of the internal data
  106. * according to the specified problem \a size.
  107. * \sa EigenSolver()
  108. */
  109. explicit EigenSolver(Index size)
  110. : m_eivec(size, size),
  111. m_eivalues(size),
  112. m_isInitialized(false),
  113. m_eigenvectorsOk(false),
  114. m_realSchur(size),
  115. m_matT(size, size),
  116. m_tmp(size)
  117. {}
  118. /** \brief Constructor; computes eigendecomposition of given matrix.
  119. *
  120. * \param[in] matrix Square matrix whose eigendecomposition is to be computed.
  121. * \param[in] computeEigenvectors If true, both the eigenvectors and the
  122. * eigenvalues are computed; if false, only the eigenvalues are
  123. * computed.
  124. *
  125. * This constructor calls compute() to compute the eigenvalues
  126. * and eigenvectors.
  127. *
  128. * Example: \include EigenSolver_EigenSolver_MatrixType.cpp
  129. * Output: \verbinclude EigenSolver_EigenSolver_MatrixType.out
  130. *
  131. * \sa compute()
  132. */
  133. template<typename InputType>
  134. explicit EigenSolver(const EigenBase<InputType>& matrix, bool computeEigenvectors = true)
  135. : m_eivec(matrix.rows(), matrix.cols()),
  136. m_eivalues(matrix.cols()),
  137. m_isInitialized(false),
  138. m_eigenvectorsOk(false),
  139. m_realSchur(matrix.cols()),
  140. m_matT(matrix.rows(), matrix.cols()),
  141. m_tmp(matrix.cols())
  142. {
  143. compute(matrix.derived(), computeEigenvectors);
  144. }
  145. /** \brief Returns the eigenvectors of given matrix.
  146. *
  147. * \returns %Matrix whose columns are the (possibly complex) eigenvectors.
  148. *
  149. * \pre Either the constructor
  150. * EigenSolver(const MatrixType&,bool) or the member function
  151. * compute(const MatrixType&, bool) has been called before, and
  152. * \p computeEigenvectors was set to true (the default).
  153. *
  154. * Column \f$ k \f$ of the returned matrix is an eigenvector corresponding
  155. * to eigenvalue number \f$ k \f$ as returned by eigenvalues(). The
  156. * eigenvectors are normalized to have (Euclidean) norm equal to one. The
  157. * matrix returned by this function is the matrix \f$ V \f$ in the
  158. * eigendecomposition \f$ A = V D V^{-1} \f$, if it exists.
  159. *
  160. * Example: \include EigenSolver_eigenvectors.cpp
  161. * Output: \verbinclude EigenSolver_eigenvectors.out
  162. *
  163. * \sa eigenvalues(), pseudoEigenvectors()
  164. */
  165. EigenvectorsType eigenvectors() const;
  166. /** \brief Returns the pseudo-eigenvectors of given matrix.
  167. *
  168. * \returns Const reference to matrix whose columns are the pseudo-eigenvectors.
  169. *
  170. * \pre Either the constructor
  171. * EigenSolver(const MatrixType&,bool) or the member function
  172. * compute(const MatrixType&, bool) has been called before, and
  173. * \p computeEigenvectors was set to true (the default).
  174. *
  175. * The real matrix \f$ V \f$ returned by this function and the
  176. * block-diagonal matrix \f$ D \f$ returned by pseudoEigenvalueMatrix()
  177. * satisfy \f$ AV = VD \f$.
  178. *
  179. * Example: \include EigenSolver_pseudoEigenvectors.cpp
  180. * Output: \verbinclude EigenSolver_pseudoEigenvectors.out
  181. *
  182. * \sa pseudoEigenvalueMatrix(), eigenvectors()
  183. */
  184. const MatrixType& pseudoEigenvectors() const
  185. {
  186. eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
  187. eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
  188. return m_eivec;
  189. }
  190. /** \brief Returns the block-diagonal matrix in the pseudo-eigendecomposition.
  191. *
  192. * \returns A block-diagonal matrix.
  193. *
  194. * \pre Either the constructor
  195. * EigenSolver(const MatrixType&,bool) or the member function
  196. * compute(const MatrixType&, bool) has been called before.
  197. *
  198. * The matrix \f$ D \f$ returned by this function is real and
  199. * block-diagonal. The blocks on the diagonal are either 1-by-1 or 2-by-2
  200. * blocks of the form
  201. * \f$ \begin{bmatrix} u & v \\ -v & u \end{bmatrix} \f$.
  202. * These blocks are not sorted in any particular order.
  203. * The matrix \f$ D \f$ and the matrix \f$ V \f$ returned by
  204. * pseudoEigenvectors() satisfy \f$ AV = VD \f$.
  205. *
  206. * \sa pseudoEigenvectors() for an example, eigenvalues()
  207. */
  208. MatrixType pseudoEigenvalueMatrix() const;
  209. /** \brief Returns the eigenvalues of given matrix.
  210. *
  211. * \returns A const reference to the column vector containing the eigenvalues.
  212. *
  213. * \pre Either the constructor
  214. * EigenSolver(const MatrixType&,bool) or the member function
  215. * compute(const MatrixType&, bool) has been called before.
  216. *
  217. * The eigenvalues are repeated according to their algebraic multiplicity,
  218. * so there are as many eigenvalues as rows in the matrix. The eigenvalues
  219. * are not sorted in any particular order.
  220. *
  221. * Example: \include EigenSolver_eigenvalues.cpp
  222. * Output: \verbinclude EigenSolver_eigenvalues.out
  223. *
  224. * \sa eigenvectors(), pseudoEigenvalueMatrix(),
  225. * MatrixBase::eigenvalues()
  226. */
  227. const EigenvalueType& eigenvalues() const
  228. {
  229. eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
  230. return m_eivalues;
  231. }
  232. /** \brief Computes eigendecomposition of given matrix.
  233. *
  234. * \param[in] matrix Square matrix whose eigendecomposition is to be computed.
  235. * \param[in] computeEigenvectors If true, both the eigenvectors and the
  236. * eigenvalues are computed; if false, only the eigenvalues are
  237. * computed.
  238. * \returns Reference to \c *this
  239. *
  240. * This function computes the eigenvalues of the real matrix \p matrix.
  241. * The eigenvalues() function can be used to retrieve them. If
  242. * \p computeEigenvectors is true, then the eigenvectors are also computed
  243. * and can be retrieved by calling eigenvectors().
  244. *
  245. * The matrix is first reduced to real Schur form using the RealSchur
  246. * class. The Schur decomposition is then used to compute the eigenvalues
  247. * and eigenvectors.
  248. *
  249. * The cost of the computation is dominated by the cost of the
  250. * Schur decomposition, which is very approximately \f$ 25n^3 \f$
  251. * (where \f$ n \f$ is the size of the matrix) if \p computeEigenvectors
  252. * is true, and \f$ 10n^3 \f$ if \p computeEigenvectors is false.
  253. *
  254. * This method reuses of the allocated data in the EigenSolver object.
  255. *
  256. * Example: \include EigenSolver_compute.cpp
  257. * Output: \verbinclude EigenSolver_compute.out
  258. */
  259. template<typename InputType>
  260. EigenSolver& compute(const EigenBase<InputType>& matrix, bool computeEigenvectors = true);
  261. /** \returns NumericalIssue if the input contains INF or NaN values or overflow occured. Returns Success otherwise. */
  262. ComputationInfo info() const
  263. {
  264. eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
  265. return m_info;
  266. }
  267. /** \brief Sets the maximum number of iterations allowed. */
  268. EigenSolver& setMaxIterations(Index maxIters)
  269. {
  270. m_realSchur.setMaxIterations(maxIters);
  271. return *this;
  272. }
  273. /** \brief Returns the maximum number of iterations. */
  274. Index getMaxIterations()
  275. {
  276. return m_realSchur.getMaxIterations();
  277. }
  278. private:
  279. void doComputeEigenvectors();
  280. protected:
  281. static void check_template_parameters()
  282. {
  283. EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
  284. EIGEN_STATIC_ASSERT(!NumTraits<Scalar>::IsComplex, NUMERIC_TYPE_MUST_BE_REAL);
  285. }
  286. MatrixType m_eivec;
  287. EigenvalueType m_eivalues;
  288. bool m_isInitialized;
  289. bool m_eigenvectorsOk;
  290. ComputationInfo m_info;
  291. RealSchur<MatrixType> m_realSchur;
  292. MatrixType m_matT;
  293. typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
  294. ColumnVectorType m_tmp;
  295. };
  296. template<typename MatrixType>
  297. MatrixType EigenSolver<MatrixType>::pseudoEigenvalueMatrix() const
  298. {
  299. eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
  300. const RealScalar precision = RealScalar(2)*NumTraits<RealScalar>::epsilon();
  301. Index n = m_eivalues.rows();
  302. MatrixType matD = MatrixType::Zero(n,n);
  303. for (Index i=0; i<n; ++i)
  304. {
  305. if (internal::isMuchSmallerThan(numext::imag(m_eivalues.coeff(i)), numext::real(m_eivalues.coeff(i)), precision))
  306. matD.coeffRef(i,i) = numext::real(m_eivalues.coeff(i));
  307. else
  308. {
  309. matD.template block<2,2>(i,i) << numext::real(m_eivalues.coeff(i)), numext::imag(m_eivalues.coeff(i)),
  310. -numext::imag(m_eivalues.coeff(i)), numext::real(m_eivalues.coeff(i));
  311. ++i;
  312. }
  313. }
  314. return matD;
  315. }
  316. template<typename MatrixType>
  317. typename EigenSolver<MatrixType>::EigenvectorsType EigenSolver<MatrixType>::eigenvectors() const
  318. {
  319. eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
  320. eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
  321. const RealScalar precision = RealScalar(2)*NumTraits<RealScalar>::epsilon();
  322. Index n = m_eivec.cols();
  323. EigenvectorsType matV(n,n);
  324. for (Index j=0; j<n; ++j)
  325. {
  326. if (internal::isMuchSmallerThan(numext::imag(m_eivalues.coeff(j)), numext::real(m_eivalues.coeff(j)), precision) || j+1==n)
  327. {
  328. // we have a real eigen value
  329. matV.col(j) = m_eivec.col(j).template cast<ComplexScalar>();
  330. matV.col(j).normalize();
  331. }
  332. else
  333. {
  334. // we have a pair of complex eigen values
  335. for (Index i=0; i<n; ++i)
  336. {
  337. matV.coeffRef(i,j) = ComplexScalar(m_eivec.coeff(i,j), m_eivec.coeff(i,j+1));
  338. matV.coeffRef(i,j+1) = ComplexScalar(m_eivec.coeff(i,j), -m_eivec.coeff(i,j+1));
  339. }
  340. matV.col(j).normalize();
  341. matV.col(j+1).normalize();
  342. ++j;
  343. }
  344. }
  345. return matV;
  346. }
  347. template<typename MatrixType>
  348. template<typename InputType>
  349. EigenSolver<MatrixType>&
  350. EigenSolver<MatrixType>::compute(const EigenBase<InputType>& matrix, bool computeEigenvectors)
  351. {
  352. check_template_parameters();
  353. using std::sqrt;
  354. using std::abs;
  355. using numext::isfinite;
  356. eigen_assert(matrix.cols() == matrix.rows());
  357. // Reduce to real Schur form.
  358. m_realSchur.compute(matrix.derived(), computeEigenvectors);
  359. m_info = m_realSchur.info();
  360. if (m_info == Success)
  361. {
  362. m_matT = m_realSchur.matrixT();
  363. if (computeEigenvectors)
  364. m_eivec = m_realSchur.matrixU();
  365. // Compute eigenvalues from matT
  366. m_eivalues.resize(matrix.cols());
  367. Index i = 0;
  368. while (i < matrix.cols())
  369. {
  370. if (i == matrix.cols() - 1 || m_matT.coeff(i+1, i) == Scalar(0))
  371. {
  372. m_eivalues.coeffRef(i) = m_matT.coeff(i, i);
  373. if(!(isfinite)(m_eivalues.coeffRef(i)))
  374. {
  375. m_isInitialized = true;
  376. m_eigenvectorsOk = false;
  377. m_info = NumericalIssue;
  378. return *this;
  379. }
  380. ++i;
  381. }
  382. else
  383. {
  384. Scalar p = Scalar(0.5) * (m_matT.coeff(i, i) - m_matT.coeff(i+1, i+1));
  385. Scalar z;
  386. // Compute z = sqrt(abs(p * p + m_matT.coeff(i+1, i) * m_matT.coeff(i, i+1)));
  387. // without overflow
  388. {
  389. Scalar t0 = m_matT.coeff(i+1, i);
  390. Scalar t1 = m_matT.coeff(i, i+1);
  391. Scalar maxval = numext::maxi<Scalar>(abs(p),numext::maxi<Scalar>(abs(t0),abs(t1)));
  392. t0 /= maxval;
  393. t1 /= maxval;
  394. Scalar p0 = p/maxval;
  395. z = maxval * sqrt(abs(p0 * p0 + t0 * t1));
  396. }
  397. m_eivalues.coeffRef(i) = ComplexScalar(m_matT.coeff(i+1, i+1) + p, z);
  398. m_eivalues.coeffRef(i+1) = ComplexScalar(m_matT.coeff(i+1, i+1) + p, -z);
  399. if(!((isfinite)(m_eivalues.coeffRef(i)) && (isfinite)(m_eivalues.coeffRef(i+1))))
  400. {
  401. m_isInitialized = true;
  402. m_eigenvectorsOk = false;
  403. m_info = NumericalIssue;
  404. return *this;
  405. }
  406. i += 2;
  407. }
  408. }
  409. // Compute eigenvectors.
  410. if (computeEigenvectors)
  411. doComputeEigenvectors();
  412. }
  413. m_isInitialized = true;
  414. m_eigenvectorsOk = computeEigenvectors;
  415. return *this;
  416. }
  417. template<typename MatrixType>
  418. void EigenSolver<MatrixType>::doComputeEigenvectors()
  419. {
  420. using std::abs;
  421. const Index size = m_eivec.cols();
  422. const Scalar eps = NumTraits<Scalar>::epsilon();
  423. // inefficient! this is already computed in RealSchur
  424. Scalar norm(0);
  425. for (Index j = 0; j < size; ++j)
  426. {
  427. norm += m_matT.row(j).segment((std::max)(j-1,Index(0)), size-(std::max)(j-1,Index(0))).cwiseAbs().sum();
  428. }
  429. // Backsubstitute to find vectors of upper triangular form
  430. if (norm == Scalar(0))
  431. {
  432. return;
  433. }
  434. for (Index n = size-1; n >= 0; n--)
  435. {
  436. Scalar p = m_eivalues.coeff(n).real();
  437. Scalar q = m_eivalues.coeff(n).imag();
  438. // Scalar vector
  439. if (q == Scalar(0))
  440. {
  441. Scalar lastr(0), lastw(0);
  442. Index l = n;
  443. m_matT.coeffRef(n,n) = Scalar(1);
  444. for (Index i = n-1; i >= 0; i--)
  445. {
  446. Scalar w = m_matT.coeff(i,i) - p;
  447. Scalar r = m_matT.row(i).segment(l,n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
  448. if (m_eivalues.coeff(i).imag() < Scalar(0))
  449. {
  450. lastw = w;
  451. lastr = r;
  452. }
  453. else
  454. {
  455. l = i;
  456. if (m_eivalues.coeff(i).imag() == Scalar(0))
  457. {
  458. if (w != Scalar(0))
  459. m_matT.coeffRef(i,n) = -r / w;
  460. else
  461. m_matT.coeffRef(i,n) = -r / (eps * norm);
  462. }
  463. else // Solve real equations
  464. {
  465. Scalar x = m_matT.coeff(i,i+1);
  466. Scalar y = m_matT.coeff(i+1,i);
  467. Scalar denom = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) + m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag();
  468. Scalar t = (x * lastr - lastw * r) / denom;
  469. m_matT.coeffRef(i,n) = t;
  470. if (abs(x) > abs(lastw))
  471. m_matT.coeffRef(i+1,n) = (-r - w * t) / x;
  472. else
  473. m_matT.coeffRef(i+1,n) = (-lastr - y * t) / lastw;
  474. }
  475. // Overflow control
  476. Scalar t = abs(m_matT.coeff(i,n));
  477. if ((eps * t) * t > Scalar(1))
  478. m_matT.col(n).tail(size-i) /= t;
  479. }
  480. }
  481. }
  482. else if (q < Scalar(0) && n > 0) // Complex vector
  483. {
  484. Scalar lastra(0), lastsa(0), lastw(0);
  485. Index l = n-1;
  486. // Last vector component imaginary so matrix is triangular
  487. if (abs(m_matT.coeff(n,n-1)) > abs(m_matT.coeff(n-1,n)))
  488. {
  489. m_matT.coeffRef(n-1,n-1) = q / m_matT.coeff(n,n-1);
  490. m_matT.coeffRef(n-1,n) = -(m_matT.coeff(n,n) - p) / m_matT.coeff(n,n-1);
  491. }
  492. else
  493. {
  494. ComplexScalar cc = ComplexScalar(Scalar(0),-m_matT.coeff(n-1,n)) / ComplexScalar(m_matT.coeff(n-1,n-1)-p,q);
  495. m_matT.coeffRef(n-1,n-1) = numext::real(cc);
  496. m_matT.coeffRef(n-1,n) = numext::imag(cc);
  497. }
  498. m_matT.coeffRef(n,n-1) = Scalar(0);
  499. m_matT.coeffRef(n,n) = Scalar(1);
  500. for (Index i = n-2; i >= 0; i--)
  501. {
  502. Scalar ra = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n-1).segment(l, n-l+1));
  503. Scalar sa = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
  504. Scalar w = m_matT.coeff(i,i) - p;
  505. if (m_eivalues.coeff(i).imag() < Scalar(0))
  506. {
  507. lastw = w;
  508. lastra = ra;
  509. lastsa = sa;
  510. }
  511. else
  512. {
  513. l = i;
  514. if (m_eivalues.coeff(i).imag() == RealScalar(0))
  515. {
  516. ComplexScalar cc = ComplexScalar(-ra,-sa) / ComplexScalar(w,q);
  517. m_matT.coeffRef(i,n-1) = numext::real(cc);
  518. m_matT.coeffRef(i,n) = numext::imag(cc);
  519. }
  520. else
  521. {
  522. // Solve complex equations
  523. Scalar x = m_matT.coeff(i,i+1);
  524. Scalar y = m_matT.coeff(i+1,i);
  525. Scalar vr = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) + m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag() - q * q;
  526. Scalar vi = (m_eivalues.coeff(i).real() - p) * Scalar(2) * q;
  527. if ((vr == Scalar(0)) && (vi == Scalar(0)))
  528. vr = eps * norm * (abs(w) + abs(q) + abs(x) + abs(y) + abs(lastw));
  529. ComplexScalar cc = ComplexScalar(x*lastra-lastw*ra+q*sa,x*lastsa-lastw*sa-q*ra) / ComplexScalar(vr,vi);
  530. m_matT.coeffRef(i,n-1) = numext::real(cc);
  531. m_matT.coeffRef(i,n) = numext::imag(cc);
  532. if (abs(x) > (abs(lastw) + abs(q)))
  533. {
  534. m_matT.coeffRef(i+1,n-1) = (-ra - w * m_matT.coeff(i,n-1) + q * m_matT.coeff(i,n)) / x;
  535. m_matT.coeffRef(i+1,n) = (-sa - w * m_matT.coeff(i,n) - q * m_matT.coeff(i,n-1)) / x;
  536. }
  537. else
  538. {
  539. cc = ComplexScalar(-lastra-y*m_matT.coeff(i,n-1),-lastsa-y*m_matT.coeff(i,n)) / ComplexScalar(lastw,q);
  540. m_matT.coeffRef(i+1,n-1) = numext::real(cc);
  541. m_matT.coeffRef(i+1,n) = numext::imag(cc);
  542. }
  543. }
  544. // Overflow control
  545. Scalar t = numext::maxi<Scalar>(abs(m_matT.coeff(i,n-1)),abs(m_matT.coeff(i,n)));
  546. if ((eps * t) * t > Scalar(1))
  547. m_matT.block(i, n-1, size-i, 2) /= t;
  548. }
  549. }
  550. // We handled a pair of complex conjugate eigenvalues, so need to skip them both
  551. n--;
  552. }
  553. else
  554. {
  555. eigen_assert(0 && "Internal bug in EigenSolver (INF or NaN has not been detected)"); // this should not happen
  556. }
  557. }
  558. // Back transformation to get eigenvectors of original matrix
  559. for (Index j = size-1; j >= 0; j--)
  560. {
  561. m_tmp.noalias() = m_eivec.leftCols(j+1) * m_matT.col(j).segment(0, j+1);
  562. m_eivec.col(j) = m_tmp;
  563. }
  564. }
  565. } // end namespace Eigen
  566. #endif // EIGEN_EIGENSOLVER_H