LDLT.h 24 KB

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