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