Eigen  3.3.0
 
Loading...
Searching...
No Matches
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
16namespace Eigen {
17
18namespace 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
50template<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;
64 typedef typename MatrixType::StorageIndex StorageIndex;
66
69
70 typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
71
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 m_info;
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 ComputationInfo m_info;
285};
286
287namespace internal {
288
289template<int UpLo> struct ldlt_inplace;
290
291template<> struct ldlt_inplace<Lower>
292{
293 template<typename MatrixType, typename TranspositionType, typename Workspace>
294 static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
295 {
296 using std::abs;
297 typedef typename MatrixType::Scalar Scalar;
298 typedef typename MatrixType::RealScalar RealScalar;
299 typedef typename TranspositionType::StorageIndex IndexType;
300 eigen_assert(mat.rows()==mat.cols());
301 const Index size = mat.rows();
302 bool found_zero_pivot = false;
303 bool ret = true;
304
305 if (size <= 1)
306 {
307 transpositions.setIdentity();
308 if (numext::real(mat.coeff(0,0)) > static_cast<RealScalar>(0) ) sign = PositiveSemiDef;
309 else if (numext::real(mat.coeff(0,0)) < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
310 else sign = ZeroSign;
311 return true;
312 }
313
314 for (Index k = 0; k < size; ++k)
315 {
316 // Find largest diagonal element
317 Index index_of_biggest_in_corner;
318 mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
319 index_of_biggest_in_corner += k;
320
321 transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
322 if(k != index_of_biggest_in_corner)
323 {
324 // apply the transposition while taking care to consider only
325 // the lower triangular part
326 Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
327 mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
328 mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
329 std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
330 for(Index i=k+1;i<index_of_biggest_in_corner;++i)
331 {
332 Scalar tmp = mat.coeffRef(i,k);
333 mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
334 mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
335 }
336 if(NumTraits<Scalar>::IsComplex)
337 mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
338 }
339
340 // partition the matrix:
341 // A00 | - | -
342 // lu = A10 | A11 | -
343 // A20 | A21 | A22
344 Index rs = size - k - 1;
345 Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
346 Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
347 Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
348
349 if(k>0)
350 {
351 temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
352 mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
353 if(rs>0)
354 A21.noalias() -= A20 * temp.head(k);
355 }
356
357 // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot
358 // was smaller than the cutoff value. However, since LDLT is not rank-revealing
359 // we should only make sure that we do not introduce INF or NaN values.
360 // Remark that LAPACK also uses 0 as the cutoff value.
361 RealScalar realAkk = numext::real(mat.coeffRef(k,k));
362 bool pivot_is_valid = (abs(realAkk) > RealScalar(0));
363
364 if(k==0 && !pivot_is_valid)
365 {
366 // The entire diagonal is zero, there is nothing more to do
367 // except filling the transpositions, and checking whether the matrix is zero.
368 sign = ZeroSign;
369 for(Index j = 0; j<size; ++j)
370 {
371 transpositions.coeffRef(j) = IndexType(j);
372 ret = ret && (mat.col(j).tail(size-j-1).array()==Scalar(0)).all();
373 }
374 return ret;
375 }
376
377 if((rs>0) && pivot_is_valid)
378 A21 /= realAkk;
379
380 if(found_zero_pivot && pivot_is_valid) ret = false; // factorization failed
381 else if(!pivot_is_valid) found_zero_pivot = true;
382
383 if (sign == PositiveSemiDef) {
384 if (realAkk < static_cast<RealScalar>(0)) sign = Indefinite;
385 } else if (sign == NegativeSemiDef) {
386 if (realAkk > static_cast<RealScalar>(0)) sign = Indefinite;
387 } else if (sign == ZeroSign) {
388 if (realAkk > static_cast<RealScalar>(0)) sign = PositiveSemiDef;
389 else if (realAkk < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
390 }
391 }
392
393 return ret;
394 }
395
396 // Reference for the algorithm: Davis and Hager, "Multiple Rank
397 // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
398 // Trivial rearrangements of their computations (Timothy E. Holy)
399 // allow their algorithm to work for rank-1 updates even if the
400 // original matrix is not of full rank.
401 // Here only rank-1 updates are implemented, to reduce the
402 // requirement for intermediate storage and improve accuracy
403 template<typename MatrixType, typename WDerived>
404 static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, const typename MatrixType::RealScalar& sigma=1)
405 {
406 using numext::isfinite;
407 typedef typename MatrixType::Scalar Scalar;
408 typedef typename MatrixType::RealScalar RealScalar;
409
410 const Index size = mat.rows();
411 eigen_assert(mat.cols() == size && w.size()==size);
412
413 RealScalar alpha = 1;
414
415 // Apply the update
416 for (Index j = 0; j < size; j++)
417 {
418 // Check for termination due to an original decomposition of low-rank
419 if (!(isfinite)(alpha))
420 break;
421
422 // Update the diagonal terms
423 RealScalar dj = numext::real(mat.coeff(j,j));
424 Scalar wj = w.coeff(j);
425 RealScalar swj2 = sigma*numext::abs2(wj);
426 RealScalar gamma = dj*alpha + swj2;
427
428 mat.coeffRef(j,j) += swj2/alpha;
429 alpha += swj2/dj;
430
431
432 // Update the terms of L
433 Index rs = size-j-1;
434 w.tail(rs) -= wj * mat.col(j).tail(rs);
435 if(gamma != 0)
436 mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
437 }
438 return true;
439 }
440
441 template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
442 static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1)
443 {
444 // Apply the permutation to the input w
445 tmp = transpositions * w;
446
447 return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
448 }
449};
450
451template<> struct ldlt_inplace<Upper>
452{
453 template<typename MatrixType, typename TranspositionType, typename Workspace>
454 static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
455 {
456 Transpose<MatrixType> matt(mat);
457 return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
458 }
459
460 template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
461 static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1)
462 {
463 Transpose<MatrixType> matt(mat);
464 return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
465 }
466};
467
468template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
469{
470 typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
471 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
472 static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
473 static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
474};
475
476template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
477{
478 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
479 typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
480 static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); }
481 static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); }
482};
483
484} // end namespace internal
485
488template<typename MatrixType, int _UpLo>
489template<typename InputType>
491{
492 check_template_parameters();
493
494 eigen_assert(a.rows()==a.cols());
495 const Index size = a.rows();
496
497 m_matrix = a.derived();
498
499 // Compute matrix L1 norm = max abs column sum.
500 m_l1_norm = RealScalar(0);
501 // TODO move this code to SelfAdjointView
502 for (Index col = 0; col < size; ++col) {
503 RealScalar abs_col_sum;
504 if (_UpLo == Lower)
505 abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
506 else
507 abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
508 if (abs_col_sum > m_l1_norm)
509 m_l1_norm = abs_col_sum;
510 }
511
512 m_transpositions.resize(size);
513 m_isInitialized = false;
514 m_temporary.resize(size);
515 m_sign = internal::ZeroSign;
516
517 m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ? Success : NumericalIssue;
518
519 m_isInitialized = true;
520 return *this;
521}
522
528template<typename MatrixType, int _UpLo>
529template<typename Derived>
530LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w, const typename LDLT<MatrixType,_UpLo>::RealScalar& sigma)
531{
532 typedef typename TranspositionType::StorageIndex IndexType;
533 const Index size = w.rows();
534 if (m_isInitialized)
535 {
536 eigen_assert(m_matrix.rows()==size);
537 }
538 else
539 {
540 m_matrix.resize(size,size);
541 m_matrix.setZero();
542 m_transpositions.resize(size);
543 for (Index i = 0; i < size; i++)
544 m_transpositions.coeffRef(i) = IndexType(i);
545 m_temporary.resize(size);
546 m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
547 m_isInitialized = true;
548 }
549
550 internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
551
552 return *this;
553}
554
555#ifndef EIGEN_PARSED_BY_DOXYGEN
556template<typename _MatrixType, int _UpLo>
557template<typename RhsType, typename DstType>
558void LDLT<_MatrixType,_UpLo>::_solve_impl(const RhsType &rhs, DstType &dst) const
559{
560 eigen_assert(rhs.rows() == rows());
561 // dst = P b
562 dst = m_transpositions * rhs;
563
564 // dst = L^-1 (P b)
565 matrixL().solveInPlace(dst);
566
567 // dst = D^-1 (L^-1 P b)
568 // more precisely, use pseudo-inverse of D (see bug 241)
569 using std::abs;
570 const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD());
571 // In some previous versions, tolerance was set to the max of 1/highest and the maximal diagonal entry * epsilon
572 // as motivated by LAPACK's xGELSS:
573 // RealScalar tolerance = numext::maxi(vecD.array().abs().maxCoeff() * NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest());
574 // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest
575 // diagonal element is not well justified and leads to numerical issues in some cases.
576 // Moreover, Lapack's xSYTRS routines use 0 for the tolerance.
577 RealScalar tolerance = RealScalar(1) / NumTraits<RealScalar>::highest();
578
579 for (Index i = 0; i < vecD.size(); ++i)
580 {
581 if(abs(vecD(i)) > tolerance)
582 dst.row(i) /= vecD(i);
583 else
584 dst.row(i).setZero();
585 }
586
587 // dst = L^-T (D^-1 L^-1 P b)
588 matrixU().solveInPlace(dst);
589
590 // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
591 dst = m_transpositions.transpose() * dst;
592}
593#endif
594
608template<typename MatrixType,int _UpLo>
609template<typename Derived>
610bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
611{
612 eigen_assert(m_isInitialized && "LDLT is not initialized.");
613 eigen_assert(m_matrix.rows() == bAndX.rows());
614
615 bAndX = this->solve(bAndX);
616
617 return true;
618}
619
623template<typename MatrixType, int _UpLo>
625{
626 eigen_assert(m_isInitialized && "LDLT is not initialized.");
627 const Index size = m_matrix.rows();
628 MatrixType res(size,size);
629
630 // P
631 res.setIdentity();
632 res = transpositionsP() * res;
633 // L^* P
634 res = matrixU() * res;
635 // D(L^*P)
636 res = vectorD().real().asDiagonal() * res;
637 // L(DL^*P)
638 res = matrixL() * res;
639 // P^T (LDL^*P)
640 res = transpositionsP().transpose() * res;
641
642 return res;
643}
644
649template<typename MatrixType, unsigned int UpLo>
652{
653 return LDLT<PlainObject,UpLo>(m_matrix);
654}
655
660template<typename Derived>
663{
664 return LDLT<PlainObject>(derived());
665}
666
667} // end namespace Eigen
668
669#endif // EIGEN_LDLT_H
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:66
Derived & derived()
Definition: EigenBase.h:44
Index rows() const
Definition: EigenBase.h:58
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:65
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:51
const Solve< LDLT, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: LDLT.h:200
LDLT(Index size)
Default Constructor with memory preallocation.
Definition: LDLT.h:90
LDLT()
Default Constructor.
Definition: LDLT.h:77
Traits::MatrixU matrixU() const
Definition: LDLT.h:141
bool isPositive() const
Definition: LDLT.h:170
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LDLT.h:253
void setZero()
Definition: LDLT.h:135
const MatrixType & matrixLDLT() const
Definition: LDLT.h:230
bool isNegative(void) const
Definition: LDLT.h:177
Diagonal< const MatrixType > vectorD() const
Definition: LDLT.h:163
LDLT(const EigenBase< InputType > &matrix)
Constructor with decomposition.
Definition: LDLT.h:105
Eigen::Index Index
Definition: LDLT.h:63
LDLT(EigenBase< InputType > &matrix)
Constructs a LDLT factorization from a given matrix.
Definition: LDLT.h:122
const LDLT & adjoint() const
Definition: LDLT.h:243
MatrixType reconstructedMatrix() const
Definition: LDLT.h:624
RealScalar rcond() const
Definition: LDLT.h:217
Traits::MatrixL matrixL() const
Definition: LDLT.h:148
const TranspositionType & transpositionsP() const
Definition: LDLT.h:156
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:50
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:180
Permutation matrix.
Definition: PermutationMatrix.h:309
Expression of a selfadjoint matrix from a triangular part of a dense matrix.
Definition: SelfAdjointView.h:51
Pseudo expression representing a solving operation.
Definition: Solve.h:63
Represents a sequence of transpositions (row/column interchange)
Definition: Transpositions.h:159
ComputationInfo
Definition: Constants.h:430
@ Lower
Definition: Constants.h:204
@ Upper
Definition: Constants.h:206
@ NumericalIssue
Definition: Constants.h:434
@ Success
Definition: Constants.h:432
Namespace containing all symbols from the Eigen library.
Definition: Core:287
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sign_op< typename Derived::Scalar >, const Derived > sign(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_isfinite_op< typename Derived::Scalar >, const Derived > isfinite(const Eigen::ArrayBase< Derived > &x)
Definition: EigenBase.h:29
Index cols() const
Definition: EigenBase.h:61
Derived & derived()
Definition: EigenBase.h:44
Index rows() const
Definition: EigenBase.h:58
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:151