Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/eigen3/Eigen/src/Eigenvalues/EigenSolver.h was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
0005 // Copyright (C) 2010,2012 Jitse Niesen <jitse@maths.leeds.ac.uk>
0006 //
0007 // This Source Code Form is subject to the terms of the Mozilla
0008 // Public License v. 2.0. If a copy of the MPL was not distributed
0009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0010 
0011 #ifndef EIGEN_EIGENSOLVER_H
0012 #define EIGEN_EIGENSOLVER_H
0013 
0014 #include "./RealSchur.h"
0015 
0016 namespace Eigen { 
0017 
0018 /** \eigenvalues_module \ingroup Eigenvalues_Module
0019   *
0020   *
0021   * \class EigenSolver
0022   *
0023   * \brief Computes eigenvalues and eigenvectors of general matrices
0024   *
0025   * \tparam _MatrixType the type of the matrix of which we are computing the
0026   * eigendecomposition; this is expected to be an instantiation of the Matrix
0027   * class template. Currently, only real matrices are supported.
0028   *
0029   * The eigenvalues and eigenvectors of a matrix \f$ A \f$ are scalars
0030   * \f$ \lambda \f$ and vectors \f$ v \f$ such that \f$ Av = \lambda v \f$.  If
0031   * \f$ D \f$ is a diagonal matrix with the eigenvalues on the diagonal, and
0032   * \f$ V \f$ is a matrix with the eigenvectors as its columns, then \f$ A V =
0033   * V D \f$. The matrix \f$ V \f$ is almost always invertible, in which case we
0034   * have \f$ A = V D V^{-1} \f$. This is called the eigendecomposition.
0035   *
0036   * The eigenvalues and eigenvectors of a matrix may be complex, even when the
0037   * matrix is real. However, we can choose real matrices \f$ V \f$ and \f$ D
0038   * \f$ satisfying \f$ A V = V D \f$, just like the eigendecomposition, if the
0039   * matrix \f$ D \f$ is not required to be diagonal, but if it is allowed to
0040   * have blocks of the form
0041   * \f[ \begin{bmatrix} u & v \\ -v & u \end{bmatrix} \f]
0042   * (where \f$ u \f$ and \f$ v \f$ are real numbers) on the diagonal.  These
0043   * blocks correspond to complex eigenvalue pairs \f$ u \pm iv \f$. We call
0044   * this variant of the eigendecomposition the pseudo-eigendecomposition.
0045   *
0046   * Call the function compute() to compute the eigenvalues and eigenvectors of
0047   * a given matrix. Alternatively, you can use the 
0048   * EigenSolver(const MatrixType&, bool) constructor which computes the
0049   * eigenvalues and eigenvectors at construction time. Once the eigenvalue and
0050   * eigenvectors are computed, they can be retrieved with the eigenvalues() and
0051   * eigenvectors() functions. The pseudoEigenvalueMatrix() and
0052   * pseudoEigenvectors() methods allow the construction of the
0053   * pseudo-eigendecomposition.
0054   *
0055   * The documentation for EigenSolver(const MatrixType&, bool) contains an
0056   * example of the typical use of this class.
0057   *
0058   * \note The implementation is adapted from
0059   * <a href="http://math.nist.gov/javanumerics/jama/">JAMA</a> (public domain).
0060   * Their code is based on EISPACK.
0061   *
0062   * \sa MatrixBase::eigenvalues(), class ComplexEigenSolver, class SelfAdjointEigenSolver
0063   */
0064 template<typename _MatrixType> class EigenSolver
0065 {
0066   public:
0067 
0068     /** \brief Synonym for the template parameter \p _MatrixType. */
0069     typedef _MatrixType MatrixType;
0070 
0071     enum {
0072       RowsAtCompileTime = MatrixType::RowsAtCompileTime,
0073       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
0074       Options = MatrixType::Options,
0075       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
0076       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
0077     };
0078 
0079     /** \brief Scalar type for matrices of type #MatrixType. */
0080     typedef typename MatrixType::Scalar Scalar;
0081     typedef typename NumTraits<Scalar>::Real RealScalar;
0082     typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
0083 
0084     /** \brief Complex scalar type for #MatrixType. 
0085       *
0086       * This is \c std::complex<Scalar> if #Scalar is real (e.g.,
0087       * \c float or \c double) and just \c Scalar if #Scalar is
0088       * complex.
0089       */
0090     typedef std::complex<RealScalar> ComplexScalar;
0091 
0092     /** \brief Type for vector of eigenvalues as returned by eigenvalues(). 
0093       *
0094       * This is a column vector with entries of type #ComplexScalar.
0095       * The length of the vector is the size of #MatrixType.
0096       */
0097     typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> EigenvalueType;
0098 
0099     /** \brief Type for matrix of eigenvectors as returned by eigenvectors(). 
0100       *
0101       * This is a square matrix with entries of type #ComplexScalar. 
0102       * The size is the same as the size of #MatrixType.
0103       */
0104     typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> EigenvectorsType;
0105 
0106     /** \brief Default constructor.
0107       *
0108       * The default constructor is useful in cases in which the user intends to
0109       * perform decompositions via EigenSolver::compute(const MatrixType&, bool).
0110       *
0111       * \sa compute() for an example.
0112       */
0113     EigenSolver() : m_eivec(), m_eivalues(), m_isInitialized(false), m_eigenvectorsOk(false), m_realSchur(), m_matT(), m_tmp() {}
0114 
0115     /** \brief Default constructor with memory preallocation
0116       *
0117       * Like the default constructor but with preallocation of the internal data
0118       * according to the specified problem \a size.
0119       * \sa EigenSolver()
0120       */
0121     explicit EigenSolver(Index size)
0122       : m_eivec(size, size),
0123         m_eivalues(size),
0124         m_isInitialized(false),
0125         m_eigenvectorsOk(false),
0126         m_realSchur(size),
0127         m_matT(size, size), 
0128         m_tmp(size)
0129     {}
0130 
0131     /** \brief Constructor; computes eigendecomposition of given matrix. 
0132       * 
0133       * \param[in]  matrix  Square matrix whose eigendecomposition is to be computed.
0134       * \param[in]  computeEigenvectors  If true, both the eigenvectors and the
0135       *    eigenvalues are computed; if false, only the eigenvalues are
0136       *    computed. 
0137       *
0138       * This constructor calls compute() to compute the eigenvalues
0139       * and eigenvectors.
0140       *
0141       * Example: \include EigenSolver_EigenSolver_MatrixType.cpp
0142       * Output: \verbinclude EigenSolver_EigenSolver_MatrixType.out
0143       *
0144       * \sa compute()
0145       */
0146     template<typename InputType>
0147     explicit EigenSolver(const EigenBase<InputType>& matrix, bool computeEigenvectors = true)
0148       : m_eivec(matrix.rows(), matrix.cols()),
0149         m_eivalues(matrix.cols()),
0150         m_isInitialized(false),
0151         m_eigenvectorsOk(false),
0152         m_realSchur(matrix.cols()),
0153         m_matT(matrix.rows(), matrix.cols()), 
0154         m_tmp(matrix.cols())
0155     {
0156       compute(matrix.derived(), computeEigenvectors);
0157     }
0158 
0159     /** \brief Returns the eigenvectors of given matrix. 
0160       *
0161       * \returns  %Matrix whose columns are the (possibly complex) eigenvectors.
0162       *
0163       * \pre Either the constructor 
0164       * EigenSolver(const MatrixType&,bool) or the member function
0165       * compute(const MatrixType&, bool) has been called before, and
0166       * \p computeEigenvectors was set to true (the default).
0167       *
0168       * Column \f$ k \f$ of the returned matrix is an eigenvector corresponding
0169       * to eigenvalue number \f$ k \f$ as returned by eigenvalues().  The
0170       * eigenvectors are normalized to have (Euclidean) norm equal to one. The
0171       * matrix returned by this function is the matrix \f$ V \f$ in the
0172       * eigendecomposition \f$ A = V D V^{-1} \f$, if it exists.
0173       *
0174       * Example: \include EigenSolver_eigenvectors.cpp
0175       * Output: \verbinclude EigenSolver_eigenvectors.out
0176       *
0177       * \sa eigenvalues(), pseudoEigenvectors()
0178       */
0179     EigenvectorsType eigenvectors() const;
0180 
0181     /** \brief Returns the pseudo-eigenvectors of given matrix. 
0182       *
0183       * \returns  Const reference to matrix whose columns are the pseudo-eigenvectors.
0184       *
0185       * \pre Either the constructor 
0186       * EigenSolver(const MatrixType&,bool) or the member function
0187       * compute(const MatrixType&, bool) has been called before, and
0188       * \p computeEigenvectors was set to true (the default).
0189       *
0190       * The real matrix \f$ V \f$ returned by this function and the
0191       * block-diagonal matrix \f$ D \f$ returned by pseudoEigenvalueMatrix()
0192       * satisfy \f$ AV = VD \f$.
0193       *
0194       * Example: \include EigenSolver_pseudoEigenvectors.cpp
0195       * Output: \verbinclude EigenSolver_pseudoEigenvectors.out
0196       *
0197       * \sa pseudoEigenvalueMatrix(), eigenvectors()
0198       */
0199     const MatrixType& pseudoEigenvectors() const
0200     {
0201       eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
0202       eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
0203       return m_eivec;
0204     }
0205 
0206     /** \brief Returns the block-diagonal matrix in the pseudo-eigendecomposition.
0207       *
0208       * \returns  A block-diagonal matrix.
0209       *
0210       * \pre Either the constructor 
0211       * EigenSolver(const MatrixType&,bool) or the member function
0212       * compute(const MatrixType&, bool) has been called before.
0213       *
0214       * The matrix \f$ D \f$ returned by this function is real and
0215       * block-diagonal. The blocks on the diagonal are either 1-by-1 or 2-by-2
0216       * blocks of the form
0217       * \f$ \begin{bmatrix} u & v \\ -v & u \end{bmatrix} \f$.
0218       * These blocks are not sorted in any particular order.
0219       * The matrix \f$ D \f$ and the matrix \f$ V \f$ returned by
0220       * pseudoEigenvectors() satisfy \f$ AV = VD \f$.
0221       *
0222       * \sa pseudoEigenvectors() for an example, eigenvalues()
0223       */
0224     MatrixType pseudoEigenvalueMatrix() const;
0225 
0226     /** \brief Returns the eigenvalues of given matrix. 
0227       *
0228       * \returns A const reference to the column vector containing the eigenvalues.
0229       *
0230       * \pre Either the constructor 
0231       * EigenSolver(const MatrixType&,bool) or the member function
0232       * compute(const MatrixType&, bool) has been called before.
0233       *
0234       * The eigenvalues are repeated according to their algebraic multiplicity,
0235       * so there are as many eigenvalues as rows in the matrix. The eigenvalues 
0236       * are not sorted in any particular order.
0237       *
0238       * Example: \include EigenSolver_eigenvalues.cpp
0239       * Output: \verbinclude EigenSolver_eigenvalues.out
0240       *
0241       * \sa eigenvectors(), pseudoEigenvalueMatrix(),
0242       *     MatrixBase::eigenvalues()
0243       */
0244     const EigenvalueType& eigenvalues() const
0245     {
0246       eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
0247       return m_eivalues;
0248     }
0249 
0250     /** \brief Computes eigendecomposition of given matrix. 
0251       * 
0252       * \param[in]  matrix  Square matrix whose eigendecomposition is to be computed.
0253       * \param[in]  computeEigenvectors  If true, both the eigenvectors and the
0254       *    eigenvalues are computed; if false, only the eigenvalues are
0255       *    computed. 
0256       * \returns    Reference to \c *this
0257       *
0258       * This function computes the eigenvalues of the real matrix \p matrix.
0259       * The eigenvalues() function can be used to retrieve them.  If 
0260       * \p computeEigenvectors is true, then the eigenvectors are also computed
0261       * and can be retrieved by calling eigenvectors().
0262       *
0263       * The matrix is first reduced to real Schur form using the RealSchur
0264       * class. The Schur decomposition is then used to compute the eigenvalues
0265       * and eigenvectors.
0266       *
0267       * The cost of the computation is dominated by the cost of the
0268       * Schur decomposition, which is very approximately \f$ 25n^3 \f$
0269       * (where \f$ n \f$ is the size of the matrix) if \p computeEigenvectors 
0270       * is true, and \f$ 10n^3 \f$ if \p computeEigenvectors is false.
0271       *
0272       * This method reuses of the allocated data in the EigenSolver object.
0273       *
0274       * Example: \include EigenSolver_compute.cpp
0275       * Output: \verbinclude EigenSolver_compute.out
0276       */
0277     template<typename InputType>
0278     EigenSolver& compute(const EigenBase<InputType>& matrix, bool computeEigenvectors = true);
0279 
0280     /** \returns NumericalIssue if the input contains INF or NaN values or overflow occurred. Returns Success otherwise. */
0281     ComputationInfo info() const
0282     {
0283       eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
0284       return m_info;
0285     }
0286 
0287     /** \brief Sets the maximum number of iterations allowed. */
0288     EigenSolver& setMaxIterations(Index maxIters)
0289     {
0290       m_realSchur.setMaxIterations(maxIters);
0291       return *this;
0292     }
0293 
0294     /** \brief Returns the maximum number of iterations. */
0295     Index getMaxIterations()
0296     {
0297       return m_realSchur.getMaxIterations();
0298     }
0299 
0300   private:
0301     void doComputeEigenvectors();
0302 
0303   protected:
0304     
0305     static void check_template_parameters()
0306     {
0307       EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
0308       EIGEN_STATIC_ASSERT(!NumTraits<Scalar>::IsComplex, NUMERIC_TYPE_MUST_BE_REAL);
0309     }
0310     
0311     MatrixType m_eivec;
0312     EigenvalueType m_eivalues;
0313     bool m_isInitialized;
0314     bool m_eigenvectorsOk;
0315     ComputationInfo m_info;
0316     RealSchur<MatrixType> m_realSchur;
0317     MatrixType m_matT;
0318 
0319     typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
0320     ColumnVectorType m_tmp;
0321 };
0322 
0323 template<typename MatrixType>
0324 MatrixType EigenSolver<MatrixType>::pseudoEigenvalueMatrix() const
0325 {
0326   eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
0327   const RealScalar precision = RealScalar(2)*NumTraits<RealScalar>::epsilon();
0328   Index n = m_eivalues.rows();
0329   MatrixType matD = MatrixType::Zero(n,n);
0330   for (Index i=0; i<n; ++i)
0331   {
0332     if (internal::isMuchSmallerThan(numext::imag(m_eivalues.coeff(i)), numext::real(m_eivalues.coeff(i)), precision))
0333       matD.coeffRef(i,i) = numext::real(m_eivalues.coeff(i));
0334     else
0335     {
0336       matD.template block<2,2>(i,i) <<  numext::real(m_eivalues.coeff(i)), numext::imag(m_eivalues.coeff(i)),
0337                                        -numext::imag(m_eivalues.coeff(i)), numext::real(m_eivalues.coeff(i));
0338       ++i;
0339     }
0340   }
0341   return matD;
0342 }
0343 
0344 template<typename MatrixType>
0345 typename EigenSolver<MatrixType>::EigenvectorsType EigenSolver<MatrixType>::eigenvectors() const
0346 {
0347   eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
0348   eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
0349   const RealScalar precision = RealScalar(2)*NumTraits<RealScalar>::epsilon();
0350   Index n = m_eivec.cols();
0351   EigenvectorsType matV(n,n);
0352   for (Index j=0; j<n; ++j)
0353   {
0354     if (internal::isMuchSmallerThan(numext::imag(m_eivalues.coeff(j)), numext::real(m_eivalues.coeff(j)), precision) || j+1==n)
0355     {
0356       // we have a real eigen value
0357       matV.col(j) = m_eivec.col(j).template cast<ComplexScalar>();
0358       matV.col(j).normalize();
0359     }
0360     else
0361     {
0362       // we have a pair of complex eigen values
0363       for (Index i=0; i<n; ++i)
0364       {
0365         matV.coeffRef(i,j)   = ComplexScalar(m_eivec.coeff(i,j),  m_eivec.coeff(i,j+1));
0366         matV.coeffRef(i,j+1) = ComplexScalar(m_eivec.coeff(i,j), -m_eivec.coeff(i,j+1));
0367       }
0368       matV.col(j).normalize();
0369       matV.col(j+1).normalize();
0370       ++j;
0371     }
0372   }
0373   return matV;
0374 }
0375 
0376 template<typename MatrixType>
0377 template<typename InputType>
0378 EigenSolver<MatrixType>& 
0379 EigenSolver<MatrixType>::compute(const EigenBase<InputType>& matrix, bool computeEigenvectors)
0380 {
0381   check_template_parameters();
0382   
0383   using std::sqrt;
0384   using std::abs;
0385   using numext::isfinite;
0386   eigen_assert(matrix.cols() == matrix.rows());
0387 
0388   // Reduce to real Schur form.
0389   m_realSchur.compute(matrix.derived(), computeEigenvectors);
0390   
0391   m_info = m_realSchur.info();
0392 
0393   if (m_info == Success)
0394   {
0395     m_matT = m_realSchur.matrixT();
0396     if (computeEigenvectors)
0397       m_eivec = m_realSchur.matrixU();
0398   
0399     // Compute eigenvalues from matT
0400     m_eivalues.resize(matrix.cols());
0401     Index i = 0;
0402     while (i < matrix.cols()) 
0403     {
0404       if (i == matrix.cols() - 1 || m_matT.coeff(i+1, i) == Scalar(0)) 
0405       {
0406         m_eivalues.coeffRef(i) = m_matT.coeff(i, i);
0407         if(!(isfinite)(m_eivalues.coeffRef(i)))
0408         {
0409           m_isInitialized = true;
0410           m_eigenvectorsOk = false;
0411           m_info = NumericalIssue;
0412           return *this;
0413         }
0414         ++i;
0415       }
0416       else
0417       {
0418         Scalar p = Scalar(0.5) * (m_matT.coeff(i, i) - m_matT.coeff(i+1, i+1));
0419         Scalar z;
0420         // Compute z = sqrt(abs(p * p + m_matT.coeff(i+1, i) * m_matT.coeff(i, i+1)));
0421         // without overflow
0422         {
0423           Scalar t0 = m_matT.coeff(i+1, i);
0424           Scalar t1 = m_matT.coeff(i, i+1);
0425           Scalar maxval = numext::maxi<Scalar>(abs(p),numext::maxi<Scalar>(abs(t0),abs(t1)));
0426           t0 /= maxval;
0427           t1 /= maxval;
0428           Scalar p0 = p/maxval;
0429           z = maxval * sqrt(abs(p0 * p0 + t0 * t1));
0430         }
0431         
0432         m_eivalues.coeffRef(i)   = ComplexScalar(m_matT.coeff(i+1, i+1) + p, z);
0433         m_eivalues.coeffRef(i+1) = ComplexScalar(m_matT.coeff(i+1, i+1) + p, -z);
0434         if(!((isfinite)(m_eivalues.coeffRef(i)) && (isfinite)(m_eivalues.coeffRef(i+1))))
0435         {
0436           m_isInitialized = true;
0437           m_eigenvectorsOk = false;
0438           m_info = NumericalIssue;
0439           return *this;
0440         }
0441         i += 2;
0442       }
0443     }
0444     
0445     // Compute eigenvectors.
0446     if (computeEigenvectors)
0447       doComputeEigenvectors();
0448   }
0449 
0450   m_isInitialized = true;
0451   m_eigenvectorsOk = computeEigenvectors;
0452 
0453   return *this;
0454 }
0455 
0456 
0457 template<typename MatrixType>
0458 void EigenSolver<MatrixType>::doComputeEigenvectors()
0459 {
0460   using std::abs;
0461   const Index size = m_eivec.cols();
0462   const Scalar eps = NumTraits<Scalar>::epsilon();
0463 
0464   // inefficient! this is already computed in RealSchur
0465   Scalar norm(0);
0466   for (Index j = 0; j < size; ++j)
0467   {
0468     norm += m_matT.row(j).segment((std::max)(j-1,Index(0)), size-(std::max)(j-1,Index(0))).cwiseAbs().sum();
0469   }
0470   
0471   // Backsubstitute to find vectors of upper triangular form
0472   if (norm == Scalar(0))
0473   {
0474     return;
0475   }
0476 
0477   for (Index n = size-1; n >= 0; n--)
0478   {
0479     Scalar p = m_eivalues.coeff(n).real();
0480     Scalar q = m_eivalues.coeff(n).imag();
0481 
0482     // Scalar vector
0483     if (q == Scalar(0))
0484     {
0485       Scalar lastr(0), lastw(0);
0486       Index l = n;
0487 
0488       m_matT.coeffRef(n,n) = Scalar(1);
0489       for (Index i = n-1; i >= 0; i--)
0490       {
0491         Scalar w = m_matT.coeff(i,i) - p;
0492         Scalar r = m_matT.row(i).segment(l,n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
0493 
0494         if (m_eivalues.coeff(i).imag() < Scalar(0))
0495         {
0496           lastw = w;
0497           lastr = r;
0498         }
0499         else
0500         {
0501           l = i;
0502           if (m_eivalues.coeff(i).imag() == Scalar(0))
0503           {
0504             if (w != Scalar(0))
0505               m_matT.coeffRef(i,n) = -r / w;
0506             else
0507               m_matT.coeffRef(i,n) = -r / (eps * norm);
0508           }
0509           else // Solve real equations
0510           {
0511             Scalar x = m_matT.coeff(i,i+1);
0512             Scalar y = m_matT.coeff(i+1,i);
0513             Scalar denom = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) + m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag();
0514             Scalar t = (x * lastr - lastw * r) / denom;
0515             m_matT.coeffRef(i,n) = t;
0516             if (abs(x) > abs(lastw))
0517               m_matT.coeffRef(i+1,n) = (-r - w * t) / x;
0518             else
0519               m_matT.coeffRef(i+1,n) = (-lastr - y * t) / lastw;
0520           }
0521 
0522           // Overflow control
0523           Scalar t = abs(m_matT.coeff(i,n));
0524           if ((eps * t) * t > Scalar(1))
0525             m_matT.col(n).tail(size-i) /= t;
0526         }
0527       }
0528     }
0529     else if (q < Scalar(0) && n > 0) // Complex vector
0530     {
0531       Scalar lastra(0), lastsa(0), lastw(0);
0532       Index l = n-1;
0533 
0534       // Last vector component imaginary so matrix is triangular
0535       if (abs(m_matT.coeff(n,n-1)) > abs(m_matT.coeff(n-1,n)))
0536       {
0537         m_matT.coeffRef(n-1,n-1) = q / m_matT.coeff(n,n-1);
0538         m_matT.coeffRef(n-1,n) = -(m_matT.coeff(n,n) - p) / m_matT.coeff(n,n-1);
0539       }
0540       else
0541       {
0542         ComplexScalar cc = ComplexScalar(Scalar(0),-m_matT.coeff(n-1,n)) / ComplexScalar(m_matT.coeff(n-1,n-1)-p,q);
0543         m_matT.coeffRef(n-1,n-1) = numext::real(cc);
0544         m_matT.coeffRef(n-1,n) = numext::imag(cc);
0545       }
0546       m_matT.coeffRef(n,n-1) = Scalar(0);
0547       m_matT.coeffRef(n,n) = Scalar(1);
0548       for (Index i = n-2; i >= 0; i--)
0549       {
0550         Scalar ra = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n-1).segment(l, n-l+1));
0551         Scalar sa = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
0552         Scalar w = m_matT.coeff(i,i) - p;
0553 
0554         if (m_eivalues.coeff(i).imag() < Scalar(0))
0555         {
0556           lastw = w;
0557           lastra = ra;
0558           lastsa = sa;
0559         }
0560         else
0561         {
0562           l = i;
0563           if (m_eivalues.coeff(i).imag() == RealScalar(0))
0564           {
0565             ComplexScalar cc = ComplexScalar(-ra,-sa) / ComplexScalar(w,q);
0566             m_matT.coeffRef(i,n-1) = numext::real(cc);
0567             m_matT.coeffRef(i,n) = numext::imag(cc);
0568           }
0569           else
0570           {
0571             // Solve complex equations
0572             Scalar x = m_matT.coeff(i,i+1);
0573             Scalar y = m_matT.coeff(i+1,i);
0574             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;
0575             Scalar vi = (m_eivalues.coeff(i).real() - p) * Scalar(2) * q;
0576             if ((vr == Scalar(0)) && (vi == Scalar(0)))
0577               vr = eps * norm * (abs(w) + abs(q) + abs(x) + abs(y) + abs(lastw));
0578 
0579             ComplexScalar cc = ComplexScalar(x*lastra-lastw*ra+q*sa,x*lastsa-lastw*sa-q*ra) / ComplexScalar(vr,vi);
0580             m_matT.coeffRef(i,n-1) = numext::real(cc);
0581             m_matT.coeffRef(i,n) = numext::imag(cc);
0582             if (abs(x) > (abs(lastw) + abs(q)))
0583             {
0584               m_matT.coeffRef(i+1,n-1) = (-ra - w * m_matT.coeff(i,n-1) + q * m_matT.coeff(i,n)) / x;
0585               m_matT.coeffRef(i+1,n) = (-sa - w * m_matT.coeff(i,n) - q * m_matT.coeff(i,n-1)) / x;
0586             }
0587             else
0588             {
0589               cc = ComplexScalar(-lastra-y*m_matT.coeff(i,n-1),-lastsa-y*m_matT.coeff(i,n)) / ComplexScalar(lastw,q);
0590               m_matT.coeffRef(i+1,n-1) = numext::real(cc);
0591               m_matT.coeffRef(i+1,n) = numext::imag(cc);
0592             }
0593           }
0594 
0595           // Overflow control
0596           Scalar t = numext::maxi<Scalar>(abs(m_matT.coeff(i,n-1)),abs(m_matT.coeff(i,n)));
0597           if ((eps * t) * t > Scalar(1))
0598             m_matT.block(i, n-1, size-i, 2) /= t;
0599 
0600         }
0601       }
0602       
0603       // We handled a pair of complex conjugate eigenvalues, so need to skip them both
0604       n--;
0605     }
0606     else
0607     {
0608       eigen_assert(0 && "Internal bug in EigenSolver (INF or NaN has not been detected)"); // this should not happen
0609     }
0610   }
0611 
0612   // Back transformation to get eigenvectors of original matrix
0613   for (Index j = size-1; j >= 0; j--)
0614   {
0615     m_tmp.noalias() = m_eivec.leftCols(j+1) * m_matT.col(j).segment(0, j+1);
0616     m_eivec.col(j) = m_tmp;
0617   }
0618 }
0619 
0620 } // end namespace Eigen
0621 
0622 #endif // EIGEN_EIGENSOLVER_H