Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/eigen3/Eigen/src/Eigenvalues/RealSchur.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_REAL_SCHUR_H
0012 #define EIGEN_REAL_SCHUR_H
0013 
0014 #include "./HessenbergDecomposition.h"
0015 
0016 namespace Eigen { 
0017 
0018 /** \eigenvalues_module \ingroup Eigenvalues_Module
0019   *
0020   *
0021   * \class RealSchur
0022   *
0023   * \brief Performs a real Schur decomposition of a square matrix
0024   *
0025   * \tparam _MatrixType the type of the matrix of which we are computing the
0026   * real Schur decomposition; this is expected to be an instantiation of the
0027   * Matrix class template.
0028   *
0029   * Given a real square matrix A, this class computes the real Schur
0030   * decomposition: \f$ A = U T U^T \f$ where U is a real orthogonal matrix and
0031   * T is a real quasi-triangular matrix. An orthogonal matrix is a matrix whose
0032   * inverse is equal to its transpose, \f$ U^{-1} = U^T \f$. A quasi-triangular
0033   * matrix is a block-triangular matrix whose diagonal consists of 1-by-1
0034   * blocks and 2-by-2 blocks with complex eigenvalues. The eigenvalues of the
0035   * blocks on the diagonal of T are the same as the eigenvalues of the matrix
0036   * A, and thus the real Schur decomposition is used in EigenSolver to compute
0037   * the eigendecomposition of a matrix.
0038   *
0039   * Call the function compute() to compute the real Schur decomposition of a
0040   * given matrix. Alternatively, you can use the RealSchur(const MatrixType&, bool)
0041   * constructor which computes the real Schur decomposition at construction
0042   * time. Once the decomposition is computed, you can use the matrixU() and
0043   * matrixT() functions to retrieve the matrices U and T in the decomposition.
0044   *
0045   * The documentation of RealSchur(const MatrixType&, bool) contains an example
0046   * of the typical use of this class.
0047   *
0048   * \note The implementation is adapted from
0049   * <a href="http://math.nist.gov/javanumerics/jama/">JAMA</a> (public domain).
0050   * Their code is based on EISPACK.
0051   *
0052   * \sa class ComplexSchur, class EigenSolver, class ComplexEigenSolver
0053   */
0054 template<typename _MatrixType> class RealSchur
0055 {
0056   public:
0057     typedef _MatrixType MatrixType;
0058     enum {
0059       RowsAtCompileTime = MatrixType::RowsAtCompileTime,
0060       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
0061       Options = MatrixType::Options,
0062       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
0063       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
0064     };
0065     typedef typename MatrixType::Scalar Scalar;
0066     typedef std::complex<typename NumTraits<Scalar>::Real> ComplexScalar;
0067     typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
0068 
0069     typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> EigenvalueType;
0070     typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
0071 
0072     /** \brief Default constructor.
0073       *
0074       * \param [in] size  Positive integer, size of the matrix whose Schur decomposition will be computed.
0075       *
0076       * The default constructor is useful in cases in which the user intends to
0077       * perform decompositions via compute().  The \p size parameter is only
0078       * used as a hint. It is not an error to give a wrong \p size, but it may
0079       * impair performance.
0080       *
0081       * \sa compute() for an example.
0082       */
0083     explicit RealSchur(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime)
0084             : m_matT(size, size),
0085               m_matU(size, size),
0086               m_workspaceVector(size),
0087               m_hess(size),
0088               m_isInitialized(false),
0089               m_matUisUptodate(false),
0090               m_maxIters(-1)
0091     { }
0092 
0093     /** \brief Constructor; computes real Schur decomposition of given matrix. 
0094       * 
0095       * \param[in]  matrix    Square matrix whose Schur decomposition is to be computed.
0096       * \param[in]  computeU  If true, both T and U are computed; if false, only T is computed.
0097       *
0098       * This constructor calls compute() to compute the Schur decomposition.
0099       *
0100       * Example: \include RealSchur_RealSchur_MatrixType.cpp
0101       * Output: \verbinclude RealSchur_RealSchur_MatrixType.out
0102       */
0103     template<typename InputType>
0104     explicit RealSchur(const EigenBase<InputType>& matrix, bool computeU = true)
0105             : m_matT(matrix.rows(),matrix.cols()),
0106               m_matU(matrix.rows(),matrix.cols()),
0107               m_workspaceVector(matrix.rows()),
0108               m_hess(matrix.rows()),
0109               m_isInitialized(false),
0110               m_matUisUptodate(false),
0111               m_maxIters(-1)
0112     {
0113       compute(matrix.derived(), computeU);
0114     }
0115 
0116     /** \brief Returns the orthogonal matrix in the Schur decomposition. 
0117       *
0118       * \returns A const reference to the matrix U.
0119       *
0120       * \pre Either the constructor RealSchur(const MatrixType&, bool) or the
0121       * member function compute(const MatrixType&, bool) has been called before
0122       * to compute the Schur decomposition of a matrix, and \p computeU was set
0123       * to true (the default value).
0124       *
0125       * \sa RealSchur(const MatrixType&, bool) for an example
0126       */
0127     const MatrixType& matrixU() const
0128     {
0129       eigen_assert(m_isInitialized && "RealSchur is not initialized.");
0130       eigen_assert(m_matUisUptodate && "The matrix U has not been computed during the RealSchur decomposition.");
0131       return m_matU;
0132     }
0133 
0134     /** \brief Returns the quasi-triangular matrix in the Schur decomposition. 
0135       *
0136       * \returns A const reference to the matrix T.
0137       *
0138       * \pre Either the constructor RealSchur(const MatrixType&, bool) or the
0139       * member function compute(const MatrixType&, bool) has been called before
0140       * to compute the Schur decomposition of a matrix.
0141       *
0142       * \sa RealSchur(const MatrixType&, bool) for an example
0143       */
0144     const MatrixType& matrixT() const
0145     {
0146       eigen_assert(m_isInitialized && "RealSchur is not initialized.");
0147       return m_matT;
0148     }
0149   
0150     /** \brief Computes Schur decomposition of given matrix. 
0151       * 
0152       * \param[in]  matrix    Square matrix whose Schur decomposition is to be computed.
0153       * \param[in]  computeU  If true, both T and U are computed; if false, only T is computed.
0154       * \returns    Reference to \c *this
0155       *
0156       * The Schur decomposition is computed by first reducing the matrix to
0157       * Hessenberg form using the class HessenbergDecomposition. The Hessenberg
0158       * matrix is then reduced to triangular form by performing Francis QR
0159       * iterations with implicit double shift. The cost of computing the Schur
0160       * decomposition depends on the number of iterations; as a rough guide, it
0161       * may be taken to be \f$25n^3\f$ flops if \a computeU is true and
0162       * \f$10n^3\f$ flops if \a computeU is false.
0163       *
0164       * Example: \include RealSchur_compute.cpp
0165       * Output: \verbinclude RealSchur_compute.out
0166       *
0167       * \sa compute(const MatrixType&, bool, Index)
0168       */
0169     template<typename InputType>
0170     RealSchur& compute(const EigenBase<InputType>& matrix, bool computeU = true);
0171 
0172     /** \brief Computes Schur decomposition of a Hessenberg matrix H = Z T Z^T
0173      *  \param[in] matrixH Matrix in Hessenberg form H
0174      *  \param[in] matrixQ orthogonal matrix Q that transform a matrix A to H : A = Q H Q^T
0175      *  \param computeU Computes the matriX U of the Schur vectors
0176      * \return Reference to \c *this
0177      * 
0178      *  This routine assumes that the matrix is already reduced in Hessenberg form matrixH
0179      *  using either the class HessenbergDecomposition or another mean. 
0180      *  It computes the upper quasi-triangular matrix T of the Schur decomposition of H
0181      *  When computeU is true, this routine computes the matrix U such that 
0182      *  A = U T U^T =  (QZ) T (QZ)^T = Q H Q^T where A is the initial matrix
0183      * 
0184      * NOTE Q is referenced if computeU is true; so, if the initial orthogonal matrix
0185      * is not available, the user should give an identity matrix (Q.setIdentity())
0186      * 
0187      * \sa compute(const MatrixType&, bool)
0188      */
0189     template<typename HessMatrixType, typename OrthMatrixType>
0190     RealSchur& computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ,  bool computeU);
0191     /** \brief Reports whether previous computation was successful.
0192       *
0193       * \returns \c Success if computation was successful, \c NoConvergence otherwise.
0194       */
0195     ComputationInfo info() const
0196     {
0197       eigen_assert(m_isInitialized && "RealSchur is not initialized.");
0198       return m_info;
0199     }
0200 
0201     /** \brief Sets the maximum number of iterations allowed. 
0202       *
0203       * If not specified by the user, the maximum number of iterations is m_maxIterationsPerRow times the size
0204       * of the matrix.
0205       */
0206     RealSchur& setMaxIterations(Index maxIters)
0207     {
0208       m_maxIters = maxIters;
0209       return *this;
0210     }
0211 
0212     /** \brief Returns the maximum number of iterations. */
0213     Index getMaxIterations()
0214     {
0215       return m_maxIters;
0216     }
0217 
0218     /** \brief Maximum number of iterations per row.
0219       *
0220       * If not otherwise specified, the maximum number of iterations is this number times the size of the
0221       * matrix. It is currently set to 40.
0222       */
0223     static const int m_maxIterationsPerRow = 40;
0224 
0225   private:
0226     
0227     MatrixType m_matT;
0228     MatrixType m_matU;
0229     ColumnVectorType m_workspaceVector;
0230     HessenbergDecomposition<MatrixType> m_hess;
0231     ComputationInfo m_info;
0232     bool m_isInitialized;
0233     bool m_matUisUptodate;
0234     Index m_maxIters;
0235 
0236     typedef Matrix<Scalar,3,1> Vector3s;
0237 
0238     Scalar computeNormOfT();
0239     Index findSmallSubdiagEntry(Index iu, const Scalar& considerAsZero);
0240     void splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift);
0241     void computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo);
0242     void initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector);
0243     void performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace);
0244 };
0245 
0246 
0247 template<typename MatrixType>
0248 template<typename InputType>
0249 RealSchur<MatrixType>& RealSchur<MatrixType>::compute(const EigenBase<InputType>& matrix, bool computeU)
0250 {
0251   const Scalar considerAsZero = (std::numeric_limits<Scalar>::min)();
0252 
0253   eigen_assert(matrix.cols() == matrix.rows());
0254   Index maxIters = m_maxIters;
0255   if (maxIters == -1)
0256     maxIters = m_maxIterationsPerRow * matrix.rows();
0257 
0258   Scalar scale = matrix.derived().cwiseAbs().maxCoeff();
0259   if(scale<considerAsZero)
0260   {
0261     m_matT.setZero(matrix.rows(),matrix.cols());
0262     if(computeU)
0263       m_matU.setIdentity(matrix.rows(),matrix.cols());
0264     m_info = Success;
0265     m_isInitialized = true;
0266     m_matUisUptodate = computeU;
0267     return *this;
0268   }
0269 
0270   // Step 1. Reduce to Hessenberg form
0271   m_hess.compute(matrix.derived()/scale);
0272 
0273   // Step 2. Reduce to real Schur form
0274   // Note: we copy m_hess.matrixQ() into m_matU here and not in computeFromHessenberg
0275   //       to be able to pass our working-space buffer for the Householder to Dense evaluation.
0276   m_workspaceVector.resize(matrix.cols());
0277   if(computeU)
0278     m_hess.matrixQ().evalTo(m_matU, m_workspaceVector);
0279   computeFromHessenberg(m_hess.matrixH(), m_matU, computeU);
0280 
0281   m_matT *= scale;
0282   
0283   return *this;
0284 }
0285 template<typename MatrixType>
0286 template<typename HessMatrixType, typename OrthMatrixType>
0287 RealSchur<MatrixType>& RealSchur<MatrixType>::computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ,  bool computeU)
0288 {
0289   using std::abs;
0290 
0291   m_matT = matrixH;
0292   m_workspaceVector.resize(m_matT.cols());
0293   if(computeU && !internal::is_same_dense(m_matU,matrixQ))
0294     m_matU = matrixQ;
0295   
0296   Index maxIters = m_maxIters;
0297   if (maxIters == -1)
0298     maxIters = m_maxIterationsPerRow * matrixH.rows();
0299   Scalar* workspace = &m_workspaceVector.coeffRef(0);
0300 
0301   // The matrix m_matT is divided in three parts. 
0302   // Rows 0,...,il-1 are decoupled from the rest because m_matT(il,il-1) is zero. 
0303   // Rows il,...,iu is the part we are working on (the active window).
0304   // Rows iu+1,...,end are already brought in triangular form.
0305   Index iu = m_matT.cols() - 1;
0306   Index iter = 0;      // iteration count for current eigenvalue
0307   Index totalIter = 0; // iteration count for whole matrix
0308   Scalar exshift(0);   // sum of exceptional shifts
0309   Scalar norm = computeNormOfT();
0310   // sub-diagonal entries smaller than considerAsZero will be treated as zero.
0311   // We use eps^2 to enable more precision in small eigenvalues.
0312   Scalar considerAsZero = numext::maxi<Scalar>( norm * numext::abs2(NumTraits<Scalar>::epsilon()),
0313                                                 (std::numeric_limits<Scalar>::min)() );
0314 
0315   if(norm!=Scalar(0))
0316   {
0317     while (iu >= 0)
0318     {
0319       Index il = findSmallSubdiagEntry(iu,considerAsZero);
0320 
0321       // Check for convergence
0322       if (il == iu) // One root found
0323       {
0324         m_matT.coeffRef(iu,iu) = m_matT.coeff(iu,iu) + exshift;
0325         if (iu > 0)
0326           m_matT.coeffRef(iu, iu-1) = Scalar(0);
0327         iu--;
0328         iter = 0;
0329       }
0330       else if (il == iu-1) // Two roots found
0331       {
0332         splitOffTwoRows(iu, computeU, exshift);
0333         iu -= 2;
0334         iter = 0;
0335       }
0336       else // No convergence yet
0337       {
0338         // The firstHouseholderVector vector has to be initialized to something to get rid of a silly GCC warning (-O1 -Wall -DNDEBUG )
0339         Vector3s firstHouseholderVector = Vector3s::Zero(), shiftInfo;
0340         computeShift(iu, iter, exshift, shiftInfo);
0341         iter = iter + 1;
0342         totalIter = totalIter + 1;
0343         if (totalIter > maxIters) break;
0344         Index im;
0345         initFrancisQRStep(il, iu, shiftInfo, im, firstHouseholderVector);
0346         performFrancisQRStep(il, im, iu, computeU, firstHouseholderVector, workspace);
0347       }
0348     }
0349   }
0350   if(totalIter <= maxIters)
0351     m_info = Success;
0352   else
0353     m_info = NoConvergence;
0354 
0355   m_isInitialized = true;
0356   m_matUisUptodate = computeU;
0357   return *this;
0358 }
0359 
0360 /** \internal Computes and returns vector L1 norm of T */
0361 template<typename MatrixType>
0362 inline typename MatrixType::Scalar RealSchur<MatrixType>::computeNormOfT()
0363 {
0364   const Index size = m_matT.cols();
0365   // FIXME to be efficient the following would requires a triangular reduxion code
0366   // Scalar norm = m_matT.upper().cwiseAbs().sum() 
0367   //               + m_matT.bottomLeftCorner(size-1,size-1).diagonal().cwiseAbs().sum();
0368   Scalar norm(0);
0369   for (Index j = 0; j < size; ++j)
0370     norm += m_matT.col(j).segment(0, (std::min)(size,j+2)).cwiseAbs().sum();
0371   return norm;
0372 }
0373 
0374 /** \internal Look for single small sub-diagonal element and returns its index */
0375 template<typename MatrixType>
0376 inline Index RealSchur<MatrixType>::findSmallSubdiagEntry(Index iu, const Scalar& considerAsZero)
0377 {
0378   using std::abs;
0379   Index res = iu;
0380   while (res > 0)
0381   {
0382     Scalar s = abs(m_matT.coeff(res-1,res-1)) + abs(m_matT.coeff(res,res));
0383 
0384     s = numext::maxi<Scalar>(s * NumTraits<Scalar>::epsilon(), considerAsZero);
0385     
0386     if (abs(m_matT.coeff(res,res-1)) <= s)
0387       break;
0388     res--;
0389   }
0390   return res;
0391 }
0392 
0393 /** \internal Update T given that rows iu-1 and iu decouple from the rest. */
0394 template<typename MatrixType>
0395 inline void RealSchur<MatrixType>::splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift)
0396 {
0397   using std::sqrt;
0398   using std::abs;
0399   const Index size = m_matT.cols();
0400 
0401   // The eigenvalues of the 2x2 matrix [a b; c d] are 
0402   // trace +/- sqrt(discr/4) where discr = tr^2 - 4*det, tr = a + d, det = ad - bc
0403   Scalar p = Scalar(0.5) * (m_matT.coeff(iu-1,iu-1) - m_matT.coeff(iu,iu));
0404   Scalar q = p * p + m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);   // q = tr^2 / 4 - det = discr/4
0405   m_matT.coeffRef(iu,iu) += exshift;
0406   m_matT.coeffRef(iu-1,iu-1) += exshift;
0407 
0408   if (q >= Scalar(0)) // Two real eigenvalues
0409   {
0410     Scalar z = sqrt(abs(q));
0411     JacobiRotation<Scalar> rot;
0412     if (p >= Scalar(0))
0413       rot.makeGivens(p + z, m_matT.coeff(iu, iu-1));
0414     else
0415       rot.makeGivens(p - z, m_matT.coeff(iu, iu-1));
0416 
0417     m_matT.rightCols(size-iu+1).applyOnTheLeft(iu-1, iu, rot.adjoint());
0418     m_matT.topRows(iu+1).applyOnTheRight(iu-1, iu, rot);
0419     m_matT.coeffRef(iu, iu-1) = Scalar(0); 
0420     if (computeU)
0421       m_matU.applyOnTheRight(iu-1, iu, rot);
0422   }
0423 
0424   if (iu > 1) 
0425     m_matT.coeffRef(iu-1, iu-2) = Scalar(0);
0426 }
0427 
0428 /** \internal Form shift in shiftInfo, and update exshift if an exceptional shift is performed. */
0429 template<typename MatrixType>
0430 inline void RealSchur<MatrixType>::computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo)
0431 {
0432   using std::sqrt;
0433   using std::abs;
0434   shiftInfo.coeffRef(0) = m_matT.coeff(iu,iu);
0435   shiftInfo.coeffRef(1) = m_matT.coeff(iu-1,iu-1);
0436   shiftInfo.coeffRef(2) = m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
0437 
0438   // Wilkinson's original ad hoc shift
0439   if (iter == 10)
0440   {
0441     exshift += shiftInfo.coeff(0);
0442     for (Index i = 0; i <= iu; ++i)
0443       m_matT.coeffRef(i,i) -= shiftInfo.coeff(0);
0444     Scalar s = abs(m_matT.coeff(iu,iu-1)) + abs(m_matT.coeff(iu-1,iu-2));
0445     shiftInfo.coeffRef(0) = Scalar(0.75) * s;
0446     shiftInfo.coeffRef(1) = Scalar(0.75) * s;
0447     shiftInfo.coeffRef(2) = Scalar(-0.4375) * s * s;
0448   }
0449 
0450   // MATLAB's new ad hoc shift
0451   if (iter == 30)
0452   {
0453     Scalar s = (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
0454     s = s * s + shiftInfo.coeff(2);
0455     if (s > Scalar(0))
0456     {
0457       s = sqrt(s);
0458       if (shiftInfo.coeff(1) < shiftInfo.coeff(0))
0459         s = -s;
0460       s = s + (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
0461       s = shiftInfo.coeff(0) - shiftInfo.coeff(2) / s;
0462       exshift += s;
0463       for (Index i = 0; i <= iu; ++i)
0464         m_matT.coeffRef(i,i) -= s;
0465       shiftInfo.setConstant(Scalar(0.964));
0466     }
0467   }
0468 }
0469 
0470 /** \internal Compute index im at which Francis QR step starts and the first Householder vector. */
0471 template<typename MatrixType>
0472 inline void RealSchur<MatrixType>::initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector)
0473 {
0474   using std::abs;
0475   Vector3s& v = firstHouseholderVector; // alias to save typing
0476 
0477   for (im = iu-2; im >= il; --im)
0478   {
0479     const Scalar Tmm = m_matT.coeff(im,im);
0480     const Scalar r = shiftInfo.coeff(0) - Tmm;
0481     const Scalar s = shiftInfo.coeff(1) - Tmm;
0482     v.coeffRef(0) = (r * s - shiftInfo.coeff(2)) / m_matT.coeff(im+1,im) + m_matT.coeff(im,im+1);
0483     v.coeffRef(1) = m_matT.coeff(im+1,im+1) - Tmm - r - s;
0484     v.coeffRef(2) = m_matT.coeff(im+2,im+1);
0485     if (im == il) {
0486       break;
0487     }
0488     const Scalar lhs = m_matT.coeff(im,im-1) * (abs(v.coeff(1)) + abs(v.coeff(2)));
0489     const Scalar rhs = v.coeff(0) * (abs(m_matT.coeff(im-1,im-1)) + abs(Tmm) + abs(m_matT.coeff(im+1,im+1)));
0490     if (abs(lhs) < NumTraits<Scalar>::epsilon() * rhs)
0491       break;
0492   }
0493 }
0494 
0495 /** \internal Perform a Francis QR step involving rows il:iu and columns im:iu. */
0496 template<typename MatrixType>
0497 inline void RealSchur<MatrixType>::performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace)
0498 {
0499   eigen_assert(im >= il);
0500   eigen_assert(im <= iu-2);
0501 
0502   const Index size = m_matT.cols();
0503 
0504   for (Index k = im; k <= iu-2; ++k)
0505   {
0506     bool firstIteration = (k == im);
0507 
0508     Vector3s v;
0509     if (firstIteration)
0510       v = firstHouseholderVector;
0511     else
0512       v = m_matT.template block<3,1>(k,k-1);
0513 
0514     Scalar tau, beta;
0515     Matrix<Scalar, 2, 1> ess;
0516     v.makeHouseholder(ess, tau, beta);
0517     
0518     if (beta != Scalar(0)) // if v is not zero
0519     {
0520       if (firstIteration && k > il)
0521         m_matT.coeffRef(k,k-1) = -m_matT.coeff(k,k-1);
0522       else if (!firstIteration)
0523         m_matT.coeffRef(k,k-1) = beta;
0524 
0525       // These Householder transformations form the O(n^3) part of the algorithm
0526       m_matT.block(k, k, 3, size-k).applyHouseholderOnTheLeft(ess, tau, workspace);
0527       m_matT.block(0, k, (std::min)(iu,k+3) + 1, 3).applyHouseholderOnTheRight(ess, tau, workspace);
0528       if (computeU)
0529         m_matU.block(0, k, size, 3).applyHouseholderOnTheRight(ess, tau, workspace);
0530     }
0531   }
0532 
0533   Matrix<Scalar, 2, 1> v = m_matT.template block<2,1>(iu-1, iu-2);
0534   Scalar tau, beta;
0535   Matrix<Scalar, 1, 1> ess;
0536   v.makeHouseholder(ess, tau, beta);
0537 
0538   if (beta != Scalar(0)) // if v is not zero
0539   {
0540     m_matT.coeffRef(iu-1, iu-2) = beta;
0541     m_matT.block(iu-1, iu-1, 2, size-iu+1).applyHouseholderOnTheLeft(ess, tau, workspace);
0542     m_matT.block(0, iu-1, iu+1, 2).applyHouseholderOnTheRight(ess, tau, workspace);
0543     if (computeU)
0544       m_matU.block(0, iu-1, size, 2).applyHouseholderOnTheRight(ess, tau, workspace);
0545   }
0546 
0547   // clean up pollution due to round-off errors
0548   for (Index i = im+2; i <= iu; ++i)
0549   {
0550     m_matT.coeffRef(i,i-2) = Scalar(0);
0551     if (i > im+2)
0552       m_matT.coeffRef(i,i-3) = Scalar(0);
0553   }
0554 }
0555 
0556 } // end namespace Eigen
0557 
0558 #endif // EIGEN_REAL_SCHUR_H