SuiteSparseQRSupport.h 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2012 Desire Nuentsa <desire.nuentsa_wakam@inria.fr>
  5. // Copyright (C) 2014 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_SUITESPARSEQRSUPPORT_H
  11. #define EIGEN_SUITESPARSEQRSUPPORT_H
  12. namespace Eigen {
  13. template<typename MatrixType> class SPQR;
  14. template<typename SPQRType> struct SPQRMatrixQReturnType;
  15. template<typename SPQRType> struct SPQRMatrixQTransposeReturnType;
  16. template <typename SPQRType, typename Derived> struct SPQR_QProduct;
  17. namespace internal {
  18. template <typename SPQRType> struct traits<SPQRMatrixQReturnType<SPQRType> >
  19. {
  20. typedef typename SPQRType::MatrixType ReturnType;
  21. };
  22. template <typename SPQRType> struct traits<SPQRMatrixQTransposeReturnType<SPQRType> >
  23. {
  24. typedef typename SPQRType::MatrixType ReturnType;
  25. };
  26. template <typename SPQRType, typename Derived> struct traits<SPQR_QProduct<SPQRType, Derived> >
  27. {
  28. typedef typename Derived::PlainObject ReturnType;
  29. };
  30. } // End namespace internal
  31. /**
  32. * \ingroup SPQRSupport_Module
  33. * \class SPQR
  34. * \brief Sparse QR factorization based on SuiteSparseQR library
  35. *
  36. * This class is used to perform a multithreaded and multifrontal rank-revealing QR decomposition
  37. * of sparse matrices. The result is then used to solve linear leasts_square systems.
  38. * Clearly, a QR factorization is returned such that A*P = Q*R where :
  39. *
  40. * P is the column permutation. Use colsPermutation() to get it.
  41. *
  42. * Q is the orthogonal matrix represented as Householder reflectors.
  43. * Use matrixQ() to get an expression and matrixQ().transpose() to get the transpose.
  44. * You can then apply it to a vector.
  45. *
  46. * R is the sparse triangular factor. Use matrixQR() to get it as SparseMatrix.
  47. * NOTE : The Index type of R is always SuiteSparse_long. You can get it with SPQR::Index
  48. *
  49. * \tparam _MatrixType The type of the sparse matrix A, must be a column-major SparseMatrix<>
  50. *
  51. * \implsparsesolverconcept
  52. *
  53. *
  54. */
  55. template<typename _MatrixType>
  56. class SPQR : public SparseSolverBase<SPQR<_MatrixType> >
  57. {
  58. protected:
  59. typedef SparseSolverBase<SPQR<_MatrixType> > Base;
  60. using Base::m_isInitialized;
  61. public:
  62. typedef typename _MatrixType::Scalar Scalar;
  63. typedef typename _MatrixType::RealScalar RealScalar;
  64. typedef SuiteSparse_long StorageIndex ;
  65. typedef SparseMatrix<Scalar, ColMajor, StorageIndex> MatrixType;
  66. typedef Map<PermutationMatrix<Dynamic, Dynamic, StorageIndex> > PermutationType;
  67. enum {
  68. ColsAtCompileTime = Dynamic,
  69. MaxColsAtCompileTime = Dynamic
  70. };
  71. public:
  72. SPQR()
  73. : m_analysisIsOk(false),
  74. m_factorizationIsOk(false),
  75. m_isRUpToDate(false),
  76. m_ordering(SPQR_ORDERING_DEFAULT),
  77. m_allow_tol(SPQR_DEFAULT_TOL),
  78. m_tolerance (NumTraits<Scalar>::epsilon()),
  79. m_cR(0),
  80. m_E(0),
  81. m_H(0),
  82. m_HPinv(0),
  83. m_HTau(0),
  84. m_useDefaultThreshold(true)
  85. {
  86. cholmod_l_start(&m_cc);
  87. }
  88. explicit SPQR(const _MatrixType& matrix)
  89. : m_analysisIsOk(false),
  90. m_factorizationIsOk(false),
  91. m_isRUpToDate(false),
  92. m_ordering(SPQR_ORDERING_DEFAULT),
  93. m_allow_tol(SPQR_DEFAULT_TOL),
  94. m_tolerance (NumTraits<Scalar>::epsilon()),
  95. m_cR(0),
  96. m_E(0),
  97. m_H(0),
  98. m_HPinv(0),
  99. m_HTau(0),
  100. m_useDefaultThreshold(true)
  101. {
  102. cholmod_l_start(&m_cc);
  103. compute(matrix);
  104. }
  105. ~SPQR()
  106. {
  107. SPQR_free();
  108. cholmod_l_finish(&m_cc);
  109. }
  110. void SPQR_free()
  111. {
  112. cholmod_l_free_sparse(&m_H, &m_cc);
  113. cholmod_l_free_sparse(&m_cR, &m_cc);
  114. cholmod_l_free_dense(&m_HTau, &m_cc);
  115. std::free(m_E);
  116. std::free(m_HPinv);
  117. }
  118. void compute(const _MatrixType& matrix)
  119. {
  120. if(m_isInitialized) SPQR_free();
  121. MatrixType mat(matrix);
  122. /* Compute the default threshold as in MatLab, see:
  123. * Tim Davis, "Algorithm 915, SuiteSparseQR: Multifrontal Multithreaded Rank-Revealing
  124. * Sparse QR Factorization, ACM Trans. on Math. Soft. 38(1), 2011, Page 8:3
  125. */
  126. RealScalar pivotThreshold = m_tolerance;
  127. if(m_useDefaultThreshold)
  128. {
  129. RealScalar max2Norm = 0.0;
  130. for (int j = 0; j < mat.cols(); j++) max2Norm = numext::maxi(max2Norm, mat.col(j).norm());
  131. if(max2Norm==RealScalar(0))
  132. max2Norm = RealScalar(1);
  133. pivotThreshold = 20 * (mat.rows() + mat.cols()) * max2Norm * NumTraits<RealScalar>::epsilon();
  134. }
  135. cholmod_sparse A;
  136. A = viewAsCholmod(mat);
  137. m_rows = matrix.rows();
  138. Index col = matrix.cols();
  139. m_rank = SuiteSparseQR<Scalar>(m_ordering, pivotThreshold, col, &A,
  140. &m_cR, &m_E, &m_H, &m_HPinv, &m_HTau, &m_cc);
  141. if (!m_cR)
  142. {
  143. m_info = NumericalIssue;
  144. m_isInitialized = false;
  145. return;
  146. }
  147. m_info = Success;
  148. m_isInitialized = true;
  149. m_isRUpToDate = false;
  150. }
  151. /**
  152. * Get the number of rows of the input matrix and the Q matrix
  153. */
  154. inline Index rows() const {return m_rows; }
  155. /**
  156. * Get the number of columns of the input matrix.
  157. */
  158. inline Index cols() const { return m_cR->ncol; }
  159. template<typename Rhs, typename Dest>
  160. void _solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest> &dest) const
  161. {
  162. eigen_assert(m_isInitialized && " The QR factorization should be computed first, call compute()");
  163. eigen_assert(b.cols()==1 && "This method is for vectors only");
  164. //Compute Q^T * b
  165. typename Dest::PlainObject y, y2;
  166. y = matrixQ().transpose() * b;
  167. // Solves with the triangular matrix R
  168. Index rk = this->rank();
  169. y2 = y;
  170. y.resize((std::max)(cols(),Index(y.rows())),y.cols());
  171. y.topRows(rk) = this->matrixR().topLeftCorner(rk, rk).template triangularView<Upper>().solve(y2.topRows(rk));
  172. // Apply the column permutation
  173. // colsPermutation() performs a copy of the permutation,
  174. // so let's apply it manually:
  175. for(Index i = 0; i < rk; ++i) dest.row(m_E[i]) = y.row(i);
  176. for(Index i = rk; i < cols(); ++i) dest.row(m_E[i]).setZero();
  177. // y.bottomRows(y.rows()-rk).setZero();
  178. // dest = colsPermutation() * y.topRows(cols());
  179. m_info = Success;
  180. }
  181. /** \returns the sparse triangular factor R. It is a sparse matrix
  182. */
  183. const MatrixType matrixR() const
  184. {
  185. eigen_assert(m_isInitialized && " The QR factorization should be computed first, call compute()");
  186. if(!m_isRUpToDate) {
  187. m_R = viewAsEigen<Scalar,ColMajor, typename MatrixType::StorageIndex>(*m_cR);
  188. m_isRUpToDate = true;
  189. }
  190. return m_R;
  191. }
  192. /// Get an expression of the matrix Q
  193. SPQRMatrixQReturnType<SPQR> matrixQ() const
  194. {
  195. return SPQRMatrixQReturnType<SPQR>(*this);
  196. }
  197. /// Get the permutation that was applied to columns of A
  198. PermutationType colsPermutation() const
  199. {
  200. eigen_assert(m_isInitialized && "Decomposition is not initialized.");
  201. return PermutationType(m_E, m_cR->ncol);
  202. }
  203. /**
  204. * Gets the rank of the matrix.
  205. * It should be equal to matrixQR().cols if the matrix is full-rank
  206. */
  207. Index rank() const
  208. {
  209. eigen_assert(m_isInitialized && "Decomposition is not initialized.");
  210. return m_cc.SPQR_istat[4];
  211. }
  212. /// Set the fill-reducing ordering method to be used
  213. void setSPQROrdering(int ord) { m_ordering = ord;}
  214. /// Set the tolerance tol to treat columns with 2-norm < =tol as zero
  215. void setPivotThreshold(const RealScalar& tol)
  216. {
  217. m_useDefaultThreshold = false;
  218. m_tolerance = tol;
  219. }
  220. /** \returns a pointer to the SPQR workspace */
  221. cholmod_common *cholmodCommon() const { return &m_cc; }
  222. /** \brief Reports whether previous computation was successful.
  223. *
  224. * \returns \c Success if computation was successful,
  225. * \c NumericalIssue if the sparse QR can not be computed
  226. */
  227. ComputationInfo info() const
  228. {
  229. eigen_assert(m_isInitialized && "Decomposition is not initialized.");
  230. return m_info;
  231. }
  232. protected:
  233. bool m_analysisIsOk;
  234. bool m_factorizationIsOk;
  235. mutable bool m_isRUpToDate;
  236. mutable ComputationInfo m_info;
  237. int m_ordering; // Ordering method to use, see SPQR's manual
  238. int m_allow_tol; // Allow to use some tolerance during numerical factorization.
  239. RealScalar m_tolerance; // treat columns with 2-norm below this tolerance as zero
  240. mutable cholmod_sparse *m_cR; // The sparse R factor in cholmod format
  241. mutable MatrixType m_R; // The sparse matrix R in Eigen format
  242. mutable StorageIndex *m_E; // The permutation applied to columns
  243. mutable cholmod_sparse *m_H; //The householder vectors
  244. mutable StorageIndex *m_HPinv; // The row permutation of H
  245. mutable cholmod_dense *m_HTau; // The Householder coefficients
  246. mutable Index m_rank; // The rank of the matrix
  247. mutable cholmod_common m_cc; // Workspace and parameters
  248. bool m_useDefaultThreshold; // Use default threshold
  249. Index m_rows;
  250. template<typename ,typename > friend struct SPQR_QProduct;
  251. };
  252. template <typename SPQRType, typename Derived>
  253. struct SPQR_QProduct : ReturnByValue<SPQR_QProduct<SPQRType,Derived> >
  254. {
  255. typedef typename SPQRType::Scalar Scalar;
  256. typedef typename SPQRType::StorageIndex StorageIndex;
  257. //Define the constructor to get reference to argument types
  258. SPQR_QProduct(const SPQRType& spqr, const Derived& other, bool transpose) : m_spqr(spqr),m_other(other),m_transpose(transpose) {}
  259. inline Index rows() const { return m_transpose ? m_spqr.rows() : m_spqr.cols(); }
  260. inline Index cols() const { return m_other.cols(); }
  261. // Assign to a vector
  262. template<typename ResType>
  263. void evalTo(ResType& res) const
  264. {
  265. cholmod_dense y_cd;
  266. cholmod_dense *x_cd;
  267. int method = m_transpose ? SPQR_QTX : SPQR_QX;
  268. cholmod_common *cc = m_spqr.cholmodCommon();
  269. y_cd = viewAsCholmod(m_other.const_cast_derived());
  270. x_cd = SuiteSparseQR_qmult<Scalar>(method, m_spqr.m_H, m_spqr.m_HTau, m_spqr.m_HPinv, &y_cd, cc);
  271. res = Matrix<Scalar,ResType::RowsAtCompileTime,ResType::ColsAtCompileTime>::Map(reinterpret_cast<Scalar*>(x_cd->x), x_cd->nrow, x_cd->ncol);
  272. cholmod_l_free_dense(&x_cd, cc);
  273. }
  274. const SPQRType& m_spqr;
  275. const Derived& m_other;
  276. bool m_transpose;
  277. };
  278. template<typename SPQRType>
  279. struct SPQRMatrixQReturnType{
  280. SPQRMatrixQReturnType(const SPQRType& spqr) : m_spqr(spqr) {}
  281. template<typename Derived>
  282. SPQR_QProduct<SPQRType, Derived> operator*(const MatrixBase<Derived>& other)
  283. {
  284. return SPQR_QProduct<SPQRType,Derived>(m_spqr,other.derived(),false);
  285. }
  286. SPQRMatrixQTransposeReturnType<SPQRType> adjoint() const
  287. {
  288. return SPQRMatrixQTransposeReturnType<SPQRType>(m_spqr);
  289. }
  290. // To use for operations with the transpose of Q
  291. SPQRMatrixQTransposeReturnType<SPQRType> transpose() const
  292. {
  293. return SPQRMatrixQTransposeReturnType<SPQRType>(m_spqr);
  294. }
  295. const SPQRType& m_spqr;
  296. };
  297. template<typename SPQRType>
  298. struct SPQRMatrixQTransposeReturnType{
  299. SPQRMatrixQTransposeReturnType(const SPQRType& spqr) : m_spqr(spqr) {}
  300. template<typename Derived>
  301. SPQR_QProduct<SPQRType,Derived> operator*(const MatrixBase<Derived>& other)
  302. {
  303. return SPQR_QProduct<SPQRType,Derived>(m_spqr,other.derived(), true);
  304. }
  305. const SPQRType& m_spqr;
  306. };
  307. }// End namespace Eigen
  308. #endif