Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
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 <[email protected]>
5 // Copyright (C) 2009 Keir Mierle <[email protected]>
6 // Copyright (C) 2009 Benoit Jacob <[email protected]>
7 // Copyright (C) 2011 Timothy E. Holy <[email protected] >
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  // PositiveSemiDef means positive semi-definite and non-zero; same for NegativeSemiDef
22  enum SignMatrix { PositiveSemiDef, NegativeSemiDef, ZeroSign, Indefinite };
23 }
24 
48 template<typename _MatrixType, int _UpLo> class LDLT
49 {
50  public:
51  typedef _MatrixType MatrixType;
52  enum {
53  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
54  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
55  Options = MatrixType::Options & ~RowMajorBit, // these are the options for the TmpMatrixType, we need a ColMajor matrix here!
56  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
57  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
58  UpLo = _UpLo
59  };
60  typedef typename MatrixType::Scalar Scalar;
61  typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
62  typedef typename MatrixType::Index Index;
64 
67 
68  typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
69 
75  LDLT()
76  : m_matrix(),
77  m_transpositions(),
78  m_sign(internal::ZeroSign),
79  m_isInitialized(false)
80  {}
81 
88  LDLT(Index size)
89  : m_matrix(size, size),
90  m_transpositions(size),
91  m_temporary(size),
92  m_sign(internal::ZeroSign),
93  m_isInitialized(false)
94  {}
95 
101  LDLT(const MatrixType& matrix)
102  : m_matrix(matrix.rows(), matrix.cols()),
103  m_transpositions(matrix.rows()),
104  m_temporary(matrix.rows()),
105  m_sign(internal::ZeroSign),
106  m_isInitialized(false)
107  {
108  compute(matrix);
109  }
110 
114  void setZero()
115  {
116  m_isInitialized = false;
117  }
118 
120  inline typename Traits::MatrixU matrixU() const
121  {
122  eigen_assert(m_isInitialized && "LDLT is not initialized.");
123  return Traits::getU(m_matrix);
124  }
125 
127  inline typename Traits::MatrixL matrixL() const
128  {
129  eigen_assert(m_isInitialized && "LDLT is not initialized.");
130  return Traits::getL(m_matrix);
131  }
132 
135  inline const TranspositionType& transpositionsP() const
136  {
137  eigen_assert(m_isInitialized && "LDLT is not initialized.");
138  return m_transpositions;
139  }
140 
143  {
144  eigen_assert(m_isInitialized && "LDLT is not initialized.");
145  return m_matrix.diagonal();
146  }
147 
149  inline bool isPositive() const
150  {
151  eigen_assert(m_isInitialized && "LDLT is not initialized.");
152  return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign;
153  }
154 
155  #ifdef EIGEN2_SUPPORT
156  inline bool isPositiveDefinite() const
157  {
158  return isPositive();
159  }
160  #endif
161 
163  inline bool isNegative(void) const
164  {
165  eigen_assert(m_isInitialized && "LDLT is not initialized.");
166  return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign;
167  }
168 
184  template<typename Rhs>
185  inline const internal::solve_retval<LDLT, Rhs>
186  solve(const MatrixBase<Rhs>& b) const
187  {
188  eigen_assert(m_isInitialized && "LDLT is not initialized.");
189  eigen_assert(m_matrix.rows()==b.rows()
190  && "LDLT::solve(): invalid number of rows of the right hand side matrix b");
191  return internal::solve_retval<LDLT, Rhs>(*this, b.derived());
192  }
193 
194  #ifdef EIGEN2_SUPPORT
195  template<typename OtherDerived, typename ResultType>
196  bool solve(const MatrixBase<OtherDerived>& b, ResultType *result) const
197  {
198  *result = this->solve(b);
199  return true;
200  }
201  #endif
202 
203  template<typename Derived>
204  bool solveInPlace(MatrixBase<Derived> &bAndX) const;
205 
206  LDLT& compute(const MatrixType& matrix);
207 
208  template <typename Derived>
209  LDLT& rankUpdate(const MatrixBase<Derived>& w, const RealScalar& alpha=1);
210 
215  inline const MatrixType& matrixLDLT() const
216  {
217  eigen_assert(m_isInitialized && "LDLT is not initialized.");
218  return m_matrix;
219  }
220 
221  MatrixType reconstructedMatrix() const;
222 
223  inline Index rows() const { return m_matrix.rows(); }
224  inline Index cols() const { return m_matrix.cols(); }
225 
232  {
233  eigen_assert(m_isInitialized && "LDLT is not initialized.");
234  return Success;
235  }
236 
237  protected:
238 
239  static void check_template_parameters()
240  {
241  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
242  }
243 
250  MatrixType m_matrix;
251  TranspositionType m_transpositions;
252  TmpMatrixType m_temporary;
253  internal::SignMatrix m_sign;
254  bool m_isInitialized;
255 };
256 
257 namespace internal {
258 
259 template<int UpLo> struct ldlt_inplace;
260 
261 template<> struct ldlt_inplace<Lower>
262 {
263  template<typename MatrixType, typename TranspositionType, typename Workspace>
264  static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
265  {
266  using std::abs;
267  typedef typename MatrixType::Scalar Scalar;
268  typedef typename MatrixType::RealScalar RealScalar;
269  typedef typename MatrixType::Index Index;
270  eigen_assert(mat.rows()==mat.cols());
271  const Index size = mat.rows();
272 
273  if (size <= 1)
274  {
275  transpositions.setIdentity();
276  if (numext::real(mat.coeff(0,0)) > 0) sign = PositiveSemiDef;
277  else if (numext::real(mat.coeff(0,0)) < 0) sign = NegativeSemiDef;
278  else sign = ZeroSign;
279  return true;
280  }
281 
282  for (Index k = 0; k < size; ++k)
283  {
284  // Find largest diagonal element
285  Index index_of_biggest_in_corner;
286  mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
287  index_of_biggest_in_corner += k;
288 
289  transpositions.coeffRef(k) = index_of_biggest_in_corner;
290  if(k != index_of_biggest_in_corner)
291  {
292  // apply the transposition while taking care to consider only
293  // the lower triangular part
294  Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
295  mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
296  mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
297  std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
298  for(int i=k+1;i<index_of_biggest_in_corner;++i)
299  {
300  Scalar tmp = mat.coeffRef(i,k);
301  mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
302  mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
303  }
304  if(NumTraits<Scalar>::IsComplex)
305  mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
306  }
307 
308  // partition the matrix:
309  // A00 | - | -
310  // lu = A10 | A11 | -
311  // A20 | A21 | A22
312  Index rs = size - k - 1;
313  Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
314  Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
315  Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
316 
317  if(k>0)
318  {
319  temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
320  mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
321  if(rs>0)
322  A21.noalias() -= A20 * temp.head(k);
323  }
324 
325  // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot
326  // was smaller than the cutoff value. However, soince LDLT is not rank-revealing
327  // we should only make sure we do not introduce INF or NaN values.
328  // LAPACK also uses 0 as the cutoff value.
329  RealScalar realAkk = numext::real(mat.coeffRef(k,k));
330  if((rs>0) && (abs(realAkk) > RealScalar(0)))
331  A21 /= realAkk;
332 
333  if (sign == PositiveSemiDef) {
334  if (realAkk < 0) sign = Indefinite;
335  } else if (sign == NegativeSemiDef) {
336  if (realAkk > 0) sign = Indefinite;
337  } else if (sign == ZeroSign) {
338  if (realAkk > 0) sign = PositiveSemiDef;
339  else if (realAkk < 0) sign = NegativeSemiDef;
340  }
341  }
342 
343  return true;
344  }
345 
346  // Reference for the algorithm: Davis and Hager, "Multiple Rank
347  // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
348  // Trivial rearrangements of their computations (Timothy E. Holy)
349  // allow their algorithm to work for rank-1 updates even if the
350  // original matrix is not of full rank.
351  // Here only rank-1 updates are implemented, to reduce the
352  // requirement for intermediate storage and improve accuracy
353  template<typename MatrixType, typename WDerived>
354  static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, const typename MatrixType::RealScalar& sigma=1)
355  {
356  using numext::isfinite;
357  typedef typename MatrixType::Scalar Scalar;
358  typedef typename MatrixType::RealScalar RealScalar;
359  typedef typename MatrixType::Index Index;
360 
361  const Index size = mat.rows();
362  eigen_assert(mat.cols() == size && w.size()==size);
363 
364  RealScalar alpha = 1;
365 
366  // Apply the update
367  for (Index j = 0; j < size; j++)
368  {
369  // Check for termination due to an original decomposition of low-rank
370  if (!(isfinite)(alpha))
371  break;
372 
373  // Update the diagonal terms
374  RealScalar dj = numext::real(mat.coeff(j,j));
375  Scalar wj = w.coeff(j);
376  RealScalar swj2 = sigma*numext::abs2(wj);
377  RealScalar gamma = dj*alpha + swj2;
378 
379  mat.coeffRef(j,j) += swj2/alpha;
380  alpha += swj2/dj;
381 
382 
383  // Update the terms of L
384  Index rs = size-j-1;
385  w.tail(rs) -= wj * mat.col(j).tail(rs);
386  if(gamma != 0)
387  mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
388  }
389  return true;
390  }
391 
392  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
393  static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1)
394  {
395  // Apply the permutation to the input w
396  tmp = transpositions * w;
397 
398  return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
399  }
400 };
401 
402 template<> struct ldlt_inplace<Upper>
403 {
404  template<typename MatrixType, typename TranspositionType, typename Workspace>
405  static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
406  {
407  Transpose<MatrixType> matt(mat);
408  return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
409  }
410 
411  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
412  static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1)
413  {
414  Transpose<MatrixType> matt(mat);
415  return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
416  }
417 };
418 
419 template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
420 {
421  typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
422  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
423  static inline MatrixL getL(const MatrixType& m) { return m; }
424  static inline MatrixU getU(const MatrixType& m) { return m.adjoint(); }
425 };
426 
427 template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
428 {
429  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
430  typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
431  static inline MatrixL getL(const MatrixType& m) { return m.adjoint(); }
432  static inline MatrixU getU(const MatrixType& m) { return m; }
433 };
434 
435 } // end namespace internal
436 
439 template<typename MatrixType, int _UpLo>
441 {
442  check_template_parameters();
443 
444  eigen_assert(a.rows()==a.cols());
445  const Index size = a.rows();
446 
447  m_matrix = a;
448 
449  m_transpositions.resize(size);
450  m_isInitialized = false;
451  m_temporary.resize(size);
452  m_sign = internal::ZeroSign;
453 
454  internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign);
455 
456  m_isInitialized = true;
457  return *this;
458 }
459 
465 template<typename MatrixType, int _UpLo>
466 template<typename Derived>
467 LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w, const typename LDLT<MatrixType,_UpLo>::RealScalar& sigma)
468 {
469  const Index size = w.rows();
470  if (m_isInitialized)
471  {
472  eigen_assert(m_matrix.rows()==size);
473  }
474  else
475  {
476  m_matrix.resize(size,size);
477  m_matrix.setZero();
478  m_transpositions.resize(size);
479  for (Index i = 0; i < size; i++)
480  m_transpositions.coeffRef(i) = i;
481  m_temporary.resize(size);
482  m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
483  m_isInitialized = true;
484  }
485 
486  internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
487 
488  return *this;
489 }
490 
491 namespace internal {
492 template<typename _MatrixType, int _UpLo, typename Rhs>
493 struct solve_retval<LDLT<_MatrixType,_UpLo>, Rhs>
494  : solve_retval_base<LDLT<_MatrixType,_UpLo>, Rhs>
495 {
496  typedef LDLT<_MatrixType,_UpLo> LDLTType;
497  EIGEN_MAKE_SOLVE_HELPERS(LDLTType,Rhs)
498 
499  template<typename Dest> void evalTo(Dest& dst) const
500  {
501  eigen_assert(rhs().rows() == dec().matrixLDLT().rows());
502  // dst = P b
503  dst = dec().transpositionsP() * rhs();
504 
505  // dst = L^-1 (P b)
506  dec().matrixL().solveInPlace(dst);
507 
508  // dst = D^-1 (L^-1 P b)
509  // more precisely, use pseudo-inverse of D (see bug 241)
510  using std::abs;
511  using std::max;
512  typedef typename LDLTType::MatrixType MatrixType;
513  typedef typename LDLTType::RealScalar RealScalar;
514  const typename Diagonal<const MatrixType>::RealReturnType vectorD(dec().vectorD());
515  // In some previous versions, tolerance was set to the max of 1/highest and the maximal diagonal entry * epsilon
516  // as motivated by LAPACK's xGELSS:
517  // RealScalar tolerance = (max)(vectorD.array().abs().maxCoeff() *NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest());
518  // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest
519  // diagonal element is not well justified and to numerical issues in some cases.
520  // Moreover, Lapack's xSYTRS routines use 0 for the tolerance.
521  RealScalar tolerance = RealScalar(1) / NumTraits<RealScalar>::highest();
522 
523  for (Index i = 0; i < vectorD.size(); ++i) {
524  if(abs(vectorD(i)) > tolerance)
525  dst.row(i) /= vectorD(i);
526  else
527  dst.row(i).setZero();
528  }
529 
530  // dst = L^-T (D^-1 L^-1 P b)
531  dec().matrixU().solveInPlace(dst);
532 
533  // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
534  dst = dec().transpositionsP().transpose() * dst;
535  }
536 };
537 }
538 
552 template<typename MatrixType,int _UpLo>
553 template<typename Derived>
554 bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
555 {
556  eigen_assert(m_isInitialized && "LDLT is not initialized.");
557  eigen_assert(m_matrix.rows() == bAndX.rows());
558 
559  bAndX = this->solve(bAndX);
560 
561  return true;
562 }
563 
567 template<typename MatrixType, int _UpLo>
569 {
570  eigen_assert(m_isInitialized && "LDLT is not initialized.");
571  const Index size = m_matrix.rows();
572  MatrixType res(size,size);
573 
574  // P
575  res.setIdentity();
576  res = transpositionsP() * res;
577  // L^* P
578  res = matrixU() * res;
579  // D(L^*P)
580  res = vectorD().real().asDiagonal() * res;
581  // L(DL^*P)
582  res = matrixL() * res;
583  // P^T (LDL^*P)
584  res = transpositionsP().transpose() * res;
585 
586  return res;
587 }
588 
592 template<typename MatrixType, unsigned int UpLo>
595 {
596  return LDLT<PlainObject,UpLo>(m_matrix);
597 }
598 
602 template<typename Derived>
605 {
606  return LDLT<PlainObject>(derived());
607 }
608 
609 } // end namespace Eigen
610 
611 #endif // EIGEN_LDLT_H
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:48
Definition: Constants.h:167
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LDLT.h:231
LDLT(const MatrixType &matrix)
Constructor with decomposition.
Definition: LDLT.h:101
MatrixType reconstructedMatrix() const
Definition: LDLT.h:568
const TranspositionType & transpositionsP() const
Definition: LDLT.h:135
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
Definition: Constants.h:169
Traits::MatrixL matrixL() const
Definition: LDLT.h:127
const internal::solve_retval< LDLT, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: LDLT.h:186
bool isPositive() const
Definition: LDLT.h:149
const LDLT< PlainObject, UpLo > ldlt() const
Definition: LDLT.h:594
LDLT & compute(const MatrixType &matrix)
Definition: LDLT.h:440
LDLT(Index size)
Default Constructor with memory preallocation.
Definition: LDLT.h:88
void setZero()
Definition: LDLT.h:114
Diagonal< const MatrixType > vectorD() const
Definition: LDLT.h:142
LDLT()
Default Constructor.
Definition: LDLT.h:75
Definition: Constants.h:376
bool isNegative(void) const
Definition: LDLT.h:163
const unsigned int RowMajorBit
Definition: Constants.h:53
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:64
ComputationInfo
Definition: Constants.h:374
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
const LDLT< PlainObject > ldlt() const
Definition: LDLT.h:604
Traits::MatrixU matrixU() const
Definition: LDLT.h:120
const MatrixType & matrixLDLT() const
Definition: LDLT.h:215