LLT.h 18 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #ifndef EIGEN_LLT_H
  10. #define EIGEN_LLT_H
  11. namespace Eigen {
  12. namespace internal{
  13. template<typename MatrixType, int UpLo> struct LLT_Traits;
  14. }
  15. /** \ingroup Cholesky_Module
  16. *
  17. * \class LLT
  18. *
  19. * \brief Standard Cholesky decomposition (LL^T) of a matrix and associated features
  20. *
  21. * \tparam _MatrixType the type of the matrix of which we are computing the LL^T Cholesky decomposition
  22. * \tparam _UpLo the triangular part that will be used for the decompositon: Lower (default) or Upper.
  23. * The other triangular part won't be read.
  24. *
  25. * This class performs a LL^T Cholesky decomposition of a symmetric, positive definite
  26. * matrix A such that A = LL^* = U^*U, where L is lower triangular.
  27. *
  28. * While the Cholesky decomposition is particularly useful to solve selfadjoint problems like D^*D x = b,
  29. * for that purpose, we recommend the Cholesky decomposition without square root which is more stable
  30. * and even faster. Nevertheless, this standard Cholesky decomposition remains useful in many other
  31. * situations like generalised eigen problems with hermitian matrices.
  32. *
  33. * Remember that Cholesky decompositions are not rank-revealing. This LLT decomposition is only stable on positive definite matrices,
  34. * use LDLT instead for the semidefinite case. Also, do not use a Cholesky decomposition to determine whether a system of equations
  35. * has a solution.
  36. *
  37. * Example: \include LLT_example.cpp
  38. * Output: \verbinclude LLT_example.out
  39. *
  40. * \b Performance: for best performance, it is recommended to use a column-major storage format
  41. * with the Lower triangular part (the default), or, equivalently, a row-major storage format
  42. * with the Upper triangular part. Otherwise, you might get a 20% slowdown for the full factorization
  43. * step, and rank-updates can be up to 3 times slower.
  44. *
  45. * This class supports the \link InplaceDecomposition inplace decomposition \endlink mechanism.
  46. *
  47. * Note that during the decomposition, only the lower (or upper, as defined by _UpLo) triangular part of A is considered.
  48. * Therefore, the strict lower part does not have to store correct values.
  49. *
  50. * \sa MatrixBase::llt(), SelfAdjointView::llt(), class LDLT
  51. */
  52. template<typename _MatrixType, int _UpLo> class LLT
  53. {
  54. public:
  55. typedef _MatrixType MatrixType;
  56. enum {
  57. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  58. ColsAtCompileTime = MatrixType::ColsAtCompileTime,
  59. MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
  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. enum {
  66. PacketSize = internal::packet_traits<Scalar>::size,
  67. AlignmentMask = int(PacketSize)-1,
  68. UpLo = _UpLo
  69. };
  70. typedef internal::LLT_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 LLT::compute(const MatrixType&).
  76. */
  77. LLT() : m_matrix(), m_isInitialized(false) {}
  78. /** \brief Default Constructor with memory preallocation
  79. *
  80. * Like the default constructor but with preallocation of the internal data
  81. * according to the specified problem \a size.
  82. * \sa LLT()
  83. */
  84. explicit LLT(Index size) : m_matrix(size, size),
  85. m_isInitialized(false) {}
  86. template<typename InputType>
  87. explicit LLT(const EigenBase<InputType>& matrix)
  88. : m_matrix(matrix.rows(), matrix.cols()),
  89. m_isInitialized(false)
  90. {
  91. compute(matrix.derived());
  92. }
  93. /** \brief Constructs a LDLT factorization from a given matrix
  94. *
  95. * This overloaded constructor is provided for \link InplaceDecomposition inplace decomposition \endlink when
  96. * \c MatrixType is a Eigen::Ref.
  97. *
  98. * \sa LLT(const EigenBase&)
  99. */
  100. template<typename InputType>
  101. explicit LLT(EigenBase<InputType>& matrix)
  102. : m_matrix(matrix.derived()),
  103. m_isInitialized(false)
  104. {
  105. compute(matrix.derived());
  106. }
  107. /** \returns a view of the upper triangular matrix U */
  108. inline typename Traits::MatrixU matrixU() const
  109. {
  110. eigen_assert(m_isInitialized && "LLT is not initialized.");
  111. return Traits::getU(m_matrix);
  112. }
  113. /** \returns a view of the lower triangular matrix L */
  114. inline typename Traits::MatrixL matrixL() const
  115. {
  116. eigen_assert(m_isInitialized && "LLT is not initialized.");
  117. return Traits::getL(m_matrix);
  118. }
  119. /** \returns the solution x of \f$ A x = b \f$ using the current decomposition of A.
  120. *
  121. * Since this LLT class assumes anyway that the matrix A is invertible, the solution
  122. * theoretically exists and is unique regardless of b.
  123. *
  124. * Example: \include LLT_solve.cpp
  125. * Output: \verbinclude LLT_solve.out
  126. *
  127. * \sa solveInPlace(), MatrixBase::llt(), SelfAdjointView::llt()
  128. */
  129. template<typename Rhs>
  130. inline const Solve<LLT, Rhs>
  131. solve(const MatrixBase<Rhs>& b) const
  132. {
  133. eigen_assert(m_isInitialized && "LLT is not initialized.");
  134. eigen_assert(m_matrix.rows()==b.rows()
  135. && "LLT::solve(): invalid number of rows of the right hand side matrix b");
  136. return Solve<LLT, Rhs>(*this, b.derived());
  137. }
  138. template<typename Derived>
  139. void solveInPlace(const MatrixBase<Derived> &bAndX) const;
  140. template<typename InputType>
  141. LLT& compute(const EigenBase<InputType>& matrix);
  142. /** \returns an estimate of the reciprocal condition number of the matrix of
  143. * which \c *this is the Cholesky decomposition.
  144. */
  145. RealScalar rcond() const
  146. {
  147. eigen_assert(m_isInitialized && "LLT is not initialized.");
  148. eigen_assert(m_info == Success && "LLT failed because matrix appears to be negative");
  149. return internal::rcond_estimate_helper(m_l1_norm, *this);
  150. }
  151. /** \returns the LLT decomposition matrix
  152. *
  153. * TODO: document the storage layout
  154. */
  155. inline const MatrixType& matrixLLT() const
  156. {
  157. eigen_assert(m_isInitialized && "LLT is not initialized.");
  158. return m_matrix;
  159. }
  160. MatrixType reconstructedMatrix() const;
  161. /** \brief Reports whether previous computation was successful.
  162. *
  163. * \returns \c Success if computation was succesful,
  164. * \c NumericalIssue if the matrix.appears not to be positive definite.
  165. */
  166. ComputationInfo info() const
  167. {
  168. eigen_assert(m_isInitialized && "LLT is not initialized.");
  169. return m_info;
  170. }
  171. /** \returns the adjoint of \c *this, that is, a const reference to the decomposition itself as the underlying matrix is self-adjoint.
  172. *
  173. * This method is provided for compatibility with other matrix decompositions, thus enabling generic code such as:
  174. * \code x = decomposition.adjoint().solve(b) \endcode
  175. */
  176. const LLT& adjoint() const { return *this; };
  177. inline Index rows() const { return m_matrix.rows(); }
  178. inline Index cols() const { return m_matrix.cols(); }
  179. template<typename VectorType>
  180. LLT rankUpdate(const VectorType& vec, const RealScalar& sigma = 1);
  181. #ifndef EIGEN_PARSED_BY_DOXYGEN
  182. template<typename RhsType, typename DstType>
  183. EIGEN_DEVICE_FUNC
  184. void _solve_impl(const RhsType &rhs, DstType &dst) const;
  185. #endif
  186. protected:
  187. static void check_template_parameters()
  188. {
  189. EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
  190. }
  191. /** \internal
  192. * Used to compute and store L
  193. * The strict upper part is not used and even not initialized.
  194. */
  195. MatrixType m_matrix;
  196. RealScalar m_l1_norm;
  197. bool m_isInitialized;
  198. ComputationInfo m_info;
  199. };
  200. namespace internal {
  201. template<typename Scalar, int UpLo> struct llt_inplace;
  202. template<typename MatrixType, typename VectorType>
  203. static Index llt_rank_update_lower(MatrixType& mat, const VectorType& vec, const typename MatrixType::RealScalar& sigma)
  204. {
  205. using std::sqrt;
  206. typedef typename MatrixType::Scalar Scalar;
  207. typedef typename MatrixType::RealScalar RealScalar;
  208. typedef typename MatrixType::ColXpr ColXpr;
  209. typedef typename internal::remove_all<ColXpr>::type ColXprCleaned;
  210. typedef typename ColXprCleaned::SegmentReturnType ColXprSegment;
  211. typedef Matrix<Scalar,Dynamic,1> TempVectorType;
  212. typedef typename TempVectorType::SegmentReturnType TempVecSegment;
  213. Index n = mat.cols();
  214. eigen_assert(mat.rows()==n && vec.size()==n);
  215. TempVectorType temp;
  216. if(sigma>0)
  217. {
  218. // This version is based on Givens rotations.
  219. // It is faster than the other one below, but only works for updates,
  220. // i.e., for sigma > 0
  221. temp = sqrt(sigma) * vec;
  222. for(Index i=0; i<n; ++i)
  223. {
  224. JacobiRotation<Scalar> g;
  225. g.makeGivens(mat(i,i), -temp(i), &mat(i,i));
  226. Index rs = n-i-1;
  227. if(rs>0)
  228. {
  229. ColXprSegment x(mat.col(i).tail(rs));
  230. TempVecSegment y(temp.tail(rs));
  231. apply_rotation_in_the_plane(x, y, g);
  232. }
  233. }
  234. }
  235. else
  236. {
  237. temp = vec;
  238. RealScalar beta = 1;
  239. for(Index j=0; j<n; ++j)
  240. {
  241. RealScalar Ljj = numext::real(mat.coeff(j,j));
  242. RealScalar dj = numext::abs2(Ljj);
  243. Scalar wj = temp.coeff(j);
  244. RealScalar swj2 = sigma*numext::abs2(wj);
  245. RealScalar gamma = dj*beta + swj2;
  246. RealScalar x = dj + swj2/beta;
  247. if (x<=RealScalar(0))
  248. return j;
  249. RealScalar nLjj = sqrt(x);
  250. mat.coeffRef(j,j) = nLjj;
  251. beta += swj2/dj;
  252. // Update the terms of L
  253. Index rs = n-j-1;
  254. if(rs)
  255. {
  256. temp.tail(rs) -= (wj/Ljj) * mat.col(j).tail(rs);
  257. if(gamma != 0)
  258. mat.col(j).tail(rs) = (nLjj/Ljj) * mat.col(j).tail(rs) + (nLjj * sigma*numext::conj(wj)/gamma)*temp.tail(rs);
  259. }
  260. }
  261. }
  262. return -1;
  263. }
  264. template<typename Scalar> struct llt_inplace<Scalar, Lower>
  265. {
  266. typedef typename NumTraits<Scalar>::Real RealScalar;
  267. template<typename MatrixType>
  268. static Index unblocked(MatrixType& mat)
  269. {
  270. using std::sqrt;
  271. eigen_assert(mat.rows()==mat.cols());
  272. const Index size = mat.rows();
  273. for(Index k = 0; k < size; ++k)
  274. {
  275. Index rs = size-k-1; // remaining size
  276. Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
  277. Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
  278. Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
  279. RealScalar x = numext::real(mat.coeff(k,k));
  280. if (k>0) x -= A10.squaredNorm();
  281. if (x<=RealScalar(0))
  282. return k;
  283. mat.coeffRef(k,k) = x = sqrt(x);
  284. if (k>0 && rs>0) A21.noalias() -= A20 * A10.adjoint();
  285. if (rs>0) A21 /= x;
  286. }
  287. return -1;
  288. }
  289. template<typename MatrixType>
  290. static Index blocked(MatrixType& m)
  291. {
  292. eigen_assert(m.rows()==m.cols());
  293. Index size = m.rows();
  294. if(size<32)
  295. return unblocked(m);
  296. Index blockSize = size/8;
  297. blockSize = (blockSize/16)*16;
  298. blockSize = (std::min)((std::max)(blockSize,Index(8)), Index(128));
  299. for (Index k=0; k<size; k+=blockSize)
  300. {
  301. // partition the matrix:
  302. // A00 | - | -
  303. // lu = A10 | A11 | -
  304. // A20 | A21 | A22
  305. Index bs = (std::min)(blockSize, size-k);
  306. Index rs = size - k - bs;
  307. Block<MatrixType,Dynamic,Dynamic> A11(m,k, k, bs,bs);
  308. Block<MatrixType,Dynamic,Dynamic> A21(m,k+bs,k, rs,bs);
  309. Block<MatrixType,Dynamic,Dynamic> A22(m,k+bs,k+bs,rs,rs);
  310. Index ret;
  311. if((ret=unblocked(A11))>=0) return k+ret;
  312. if(rs>0) A11.adjoint().template triangularView<Upper>().template solveInPlace<OnTheRight>(A21);
  313. if(rs>0) A22.template selfadjointView<Lower>().rankUpdate(A21,typename NumTraits<RealScalar>::Literal(-1)); // bottleneck
  314. }
  315. return -1;
  316. }
  317. template<typename MatrixType, typename VectorType>
  318. static Index rankUpdate(MatrixType& mat, const VectorType& vec, const RealScalar& sigma)
  319. {
  320. return Eigen::internal::llt_rank_update_lower(mat, vec, sigma);
  321. }
  322. };
  323. template<typename Scalar> struct llt_inplace<Scalar, Upper>
  324. {
  325. typedef typename NumTraits<Scalar>::Real RealScalar;
  326. template<typename MatrixType>
  327. static EIGEN_STRONG_INLINE Index unblocked(MatrixType& mat)
  328. {
  329. Transpose<MatrixType> matt(mat);
  330. return llt_inplace<Scalar, Lower>::unblocked(matt);
  331. }
  332. template<typename MatrixType>
  333. static EIGEN_STRONG_INLINE Index blocked(MatrixType& mat)
  334. {
  335. Transpose<MatrixType> matt(mat);
  336. return llt_inplace<Scalar, Lower>::blocked(matt);
  337. }
  338. template<typename MatrixType, typename VectorType>
  339. static Index rankUpdate(MatrixType& mat, const VectorType& vec, const RealScalar& sigma)
  340. {
  341. Transpose<MatrixType> matt(mat);
  342. return llt_inplace<Scalar, Lower>::rankUpdate(matt, vec.conjugate(), sigma);
  343. }
  344. };
  345. template<typename MatrixType> struct LLT_Traits<MatrixType,Lower>
  346. {
  347. typedef const TriangularView<const MatrixType, Lower> MatrixL;
  348. typedef const TriangularView<const typename MatrixType::AdjointReturnType, Upper> MatrixU;
  349. static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
  350. static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
  351. static bool inplace_decomposition(MatrixType& m)
  352. { return llt_inplace<typename MatrixType::Scalar, Lower>::blocked(m)==-1; }
  353. };
  354. template<typename MatrixType> struct LLT_Traits<MatrixType,Upper>
  355. {
  356. typedef const TriangularView<const typename MatrixType::AdjointReturnType, Lower> MatrixL;
  357. typedef const TriangularView<const MatrixType, Upper> MatrixU;
  358. static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); }
  359. static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); }
  360. static bool inplace_decomposition(MatrixType& m)
  361. { return llt_inplace<typename MatrixType::Scalar, Upper>::blocked(m)==-1; }
  362. };
  363. } // end namespace internal
  364. /** Computes / recomputes the Cholesky decomposition A = LL^* = U^*U of \a matrix
  365. *
  366. * \returns a reference to *this
  367. *
  368. * Example: \include TutorialLinAlgComputeTwice.cpp
  369. * Output: \verbinclude TutorialLinAlgComputeTwice.out
  370. */
  371. template<typename MatrixType, int _UpLo>
  372. template<typename InputType>
  373. LLT<MatrixType,_UpLo>& LLT<MatrixType,_UpLo>::compute(const EigenBase<InputType>& a)
  374. {
  375. check_template_parameters();
  376. eigen_assert(a.rows()==a.cols());
  377. const Index size = a.rows();
  378. m_matrix.resize(size, size);
  379. if (!internal::is_same_dense(m_matrix, a.derived()))
  380. m_matrix = a.derived();
  381. // Compute matrix L1 norm = max abs column sum.
  382. m_l1_norm = RealScalar(0);
  383. // TODO move this code to SelfAdjointView
  384. for (Index col = 0; col < size; ++col) {
  385. RealScalar abs_col_sum;
  386. if (_UpLo == Lower)
  387. abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
  388. else
  389. abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
  390. if (abs_col_sum > m_l1_norm)
  391. m_l1_norm = abs_col_sum;
  392. }
  393. m_isInitialized = true;
  394. bool ok = Traits::inplace_decomposition(m_matrix);
  395. m_info = ok ? Success : NumericalIssue;
  396. return *this;
  397. }
  398. /** Performs a rank one update (or dowdate) of the current decomposition.
  399. * If A = LL^* before the rank one update,
  400. * then after it we have LL^* = A + sigma * v v^* where \a v must be a vector
  401. * of same dimension.
  402. */
  403. template<typename _MatrixType, int _UpLo>
  404. template<typename VectorType>
  405. LLT<_MatrixType,_UpLo> LLT<_MatrixType,_UpLo>::rankUpdate(const VectorType& v, const RealScalar& sigma)
  406. {
  407. EIGEN_STATIC_ASSERT_VECTOR_ONLY(VectorType);
  408. eigen_assert(v.size()==m_matrix.cols());
  409. eigen_assert(m_isInitialized);
  410. if(internal::llt_inplace<typename MatrixType::Scalar, UpLo>::rankUpdate(m_matrix,v,sigma)>=0)
  411. m_info = NumericalIssue;
  412. else
  413. m_info = Success;
  414. return *this;
  415. }
  416. #ifndef EIGEN_PARSED_BY_DOXYGEN
  417. template<typename _MatrixType,int _UpLo>
  418. template<typename RhsType, typename DstType>
  419. void LLT<_MatrixType,_UpLo>::_solve_impl(const RhsType &rhs, DstType &dst) const
  420. {
  421. dst = rhs;
  422. solveInPlace(dst);
  423. }
  424. #endif
  425. /** \internal use x = llt_object.solve(x);
  426. *
  427. * This is the \em in-place version of solve().
  428. *
  429. * \param bAndX represents both the right-hand side matrix b and result x.
  430. *
  431. * This version avoids a copy when the right hand side matrix b is not needed anymore.
  432. *
  433. * \warning The parameter is only marked 'const' to make the C++ compiler accept a temporary expression here.
  434. * This function will const_cast it, so constness isn't honored here.
  435. *
  436. * \sa LLT::solve(), MatrixBase::llt()
  437. */
  438. template<typename MatrixType, int _UpLo>
  439. template<typename Derived>
  440. void LLT<MatrixType,_UpLo>::solveInPlace(const MatrixBase<Derived> &bAndX) const
  441. {
  442. eigen_assert(m_isInitialized && "LLT is not initialized.");
  443. eigen_assert(m_matrix.rows()==bAndX.rows());
  444. matrixL().solveInPlace(bAndX);
  445. matrixU().solveInPlace(bAndX);
  446. }
  447. /** \returns the matrix represented by the decomposition,
  448. * i.e., it returns the product: L L^*.
  449. * This function is provided for debug purpose. */
  450. template<typename MatrixType, int _UpLo>
  451. MatrixType LLT<MatrixType,_UpLo>::reconstructedMatrix() const
  452. {
  453. eigen_assert(m_isInitialized && "LLT is not initialized.");
  454. return matrixL() * matrixL().adjoint().toDenseMatrix();
  455. }
  456. /** \cholesky_module
  457. * \returns the LLT decomposition of \c *this
  458. * \sa SelfAdjointView::llt()
  459. */
  460. template<typename Derived>
  461. inline const LLT<typename MatrixBase<Derived>::PlainObject>
  462. MatrixBase<Derived>::llt() const
  463. {
  464. return LLT<PlainObject>(derived());
  465. }
  466. /** \cholesky_module
  467. * \returns the LLT decomposition of \c *this
  468. * \sa SelfAdjointView::llt()
  469. */
  470. template<typename MatrixType, unsigned int UpLo>
  471. inline const LLT<typename SelfAdjointView<MatrixType, UpLo>::PlainObject, UpLo>
  472. SelfAdjointView<MatrixType, UpLo>::llt() const
  473. {
  474. return LLT<PlainObject,UpLo>(m_matrix);
  475. }
  476. } // end namespace Eigen
  477. #endif // EIGEN_LLT_H