MOAB  4.9.3pre
SelfAdjointEigenSolver.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-2010 Gael Guennebaud <[email protected]>
00005 // Copyright (C) 2010 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_SELFADJOINTEIGENSOLVER_H
00012 #define EIGEN_SELFADJOINTEIGENSOLVER_H
00013 
00014 #include "./Tridiagonalization.h"
00015 
00016 namespace Eigen { 
00017 
00018 template<typename _MatrixType>
00019 class GeneralizedSelfAdjointEigenSolver;
00020 
00021 namespace internal {
00022 template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues;
00023 template<typename MatrixType, typename DiagType, typename SubDiagType>
00024 ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType& eivec);
00025 }
00026 
00070 template<typename _MatrixType> class SelfAdjointEigenSolver
00071 {
00072   public:
00073 
00074     typedef _MatrixType MatrixType;
00075     enum {
00076       Size = MatrixType::RowsAtCompileTime,
00077       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00078       Options = MatrixType::Options,
00079       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
00080     };
00081     
00083     typedef typename MatrixType::Scalar Scalar;
00084     typedef Eigen::Index Index; 
00085     
00086     typedef Matrix<Scalar,Size,Size,ColMajor,MaxColsAtCompileTime,MaxColsAtCompileTime> EigenvectorsType;
00087 
00094     typedef typename NumTraits<Scalar>::Real RealScalar;
00095     
00096     friend struct internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>;
00097 
00103     typedef typename internal::plain_col_type<MatrixType, RealScalar>::type RealVectorType;
00104     typedef Tridiagonalization<MatrixType> TridiagonalizationType;
00105     typedef typename TridiagonalizationType::SubDiagonalType SubDiagonalType;
00106 
00117     EIGEN_DEVICE_FUNC
00118     SelfAdjointEigenSolver()
00119         : m_eivec(),
00120           m_eivalues(),
00121           m_subdiag(),
00122           m_isInitialized(false)
00123     { }
00124 
00137     EIGEN_DEVICE_FUNC
00138     explicit SelfAdjointEigenSolver(Index size)
00139         : m_eivec(size, size),
00140           m_eivalues(size),
00141           m_subdiag(size > 1 ? size - 1 : 1),
00142           m_isInitialized(false)
00143     {}
00144 
00160     template<typename InputType>
00161     EIGEN_DEVICE_FUNC
00162     explicit SelfAdjointEigenSolver(const EigenBase<InputType>& matrix, int options = ComputeEigenvectors)
00163       : m_eivec(matrix.rows(), matrix.cols()),
00164         m_eivalues(matrix.cols()),
00165         m_subdiag(matrix.rows() > 1 ? matrix.rows() - 1 : 1),
00166         m_isInitialized(false)
00167     {
00168       compute(matrix.derived(), options);
00169     }
00170 
00201     template<typename InputType>
00202     EIGEN_DEVICE_FUNC
00203     SelfAdjointEigenSolver& compute(const EigenBase<InputType>& matrix, int options = ComputeEigenvectors);
00204     
00223     EIGEN_DEVICE_FUNC
00224     SelfAdjointEigenSolver& computeDirect(const MatrixType& matrix, int options = ComputeEigenvectors);
00225 
00238     SelfAdjointEigenSolver& computeFromTridiagonal(const RealVectorType& diag, const SubDiagonalType& subdiag , int options=ComputeEigenvectors);
00239 
00258     EIGEN_DEVICE_FUNC
00259     const EigenvectorsType& eigenvectors() const
00260     {
00261       eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
00262       eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
00263       return m_eivec;
00264     }
00265 
00281     EIGEN_DEVICE_FUNC
00282     const RealVectorType& eigenvalues() const
00283     {
00284       eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
00285       return m_eivalues;
00286     }
00287 
00305     EIGEN_DEVICE_FUNC
00306     MatrixType operatorSqrt() const
00307     {
00308       eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
00309       eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
00310       return m_eivec * m_eivalues.cwiseSqrt().asDiagonal() * m_eivec.adjoint();
00311     }
00312 
00330     EIGEN_DEVICE_FUNC
00331     MatrixType operatorInverseSqrt() const
00332     {
00333       eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
00334       eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
00335       return m_eivec * m_eivalues.cwiseInverse().cwiseSqrt().asDiagonal() * m_eivec.adjoint();
00336     }
00337 
00342     EIGEN_DEVICE_FUNC
00343     ComputationInfo info() const
00344     {
00345       eigen_assert(m_isInitialized && "SelfAdjointEigenSolver is not initialized.");
00346       return m_info;
00347     }
00348 
00354     static const int m_maxIterations = 30;
00355 
00356   protected:
00357     static void check_template_parameters()
00358     {
00359       EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
00360     }
00361     
00362     EigenvectorsType m_eivec;
00363     RealVectorType m_eivalues;
00364     typename TridiagonalizationType::SubDiagonalType m_subdiag;
00365     ComputationInfo m_info;
00366     bool m_isInitialized;
00367     bool m_eigenvectorsOk;
00368 };
00369 
00370 namespace internal {
00391 template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
00392 EIGEN_DEVICE_FUNC
00393 static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n);
00394 }
00395 
00396 template<typename MatrixType>
00397 template<typename InputType>
00398 EIGEN_DEVICE_FUNC
00399 SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
00400 ::compute(const EigenBase<InputType>& a_matrix, int options)
00401 {
00402   check_template_parameters();
00403   
00404   const InputType &matrix(a_matrix.derived());
00405   
00406   using std::abs;
00407   eigen_assert(matrix.cols() == matrix.rows());
00408   eigen_assert((options&~(EigVecMask|GenEigMask))==0
00409           && (options&EigVecMask)!=EigVecMask
00410           && "invalid option parameter");
00411   bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
00412   Index n = matrix.cols();
00413   m_eivalues.resize(n,1);
00414 
00415   if(n==1)
00416   {
00417     m_eivalues.coeffRef(0,0) = numext::real(matrix(0,0));
00418     if(computeEigenvectors)
00419       m_eivec.setOnes(n,n);
00420     m_info = Success;
00421     m_isInitialized = true;
00422     m_eigenvectorsOk = computeEigenvectors;
00423     return *this;
00424   }
00425 
00426   // declare some aliases
00427   RealVectorType& diag = m_eivalues;
00428   EigenvectorsType& mat = m_eivec;
00429 
00430   // map the matrix coefficients to [-1:1] to avoid over- and underflow.
00431   mat = matrix.template triangularView<Lower>();
00432   RealScalar scale = mat.cwiseAbs().maxCoeff();
00433   if(scale==RealScalar(0)) scale = RealScalar(1);
00434   mat.template triangularView<Lower>() /= scale;
00435   m_subdiag.resize(n-1);
00436   internal::tridiagonalization_inplace(mat, diag, m_subdiag, computeEigenvectors);
00437 
00438   m_info = internal::computeFromTridiagonal_impl(diag, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
00439   
00440   // scale back the eigen values
00441   m_eivalues *= scale;
00442 
00443   m_isInitialized = true;
00444   m_eigenvectorsOk = computeEigenvectors;
00445   return *this;
00446 }
00447 
00448 template<typename MatrixType>
00449 SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
00450 ::computeFromTridiagonal(const RealVectorType& diag, const SubDiagonalType& subdiag , int options)
00451 {
00452   //TODO : Add an option to scale the values beforehand
00453   bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
00454 
00455   m_eivalues = diag;
00456   m_subdiag = subdiag;
00457   if (computeEigenvectors)
00458   {
00459     m_eivec.setIdentity(diag.size(), diag.size());
00460   }
00461   m_info = computeFromTridiagonal_impl(m_eivalues, m_subdiag, m_maxIterations, computeEigenvectors, m_eivec);
00462 
00463   m_isInitialized = true;
00464   m_eigenvectorsOk = computeEigenvectors;
00465   return *this;
00466 }
00467 
00468 namespace internal {
00480 template<typename MatrixType, typename DiagType, typename SubDiagType>
00481 ComputationInfo computeFromTridiagonal_impl(DiagType& diag, SubDiagType& subdiag, const Index maxIterations, bool computeEigenvectors, MatrixType& eivec)
00482 {
00483   using std::abs;
00484 
00485   ComputationInfo info;
00486   typedef typename MatrixType::Scalar Scalar;
00487 
00488   Index n = diag.size();
00489   Index end = n-1;
00490   Index start = 0;
00491   Index iter = 0; // total number of iterations
00492   
00493   typedef typename DiagType::RealScalar RealScalar;
00494   const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
00495   
00496   while (end>0)
00497   {
00498     for (Index i = start; i<end; ++i)
00499       if (internal::isMuchSmallerThan(abs(subdiag[i]),(abs(diag[i])+abs(diag[i+1]))) || abs(subdiag[i]) <= considerAsZero)
00500         subdiag[i] = 0;
00501 
00502     // find the largest unreduced block
00503     while (end>0 && subdiag[end-1]==0)
00504     {
00505       end--;
00506     }
00507     if (end<=0)
00508       break;
00509 
00510     // if we spent too many iterations, we give up
00511     iter++;
00512     if(iter > maxIterations * n) break;
00513 
00514     start = end - 1;
00515     while (start>0 && subdiag[start-1]!=0)
00516       start--;
00517 
00518     internal::tridiagonal_qr_step<MatrixType::Flags&RowMajorBit ? RowMajor : ColMajor>(diag.data(), subdiag.data(), start, end, computeEigenvectors ? eivec.data() : (Scalar*)0, n);
00519   }
00520   if (iter <= maxIterations * n)
00521     info = Success;
00522   else
00523     info = NoConvergence;
00524 
00525   // Sort eigenvalues and corresponding vectors.
00526   // TODO make the sort optional ?
00527   // TODO use a better sort algorithm !!
00528   if (info == Success)
00529   {
00530     for (Index i = 0; i < n-1; ++i)
00531     {
00532       Index k;
00533       diag.segment(i,n-i).minCoeff(&k);
00534       if (k > 0)
00535       {
00536         std::swap(diag[i], diag[k+i]);
00537         if(computeEigenvectors)
00538           eivec.col(i).swap(eivec.col(k+i));
00539       }
00540     }
00541   }
00542   return info;
00543 }
00544   
00545 template<typename SolverType,int Size,bool IsComplex> struct direct_selfadjoint_eigenvalues
00546 {
00547   EIGEN_DEVICE_FUNC
00548   static inline void run(SolverType& eig, const typename SolverType::MatrixType& A, int options)
00549   { eig.compute(A,options); }
00550 };
00551 
00552 template<typename SolverType> struct direct_selfadjoint_eigenvalues<SolverType,3,false>
00553 {
00554   typedef typename SolverType::MatrixType MatrixType;
00555   typedef typename SolverType::RealVectorType VectorType;
00556   typedef typename SolverType::Scalar Scalar;
00557   typedef typename SolverType::EigenvectorsType EigenvectorsType;
00558   
00559 
00564   EIGEN_DEVICE_FUNC
00565   static inline void computeRoots(const MatrixType& m, VectorType& roots)
00566   {
00567     EIGEN_USING_STD_MATH(sqrt)
00568     EIGEN_USING_STD_MATH(atan2)
00569     EIGEN_USING_STD_MATH(cos)
00570     EIGEN_USING_STD_MATH(sin)
00571     const Scalar s_inv3 = Scalar(1.0)/Scalar(3.0);
00572     const Scalar s_sqrt3 = sqrt(Scalar(3.0));
00573 
00574     // The characteristic equation is x^3 - c2*x^2 + c1*x - c0 = 0.  The
00575     // eigenvalues are the roots to this equation, all guaranteed to be
00576     // real-valued, because the matrix is symmetric.
00577     Scalar c0 = m(0,0)*m(1,1)*m(2,2) + Scalar(2)*m(1,0)*m(2,0)*m(2,1) - m(0,0)*m(2,1)*m(2,1) - m(1,1)*m(2,0)*m(2,0) - m(2,2)*m(1,0)*m(1,0);
00578     Scalar c1 = m(0,0)*m(1,1) - m(1,0)*m(1,0) + m(0,0)*m(2,2) - m(2,0)*m(2,0) + m(1,1)*m(2,2) - m(2,1)*m(2,1);
00579     Scalar c2 = m(0,0) + m(1,1) + m(2,2);
00580 
00581     // Construct the parameters used in classifying the roots of the equation
00582     // and in solving the equation for the roots in closed form.
00583     Scalar c2_over_3 = c2*s_inv3;
00584     Scalar a_over_3 = (c2*c2_over_3 - c1)*s_inv3;
00585     a_over_3 = numext::maxi(a_over_3, Scalar(0));
00586 
00587     Scalar half_b = Scalar(0.5)*(c0 + c2_over_3*(Scalar(2)*c2_over_3*c2_over_3 - c1));
00588 
00589     Scalar q = a_over_3*a_over_3*a_over_3 - half_b*half_b;
00590     q = numext::maxi(q, Scalar(0));
00591 
00592     // Compute the eigenvalues by solving for the roots of the polynomial.
00593     Scalar rho = sqrt(a_over_3);
00594     Scalar theta = atan2(sqrt(q),half_b)*s_inv3;  // since sqrt(q) > 0, atan2 is in [0, pi] and theta is in [0, pi/3]
00595     Scalar cos_theta = cos(theta);
00596     Scalar sin_theta = sin(theta);
00597     // roots are already sorted, since cos is monotonically decreasing on [0, pi]
00598     roots(0) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta); // == 2*rho*cos(theta+2pi/3)
00599     roots(1) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta); // == 2*rho*cos(theta+ pi/3)
00600     roots(2) = c2_over_3 + Scalar(2)*rho*cos_theta;
00601   }
00602 
00603   EIGEN_DEVICE_FUNC
00604   static inline bool extract_kernel(MatrixType& mat, Ref<VectorType> res, Ref<VectorType> representative)
00605   {
00606     using std::abs;
00607     Index i0;
00608     // Find non-zero column i0 (by construction, there must exist a non zero coefficient on the diagonal):
00609     mat.diagonal().cwiseAbs().maxCoeff(&i0);
00610     // mat.col(i0) is a good candidate for an orthogonal vector to the current eigenvector,
00611     // so let's save it:
00612     representative = mat.col(i0);
00613     Scalar n0, n1;
00614     VectorType c0, c1;
00615     n0 = (c0 = representative.cross(mat.col((i0+1)%3))).squaredNorm();
00616     n1 = (c1 = representative.cross(mat.col((i0+2)%3))).squaredNorm();
00617     if(n0>n1) res = c0/std::sqrt(n0);
00618     else      res = c1/std::sqrt(n1);
00619 
00620     return true;
00621   }
00622 
00623   EIGEN_DEVICE_FUNC
00624   static inline void run(SolverType& solver, const MatrixType& mat, int options)
00625   {
00626     eigen_assert(mat.cols() == 3 && mat.cols() == mat.rows());
00627     eigen_assert((options&~(EigVecMask|GenEigMask))==0
00628             && (options&EigVecMask)!=EigVecMask
00629             && "invalid option parameter");
00630     bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
00631     
00632     EigenvectorsType& eivecs = solver.m_eivec;
00633     VectorType& eivals = solver.m_eivalues;
00634   
00635     // Shift the matrix to the mean eigenvalue and map the matrix coefficients to [-1:1] to avoid over- and underflow.
00636     Scalar shift = mat.trace() / Scalar(3);
00637     // TODO Avoid this copy. Currently it is necessary to suppress bogus values when determining maxCoeff and for computing the eigenvectors later
00638     MatrixType scaledMat = mat.template selfadjointView<Lower>();
00639     scaledMat.diagonal().array() -= shift;
00640     Scalar scale = scaledMat.cwiseAbs().maxCoeff();
00641     if(scale > 0) scaledMat /= scale;   // TODO for scale==0 we could save the remaining operations
00642 
00643     // compute the eigenvalues
00644     computeRoots(scaledMat,eivals);
00645 
00646     // compute the eigenvectors
00647     if(computeEigenvectors)
00648     {
00649       if((eivals(2)-eivals(0))<=Eigen::NumTraits<Scalar>::epsilon())
00650       {
00651         // All three eigenvalues are numerically the same
00652         eivecs.setIdentity();
00653       }
00654       else
00655       {
00656         MatrixType tmp;
00657         tmp = scaledMat;
00658 
00659         // Compute the eigenvector of the most distinct eigenvalue
00660         Scalar d0 = eivals(2) - eivals(1);
00661         Scalar d1 = eivals(1) - eivals(0);
00662         Index k(0), l(2);
00663         if(d0 > d1)
00664         {
00665           numext::swap(k,l);
00666           d0 = d1;
00667         }
00668 
00669         // Compute the eigenvector of index k
00670         {
00671           tmp.diagonal().array () -= eivals(k);
00672           // By construction, 'tmp' is of rank 2, and its kernel corresponds to the respective eigenvector.
00673           extract_kernel(tmp, eivecs.col(k), eivecs.col(l));
00674         }
00675 
00676         // Compute eigenvector of index l
00677         if(d0<=2*Eigen::NumTraits<Scalar>::epsilon()*d1)
00678         {
00679           // If d0 is too small, then the two other eigenvalues are numerically the same,
00680           // and thus we only have to ortho-normalize the near orthogonal vector we saved above.
00681           eivecs.col(l) -= eivecs.col(k).dot(eivecs.col(l))*eivecs.col(l);
00682           eivecs.col(l).normalize();
00683         }
00684         else
00685         {
00686           tmp = scaledMat;
00687           tmp.diagonal().array () -= eivals(l);
00688 
00689           VectorType dummy;
00690           extract_kernel(tmp, eivecs.col(l), dummy);
00691         }
00692 
00693         // Compute last eigenvector from the other two
00694         eivecs.col(1) = eivecs.col(2).cross(eivecs.col(0)).normalized();
00695       }
00696     }
00697 
00698     // Rescale back to the original size.
00699     eivals *= scale;
00700     eivals.array() += shift;
00701     
00702     solver.m_info = Success;
00703     solver.m_isInitialized = true;
00704     solver.m_eigenvectorsOk = computeEigenvectors;
00705   }
00706 };
00707 
00708 // 2x2 direct eigenvalues decomposition, code from Hauke Heibel
00709 template<typename SolverType> 
00710 struct direct_selfadjoint_eigenvalues<SolverType,2,false>
00711 {
00712   typedef typename SolverType::MatrixType MatrixType;
00713   typedef typename SolverType::RealVectorType VectorType;
00714   typedef typename SolverType::Scalar Scalar;
00715   typedef typename SolverType::EigenvectorsType EigenvectorsType;
00716   
00717   EIGEN_DEVICE_FUNC
00718   static inline void computeRoots(const MatrixType& m, VectorType& roots)
00719   {
00720     using std::sqrt;
00721     const Scalar t0 = Scalar(0.5) * sqrt( numext::abs2(m(0,0)-m(1,1)) + Scalar(4)*numext::abs2(m(1,0)));
00722     const Scalar t1 = Scalar(0.5) * (m(0,0) + m(1,1));
00723     roots(0) = t1 - t0;
00724     roots(1) = t1 + t0;
00725   }
00726   
00727   EIGEN_DEVICE_FUNC
00728   static inline void run(SolverType& solver, const MatrixType& mat, int options)
00729   {
00730     EIGEN_USING_STD_MATH(sqrt);
00731     EIGEN_USING_STD_MATH(abs);
00732     
00733     eigen_assert(mat.cols() == 2 && mat.cols() == mat.rows());
00734     eigen_assert((options&~(EigVecMask|GenEigMask))==0
00735             && (options&EigVecMask)!=EigVecMask
00736             && "invalid option parameter");
00737     bool computeEigenvectors = (options&ComputeEigenvectors)==ComputeEigenvectors;
00738     
00739     EigenvectorsType& eivecs = solver.m_eivec;
00740     VectorType& eivals = solver.m_eivalues;
00741   
00742     // map the matrix coefficients to [-1:1] to avoid over- and underflow.
00743     Scalar scale = mat.cwiseAbs().maxCoeff();
00744     scale = numext::maxi(scale,Scalar(1));
00745     MatrixType scaledMat = mat / scale;
00746     
00747     // Compute the eigenvalues
00748     computeRoots(scaledMat,eivals);
00749     
00750     // compute the eigen vectors
00751     if(computeEigenvectors)
00752     {
00753       if((eivals(1)-eivals(0))<=abs(eivals(1))*Eigen::NumTraits<Scalar>::epsilon())
00754       {
00755         eivecs.setIdentity();
00756       }
00757       else
00758       {
00759         scaledMat.diagonal().array () -= eivals(1);
00760         Scalar a2 = numext::abs2(scaledMat(0,0));
00761         Scalar c2 = numext::abs2(scaledMat(1,1));
00762         Scalar b2 = numext::abs2(scaledMat(1,0));
00763         if(a2>c2)
00764         {
00765           eivecs.col(1) << -scaledMat(1,0), scaledMat(0,0);
00766           eivecs.col(1) /= sqrt(a2+b2);
00767         }
00768         else
00769         {
00770           eivecs.col(1) << -scaledMat(1,1), scaledMat(1,0);
00771           eivecs.col(1) /= sqrt(c2+b2);
00772         }
00773 
00774         eivecs.col(0) << eivecs.col(1).unitOrthogonal();
00775       }
00776     }
00777     
00778     // Rescale back to the original size.
00779     eivals *= scale;
00780     
00781     solver.m_info = Success;
00782     solver.m_isInitialized = true;
00783     solver.m_eigenvectorsOk = computeEigenvectors;
00784   }
00785 };
00786 
00787 }
00788 
00789 template<typename MatrixType>
00790 EIGEN_DEVICE_FUNC
00791 SelfAdjointEigenSolver<MatrixType>& SelfAdjointEigenSolver<MatrixType>
00792 ::computeDirect(const MatrixType& matrix, int options)
00793 {
00794   internal::direct_selfadjoint_eigenvalues<SelfAdjointEigenSolver,Size,NumTraits<Scalar>::IsComplex>::run(*this,matrix,options);
00795   return *this;
00796 }
00797 
00798 namespace internal {
00799 template<int StorageOrder,typename RealScalar, typename Scalar, typename Index>
00800 EIGEN_DEVICE_FUNC
00801 static void tridiagonal_qr_step(RealScalar* diag, RealScalar* subdiag, Index start, Index end, Scalar* matrixQ, Index n)
00802 {
00803   using std::abs;
00804   RealScalar td = (diag[end-1] - diag[end])*RealScalar(0.5);
00805   RealScalar e = subdiag[end-1];
00806   // Note that thanks to scaling, e^2 or td^2 cannot overflow, however they can still
00807   // underflow thus leading to inf/NaN values when using the following commented code:
00808 //   RealScalar e2 = numext::abs2(subdiag[end-1]);
00809 //   RealScalar mu = diag[end] - e2 / (td + (td>0 ? 1 : -1) * sqrt(td*td + e2));
00810   // This explain the following, somewhat more complicated, version:
00811   RealScalar mu = diag[end];
00812   if(td==0)
00813     mu -= abs(e);
00814   else
00815   {
00816     RealScalar e2 = numext::abs2(subdiag[end-1]);
00817     RealScalar h = numext::hypot(td,e);
00818     if(e2==0)  mu -= (e / (td + (td>0 ? 1 : -1))) * (e / h);
00819     else       mu -= e2 / (td + (td>0 ? h : -h));
00820   }
00821   
00822   RealScalar x = diag[start] - mu;
00823   RealScalar z = subdiag[start];
00824   for (Index k = start; k < end; ++k)
00825   {
00826     JacobiRotation<RealScalar> rot;
00827     rot.makeGivens(x, z);
00828 
00829     // do T = G' T G
00830     RealScalar sdk = rot.s() * diag[k] + rot.c() * subdiag[k];
00831     RealScalar dkp1 = rot.s() * subdiag[k] + rot.c() * diag[k+1];
00832 
00833     diag[k] = rot.c() * (rot.c() * diag[k] - rot.s() * subdiag[k]) - rot.s() * (rot.c() * subdiag[k] - rot.s() * diag[k+1]);
00834     diag[k+1] = rot.s() * sdk + rot.c() * dkp1;
00835     subdiag[k] = rot.c() * sdk - rot.s() * dkp1;
00836     
00837 
00838     if (k > start)
00839       subdiag[k - 1] = rot.c() * subdiag[k-1] - rot.s() * z;
00840 
00841     x = subdiag[k];
00842 
00843     if (k < end - 1)
00844     {
00845       z = -rot.s() * subdiag[k+1];
00846       subdiag[k + 1] = rot.c() * subdiag[k+1];
00847     }
00848     
00849     // apply the givens rotation to the unit matrix Q = Q * G
00850     if (matrixQ)
00851     {
00852       // FIXME if StorageOrder == RowMajor this operation is not very efficient
00853       Map<Matrix<Scalar,Dynamic,Dynamic,StorageOrder> > q(matrixQ,n,n);
00854       q.applyOnTheRight(k,k+1,rot);
00855     }
00856   }
00857 }
00858 
00859 } // end namespace internal
00860 
00861 } // end namespace Eigen
00862 
00863 #endif // EIGEN_SELFADJOINTEIGENSOLVER_H
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines