GeneralizedSelfAdjointEigenSolver.h 9.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226
  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_GENERALIZEDSELFADJOINTEIGENSOLVER_H
  11. #define EIGEN_GENERALIZEDSELFADJOINTEIGENSOLVER_H
  12. #include "./Tridiagonalization.h"
  13. namespace Eigen {
  14. /** \eigenvalues_module \ingroup Eigenvalues_Module
  15. *
  16. *
  17. * \class GeneralizedSelfAdjointEigenSolver
  18. *
  19. * \brief Computes eigenvalues and eigenvectors of the generalized selfadjoint eigen problem
  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.
  24. *
  25. * This class solves the generalized eigenvalue problem
  26. * \f$ Av = \lambda Bv \f$. In this case, the matrix \f$ A \f$ should be
  27. * selfadjoint and the matrix \f$ B \f$ should be positive definite.
  28. *
  29. * Only the \b lower \b triangular \b part of the input matrix is referenced.
  30. *
  31. * Call the function compute() to compute the eigenvalues and eigenvectors of
  32. * a given matrix. Alternatively, you can use the
  33. * GeneralizedSelfAdjointEigenSolver(const MatrixType&, const MatrixType&, int)
  34. * constructor which computes the eigenvalues and eigenvectors at construction time.
  35. * Once the eigenvalue and eigenvectors are computed, they can be retrieved with the eigenvalues()
  36. * and eigenvectors() functions.
  37. *
  38. * The documentation for GeneralizedSelfAdjointEigenSolver(const MatrixType&, const MatrixType&, int)
  39. * contains an example of the typical use of this class.
  40. *
  41. * \sa class SelfAdjointEigenSolver, class EigenSolver, class ComplexEigenSolver
  42. */
  43. template<typename _MatrixType>
  44. class GeneralizedSelfAdjointEigenSolver : public SelfAdjointEigenSolver<_MatrixType>
  45. {
  46. typedef SelfAdjointEigenSolver<_MatrixType> Base;
  47. public:
  48. typedef _MatrixType MatrixType;
  49. /** \brief Default constructor for fixed-size matrices.
  50. *
  51. * The default constructor is useful in cases in which the user intends to
  52. * perform decompositions via compute(). This constructor
  53. * can only be used if \p _MatrixType is a fixed-size matrix; use
  54. * GeneralizedSelfAdjointEigenSolver(Index) for dynamic-size matrices.
  55. */
  56. GeneralizedSelfAdjointEigenSolver() : Base() {}
  57. /** \brief Constructor, pre-allocates memory for dynamic-size matrices.
  58. *
  59. * \param [in] size Positive integer, size of the matrix whose
  60. * eigenvalues and eigenvectors will be computed.
  61. *
  62. * This constructor is useful for dynamic-size matrices, when the user
  63. * intends to perform decompositions via compute(). The \p size
  64. * parameter is only used as a hint. It is not an error to give a wrong
  65. * \p size, but it may impair performance.
  66. *
  67. * \sa compute() for an example
  68. */
  69. explicit GeneralizedSelfAdjointEigenSolver(Index size)
  70. : Base(size)
  71. {}
  72. /** \brief Constructor; computes generalized eigendecomposition of given matrix pencil.
  73. *
  74. * \param[in] matA Selfadjoint matrix in matrix pencil.
  75. * Only the lower triangular part of the matrix is referenced.
  76. * \param[in] matB Positive-definite matrix in matrix pencil.
  77. * Only the lower triangular part of the matrix is referenced.
  78. * \param[in] options A or-ed set of flags {#ComputeEigenvectors,#EigenvaluesOnly} | {#Ax_lBx,#ABx_lx,#BAx_lx}.
  79. * Default is #ComputeEigenvectors|#Ax_lBx.
  80. *
  81. * This constructor calls compute(const MatrixType&, const MatrixType&, int)
  82. * to compute the eigenvalues and (if requested) the eigenvectors of the
  83. * generalized eigenproblem \f$ Ax = \lambda B x \f$ with \a matA the
  84. * selfadjoint matrix \f$ A \f$ and \a matB the positive definite matrix
  85. * \f$ B \f$. Each eigenvector \f$ x \f$ satisfies the property
  86. * \f$ x^* B x = 1 \f$. The eigenvectors are computed if
  87. * \a options contains ComputeEigenvectors.
  88. *
  89. * In addition, the two following variants can be solved via \p options:
  90. * - \c ABx_lx: \f$ ABx = \lambda x \f$
  91. * - \c BAx_lx: \f$ BAx = \lambda x \f$
  92. *
  93. * Example: \include SelfAdjointEigenSolver_SelfAdjointEigenSolver_MatrixType2.cpp
  94. * Output: \verbinclude SelfAdjointEigenSolver_SelfAdjointEigenSolver_MatrixType2.out
  95. *
  96. * \sa compute(const MatrixType&, const MatrixType&, int)
  97. */
  98. GeneralizedSelfAdjointEigenSolver(const MatrixType& matA, const MatrixType& matB,
  99. int options = ComputeEigenvectors|Ax_lBx)
  100. : Base(matA.cols())
  101. {
  102. compute(matA, matB, options);
  103. }
  104. /** \brief Computes generalized eigendecomposition of given matrix pencil.
  105. *
  106. * \param[in] matA Selfadjoint matrix in matrix pencil.
  107. * Only the lower triangular part of the matrix is referenced.
  108. * \param[in] matB Positive-definite matrix in matrix pencil.
  109. * Only the lower triangular part of the matrix is referenced.
  110. * \param[in] options A or-ed set of flags {#ComputeEigenvectors,#EigenvaluesOnly} | {#Ax_lBx,#ABx_lx,#BAx_lx}.
  111. * Default is #ComputeEigenvectors|#Ax_lBx.
  112. *
  113. * \returns Reference to \c *this
  114. *
  115. * Accoring to \p options, this function computes eigenvalues and (if requested)
  116. * the eigenvectors of one of the following three generalized eigenproblems:
  117. * - \c Ax_lBx: \f$ Ax = \lambda B x \f$
  118. * - \c ABx_lx: \f$ ABx = \lambda x \f$
  119. * - \c BAx_lx: \f$ BAx = \lambda x \f$
  120. * with \a matA the selfadjoint matrix \f$ A \f$ and \a matB the positive definite
  121. * matrix \f$ B \f$.
  122. * In addition, each eigenvector \f$ x \f$ satisfies the property \f$ x^* B x = 1 \f$.
  123. *
  124. * The eigenvalues() function can be used to retrieve
  125. * the eigenvalues. If \p options contains ComputeEigenvectors, then the
  126. * eigenvectors are also computed and can be retrieved by calling
  127. * eigenvectors().
  128. *
  129. * The implementation uses LLT to compute the Cholesky decomposition
  130. * \f$ B = LL^* \f$ and computes the classical eigendecomposition
  131. * of the selfadjoint matrix \f$ L^{-1} A (L^*)^{-1} \f$ if \p options contains Ax_lBx
  132. * and of \f$ L^{*} A L \f$ otherwise. This solves the
  133. * generalized eigenproblem, because any solution of the generalized
  134. * eigenproblem \f$ Ax = \lambda B x \f$ corresponds to a solution
  135. * \f$ L^{-1} A (L^*)^{-1} (L^* x) = \lambda (L^* x) \f$ of the
  136. * eigenproblem for \f$ L^{-1} A (L^*)^{-1} \f$. Similar statements
  137. * can be made for the two other variants.
  138. *
  139. * Example: \include SelfAdjointEigenSolver_compute_MatrixType2.cpp
  140. * Output: \verbinclude SelfAdjointEigenSolver_compute_MatrixType2.out
  141. *
  142. * \sa GeneralizedSelfAdjointEigenSolver(const MatrixType&, const MatrixType&, int)
  143. */
  144. GeneralizedSelfAdjointEigenSolver& compute(const MatrixType& matA, const MatrixType& matB,
  145. int options = ComputeEigenvectors|Ax_lBx);
  146. protected:
  147. };
  148. template<typename MatrixType>
  149. GeneralizedSelfAdjointEigenSolver<MatrixType>& GeneralizedSelfAdjointEigenSolver<MatrixType>::
  150. compute(const MatrixType& matA, const MatrixType& matB, int options)
  151. {
  152. eigen_assert(matA.cols()==matA.rows() && matB.rows()==matA.rows() && matB.cols()==matB.rows());
  153. eigen_assert((options&~(EigVecMask|GenEigMask))==0
  154. && (options&EigVecMask)!=EigVecMask
  155. && ((options&GenEigMask)==0 || (options&GenEigMask)==Ax_lBx
  156. || (options&GenEigMask)==ABx_lx || (options&GenEigMask)==BAx_lx)
  157. && "invalid option parameter");
  158. bool computeEigVecs = ((options&EigVecMask)==0) || ((options&EigVecMask)==ComputeEigenvectors);
  159. // Compute the cholesky decomposition of matB = L L' = U'U
  160. LLT<MatrixType> cholB(matB);
  161. int type = (options&GenEigMask);
  162. if(type==0)
  163. type = Ax_lBx;
  164. if(type==Ax_lBx)
  165. {
  166. // compute C = inv(L) A inv(L')
  167. MatrixType matC = matA.template selfadjointView<Lower>();
  168. cholB.matrixL().template solveInPlace<OnTheLeft>(matC);
  169. cholB.matrixU().template solveInPlace<OnTheRight>(matC);
  170. Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly );
  171. // transform back the eigen vectors: evecs = inv(U) * evecs
  172. if(computeEigVecs)
  173. cholB.matrixU().solveInPlace(Base::m_eivec);
  174. }
  175. else if(type==ABx_lx)
  176. {
  177. // compute C = L' A L
  178. MatrixType matC = matA.template selfadjointView<Lower>();
  179. matC = matC * cholB.matrixL();
  180. matC = cholB.matrixU() * matC;
  181. Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly);
  182. // transform back the eigen vectors: evecs = inv(U) * evecs
  183. if(computeEigVecs)
  184. cholB.matrixU().solveInPlace(Base::m_eivec);
  185. }
  186. else if(type==BAx_lx)
  187. {
  188. // compute C = L' A L
  189. MatrixType matC = matA.template selfadjointView<Lower>();
  190. matC = matC * cholB.matrixL();
  191. matC = cholB.matrixU() * matC;
  192. Base::compute(matC, computeEigVecs ? ComputeEigenvectors : EigenvaluesOnly);
  193. // transform back the eigen vectors: evecs = L * evecs
  194. if(computeEigVecs)
  195. Base::m_eivec = cholB.matrixL() * Base::m_eivec;
  196. }
  197. return *this;
  198. }
  199. } // end namespace Eigen
  200. #endif // EIGEN_GENERALIZEDSELFADJOINTEIGENSOLVER_H