Eigen  3.2.91
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  // 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 Eigen::Index Index;
63  typedef typename MatrixType::StorageIndex StorageIndex;
65 
68 
69  typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
70 
76  LDLT()
77  : m_matrix(),
78  m_transpositions(),
79  m_sign(internal::ZeroSign),
80  m_isInitialized(false)
81  {}
82 
89  explicit LDLT(Index size)
90  : m_matrix(size, size),
91  m_transpositions(size),
92  m_temporary(size),
93  m_sign(internal::ZeroSign),
94  m_isInitialized(false)
95  {}
96 
102  explicit LDLT(const MatrixType& matrix)
103  : m_matrix(matrix.rows(), matrix.cols()),
104  m_transpositions(matrix.rows()),
105  m_temporary(matrix.rows()),
106  m_sign(internal::ZeroSign),
107  m_isInitialized(false)
108  {
109  compute(matrix);
110  }
111 
115  void setZero()
116  {
117  m_isInitialized = false;
118  }
119 
121  inline typename Traits::MatrixU matrixU() const
122  {
123  eigen_assert(m_isInitialized && "LDLT is not initialized.");
124  return Traits::getU(m_matrix);
125  }
126 
128  inline typename Traits::MatrixL matrixL() const
129  {
130  eigen_assert(m_isInitialized && "LDLT is not initialized.");
131  return Traits::getL(m_matrix);
132  }
133 
136  inline const TranspositionType& transpositionsP() const
137  {
138  eigen_assert(m_isInitialized && "LDLT is not initialized.");
139  return m_transpositions;
140  }
141 
144  {
145  eigen_assert(m_isInitialized && "LDLT is not initialized.");
146  return m_matrix.diagonal();
147  }
148 
150  inline bool isPositive() const
151  {
152  eigen_assert(m_isInitialized && "LDLT is not initialized.");
153  return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign;
154  }
155 
157  inline bool isNegative(void) const
158  {
159  eigen_assert(m_isInitialized && "LDLT is not initialized.");
160  return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign;
161  }
162 
178  template<typename Rhs>
179  inline const Solve<LDLT, Rhs>
180  solve(const MatrixBase<Rhs>& b) const
181  {
182  eigen_assert(m_isInitialized && "LDLT is not initialized.");
183  eigen_assert(m_matrix.rows()==b.rows()
184  && "LDLT::solve(): invalid number of rows of the right hand side matrix b");
185  return Solve<LDLT, Rhs>(*this, b.derived());
186  }
187 
188  template<typename Derived>
189  bool solveInPlace(MatrixBase<Derived> &bAndX) const;
190 
191  LDLT& compute(const MatrixType& matrix);
192 
193  template <typename Derived>
194  LDLT& rankUpdate(const MatrixBase<Derived>& w, const RealScalar& alpha=1);
195 
200  inline const MatrixType& matrixLDLT() const
201  {
202  eigen_assert(m_isInitialized && "LDLT is not initialized.");
203  return m_matrix;
204  }
205 
206  MatrixType reconstructedMatrix() const;
207 
208  inline Index rows() const { return m_matrix.rows(); }
209  inline Index cols() const { return m_matrix.cols(); }
210 
217  {
218  eigen_assert(m_isInitialized && "LDLT is not initialized.");
219  return Success;
220  }
221 
222  #ifndef EIGEN_PARSED_BY_DOXYGEN
223  template<typename RhsType, typename DstType>
224  EIGEN_DEVICE_FUNC
225  void _solve_impl(const RhsType &rhs, DstType &dst) const;
226  #endif
227 
228  protected:
229 
230  static void check_template_parameters()
231  {
232  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
233  }
234 
241  MatrixType m_matrix;
242  TranspositionType m_transpositions;
243  TmpMatrixType m_temporary;
244  internal::SignMatrix m_sign;
245  bool m_isInitialized;
246 };
247 
248 namespace internal {
249 
250 template<int UpLo> struct ldlt_inplace;
251 
252 template<> struct ldlt_inplace<Lower>
253 {
254  template<typename MatrixType, typename TranspositionType, typename Workspace>
255  static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
256  {
257  using std::abs;
258  typedef typename MatrixType::Scalar Scalar;
259  typedef typename MatrixType::RealScalar RealScalar;
260  typedef typename TranspositionType::StorageIndex IndexType;
261  eigen_assert(mat.rows()==mat.cols());
262  const Index size = mat.rows();
263 
264  if (size <= 1)
265  {
266  transpositions.setIdentity();
267  if (numext::real(mat.coeff(0,0)) > 0) sign = PositiveSemiDef;
268  else if (numext::real(mat.coeff(0,0)) < 0) sign = NegativeSemiDef;
269  else sign = ZeroSign;
270  return true;
271  }
272 
273  for (Index k = 0; k < size; ++k)
274  {
275  // Find largest diagonal element
276  Index index_of_biggest_in_corner;
277  mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
278  index_of_biggest_in_corner += k;
279 
280  transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
281  if(k != index_of_biggest_in_corner)
282  {
283  // apply the transposition while taking care to consider only
284  // the lower triangular part
285  Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
286  mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
287  mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
288  std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
289  for(Index i=k+1;i<index_of_biggest_in_corner;++i)
290  {
291  Scalar tmp = mat.coeffRef(i,k);
292  mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
293  mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
294  }
295  if(NumTraits<Scalar>::IsComplex)
296  mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
297  }
298 
299  // partition the matrix:
300  // A00 | - | -
301  // lu = A10 | A11 | -
302  // A20 | A21 | A22
303  Index rs = size - k - 1;
304  Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
305  Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
306  Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
307 
308  if(k>0)
309  {
310  temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
311  mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
312  if(rs>0)
313  A21.noalias() -= A20 * temp.head(k);
314  }
315 
316  // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot
317  // was smaller than the cutoff value. However, since LDLT is not rank-revealing
318  // we should only make sure that we do not introduce INF or NaN values.
319  // Remark that LAPACK also uses 0 as the cutoff value.
320  RealScalar realAkk = numext::real(mat.coeffRef(k,k));
321  if((rs>0) && (abs(realAkk) > RealScalar(0)))
322  A21 /= realAkk;
323 
324  if (sign == PositiveSemiDef) {
325  if (realAkk < 0) sign = Indefinite;
326  } else if (sign == NegativeSemiDef) {
327  if (realAkk > 0) sign = Indefinite;
328  } else if (sign == ZeroSign) {
329  if (realAkk > 0) sign = PositiveSemiDef;
330  else if (realAkk < 0) sign = NegativeSemiDef;
331  }
332  }
333 
334  return true;
335  }
336 
337  // Reference for the algorithm: Davis and Hager, "Multiple Rank
338  // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
339  // Trivial rearrangements of their computations (Timothy E. Holy)
340  // allow their algorithm to work for rank-1 updates even if the
341  // original matrix is not of full rank.
342  // Here only rank-1 updates are implemented, to reduce the
343  // requirement for intermediate storage and improve accuracy
344  template<typename MatrixType, typename WDerived>
345  static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, const typename MatrixType::RealScalar& sigma=1)
346  {
347  using numext::isfinite;
348  typedef typename MatrixType::Scalar Scalar;
349  typedef typename MatrixType::RealScalar RealScalar;
350 
351  const Index size = mat.rows();
352  eigen_assert(mat.cols() == size && w.size()==size);
353 
354  RealScalar alpha = 1;
355 
356  // Apply the update
357  for (Index j = 0; j < size; j++)
358  {
359  // Check for termination due to an original decomposition of low-rank
360  if (!(isfinite)(alpha))
361  break;
362 
363  // Update the diagonal terms
364  RealScalar dj = numext::real(mat.coeff(j,j));
365  Scalar wj = w.coeff(j);
366  RealScalar swj2 = sigma*numext::abs2(wj);
367  RealScalar gamma = dj*alpha + swj2;
368 
369  mat.coeffRef(j,j) += swj2/alpha;
370  alpha += swj2/dj;
371 
372 
373  // Update the terms of L
374  Index rs = size-j-1;
375  w.tail(rs) -= wj * mat.col(j).tail(rs);
376  if(gamma != 0)
377  mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
378  }
379  return true;
380  }
381 
382  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
383  static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1)
384  {
385  // Apply the permutation to the input w
386  tmp = transpositions * w;
387 
388  return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
389  }
390 };
391 
392 template<> struct ldlt_inplace<Upper>
393 {
394  template<typename MatrixType, typename TranspositionType, typename Workspace>
395  static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
396  {
397  Transpose<MatrixType> matt(mat);
398  return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
399  }
400 
401  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
402  static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1)
403  {
404  Transpose<MatrixType> matt(mat);
405  return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
406  }
407 };
408 
409 template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
410 {
411  typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
412  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
413  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
414  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
415 };
416 
417 template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
418 {
419  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
420  typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
421  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); }
422  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); }
423 };
424 
425 } // end namespace internal
426 
429 template<typename MatrixType, int _UpLo>
431 {
432  check_template_parameters();
433 
434  eigen_assert(a.rows()==a.cols());
435  const Index size = a.rows();
436 
437  m_matrix = a;
438 
439  m_transpositions.resize(size);
440  m_isInitialized = false;
441  m_temporary.resize(size);
442  m_sign = internal::ZeroSign;
443 
444  internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign);
445 
446  m_isInitialized = true;
447  return *this;
448 }
449 
455 template<typename MatrixType, int _UpLo>
456 template<typename Derived>
457 LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w, const typename LDLT<MatrixType,_UpLo>::RealScalar& sigma)
458 {
459  typedef typename TranspositionType::StorageIndex IndexType;
460  const Index size = w.rows();
461  if (m_isInitialized)
462  {
463  eigen_assert(m_matrix.rows()==size);
464  }
465  else
466  {
467  m_matrix.resize(size,size);
468  m_matrix.setZero();
469  m_transpositions.resize(size);
470  for (Index i = 0; i < size; i++)
471  m_transpositions.coeffRef(i) = IndexType(i);
472  m_temporary.resize(size);
473  m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
474  m_isInitialized = true;
475  }
476 
477  internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
478 
479  return *this;
480 }
481 
482 #ifndef EIGEN_PARSED_BY_DOXYGEN
483 template<typename _MatrixType, int _UpLo>
484 template<typename RhsType, typename DstType>
485 void LDLT<_MatrixType,_UpLo>::_solve_impl(const RhsType &rhs, DstType &dst) const
486 {
487  eigen_assert(rhs.rows() == rows());
488  // dst = P b
489  dst = m_transpositions * rhs;
490 
491  // dst = L^-1 (P b)
492  matrixL().solveInPlace(dst);
493 
494  // dst = D^-1 (L^-1 P b)
495  // more precisely, use pseudo-inverse of D (see bug 241)
496  using std::abs;
497  const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD());
498  // In some previous versions, tolerance was set to the max of 1/highest and the maximal diagonal entry * epsilon
499  // as motivated by LAPACK's xGELSS:
500  // RealScalar tolerance = numext::maxi(vecD.array().abs().maxCoeff() * NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest());
501  // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest
502  // diagonal element is not well justified and leads to numerical issues in some cases.
503  // Moreover, Lapack's xSYTRS routines use 0 for the tolerance.
504  RealScalar tolerance = RealScalar(1) / NumTraits<RealScalar>::highest();
505 
506  for (Index i = 0; i < vecD.size(); ++i)
507  {
508  if(abs(vecD(i)) > tolerance)
509  dst.row(i) /= vecD(i);
510  else
511  dst.row(i).setZero();
512  }
513 
514  // dst = L^-T (D^-1 L^-1 P b)
515  matrixU().solveInPlace(dst);
516 
517  // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
518  dst = m_transpositions.transpose() * dst;
519 }
520 #endif
521 
535 template<typename MatrixType,int _UpLo>
536 template<typename Derived>
537 bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
538 {
539  eigen_assert(m_isInitialized && "LDLT is not initialized.");
540  eigen_assert(m_matrix.rows() == bAndX.rows());
541 
542  bAndX = this->solve(bAndX);
543 
544  return true;
545 }
546 
550 template<typename MatrixType, int _UpLo>
552 {
553  eigen_assert(m_isInitialized && "LDLT is not initialized.");
554  const Index size = m_matrix.rows();
555  MatrixType res(size,size);
556 
557  // P
558  res.setIdentity();
559  res = transpositionsP() * res;
560  // L^* P
561  res = matrixU() * res;
562  // D(L^*P)
563  res = vectorD().real().asDiagonal() * res;
564  // L(DL^*P)
565  res = matrixL() * res;
566  // P^T (LDL^*P)
567  res = transpositionsP().transpose() * res;
568 
569  return res;
570 }
571 
572 #ifndef __CUDACC__
573 
577 template<typename MatrixType, unsigned int UpLo>
580 {
581  return LDLT<PlainObject,UpLo>(m_matrix);
582 }
583 
588 template<typename Derived>
591 {
592  return LDLT<PlainObject>(derived());
593 }
594 #endif // __CUDACC__
595 
596 } // end namespace Eigen
597 
598 #endif // EIGEN_LDLT_H
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:48
const LDLT< PlainObject > ldlt() const
Definition: LDLT.h:590
Definition: Constants.h:196
LDLT(Index size)
Default Constructor with memory preallocation.
Definition: LDLT.h:89
Definition: LDLT.h:16
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:107
void setZero()
Definition: LDLT.h:115
const unsigned int RowMajorBit
Definition: Constants.h:53
const MatrixType & matrixLDLT() const
Definition: LDLT.h:200
Definition: Constants.h:198
LDLT()
Default Constructor.
Definition: LDLT.h:76
bool isNegative(void) const
Definition: LDLT.h:157
MatrixType reconstructedMatrix() const
Definition: LDLT.h:551
const TranspositionType & transpositionsP() const
Definition: LDLT.h:136
Traits::MatrixU matrixU() const
Definition: LDLT.h:121
Definition: Constants.h:424
Eigen::Index Index
Definition: LDLT.h:62
const Solve< LDLT, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: LDLT.h:180
Definition: Eigen_Colamd.h:54
bool isPositive() const
Definition: LDLT.h:150
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LDLT.h:216
Diagonal< const MatrixType > vectorD() const
Definition: LDLT.h:143
Traits::MatrixL matrixL() const
Definition: LDLT.h:128
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:63
Pseudo expression representing a solving operation.
Definition: Solve.h:63
ComputationInfo
Definition: Constants.h:422
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
LDLT & compute(const MatrixType &matrix)
Definition: LDLT.h:430
LDLT(const MatrixType &matrix)
Constructor with decomposition.
Definition: LDLT.h:102
const LDLT< PlainObject, UpLo > ldlt() const
Definition: LDLT.h:579