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