LDLT.h revision c981c48f5bc9aefeffc0bcb0cc3934c2fae179dd
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
5// Copyright (C) 2009 Keir Mierle <mierle@gmail.com>
6// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
7// Copyright (C) 2011 Timothy E. Holy <tim.holy@gmail.com >
8//
9// This Source Code Form is subject to the terms of the Mozilla
10// Public License v. 2.0. If a copy of the MPL was not distributed
11// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
12
13#ifndef EIGEN_LDLT_H
14#define EIGEN_LDLT_H
15
16namespace Eigen {
17
18namespace internal {
19template<typename MatrixType, int UpLo> struct LDLT_Traits;
20}
21
22/** \ingroup Cholesky_Module
23  *
24  * \class LDLT
25  *
26  * \brief Robust Cholesky decomposition of a matrix with pivoting
27  *
28  * \param MatrixType the type of the matrix of which to compute the LDL^T Cholesky decomposition
29  * \param UpLo the triangular part that will be used for the decompositon: Lower (default) or Upper.
30  *             The other triangular part won't be read.
31  *
32  * Perform a robust Cholesky decomposition of a positive semidefinite or negative semidefinite
33  * matrix \f$ A \f$ such that \f$ A =  P^TLDL^*P \f$, where P is a permutation matrix, L
34  * is lower triangular with a unit diagonal and D is a diagonal matrix.
35  *
36  * The decomposition uses pivoting to ensure stability, so that L will have
37  * zeros in the bottom right rank(A) - n submatrix. Avoiding the square root
38  * on D also stabilizes the computation.
39  *
40  * Remember that Cholesky decompositions are not rank-revealing. Also, do not use a Cholesky
41  * decomposition to determine whether a system of equations has a solution.
42  *
43  * \sa MatrixBase::ldlt(), class LLT
44  */
45template<typename _MatrixType, int _UpLo> class LDLT
46{
47  public:
48    typedef _MatrixType MatrixType;
49    enum {
50      RowsAtCompileTime = MatrixType::RowsAtCompileTime,
51      ColsAtCompileTime = MatrixType::ColsAtCompileTime,
52      Options = MatrixType::Options & ~RowMajorBit, // these are the options for the TmpMatrixType, we need a ColMajor matrix here!
53      MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
54      MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
55      UpLo = _UpLo
56    };
57    typedef typename MatrixType::Scalar Scalar;
58    typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
59    typedef typename MatrixType::Index Index;
60    typedef Matrix<Scalar, RowsAtCompileTime, 1, Options, MaxRowsAtCompileTime, 1> TmpMatrixType;
61
62    typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType;
63    typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType;
64
65    typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
66
67    /** \brief Default Constructor.
68      *
69      * The default constructor is useful in cases in which the user intends to
70      * perform decompositions via LDLT::compute(const MatrixType&).
71      */
72    LDLT() : m_matrix(), m_transpositions(), m_isInitialized(false) {}
73
74    /** \brief Default Constructor with memory preallocation
75      *
76      * Like the default constructor but with preallocation of the internal data
77      * according to the specified problem \a size.
78      * \sa LDLT()
79      */
80    LDLT(Index size)
81      : m_matrix(size, size),
82        m_transpositions(size),
83        m_temporary(size),
84        m_isInitialized(false)
85    {}
86
87    /** \brief Constructor with decomposition
88      *
89      * This calculates the decomposition for the input \a matrix.
90      * \sa LDLT(Index size)
91      */
92    LDLT(const MatrixType& matrix)
93      : m_matrix(matrix.rows(), matrix.cols()),
94        m_transpositions(matrix.rows()),
95        m_temporary(matrix.rows()),
96        m_isInitialized(false)
97    {
98      compute(matrix);
99    }
100
101    /** Clear any existing decomposition
102     * \sa rankUpdate(w,sigma)
103     */
104    void setZero()
105    {
106      m_isInitialized = false;
107    }
108
109    /** \returns a view of the upper triangular matrix U */
110    inline typename Traits::MatrixU matrixU() const
111    {
112      eigen_assert(m_isInitialized && "LDLT is not initialized.");
113      return Traits::getU(m_matrix);
114    }
115
116    /** \returns a view of the lower triangular matrix L */
117    inline typename Traits::MatrixL matrixL() const
118    {
119      eigen_assert(m_isInitialized && "LDLT is not initialized.");
120      return Traits::getL(m_matrix);
121    }
122
123    /** \returns the permutation matrix P as a transposition sequence.
124      */
125    inline const TranspositionType& transpositionsP() const
126    {
127      eigen_assert(m_isInitialized && "LDLT is not initialized.");
128      return m_transpositions;
129    }
130
131    /** \returns the coefficients of the diagonal matrix D */
132    inline Diagonal<const MatrixType> vectorD() const
133    {
134      eigen_assert(m_isInitialized && "LDLT is not initialized.");
135      return m_matrix.diagonal();
136    }
137
138    /** \returns true if the matrix is positive (semidefinite) */
139    inline bool isPositive() const
140    {
141      eigen_assert(m_isInitialized && "LDLT is not initialized.");
142      return m_sign == 1;
143    }
144
145    #ifdef EIGEN2_SUPPORT
146    inline bool isPositiveDefinite() const
147    {
148      return isPositive();
149    }
150    #endif
151
152    /** \returns true if the matrix is negative (semidefinite) */
153    inline bool isNegative(void) const
154    {
155      eigen_assert(m_isInitialized && "LDLT is not initialized.");
156      return m_sign == -1;
157    }
158
159    /** \returns a solution x of \f$ A x = b \f$ using the current decomposition of A.
160      *
161      * This function also supports in-place solves using the syntax <tt>x = decompositionObject.solve(x)</tt> .
162      *
163      * \note_about_checking_solutions
164      *
165      * More precisely, this method solves \f$ A x = b \f$ using the decomposition \f$ A = P^T L D L^* P \f$
166      * by solving the systems \f$ P^T y_1 = b \f$, \f$ L y_2 = y_1 \f$, \f$ D y_3 = y_2 \f$,
167      * \f$ L^* y_4 = y_3 \f$ and \f$ P x = y_4 \f$ in succession. If the matrix \f$ A \f$ is singular, then
168      * \f$ D \f$ will also be singular (all the other matrices are invertible). In that case, the
169      * least-square solution of \f$ D y_3 = y_2 \f$ is computed. This does not mean that this function
170      * computes the least-square solution of \f$ A x = b \f$ is \f$ A \f$ is singular.
171      *
172      * \sa MatrixBase::ldlt()
173      */
174    template<typename Rhs>
175    inline const internal::solve_retval<LDLT, Rhs>
176    solve(const MatrixBase<Rhs>& b) const
177    {
178      eigen_assert(m_isInitialized && "LDLT is not initialized.");
179      eigen_assert(m_matrix.rows()==b.rows()
180                && "LDLT::solve(): invalid number of rows of the right hand side matrix b");
181      return internal::solve_retval<LDLT, Rhs>(*this, b.derived());
182    }
183
184    #ifdef EIGEN2_SUPPORT
185    template<typename OtherDerived, typename ResultType>
186    bool solve(const MatrixBase<OtherDerived>& b, ResultType *result) const
187    {
188      *result = this->solve(b);
189      return true;
190    }
191    #endif
192
193    template<typename Derived>
194    bool solveInPlace(MatrixBase<Derived> &bAndX) const;
195
196    LDLT& compute(const MatrixType& matrix);
197
198    template <typename Derived>
199    LDLT& rankUpdate(const MatrixBase<Derived>& w,RealScalar alpha=1);
200
201    /** \returns the internal LDLT decomposition matrix
202      *
203      * TODO: document the storage layout
204      */
205    inline const MatrixType& matrixLDLT() const
206    {
207      eigen_assert(m_isInitialized && "LDLT is not initialized.");
208      return m_matrix;
209    }
210
211    MatrixType reconstructedMatrix() const;
212
213    inline Index rows() const { return m_matrix.rows(); }
214    inline Index cols() const { return m_matrix.cols(); }
215
216    /** \brief Reports whether previous computation was successful.
217      *
218      * \returns \c Success if computation was succesful,
219      *          \c NumericalIssue if the matrix.appears to be negative.
220      */
221    ComputationInfo info() const
222    {
223      eigen_assert(m_isInitialized && "LDLT is not initialized.");
224      return Success;
225    }
226
227  protected:
228
229    /** \internal
230      * Used to compute and store the Cholesky decomposition A = L D L^* = U^* D U.
231      * The strict upper part is used during the decomposition, the strict lower
232      * part correspond to the coefficients of L (its diagonal is equal to 1 and
233      * is not stored), and the diagonal entries correspond to D.
234      */
235    MatrixType m_matrix;
236    TranspositionType m_transpositions;
237    TmpMatrixType m_temporary;
238    int m_sign;
239    bool m_isInitialized;
240};
241
242namespace internal {
243
244template<int UpLo> struct ldlt_inplace;
245
246template<> struct ldlt_inplace<Lower>
247{
248  template<typename MatrixType, typename TranspositionType, typename Workspace>
249  static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, int* sign=0)
250  {
251    typedef typename MatrixType::Scalar Scalar;
252    typedef typename MatrixType::RealScalar RealScalar;
253    typedef typename MatrixType::Index Index;
254    eigen_assert(mat.rows()==mat.cols());
255    const Index size = mat.rows();
256
257    if (size <= 1)
258    {
259      transpositions.setIdentity();
260      if(sign)
261        *sign = real(mat.coeff(0,0))>0 ? 1:-1;
262      return true;
263    }
264
265    RealScalar cutoff(0), biggest_in_corner;
266
267    for (Index k = 0; k < size; ++k)
268    {
269      // Find largest diagonal element
270      Index index_of_biggest_in_corner;
271      biggest_in_corner = mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
272      index_of_biggest_in_corner += k;
273
274      if(k == 0)
275      {
276        // The biggest overall is the point of reference to which further diagonals
277        // are compared; if any diagonal is negligible compared
278        // to the largest overall, the algorithm bails.
279        cutoff = abs(NumTraits<Scalar>::epsilon() * biggest_in_corner);
280
281        if(sign)
282          *sign = real(mat.diagonal().coeff(index_of_biggest_in_corner)) > 0 ? 1 : -1;
283      }
284
285      // Finish early if the matrix is not full rank.
286      if(biggest_in_corner < cutoff)
287      {
288        for(Index i = k; i < size; i++) transpositions.coeffRef(i) = i;
289        break;
290      }
291
292      transpositions.coeffRef(k) = index_of_biggest_in_corner;
293      if(k != index_of_biggest_in_corner)
294      {
295        // apply the transposition while taking care to consider only
296        // the lower triangular part
297        Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
298        mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
299        mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
300        std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
301        for(int i=k+1;i<index_of_biggest_in_corner;++i)
302        {
303          Scalar tmp = mat.coeffRef(i,k);
304          mat.coeffRef(i,k) = conj(mat.coeffRef(index_of_biggest_in_corner,i));
305          mat.coeffRef(index_of_biggest_in_corner,i) = conj(tmp);
306        }
307        if(NumTraits<Scalar>::IsComplex)
308          mat.coeffRef(index_of_biggest_in_corner,k) = conj(mat.coeff(index_of_biggest_in_corner,k));
309      }
310
311      // partition the matrix:
312      //       A00 |  -  |  -
313      // lu  = A10 | A11 |  -
314      //       A20 | A21 | A22
315      Index rs = size - k - 1;
316      Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
317      Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
318      Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
319
320      if(k>0)
321      {
322        temp.head(k) = mat.diagonal().head(k).asDiagonal() * A10.adjoint();
323        mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
324        if(rs>0)
325          A21.noalias() -= A20 * temp.head(k);
326      }
327      if((rs>0) && (abs(mat.coeffRef(k,k)) > cutoff))
328        A21 /= mat.coeffRef(k,k);
329    }
330
331    return true;
332  }
333
334  // Reference for the algorithm: Davis and Hager, "Multiple Rank
335  // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
336  // Trivial rearrangements of their computations (Timothy E. Holy)
337  // allow their algorithm to work for rank-1 updates even if the
338  // original matrix is not of full rank.
339  // Here only rank-1 updates are implemented, to reduce the
340  // requirement for intermediate storage and improve accuracy
341  template<typename MatrixType, typename WDerived>
342  static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, typename MatrixType::RealScalar sigma=1)
343  {
344    using internal::isfinite;
345    typedef typename MatrixType::Scalar Scalar;
346    typedef typename MatrixType::RealScalar RealScalar;
347    typedef typename MatrixType::Index Index;
348
349    const Index size = mat.rows();
350    eigen_assert(mat.cols() == size && w.size()==size);
351
352    RealScalar alpha = 1;
353
354    // Apply the update
355    for (Index j = 0; j < size; j++)
356    {
357      // Check for termination due to an original decomposition of low-rank
358      if (!(isfinite)(alpha))
359        break;
360
361      // Update the diagonal terms
362      RealScalar dj = real(mat.coeff(j,j));
363      Scalar wj = w.coeff(j);
364      RealScalar swj2 = sigma*abs2(wj);
365      RealScalar gamma = dj*alpha + swj2;
366
367      mat.coeffRef(j,j) += swj2/alpha;
368      alpha += swj2/dj;
369
370
371      // Update the terms of L
372      Index rs = size-j-1;
373      w.tail(rs) -= wj * mat.col(j).tail(rs);
374      if(gamma != 0)
375        mat.col(j).tail(rs) += (sigma*conj(wj)/gamma)*w.tail(rs);
376    }
377    return true;
378  }
379
380  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
381  static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, typename MatrixType::RealScalar sigma=1)
382  {
383    // Apply the permutation to the input w
384    tmp = transpositions * w;
385
386    return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
387  }
388};
389
390template<> struct ldlt_inplace<Upper>
391{
392  template<typename MatrixType, typename TranspositionType, typename Workspace>
393  static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, int* sign=0)
394  {
395    Transpose<MatrixType> matt(mat);
396    return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
397  }
398
399  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
400  static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, typename MatrixType::RealScalar sigma=1)
401  {
402    Transpose<MatrixType> matt(mat);
403    return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
404  }
405};
406
407template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
408{
409  typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
410  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
411  static inline MatrixL getL(const MatrixType& m) { return m; }
412  static inline MatrixU getU(const MatrixType& m) { return m.adjoint(); }
413};
414
415template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
416{
417  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
418  typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
419  static inline MatrixL getL(const MatrixType& m) { return m.adjoint(); }
420  static inline MatrixU getU(const MatrixType& m) { return m; }
421};
422
423} // end namespace internal
424
425/** Compute / recompute the LDLT decomposition A = L D L^* = U^* D U of \a matrix
426  */
427template<typename MatrixType, int _UpLo>
428LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::compute(const MatrixType& a)
429{
430  eigen_assert(a.rows()==a.cols());
431  const Index size = a.rows();
432
433  m_matrix = a;
434
435  m_transpositions.resize(size);
436  m_isInitialized = false;
437  m_temporary.resize(size);
438
439  internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, &m_sign);
440
441  m_isInitialized = true;
442  return *this;
443}
444
445/** Update the LDLT decomposition:  given A = L D L^T, efficiently compute the decomposition of A + sigma w w^T.
446 * \param w a vector to be incorporated into the decomposition.
447 * \param sigma a scalar, +1 for updates and -1 for "downdates," which correspond to removing previously-added column vectors. Optional; default value is +1.
448 * \sa setZero()
449  */
450template<typename MatrixType, int _UpLo>
451template<typename Derived>
452LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w,typename NumTraits<typename MatrixType::Scalar>::Real sigma)
453{
454  const Index size = w.rows();
455  if (m_isInitialized)
456  {
457    eigen_assert(m_matrix.rows()==size);
458  }
459  else
460  {
461    m_matrix.resize(size,size);
462    m_matrix.setZero();
463    m_transpositions.resize(size);
464    for (Index i = 0; i < size; i++)
465      m_transpositions.coeffRef(i) = i;
466    m_temporary.resize(size);
467    m_sign = sigma>=0 ? 1 : -1;
468    m_isInitialized = true;
469  }
470
471  internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
472
473  return *this;
474}
475
476namespace internal {
477template<typename _MatrixType, int _UpLo, typename Rhs>
478struct solve_retval<LDLT<_MatrixType,_UpLo>, Rhs>
479  : solve_retval_base<LDLT<_MatrixType,_UpLo>, Rhs>
480{
481  typedef LDLT<_MatrixType,_UpLo> LDLTType;
482  EIGEN_MAKE_SOLVE_HELPERS(LDLTType,Rhs)
483
484  template<typename Dest> void evalTo(Dest& dst) const
485  {
486    eigen_assert(rhs().rows() == dec().matrixLDLT().rows());
487    // dst = P b
488    dst = dec().transpositionsP() * rhs();
489
490    // dst = L^-1 (P b)
491    dec().matrixL().solveInPlace(dst);
492
493    // dst = D^-1 (L^-1 P b)
494    // more precisely, use pseudo-inverse of D (see bug 241)
495    using std::abs;
496    using std::max;
497    typedef typename LDLTType::MatrixType MatrixType;
498    typedef typename LDLTType::Scalar Scalar;
499    typedef typename LDLTType::RealScalar RealScalar;
500    const Diagonal<const MatrixType> vectorD = dec().vectorD();
501    RealScalar tolerance = (max)(vectorD.array().abs().maxCoeff() * NumTraits<Scalar>::epsilon(),
502				 RealScalar(1) / NumTraits<RealScalar>::highest()); // motivated by LAPACK's xGELSS
503    for (Index i = 0; i < vectorD.size(); ++i) {
504      if(abs(vectorD(i)) > tolerance)
505	dst.row(i) /= vectorD(i);
506      else
507	dst.row(i).setZero();
508    }
509
510    // dst = L^-T (D^-1 L^-1 P b)
511    dec().matrixU().solveInPlace(dst);
512
513    // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
514    dst = dec().transpositionsP().transpose() * dst;
515  }
516};
517}
518
519/** \internal use x = ldlt_object.solve(x);
520  *
521  * This is the \em in-place version of solve().
522  *
523  * \param bAndX represents both the right-hand side matrix b and result x.
524  *
525  * \returns true always! If you need to check for existence of solutions, use another decomposition like LU, QR, or SVD.
526  *
527  * This version avoids a copy when the right hand side matrix b is not
528  * needed anymore.
529  *
530  * \sa LDLT::solve(), MatrixBase::ldlt()
531  */
532template<typename MatrixType,int _UpLo>
533template<typename Derived>
534bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
535{
536  eigen_assert(m_isInitialized && "LDLT is not initialized.");
537  const Index size = m_matrix.rows();
538  eigen_assert(size == bAndX.rows());
539
540  bAndX = this->solve(bAndX);
541
542  return true;
543}
544
545/** \returns the matrix represented by the decomposition,
546 * i.e., it returns the product: P^T L D L^* P.
547 * This function is provided for debug purpose. */
548template<typename MatrixType, int _UpLo>
549MatrixType LDLT<MatrixType,_UpLo>::reconstructedMatrix() const
550{
551  eigen_assert(m_isInitialized && "LDLT is not initialized.");
552  const Index size = m_matrix.rows();
553  MatrixType res(size,size);
554
555  // P
556  res.setIdentity();
557  res = transpositionsP() * res;
558  // L^* P
559  res = matrixU() * res;
560  // D(L^*P)
561  res = vectorD().asDiagonal() * res;
562  // L(DL^*P)
563  res = matrixL() * res;
564  // P^T (LDL^*P)
565  res = transpositionsP().transpose() * res;
566
567  return res;
568}
569
570/** \cholesky_module
571  * \returns the Cholesky decomposition with full pivoting without square root of \c *this
572  */
573template<typename MatrixType, unsigned int UpLo>
574inline const LDLT<typename SelfAdjointView<MatrixType, UpLo>::PlainObject, UpLo>
575SelfAdjointView<MatrixType, UpLo>::ldlt() const
576{
577  return LDLT<PlainObject,UpLo>(m_matrix);
578}
579
580/** \cholesky_module
581  * \returns the Cholesky decomposition with full pivoting without square root of \c *this
582  */
583template<typename Derived>
584inline const LDLT<typename MatrixBase<Derived>::PlainObject>
585MatrixBase<Derived>::ldlt() const
586{
587  return LDLT<PlainObject>(derived());
588}
589
590} // end namespace Eigen
591
592#endif // EIGEN_LDLT_H
593