Eigen  3.2.93
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 
50 template<typename _MatrixType, int _UpLo> class LDLT
51 {
52  public:
53  typedef _MatrixType MatrixType;
54  enum {
55  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
56  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
57  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
58  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
59  UpLo = _UpLo
60  };
61  typedef typename MatrixType::Scalar Scalar;
62  typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
63  typedef Eigen::Index Index;
64  typedef typename MatrixType::StorageIndex StorageIndex;
66 
69 
70  typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
71 
77  LDLT()
78  : m_matrix(),
79  m_transpositions(),
80  m_sign(internal::ZeroSign),
81  m_isInitialized(false)
82  {}
83 
90  explicit LDLT(Index size)
91  : m_matrix(size, size),
92  m_transpositions(size),
93  m_temporary(size),
94  m_sign(internal::ZeroSign),
95  m_isInitialized(false)
96  {}
97 
104  template<typename InputType>
105  explicit LDLT(const EigenBase<InputType>& matrix)
106  : m_matrix(matrix.rows(), matrix.cols()),
107  m_transpositions(matrix.rows()),
108  m_temporary(matrix.rows()),
109  m_sign(internal::ZeroSign),
110  m_isInitialized(false)
111  {
112  compute(matrix.derived());
113  }
114 
121  template<typename InputType>
122  explicit LDLT(EigenBase<InputType>& matrix)
123  : m_matrix(matrix.derived()),
124  m_transpositions(matrix.rows()),
125  m_temporary(matrix.rows()),
126  m_sign(internal::ZeroSign),
127  m_isInitialized(false)
128  {
129  compute(matrix.derived());
130  }
131 
135  void setZero()
136  {
137  m_isInitialized = false;
138  }
139 
141  inline typename Traits::MatrixU matrixU() const
142  {
143  eigen_assert(m_isInitialized && "LDLT is not initialized.");
144  return Traits::getU(m_matrix);
145  }
146 
148  inline typename Traits::MatrixL matrixL() const
149  {
150  eigen_assert(m_isInitialized && "LDLT is not initialized.");
151  return Traits::getL(m_matrix);
152  }
153 
156  inline const TranspositionType& transpositionsP() const
157  {
158  eigen_assert(m_isInitialized && "LDLT is not initialized.");
159  return m_transpositions;
160  }
161 
164  {
165  eigen_assert(m_isInitialized && "LDLT is not initialized.");
166  return m_matrix.diagonal();
167  }
168 
170  inline bool isPositive() const
171  {
172  eigen_assert(m_isInitialized && "LDLT is not initialized.");
173  return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign;
174  }
175 
177  inline bool isNegative(void) const
178  {
179  eigen_assert(m_isInitialized && "LDLT is not initialized.");
180  return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign;
181  }
182 
198  template<typename Rhs>
199  inline const Solve<LDLT, Rhs>
200  solve(const MatrixBase<Rhs>& b) const
201  {
202  eigen_assert(m_isInitialized && "LDLT is not initialized.");
203  eigen_assert(m_matrix.rows()==b.rows()
204  && "LDLT::solve(): invalid number of rows of the right hand side matrix b");
205  return Solve<LDLT, Rhs>(*this, b.derived());
206  }
207 
208  template<typename Derived>
209  bool solveInPlace(MatrixBase<Derived> &bAndX) const;
210 
211  template<typename InputType>
212  LDLT& compute(const EigenBase<InputType>& matrix);
213 
217  RealScalar rcond() const
218  {
219  eigen_assert(m_isInitialized && "LDLT is not initialized.");
220  return internal::rcond_estimate_helper(m_l1_norm, *this);
221  }
222 
223  template <typename Derived>
224  LDLT& rankUpdate(const MatrixBase<Derived>& w, const RealScalar& alpha=1);
225 
230  inline const MatrixType& matrixLDLT() const
231  {
232  eigen_assert(m_isInitialized && "LDLT is not initialized.");
233  return m_matrix;
234  }
235 
236  MatrixType reconstructedMatrix() const;
237 
243  const LDLT& adjoint() const { return *this; };
244 
245  inline Index rows() const { return m_matrix.rows(); }
246  inline Index cols() const { return m_matrix.cols(); }
247 
254  {
255  eigen_assert(m_isInitialized && "LDLT is not initialized.");
256  return Success;
257  }
258 
259  #ifndef EIGEN_PARSED_BY_DOXYGEN
260  template<typename RhsType, typename DstType>
261  EIGEN_DEVICE_FUNC
262  void _solve_impl(const RhsType &rhs, DstType &dst) const;
263  #endif
264 
265  protected:
266 
267  static void check_template_parameters()
268  {
269  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
270  }
271 
278  MatrixType m_matrix;
279  RealScalar m_l1_norm;
280  TranspositionType m_transpositions;
281  TmpMatrixType m_temporary;
282  internal::SignMatrix m_sign;
283  bool m_isInitialized;
284 };
285 
286 namespace internal {
287 
288 template<int UpLo> struct ldlt_inplace;
289 
290 template<> struct ldlt_inplace<Lower>
291 {
292  template<typename MatrixType, typename TranspositionType, typename Workspace>
293  static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
294  {
295  using std::abs;
296  typedef typename MatrixType::Scalar Scalar;
297  typedef typename MatrixType::RealScalar RealScalar;
298  typedef typename TranspositionType::StorageIndex IndexType;
299  eigen_assert(mat.rows()==mat.cols());
300  const Index size = mat.rows();
301 
302  if (size <= 1)
303  {
304  transpositions.setIdentity();
305  if (numext::real(mat.coeff(0,0)) > static_cast<RealScalar>(0) ) sign = PositiveSemiDef;
306  else if (numext::real(mat.coeff(0,0)) < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
307  else sign = ZeroSign;
308  return true;
309  }
310 
311  for (Index k = 0; k < size; ++k)
312  {
313  // Find largest diagonal element
314  Index index_of_biggest_in_corner;
315  mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
316  index_of_biggest_in_corner += k;
317 
318  transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
319  if(k != index_of_biggest_in_corner)
320  {
321  // apply the transposition while taking care to consider only
322  // the lower triangular part
323  Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
324  mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
325  mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
326  std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
327  for(Index i=k+1;i<index_of_biggest_in_corner;++i)
328  {
329  Scalar tmp = mat.coeffRef(i,k);
330  mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
331  mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
332  }
334  mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
335  }
336 
337  // partition the matrix:
338  // A00 | - | -
339  // lu = A10 | A11 | -
340  // A20 | A21 | A22
341  Index rs = size - k - 1;
342  Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
343  Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
344  Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
345 
346  if(k>0)
347  {
348  temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
349  mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
350  if(rs>0)
351  A21.noalias() -= A20 * temp.head(k);
352  }
353 
354  // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot
355  // was smaller than the cutoff value. However, since LDLT is not rank-revealing
356  // we should only make sure that we do not introduce INF or NaN values.
357  // Remark that LAPACK also uses 0 as the cutoff value.
358  RealScalar realAkk = numext::real(mat.coeffRef(k,k));
359  if((rs>0) && (abs(realAkk) > RealScalar(0)))
360  A21 /= realAkk;
361 
362  if (sign == PositiveSemiDef) {
363  if (realAkk < static_cast<RealScalar>(0)) sign = Indefinite;
364  } else if (sign == NegativeSemiDef) {
365  if (realAkk > static_cast<RealScalar>(0)) sign = Indefinite;
366  } else if (sign == ZeroSign) {
367  if (realAkk > static_cast<RealScalar>(0)) sign = PositiveSemiDef;
368  else if (realAkk < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
369  }
370  }
371 
372  return true;
373  }
374 
375  // Reference for the algorithm: Davis and Hager, "Multiple Rank
376  // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
377  // Trivial rearrangements of their computations (Timothy E. Holy)
378  // allow their algorithm to work for rank-1 updates even if the
379  // original matrix is not of full rank.
380  // Here only rank-1 updates are implemented, to reduce the
381  // requirement for intermediate storage and improve accuracy
382  template<typename MatrixType, typename WDerived>
383  static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, const typename MatrixType::RealScalar& sigma=1)
384  {
385  using numext::isfinite;
386  typedef typename MatrixType::Scalar Scalar;
387  typedef typename MatrixType::RealScalar RealScalar;
388 
389  const Index size = mat.rows();
390  eigen_assert(mat.cols() == size && w.size()==size);
391 
392  RealScalar alpha = 1;
393 
394  // Apply the update
395  for (Index j = 0; j < size; j++)
396  {
397  // Check for termination due to an original decomposition of low-rank
398  if (!(isfinite)(alpha))
399  break;
400 
401  // Update the diagonal terms
402  RealScalar dj = numext::real(mat.coeff(j,j));
403  Scalar wj = w.coeff(j);
404  RealScalar swj2 = sigma*numext::abs2(wj);
405  RealScalar gamma = dj*alpha + swj2;
406 
407  mat.coeffRef(j,j) += swj2/alpha;
408  alpha += swj2/dj;
409 
410 
411  // Update the terms of L
412  Index rs = size-j-1;
413  w.tail(rs) -= wj * mat.col(j).tail(rs);
414  if(gamma != 0)
415  mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
416  }
417  return true;
418  }
419 
420  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
421  static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1)
422  {
423  // Apply the permutation to the input w
424  tmp = transpositions * w;
425 
426  return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
427  }
428 };
429 
430 template<> struct ldlt_inplace<Upper>
431 {
432  template<typename MatrixType, typename TranspositionType, typename Workspace>
433  static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
434  {
435  Transpose<MatrixType> matt(mat);
436  return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
437  }
438 
439  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
440  static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1)
441  {
442  Transpose<MatrixType> matt(mat);
443  return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
444  }
445 };
446 
447 template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
448 {
449  typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
451  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
452  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
453 };
454 
455 template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
456 {
458  typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
459  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); }
460  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); }
461 };
462 
463 } // end namespace internal
464 
467 template<typename MatrixType, int _UpLo>
468 template<typename InputType>
470 {
471  check_template_parameters();
472 
473  eigen_assert(a.rows()==a.cols());
474  const Index size = a.rows();
475 
476  m_matrix = a.derived();
477 
478  // Compute matrix L1 norm = max abs column sum.
479  m_l1_norm = RealScalar(0);
480  // TODO move this code to SelfAdjointView
481  for (Index col = 0; col < size; ++col) {
482  RealScalar abs_col_sum;
483  if (_UpLo == Lower)
484  abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
485  else
486  abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
487  if (abs_col_sum > m_l1_norm)
488  m_l1_norm = abs_col_sum;
489  }
490 
491  m_transpositions.resize(size);
492  m_isInitialized = false;
493  m_temporary.resize(size);
494  m_sign = internal::ZeroSign;
495 
496  internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign);
497 
498  m_isInitialized = true;
499  return *this;
500 }
501 
507 template<typename MatrixType, int _UpLo>
508 template<typename Derived>
509 LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w, const typename LDLT<MatrixType,_UpLo>::RealScalar& sigma)
510 {
511  typedef typename TranspositionType::StorageIndex IndexType;
512  const Index size = w.rows();
513  if (m_isInitialized)
514  {
515  eigen_assert(m_matrix.rows()==size);
516  }
517  else
518  {
519  m_matrix.resize(size,size);
520  m_matrix.setZero();
521  m_transpositions.resize(size);
522  for (Index i = 0; i < size; i++)
523  m_transpositions.coeffRef(i) = IndexType(i);
524  m_temporary.resize(size);
525  m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
526  m_isInitialized = true;
527  }
528 
529  internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
530 
531  return *this;
532 }
533 
534 #ifndef EIGEN_PARSED_BY_DOXYGEN
535 template<typename _MatrixType, int _UpLo>
536 template<typename RhsType, typename DstType>
537 void LDLT<_MatrixType,_UpLo>::_solve_impl(const RhsType &rhs, DstType &dst) const
538 {
539  eigen_assert(rhs.rows() == rows());
540  // dst = P b
541  dst = m_transpositions * rhs;
542 
543  // dst = L^-1 (P b)
544  matrixL().solveInPlace(dst);
545 
546  // dst = D^-1 (L^-1 P b)
547  // more precisely, use pseudo-inverse of D (see bug 241)
548  using std::abs;
549  const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD());
550  // In some previous versions, tolerance was set to the max of 1/highest and the maximal diagonal entry * epsilon
551  // as motivated by LAPACK's xGELSS:
552  // RealScalar tolerance = numext::maxi(vecD.array().abs().maxCoeff() * NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest());
553  // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest
554  // diagonal element is not well justified and leads to numerical issues in some cases.
555  // Moreover, Lapack's xSYTRS routines use 0 for the tolerance.
556  RealScalar tolerance = RealScalar(1) / NumTraits<RealScalar>::highest();
557 
558  for (Index i = 0; i < vecD.size(); ++i)
559  {
560  if(abs(vecD(i)) > tolerance)
561  dst.row(i) /= vecD(i);
562  else
563  dst.row(i).setZero();
564  }
565 
566  // dst = L^-T (D^-1 L^-1 P b)
567  matrixU().solveInPlace(dst);
568 
569  // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
570  dst = m_transpositions.transpose() * dst;
571 }
572 #endif
573 
587 template<typename MatrixType,int _UpLo>
588 template<typename Derived>
590 {
591  eigen_assert(m_isInitialized && "LDLT is not initialized.");
592  eigen_assert(m_matrix.rows() == bAndX.rows());
593 
594  bAndX = this->solve(bAndX);
595 
596  return true;
597 }
598 
602 template<typename MatrixType, int _UpLo>
604 {
605  eigen_assert(m_isInitialized && "LDLT is not initialized.");
606  const Index size = m_matrix.rows();
607  MatrixType res(size,size);
608 
609  // P
610  res.setIdentity();
611  res = transpositionsP() * res;
612  // L^* P
613  res = matrixU() * res;
614  // D(L^*P)
615  res = vectorD().real().asDiagonal() * res;
616  // L(DL^*P)
617  res = matrixL() * res;
618  // P^T (LDL^*P)
619  res = transpositionsP().transpose() * res;
620 
621  return res;
622 }
623 
624 #ifndef __CUDACC__
625 
629 template<typename MatrixType, unsigned int UpLo>
632 {
633  return LDLT<PlainObject,UpLo>(m_matrix);
634 }
635 
640 template<typename Derived>
643 {
644  return LDLT<PlainObject>(derived());
645 }
646 #endif // __CUDACC__
647 
648 } // end namespace Eigen
649 
650 #endif // EIGEN_LDLT_H
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:50
CoeffReturnType coeff(Index row, Index col) const
Definition: DenseCoeffsBase.h:96
const LDLT< PlainObject > ldlt() const
Definition: LDLT.h:642
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_isfinite_op< typename Derived::Scalar >, const Derived > isfinite(const Eigen::ArrayBase< Derived > &x)
Expression of the transpose of a matrix.
Definition: Transpose.h:52
LDLT(Index size)
Default Constructor with memory preallocation.
Definition: LDLT.h:90
Namespace containing all symbols from the Eigen library.
Definition: Core:271
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:167
Derived & derived()
Definition: EigenBase.h:44
void setZero()
Definition: LDLT.h:135
Index rows() const
Definition: EigenBase.h:58
Definition: EigenBase.h:28
const MatrixType & matrixLDLT() const
Definition: LDLT.h:230
LDLT()
Default Constructor.
Definition: LDLT.h:77
Definition: Constants.h:204
LDLT(const EigenBase< InputType > &matrix)
Constructor with decomposition.
Definition: LDLT.h:105
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: XprHelper.h:35
bool isNegative(void) const
Definition: LDLT.h:177
RealScalar rcond() const
Definition: LDLT.h:217
MatrixType reconstructedMatrix() const
Definition: LDLT.h:603
const LDLT & adjoint() const
Definition: LDLT.h:243
const TranspositionType & transpositionsP() const
Definition: LDLT.h:156
Traits::MatrixU matrixU() const
Definition: LDLT.h:141
Definition: Constants.h:432
Eigen::Index Index
Definition: LDLT.h:63
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
const Solve< LDLT, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: LDLT.h:200
Definition: Eigen_Colamd.h:50
bool isPositive() const
Definition: LDLT.h:170
Index cols() const
Definition: EigenBase.h:61
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LDLT.h:253
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:103
Diagonal< const MatrixType > vectorD() const
Definition: LDLT.h:163
Traits::MatrixL matrixL() const
Definition: LDLT.h:148
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:892
LDLT(EigenBase< InputType > &matrix)
Constructs a LDLT factorization from a given matrix.
Definition: LDLT.h:122
ComputationInfo
Definition: Constants.h:430
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sign_op< typename Derived::Scalar >, const Derived > sign(const Eigen::ArrayBase< Derived > &x)
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
const LDLT< PlainObject, UpLo > ldlt() const
Definition: LDLT.h:631
Index size() const
Definition: EigenBase.h:65