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