MOAB  4.9.3pre
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 <[email protected]>
00005 // Copyright (C) 2010,2012 Jitse Niesen <[email protected]>
00006 //
00007 // This Source Code Form is subject to the terms of the Mozilla
00008 // Public License v. 2.0. If a copy of the MPL was not distributed
00009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00010 
00011 #ifndef EIGEN_REAL_SCHUR_H
00012 #define EIGEN_REAL_SCHUR_H
00013 
00014 #include "./HessenbergDecomposition.h"
00015 
00016 namespace Eigen { 
00017 
00054 template<typename _MatrixType> class RealSchur
00055 {
00056   public:
00057     typedef _MatrixType MatrixType;
00058     enum {
00059       RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00060       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00061       Options = MatrixType::Options,
00062       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00063       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
00064     };
00065     typedef typename MatrixType::Scalar Scalar;
00066     typedef std::complex<typename NumTraits<Scalar>::Real> ComplexScalar;
00067     typedef Eigen::Index Index; 
00068 
00069     typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> EigenvalueType;
00070     typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
00071 
00083     explicit RealSchur(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime)
00084             : m_matT(size, size),
00085               m_matU(size, size),
00086               m_workspaceVector(size),
00087               m_hess(size),
00088               m_isInitialized(false),
00089               m_matUisUptodate(false),
00090               m_maxIters(-1)
00091     { }
00092 
00103     template<typename InputType>
00104     explicit RealSchur(const EigenBase<InputType>& matrix, bool computeU = true)
00105             : m_matT(matrix.rows(),matrix.cols()),
00106               m_matU(matrix.rows(),matrix.cols()),
00107               m_workspaceVector(matrix.rows()),
00108               m_hess(matrix.rows()),
00109               m_isInitialized(false),
00110               m_matUisUptodate(false),
00111               m_maxIters(-1)
00112     {
00113       compute(matrix.derived(), computeU);
00114     }
00115 
00127     const MatrixType& matrixU() const
00128     {
00129       eigen_assert(m_isInitialized && "RealSchur is not initialized.");
00130       eigen_assert(m_matUisUptodate && "The matrix U has not been computed during the RealSchur decomposition.");
00131       return m_matU;
00132     }
00133 
00144     const MatrixType& matrixT() const
00145     {
00146       eigen_assert(m_isInitialized && "RealSchur is not initialized.");
00147       return m_matT;
00148     }
00149   
00169     template<typename InputType>
00170     RealSchur& compute(const EigenBase<InputType>& matrix, bool computeU = true);
00171 
00189     template<typename HessMatrixType, typename OrthMatrixType>
00190     RealSchur& computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ,  bool computeU);
00195     ComputationInfo info() const
00196     {
00197       eigen_assert(m_isInitialized && "RealSchur is not initialized.");
00198       return m_info;
00199     }
00200 
00206     RealSchur& setMaxIterations(Index maxIters)
00207     {
00208       m_maxIters = maxIters;
00209       return *this;
00210     }
00211 
00213     Index getMaxIterations()
00214     {
00215       return m_maxIters;
00216     }
00217 
00223     static const int m_maxIterationsPerRow = 40;
00224 
00225   private:
00226     
00227     MatrixType m_matT;
00228     MatrixType m_matU;
00229     ColumnVectorType m_workspaceVector;
00230     HessenbergDecomposition<MatrixType> m_hess;
00231     ComputationInfo m_info;
00232     bool m_isInitialized;
00233     bool m_matUisUptodate;
00234     Index m_maxIters;
00235 
00236     typedef Matrix<Scalar,3,1> Vector3s;
00237 
00238     Scalar computeNormOfT();
00239     Index findSmallSubdiagEntry(Index iu);
00240     void splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift);
00241     void computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo);
00242     void initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector);
00243     void performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace);
00244 };
00245 
00246 
00247 template<typename MatrixType>
00248 template<typename InputType>
00249 RealSchur<MatrixType>& RealSchur<MatrixType>::compute(const EigenBase<InputType>& matrix, bool computeU)
00250 {
00251   eigen_assert(matrix.cols() == matrix.rows());
00252   Index maxIters = m_maxIters;
00253   if (maxIters == -1)
00254     maxIters = m_maxIterationsPerRow * matrix.rows();
00255 
00256   // Step 1. Reduce to Hessenberg form
00257   m_hess.compute(matrix.derived());
00258 
00259   // Step 2. Reduce to real Schur form  
00260   computeFromHessenberg(m_hess.matrixH(), m_hess.matrixQ(), computeU);
00261   
00262   return *this;
00263 }
00264 template<typename MatrixType>
00265 template<typename HessMatrixType, typename OrthMatrixType>
00266 RealSchur<MatrixType>& RealSchur<MatrixType>::computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ,  bool computeU)
00267 {  
00268   m_matT = matrixH; 
00269   if(computeU)
00270     m_matU = matrixQ;
00271   
00272   Index maxIters = m_maxIters;
00273   if (maxIters == -1)
00274     maxIters = m_maxIterationsPerRow * matrixH.rows();
00275   m_workspaceVector.resize(m_matT.cols());
00276   Scalar* workspace = &m_workspaceVector.coeffRef(0);
00277 
00278   // The matrix m_matT is divided in three parts. 
00279   // Rows 0,...,il-1 are decoupled from the rest because m_matT(il,il-1) is zero. 
00280   // Rows il,...,iu is the part we are working on (the active window).
00281   // Rows iu+1,...,end are already brought in triangular form.
00282   Index iu = m_matT.cols() - 1;
00283   Index iter = 0;      // iteration count for current eigenvalue
00284   Index totalIter = 0; // iteration count for whole matrix
00285   Scalar exshift(0);   // sum of exceptional shifts
00286   Scalar norm = computeNormOfT();
00287 
00288   if(norm!=0)
00289   {
00290     while (iu >= 0)
00291     {
00292       Index il = findSmallSubdiagEntry(iu);
00293 
00294       // Check for convergence
00295       if (il == iu) // One root found
00296       {
00297         m_matT.coeffRef(iu,iu) = m_matT.coeff(iu,iu) + exshift;
00298         if (iu > 0)
00299           m_matT.coeffRef(iu, iu-1) = Scalar(0);
00300         iu--;
00301         iter = 0;
00302       }
00303       else if (il == iu-1) // Two roots found
00304       {
00305         splitOffTwoRows(iu, computeU, exshift);
00306         iu -= 2;
00307         iter = 0;
00308       }
00309       else // No convergence yet
00310       {
00311         // The firstHouseholderVector vector has to be initialized to something to get rid of a silly GCC warning (-O1 -Wall -DNDEBUG )
00312         Vector3s firstHouseholderVector(0,0,0), shiftInfo;
00313         computeShift(iu, iter, exshift, shiftInfo);
00314         iter = iter + 1;
00315         totalIter = totalIter + 1;
00316         if (totalIter > maxIters) break;
00317         Index im;
00318         initFrancisQRStep(il, iu, shiftInfo, im, firstHouseholderVector);
00319         performFrancisQRStep(il, im, iu, computeU, firstHouseholderVector, workspace);
00320       }
00321     }
00322   }
00323   if(totalIter <= maxIters)
00324     m_info = Success;
00325   else
00326     m_info = NoConvergence;
00327 
00328   m_isInitialized = true;
00329   m_matUisUptodate = computeU;
00330   return *this;
00331 }
00332 
00334 template<typename MatrixType>
00335 inline typename MatrixType::Scalar RealSchur<MatrixType>::computeNormOfT()
00336 {
00337   const Index size = m_matT.cols();
00338   // FIXME to be efficient the following would requires a triangular reduxion code
00339   // Scalar norm = m_matT.upper().cwiseAbs().sum() 
00340   //               + m_matT.bottomLeftCorner(size-1,size-1).diagonal().cwiseAbs().sum();
00341   Scalar norm(0);
00342   for (Index j = 0; j < size; ++j)
00343     norm += m_matT.col(j).segment(0, (std::min)(size,j+2)).cwiseAbs().sum();
00344   return norm;
00345 }
00346 
00348 template<typename MatrixType>
00349 inline Index RealSchur<MatrixType>::findSmallSubdiagEntry(Index iu)
00350 {
00351   using std::abs;
00352   Index res = iu;
00353   while (res > 0)
00354   {
00355     Scalar s = abs(m_matT.coeff(res-1,res-1)) + abs(m_matT.coeff(res,res));
00356     if (abs(m_matT.coeff(res,res-1)) <= NumTraits<Scalar>::epsilon() * s)
00357       break;
00358     res--;
00359   }
00360   return res;
00361 }
00362 
00364 template<typename MatrixType>
00365 inline void RealSchur<MatrixType>::splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift)
00366 {
00367   using std::sqrt;
00368   using std::abs;
00369   const Index size = m_matT.cols();
00370 
00371   // The eigenvalues of the 2x2 matrix [a b; c d] are 
00372   // trace +/- sqrt(discr/4) where discr = tr^2 - 4*det, tr = a + d, det = ad - bc
00373   Scalar p = Scalar(0.5) * (m_matT.coeff(iu-1,iu-1) - m_matT.coeff(iu,iu));
00374   Scalar q = p * p + m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);   // q = tr^2 / 4 - det = discr/4
00375   m_matT.coeffRef(iu,iu) += exshift;
00376   m_matT.coeffRef(iu-1,iu-1) += exshift;
00377 
00378   if (q >= Scalar(0)) // Two real eigenvalues
00379   {
00380     Scalar z = sqrt(abs(q));
00381     JacobiRotation<Scalar> rot;
00382     if (p >= Scalar(0))
00383       rot.makeGivens(p + z, m_matT.coeff(iu, iu-1));
00384     else
00385       rot.makeGivens(p - z, m_matT.coeff(iu, iu-1));
00386 
00387     m_matT.rightCols(size-iu+1).applyOnTheLeft(iu-1, iu, rot.adjoint());
00388     m_matT.topRows(iu+1).applyOnTheRight(iu-1, iu, rot);
00389     m_matT.coeffRef(iu, iu-1) = Scalar(0); 
00390     if (computeU)
00391       m_matU.applyOnTheRight(iu-1, iu, rot);
00392   }
00393 
00394   if (iu > 1) 
00395     m_matT.coeffRef(iu-1, iu-2) = Scalar(0);
00396 }
00397 
00399 template<typename MatrixType>
00400 inline void RealSchur<MatrixType>::computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo)
00401 {
00402   using std::sqrt;
00403   using std::abs;
00404   shiftInfo.coeffRef(0) = m_matT.coeff(iu,iu);
00405   shiftInfo.coeffRef(1) = m_matT.coeff(iu-1,iu-1);
00406   shiftInfo.coeffRef(2) = m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
00407 
00408   // Wilkinson's original ad hoc shift
00409   if (iter == 10)
00410   {
00411     exshift += shiftInfo.coeff(0);
00412     for (Index i = 0; i <= iu; ++i)
00413       m_matT.coeffRef(i,i) -= shiftInfo.coeff(0);
00414     Scalar s = abs(m_matT.coeff(iu,iu-1)) + abs(m_matT.coeff(iu-1,iu-2));
00415     shiftInfo.coeffRef(0) = Scalar(0.75) * s;
00416     shiftInfo.coeffRef(1) = Scalar(0.75) * s;
00417     shiftInfo.coeffRef(2) = Scalar(-0.4375) * s * s;
00418   }
00419 
00420   // MATLAB's new ad hoc shift
00421   if (iter == 30)
00422   {
00423     Scalar s = (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
00424     s = s * s + shiftInfo.coeff(2);
00425     if (s > Scalar(0))
00426     {
00427       s = sqrt(s);
00428       if (shiftInfo.coeff(1) < shiftInfo.coeff(0))
00429         s = -s;
00430       s = s + (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
00431       s = shiftInfo.coeff(0) - shiftInfo.coeff(2) / s;
00432       exshift += s;
00433       for (Index i = 0; i <= iu; ++i)
00434         m_matT.coeffRef(i,i) -= s;
00435       shiftInfo.setConstant(Scalar(0.964));
00436     }
00437   }
00438 }
00439 
00441 template<typename MatrixType>
00442 inline void RealSchur<MatrixType>::initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector)
00443 {
00444   using std::abs;
00445   Vector3s& v = firstHouseholderVector; // alias to save typing
00446 
00447   for (im = iu-2; im >= il; --im)
00448   {
00449     const Scalar Tmm = m_matT.coeff(im,im);
00450     const Scalar r = shiftInfo.coeff(0) - Tmm;
00451     const Scalar s = shiftInfo.coeff(1) - Tmm;
00452     v.coeffRef(0) = (r * s - shiftInfo.coeff(2)) / m_matT.coeff(im+1,im) + m_matT.coeff(im,im+1);
00453     v.coeffRef(1) = m_matT.coeff(im+1,im+1) - Tmm - r - s;
00454     v.coeffRef(2) = m_matT.coeff(im+2,im+1);
00455     if (im == il) {
00456       break;
00457     }
00458     const Scalar lhs = m_matT.coeff(im,im-1) * (abs(v.coeff(1)) + abs(v.coeff(2)));
00459     const Scalar rhs = v.coeff(0) * (abs(m_matT.coeff(im-1,im-1)) + abs(Tmm) + abs(m_matT.coeff(im+1,im+1)));
00460     if (abs(lhs) < NumTraits<Scalar>::epsilon() * rhs)
00461       break;
00462   }
00463 }
00464 
00466 template<typename MatrixType>
00467 inline void RealSchur<MatrixType>::performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace)
00468 {
00469   eigen_assert(im >= il);
00470   eigen_assert(im <= iu-2);
00471 
00472   const Index size = m_matT.cols();
00473 
00474   for (Index k = im; k <= iu-2; ++k)
00475   {
00476     bool firstIteration = (k == im);
00477 
00478     Vector3s v;
00479     if (firstIteration)
00480       v = firstHouseholderVector;
00481     else
00482       v = m_matT.template block<3,1>(k,k-1);
00483 
00484     Scalar tau, beta;
00485     Matrix<Scalar, 2, 1> ess;
00486     v.makeHouseholder(ess, tau, beta);
00487     
00488     if (beta != Scalar(0)) // if v is not zero
00489     {
00490       if (firstIteration && k > il)
00491         m_matT.coeffRef(k,k-1) = -m_matT.coeff(k,k-1);
00492       else if (!firstIteration)
00493         m_matT.coeffRef(k,k-1) = beta;
00494 
00495       // These Householder transformations form the O(n^3) part of the algorithm
00496       m_matT.block(k, k, 3, size-k).applyHouseholderOnTheLeft(ess, tau, workspace);
00497       m_matT.block(0, k, (std::min)(iu,k+3) + 1, 3).applyHouseholderOnTheRight(ess, tau, workspace);
00498       if (computeU)
00499         m_matU.block(0, k, size, 3).applyHouseholderOnTheRight(ess, tau, workspace);
00500     }
00501   }
00502 
00503   Matrix<Scalar, 2, 1> v = m_matT.template block<2,1>(iu-1, iu-2);
00504   Scalar tau, beta;
00505   Matrix<Scalar, 1, 1> ess;
00506   v.makeHouseholder(ess, tau, beta);
00507 
00508   if (beta != Scalar(0)) // if v is not zero
00509   {
00510     m_matT.coeffRef(iu-1, iu-2) = beta;
00511     m_matT.block(iu-1, iu-1, 2, size-iu+1).applyHouseholderOnTheLeft(ess, tau, workspace);
00512     m_matT.block(0, iu-1, iu+1, 2).applyHouseholderOnTheRight(ess, tau, workspace);
00513     if (computeU)
00514       m_matU.block(0, iu-1, size, 2).applyHouseholderOnTheRight(ess, tau, workspace);
00515   }
00516 
00517   // clean up pollution due to round-off errors
00518   for (Index i = im+2; i <= iu; ++i)
00519   {
00520     m_matT.coeffRef(i,i-2) = Scalar(0);
00521     if (i > im+2)
00522       m_matT.coeffRef(i,i-3) = Scalar(0);
00523   }
00524 }
00525 
00526 } // end namespace Eigen
00527 
00528 #endif // EIGEN_REAL_SCHUR_H
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines