LDLT.h 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673
  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. #ifndef EIGEN_LDLT_H
  13. #define EIGEN_LDLT_H
  14. namespace Eigen {
  15. namespace internal {
  16. template<typename MatrixType, int UpLo> struct LDLT_Traits;
  17. // PositiveSemiDef means positive semi-definite and non-zero; same for NegativeSemiDef
  18. enum SignMatrix { PositiveSemiDef, NegativeSemiDef, ZeroSign, Indefinite };
  19. }
  20. /** \ingroup Cholesky_Module
  21. *
  22. * \class LDLT
  23. *
  24. * \brief Robust Cholesky decomposition of a matrix with pivoting
  25. *
  26. * \tparam _MatrixType the type of the matrix of which to compute the LDL^T Cholesky decomposition
  27. * \tparam _UpLo the triangular part that will be used for the decompositon: Lower (default) or Upper.
  28. * The other triangular part won't be read.
  29. *
  30. * Perform a robust Cholesky decomposition of a positive semidefinite or negative semidefinite
  31. * matrix \f$ A \f$ such that \f$ A = P^TLDL^*P \f$, where P is a permutation matrix, L
  32. * is lower triangular with a unit diagonal and D is a diagonal matrix.
  33. *
  34. * The decomposition uses pivoting to ensure stability, so that L will have
  35. * zeros in the bottom right rank(A) - n submatrix. Avoiding the square root
  36. * on D also stabilizes the computation.
  37. *
  38. * Remember that Cholesky decompositions are not rank-revealing. Also, do not use a Cholesky
  39. * decomposition to determine whether a system of equations has a solution.
  40. *
  41. * This class supports the \link InplaceDecomposition inplace decomposition \endlink mechanism.
  42. *
  43. * \sa MatrixBase::ldlt(), SelfAdjointView::ldlt(), class LLT
  44. */
  45. template<typename _MatrixType, int _UpLo> class LDLT
  46. {
  47. public:
  48. typedef _MatrixType MatrixType;
  49. enum {
  50. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  51. ColsAtCompileTime = MatrixType::ColsAtCompileTime,
  52. MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  53. MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
  54. UpLo = _UpLo
  55. };
  56. typedef typename MatrixType::Scalar Scalar;
  57. typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
  58. typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
  59. typedef typename MatrixType::StorageIndex StorageIndex;
  60. typedef Matrix<Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1> TmpMatrixType;
  61. typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType;
  62. typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType;
  63. typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
  64. /** \brief Default Constructor.
  65. *
  66. * The default constructor is useful in cases in which the user intends to
  67. * perform decompositions via LDLT::compute(const MatrixType&).
  68. */
  69. LDLT()
  70. : m_matrix(),
  71. m_transpositions(),
  72. m_sign(internal::ZeroSign),
  73. m_isInitialized(false)
  74. {}
  75. /** \brief Default Constructor with memory preallocation
  76. *
  77. * Like the default constructor but with preallocation of the internal data
  78. * according to the specified problem \a size.
  79. * \sa LDLT()
  80. */
  81. explicit LDLT(Index size)
  82. : m_matrix(size, size),
  83. m_transpositions(size),
  84. m_temporary(size),
  85. m_sign(internal::ZeroSign),
  86. m_isInitialized(false)
  87. {}
  88. /** \brief Constructor with decomposition
  89. *
  90. * This calculates the decomposition for the input \a matrix.
  91. *
  92. * \sa LDLT(Index size)
  93. */
  94. template<typename InputType>
  95. explicit LDLT(const EigenBase<InputType>& matrix)
  96. : m_matrix(matrix.rows(), matrix.cols()),
  97. m_transpositions(matrix.rows()),
  98. m_temporary(matrix.rows()),
  99. m_sign(internal::ZeroSign),
  100. m_isInitialized(false)
  101. {
  102. compute(matrix.derived());
  103. }
  104. /** \brief Constructs a LDLT factorization from a given matrix
  105. *
  106. * This overloaded constructor is provided for \link InplaceDecomposition inplace decomposition \endlink when \c MatrixType is a Eigen::Ref.
  107. *
  108. * \sa LDLT(const EigenBase&)
  109. */
  110. template<typename InputType>
  111. explicit LDLT(EigenBase<InputType>& matrix)
  112. : m_matrix(matrix.derived()),
  113. m_transpositions(matrix.rows()),
  114. m_temporary(matrix.rows()),
  115. m_sign(internal::ZeroSign),
  116. m_isInitialized(false)
  117. {
  118. compute(matrix.derived());
  119. }
  120. /** Clear any existing decomposition
  121. * \sa rankUpdate(w,sigma)
  122. */
  123. void setZero()
  124. {
  125. m_isInitialized = false;
  126. }
  127. /** \returns a view of the upper triangular matrix U */
  128. inline typename Traits::MatrixU matrixU() const
  129. {
  130. eigen_assert(m_isInitialized && "LDLT is not initialized.");
  131. return Traits::getU(m_matrix);
  132. }
  133. /** \returns a view of the lower triangular matrix L */
  134. inline typename Traits::MatrixL matrixL() const
  135. {
  136. eigen_assert(m_isInitialized && "LDLT is not initialized.");
  137. return Traits::getL(m_matrix);
  138. }
  139. /** \returns the permutation matrix P as a transposition sequence.
  140. */
  141. inline const TranspositionType& transpositionsP() const
  142. {
  143. eigen_assert(m_isInitialized && "LDLT is not initialized.");
  144. return m_transpositions;
  145. }
  146. /** \returns the coefficients of the diagonal matrix D */
  147. inline Diagonal<const MatrixType> vectorD() const
  148. {
  149. eigen_assert(m_isInitialized && "LDLT is not initialized.");
  150. return m_matrix.diagonal();
  151. }
  152. /** \returns true if the matrix is positive (semidefinite) */
  153. inline bool isPositive() const
  154. {
  155. eigen_assert(m_isInitialized && "LDLT is not initialized.");
  156. return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign;
  157. }
  158. /** \returns true if the matrix is negative (semidefinite) */
  159. inline bool isNegative(void) const
  160. {
  161. eigen_assert(m_isInitialized && "LDLT is not initialized.");
  162. return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign;
  163. }
  164. /** \returns a solution x of \f$ A x = b \f$ using the current decomposition of A.
  165. *
  166. * This function also supports in-place solves using the syntax <tt>x = decompositionObject.solve(x)</tt> .
  167. *
  168. * \note_about_checking_solutions
  169. *
  170. * More precisely, this method solves \f$ A x = b \f$ using the decomposition \f$ A = P^T L D L^* P \f$
  171. * 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$,
  172. * \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
  173. * \f$ D \f$ will also be singular (all the other matrices are invertible). In that case, the
  174. * least-square solution of \f$ D y_3 = y_2 \f$ is computed. This does not mean that this function
  175. * computes the least-square solution of \f$ A x = b \f$ is \f$ A \f$ is singular.
  176. *
  177. * \sa MatrixBase::ldlt(), SelfAdjointView::ldlt()
  178. */
  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. template<typename Derived>
  189. bool solveInPlace(MatrixBase<Derived> &bAndX) const;
  190. template<typename InputType>
  191. LDLT& compute(const EigenBase<InputType>& matrix);
  192. /** \returns an estimate of the reciprocal condition number of the matrix of
  193. * which \c *this is the LDLT decomposition.
  194. */
  195. RealScalar rcond() const
  196. {
  197. eigen_assert(m_isInitialized && "LDLT is not initialized.");
  198. return internal::rcond_estimate_helper(m_l1_norm, *this);
  199. }
  200. template <typename Derived>
  201. LDLT& rankUpdate(const MatrixBase<Derived>& w, const RealScalar& alpha=1);
  202. /** \returns the internal LDLT decomposition matrix
  203. *
  204. * TODO: document the storage layout
  205. */
  206. inline const MatrixType& matrixLDLT() const
  207. {
  208. eigen_assert(m_isInitialized && "LDLT is not initialized.");
  209. return m_matrix;
  210. }
  211. MatrixType reconstructedMatrix() const;
  212. /** \returns the adjoint of \c *this, that is, a const reference to the decomposition itself as the underlying matrix is self-adjoint.
  213. *
  214. * This method is provided for compatibility with other matrix decompositions, thus enabling generic code such as:
  215. * \code x = decomposition.adjoint().solve(b) \endcode
  216. */
  217. const LDLT& adjoint() const { return *this; };
  218. inline Index rows() const { return m_matrix.rows(); }
  219. inline Index cols() const { return m_matrix.cols(); }
  220. /** \brief Reports whether previous computation was successful.
  221. *
  222. * \returns \c Success if computation was succesful,
  223. * \c NumericalIssue if the factorization failed because of a zero pivot.
  224. */
  225. ComputationInfo info() const
  226. {
  227. eigen_assert(m_isInitialized && "LDLT is not initialized.");
  228. return m_info;
  229. }
  230. #ifndef EIGEN_PARSED_BY_DOXYGEN
  231. template<typename RhsType, typename DstType>
  232. EIGEN_DEVICE_FUNC
  233. void _solve_impl(const RhsType &rhs, DstType &dst) const;
  234. #endif
  235. protected:
  236. static void check_template_parameters()
  237. {
  238. EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
  239. }
  240. /** \internal
  241. * Used to compute and store the Cholesky decomposition A = L D L^* = U^* D U.
  242. * The strict upper part is used during the decomposition, the strict lower
  243. * part correspond to the coefficients of L (its diagonal is equal to 1 and
  244. * is not stored), and the diagonal entries correspond to D.
  245. */
  246. MatrixType m_matrix;
  247. RealScalar m_l1_norm;
  248. TranspositionType m_transpositions;
  249. TmpMatrixType m_temporary;
  250. internal::SignMatrix m_sign;
  251. bool m_isInitialized;
  252. ComputationInfo m_info;
  253. };
  254. namespace internal {
  255. template<int UpLo> struct ldlt_inplace;
  256. template<> struct ldlt_inplace<Lower>
  257. {
  258. template<typename MatrixType, typename TranspositionType, typename Workspace>
  259. static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
  260. {
  261. using std::abs;
  262. typedef typename MatrixType::Scalar Scalar;
  263. typedef typename MatrixType::RealScalar RealScalar;
  264. typedef typename TranspositionType::StorageIndex IndexType;
  265. eigen_assert(mat.rows()==mat.cols());
  266. const Index size = mat.rows();
  267. bool found_zero_pivot = false;
  268. bool ret = true;
  269. if (size <= 1)
  270. {
  271. transpositions.setIdentity();
  272. if(size==0) sign = ZeroSign;
  273. else if (numext::real(mat.coeff(0,0)) > static_cast<RealScalar>(0) ) sign = PositiveSemiDef;
  274. else if (numext::real(mat.coeff(0,0)) < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
  275. else sign = ZeroSign;
  276. return true;
  277. }
  278. for (Index k = 0; k < size; ++k)
  279. {
  280. // Find largest diagonal element
  281. Index index_of_biggest_in_corner;
  282. mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
  283. index_of_biggest_in_corner += k;
  284. transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
  285. if(k != index_of_biggest_in_corner)
  286. {
  287. // apply the transposition while taking care to consider only
  288. // the lower triangular part
  289. Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
  290. mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
  291. mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
  292. std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
  293. for(Index i=k+1;i<index_of_biggest_in_corner;++i)
  294. {
  295. Scalar tmp = mat.coeffRef(i,k);
  296. mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
  297. mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
  298. }
  299. if(NumTraits<Scalar>::IsComplex)
  300. mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
  301. }
  302. // partition the matrix:
  303. // A00 | - | -
  304. // lu = A10 | A11 | -
  305. // A20 | A21 | A22
  306. Index rs = size - k - 1;
  307. Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
  308. Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
  309. Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
  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. // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot
  318. // was smaller than the cutoff value. However, since LDLT is not rank-revealing
  319. // we should only make sure that we do not introduce INF or NaN values.
  320. // Remark that LAPACK also uses 0 as the cutoff value.
  321. RealScalar realAkk = numext::real(mat.coeffRef(k,k));
  322. bool pivot_is_valid = (abs(realAkk) > RealScalar(0));
  323. if(k==0 && !pivot_is_valid)
  324. {
  325. // The entire diagonal is zero, there is nothing more to do
  326. // except filling the transpositions, and checking whether the matrix is zero.
  327. sign = ZeroSign;
  328. for(Index j = 0; j<size; ++j)
  329. {
  330. transpositions.coeffRef(j) = IndexType(j);
  331. ret = ret && (mat.col(j).tail(size-j-1).array()==Scalar(0)).all();
  332. }
  333. return ret;
  334. }
  335. if((rs>0) && pivot_is_valid)
  336. A21 /= realAkk;
  337. else if(rs>0)
  338. ret = ret && (A21.array()==Scalar(0)).all();
  339. if(found_zero_pivot && pivot_is_valid) ret = false; // factorization failed
  340. else if(!pivot_is_valid) found_zero_pivot = true;
  341. if (sign == PositiveSemiDef) {
  342. if (realAkk < static_cast<RealScalar>(0)) sign = Indefinite;
  343. } else if (sign == NegativeSemiDef) {
  344. if (realAkk > static_cast<RealScalar>(0)) sign = Indefinite;
  345. } else if (sign == ZeroSign) {
  346. if (realAkk > static_cast<RealScalar>(0)) sign = PositiveSemiDef;
  347. else if (realAkk < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
  348. }
  349. }
  350. return ret;
  351. }
  352. // Reference for the algorithm: Davis and Hager, "Multiple Rank
  353. // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
  354. // Trivial rearrangements of their computations (Timothy E. Holy)
  355. // allow their algorithm to work for rank-1 updates even if the
  356. // original matrix is not of full rank.
  357. // Here only rank-1 updates are implemented, to reduce the
  358. // requirement for intermediate storage and improve accuracy
  359. template<typename MatrixType, typename WDerived>
  360. static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, const typename MatrixType::RealScalar& sigma=1)
  361. {
  362. using numext::isfinite;
  363. typedef typename MatrixType::Scalar Scalar;
  364. typedef typename MatrixType::RealScalar RealScalar;
  365. const Index size = mat.rows();
  366. eigen_assert(mat.cols() == size && w.size()==size);
  367. RealScalar alpha = 1;
  368. // Apply the update
  369. for (Index j = 0; j < size; j++)
  370. {
  371. // Check for termination due to an original decomposition of low-rank
  372. if (!(isfinite)(alpha))
  373. break;
  374. // Update the diagonal terms
  375. RealScalar dj = numext::real(mat.coeff(j,j));
  376. Scalar wj = w.coeff(j);
  377. RealScalar swj2 = sigma*numext::abs2(wj);
  378. RealScalar gamma = dj*alpha + swj2;
  379. mat.coeffRef(j,j) += swj2/alpha;
  380. alpha += swj2/dj;
  381. // Update the terms of L
  382. Index rs = size-j-1;
  383. w.tail(rs) -= wj * mat.col(j).tail(rs);
  384. if(gamma != 0)
  385. mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
  386. }
  387. return true;
  388. }
  389. template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
  390. static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1)
  391. {
  392. // Apply the permutation to the input w
  393. tmp = transpositions * w;
  394. return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
  395. }
  396. };
  397. template<> struct ldlt_inplace<Upper>
  398. {
  399. template<typename MatrixType, typename TranspositionType, typename Workspace>
  400. static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
  401. {
  402. Transpose<MatrixType> matt(mat);
  403. return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
  404. }
  405. template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
  406. static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1)
  407. {
  408. Transpose<MatrixType> matt(mat);
  409. return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
  410. }
  411. };
  412. template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
  413. {
  414. typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
  415. typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
  416. static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
  417. static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
  418. };
  419. template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
  420. {
  421. typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
  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. } // end namespace internal
  427. /** Compute / recompute the LDLT decomposition A = L D L^* = U^* D U of \a matrix
  428. */
  429. template<typename MatrixType, int _UpLo>
  430. template<typename InputType>
  431. LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::compute(const EigenBase<InputType>& a)
  432. {
  433. check_template_parameters();
  434. eigen_assert(a.rows()==a.cols());
  435. const Index size = a.rows();
  436. m_matrix = a.derived();
  437. // Compute matrix L1 norm = max abs column sum.
  438. m_l1_norm = RealScalar(0);
  439. // TODO move this code to SelfAdjointView
  440. for (Index col = 0; col < size; ++col) {
  441. RealScalar abs_col_sum;
  442. if (_UpLo == Lower)
  443. abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
  444. else
  445. abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
  446. if (abs_col_sum > m_l1_norm)
  447. m_l1_norm = abs_col_sum;
  448. }
  449. m_transpositions.resize(size);
  450. m_isInitialized = false;
  451. m_temporary.resize(size);
  452. m_sign = internal::ZeroSign;
  453. m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ? Success : NumericalIssue;
  454. m_isInitialized = true;
  455. return *this;
  456. }
  457. /** Update the LDLT decomposition: given A = L D L^T, efficiently compute the decomposition of A + sigma w w^T.
  458. * \param w a vector to be incorporated into the decomposition.
  459. * \param sigma a scalar, +1 for updates and -1 for "downdates," which correspond to removing previously-added column vectors. Optional; default value is +1.
  460. * \sa setZero()
  461. */
  462. template<typename MatrixType, int _UpLo>
  463. template<typename Derived>
  464. LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w, const typename LDLT<MatrixType,_UpLo>::RealScalar& sigma)
  465. {
  466. typedef typename TranspositionType::StorageIndex IndexType;
  467. const Index size = w.rows();
  468. if (m_isInitialized)
  469. {
  470. eigen_assert(m_matrix.rows()==size);
  471. }
  472. else
  473. {
  474. m_matrix.resize(size,size);
  475. m_matrix.setZero();
  476. m_transpositions.resize(size);
  477. for (Index i = 0; i < size; i++)
  478. m_transpositions.coeffRef(i) = IndexType(i);
  479. m_temporary.resize(size);
  480. m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
  481. m_isInitialized = true;
  482. }
  483. internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
  484. return *this;
  485. }
  486. #ifndef EIGEN_PARSED_BY_DOXYGEN
  487. template<typename _MatrixType, int _UpLo>
  488. template<typename RhsType, typename DstType>
  489. void LDLT<_MatrixType,_UpLo>::_solve_impl(const RhsType &rhs, DstType &dst) const
  490. {
  491. eigen_assert(rhs.rows() == rows());
  492. // dst = P b
  493. dst = m_transpositions * rhs;
  494. // dst = L^-1 (P b)
  495. matrixL().solveInPlace(dst);
  496. // dst = D^-1 (L^-1 P b)
  497. // more precisely, use pseudo-inverse of D (see bug 241)
  498. using std::abs;
  499. const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD());
  500. // In some previous versions, tolerance was set to the max of 1/highest (or rather numeric_limits::min())
  501. // and the maximal diagonal entry * epsilon as motivated by LAPACK's xGELSS:
  502. // RealScalar tolerance = numext::maxi(vecD.array().abs().maxCoeff() * NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest());
  503. // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest
  504. // diagonal element is not well justified and leads to numerical issues in some cases.
  505. // Moreover, Lapack's xSYTRS routines use 0 for the tolerance.
  506. // Using numeric_limits::min() gives us more robustness to denormals.
  507. RealScalar tolerance = (std::numeric_limits<RealScalar>::min)();
  508. for (Index i = 0; i < vecD.size(); ++i)
  509. {
  510. if(abs(vecD(i)) > tolerance)
  511. dst.row(i) /= vecD(i);
  512. else
  513. dst.row(i).setZero();
  514. }
  515. // dst = L^-T (D^-1 L^-1 P b)
  516. matrixU().solveInPlace(dst);
  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. /** \internal use x = ldlt_object.solve(x);
  522. *
  523. * This is the \em in-place version of solve().
  524. *
  525. * \param bAndX represents both the right-hand side matrix b and result x.
  526. *
  527. * \returns true always! If you need to check for existence of solutions, use another decomposition like LU, QR, or SVD.
  528. *
  529. * This version avoids a copy when the right hand side matrix b is not
  530. * needed anymore.
  531. *
  532. * \sa LDLT::solve(), MatrixBase::ldlt()
  533. */
  534. template<typename MatrixType,int _UpLo>
  535. template<typename Derived>
  536. bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
  537. {
  538. eigen_assert(m_isInitialized && "LDLT is not initialized.");
  539. eigen_assert(m_matrix.rows() == bAndX.rows());
  540. bAndX = this->solve(bAndX);
  541. return true;
  542. }
  543. /** \returns the matrix represented by the decomposition,
  544. * i.e., it returns the product: P^T L D L^* P.
  545. * This function is provided for debug purpose. */
  546. template<typename MatrixType, int _UpLo>
  547. MatrixType LDLT<MatrixType,_UpLo>::reconstructedMatrix() const
  548. {
  549. eigen_assert(m_isInitialized && "LDLT is not initialized.");
  550. const Index size = m_matrix.rows();
  551. MatrixType res(size,size);
  552. // P
  553. res.setIdentity();
  554. res = transpositionsP() * res;
  555. // L^* P
  556. res = matrixU() * res;
  557. // D(L^*P)
  558. res = vectorD().real().asDiagonal() * res;
  559. // L(DL^*P)
  560. res = matrixL() * res;
  561. // P^T (LDL^*P)
  562. res = transpositionsP().transpose() * res;
  563. return res;
  564. }
  565. /** \cholesky_module
  566. * \returns the Cholesky decomposition with full pivoting without square root of \c *this
  567. * \sa MatrixBase::ldlt()
  568. */
  569. template<typename MatrixType, unsigned int UpLo>
  570. inline const LDLT<typename SelfAdjointView<MatrixType, UpLo>::PlainObject, UpLo>
  571. SelfAdjointView<MatrixType, UpLo>::ldlt() const
  572. {
  573. return LDLT<PlainObject,UpLo>(m_matrix);
  574. }
  575. /** \cholesky_module
  576. * \returns the Cholesky decomposition with full pivoting without square root of \c *this
  577. * \sa SelfAdjointView::ldlt()
  578. */
  579. template<typename Derived>
  580. inline const LDLT<typename MatrixBase<Derived>::PlainObject>
  581. MatrixBase<Derived>::ldlt() const
  582. {
  583. return LDLT<PlainObject>(derived());
  584. }
  585. } // end namespace Eigen
  586. #endif // EIGEN_LDLT_H