You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

248 lines
7.1 KiB

  1. namespace Eigen {
  2. /** \eigenManualPage SparseQuickRefPage Quick reference guide for sparse matrices
  3. \eigenAutoToc
  4. <hr>
  5. In this page, we give a quick summary of the main operations available for sparse matrices in the class SparseMatrix. First, it is recommended to read the introductory tutorial at \ref TutorialSparse. The important point to have in mind when working on sparse matrices is how they are stored :
  6. i.e either row major or column major. The default is column major. Most arithmetic operations on sparse matrices will assert that they have the same storage order.
  7. \section SparseMatrixInit Sparse Matrix Initialization
  8. <table class="manual">
  9. <tr><th> Category </th> <th> Operations</th> <th>Notes</th></tr>
  10. <tr><td>Constructor</td>
  11. <td>
  12. \code
  13. SparseMatrix<double> sm1(1000,1000);
  14. SparseMatrix<std::complex<double>,RowMajor> sm2;
  15. \endcode
  16. </td> <td> Default is ColMajor</td> </tr>
  17. <tr class="alt">
  18. <td> Resize/Reserve</td>
  19. <td>
  20. \code
  21. sm1.resize(m,n); //Change sm1 to a m x n matrix.
  22. sm1.reserve(nnz); // Allocate room for nnz nonzeros elements.
  23. \endcode
  24. </td>
  25. <td> Note that when calling reserve(), it is not required that nnz is the exact number of nonzero elements in the final matrix. However, an exact estimation will avoid multiple reallocations during the insertion phase. </td>
  26. </tr>
  27. <tr>
  28. <td> Assignment </td>
  29. <td>
  30. \code
  31. SparseMatrix<double,Colmajor> sm1;
  32. // Initialize sm2 with sm1.
  33. SparseMatrix<double,Rowmajor> sm2(sm1), sm3;
  34. // Assignment and evaluations modify the storage order.
  35. sm3 = sm1;
  36. \endcode
  37. </td>
  38. <td> The copy constructor can be used to convert from a storage order to another</td>
  39. </tr>
  40. <tr class="alt">
  41. <td> Element-wise Insertion</td>
  42. <td>
  43. \code
  44. // Insert a new element;
  45. sm1.insert(i, j) = v_ij;
  46. // Update the value v_ij
  47. sm1.coeffRef(i,j) = v_ij;
  48. sm1.coeffRef(i,j) += v_ij;
  49. sm1.coeffRef(i,j) -= v_ij;
  50. \endcode
  51. </td>
  52. <td> insert() assumes that the element does not already exist; otherwise, use coeffRef()</td>
  53. </tr>
  54. <tr>
  55. <td> Batch insertion</td>
  56. <td>
  57. \code
  58. std::vector< Eigen::Triplet<double> > tripletList;
  59. tripletList.reserve(estimation_of_entries);
  60. // -- Fill tripletList with nonzero elements...
  61. sm1.setFromTriplets(TripletList.begin(), TripletList.end());
  62. \endcode
  63. </td>
  64. <td>A complete example is available at \link TutorialSparseFilling Triplet Insertion \endlink.</td>
  65. </tr>
  66. <tr class="alt">
  67. <td> Constant or Random Insertion</td>
  68. <td>
  69. \code
  70. sm1.setZero();
  71. \endcode
  72. </td>
  73. <td>Remove all non-zero coefficients</td>
  74. </tr>
  75. </table>
  76. \section SparseBasicInfos Matrix properties
  77. Beyond the basic functions rows() and cols(), there are some useful functions that are available to easily get some informations from the matrix.
  78. <table class="manual">
  79. <tr>
  80. <td> \code
  81. sm1.rows(); // Number of rows
  82. sm1.cols(); // Number of columns
  83. sm1.nonZeros(); // Number of non zero values
  84. sm1.outerSize(); // Number of columns (resp. rows) for a column major (resp. row major )
  85. sm1.innerSize(); // Number of rows (resp. columns) for a row major (resp. column major)
  86. sm1.norm(); // Euclidian norm of the matrix
  87. sm1.squaredNorm(); // Squared norm of the matrix
  88. sm1.blueNorm();
  89. sm1.isVector(); // Check if sm1 is a sparse vector or a sparse matrix
  90. sm1.isCompressed(); // Check if sm1 is in compressed form
  91. ...
  92. \endcode </td>
  93. </tr>
  94. </table>
  95. \section SparseBasicOps Arithmetic operations
  96. It is easy to perform arithmetic operations on sparse matrices provided that the dimensions are adequate and that the matrices have the same storage order. Note that the evaluation can always be done in a matrix with a different storage order. In the following, \b sm denotes a sparse matrix, \b dm a dense matrix and \b dv a dense vector.
  97. <table class="manual">
  98. <tr><th> Operations </th> <th> Code </th> <th> Notes </th></tr>
  99. <tr>
  100. <td> add subtract </td>
  101. <td> \code
  102. sm3 = sm1 + sm2;
  103. sm3 = sm1 - sm2;
  104. sm2 += sm1;
  105. sm2 -= sm1; \endcode
  106. </td>
  107. <td>
  108. sm1 and sm2 should have the same storage order
  109. </td>
  110. </tr>
  111. <tr class="alt"><td>
  112. scalar product</td><td>\code
  113. sm3 = sm1 * s1; sm3 *= s1;
  114. sm3 = s1 * sm1 + s2 * sm2; sm3 /= s1;\endcode
  115. </td>
  116. <td>
  117. Many combinations are possible if the dimensions and the storage order agree.
  118. </tr>
  119. <tr>
  120. <td> %Sparse %Product </td>
  121. <td> \code
  122. sm3 = sm1 * sm2;
  123. dm2 = sm1 * dm1;
  124. dv2 = sm1 * dv1;
  125. \endcode </td>
  126. <td>
  127. </td>
  128. </tr>
  129. <tr class='alt'>
  130. <td> transposition, adjoint</td>
  131. <td> \code
  132. sm2 = sm1.transpose();
  133. sm2 = sm1.adjoint();
  134. \endcode </td>
  135. <td>
  136. Note that the transposition change the storage order. There is no support for transposeInPlace().
  137. </td>
  138. </tr>
  139. <tr>
  140. <td> Permutation </td>
  141. <td>
  142. \code
  143. perm.indices(); // Reference to the vector of indices
  144. sm1.twistedBy(perm); // Permute rows and columns
  145. sm2 = sm1 * perm; //Permute the columns
  146. sm2 = perm * sm1; // Permute the columns
  147. \endcode
  148. </td>
  149. <td>
  150. </td>
  151. </tr>
  152. <tr>
  153. <td>
  154. Component-wise ops
  155. </td>
  156. <td>\code
  157. sm1.cwiseProduct(sm2);
  158. sm1.cwiseQuotient(sm2);
  159. sm1.cwiseMin(sm2);
  160. sm1.cwiseMax(sm2);
  161. sm1.cwiseAbs();
  162. sm1.cwiseSqrt();
  163. \endcode</td>
  164. <td>
  165. sm1 and sm2 should have the same storage order
  166. </td>
  167. </tr>
  168. </table>
  169. \section sparseotherops Other supported operations
  170. <table class="manual">
  171. <tr><th>Operations</th> <th> Code </th> <th> Notes</th> </tr>
  172. <tr>
  173. <td>Sub-matrices</td>
  174. <td>
  175. \code
  176. sm1.block(startRow, startCol, rows, cols);
  177. sm1.block(startRow, startCol);
  178. sm1.topLeftCorner(rows, cols);
  179. sm1.topRightCorner(rows, cols);
  180. sm1.bottomLeftCorner( rows, cols);
  181. sm1.bottomRightCorner( rows, cols);
  182. \endcode
  183. </td> <td> </td>
  184. </tr>
  185. <tr>
  186. <td> Range </td>
  187. <td>
  188. \code
  189. sm1.innerVector(outer);
  190. sm1.innerVectors(start, size);
  191. sm1.leftCols(size);
  192. sm2.rightCols(size);
  193. sm1.middleRows(start, numRows);
  194. sm1.middleCols(start, numCols);
  195. sm1.col(j);
  196. \endcode
  197. </td>
  198. <td>A inner vector is either a row (for row-major) or a column (for column-major). As stated earlier, the evaluation can be done in a matrix with different storage order </td>
  199. </tr>
  200. <tr>
  201. <td> Triangular and selfadjoint views</td>
  202. <td>
  203. \code
  204. sm2 = sm1.triangularview<Lower>();
  205. sm2 = sm1.selfadjointview<Lower>();
  206. \endcode
  207. </td>
  208. <td> Several combination between triangular views and blocks views are possible
  209. \code
  210. \endcode </td>
  211. </tr>
  212. <tr>
  213. <td>Triangular solve </td>
  214. <td>
  215. \code
  216. dv2 = sm1.triangularView<Upper>().solve(dv1);
  217. dv2 = sm1.topLeftCorner(size, size).triangularView<Lower>().solve(dv1);
  218. \endcode
  219. </td>
  220. <td> For general sparse solve, Use any suitable module described at \ref TopicSparseSystems </td>
  221. </tr>
  222. <tr>
  223. <td> Low-level API</td>
  224. <td>
  225. \code
  226. sm1.valuePtr(); // Pointer to the values
  227. sm1.innerIndextr(); // Pointer to the indices.
  228. sm1.outerIndexPtr(); //Pointer to the beginning of each inner vector
  229. \endcode
  230. </td>
  231. <td> If the matrix is not in compressed form, makeCompressed() should be called before. Note that these functions are mostly provided for interoperability purposes with external libraries. A better access to the values of the matrix is done by using the InnerIterator class as described in \link TutorialSparse the Tutorial Sparse \endlink section</td>
  232. </tr>
  233. </table>
  234. */
  235. }