IncompleteCholesky.h 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
  5. // Copyright (C) 2015 Gael Guennebaud <gael.guennebaud@inria.fr>
  6. //
  7. // This Source Code Form is subject to the terms of the Mozilla
  8. // Public License v. 2.0. If a copy of the MPL was not distributed
  9. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  10. #ifndef EIGEN_INCOMPLETE_CHOlESKY_H
  11. #define EIGEN_INCOMPLETE_CHOlESKY_H
  12. #include <vector>
  13. #include <list>
  14. namespace Eigen {
  15. /**
  16. * \brief Modified Incomplete Cholesky with dual threshold
  17. *
  18. * References : C-J. Lin and J. J. Moré, Incomplete Cholesky Factorizations with
  19. * Limited memory, SIAM J. Sci. Comput. 21(1), pp. 24-45, 1999
  20. *
  21. * \tparam Scalar the scalar type of the input matrices
  22. * \tparam _UpLo The triangular part that will be used for the computations. It can be Lower
  23. * or Upper. Default is Lower.
  24. * \tparam _OrderingType The ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<int>,
  25. * unless EIGEN_MPL2_ONLY is defined, in which case the default is NaturalOrdering<int>.
  26. *
  27. * \implsparsesolverconcept
  28. *
  29. * It performs the following incomplete factorization: \f$ S P A P' S \approx L L' \f$
  30. * where L is a lower triangular factor, S is a diagonal scaling matrix, and P is a
  31. * fill-in reducing permutation as computed by the ordering method.
  32. *
  33. * \b Shifting \b strategy: Let \f$ B = S P A P' S \f$ be the scaled matrix on which the factorization is carried out,
  34. * and \f$ \beta \f$ be the minimum value of the diagonal. If \f$ \beta > 0 \f$ then, the factorization is directly performed
  35. * on the matrix B. Otherwise, the factorization is performed on the shifted matrix \f$ B + (\sigma+|\beta| I \f$ where
  36. * \f$ \sigma \f$ is the initial shift value as returned and set by setInitialShift() method. The default value is \f$ \sigma = 10^{-3} \f$.
  37. * If the factorization fails, then the shift in doubled until it succeed or a maximum of ten attempts. If it still fails, as returned by
  38. * the info() method, then you can either increase the initial shift, or better use another preconditioning technique.
  39. *
  40. */
  41. template <typename Scalar, int _UpLo = Lower, typename _OrderingType = AMDOrdering<int> >
  42. class IncompleteCholesky : public SparseSolverBase<IncompleteCholesky<Scalar,_UpLo,_OrderingType> >
  43. {
  44. protected:
  45. typedef SparseSolverBase<IncompleteCholesky<Scalar,_UpLo,_OrderingType> > Base;
  46. using Base::m_isInitialized;
  47. public:
  48. typedef typename NumTraits<Scalar>::Real RealScalar;
  49. typedef _OrderingType OrderingType;
  50. typedef typename OrderingType::PermutationType PermutationType;
  51. typedef typename PermutationType::StorageIndex StorageIndex;
  52. typedef SparseMatrix<Scalar,ColMajor,StorageIndex> FactorType;
  53. typedef Matrix<Scalar,Dynamic,1> VectorSx;
  54. typedef Matrix<RealScalar,Dynamic,1> VectorRx;
  55. typedef Matrix<StorageIndex,Dynamic, 1> VectorIx;
  56. typedef std::vector<std::list<StorageIndex> > VectorList;
  57. enum { UpLo = _UpLo };
  58. enum {
  59. ColsAtCompileTime = Dynamic,
  60. MaxColsAtCompileTime = Dynamic
  61. };
  62. public:
  63. /** Default constructor leaving the object in a partly non-initialized stage.
  64. *
  65. * You must call compute() or the pair analyzePattern()/factorize() to make it valid.
  66. *
  67. * \sa IncompleteCholesky(const MatrixType&)
  68. */
  69. IncompleteCholesky() : m_initialShift(1e-3),m_analysisIsOk(false),m_factorizationIsOk(false) {}
  70. /** Constructor computing the incomplete factorization for the given matrix \a matrix.
  71. */
  72. template<typename MatrixType>
  73. IncompleteCholesky(const MatrixType& matrix) : m_initialShift(1e-3),m_analysisIsOk(false),m_factorizationIsOk(false)
  74. {
  75. compute(matrix);
  76. }
  77. /** \returns number of rows of the factored matrix */
  78. EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT { return m_L.rows(); }
  79. /** \returns number of columns of the factored matrix */
  80. EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT { return m_L.cols(); }
  81. /** \brief Reports whether previous computation was successful.
  82. *
  83. * It triggers an assertion if \c *this has not been initialized through the respective constructor,
  84. * or a call to compute() or analyzePattern().
  85. *
  86. * \returns \c Success if computation was successful,
  87. * \c NumericalIssue if the matrix appears to be negative.
  88. */
  89. ComputationInfo info() const
  90. {
  91. eigen_assert(m_isInitialized && "IncompleteCholesky is not initialized.");
  92. return m_info;
  93. }
  94. /** \brief Set the initial shift parameter \f$ \sigma \f$.
  95. */
  96. void setInitialShift(RealScalar shift) { m_initialShift = shift; }
  97. /** \brief Computes the fill reducing permutation vector using the sparsity pattern of \a mat
  98. */
  99. template<typename MatrixType>
  100. void analyzePattern(const MatrixType& mat)
  101. {
  102. OrderingType ord;
  103. PermutationType pinv;
  104. ord(mat.template selfadjointView<UpLo>(), pinv);
  105. if(pinv.size()>0) m_perm = pinv.inverse();
  106. else m_perm.resize(0);
  107. m_L.resize(mat.rows(), mat.cols());
  108. m_analysisIsOk = true;
  109. m_isInitialized = true;
  110. m_info = Success;
  111. }
  112. /** \brief Performs the numerical factorization of the input matrix \a mat
  113. *
  114. * The method analyzePattern() or compute() must have been called beforehand
  115. * with a matrix having the same pattern.
  116. *
  117. * \sa compute(), analyzePattern()
  118. */
  119. template<typename MatrixType>
  120. void factorize(const MatrixType& mat);
  121. /** Computes or re-computes the incomplete Cholesky factorization of the input matrix \a mat
  122. *
  123. * It is a shortcut for a sequential call to the analyzePattern() and factorize() methods.
  124. *
  125. * \sa analyzePattern(), factorize()
  126. */
  127. template<typename MatrixType>
  128. void compute(const MatrixType& mat)
  129. {
  130. analyzePattern(mat);
  131. factorize(mat);
  132. }
  133. // internal
  134. template<typename Rhs, typename Dest>
  135. void _solve_impl(const Rhs& b, Dest& x) const
  136. {
  137. eigen_assert(m_factorizationIsOk && "factorize() should be called first");
  138. if (m_perm.rows() == b.rows()) x = m_perm * b;
  139. else x = b;
  140. x = m_scale.asDiagonal() * x;
  141. x = m_L.template triangularView<Lower>().solve(x);
  142. x = m_L.adjoint().template triangularView<Upper>().solve(x);
  143. x = m_scale.asDiagonal() * x;
  144. if (m_perm.rows() == b.rows())
  145. x = m_perm.inverse() * x;
  146. }
  147. /** \returns the sparse lower triangular factor L */
  148. const FactorType& matrixL() const { eigen_assert("m_factorizationIsOk"); return m_L; }
  149. /** \returns a vector representing the scaling factor S */
  150. const VectorRx& scalingS() const { eigen_assert("m_factorizationIsOk"); return m_scale; }
  151. /** \returns the fill-in reducing permutation P (can be empty for a natural ordering) */
  152. const PermutationType& permutationP() const { eigen_assert("m_analysisIsOk"); return m_perm; }
  153. protected:
  154. FactorType m_L; // The lower part stored in CSC
  155. VectorRx m_scale; // The vector for scaling the matrix
  156. RealScalar m_initialShift; // The initial shift parameter
  157. bool m_analysisIsOk;
  158. bool m_factorizationIsOk;
  159. ComputationInfo m_info;
  160. PermutationType m_perm;
  161. private:
  162. inline void updateList(Ref<const VectorIx> colPtr, Ref<VectorIx> rowIdx, Ref<VectorSx> vals, const Index& col, const Index& jk, VectorIx& firstElt, VectorList& listCol);
  163. };
  164. // Based on the following paper:
  165. // C-J. Lin and J. J. Moré, Incomplete Cholesky Factorizations with
  166. // Limited memory, SIAM J. Sci. Comput. 21(1), pp. 24-45, 1999
  167. // http://ftp.mcs.anl.gov/pub/tech_reports/reports/P682.pdf
  168. template<typename Scalar, int _UpLo, typename OrderingType>
  169. template<typename _MatrixType>
  170. void IncompleteCholesky<Scalar,_UpLo, OrderingType>::factorize(const _MatrixType& mat)
  171. {
  172. using std::sqrt;
  173. eigen_assert(m_analysisIsOk && "analyzePattern() should be called first");
  174. // Dropping strategy : Keep only the p largest elements per column, where p is the number of elements in the column of the original matrix. Other strategies will be added
  175. // Apply the fill-reducing permutation computed in analyzePattern()
  176. if (m_perm.rows() == mat.rows() ) // To detect the null permutation
  177. {
  178. // The temporary is needed to make sure that the diagonal entry is properly sorted
  179. FactorType tmp(mat.rows(), mat.cols());
  180. tmp = mat.template selfadjointView<_UpLo>().twistedBy(m_perm);
  181. m_L.template selfadjointView<Lower>() = tmp.template selfadjointView<Lower>();
  182. }
  183. else
  184. {
  185. m_L.template selfadjointView<Lower>() = mat.template selfadjointView<_UpLo>();
  186. }
  187. Index n = m_L.cols();
  188. Index nnz = m_L.nonZeros();
  189. Map<VectorSx> vals(m_L.valuePtr(), nnz); //values
  190. Map<VectorIx> rowIdx(m_L.innerIndexPtr(), nnz); //Row indices
  191. Map<VectorIx> colPtr( m_L.outerIndexPtr(), n+1); // Pointer to the beginning of each row
  192. VectorIx firstElt(n-1); // for each j, points to the next entry in vals that will be used in the factorization
  193. VectorList listCol(n); // listCol(j) is a linked list of columns to update column j
  194. VectorSx col_vals(n); // Store a nonzero values in each column
  195. VectorIx col_irow(n); // Row indices of nonzero elements in each column
  196. VectorIx col_pattern(n);
  197. col_pattern.fill(-1);
  198. StorageIndex col_nnz;
  199. // Computes the scaling factors
  200. m_scale.resize(n);
  201. m_scale.setZero();
  202. for (Index j = 0; j < n; j++)
  203. for (Index k = colPtr[j]; k < colPtr[j+1]; k++)
  204. {
  205. m_scale(j) += numext::abs2(vals(k));
  206. if(rowIdx[k]!=j)
  207. m_scale(rowIdx[k]) += numext::abs2(vals(k));
  208. }
  209. m_scale = m_scale.cwiseSqrt().cwiseSqrt();
  210. for (Index j = 0; j < n; ++j)
  211. if(m_scale(j)>(std::numeric_limits<RealScalar>::min)())
  212. m_scale(j) = RealScalar(1)/m_scale(j);
  213. else
  214. m_scale(j) = 1;
  215. // TODO disable scaling if not needed, i.e., if it is roughly uniform? (this will make solve() faster)
  216. // Scale and compute the shift for the matrix
  217. RealScalar mindiag = NumTraits<RealScalar>::highest();
  218. for (Index j = 0; j < n; j++)
  219. {
  220. for (Index k = colPtr[j]; k < colPtr[j+1]; k++)
  221. vals[k] *= (m_scale(j)*m_scale(rowIdx[k]));
  222. eigen_internal_assert(rowIdx[colPtr[j]]==j && "IncompleteCholesky: only the lower triangular part must be stored");
  223. mindiag = numext::mini(numext::real(vals[colPtr[j]]), mindiag);
  224. }
  225. FactorType L_save = m_L;
  226. RealScalar shift = 0;
  227. if(mindiag <= RealScalar(0.))
  228. shift = m_initialShift - mindiag;
  229. m_info = NumericalIssue;
  230. // Try to perform the incomplete factorization using the current shift
  231. int iter = 0;
  232. do
  233. {
  234. // Apply the shift to the diagonal elements of the matrix
  235. for (Index j = 0; j < n; j++)
  236. vals[colPtr[j]] += shift;
  237. // jki version of the Cholesky factorization
  238. Index j=0;
  239. for (; j < n; ++j)
  240. {
  241. // Left-looking factorization of the j-th column
  242. // First, load the j-th column into col_vals
  243. Scalar diag = vals[colPtr[j]]; // It is assumed that only the lower part is stored
  244. col_nnz = 0;
  245. for (Index i = colPtr[j] + 1; i < colPtr[j+1]; i++)
  246. {
  247. StorageIndex l = rowIdx[i];
  248. col_vals(col_nnz) = vals[i];
  249. col_irow(col_nnz) = l;
  250. col_pattern(l) = col_nnz;
  251. col_nnz++;
  252. }
  253. {
  254. typename std::list<StorageIndex>::iterator k;
  255. // Browse all previous columns that will update column j
  256. for(k = listCol[j].begin(); k != listCol[j].end(); k++)
  257. {
  258. Index jk = firstElt(*k); // First element to use in the column
  259. eigen_internal_assert(rowIdx[jk]==j);
  260. Scalar v_j_jk = numext::conj(vals[jk]);
  261. jk += 1;
  262. for (Index i = jk; i < colPtr[*k+1]; i++)
  263. {
  264. StorageIndex l = rowIdx[i];
  265. if(col_pattern[l]<0)
  266. {
  267. col_vals(col_nnz) = vals[i] * v_j_jk;
  268. col_irow[col_nnz] = l;
  269. col_pattern(l) = col_nnz;
  270. col_nnz++;
  271. }
  272. else
  273. col_vals(col_pattern[l]) -= vals[i] * v_j_jk;
  274. }
  275. updateList(colPtr,rowIdx,vals, *k, jk, firstElt, listCol);
  276. }
  277. }
  278. // Scale the current column
  279. if(numext::real(diag) <= 0)
  280. {
  281. if(++iter>=10)
  282. return;
  283. // increase shift
  284. shift = numext::maxi(m_initialShift,RealScalar(2)*shift);
  285. // restore m_L, col_pattern, and listCol
  286. vals = Map<const VectorSx>(L_save.valuePtr(), nnz);
  287. rowIdx = Map<const VectorIx>(L_save.innerIndexPtr(), nnz);
  288. colPtr = Map<const VectorIx>(L_save.outerIndexPtr(), n+1);
  289. col_pattern.fill(-1);
  290. for(Index i=0; i<n; ++i)
  291. listCol[i].clear();
  292. break;
  293. }
  294. RealScalar rdiag = sqrt(numext::real(diag));
  295. vals[colPtr[j]] = rdiag;
  296. for (Index k = 0; k<col_nnz; ++k)
  297. {
  298. Index i = col_irow[k];
  299. //Scale
  300. col_vals(k) /= rdiag;
  301. //Update the remaining diagonals with col_vals
  302. vals[colPtr[i]] -= numext::abs2(col_vals(k));
  303. }
  304. // Select the largest p elements
  305. // p is the original number of elements in the column (without the diagonal)
  306. Index p = colPtr[j+1] - colPtr[j] - 1 ;
  307. Ref<VectorSx> cvals = col_vals.head(col_nnz);
  308. Ref<VectorIx> cirow = col_irow.head(col_nnz);
  309. internal::QuickSplit(cvals,cirow, p);
  310. // Insert the largest p elements in the matrix
  311. Index cpt = 0;
  312. for (Index i = colPtr[j]+1; i < colPtr[j+1]; i++)
  313. {
  314. vals[i] = col_vals(cpt);
  315. rowIdx[i] = col_irow(cpt);
  316. // restore col_pattern:
  317. col_pattern(col_irow(cpt)) = -1;
  318. cpt++;
  319. }
  320. // Get the first smallest row index and put it after the diagonal element
  321. Index jk = colPtr(j)+1;
  322. updateList(colPtr,rowIdx,vals,j,jk,firstElt,listCol);
  323. }
  324. if(j==n)
  325. {
  326. m_factorizationIsOk = true;
  327. m_info = Success;
  328. }
  329. } while(m_info!=Success);
  330. }
  331. template<typename Scalar, int _UpLo, typename OrderingType>
  332. inline void IncompleteCholesky<Scalar,_UpLo, OrderingType>::updateList(Ref<const VectorIx> colPtr, Ref<VectorIx> rowIdx, Ref<VectorSx> vals, const Index& col, const Index& jk, VectorIx& firstElt, VectorList& listCol)
  333. {
  334. if (jk < colPtr(col+1) )
  335. {
  336. Index p = colPtr(col+1) - jk;
  337. Index minpos;
  338. rowIdx.segment(jk,p).minCoeff(&minpos);
  339. minpos += jk;
  340. if (rowIdx(minpos) != rowIdx(jk))
  341. {
  342. //Swap
  343. std::swap(rowIdx(jk),rowIdx(minpos));
  344. std::swap(vals(jk),vals(minpos));
  345. }
  346. firstElt(col) = internal::convert_index<StorageIndex,Index>(jk);
  347. listCol[rowIdx(jk)].push_back(internal::convert_index<StorageIndex,Index>(col));
  348. }
  349. }
  350. } // end namespace Eigen
  351. #endif