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
25/** \ingroup Cholesky_Module
26 *
27 * \class LDLT
28 *
29 * \brief Robust Cholesky decomposition of a matrix with pivoting
30 *
31 * \tparam _MatrixType the type of the matrix of which to compute the LDL^T Cholesky decomposition
32 * \tparam _UpLo the triangular part that will be used for the decompositon: Lower (default) or Upper.
33 * The other triangular part won't be read.
34 *
35 * Perform a robust Cholesky decomposition of a positive semidefinite or negative semidefinite
36 * matrix \f$ A \f$ such that \f$ A = P^TLDL^*P \f$, where P is a permutation matrix, L
37 * is lower triangular with a unit diagonal and D is a diagonal matrix.
38 *
39 * The decomposition uses pivoting to ensure stability, so that L will have
40 * zeros in the bottom right rank(A) - n submatrix. Avoiding the square root
41 * on D also stabilizes the computation.
42 *
43 * Remember that Cholesky decompositions are not rank-revealing. Also, do not use a Cholesky
44 * decomposition to determine whether a system of equations has a solution.
45 *
46 * This class supports the \link InplaceDecomposition inplace decomposition \endlink mechanism.
47 *
48 * \sa MatrixBase::ldlt(), SelfAdjointView::ldlt(), class LLT
49 */
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;
63 typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
64 typedef typename MatrixType::StorageIndex StorageIndex;
65 typedef Matrix<Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1> TmpMatrixType;
66
67 typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType;
68 typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType;
69
70 typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
71
72 /** \brief Default Constructor.
73 *
74 * The default constructor is useful in cases in which the user intends to
75 * perform decompositions via LDLT::compute(const MatrixType&).
76 */
77 LDLT()
78 : m_matrix(),
79 m_transpositions(),
80 m_sign(internal::ZeroSign),
81 m_isInitialized(false)
82 {}
83
84 /** \brief Default Constructor with memory preallocation
85 *
86 * Like the default constructor but with preallocation of the internal data
87 * according to the specified problem \a size.
88 * \sa LDLT()
89 */
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
98 /** \brief Constructor with decomposition
99 *
100 * This calculates the decomposition for the input \a matrix.
101 *
102 * \sa LDLT(Index size)
103 */
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
115 /** \brief Constructs a LDLT factorization from a given matrix
116 *
117 * This overloaded constructor is provided for \link InplaceDecomposition inplace decomposition \endlink when \c MatrixType is a Eigen::Ref.
118 *
119 * \sa LDLT(const EigenBase&)
120 */
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
132 /** Clear any existing decomposition
133 * \sa rankUpdate(w,sigma)
134 */
135 void setZero()
136 {
137 m_isInitialized = false;
138 }
139
140 /** \returns a view of the upper triangular matrix U */
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
147 /** \returns a view of the lower triangular matrix L */
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
154 /** \returns the permutation matrix P as a transposition sequence.
155 */
156 inline const TranspositionType& transpositionsP() const
157 {
158 eigen_assert(m_isInitialized && "LDLT is not initialized.");
159 return m_transpositions;
160 }
161
162 /** \returns the coefficients of the diagonal matrix D */
163 inline Diagonal<const MatrixType> vectorD() const
164 {
165 eigen_assert(m_isInitialized && "LDLT is not initialized.");
166 return m_matrix.diagonal();
167 }
168
169 /** \returns true if the matrix is positive (semidefinite) */
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
176 /** \returns true if the matrix is negative (semidefinite) */
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
183 /** \returns a solution x of \f$ A x = b \f$ using the current decomposition of A.
184 *
185 * This function also supports in-place solves using the syntax <tt>x = decompositionObject.solve(x)</tt> .
186 *
187 * \note_about_checking_solutions
188 *
189 * More precisely, this method solves \f$ A x = b \f$ using the decomposition \f$ A = P^T L D L^* P \f$
190 * by solving the systems \f$ P^T y_1 = b \f$, \f$ L y_2 = y_1 \f$, \f$ D y_3 = y_2 \f$,
191 * \f$ L^* y_4 = y_3 \f$ and \f$ P x = y_4 \f$ in succession. If the matrix \f$ A \f$ is singular, then
192 * \f$ D \f$ will also be singular (all the other matrices are invertible). In that case, the
193 * least-square solution of \f$ D y_3 = y_2 \f$ is computed. This does not mean that this function
194 * computes the least-square solution of \f$ A x = b \f$ is \f$ A \f$ is singular.
195 *
196 * \sa MatrixBase::ldlt(), SelfAdjointView::ldlt()
197 */
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
214 /** \returns an estimate of the reciprocal condition number of the matrix of
215 * which \c *this is the LDLT decomposition.
216 */
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
226 /** \returns the internal LDLT decomposition matrix
227 *
228 * TODO: document the storage layout
229 */
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
238 /** \returns the adjoint of \c *this, that is, a const reference to the decomposition itself as the underlying matrix is self-adjoint.
239 *
240 * This method is provided for compatibility with other matrix decompositions, thus enabling generic code such as:
241 * \code x = decomposition.adjoint().solve(b) \endcode
242 */
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
248 /** \brief Reports whether previous computation was successful.
249 *
250 * \returns \c Success if computation was succesful,
251 * \c NumericalIssue if the factorization failed because of a zero pivot.
252 */
253 ComputationInfo info() const
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
272 /** \internal
273 * Used to compute and store the Cholesky decomposition A = L D L^* = U^* D U.
274 * The strict upper part is used during the decomposition, the strict lower
275 * part correspond to the coefficients of L (its diagonal is equal to 1 and
276 * is not stored), and the diagonal entries correspond to D.
277 */
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(size==0) sign = ZeroSign;
309 else if (numext::real(mat.coeff(0,0)) > static_cast<RealScalar>(0) ) sign = PositiveSemiDef;
310 else if (numext::real(mat.coeff(0,0)) < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
311 else sign = ZeroSign;
312 return true;
313 }
314
315 for (Index k = 0; k < size; ++k)
316 {
317 // Find largest diagonal element
318 Index index_of_biggest_in_corner;
319 mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
320 index_of_biggest_in_corner += k;
321
322 transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
323 if(k != index_of_biggest_in_corner)
324 {
325 // apply the transposition while taking care to consider only
326 // the lower triangular part
327 Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
328 mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
329 mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
330 std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
331 for(Index i=k+1;i<index_of_biggest_in_corner;++i)
332 {
333 Scalar tmp = mat.coeffRef(i,k);
334 mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
335 mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
336 }
337 if(NumTraits<Scalar>::IsComplex)
338 mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
339 }
340
341 // partition the matrix:
342 // A00 | - | -
343 // lu = A10 | A11 | -
344 // A20 | A21 | A22
345 Index rs = size - k - 1;
346 Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
347 Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
348 Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
349
350 if(k>0)
351 {
352 temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
353 mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
354 if(rs>0)
355 A21.noalias() -= A20 * temp.head(k);
356 }
357
358 // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot
359 // was smaller than the cutoff value. However, since LDLT is not rank-revealing
360 // we should only make sure that we do not introduce INF or NaN values.
361 // Remark that LAPACK also uses 0 as the cutoff value.
362 RealScalar realAkk = numext::real(mat.coeffRef(k,k));
363 bool pivot_is_valid = (abs(realAkk) > RealScalar(0));
364
365 if(k==0 && !pivot_is_valid)
366 {
367 // The entire diagonal is zero, there is nothing more to do
368 // except filling the transpositions, and checking whether the matrix is zero.
369 sign = ZeroSign;
370 for(Index j = 0; j<size; ++j)
371 {
372 transpositions.coeffRef(j) = IndexType(j);
373 ret = ret && (mat.col(j).tail(size-j-1).array()==Scalar(0)).all();
374 }
375 return ret;
376 }
377
378 if((rs>0) && pivot_is_valid)
379 A21 /= realAkk;
380 else if(rs>0)
381 ret = ret && (A21.array()==Scalar(0)).all();
382
383 if(found_zero_pivot && pivot_is_valid) ret = false; // factorization failed
384 else if(!pivot_is_valid) found_zero_pivot = true;
385
386 if (sign == PositiveSemiDef) {
387 if (realAkk < static_cast<RealScalar>(0)) sign = Indefinite;
388 } else if (sign == NegativeSemiDef) {
389 if (realAkk > static_cast<RealScalar>(0)) sign = Indefinite;
390 } else if (sign == ZeroSign) {
391 if (realAkk > static_cast<RealScalar>(0)) sign = PositiveSemiDef;
392 else if (realAkk < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
393 }
394 }
395
396 return ret;
397 }
398
399 // Reference for the algorithm: Davis and Hager, "Multiple Rank
400 // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
401 // Trivial rearrangements of their computations (Timothy E. Holy)
402 // allow their algorithm to work for rank-1 updates even if the
403 // original matrix is not of full rank.
404 // Here only rank-1 updates are implemented, to reduce the
405 // requirement for intermediate storage and improve accuracy
406 template<typename MatrixType, typename WDerived>
407 static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, const typename MatrixType::RealScalar& sigma=1)
408 {
409 using numext::isfinite;
410 typedef typename MatrixType::Scalar Scalar;
411 typedef typename MatrixType::RealScalar RealScalar;
412
413 const Index size = mat.rows();
414 eigen_assert(mat.cols() == size && w.size()==size);
415
416 RealScalar alpha = 1;
417
418 // Apply the update
419 for (Index j = 0; j < size; j++)
420 {
421 // Check for termination due to an original decomposition of low-rank
422 if (!(isfinite)(alpha))
423 break;
424
425 // Update the diagonal terms
426 RealScalar dj = numext::real(mat.coeff(j,j));
427 Scalar wj = w.coeff(j);
428 RealScalar swj2 = sigma*numext::abs2(wj);
429 RealScalar gamma = dj*alpha + swj2;
430
431 mat.coeffRef(j,j) += swj2/alpha;
432 alpha += swj2/dj;
433
434
435 // Update the terms of L
436 Index rs = size-j-1;
437 w.tail(rs) -= wj * mat.col(j).tail(rs);
438 if(gamma != 0)
439 mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
440 }
441 return true;
442 }
443
444 template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
445 static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1)
446 {
447 // Apply the permutation to the input w
448 tmp = transpositions * w;
449
450 return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
451 }
452};
453
454template<> struct ldlt_inplace<Upper>
455{
456 template<typename MatrixType, typename TranspositionType, typename Workspace>
457 static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
458 {
459 Transpose<MatrixType> matt(mat);
460 return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
461 }
462
463 template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
464 static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1)
465 {
466 Transpose<MatrixType> matt(mat);
467 return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
468 }
469};
470
471template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
472{
473 typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
474 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
475 static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
476 static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
477};
478
479template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
480{
481 typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
482 typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
483 static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); }
484 static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); }
485};
486
487} // end namespace internal
488
489/** Compute / recompute the LDLT decomposition A = L D L^* = U^* D U of \a matrix
490 */
491template<typename MatrixType, int _UpLo>
492template<typename InputType>
493LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::compute(const EigenBase<InputType>& a)
494{
495 check_template_parameters();
496
497 eigen_assert(a.rows()==a.cols());
498 const Index size = a.rows();
499
500 m_matrix = a.derived();
501
502 // Compute matrix L1 norm = max abs column sum.
503 m_l1_norm = RealScalar(0);
504 // TODO move this code to SelfAdjointView
505 for (Index col = 0; col < size; ++col) {
506 RealScalar abs_col_sum;
507 if (_UpLo == Lower)
508 abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
509 else
510 abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
511 if (abs_col_sum > m_l1_norm)
512 m_l1_norm = abs_col_sum;
513 }
514
515 m_transpositions.resize(size);
516 m_isInitialized = false;
517 m_temporary.resize(size);
518 m_sign = internal::ZeroSign;
519
520 m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ? Success : NumericalIssue;
521
522 m_isInitialized = true;
523 return *this;
524}
525
526/** Update the LDLT decomposition: given A = L D L^T, efficiently compute the decomposition of A + sigma w w^T.
527 * \param w a vector to be incorporated into the decomposition.
528 * \param sigma a scalar, +1 for updates and -1 for "downdates," which correspond to removing previously-added column vectors. Optional; default value is +1.
529 * \sa setZero()
530 */
531template<typename MatrixType, int _UpLo>
532template<typename Derived>
533LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w, const typename LDLT<MatrixType,_UpLo>::RealScalar& sigma)
534{
535 typedef typename TranspositionType::StorageIndex IndexType;
536 const Index size = w.rows();
537 if (m_isInitialized)
538 {
539 eigen_assert(m_matrix.rows()==size);
540 }
541 else
542 {
543 m_matrix.resize(size,size);
544 m_matrix.setZero();
545 m_transpositions.resize(size);
546 for (Index i = 0; i < size; i++)
547 m_transpositions.coeffRef(i) = IndexType(i);
548 m_temporary.resize(size);
549 m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
550 m_isInitialized = true;
551 }
552
553 internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
554
555 return *this;
556}
557
558#ifndef EIGEN_PARSED_BY_DOXYGEN
559template<typename _MatrixType, int _UpLo>
560template<typename RhsType, typename DstType>
561void LDLT<_MatrixType,_UpLo>::_solve_impl(const RhsType &rhs, DstType &dst) const
562{
563 eigen_assert(rhs.rows() == rows());
564 // dst = P b
565 dst = m_transpositions * rhs;
566
567 // dst = L^-1 (P b)
568 matrixL().solveInPlace(dst);
569
570 // dst = D^-1 (L^-1 P b)
571 // more precisely, use pseudo-inverse of D (see bug 241)
572 using std::abs;
573 const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD());
574 // In some previous versions, tolerance was set to the max of 1/highest (or rather numeric_limits::min())
575 // and the maximal diagonal entry * epsilon as motivated by LAPACK's xGELSS:
576 // RealScalar tolerance = numext::maxi(vecD.array().abs().maxCoeff() * NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest());
577 // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest
578 // diagonal element is not well justified and leads to numerical issues in some cases.
579 // Moreover, Lapack's xSYTRS routines use 0 for the tolerance.
580 // Using numeric_limits::min() gives us more robustness to denormals.
581 RealScalar tolerance = (std::numeric_limits<RealScalar>::min)();
582
583 for (Index i = 0; i < vecD.size(); ++i)
584 {
585 if(abs(vecD(i)) > tolerance)
586 dst.row(i) /= vecD(i);
587 else
588 dst.row(i).setZero();
589 }
590
591 // dst = L^-T (D^-1 L^-1 P b)
592 matrixU().solveInPlace(dst);
593
594 // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
595 dst = m_transpositions.transpose() * dst;
596}
597#endif
598
599/** \internal use x = ldlt_object.solve(x);
600 *
601 * This is the \em in-place version of solve().
602 *
603 * \param bAndX represents both the right-hand side matrix b and result x.
604 *
605 * \returns true always! If you need to check for existence of solutions, use another decomposition like LU, QR, or SVD.
606 *
607 * This version avoids a copy when the right hand side matrix b is not
608 * needed anymore.
609 *
610 * \sa LDLT::solve(), MatrixBase::ldlt()
611 */
612template<typename MatrixType,int _UpLo>
613template<typename Derived>
614bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
615{
616 eigen_assert(m_isInitialized && "LDLT is not initialized.");
617 eigen_assert(m_matrix.rows() == bAndX.rows());
618
619 bAndX = this->solve(bAndX);
620
621 return true;
622}
623
624/** \returns the matrix represented by the decomposition,
625 * i.e., it returns the product: P^T L D L^* P.
626 * This function is provided for debug purpose. */
627template<typename MatrixType, int _UpLo>
628MatrixType LDLT<MatrixType,_UpLo>::reconstructedMatrix() const
629{
630 eigen_assert(m_isInitialized && "LDLT is not initialized.");
631 const Index size = m_matrix.rows();
632 MatrixType res(size,size);
633
634 // P
635 res.setIdentity();
636 res = transpositionsP() * res;
637 // L^* P
638 res = matrixU() * res;
639 // D(L^*P)
640 res = vectorD().real().asDiagonal() * res;
641 // L(DL^*P)
642 res = matrixL() * res;
643 // P^T (LDL^*P)
644 res = transpositionsP().transpose() * res;
645
646 return res;
647}
648
649/** \cholesky_module
650 * \returns the Cholesky decomposition with full pivoting without square root of \c *this
651 * \sa MatrixBase::ldlt()
652 */
653template<typename MatrixType, unsigned int UpLo>
654inline const LDLT<typename SelfAdjointView<MatrixType, UpLo>::PlainObject, UpLo>
655SelfAdjointView<MatrixType, UpLo>::ldlt() const
656{
657 return LDLT<PlainObject,UpLo>(m_matrix);
658}
659
660/** \cholesky_module
661 * \returns the Cholesky decomposition with full pivoting without square root of \c *this
662 * \sa SelfAdjointView::ldlt()
663 */
664template<typename Derived>
665inline const LDLT<typename MatrixBase<Derived>::PlainObject>
666MatrixBase<Derived>::ldlt() const
667{
668 return LDLT<PlainObject>(derived());
669}
670
671} // end namespace Eigen
672
673#endif // EIGEN_LDLT_H
674