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,2012 Jitse Niesen <jitse@maths.leeds.ac.uk>
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 typename MatrixType::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     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     RealSchur(const MatrixType& matrix, bool computeU = true)
00104             : m_matT(matrix.rows(),matrix.cols()),
00105               m_matU(matrix.rows(),matrix.cols()),
00106               m_workspaceVector(matrix.rows()),
00107               m_hess(matrix.rows()),
00108               m_isInitialized(false),
00109               m_matUisUptodate(false),
00110               m_maxIters(-1)
00111     {
00112       compute(matrix, computeU);
00113     }
00114 
00126     const MatrixType& matrixU() const
00127     {
00128       eigen_assert(m_isInitialized && "RealSchur is not initialized.");
00129       eigen_assert(m_matUisUptodate && "The matrix U has not been computed during the RealSchur decomposition.");
00130       return m_matU;
00131     }
00132 
00143     const MatrixType& matrixT() const
00144     {
00145       eigen_assert(m_isInitialized && "RealSchur is not initialized.");
00146       return m_matT;
00147     }
00148   
00168     RealSchur& compute(const MatrixType& matrix, bool computeU = true);
00169 
00187     template<typename HessMatrixType, typename OrthMatrixType>
00188     RealSchur& computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ,  bool computeU);
00193     ComputationInfo info() const
00194     {
00195       eigen_assert(m_isInitialized && "RealSchur is not initialized.");
00196       return m_info;
00197     }
00198 
00204     RealSchur& setMaxIterations(Index maxIters)
00205     {
00206       m_maxIters = maxIters;
00207       return *this;
00208     }
00209 
00211     Index getMaxIterations()
00212     {
00213       return m_maxIters;
00214     }
00215 
00221     static const int m_maxIterationsPerRow = 40;
00222 
00223   private:
00224     
00225     MatrixType m_matT;
00226     MatrixType m_matU;
00227     ColumnVectorType m_workspaceVector;
00228     HessenbergDecomposition<MatrixType> m_hess;
00229     ComputationInfo m_info;
00230     bool m_isInitialized;
00231     bool m_matUisUptodate;
00232     Index m_maxIters;
00233 
00234     typedef Matrix<Scalar,3,1> Vector3s;
00235 
00236     Scalar computeNormOfT();
00237     Index findSmallSubdiagEntry(Index iu, const Scalar& norm);
00238     void splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift);
00239     void computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo);
00240     void initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector);
00241     void performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace);
00242 };
00243 
00244 
00245 template<typename MatrixType>
00246 RealSchur<MatrixType>& RealSchur<MatrixType>::compute(const MatrixType& matrix, bool computeU)
00247 {
00248   eigen_assert(matrix.cols() == matrix.rows());
00249   Index maxIters = m_maxIters;
00250   if (maxIters == -1)
00251     maxIters = m_maxIterationsPerRow * matrix.rows();
00252 
00253   // Step 1. Reduce to Hessenberg form
00254   m_hess.compute(matrix);
00255 
00256   // Step 2. Reduce to real Schur form  
00257   computeFromHessenberg(m_hess.matrixH(), m_hess.matrixQ(), computeU);
00258   
00259   return *this;
00260 }
00261 template<typename MatrixType>
00262 template<typename HessMatrixType, typename OrthMatrixType>
00263 RealSchur<MatrixType>& RealSchur<MatrixType>::computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ,  bool computeU)
00264 {  
00265   m_matT = matrixH; 
00266   if(computeU)
00267     m_matU = matrixQ;
00268   
00269   Index maxIters = m_maxIters;
00270   if (maxIters == -1)
00271     maxIters = m_maxIterationsPerRow * matrixH.rows();
00272   m_workspaceVector.resize(m_matT.cols());
00273   Scalar* workspace = &m_workspaceVector.coeffRef(0);
00274 
00275   // The matrix m_matT is divided in three parts. 
00276   // Rows 0,...,il-1 are decoupled from the rest because m_matT(il,il-1) is zero. 
00277   // Rows il,...,iu is the part we are working on (the active window).
00278   // Rows iu+1,...,end are already brought in triangular form.
00279   Index iu = m_matT.cols() - 1;
00280   Index iter = 0;      // iteration count for current eigenvalue
00281   Index totalIter = 0; // iteration count for whole matrix
00282   Scalar exshift(0);   // sum of exceptional shifts
00283   Scalar norm = computeNormOfT();
00284 
00285   if(norm!=0)
00286   {
00287     while (iu >= 0)
00288     {
00289       Index il = findSmallSubdiagEntry(iu, norm);
00290 
00291       // Check for convergence
00292       if (il == iu) // One root found
00293       {
00294         m_matT.coeffRef(iu,iu) = m_matT.coeff(iu,iu) + exshift;
00295         if (iu > 0)
00296           m_matT.coeffRef(iu, iu-1) = Scalar(0);
00297         iu--;
00298         iter = 0;
00299       }
00300       else if (il == iu-1) // Two roots found
00301       {
00302         splitOffTwoRows(iu, computeU, exshift);
00303         iu -= 2;
00304         iter = 0;
00305       }
00306       else // No convergence yet
00307       {
00308         // The firstHouseholderVector vector has to be initialized to something to get rid of a silly GCC warning (-O1 -Wall -DNDEBUG )
00309         Vector3s firstHouseholderVector(0,0,0), shiftInfo;
00310         computeShift(iu, iter, exshift, shiftInfo);
00311         iter = iter + 1;
00312         totalIter = totalIter + 1;
00313         if (totalIter > maxIters) break;
00314         Index im;
00315         initFrancisQRStep(il, iu, shiftInfo, im, firstHouseholderVector);
00316         performFrancisQRStep(il, im, iu, computeU, firstHouseholderVector, workspace);
00317       }
00318     }
00319   }
00320   if(totalIter <= maxIters)
00321     m_info = Success;
00322   else
00323     m_info = NoConvergence;
00324 
00325   m_isInitialized = true;
00326   m_matUisUptodate = computeU;
00327   return *this;
00328 }
00329 
00331 template<typename MatrixType>
00332 inline typename MatrixType::Scalar RealSchur<MatrixType>::computeNormOfT()
00333 {
00334   const Index size = m_matT.cols();
00335   // FIXME to be efficient the following would requires a triangular reduxion code
00336   // Scalar norm = m_matT.upper().cwiseAbs().sum() 
00337   //               + m_matT.bottomLeftCorner(size-1,size-1).diagonal().cwiseAbs().sum();
00338   Scalar norm(0);
00339   for (Index j = 0; j < size; ++j)
00340     norm += m_matT.col(j).segment(0, (std::min)(size,j+2)).cwiseAbs().sum();
00341   return norm;
00342 }
00343 
00345 template<typename MatrixType>
00346 inline typename MatrixType::Index RealSchur<MatrixType>::findSmallSubdiagEntry(Index iu, const Scalar& norm)
00347 {
00348   using std::abs;
00349   Index res = iu;
00350   while (res > 0)
00351   {
00352     Scalar s = abs(m_matT.coeff(res-1,res-1)) + abs(m_matT.coeff(res,res));
00353     if (s == 0.0)
00354       s = norm;
00355     if (abs(m_matT.coeff(res,res-1)) < NumTraits<Scalar>::epsilon() * s)
00356       break;
00357     res--;
00358   }
00359   return res;
00360 }
00361 
00363 template<typename MatrixType>
00364 inline void RealSchur<MatrixType>::splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift)
00365 {
00366   using std::sqrt;
00367   using std::abs;
00368   const Index size = m_matT.cols();
00369 
00370   // The eigenvalues of the 2x2 matrix [a b; c d] are 
00371   // trace +/- sqrt(discr/4) where discr = tr^2 - 4*det, tr = a + d, det = ad - bc
00372   Scalar p = Scalar(0.5) * (m_matT.coeff(iu-1,iu-1) - m_matT.coeff(iu,iu));
00373   Scalar q = p * p + m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);   // q = tr^2 / 4 - det = discr/4
00374   m_matT.coeffRef(iu,iu) += exshift;
00375   m_matT.coeffRef(iu-1,iu-1) += exshift;
00376 
00377   if (q >= Scalar(0)) // Two real eigenvalues
00378   {
00379     Scalar z = sqrt(abs(q));
00380     JacobiRotation<Scalar> rot;
00381     if (p >= Scalar(0))
00382       rot.makeGivens(p + z, m_matT.coeff(iu, iu-1));
00383     else
00384       rot.makeGivens(p - z, m_matT.coeff(iu, iu-1));
00385 
00386     m_matT.rightCols(size-iu+1).applyOnTheLeft(iu-1, iu, rot.adjoint());
00387     m_matT.topRows(iu+1).applyOnTheRight(iu-1, iu, rot);
00388     m_matT.coeffRef(iu, iu-1) = Scalar(0); 
00389     if (computeU)
00390       m_matU.applyOnTheRight(iu-1, iu, rot);
00391   }
00392 
00393   if (iu > 1) 
00394     m_matT.coeffRef(iu-1, iu-2) = Scalar(0);
00395 }
00396 
00398 template<typename MatrixType>
00399 inline void RealSchur<MatrixType>::computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo)
00400 {
00401   using std::sqrt;
00402   using std::abs;
00403   shiftInfo.coeffRef(0) = m_matT.coeff(iu,iu);
00404   shiftInfo.coeffRef(1) = m_matT.coeff(iu-1,iu-1);
00405   shiftInfo.coeffRef(2) = m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
00406 
00407   // Wilkinson's original ad hoc shift
00408   if (iter == 10)
00409   {
00410     exshift += shiftInfo.coeff(0);
00411     for (Index i = 0; i <= iu; ++i)
00412       m_matT.coeffRef(i,i) -= shiftInfo.coeff(0);
00413     Scalar s = abs(m_matT.coeff(iu,iu-1)) + abs(m_matT.coeff(iu-1,iu-2));
00414     shiftInfo.coeffRef(0) = Scalar(0.75) * s;
00415     shiftInfo.coeffRef(1) = Scalar(0.75) * s;
00416     shiftInfo.coeffRef(2) = Scalar(-0.4375) * s * s;
00417   }
00418 
00419   // MATLAB's new ad hoc shift
00420   if (iter == 30)
00421   {
00422     Scalar s = (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
00423     s = s * s + shiftInfo.coeff(2);
00424     if (s > Scalar(0))
00425     {
00426       s = sqrt(s);
00427       if (shiftInfo.coeff(1) < shiftInfo.coeff(0))
00428         s = -s;
00429       s = s + (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
00430       s = shiftInfo.coeff(0) - shiftInfo.coeff(2) / s;
00431       exshift += s;
00432       for (Index i = 0; i <= iu; ++i)
00433         m_matT.coeffRef(i,i) -= s;
00434       shiftInfo.setConstant(Scalar(0.964));
00435     }
00436   }
00437 }
00438 
00440 template<typename MatrixType>
00441 inline void RealSchur<MatrixType>::initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector)
00442 {
00443   using std::abs;
00444   Vector3s& v = firstHouseholderVector; // alias to save typing
00445 
00446   for (im = iu-2; im >= il; --im)
00447   {
00448     const Scalar Tmm = m_matT.coeff(im,im);
00449     const Scalar r = shiftInfo.coeff(0) - Tmm;
00450     const Scalar s = shiftInfo.coeff(1) - Tmm;
00451     v.coeffRef(0) = (r * s - shiftInfo.coeff(2)) / m_matT.coeff(im+1,im) + m_matT.coeff(im,im+1);
00452     v.coeffRef(1) = m_matT.coeff(im+1,im+1) - Tmm - r - s;
00453     v.coeffRef(2) = m_matT.coeff(im+2,im+1);
00454     if (im == il) {
00455       break;
00456     }
00457     const Scalar lhs = m_matT.coeff(im,im-1) * (abs(v.coeff(1)) + abs(v.coeff(2)));
00458     const Scalar rhs = v.coeff(0) * (abs(m_matT.coeff(im-1,im-1)) + abs(Tmm) + abs(m_matT.coeff(im+1,im+1)));
00459     if (abs(lhs) < NumTraits<Scalar>::epsilon() * rhs)
00460     {
00461       break;
00462     }
00463   }
00464 }
00465 
00467 template<typename MatrixType>
00468 inline void RealSchur<MatrixType>::performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace)
00469 {
00470   eigen_assert(im >= il);
00471   eigen_assert(im <= iu-2);
00472 
00473   const Index size = m_matT.cols();
00474 
00475   for (Index k = im; k <= iu-2; ++k)
00476   {
00477     bool firstIteration = (k == im);
00478 
00479     Vector3s v;
00480     if (firstIteration)
00481       v = firstHouseholderVector;
00482     else
00483       v = m_matT.template block<3,1>(k,k-1);
00484 
00485     Scalar tau, beta;
00486     Matrix<Scalar, 2, 1> ess;
00487     v.makeHouseholder(ess, tau, beta);
00488     
00489     if (beta != Scalar(0)) // if v is not zero
00490     {
00491       if (firstIteration && k > il)
00492         m_matT.coeffRef(k,k-1) = -m_matT.coeff(k,k-1);
00493       else if (!firstIteration)
00494         m_matT.coeffRef(k,k-1) = beta;
00495 
00496       // These Householder transformations form the O(n^3) part of the algorithm
00497       m_matT.block(k, k, 3, size-k).applyHouseholderOnTheLeft(ess, tau, workspace);
00498       m_matT.block(0, k, (std::min)(iu,k+3) + 1, 3).applyHouseholderOnTheRight(ess, tau, workspace);
00499       if (computeU)
00500         m_matU.block(0, k, size, 3).applyHouseholderOnTheRight(ess, tau, workspace);
00501     }
00502   }
00503 
00504   Matrix<Scalar, 2, 1> v = m_matT.template block<2,1>(iu-1, iu-2);
00505   Scalar tau, beta;
00506   Matrix<Scalar, 1, 1> ess;
00507   v.makeHouseholder(ess, tau, beta);
00508 
00509   if (beta != Scalar(0)) // if v is not zero
00510   {
00511     m_matT.coeffRef(iu-1, iu-2) = beta;
00512     m_matT.block(iu-1, iu-1, 2, size-iu+1).applyHouseholderOnTheLeft(ess, tau, workspace);
00513     m_matT.block(0, iu-1, iu+1, 2).applyHouseholderOnTheRight(ess, tau, workspace);
00514     if (computeU)
00515       m_matU.block(0, iu-1, size, 2).applyHouseholderOnTheRight(ess, tau, workspace);
00516   }
00517 
00518   // clean up pollution due to round-off errors
00519   for (Index i = im+2; i <= iu; ++i)
00520   {
00521     m_matT.coeffRef(i,i-2) = Scalar(0);
00522     if (i > im+2)
00523       m_matT.coeffRef(i,i-3) = Scalar(0);
00524   }
00525 }
00526 
00527 } // end namespace Eigen
00528 
00529 #endif // EIGEN_REAL_SCHUR_H


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Thu Aug 27 2015 11:59:53