ComplexSchur.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) 2009 Claire Maurice
00005 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
00006 // Copyright (C) 2010 Jitse Niesen <jitse@maths.leeds.ac.uk>
00007 //
00008 // This Source Code Form is subject to the terms of the Mozilla
00009 // Public License v. 2.0. If a copy of the MPL was not distributed
00010 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00011 
00012 #ifndef EIGEN_COMPLEX_SCHUR_H
00013 #define EIGEN_COMPLEX_SCHUR_H
00014 
00015 #include "./HessenbergDecomposition.h"
00016 
00017 namespace Eigen { 
00018 
00019 namespace internal {
00020 template<typename MatrixType, bool IsComplex> struct complex_schur_reduce_to_hessenberg;
00021 }
00022 
00051 template<typename _MatrixType> class ComplexSchur
00052 {
00053   public:
00054     typedef _MatrixType MatrixType;
00055     enum {
00056       RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00057       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00058       Options = MatrixType::Options,
00059       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00060       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
00061     };
00062 
00064     typedef typename MatrixType::Scalar Scalar;
00065     typedef typename NumTraits<Scalar>::Real RealScalar;
00066     typedef typename MatrixType::Index Index;
00067 
00074     typedef std::complex<RealScalar> ComplexScalar;
00075 
00081     typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> ComplexMatrixType;
00082 
00094     ComplexSchur(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime)
00095       : m_matT(size,size),
00096         m_matU(size,size),
00097         m_hess(size),
00098         m_isInitialized(false),
00099         m_matUisUptodate(false)
00100     {}
00101 
00111     ComplexSchur(const MatrixType& matrix, bool computeU = true)
00112             : m_matT(matrix.rows(),matrix.cols()),
00113               m_matU(matrix.rows(),matrix.cols()),
00114               m_hess(matrix.rows()),
00115               m_isInitialized(false),
00116               m_matUisUptodate(false)
00117     {
00118       compute(matrix, computeU);
00119     }
00120 
00135     const ComplexMatrixType& matrixU() const
00136     {
00137       eigen_assert(m_isInitialized && "ComplexSchur is not initialized.");
00138       eigen_assert(m_matUisUptodate && "The matrix U has not been computed during the ComplexSchur decomposition.");
00139       return m_matU;
00140     }
00141 
00159     const ComplexMatrixType& matrixT() const
00160     {
00161       eigen_assert(m_isInitialized && "ComplexSchur is not initialized.");
00162       return m_matT;
00163     }
00164 
00184     ComplexSchur& compute(const MatrixType& matrix, bool computeU = true);
00185 
00190     ComputationInfo info() const
00191     {
00192       eigen_assert(m_isInitialized && "RealSchur is not initialized.");
00193       return m_info;
00194     }
00195 
00200     static const int m_maxIterations = 30;
00201 
00202   protected:
00203     ComplexMatrixType m_matT, m_matU;
00204     HessenbergDecomposition<MatrixType> m_hess;
00205     ComputationInfo m_info;
00206     bool m_isInitialized;
00207     bool m_matUisUptodate;
00208 
00209   private:  
00210     bool subdiagonalEntryIsNeglegible(Index i);
00211     ComplexScalar computeShift(Index iu, Index iter);
00212     void reduceToTriangularForm(bool computeU);
00213     friend struct internal::complex_schur_reduce_to_hessenberg<MatrixType, NumTraits<Scalar>::IsComplex>;
00214 };
00215 
00219 template<typename MatrixType>
00220 inline bool ComplexSchur<MatrixType>::subdiagonalEntryIsNeglegible(Index i)
00221 {
00222   RealScalar d = internal::norm1(m_matT.coeff(i,i)) + internal::norm1(m_matT.coeff(i+1,i+1));
00223   RealScalar sd = internal::norm1(m_matT.coeff(i+1,i));
00224   if (internal::isMuchSmallerThan(sd, d, NumTraits<RealScalar>::epsilon()))
00225   {
00226     m_matT.coeffRef(i+1,i) = ComplexScalar(0);
00227     return true;
00228   }
00229   return false;
00230 }
00231 
00232 
00234 template<typename MatrixType>
00235 typename ComplexSchur<MatrixType>::ComplexScalar ComplexSchur<MatrixType>::computeShift(Index iu, Index iter)
00236 {
00237   if (iter == 10 || iter == 20) 
00238   {
00239     // exceptional shift, taken from http://www.netlib.org/eispack/comqr.f
00240     return internal::abs(internal::real(m_matT.coeff(iu,iu-1))) + internal::abs(internal::real(m_matT.coeff(iu-1,iu-2)));
00241   }
00242 
00243   // compute the shift as one of the eigenvalues of t, the 2x2
00244   // diagonal block on the bottom of the active submatrix
00245   Matrix<ComplexScalar,2,2> t = m_matT.template block<2,2>(iu-1,iu-1);
00246   RealScalar normt = t.cwiseAbs().sum();
00247   t /= normt;     // the normalization by sf is to avoid under/overflow
00248 
00249   ComplexScalar b = t.coeff(0,1) * t.coeff(1,0);
00250   ComplexScalar c = t.coeff(0,0) - t.coeff(1,1);
00251   ComplexScalar disc = sqrt(c*c + RealScalar(4)*b);
00252   ComplexScalar det = t.coeff(0,0) * t.coeff(1,1) - b;
00253   ComplexScalar trace = t.coeff(0,0) + t.coeff(1,1);
00254   ComplexScalar eival1 = (trace + disc) / RealScalar(2);
00255   ComplexScalar eival2 = (trace - disc) / RealScalar(2);
00256 
00257   if(internal::norm1(eival1) > internal::norm1(eival2))
00258     eival2 = det / eival1;
00259   else
00260     eival1 = det / eival2;
00261 
00262   // choose the eigenvalue closest to the bottom entry of the diagonal
00263   if(internal::norm1(eival1-t.coeff(1,1)) < internal::norm1(eival2-t.coeff(1,1)))
00264     return normt * eival1;
00265   else
00266     return normt * eival2;
00267 }
00268 
00269 
00270 template<typename MatrixType>
00271 ComplexSchur<MatrixType>& ComplexSchur<MatrixType>::compute(const MatrixType& matrix, bool computeU)
00272 {
00273   m_matUisUptodate = false;
00274   eigen_assert(matrix.cols() == matrix.rows());
00275 
00276   if(matrix.cols() == 1)
00277   {
00278     m_matT = matrix.template cast<ComplexScalar>();
00279     if(computeU)  m_matU = ComplexMatrixType::Identity(1,1);
00280     m_info = Success;
00281     m_isInitialized = true;
00282     m_matUisUptodate = computeU;
00283     return *this;
00284   }
00285 
00286   internal::complex_schur_reduce_to_hessenberg<MatrixType, NumTraits<Scalar>::IsComplex>::run(*this, matrix, computeU);
00287   reduceToTriangularForm(computeU);
00288   return *this;
00289 }
00290 
00291 namespace internal {
00292 
00293 /* Reduce given matrix to Hessenberg form */
00294 template<typename MatrixType, bool IsComplex>
00295 struct complex_schur_reduce_to_hessenberg
00296 {
00297   // this is the implementation for the case IsComplex = true
00298   static void run(ComplexSchur<MatrixType>& _this, const MatrixType& matrix, bool computeU)
00299   {
00300     _this.m_hess.compute(matrix);
00301     _this.m_matT = _this.m_hess.matrixH();
00302     if(computeU)  _this.m_matU = _this.m_hess.matrixQ();
00303   }
00304 };
00305 
00306 template<typename MatrixType>
00307 struct complex_schur_reduce_to_hessenberg<MatrixType, false>
00308 {
00309   static void run(ComplexSchur<MatrixType>& _this, const MatrixType& matrix, bool computeU)
00310   {
00311     typedef typename ComplexSchur<MatrixType>::ComplexScalar ComplexScalar;
00312     typedef typename ComplexSchur<MatrixType>::ComplexMatrixType ComplexMatrixType;
00313 
00314     // Note: m_hess is over RealScalar; m_matT and m_matU is over ComplexScalar
00315     _this.m_hess.compute(matrix);
00316     _this.m_matT = _this.m_hess.matrixH().template cast<ComplexScalar>();
00317     if(computeU)  
00318     {
00319       // This may cause an allocation which seems to be avoidable
00320       MatrixType Q = _this.m_hess.matrixQ(); 
00321       _this.m_matU = Q.template cast<ComplexScalar>();
00322     }
00323   }
00324 };
00325 
00326 } // end namespace internal
00327 
00328 // Reduce the Hessenberg matrix m_matT to triangular form by QR iteration.
00329 template<typename MatrixType>
00330 void ComplexSchur<MatrixType>::reduceToTriangularForm(bool computeU)
00331 {  
00332   // The matrix m_matT is divided in three parts. 
00333   // Rows 0,...,il-1 are decoupled from the rest because m_matT(il,il-1) is zero. 
00334   // Rows il,...,iu is the part we are working on (the active submatrix).
00335   // Rows iu+1,...,end are already brought in triangular form.
00336   Index iu = m_matT.cols() - 1;
00337   Index il;
00338   Index iter = 0; // number of iterations we are working on the (iu,iu) element
00339   Index totalIter = 0; // number of iterations for whole matrix
00340 
00341   while(true)
00342   {
00343     // find iu, the bottom row of the active submatrix
00344     while(iu > 0)
00345     {
00346       if(!subdiagonalEntryIsNeglegible(iu-1)) break;
00347       iter = 0;
00348       --iu;
00349     }
00350 
00351     // if iu is zero then we are done; the whole matrix is triangularized
00352     if(iu==0) break;
00353 
00354     // if we spent too many iterations, we give up
00355     iter++;
00356     totalIter++;
00357     if(totalIter > m_maxIterations * m_matT.cols()) break;
00358 
00359     // find il, the top row of the active submatrix
00360     il = iu-1;
00361     while(il > 0 && !subdiagonalEntryIsNeglegible(il-1))
00362     {
00363       --il;
00364     }
00365 
00366     /* perform the QR step using Givens rotations. The first rotation
00367        creates a bulge; the (il+2,il) element becomes nonzero. This
00368        bulge is chased down to the bottom of the active submatrix. */
00369 
00370     ComplexScalar shift = computeShift(iu, iter);
00371     JacobiRotation<ComplexScalar> rot;
00372     rot.makeGivens(m_matT.coeff(il,il) - shift, m_matT.coeff(il+1,il));
00373     m_matT.rightCols(m_matT.cols()-il).applyOnTheLeft(il, il+1, rot.adjoint());
00374     m_matT.topRows((std::min)(il+2,iu)+1).applyOnTheRight(il, il+1, rot);
00375     if(computeU) m_matU.applyOnTheRight(il, il+1, rot);
00376 
00377     for(Index i=il+1 ; i<iu ; i++)
00378     {
00379       rot.makeGivens(m_matT.coeffRef(i,i-1), m_matT.coeffRef(i+1,i-1), &m_matT.coeffRef(i,i-1));
00380       m_matT.coeffRef(i+1,i-1) = ComplexScalar(0);
00381       m_matT.rightCols(m_matT.cols()-i).applyOnTheLeft(i, i+1, rot.adjoint());
00382       m_matT.topRows((std::min)(i+2,iu)+1).applyOnTheRight(i, i+1, rot);
00383       if(computeU) m_matU.applyOnTheRight(i, i+1, rot);
00384     }
00385   }
00386 
00387   if(totalIter <= m_maxIterations * m_matT.cols())
00388     m_info = Success;
00389   else
00390     m_info = NoConvergence;
00391 
00392   m_isInitialized = true;
00393   m_matUisUptodate = computeU;
00394 }
00395 
00396 } // end namespace Eigen
00397 
00398 #endif // EIGEN_COMPLEX_SCHUR_H


win_eigen
Author(s): Daniel Stonier
autogenerated on Wed Sep 16 2015 07:10:27