RealSchur.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 // Copyright (C) 2010 Jitse Niesen <jitse@maths.leeds.ac.uk>
00006 //
00007 // Eigen is free software; you can redistribute it and/or
00008 // modify it under the terms of the GNU Lesser General Public
00009 // License as published by the Free Software Foundation; either
00010 // version 3 of the License, or (at your option) any later version.
00011 //
00012 // Alternatively, you can redistribute it and/or
00013 // modify it under the terms of the GNU General Public License as
00014 // published by the Free Software Foundation; either version 2 of
00015 // the License, or (at your option) any later version.
00016 //
00017 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00018 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00019 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00020 // GNU General Public License for more details.
00021 //
00022 // You should have received a copy of the GNU Lesser General Public
00023 // License and a copy of the GNU General Public License along with
00024 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00025 
00026 #ifndef EIGEN_REAL_SCHUR_H
00027 #define EIGEN_REAL_SCHUR_H
00028 
00029 #include "./HessenbergDecomposition.h"
00030 
00031 namespace Eigen { 
00032 
00069 template<typename _MatrixType> class RealSchur
00070 {
00071   public:
00072     typedef _MatrixType MatrixType;
00073     enum {
00074       RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00075       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00076       Options = MatrixType::Options,
00077       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00078       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
00079     };
00080     typedef typename MatrixType::Scalar Scalar;
00081     typedef std::complex<typename NumTraits<Scalar>::Real> ComplexScalar;
00082     typedef typename MatrixType::Index Index;
00083 
00084     typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> EigenvalueType;
00085     typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
00086 
00098     RealSchur(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime)
00099             : m_matT(size, size),
00100               m_matU(size, size),
00101               m_workspaceVector(size),
00102               m_hess(size),
00103               m_isInitialized(false),
00104               m_matUisUptodate(false)
00105     { }
00106 
00117     RealSchur(const MatrixType& matrix, bool computeU = true)
00118             : m_matT(matrix.rows(),matrix.cols()),
00119               m_matU(matrix.rows(),matrix.cols()),
00120               m_workspaceVector(matrix.rows()),
00121               m_hess(matrix.rows()),
00122               m_isInitialized(false),
00123               m_matUisUptodate(false)
00124     {
00125       compute(matrix, computeU);
00126     }
00127 
00139     const MatrixType& matrixU() const
00140     {
00141       eigen_assert(m_isInitialized && "RealSchur is not initialized.");
00142       eigen_assert(m_matUisUptodate && "The matrix U has not been computed during the RealSchur decomposition.");
00143       return m_matU;
00144     }
00145 
00156     const MatrixType& matrixT() const
00157     {
00158       eigen_assert(m_isInitialized && "RealSchur is not initialized.");
00159       return m_matT;
00160     }
00161   
00179     RealSchur& compute(const MatrixType& matrix, bool computeU = true);
00180 
00185     ComputationInfo info() const
00186     {
00187       eigen_assert(m_isInitialized && "RealSchur is not initialized.");
00188       return m_info;
00189     }
00190 
00195     static const int m_maxIterations = 40;
00196 
00197   private:
00198     
00199     MatrixType m_matT;
00200     MatrixType m_matU;
00201     ColumnVectorType m_workspaceVector;
00202     HessenbergDecomposition<MatrixType> m_hess;
00203     ComputationInfo m_info;
00204     bool m_isInitialized;
00205     bool m_matUisUptodate;
00206 
00207     typedef Matrix<Scalar,3,1> Vector3s;
00208 
00209     Scalar computeNormOfT();
00210     Index findSmallSubdiagEntry(Index iu, Scalar norm);
00211     void splitOffTwoRows(Index iu, bool computeU, Scalar exshift);
00212     void computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo);
00213     void initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector);
00214     void performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace);
00215 };
00216 
00217 
00218 template<typename MatrixType>
00219 RealSchur<MatrixType>& RealSchur<MatrixType>::compute(const MatrixType& matrix, bool computeU)
00220 {
00221   assert(matrix.cols() == matrix.rows());
00222 
00223   // Step 1. Reduce to Hessenberg form
00224   m_hess.compute(matrix);
00225   m_matT = m_hess.matrixH();
00226   if (computeU)
00227     m_matU = m_hess.matrixQ();
00228 
00229   // Step 2. Reduce to real Schur form  
00230   m_workspaceVector.resize(m_matT.cols());
00231   Scalar* workspace = &m_workspaceVector.coeffRef(0);
00232 
00233   // The matrix m_matT is divided in three parts. 
00234   // Rows 0,...,il-1 are decoupled from the rest because m_matT(il,il-1) is zero. 
00235   // Rows il,...,iu is the part we are working on (the active window).
00236   // Rows iu+1,...,end are already brought in triangular form.
00237   Index iu = m_matT.cols() - 1;
00238   Index iter = 0; // iteration count
00239   Scalar exshift(0); // sum of exceptional shifts
00240   Scalar norm = computeNormOfT();
00241 
00242   while (iu >= 0)
00243   {
00244     Index il = findSmallSubdiagEntry(iu, norm);
00245 
00246     // Check for convergence
00247     if (il == iu) // One root found
00248     {
00249       m_matT.coeffRef(iu,iu) = m_matT.coeff(iu,iu) + exshift;
00250       if (iu > 0) 
00251         m_matT.coeffRef(iu, iu-1) = Scalar(0);
00252       iu--;
00253       iter = 0;
00254     }
00255     else if (il == iu-1) // Two roots found
00256     {
00257       splitOffTwoRows(iu, computeU, exshift);
00258       iu -= 2;
00259       iter = 0;
00260     }
00261     else // No convergence yet
00262     {
00263       // The firstHouseholderVector vector has to be initialized to something to get rid of a silly GCC warning (-O1 -Wall -DNDEBUG )
00264       Vector3s firstHouseholderVector(0,0,0), shiftInfo;
00265       computeShift(iu, iter, exshift, shiftInfo);
00266       iter = iter + 1; 
00267       if (iter > m_maxIterations) break;
00268       Index im;
00269       initFrancisQRStep(il, iu, shiftInfo, im, firstHouseholderVector);
00270       performFrancisQRStep(il, im, iu, computeU, firstHouseholderVector, workspace);
00271     }
00272   } 
00273 
00274   if(iter <= m_maxIterations) 
00275     m_info = Success;
00276   else
00277     m_info = NoConvergence;
00278 
00279   m_isInitialized = true;
00280   m_matUisUptodate = computeU;
00281   return *this;
00282 }
00283 
00285 template<typename MatrixType>
00286 inline typename MatrixType::Scalar RealSchur<MatrixType>::computeNormOfT()
00287 {
00288   const Index size = m_matT.cols();
00289   // FIXME to be efficient the following would requires a triangular reduxion code
00290   // Scalar norm = m_matT.upper().cwiseAbs().sum() 
00291   //               + m_matT.bottomLeftCorner(size-1,size-1).diagonal().cwiseAbs().sum();
00292   Scalar norm(0);
00293   for (Index j = 0; j < size; ++j)
00294     norm += m_matT.row(j).segment((std::max)(j-1,Index(0)), size-(std::max)(j-1,Index(0))).cwiseAbs().sum();
00295   return norm;
00296 }
00297 
00299 template<typename MatrixType>
00300 inline typename MatrixType::Index RealSchur<MatrixType>::findSmallSubdiagEntry(Index iu, Scalar norm)
00301 {
00302   Index res = iu;
00303   while (res > 0)
00304   {
00305     Scalar s = internal::abs(m_matT.coeff(res-1,res-1)) + internal::abs(m_matT.coeff(res,res));
00306     if (s == 0.0)
00307       s = norm;
00308     if (internal::abs(m_matT.coeff(res,res-1)) < NumTraits<Scalar>::epsilon() * s)
00309       break;
00310     res--;
00311   }
00312   return res;
00313 }
00314 
00316 template<typename MatrixType>
00317 inline void RealSchur<MatrixType>::splitOffTwoRows(Index iu, bool computeU, Scalar exshift)
00318 {
00319   const Index size = m_matT.cols();
00320 
00321   // The eigenvalues of the 2x2 matrix [a b; c d] are 
00322   // trace +/- sqrt(discr/4) where discr = tr^2 - 4*det, tr = a + d, det = ad - bc
00323   Scalar p = Scalar(0.5) * (m_matT.coeff(iu-1,iu-1) - m_matT.coeff(iu,iu));
00324   Scalar q = p * p + m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);   // q = tr^2 / 4 - det = discr/4
00325   m_matT.coeffRef(iu,iu) += exshift;
00326   m_matT.coeffRef(iu-1,iu-1) += exshift;
00327 
00328   if (q >= Scalar(0)) // Two real eigenvalues
00329   {
00330     Scalar z = internal::sqrt(internal::abs(q));
00331     JacobiRotation<Scalar> rot;
00332     if (p >= Scalar(0))
00333       rot.makeGivens(p + z, m_matT.coeff(iu, iu-1));
00334     else
00335       rot.makeGivens(p - z, m_matT.coeff(iu, iu-1));
00336 
00337     m_matT.rightCols(size-iu+1).applyOnTheLeft(iu-1, iu, rot.adjoint());
00338     m_matT.topRows(iu+1).applyOnTheRight(iu-1, iu, rot);
00339     m_matT.coeffRef(iu, iu-1) = Scalar(0); 
00340     if (computeU)
00341       m_matU.applyOnTheRight(iu-1, iu, rot);
00342   }
00343 
00344   if (iu > 1) 
00345     m_matT.coeffRef(iu-1, iu-2) = Scalar(0);
00346 }
00347 
00349 template<typename MatrixType>
00350 inline void RealSchur<MatrixType>::computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo)
00351 {
00352   shiftInfo.coeffRef(0) = m_matT.coeff(iu,iu);
00353   shiftInfo.coeffRef(1) = m_matT.coeff(iu-1,iu-1);
00354   shiftInfo.coeffRef(2) = m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
00355 
00356   // Wilkinson's original ad hoc shift
00357   if (iter == 10)
00358   {
00359     exshift += shiftInfo.coeff(0);
00360     for (Index i = 0; i <= iu; ++i)
00361       m_matT.coeffRef(i,i) -= shiftInfo.coeff(0);
00362     Scalar s = internal::abs(m_matT.coeff(iu,iu-1)) + internal::abs(m_matT.coeff(iu-1,iu-2));
00363     shiftInfo.coeffRef(0) = Scalar(0.75) * s;
00364     shiftInfo.coeffRef(1) = Scalar(0.75) * s;
00365     shiftInfo.coeffRef(2) = Scalar(-0.4375) * s * s;
00366   }
00367 
00368   // MATLAB's new ad hoc shift
00369   if (iter == 30)
00370   {
00371     Scalar s = (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
00372     s = s * s + shiftInfo.coeff(2);
00373     if (s > Scalar(0))
00374     {
00375       s = internal::sqrt(s);
00376       if (shiftInfo.coeff(1) < shiftInfo.coeff(0))
00377         s = -s;
00378       s = s + (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
00379       s = shiftInfo.coeff(0) - shiftInfo.coeff(2) / s;
00380       exshift += s;
00381       for (Index i = 0; i <= iu; ++i)
00382         m_matT.coeffRef(i,i) -= s;
00383       shiftInfo.setConstant(Scalar(0.964));
00384     }
00385   }
00386 }
00387 
00389 template<typename MatrixType>
00390 inline void RealSchur<MatrixType>::initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector)
00391 {
00392   Vector3s& v = firstHouseholderVector; // alias to save typing
00393 
00394   for (im = iu-2; im >= il; --im)
00395   {
00396     const Scalar Tmm = m_matT.coeff(im,im);
00397     const Scalar r = shiftInfo.coeff(0) - Tmm;
00398     const Scalar s = shiftInfo.coeff(1) - Tmm;
00399     v.coeffRef(0) = (r * s - shiftInfo.coeff(2)) / m_matT.coeff(im+1,im) + m_matT.coeff(im,im+1);
00400     v.coeffRef(1) = m_matT.coeff(im+1,im+1) - Tmm - r - s;
00401     v.coeffRef(2) = m_matT.coeff(im+2,im+1);
00402     if (im == il) {
00403       break;
00404     }
00405     const Scalar lhs = m_matT.coeff(im,im-1) * (internal::abs(v.coeff(1)) + internal::abs(v.coeff(2)));
00406     const Scalar rhs = v.coeff(0) * (internal::abs(m_matT.coeff(im-1,im-1)) + internal::abs(Tmm) + internal::abs(m_matT.coeff(im+1,im+1)));
00407     if (internal::abs(lhs) < NumTraits<Scalar>::epsilon() * rhs)
00408     {
00409       break;
00410     }
00411   }
00412 }
00413 
00415 template<typename MatrixType>
00416 inline void RealSchur<MatrixType>::performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace)
00417 {
00418   assert(im >= il);
00419   assert(im <= iu-2);
00420 
00421   const Index size = m_matT.cols();
00422 
00423   for (Index k = im; k <= iu-2; ++k)
00424   {
00425     bool firstIteration = (k == im);
00426 
00427     Vector3s v;
00428     if (firstIteration)
00429       v = firstHouseholderVector;
00430     else
00431       v = m_matT.template block<3,1>(k,k-1);
00432 
00433     Scalar tau, beta;
00434     Matrix<Scalar, 2, 1> ess;
00435     v.makeHouseholder(ess, tau, beta);
00436     
00437     if (beta != Scalar(0)) // if v is not zero
00438     {
00439       if (firstIteration && k > il)
00440         m_matT.coeffRef(k,k-1) = -m_matT.coeff(k,k-1);
00441       else if (!firstIteration)
00442         m_matT.coeffRef(k,k-1) = beta;
00443 
00444       // These Householder transformations form the O(n^3) part of the algorithm
00445       m_matT.block(k, k, 3, size-k).applyHouseholderOnTheLeft(ess, tau, workspace);
00446       m_matT.block(0, k, (std::min)(iu,k+3) + 1, 3).applyHouseholderOnTheRight(ess, tau, workspace);
00447       if (computeU)
00448         m_matU.block(0, k, size, 3).applyHouseholderOnTheRight(ess, tau, workspace);
00449     }
00450   }
00451 
00452   Matrix<Scalar, 2, 1> v = m_matT.template block<2,1>(iu-1, iu-2);
00453   Scalar tau, beta;
00454   Matrix<Scalar, 1, 1> ess;
00455   v.makeHouseholder(ess, tau, beta);
00456 
00457   if (beta != Scalar(0)) // if v is not zero
00458   {
00459     m_matT.coeffRef(iu-1, iu-2) = beta;
00460     m_matT.block(iu-1, iu-1, 2, size-iu+1).applyHouseholderOnTheLeft(ess, tau, workspace);
00461     m_matT.block(0, iu-1, iu+1, 2).applyHouseholderOnTheRight(ess, tau, workspace);
00462     if (computeU)
00463       m_matU.block(0, iu-1, size, 2).applyHouseholderOnTheRight(ess, tau, workspace);
00464   }
00465 
00466   // clean up pollution due to round-off errors
00467   for (Index i = im+2; i <= iu; ++i)
00468   {
00469     m_matT.coeffRef(i,i-2) = Scalar(0);
00470     if (i > im+2)
00471       m_matT.coeffRef(i,i-3) = Scalar(0);
00472   }
00473 }
00474 
00475 } // end namespace Eigen
00476 
00477 #endif // EIGEN_REAL_SCHUR_H