LDLT.h
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 
16 namespace Eigen {
17 
18 namespace internal {
19 template<typename MatrixType, int UpLo> struct LDLT_Traits;
20 }
21 
45 template<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;
61 
64 
65  typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
66 
72  LDLT() : m_matrix(), m_transpositions(), m_isInitialized(false) {}
73 
80  LDLT(Index size)
81  : m_matrix(size, size),
82  m_transpositions(size),
83  m_temporary(size),
84  m_isInitialized(false)
85  {}
86 
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 
104  void setZero()
105  {
106  m_isInitialized = false;
107  }
108 
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 
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 
125  inline const TranspositionType& transpositionsP() const
126  {
127  eigen_assert(m_isInitialized && "LDLT is not initialized.");
128  return m_transpositions;
129  }
130 
133  {
134  eigen_assert(m_isInitialized && "LDLT is not initialized.");
135  return m_matrix.diagonal();
136  }
137 
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 
153  inline bool isNegative(void) const
154  {
155  eigen_assert(m_isInitialized && "LDLT is not initialized.");
156  return m_sign == -1;
157  }
158 
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 
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 
222  {
223  eigen_assert(m_isInitialized && "LDLT is not initialized.");
224  return Success;
225  }
226 
227  protected:
228 
235  MatrixType m_matrix;
236  TranspositionType m_transpositions;
237  TmpMatrixType m_temporary;
238  int m_sign;
239  bool m_isInitialized;
240 };
241 
242 namespace internal {
243 
244 template<int UpLo> struct ldlt_inplace;
245 
246 template<> 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 
282  // Finish early if the matrix is not full rank.
283  if(biggest_in_corner < cutoff)
284  {
285  for(Index i = k; i < size; i++) transpositions.coeffRef(i) = i;
286  if(sign) *sign = 0;
287  break;
288  }
289 
290  transpositions.coeffRef(k) = index_of_biggest_in_corner;
291  if(k != index_of_biggest_in_corner)
292  {
293  // apply the transposition while taking care to consider only
294  // the lower triangular part
295  Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
296  mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
297  mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
298  std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
299  for(int i=k+1;i<index_of_biggest_in_corner;++i)
300  {
301  Scalar tmp = mat.coeffRef(i,k);
302  mat.coeffRef(i,k) = conj(mat.coeffRef(index_of_biggest_in_corner,i));
303  mat.coeffRef(index_of_biggest_in_corner,i) = conj(tmp);
304  }
305  if(NumTraits<Scalar>::IsComplex)
306  mat.coeffRef(index_of_biggest_in_corner,k) = conj(mat.coeff(index_of_biggest_in_corner,k));
307  }
308 
309  // partition the matrix:
310  // A00 | - | -
311  // lu = A10 | A11 | -
312  // A20 | A21 | A22
313  Index rs = size - k - 1;
314  Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
315  Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
316  Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
317 
318  if(k>0)
319  {
320  temp.head(k) = mat.diagonal().head(k).asDiagonal() * A10.adjoint();
321  mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
322  if(rs>0)
323  A21.noalias() -= A20 * temp.head(k);
324  }
325  if((rs>0) && (abs(mat.coeffRef(k,k)) > cutoff))
326  A21 /= mat.coeffRef(k,k);
327 
328  if(sign)
329  {
330  // LDLT is not guaranteed to work for indefinite matrices, but let's try to get the sign right
331  int newSign = real(mat.diagonal().coeff(index_of_biggest_in_corner)) > 0;
332  if(k == 0)
333  *sign = newSign;
334  else if(*sign != newSign)
335  *sign = 0;
336  }
337  }
338 
339  return true;
340  }
341 
342  // Reference for the algorithm: Davis and Hager, "Multiple Rank
343  // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
344  // Trivial rearrangements of their computations (Timothy E. Holy)
345  // allow their algorithm to work for rank-1 updates even if the
346  // original matrix is not of full rank.
347  // Here only rank-1 updates are implemented, to reduce the
348  // requirement for intermediate storage and improve accuracy
349  template<typename MatrixType, typename WDerived>
350  static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, typename MatrixType::RealScalar sigma=1)
351  {
352  using internal::isfinite;
353  typedef typename MatrixType::Scalar Scalar;
354  typedef typename MatrixType::RealScalar RealScalar;
355  typedef typename MatrixType::Index Index;
356 
357  const Index size = mat.rows();
358  eigen_assert(mat.cols() == size && w.size()==size);
359 
360  RealScalar alpha = 1;
361 
362  // Apply the update
363  for (Index j = 0; j < size; j++)
364  {
365  // Check for termination due to an original decomposition of low-rank
366  if (!(isfinite)(alpha))
367  break;
368 
369  // Update the diagonal terms
370  RealScalar dj = real(mat.coeff(j,j));
371  Scalar wj = w.coeff(j);
372  RealScalar swj2 = sigma*abs2(wj);
373  RealScalar gamma = dj*alpha + swj2;
374 
375  mat.coeffRef(j,j) += swj2/alpha;
376  alpha += swj2/dj;
377 
378 
379  // Update the terms of L
380  Index rs = size-j-1;
381  w.tail(rs) -= wj * mat.col(j).tail(rs);
382  if(gamma != 0)
383  mat.col(j).tail(rs) += (sigma*conj(wj)/gamma)*w.tail(rs);
384  }
385  return true;
386  }
387 
388  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
389  static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, typename MatrixType::RealScalar sigma=1)
390  {
391  // Apply the permutation to the input w
392  tmp = transpositions * w;
393 
394  return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
395  }
396 };
397 
398 template<> struct ldlt_inplace<Upper>
399 {
400  template<typename MatrixType, typename TranspositionType, typename Workspace>
401  static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, int* sign=0)
402  {
403  Transpose<MatrixType> matt(mat);
404  return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
405  }
406 
407  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
408  static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, typename MatrixType::RealScalar sigma=1)
409  {
410  Transpose<MatrixType> matt(mat);
411  return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
412  }
413 };
414 
415 template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
416 {
417  typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
418  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
419  static inline MatrixL getL(const MatrixType& m) { return m; }
420  static inline MatrixU getU(const MatrixType& m) { return m.adjoint(); }
421 };
422 
423 template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
424 {
425  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
426  typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
427  static inline MatrixL getL(const MatrixType& m) { return m.adjoint(); }
428  static inline MatrixU getU(const MatrixType& m) { return m; }
429 };
430 
431 } // end namespace internal
432 
435 template<typename MatrixType, int _UpLo>
437 {
438  eigen_assert(a.rows()==a.cols());
439  const Index size = a.rows();
440 
441  m_matrix = a;
442 
443  m_transpositions.resize(size);
444  m_isInitialized = false;
445  m_temporary.resize(size);
446 
447  internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, &m_sign);
448 
449  m_isInitialized = true;
450  return *this;
451 }
452 
458 template<typename MatrixType, int _UpLo>
459 template<typename Derived>
461 {
462  const Index size = w.rows();
463  if (m_isInitialized)
464  {
465  eigen_assert(m_matrix.rows()==size);
466  }
467  else
468  {
469  m_matrix.resize(size,size);
470  m_matrix.setZero();
471  m_transpositions.resize(size);
472  for (Index i = 0; i < size; i++)
473  m_transpositions.coeffRef(i) = i;
474  m_temporary.resize(size);
475  m_sign = sigma>=0 ? 1 : -1;
476  m_isInitialized = true;
477  }
478 
479  internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
480 
481  return *this;
482 }
483 
484 namespace internal {
485 template<typename _MatrixType, int _UpLo, typename Rhs>
486 struct solve_retval<LDLT<_MatrixType,_UpLo>, Rhs>
487  : solve_retval_base<LDLT<_MatrixType,_UpLo>, Rhs>
488 {
489  typedef LDLT<_MatrixType,_UpLo> LDLTType;
490  EIGEN_MAKE_SOLVE_HELPERS(LDLTType,Rhs)
491 
492  template<typename Dest> void evalTo(Dest& dst) const
493  {
494  eigen_assert(rhs().rows() == dec().matrixLDLT().rows());
495  // dst = P b
496  dst = dec().transpositionsP() * rhs();
497 
498  // dst = L^-1 (P b)
499  dec().matrixL().solveInPlace(dst);
500 
501  // dst = D^-1 (L^-1 P b)
502  // more precisely, use pseudo-inverse of D (see bug 241)
503  using std::abs;
504  using std::max;
505  typedef typename LDLTType::MatrixType MatrixType;
506  typedef typename LDLTType::Scalar Scalar;
507  typedef typename LDLTType::RealScalar RealScalar;
508  const Diagonal<const MatrixType> vectorD = dec().vectorD();
509  RealScalar tolerance = (max)(vectorD.array().abs().maxCoeff() * NumTraits<Scalar>::epsilon(),
510  RealScalar(1) / NumTraits<RealScalar>::highest()); // motivated by LAPACK's xGELSS
511  for (Index i = 0; i < vectorD.size(); ++i) {
512  if(abs(vectorD(i)) > tolerance)
513  dst.row(i) /= vectorD(i);
514  else
515  dst.row(i).setZero();
516  }
517 
518  // dst = L^-T (D^-1 L^-1 P b)
519  dec().matrixU().solveInPlace(dst);
520 
521  // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
522  dst = dec().transpositionsP().transpose() * dst;
523  }
524 };
525 }
526 
540 template<typename MatrixType,int _UpLo>
541 template<typename Derived>
542 bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
543 {
544  eigen_assert(m_isInitialized && "LDLT is not initialized.");
545  eigen_assert(m_matrix.rows() == bAndX.rows());
546 
547  bAndX = this->solve(bAndX);
548 
549  return true;
550 }
551 
555 template<typename MatrixType, int _UpLo>
557 {
558  eigen_assert(m_isInitialized && "LDLT is not initialized.");
559  const Index size = m_matrix.rows();
560  MatrixType res(size,size);
561 
562  // P
563  res.setIdentity();
564  res = transpositionsP() * res;
565  // L^* P
566  res = matrixU() * res;
567  // D(L^*P)
568  res = vectorD().asDiagonal() * res;
569  // L(DL^*P)
570  res = matrixL() * res;
571  // P^T (LDL^*P)
572  res = transpositionsP().transpose() * res;
573 
574  return res;
575 }
576 
580 template<typename MatrixType, unsigned int UpLo>
583 {
584  return LDLT<PlainObject,UpLo>(m_matrix);
585 }
586 
590 template<typename Derived>
593 {
594  return LDLT<PlainObject>(derived());
595 }
596 
597 } // end namespace Eigen
598 
599 #endif // EIGEN_LDLT_H