LDLT.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-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 // Copyright (C) 2009 Keir Mierle <mierle@gmail.com>
00006 // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
00007 // Copyright (C) 2011 Timothy E. Holy <tim.holy@gmail.com >
00008 //
00009 // Eigen is free software; you can redistribute it and/or
00010 // modify it under the terms of the GNU Lesser General Public
00011 // License as published by the Free Software Foundation; either
00012 // version 3 of the License, or (at your option) any later version.
00013 //
00014 // Alternatively, you can redistribute it and/or
00015 // modify it under the terms of the GNU General Public License as
00016 // published by the Free Software Foundation; either version 2 of
00017 // the License, or (at your option) any later version.
00018 //
00019 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00020 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00021 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00022 // GNU General Public License for more details.
00023 //
00024 // You should have received a copy of the GNU Lesser General Public
00025 // License and a copy of the GNU General Public License along with
00026 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00027 
00028 #ifndef EIGEN_LDLT_H
00029 #define EIGEN_LDLT_H
00030 
00031 namespace Eigen { 
00032 
00033 namespace internal {
00034 template<typename MatrixType, int UpLo> struct LDLT_Traits;
00035 }
00036 
00060 template<typename _MatrixType, int _UpLo> class LDLT
00061 {
00062   public:
00063     typedef _MatrixType MatrixType;
00064     enum {
00065       RowsAtCompileTime = MatrixType::RowsAtCompileTime,
00066       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
00067       Options = MatrixType::Options & ~RowMajorBit, // these are the options for the TmpMatrixType, we need a ColMajor matrix here!
00068       MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
00069       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
00070       UpLo = _UpLo
00071     };
00072     typedef typename MatrixType::Scalar Scalar;
00073     typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
00074     typedef typename MatrixType::Index Index;
00075     typedef Matrix<Scalar, RowsAtCompileTime, 1, Options, MaxRowsAtCompileTime, 1> TmpMatrixType;
00076 
00077     typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType;
00078     typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType;
00079 
00080     typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
00081 
00087     LDLT() : m_matrix(), m_transpositions(), m_isInitialized(false) {}
00088 
00095     LDLT(Index size)
00096       : m_matrix(size, size),
00097         m_transpositions(size),
00098         m_temporary(size),
00099         m_isInitialized(false)
00100     {}
00101 
00107     LDLT(const MatrixType& matrix)
00108       : m_matrix(matrix.rows(), matrix.cols()),
00109         m_transpositions(matrix.rows()),
00110         m_temporary(matrix.rows()),
00111         m_isInitialized(false)
00112     {
00113       compute(matrix);
00114     }
00115 
00119     void setZero()
00120     {
00121       m_isInitialized = false;
00122     }
00123 
00125     inline typename Traits::MatrixU matrixU() const
00126     {
00127       eigen_assert(m_isInitialized && "LDLT is not initialized.");
00128       return Traits::getU(m_matrix);
00129     }
00130 
00132     inline typename Traits::MatrixL matrixL() const
00133     {
00134       eigen_assert(m_isInitialized && "LDLT is not initialized.");
00135       return Traits::getL(m_matrix);
00136     }
00137 
00140     inline const TranspositionType& transpositionsP() const
00141     {
00142       eigen_assert(m_isInitialized && "LDLT is not initialized.");
00143       return m_transpositions;
00144     }
00145 
00147     inline Diagonal<const MatrixType> vectorD() const
00148     {
00149       eigen_assert(m_isInitialized && "LDLT is not initialized.");
00150       return m_matrix.diagonal();
00151     }
00152 
00154     inline bool isPositive() const
00155     {
00156       eigen_assert(m_isInitialized && "LDLT is not initialized.");
00157       return m_sign == 1;
00158     }
00159     
00160     #ifdef EIGEN2_SUPPORT
00161     inline bool isPositiveDefinite() const
00162     {
00163       return isPositive();
00164     }
00165     #endif
00166 
00168     inline bool isNegative(void) const
00169     {
00170       eigen_assert(m_isInitialized && "LDLT is not initialized.");
00171       return m_sign == -1;
00172     }
00173 
00189     template<typename Rhs>
00190     inline const internal::solve_retval<LDLT, Rhs>
00191     solve(const MatrixBase<Rhs>& b) const
00192     {
00193       eigen_assert(m_isInitialized && "LDLT is not initialized.");
00194       eigen_assert(m_matrix.rows()==b.rows()
00195                 && "LDLT::solve(): invalid number of rows of the right hand side matrix b");
00196       return internal::solve_retval<LDLT, Rhs>(*this, b.derived());
00197     }
00198 
00199     #ifdef EIGEN2_SUPPORT
00200     template<typename OtherDerived, typename ResultType>
00201     bool solve(const MatrixBase<OtherDerived>& b, ResultType *result) const
00202     {
00203       *result = this->solve(b);
00204       return true;
00205     }
00206     #endif
00207 
00208     template<typename Derived>
00209     bool solveInPlace(MatrixBase<Derived> &bAndX) const;
00210 
00211     LDLT& compute(const MatrixType& matrix);
00212 
00213     template <typename Derived>
00214     LDLT& rankUpdate(const MatrixBase<Derived>& w,RealScalar alpha=1);
00215 
00220     inline const MatrixType& matrixLDLT() const
00221     {
00222       eigen_assert(m_isInitialized && "LDLT is not initialized.");
00223       return m_matrix;
00224     }
00225 
00226     MatrixType reconstructedMatrix() const;
00227 
00228     inline Index rows() const { return m_matrix.rows(); }
00229     inline Index cols() const { return m_matrix.cols(); }
00230 
00236     ComputationInfo info() const
00237     {
00238       eigen_assert(m_isInitialized && "LDLT is not initialized.");
00239       return Success;
00240     }
00241 
00242   protected:
00243 
00250     MatrixType m_matrix;
00251     TranspositionType m_transpositions;
00252     TmpMatrixType m_temporary;
00253     int m_sign;
00254     bool m_isInitialized;
00255 };
00256 
00257 namespace internal {
00258 
00259 template<int UpLo> struct ldlt_inplace;
00260 
00261 template<> struct ldlt_inplace<Lower>
00262 {
00263   template<typename MatrixType, typename TranspositionType, typename Workspace>
00264   static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, int* sign=0)
00265   {
00266     typedef typename MatrixType::Scalar Scalar;
00267     typedef typename MatrixType::RealScalar RealScalar;
00268     typedef typename MatrixType::Index Index;
00269     eigen_assert(mat.rows()==mat.cols());
00270     const Index size = mat.rows();
00271 
00272     if (size <= 1)
00273     {
00274       transpositions.setIdentity();
00275       if(sign)
00276         *sign = real(mat.coeff(0,0))>0 ? 1:-1;
00277       return true;
00278     }
00279 
00280     RealScalar cutoff(0), biggest_in_corner;
00281 
00282     for (Index k = 0; k < size; ++k)
00283     {
00284       // Find largest diagonal element
00285       Index index_of_biggest_in_corner;
00286       biggest_in_corner = mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
00287       index_of_biggest_in_corner += k;
00288 
00289       if(k == 0)
00290       {
00291         // The biggest overall is the point of reference to which further diagonals
00292         // are compared; if any diagonal is negligible compared
00293         // to the largest overall, the algorithm bails.
00294         cutoff = abs(NumTraits<Scalar>::epsilon() * biggest_in_corner);
00295 
00296         if(sign)
00297           *sign = real(mat.diagonal().coeff(index_of_biggest_in_corner)) > 0 ? 1 : -1;
00298       }
00299 
00300       // Finish early if the matrix is not full rank.
00301       if(biggest_in_corner < cutoff)
00302       {
00303         for(Index i = k; i < size; i++) transpositions.coeffRef(i) = i;
00304         break;
00305       }
00306 
00307       transpositions.coeffRef(k) = index_of_biggest_in_corner;
00308       if(k != index_of_biggest_in_corner)
00309       {
00310         // apply the transposition while taking care to consider only
00311         // the lower triangular part
00312         Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
00313         mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
00314         mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
00315         std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
00316         for(int i=k+1;i<index_of_biggest_in_corner;++i)
00317         {
00318           Scalar tmp = mat.coeffRef(i,k);
00319           mat.coeffRef(i,k) = conj(mat.coeffRef(index_of_biggest_in_corner,i));
00320           mat.coeffRef(index_of_biggest_in_corner,i) = conj(tmp);
00321         }
00322         if(NumTraits<Scalar>::IsComplex)
00323           mat.coeffRef(index_of_biggest_in_corner,k) = conj(mat.coeff(index_of_biggest_in_corner,k));
00324       }
00325 
00326       // partition the matrix:
00327       //       A00 |  -  |  -
00328       // lu  = A10 | A11 |  -
00329       //       A20 | A21 | A22
00330       Index rs = size - k - 1;
00331       Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
00332       Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
00333       Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
00334 
00335       if(k>0)
00336       {
00337         temp.head(k) = mat.diagonal().head(k).asDiagonal() * A10.adjoint();
00338         mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
00339         if(rs>0)
00340           A21.noalias() -= A20 * temp.head(k);
00341       }
00342       if((rs>0) && (abs(mat.coeffRef(k,k)) > cutoff))
00343         A21 /= mat.coeffRef(k,k);
00344     }
00345 
00346     return true;
00347   }
00348 
00349   // Reference for the algorithm: Davis and Hager, "Multiple Rank
00350   // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
00351   // Trivial rearrangements of their computations (Timothy E. Holy)
00352   // allow their algorithm to work for rank-1 updates even if the
00353   // original matrix is not of full rank.
00354   // Here only rank-1 updates are implemented, to reduce the
00355   // requirement for intermediate storage and improve accuracy
00356   template<typename MatrixType, typename WDerived>
00357   static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, typename MatrixType::RealScalar sigma=1)
00358   {
00359     using internal::isfinite;
00360     typedef typename MatrixType::Scalar Scalar;
00361     typedef typename MatrixType::RealScalar RealScalar;
00362     typedef typename MatrixType::Index Index;
00363 
00364     const Index size = mat.rows();
00365     eigen_assert(mat.cols() == size && w.size()==size);
00366 
00367     RealScalar alpha = 1;
00368 
00369     // Apply the update
00370     for (Index j = 0; j < size; j++)
00371     {
00372       // Check for termination due to an original decomposition of low-rank
00373       if (!isfinite(alpha))
00374         break;
00375 
00376       // Update the diagonal terms
00377       RealScalar dj = real(mat.coeff(j,j));
00378       Scalar wj = w.coeff(j);
00379       RealScalar swj2 = sigma*abs2(wj);
00380       RealScalar gamma = dj*alpha + swj2;
00381 
00382       mat.coeffRef(j,j) += swj2/alpha;
00383       alpha += swj2/dj;
00384 
00385 
00386       // Update the terms of L
00387       Index rs = size-j-1;
00388       w.tail(rs) -= wj * mat.col(j).tail(rs);
00389       if(gamma != 0)
00390         mat.col(j).tail(rs) += (sigma*conj(wj)/gamma)*w.tail(rs);
00391     }
00392     return true;
00393   }
00394 
00395   template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
00396   static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, typename MatrixType::RealScalar sigma=1)
00397   {
00398     // Apply the permutation to the input w
00399     tmp = transpositions * w;
00400 
00401     return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
00402   }
00403 };
00404 
00405 template<> struct ldlt_inplace<Upper>
00406 {
00407   template<typename MatrixType, typename TranspositionType, typename Workspace>
00408   static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, int* sign=0)
00409   {
00410     Transpose<MatrixType> matt(mat);
00411     return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
00412   }
00413 
00414   template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
00415   static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, typename MatrixType::RealScalar sigma=1)
00416   {
00417     Transpose<MatrixType> matt(mat);
00418     return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
00419   }
00420 };
00421 
00422 template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
00423 {
00424   typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
00425   typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
00426   static inline MatrixL getL(const MatrixType& m) { return m; }
00427   static inline MatrixU getU(const MatrixType& m) { return m.adjoint(); }
00428 };
00429 
00430 template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
00431 {
00432   typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
00433   typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
00434   static inline MatrixL getL(const MatrixType& m) { return m.adjoint(); }
00435   static inline MatrixU getU(const MatrixType& m) { return m; }
00436 };
00437 
00438 } // end namespace internal
00439 
00442 template<typename MatrixType, int _UpLo>
00443 LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::compute(const MatrixType& a)
00444 {
00445   eigen_assert(a.rows()==a.cols());
00446   const Index size = a.rows();
00447 
00448   m_matrix = a;
00449 
00450   m_transpositions.resize(size);
00451   m_isInitialized = false;
00452   m_temporary.resize(size);
00453 
00454   internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, &m_sign);
00455 
00456   m_isInitialized = true;
00457   return *this;
00458 }
00459 
00465 template<typename MatrixType, int _UpLo>
00466 template<typename Derived>
00467 LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w,typename NumTraits<typename MatrixType::Scalar>::Real sigma)
00468 {
00469   const Index size = w.rows();
00470   if (m_isInitialized)
00471   {
00472     eigen_assert(m_matrix.rows()==size);
00473   }
00474   else
00475   {    
00476     m_matrix.resize(size,size);
00477     m_matrix.setZero();
00478     m_transpositions.resize(size);
00479     for (Index i = 0; i < size; i++)
00480       m_transpositions.coeffRef(i) = i;
00481     m_temporary.resize(size);
00482     m_sign = sigma;
00483     m_isInitialized = true;
00484   }
00485 
00486   internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
00487 
00488   return *this;
00489 }
00490 
00491 namespace internal {
00492 template<typename _MatrixType, int _UpLo, typename Rhs>
00493 struct solve_retval<LDLT<_MatrixType,_UpLo>, Rhs>
00494   : solve_retval_base<LDLT<_MatrixType,_UpLo>, Rhs>
00495 {
00496   typedef LDLT<_MatrixType,_UpLo> LDLTType;
00497   EIGEN_MAKE_SOLVE_HELPERS(LDLTType,Rhs)
00498 
00499   template<typename Dest> void evalTo(Dest& dst) const
00500   {
00501     eigen_assert(rhs().rows() == dec().matrixLDLT().rows());
00502     // dst = P b
00503     dst = dec().transpositionsP() * rhs();
00504 
00505     // dst = L^-1 (P b)
00506     dec().matrixL().solveInPlace(dst);
00507 
00508     // dst = D^-1 (L^-1 P b)
00509     // more precisely, use pseudo-inverse of D (see bug 241)
00510     using std::abs;
00511     using std::max;
00512     typedef typename LDLTType::MatrixType MatrixType;
00513     typedef typename LDLTType::Scalar Scalar;
00514     typedef typename LDLTType::RealScalar RealScalar;
00515     const Diagonal<const MatrixType> vectorD = dec().vectorD();
00516     RealScalar tolerance = (max)(vectorD.array().abs().maxCoeff() * NumTraits<Scalar>::epsilon(),
00517                                  RealScalar(1) / NumTraits<RealScalar>::highest()); // motivated by LAPACK's xGELSS
00518     for (Index i = 0; i < vectorD.size(); ++i) {
00519       if(abs(vectorD(i)) > tolerance)
00520         dst.row(i) /= vectorD(i);
00521       else
00522         dst.row(i).setZero();
00523     }
00524 
00525     // dst = L^-T (D^-1 L^-1 P b)
00526     dec().matrixU().solveInPlace(dst);
00527 
00528     // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
00529     dst = dec().transpositionsP().transpose() * dst;
00530   }
00531 };
00532 }
00533 
00547 template<typename MatrixType,int _UpLo>
00548 template<typename Derived>
00549 bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
00550 {
00551   eigen_assert(m_isInitialized && "LDLT is not initialized.");
00552   const Index size = m_matrix.rows();
00553   eigen_assert(size == bAndX.rows());
00554 
00555   bAndX = this->solve(bAndX);
00556 
00557   return true;
00558 }
00559 
00563 template<typename MatrixType, int _UpLo>
00564 MatrixType LDLT<MatrixType,_UpLo>::reconstructedMatrix() const
00565 {
00566   eigen_assert(m_isInitialized && "LDLT is not initialized.");
00567   const Index size = m_matrix.rows();
00568   MatrixType res(size,size);
00569 
00570   // P
00571   res.setIdentity();
00572   res = transpositionsP() * res;
00573   // L^* P
00574   res = matrixU() * res;
00575   // D(L^*P)
00576   res = vectorD().asDiagonal() * res;
00577   // L(DL^*P)
00578   res = matrixL() * res;
00579   // P^T (LDL^*P)
00580   res = transpositionsP().transpose() * res;
00581 
00582   return res;
00583 }
00584 
00588 template<typename MatrixType, unsigned int UpLo>
00589 inline const LDLT<typename SelfAdjointView<MatrixType, UpLo>::PlainObject, UpLo>
00590 SelfAdjointView<MatrixType, UpLo>::ldlt() const
00591 {
00592   return LDLT<PlainObject,UpLo>(m_matrix);
00593 }
00594 
00598 template<typename Derived>
00599 inline const LDLT<typename MatrixBase<Derived>::PlainObject>
00600 MatrixBase<Derived>::ldlt() const
00601 {
00602   return LDLT<PlainObject>(derived());
00603 }
00604 
00605 } // end namespace Eigen
00606 
00607 #endif // EIGEN_LDLT_H