EigenSolver.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_EIGENSOLVER_H
00027 #define EIGEN_EIGENSOLVER_H
00028 
00029 #include "./EigenvaluesCommon.h"
00030 #include "./RealSchur.h"
00031 
00078 template<typename _MatrixType> class EigenSolver
00079 {
00080   public:
00081 
00083     typedef _MatrixType MatrixType;
00084 
00085     enum {
00086       RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00087       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00088       Options = MatrixType::Options,
00089       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00090       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
00091     };
00092 
00094     typedef typename MatrixType::Scalar Scalar;
00095     typedef typename NumTraits<Scalar>::Real RealScalar;
00096     typedef typename MatrixType::Index Index;
00097 
00104     typedef std::complex<RealScalar> ComplexScalar;
00105 
00111     typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> EigenvalueType;
00112 
00118     typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> EigenvectorsType;
00119 
00127  EigenSolver() : m_eivec(), m_eivalues(), m_isInitialized(false), m_realSchur(), m_matT(), m_tmp() {}
00128 
00135     EigenSolver(Index size)
00136       : m_eivec(size, size),
00137         m_eivalues(size),
00138         m_isInitialized(false),
00139         m_eigenvectorsOk(false),
00140         m_realSchur(size),
00141         m_matT(size, size), 
00142         m_tmp(size)
00143     {}
00144 
00160     EigenSolver(const MatrixType& matrix, bool computeEigenvectors = true)
00161       : m_eivec(matrix.rows(), matrix.cols()),
00162         m_eivalues(matrix.cols()),
00163         m_isInitialized(false),
00164         m_eigenvectorsOk(false),
00165         m_realSchur(matrix.cols()),
00166         m_matT(matrix.rows(), matrix.cols()), 
00167         m_tmp(matrix.cols())
00168     {
00169       compute(matrix, computeEigenvectors);
00170     }
00171 
00192     EigenvectorsType eigenvectors() const;
00193 
00212     const MatrixType& pseudoEigenvectors() const
00213     {
00214       eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
00215       eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
00216       return m_eivec;
00217     }
00218 
00237     MatrixType pseudoEigenvalueMatrix() const;
00238 
00257     const EigenvalueType& eigenvalues() const
00258     {
00259       eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
00260       return m_eivalues;
00261     }
00262 
00290     EigenSolver& compute(const MatrixType& matrix, bool computeEigenvectors = true);
00291 
00292     ComputationInfo info() const
00293     {
00294       eigen_assert(m_isInitialized && "ComplexEigenSolver is not initialized.");
00295       return m_realSchur.info();
00296     }
00297 
00298   private:
00299     void doComputeEigenvectors();
00300 
00301   protected:
00302     MatrixType m_eivec;
00303     EigenvalueType m_eivalues;
00304     bool m_isInitialized;
00305     bool m_eigenvectorsOk;
00306     RealSchur<MatrixType> m_realSchur;
00307     MatrixType m_matT;
00308 
00309     typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
00310     ColumnVectorType m_tmp;
00311 };
00312 
00313 template<typename MatrixType>
00314 MatrixType EigenSolver<MatrixType>::pseudoEigenvalueMatrix() const
00315 {
00316   eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
00317   Index n = m_eivalues.rows();
00318   MatrixType matD = MatrixType::Zero(n,n);
00319   for (Index i=0; i<n; ++i)
00320   {
00321     if (internal::isMuchSmallerThan(internal::imag(m_eivalues.coeff(i)), internal::real(m_eivalues.coeff(i))))
00322       matD.coeffRef(i,i) = internal::real(m_eivalues.coeff(i));
00323     else
00324     {
00325       matD.template block<2,2>(i,i) <<  internal::real(m_eivalues.coeff(i)), internal::imag(m_eivalues.coeff(i)),
00326                                        -internal::imag(m_eivalues.coeff(i)), internal::real(m_eivalues.coeff(i));
00327       ++i;
00328     }
00329   }
00330   return matD;
00331 }
00332 
00333 template<typename MatrixType>
00334 typename EigenSolver<MatrixType>::EigenvectorsType EigenSolver<MatrixType>::eigenvectors() const
00335 {
00336   eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
00337   eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
00338   Index n = m_eivec.cols();
00339   EigenvectorsType matV(n,n);
00340   for (Index j=0; j<n; ++j)
00341   {
00342     if (internal::isMuchSmallerThan(internal::imag(m_eivalues.coeff(j)), internal::real(m_eivalues.coeff(j))))
00343     {
00344       // we have a real eigen value
00345       matV.col(j) = m_eivec.col(j).template cast<ComplexScalar>();
00346       matV.col(j).normalize();
00347     }
00348     else
00349     {
00350       // we have a pair of complex eigen values
00351       for (Index i=0; i<n; ++i)
00352       {
00353         matV.coeffRef(i,j)   = ComplexScalar(m_eivec.coeff(i,j),  m_eivec.coeff(i,j+1));
00354         matV.coeffRef(i,j+1) = ComplexScalar(m_eivec.coeff(i,j), -m_eivec.coeff(i,j+1));
00355       }
00356       matV.col(j).normalize();
00357       matV.col(j+1).normalize();
00358       ++j;
00359     }
00360   }
00361   return matV;
00362 }
00363 
00364 template<typename MatrixType>
00365 EigenSolver<MatrixType>& EigenSolver<MatrixType>::compute(const MatrixType& matrix, bool computeEigenvectors)
00366 {
00367   assert(matrix.cols() == matrix.rows());
00368 
00369   // Reduce to real Schur form.
00370   m_realSchur.compute(matrix, computeEigenvectors);
00371   if (m_realSchur.info() == Success)
00372   {
00373     m_matT = m_realSchur.matrixT();
00374     if (computeEigenvectors)
00375       m_eivec = m_realSchur.matrixU();
00376   
00377     // Compute eigenvalues from matT
00378     m_eivalues.resize(matrix.cols());
00379     Index i = 0;
00380     while (i < matrix.cols()) 
00381     {
00382       if (i == matrix.cols() - 1 || m_matT.coeff(i+1, i) == Scalar(0)) 
00383       {
00384         m_eivalues.coeffRef(i) = m_matT.coeff(i, i);
00385         ++i;
00386       }
00387       else
00388       {
00389         Scalar p = Scalar(0.5) * (m_matT.coeff(i, i) - m_matT.coeff(i+1, i+1));
00390         Scalar z = internal::sqrt(internal::abs(p * p + m_matT.coeff(i+1, i) * m_matT.coeff(i, i+1)));
00391         m_eivalues.coeffRef(i)   = ComplexScalar(m_matT.coeff(i+1, i+1) + p, z);
00392         m_eivalues.coeffRef(i+1) = ComplexScalar(m_matT.coeff(i+1, i+1) + p, -z);
00393         i += 2;
00394       }
00395     }
00396     
00397     // Compute eigenvectors.
00398     if (computeEigenvectors)
00399       doComputeEigenvectors();
00400   }
00401 
00402   m_isInitialized = true;
00403   m_eigenvectorsOk = computeEigenvectors;
00404 
00405   return *this;
00406 }
00407 
00408 // Complex scalar division.
00409 template<typename Scalar>
00410 std::complex<Scalar> cdiv(Scalar xr, Scalar xi, Scalar yr, Scalar yi)
00411 {
00412   Scalar r,d;
00413   if (internal::abs(yr) > internal::abs(yi))
00414   {
00415       r = yi/yr;
00416       d = yr + r*yi;
00417       return std::complex<Scalar>((xr + r*xi)/d, (xi - r*xr)/d);
00418   }
00419   else
00420   {
00421       r = yr/yi;
00422       d = yi + r*yr;
00423       return std::complex<Scalar>((r*xr + xi)/d, (r*xi - xr)/d);
00424   }
00425 }
00426 
00427 
00428 template<typename MatrixType>
00429 void EigenSolver<MatrixType>::doComputeEigenvectors()
00430 {
00431   const Index size = m_eivec.cols();
00432   const Scalar eps = NumTraits<Scalar>::epsilon();
00433 
00434   // inefficient! this is already computed in RealSchur
00435   Scalar norm = 0.0;
00436   for (Index j = 0; j < size; ++j)
00437   {
00438     norm += m_matT.row(j).segment(std::max(j-1,Index(0)), size-std::max(j-1,Index(0))).cwiseAbs().sum();
00439   }
00440   
00441   // Backsubstitute to find vectors of upper triangular form
00442   if (norm == 0.0)
00443   {
00444     return;
00445   }
00446 
00447   for (Index n = size-1; n >= 0; n--)
00448   {
00449     Scalar p = m_eivalues.coeff(n).real();
00450     Scalar q = m_eivalues.coeff(n).imag();
00451 
00452     // Scalar vector
00453     if (q == Scalar(0))
00454     {
00455       Scalar lastr=0, lastw=0;
00456       Index l = n;
00457 
00458       m_matT.coeffRef(n,n) = 1.0;
00459       for (Index i = n-1; i >= 0; i--)
00460       {
00461         Scalar w = m_matT.coeff(i,i) - p;
00462         Scalar r = m_matT.row(i).segment(l,n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
00463 
00464         if (m_eivalues.coeff(i).imag() < 0.0)
00465         {
00466           lastw = w;
00467           lastr = r;
00468         }
00469         else
00470         {
00471           l = i;
00472           if (m_eivalues.coeff(i).imag() == 0.0)
00473           {
00474             if (w != 0.0)
00475               m_matT.coeffRef(i,n) = -r / w;
00476             else
00477               m_matT.coeffRef(i,n) = -r / (eps * norm);
00478           }
00479           else // Solve real equations
00480           {
00481             Scalar x = m_matT.coeff(i,i+1);
00482             Scalar y = m_matT.coeff(i+1,i);
00483             Scalar denom = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) + m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag();
00484             Scalar t = (x * lastr - lastw * r) / denom;
00485             m_matT.coeffRef(i,n) = t;
00486             if (internal::abs(x) > internal::abs(lastw))
00487               m_matT.coeffRef(i+1,n) = (-r - w * t) / x;
00488             else
00489               m_matT.coeffRef(i+1,n) = (-lastr - y * t) / lastw;
00490           }
00491 
00492           // Overflow control
00493           Scalar t = internal::abs(m_matT.coeff(i,n));
00494           if ((eps * t) * t > Scalar(1))
00495             m_matT.col(n).tail(size-i) /= t;
00496         }
00497       }
00498     }
00499     else if (q < Scalar(0) && n > 0) // Complex vector
00500     {
00501       Scalar lastra=0, lastsa=0, lastw=0;
00502       Index l = n-1;
00503 
00504       // Last vector component imaginary so matrix is triangular
00505       if (internal::abs(m_matT.coeff(n,n-1)) > internal::abs(m_matT.coeff(n-1,n)))
00506       {
00507         m_matT.coeffRef(n-1,n-1) = q / m_matT.coeff(n,n-1);
00508         m_matT.coeffRef(n-1,n) = -(m_matT.coeff(n,n) - p) / m_matT.coeff(n,n-1);
00509       }
00510       else
00511       {
00512         std::complex<Scalar> cc = cdiv<Scalar>(0.0,-m_matT.coeff(n-1,n),m_matT.coeff(n-1,n-1)-p,q);
00513         m_matT.coeffRef(n-1,n-1) = internal::real(cc);
00514         m_matT.coeffRef(n-1,n) = internal::imag(cc);
00515       }
00516       m_matT.coeffRef(n,n-1) = 0.0;
00517       m_matT.coeffRef(n,n) = 1.0;
00518       for (Index i = n-2; i >= 0; i--)
00519       {
00520         Scalar ra = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n-1).segment(l, n-l+1));
00521         Scalar sa = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
00522         Scalar w = m_matT.coeff(i,i) - p;
00523 
00524         if (m_eivalues.coeff(i).imag() < 0.0)
00525         {
00526           lastw = w;
00527           lastra = ra;
00528           lastsa = sa;
00529         }
00530         else
00531         {
00532           l = i;
00533           if (m_eivalues.coeff(i).imag() == RealScalar(0))
00534           {
00535             std::complex<Scalar> cc = cdiv(-ra,-sa,w,q);
00536             m_matT.coeffRef(i,n-1) = internal::real(cc);
00537             m_matT.coeffRef(i,n) = internal::imag(cc);
00538           }
00539           else
00540           {
00541             // Solve complex equations
00542             Scalar x = m_matT.coeff(i,i+1);
00543             Scalar y = m_matT.coeff(i+1,i);
00544             Scalar vr = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) + m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag() - q * q;
00545             Scalar vi = (m_eivalues.coeff(i).real() - p) * Scalar(2) * q;
00546             if ((vr == 0.0) && (vi == 0.0))
00547               vr = eps * norm * (internal::abs(w) + internal::abs(q) + internal::abs(x) + internal::abs(y) + internal::abs(lastw));
00548 
00549             std::complex<Scalar> cc = cdiv(x*lastra-lastw*ra+q*sa,x*lastsa-lastw*sa-q*ra,vr,vi);
00550             m_matT.coeffRef(i,n-1) = internal::real(cc);
00551             m_matT.coeffRef(i,n) = internal::imag(cc);
00552             if (internal::abs(x) > (internal::abs(lastw) + internal::abs(q)))
00553             {
00554               m_matT.coeffRef(i+1,n-1) = (-ra - w * m_matT.coeff(i,n-1) + q * m_matT.coeff(i,n)) / x;
00555               m_matT.coeffRef(i+1,n) = (-sa - w * m_matT.coeff(i,n) - q * m_matT.coeff(i,n-1)) / x;
00556             }
00557             else
00558             {
00559               cc = cdiv(-lastra-y*m_matT.coeff(i,n-1),-lastsa-y*m_matT.coeff(i,n),lastw,q);
00560               m_matT.coeffRef(i+1,n-1) = internal::real(cc);
00561               m_matT.coeffRef(i+1,n) = internal::imag(cc);
00562             }
00563           }
00564 
00565           // Overflow control
00566           using std::max;
00567           Scalar t = max(internal::abs(m_matT.coeff(i,n-1)),internal::abs(m_matT.coeff(i,n)));
00568           if ((eps * t) * t > Scalar(1))
00569             m_matT.block(i, n-1, size-i, 2) /= t;
00570 
00571         }
00572       }
00573     }
00574     else
00575     {
00576       eigen_assert("Internal bug in EigenSolver"); // this should not happen
00577     }
00578   }
00579 
00580   // Back transformation to get eigenvectors of original matrix
00581   for (Index j = size-1; j >= 0; j--)
00582   {
00583     m_tmp.noalias() = m_eivec.leftCols(j+1) * m_matT.col(j).segment(0, j+1);
00584     m_eivec.col(j) = m_tmp;
00585   }
00586 }
00587 
00588 #endif // EIGEN_EIGENSOLVER_H


re_vision
Author(s): Dorian Galvez-Lopez
autogenerated on Sun Jan 5 2014 11:31:07