GeneralizedEigenSolver.h 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2012-2016 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. // Copyright (C) 2010,2012 Jitse Niesen <jitse@maths.leeds.ac.uk>
  6. // Copyright (C) 2016 Tobias Wood <tobias@spinicist.org.uk>
  7. //
  8. // This Source Code Form is subject to the terms of the Mozilla
  9. // Public License v. 2.0. If a copy of the MPL was not distributed
  10. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  11. #ifndef EIGEN_GENERALIZEDEIGENSOLVER_H
  12. #define EIGEN_GENERALIZEDEIGENSOLVER_H
  13. #include "./RealQZ.h"
  14. namespace Eigen {
  15. /** \eigenvalues_module \ingroup Eigenvalues_Module
  16. *
  17. *
  18. * \class GeneralizedEigenSolver
  19. *
  20. * \brief Computes the generalized eigenvalues and eigenvectors of a pair of general matrices
  21. *
  22. * \tparam _MatrixType the type of the matrices of which we are computing the
  23. * eigen-decomposition; this is expected to be an instantiation of the Matrix
  24. * class template. Currently, only real matrices are supported.
  25. *
  26. * The generalized eigenvalues and eigenvectors of a matrix pair \f$ A \f$ and \f$ B \f$ are scalars
  27. * \f$ \lambda \f$ and vectors \f$ v \f$ such that \f$ Av = \lambda Bv \f$. If
  28. * \f$ D \f$ is a diagonal matrix with the eigenvalues on the diagonal, and
  29. * \f$ V \f$ is a matrix with the eigenvectors as its columns, then \f$ A V =
  30. * B V D \f$. The matrix \f$ V \f$ is almost always invertible, in which case we
  31. * have \f$ A = B V D V^{-1} \f$. This is called the generalized eigen-decomposition.
  32. *
  33. * The generalized eigenvalues and eigenvectors of a matrix pair may be complex, even when the
  34. * matrices are real. Moreover, the generalized eigenvalue might be infinite if the matrix B is
  35. * singular. To workaround this difficulty, the eigenvalues are provided as a pair of complex \f$ \alpha \f$
  36. * and real \f$ \beta \f$ such that: \f$ \lambda_i = \alpha_i / \beta_i \f$. If \f$ \beta_i \f$ is (nearly) zero,
  37. * then one can consider the well defined left eigenvalue \f$ \mu = \beta_i / \alpha_i\f$ such that:
  38. * \f$ \mu_i A v_i = B v_i \f$, or even \f$ \mu_i u_i^T A = u_i^T B \f$ where \f$ u_i \f$ is
  39. * called the left eigenvector.
  40. *
  41. * Call the function compute() to compute the generalized eigenvalues and eigenvectors of
  42. * a given matrix pair. Alternatively, you can use the
  43. * GeneralizedEigenSolver(const MatrixType&, const MatrixType&, bool) constructor which computes the
  44. * eigenvalues and eigenvectors at construction time. Once the eigenvalue and
  45. * eigenvectors are computed, they can be retrieved with the eigenvalues() and
  46. * eigenvectors() functions.
  47. *
  48. * Here is an usage example of this class:
  49. * Example: \include GeneralizedEigenSolver.cpp
  50. * Output: \verbinclude GeneralizedEigenSolver.out
  51. *
  52. * \sa MatrixBase::eigenvalues(), class ComplexEigenSolver, class SelfAdjointEigenSolver
  53. */
  54. template<typename _MatrixType> class GeneralizedEigenSolver
  55. {
  56. public:
  57. /** \brief Synonym for the template parameter \p _MatrixType. */
  58. typedef _MatrixType MatrixType;
  59. enum {
  60. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  61. ColsAtCompileTime = MatrixType::ColsAtCompileTime,
  62. Options = MatrixType::Options,
  63. MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  64. MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
  65. };
  66. /** \brief Scalar type for matrices of type #MatrixType. */
  67. typedef typename MatrixType::Scalar Scalar;
  68. typedef typename NumTraits<Scalar>::Real RealScalar;
  69. typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
  70. /** \brief Complex scalar type for #MatrixType.
  71. *
  72. * This is \c std::complex<Scalar> if #Scalar is real (e.g.,
  73. * \c float or \c double) and just \c Scalar if #Scalar is
  74. * complex.
  75. */
  76. typedef std::complex<RealScalar> ComplexScalar;
  77. /** \brief Type for vector of real scalar values eigenvalues as returned by betas().
  78. *
  79. * This is a column vector with entries of type #Scalar.
  80. * The length of the vector is the size of #MatrixType.
  81. */
  82. typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> VectorType;
  83. /** \brief Type for vector of complex scalar values eigenvalues as returned by alphas().
  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> ComplexVectorType;
  89. /** \brief Expression type for the eigenvalues as returned by eigenvalues().
  90. */
  91. typedef CwiseBinaryOp<internal::scalar_quotient_op<ComplexScalar,Scalar>,ComplexVectorType,VectorType> EigenvalueType;
  92. /** \brief Type for matrix of eigenvectors as returned by eigenvectors().
  93. *
  94. * This is a square matrix with entries of type #ComplexScalar.
  95. * The size is the same as the size of #MatrixType.
  96. */
  97. typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> EigenvectorsType;
  98. /** \brief Default constructor.
  99. *
  100. * The default constructor is useful in cases in which the user intends to
  101. * perform decompositions via EigenSolver::compute(const MatrixType&, bool).
  102. *
  103. * \sa compute() for an example.
  104. */
  105. GeneralizedEigenSolver()
  106. : m_eivec(),
  107. m_alphas(),
  108. m_betas(),
  109. m_valuesOkay(false),
  110. m_vectorsOkay(false),
  111. m_realQZ()
  112. {}
  113. /** \brief Default constructor with memory preallocation
  114. *
  115. * Like the default constructor but with preallocation of the internal data
  116. * according to the specified problem \a size.
  117. * \sa GeneralizedEigenSolver()
  118. */
  119. explicit GeneralizedEigenSolver(Index size)
  120. : m_eivec(size, size),
  121. m_alphas(size),
  122. m_betas(size),
  123. m_valuesOkay(false),
  124. m_vectorsOkay(false),
  125. m_realQZ(size),
  126. m_tmp(size)
  127. {}
  128. /** \brief Constructor; computes the generalized eigendecomposition of given matrix pair.
  129. *
  130. * \param[in] A Square matrix whose eigendecomposition is to be computed.
  131. * \param[in] B Square matrix whose eigendecomposition is to be computed.
  132. * \param[in] computeEigenvectors If true, both the eigenvectors and the
  133. * eigenvalues are computed; if false, only the eigenvalues are computed.
  134. *
  135. * This constructor calls compute() to compute the generalized eigenvalues
  136. * and eigenvectors.
  137. *
  138. * \sa compute()
  139. */
  140. GeneralizedEigenSolver(const MatrixType& A, const MatrixType& B, bool computeEigenvectors = true)
  141. : m_eivec(A.rows(), A.cols()),
  142. m_alphas(A.cols()),
  143. m_betas(A.cols()),
  144. m_valuesOkay(false),
  145. m_vectorsOkay(false),
  146. m_realQZ(A.cols()),
  147. m_tmp(A.cols())
  148. {
  149. compute(A, B, computeEigenvectors);
  150. }
  151. /* \brief Returns the computed generalized eigenvectors.
  152. *
  153. * \returns %Matrix whose columns are the (possibly complex) right eigenvectors.
  154. * i.e. the eigenvectors that solve (A - l*B)x = 0. The ordering matches the eigenvalues.
  155. *
  156. * \pre Either the constructor
  157. * GeneralizedEigenSolver(const MatrixType&,const MatrixType&, bool) or the member function
  158. * compute(const MatrixType&, const MatrixType& bool) has been called before, and
  159. * \p computeEigenvectors was set to true (the default).
  160. *
  161. * \sa eigenvalues()
  162. */
  163. EigenvectorsType eigenvectors() const {
  164. eigen_assert(m_vectorsOkay && "Eigenvectors for GeneralizedEigenSolver were not calculated.");
  165. return m_eivec;
  166. }
  167. /** \brief Returns an expression of the computed generalized eigenvalues.
  168. *
  169. * \returns An expression of the column vector containing the eigenvalues.
  170. *
  171. * It is a shortcut for \code this->alphas().cwiseQuotient(this->betas()); \endcode
  172. * Not that betas might contain zeros. It is therefore not recommended to use this function,
  173. * but rather directly deal with the alphas and betas vectors.
  174. *
  175. * \pre Either the constructor
  176. * GeneralizedEigenSolver(const MatrixType&,const MatrixType&,bool) or the member function
  177. * compute(const MatrixType&,const MatrixType&,bool) has been called before.
  178. *
  179. * The eigenvalues are repeated according to their algebraic multiplicity,
  180. * so there are as many eigenvalues as rows in the matrix. The eigenvalues
  181. * are not sorted in any particular order.
  182. *
  183. * \sa alphas(), betas(), eigenvectors()
  184. */
  185. EigenvalueType eigenvalues() const
  186. {
  187. eigen_assert(m_valuesOkay && "GeneralizedEigenSolver is not initialized.");
  188. return EigenvalueType(m_alphas,m_betas);
  189. }
  190. /** \returns A const reference to the vectors containing the alpha values
  191. *
  192. * This vector permits to reconstruct the j-th eigenvalues as alphas(i)/betas(j).
  193. *
  194. * \sa betas(), eigenvalues() */
  195. ComplexVectorType alphas() const
  196. {
  197. eigen_assert(m_valuesOkay && "GeneralizedEigenSolver is not initialized.");
  198. return m_alphas;
  199. }
  200. /** \returns A const reference to the vectors containing the beta values
  201. *
  202. * This vector permits to reconstruct the j-th eigenvalues as alphas(i)/betas(j).
  203. *
  204. * \sa alphas(), eigenvalues() */
  205. VectorType betas() const
  206. {
  207. eigen_assert(m_valuesOkay && "GeneralizedEigenSolver is not initialized.");
  208. return m_betas;
  209. }
  210. /** \brief Computes generalized eigendecomposition of given matrix.
  211. *
  212. * \param[in] A Square matrix whose eigendecomposition is to be computed.
  213. * \param[in] B Square matrix whose eigendecomposition is to be computed.
  214. * \param[in] computeEigenvectors If true, both the eigenvectors and the
  215. * eigenvalues are computed; if false, only the eigenvalues are
  216. * computed.
  217. * \returns Reference to \c *this
  218. *
  219. * This function computes the eigenvalues of the real matrix \p matrix.
  220. * The eigenvalues() function can be used to retrieve them. If
  221. * \p computeEigenvectors is true, then the eigenvectors are also computed
  222. * and can be retrieved by calling eigenvectors().
  223. *
  224. * The matrix is first reduced to real generalized Schur form using the RealQZ
  225. * class. The generalized Schur decomposition is then used to compute the eigenvalues
  226. * and eigenvectors.
  227. *
  228. * The cost of the computation is dominated by the cost of the
  229. * generalized Schur decomposition.
  230. *
  231. * This method reuses of the allocated data in the GeneralizedEigenSolver object.
  232. */
  233. GeneralizedEigenSolver& compute(const MatrixType& A, const MatrixType& B, bool computeEigenvectors = true);
  234. ComputationInfo info() const
  235. {
  236. eigen_assert(m_valuesOkay && "EigenSolver is not initialized.");
  237. return m_realQZ.info();
  238. }
  239. /** Sets the maximal number of iterations allowed.
  240. */
  241. GeneralizedEigenSolver& setMaxIterations(Index maxIters)
  242. {
  243. m_realQZ.setMaxIterations(maxIters);
  244. return *this;
  245. }
  246. protected:
  247. static void check_template_parameters()
  248. {
  249. EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
  250. EIGEN_STATIC_ASSERT(!NumTraits<Scalar>::IsComplex, NUMERIC_TYPE_MUST_BE_REAL);
  251. }
  252. EigenvectorsType m_eivec;
  253. ComplexVectorType m_alphas;
  254. VectorType m_betas;
  255. bool m_valuesOkay, m_vectorsOkay;
  256. RealQZ<MatrixType> m_realQZ;
  257. ComplexVectorType m_tmp;
  258. };
  259. template<typename MatrixType>
  260. GeneralizedEigenSolver<MatrixType>&
  261. GeneralizedEigenSolver<MatrixType>::compute(const MatrixType& A, const MatrixType& B, bool computeEigenvectors)
  262. {
  263. check_template_parameters();
  264. using std::sqrt;
  265. using std::abs;
  266. eigen_assert(A.cols() == A.rows() && B.cols() == A.rows() && B.cols() == B.rows());
  267. Index size = A.cols();
  268. m_valuesOkay = false;
  269. m_vectorsOkay = false;
  270. // Reduce to generalized real Schur form:
  271. // A = Q S Z and B = Q T Z
  272. m_realQZ.compute(A, B, computeEigenvectors);
  273. if (m_realQZ.info() == Success)
  274. {
  275. // Resize storage
  276. m_alphas.resize(size);
  277. m_betas.resize(size);
  278. if (computeEigenvectors)
  279. {
  280. m_eivec.resize(size,size);
  281. m_tmp.resize(size);
  282. }
  283. // Aliases:
  284. Map<VectorType> v(reinterpret_cast<Scalar*>(m_tmp.data()), size);
  285. ComplexVectorType &cv = m_tmp;
  286. const MatrixType &mS = m_realQZ.matrixS();
  287. const MatrixType &mT = m_realQZ.matrixT();
  288. Index i = 0;
  289. while (i < size)
  290. {
  291. if (i == size - 1 || mS.coeff(i+1, i) == Scalar(0))
  292. {
  293. // Real eigenvalue
  294. m_alphas.coeffRef(i) = mS.diagonal().coeff(i);
  295. m_betas.coeffRef(i) = mT.diagonal().coeff(i);
  296. if (computeEigenvectors)
  297. {
  298. v.setConstant(Scalar(0.0));
  299. v.coeffRef(i) = Scalar(1.0);
  300. // For singular eigenvalues do nothing more
  301. if(abs(m_betas.coeffRef(i)) >= (std::numeric_limits<RealScalar>::min)())
  302. {
  303. // Non-singular eigenvalue
  304. const Scalar alpha = real(m_alphas.coeffRef(i));
  305. const Scalar beta = m_betas.coeffRef(i);
  306. for (Index j = i-1; j >= 0; j--)
  307. {
  308. const Index st = j+1;
  309. const Index sz = i-j;
  310. if (j > 0 && mS.coeff(j, j-1) != Scalar(0))
  311. {
  312. // 2x2 block
  313. Matrix<Scalar, 2, 1> rhs = (alpha*mT.template block<2,Dynamic>(j-1,st,2,sz) - beta*mS.template block<2,Dynamic>(j-1,st,2,sz)) .lazyProduct( v.segment(st,sz) );
  314. Matrix<Scalar, 2, 2> lhs = beta * mS.template block<2,2>(j-1,j-1) - alpha * mT.template block<2,2>(j-1,j-1);
  315. v.template segment<2>(j-1) = lhs.partialPivLu().solve(rhs);
  316. j--;
  317. }
  318. else
  319. {
  320. v.coeffRef(j) = -v.segment(st,sz).transpose().cwiseProduct(beta*mS.block(j,st,1,sz) - alpha*mT.block(j,st,1,sz)).sum() / (beta*mS.coeffRef(j,j) - alpha*mT.coeffRef(j,j));
  321. }
  322. }
  323. }
  324. m_eivec.col(i).real().noalias() = m_realQZ.matrixZ().transpose() * v;
  325. m_eivec.col(i).real().normalize();
  326. m_eivec.col(i).imag().setConstant(0);
  327. }
  328. ++i;
  329. }
  330. else
  331. {
  332. // We need to extract the generalized eigenvalues of the pair of a general 2x2 block S and a positive diagonal 2x2 block T
  333. // Then taking beta=T_00*T_11, we can avoid any division, and alpha is the eigenvalues of A = (U^-1 * S * U) * diag(T_11,T_00):
  334. // T = [a 0]
  335. // [0 b]
  336. RealScalar a = mT.diagonal().coeff(i),
  337. b = mT.diagonal().coeff(i+1);
  338. const RealScalar beta = m_betas.coeffRef(i) = m_betas.coeffRef(i+1) = a*b;
  339. // ^^ NOTE: using diagonal()(i) instead of coeff(i,i) workarounds a MSVC bug.
  340. Matrix<RealScalar,2,2> S2 = mS.template block<2,2>(i,i) * Matrix<Scalar,2,1>(b,a).asDiagonal();
  341. Scalar p = Scalar(0.5) * (S2.coeff(0,0) - S2.coeff(1,1));
  342. Scalar z = sqrt(abs(p * p + S2.coeff(1,0) * S2.coeff(0,1)));
  343. const ComplexScalar alpha = ComplexScalar(S2.coeff(1,1) + p, (beta > 0) ? z : -z);
  344. m_alphas.coeffRef(i) = conj(alpha);
  345. m_alphas.coeffRef(i+1) = alpha;
  346. if (computeEigenvectors) {
  347. // Compute eigenvector in position (i+1) and then position (i) is just the conjugate
  348. cv.setZero();
  349. cv.coeffRef(i+1) = Scalar(1.0);
  350. // here, the "static_cast" workaound expression template issues.
  351. cv.coeffRef(i) = -(static_cast<Scalar>(beta*mS.coeffRef(i,i+1)) - alpha*mT.coeffRef(i,i+1))
  352. / (static_cast<Scalar>(beta*mS.coeffRef(i,i)) - alpha*mT.coeffRef(i,i));
  353. for (Index j = i-1; j >= 0; j--)
  354. {
  355. const Index st = j+1;
  356. const Index sz = i+1-j;
  357. if (j > 0 && mS.coeff(j, j-1) != Scalar(0))
  358. {
  359. // 2x2 block
  360. Matrix<ComplexScalar, 2, 1> rhs = (alpha*mT.template block<2,Dynamic>(j-1,st,2,sz) - beta*mS.template block<2,Dynamic>(j-1,st,2,sz)) .lazyProduct( cv.segment(st,sz) );
  361. Matrix<ComplexScalar, 2, 2> lhs = beta * mS.template block<2,2>(j-1,j-1) - alpha * mT.template block<2,2>(j-1,j-1);
  362. cv.template segment<2>(j-1) = lhs.partialPivLu().solve(rhs);
  363. j--;
  364. } else {
  365. cv.coeffRef(j) = cv.segment(st,sz).transpose().cwiseProduct(beta*mS.block(j,st,1,sz) - alpha*mT.block(j,st,1,sz)).sum()
  366. / (alpha*mT.coeffRef(j,j) - static_cast<Scalar>(beta*mS.coeffRef(j,j)));
  367. }
  368. }
  369. m_eivec.col(i+1).noalias() = (m_realQZ.matrixZ().transpose() * cv);
  370. m_eivec.col(i+1).normalize();
  371. m_eivec.col(i) = m_eivec.col(i+1).conjugate();
  372. }
  373. i += 2;
  374. }
  375. }
  376. m_valuesOkay = true;
  377. m_vectorsOkay = computeEigenvectors;
  378. }
  379. return *this;
  380. }
  381. } // end namespace Eigen
  382. #endif // EIGEN_GENERALIZEDEIGENSOLVER_H