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.

455 lines
24 KiB

  1. namespace Eigen {
  2. /** \page TutorialSparse Tutorial page 9 - Sparse Matrix
  3. \ingroup Tutorial
  4. \li \b Previous: \ref TutorialGeometry
  5. \li \b Next: \ref TutorialMapClass
  6. \b Table \b of \b contents \n
  7. - \ref TutorialSparseIntro
  8. - \ref TutorialSparseExample "Example"
  9. - \ref TutorialSparseSparseMatrix
  10. - \ref TutorialSparseFilling
  11. - \ref TutorialSparseDirectSolvers
  12. - \ref TutorialSparseFeatureSet
  13. - \ref TutorialSparse_BasicOps
  14. - \ref TutorialSparse_Products
  15. - \ref TutorialSparse_TriangularSelfadjoint
  16. - \ref TutorialSparse_Submat
  17. <hr>
  18. Manipulating and solving sparse problems involves various modules which are summarized below:
  19. <table class="manual">
  20. <tr><th>Module</th><th>Header file</th><th>Contents</th></tr>
  21. <tr><td>\link Sparse_Module SparseCore \endlink</td><td>\code#include <Eigen/SparseCore>\endcode</td><td>SparseMatrix and SparseVector classes, matrix assembly, basic sparse linear algebra (including sparse triangular solvers)</td></tr>
  22. <tr><td>\link SparseCholesky_Module SparseCholesky \endlink</td><td>\code#include <Eigen/SparseCholesky>\endcode</td><td>Direct sparse LLT and LDLT Cholesky factorization to solve sparse self-adjoint positive definite problems</td></tr>
  23. <tr><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>\code#include <Eigen/IterativeLinearSolvers>\endcode</td><td>Iterative solvers to solve large general linear square problems (including self-adjoint positive definite problems)</td></tr>
  24. <tr><td></td><td>\code#include <Eigen/Sparse>\endcode</td><td>Includes all the above modules</td></tr>
  25. </table>
  26. \section TutorialSparseIntro Sparse matrix representation
  27. In many applications (e.g., finite element methods) it is common to deal with very large matrices where only a few coefficients are different from zero. In such cases, memory consumption can be reduced and performance increased by using a specialized representation storing only the nonzero coefficients. Such a matrix is called a sparse matrix.
  28. \b The \b %SparseMatrix \b class
  29. The class SparseMatrix is the main sparse matrix representation of Eigen's sparse module; it offers high performance and low memory usage.
  30. It implements a more versatile variant of the widely-used Compressed Column (or Row) Storage scheme.
  31. It consists of four compact arrays:
  32. - \c Values: stores the coefficient values of the non-zeros.
  33. - \c InnerIndices: stores the row (resp. column) indices of the non-zeros.
  34. - \c OuterStarts: stores for each column (resp. row) the index of the first non-zero in the previous two arrays.
  35. - \c InnerNNZs: stores the number of non-zeros of each column (resp. row).
  36. The word \c inner refers to an \em inner \em vector that is a column for a column-major matrix, or a row for a row-major matrix.
  37. The word \c outer refers to the other direction.
  38. This storage scheme is better explained on an example. The following matrix
  39. <table class="manual">
  40. <tr><td> 0</td><td>3</td><td> 0</td><td>0</td><td> 0</td></tr>
  41. <tr><td>22</td><td>0</td><td> 0</td><td>0</td><td>17</td></tr>
  42. <tr><td> 7</td><td>5</td><td> 0</td><td>1</td><td> 0</td></tr>
  43. <tr><td> 0</td><td>0</td><td> 0</td><td>0</td><td> 0</td></tr>
  44. <tr><td> 0</td><td>0</td><td>14</td><td>0</td><td> 8</td></tr>
  45. </table>
  46. and one of its possible sparse, \b column \b major representation:
  47. <table class="manual">
  48. <tr><td>Values:</td> <td>22</td><td>7</td><td>_</td><td>3</td><td>5</td><td>14</td><td>_</td><td>_</td><td>1</td><td>_</td><td>17</td><td>8</td></tr>
  49. <tr><td>InnerIndices:</td> <td> 1</td><td>2</td><td>_</td><td>0</td><td>2</td><td> 4</td><td>_</td><td>_</td><td>2</td><td>_</td><td> 1</td><td>4</td></tr>
  50. </table>
  51. <table class="manual">
  52. <tr><td>OuterStarts:</td><td>0</td><td>3</td><td>5</td><td>8</td><td>10</td><td>\em 12 </td></tr>
  53. <tr><td>InnerNNZs:</td> <td>2</td><td>2</td><td>1</td><td>1</td><td> 2</td><td></td></tr>
  54. </table>
  55. Currently the elements of a given inner vector are guaranteed to be always sorted by increasing inner indices.
  56. The \c "_" indicates available free space to quickly insert new elements.
  57. Assuming no reallocation is needed, the insertion of a random element is therefore in O(nnz_j) where nnz_j is the number of nonzeros of the respective inner vector.
  58. On the other hand, inserting elements with increasing inner indices in a given inner vector is much more efficient since this only requires to increase the respective \c InnerNNZs entry that is a O(1) operation.
  59. The case where no empty space is available is a special case, and is refered as the \em compressed mode.
  60. It corresponds to the widely used Compressed Column (or Row) Storage schemes (CCS or CRS).
  61. Any SparseMatrix can be turned to this form by calling the SparseMatrix::makeCompressed() function.
  62. In this case, one can remark that the \c InnerNNZs array is redundant with \c OuterStarts because we the equality: \c InnerNNZs[j] = \c OuterStarts[j+1]-\c OuterStarts[j].
  63. Therefore, in practice a call to SparseMatrix::makeCompressed() frees this buffer.
  64. It is worth noting that most of our wrappers to external libraries requires compressed matrices as inputs.
  65. The results of %Eigen's operations always produces \b compressed sparse matrices.
  66. On the other hand, the insertion of a new element into a SparseMatrix converts this later to the \b uncompressed mode.
  67. Here is the previous matrix represented in compressed mode:
  68. <table class="manual">
  69. <tr><td>Values:</td> <td>22</td><td>7</td><td>3</td><td>5</td><td>14</td><td>1</td><td>17</td><td>8</td></tr>
  70. <tr><td>InnerIndices:</td> <td> 1</td><td>2</td><td>0</td><td>2</td><td> 4</td><td>2</td><td> 1</td><td>4</td></tr>
  71. </table>
  72. <table class="manual">
  73. <tr><td>OuterStarts:</td><td>0</td><td>2</td><td>4</td><td>5</td><td>6</td><td>\em 8 </td></tr>
  74. </table>
  75. A SparseVector is a special case of a SparseMatrix where only the \c Values and \c InnerIndices arrays are stored.
  76. There is no notion of compressed/uncompressed mode for a SparseVector.
  77. \section TutorialSparseExample First example
  78. Before describing each individual class, let's start with the following typical example: solving the Lapace equation \f$ \nabla u = 0 \f$ on a regular 2D grid using a finite difference scheme and Dirichlet boundary conditions.
  79. Such problem can be mathematically expressed as a linear problem of the form \f$ Ax=b \f$ where \f$ x \f$ is the vector of \c m unknowns (in our case, the values of the pixels), \f$ b \f$ is the right hand side vector resulting from the boundary conditions, and \f$ A \f$ is an \f$ m \times m \f$ matrix containing only a few non-zero elements resulting from the discretization of the Laplacian operator.
  80. <table class="manual">
  81. <tr><td>
  82. \include Tutorial_sparse_example.cpp
  83. </td>
  84. <td>
  85. \image html Tutorial_sparse_example.jpeg
  86. </td></tr></table>
  87. In this example, we start by defining a column-major sparse matrix type of double \c SparseMatrix<double>, and a triplet list of the same scalar type \c Triplet<double>. A triplet is a simple object representing a non-zero entry as the triplet: \c row index, \c column index, \c value.
  88. In the main function, we declare a list \c coefficients of triplets (as a std vector) and the right hand side vector \f$ b \f$ which are filled by the \a buildProblem function.
  89. The raw and flat list of non-zero entries is then converted to a true SparseMatrix object \c A.
  90. Note that the elements of the list do not have to be sorted, and possible duplicate entries will be summed up.
  91. The last step consists of effectively solving the assembled problem.
  92. Since the resulting matrix \c A is symmetric by construction, we can perform a direct Cholesky factorization via the SimplicialLDLT class which behaves like its LDLT counterpart for dense objects.
  93. The resulting vector \c x contains the pixel values as a 1D array which is saved to a jpeg file shown on the right of the code above.
  94. Describing the \a buildProblem and \a save functions is out of the scope of this tutorial. They are given \ref TutorialSparse_example_details "here" for the curious and reproducibility purpose.
  95. \section TutorialSparseSparseMatrix The SparseMatrix class
  96. \b %Matrix \b and \b vector \b properties \n
  97. The SparseMatrix and SparseVector classes take three template arguments:
  98. * the scalar type (e.g., double)
  99. * the storage order (ColMajor or RowMajor, the default is RowMajor)
  100. * the inner index type (default is \c int).
  101. As for dense Matrix objects, constructors takes the size of the object.
  102. Here are some examples:
  103. \code
  104. SparseMatrix<std::complex<float> > mat(1000,2000); // declares a 1000x2000 column-major compressed sparse matrix of complex<float>
  105. SparseMatrix<double,RowMajor> mat(1000,2000); // declares a 1000x2000 row-major compressed sparse matrix of double
  106. SparseVector<std::complex<float> > vec(1000); // declares a column sparse vector of complex<float> of size 1000
  107. SparseVector<double,RowMajor> vec(1000); // declares a row sparse vector of double of size 1000
  108. \endcode
  109. In the rest of the tutorial, \c mat and \c vec represent any sparse-matrix and sparse-vector objects, respectively.
  110. The dimensions of a matrix can be queried using the following functions:
  111. <table class="manual">
  112. <tr><td>Standard \n dimensions</td><td>\code
  113. mat.rows()
  114. mat.cols()\endcode</td>
  115. <td>\code
  116. vec.size() \endcode</td>
  117. </tr>
  118. <tr><td>Sizes along the \n inner/outer dimensions</td><td>\code
  119. mat.innerSize()
  120. mat.outerSize()\endcode</td>
  121. <td></td>
  122. </tr>
  123. <tr><td>Number of non \n zero coefficients</td><td>\code
  124. mat.nonZeros() \endcode</td>
  125. <td>\code
  126. vec.nonZeros() \endcode</td></tr>
  127. </table>
  128. \b Iterating \b over \b the \b nonzero \b coefficients \n
  129. Random access to the elements of a sparse object can be done through the \c coeffRef(i,j) function.
  130. However, this function involves a quite expensive binary search.
  131. In most cases, one only wants to iterate over the non-zeros elements. This is achieved by a standard loop over the outer dimension, and then by iterating over the non-zeros of the current inner vector via an InnerIterator. Thus, the non-zero entries have to be visited in the same order than the storage order.
  132. Here is an example:
  133. <table class="manual">
  134. <tr><td>
  135. \code
  136. SparseMatrix<double> mat(rows,cols);
  137. for (int k=0; k<mat.outerSize(); ++k)
  138. for (SparseMatrix<double>::InnerIterator it(mat,k); it; ++it)
  139. {
  140. it.value();
  141. it.row(); // row index
  142. it.col(); // col index (here it is equal to k)
  143. it.index(); // inner index, here it is equal to it.row()
  144. }
  145. \endcode
  146. </td><td>
  147. \code
  148. SparseVector<double> vec(size);
  149. for (SparseVector<double>::InnerIterator it(vec); it; ++it)
  150. {
  151. it.value(); // == vec[ it.index() ]
  152. it.index();
  153. }
  154. \endcode
  155. </td></tr>
  156. </table>
  157. For a writable expression, the referenced value can be modified using the valueRef() function.
  158. If the type of the sparse matrix or vector depends on a template parameter, then the \c typename keyword is
  159. required to indicate that \c InnerIterator denotes a type; see \ref TopicTemplateKeyword for details.
  160. \section TutorialSparseFilling Filling a sparse matrix
  161. Because of the special storage scheme of a SparseMatrix, special care has to be taken when adding new nonzero entries.
  162. For instance, the cost of a single purely random insertion into a SparseMatrix is \c O(nnz), where \c nnz is the current number of non-zero coefficients.
  163. The simplest way to create a sparse matrix while guaranteeing good performance is thus to first build a list of so-called \em triplets, and then convert it to a SparseMatrix.
  164. Here is a typical usage example:
  165. \code
  166. typedef Eigen::Triplet<double> T;
  167. std::vector<T> tripletList;
  168. tripletList.reserve(estimation_of_entries);
  169. for(...)
  170. {
  171. // ...
  172. tripletList.push_back(T(i,j,v_ij));
  173. }
  174. SparseMatrixType mat(rows,cols);
  175. mat.setFromTriplets(tripletList.begin(), tripletList.end());
  176. // mat is ready to go!
  177. \endcode
  178. The \c std::vector of triplets might contain the elements in arbitrary order, and might even contain duplicated elements that will be summed up by setFromTriplets().
  179. See the SparseMatrix::setFromTriplets() function and class Triplet for more details.
  180. In some cases, however, slightly higher performance, and lower memory consumption can be reached by directly inserting the non-zeros into the destination matrix.
  181. A typical scenario of this approach is illustrated bellow:
  182. \code
  183. 1: SparseMatrix<double> mat(rows,cols); // default is column major
  184. 2: mat.reserve(VectorXi::Constant(cols,6));
  185. 3: for each i,j such that v_ij != 0
  186. 4: mat.insert(i,j) = v_ij; // alternative: mat.coeffRef(i,j) += v_ij;
  187. 5: mat.makeCompressed(); // optional
  188. \endcode
  189. - The key ingredient here is the line 2 where we reserve room for 6 non-zeros per column. In many cases, the number of non-zeros per column or row can easily be known in advance. If it varies significantly for each inner vector, then it is possible to specify a reserve size for each inner vector by providing a vector object with an operator[](int j) returning the reserve size of the \c j-th inner vector (e.g., via a VectorXi or std::vector<int>). If only a rought estimate of the number of nonzeros per inner-vector can be obtained, it is highly recommended to overestimate it rather than the opposite. If this line is omitted, then the first insertion of a new element will reserve room for 2 elements per inner vector.
  190. - The line 4 performs a sorted insertion. In this example, the ideal case is when the \c j-th column is not full and contains non-zeros whose inner-indices are smaller than \c i. In this case, this operation boils down to trivial O(1) operation.
  191. - When calling insert(i,j) the element \c i \c ,j must not already exists, otherwise use the coeffRef(i,j) method that will allow to, e.g., accumulate values. This method first performs a binary search and finally calls insert(i,j) if the element does not already exist. It is more flexible than insert() but also more costly.
  192. - The line 5 suppresses the remaining empty space and transforms the matrix into a compressed column storage.
  193. \section TutorialSparseDirectSolvers Solving linear problems
  194. %Eigen currently provides a limited set of built-in solvers, as well as wrappers to external solver libraries.
  195. They are summarized in the following table:
  196. <table class="manual">
  197. <tr><th>Class</th><th>Module</th><th>Solver kind</th><th>Matrix kind</th><th>Features related to performance</th>
  198. <th>Dependencies,License</th><th class="width20em"><p>Notes</p></th></tr>
  199. <tr><td>SimplicialLLT </td><td>\link SparseCholesky_Module SparseCholesky \endlink</td><td>Direct LLt factorization</td><td>SPD</td><td>Fill-in reducing</td>
  200. <td>built-in, LGPL</td>
  201. <td>SimplicialLDLT is often preferable</td></tr>
  202. <tr><td>SimplicialLDLT </td><td>\link SparseCholesky_Module SparseCholesky \endlink</td><td>Direct LDLt factorization</td><td>SPD</td><td>Fill-in reducing</td>
  203. <td>built-in, LGPL</td>
  204. <td>Recommended for very sparse and not too large problems (e.g., 2D Poisson eq.)</td></tr>
  205. <tr><td>ConjugateGradient</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>Classic iterative CG</td><td>SPD</td><td>Preconditionning</td>
  206. <td>built-in, LGPL</td>
  207. <td>Recommended for large symmetric problems (e.g., 3D Poisson eq.)</td></tr>
  208. <tr><td>BiCGSTAB</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>Iterative stabilized bi-conjugate gradient</td><td>Square</td><td>Preconditionning</td>
  209. <td>built-in, LGPL</td>
  210. <td>Might not always converge</td></tr>
  211. <tr><td>PastixLLT \n PastixLDLT \n PastixLU</td><td>\link PaStiXSupport_Module PaStiXSupport \endlink</td><td>Direct LLt, LDLt, LU factorizations</td><td>SPD \n SPD \n Square</td><td>Fill-in reducing, Leverage fast dense algebra, Multithreading</td>
  212. <td>Requires the <a href="http://pastix.gforge.inria.fr">PaStiX</a> package, \b CeCILL-C </td>
  213. <td>optimized for tough problems and symmetric patterns</td></tr>
  214. <tr><td>CholmodSupernodalLLT</td><td>\link CholmodSupport_Module CholmodSupport \endlink</td><td>Direct LLt factorization</td><td>SPD</td><td>Fill-in reducing, Leverage fast dense algebra</td>
  215. <td>Requires the <a href="http://www.cise.ufl.edu/research/sparse/SuiteSparse/">SuiteSparse</a> package, \b GPL </td>
  216. <td></td></tr>
  217. <tr><td>UmfPackLU</td><td>\link UmfPackSupport_Module UmfPackSupport \endlink</td><td>Direct LU factorization</td><td>Square</td><td>Fill-in reducing, Leverage fast dense algebra</td>
  218. <td>Requires the <a href="http://www.cise.ufl.edu/research/sparse/SuiteSparse/">SuiteSparse</a> package, \b GPL </td>
  219. <td></td></tr>
  220. <tr><td>SuperLU</td><td>\link SuperLUSupport_Module SuperLUSupport \endlink</td><td>Direct LU factorization</td><td>Square</td><td>Fill-in reducing, Leverage fast dense algebra</td>
  221. <td>Requires the <a href="http://crd-legacy.lbl.gov/~xiaoye/SuperLU/">SuperLU</a> library, (BSD-like)</td>
  222. <td></td></tr>
  223. </table>
  224. Here \c SPD means symmetric positive definite.
  225. All these solvers follow the same general concept.
  226. Here is a typical and general example:
  227. \code
  228. #include <Eigen/RequiredModuleName>
  229. // ...
  230. SparseMatrix<double> A;
  231. // fill A
  232. VectorXd b, x;
  233. // fill b
  234. // solve Ax = b
  235. SolverClassName<SparseMatrix<double> > solver;
  236. solver.compute(A);
  237. if(solver.info()!=Success) {
  238. // decomposition failed
  239. return;
  240. }
  241. x = solver.solve(b);
  242. if(solver.info()!=Success) {
  243. // solving failed
  244. return;
  245. }
  246. // solve for another right hand side:
  247. x1 = solver.solve(b1);
  248. \endcode
  249. For \c SPD solvers, a second optional template argument allows to specify which triangular part have to be used, e.g.:
  250. \code
  251. #include <Eigen/IterativeLinearSolvers>
  252. ConjugateGradient<SparseMatrix<double>, Eigen::Upper> solver;
  253. x = solver.compute(A).solve(b);
  254. \endcode
  255. In the above example, only the upper triangular part of the input matrix A is considered for solving. The opposite triangle might either be empty or contain arbitrary values.
  256. In the case where multiple problems with the same sparcity pattern have to be solved, then the "compute" step can be decomposed as follow:
  257. \code
  258. SolverClassName<SparseMatrix<double> > solver;
  259. solver.analyzePattern(A); // for this step the numerical values of A are not used
  260. solver.factorize(A);
  261. x1 = solver.solve(b1);
  262. x2 = solver.solve(b2);
  263. ...
  264. A = ...; // modify the values of the nonzeros of A, the nonzeros pattern must stay unchanged
  265. solver.factorize(A);
  266. x1 = solver.solve(b1);
  267. x2 = solver.solve(b2);
  268. ...
  269. \endcode
  270. The compute() method is equivalent to calling both analyzePattern() and factorize().
  271. Finally, each solver provides some specific features, such as determinant, access to the factors, controls of the iterations, and so on.
  272. More details are availble in the documentations of the respective classes.
  273. \section TutorialSparseFeatureSet Supported operators and functions
  274. Because of their special storage format, sparse matrices cannot offer the same level of flexbility than dense matrices.
  275. In Eigen's sparse module we chose to expose only the subset of the dense matrix API which can be efficiently implemented.
  276. In the following \em sm denotes a sparse matrix, \em sv a sparse vector, \em dm a dense matrix, and \em dv a dense vector.
  277. \subsection TutorialSparse_BasicOps Basic operations
  278. %Sparse expressions support most of the unary and binary coefficient wise operations:
  279. \code
  280. sm1.real() sm1.imag() -sm1 0.5*sm1
  281. sm1+sm2 sm1-sm2 sm1.cwiseProduct(sm2)
  282. \endcode
  283. However, a strong restriction is that the storage orders must match. For instance, in the following example:
  284. \code
  285. sm4 = sm1 + sm2 + sm3;
  286. \endcode
  287. sm1, sm2, and sm3 must all be row-major or all column major.
  288. On the other hand, there is no restriction on the target matrix sm4.
  289. For instance, this means that for computing \f$ A^T + A \f$, the matrix \f$ A^T \f$ must be evaluated into a temporary matrix of compatible storage order:
  290. \code
  291. SparseMatrix<double> A, B;
  292. B = SparseMatrix<double>(A.transpose()) + A;
  293. \endcode
  294. Binary coefficient wise operators can also mix sparse and dense expressions:
  295. \code
  296. sm2 = sm1.cwiseProduct(dm1);
  297. dm2 = sm1 + dm1;
  298. \endcode
  299. %Sparse expressions also support transposition:
  300. \code
  301. sm1 = sm2.transpose();
  302. sm1 = sm2.adjoint();
  303. \endcode
  304. However, there is no transposeInPlace() method.
  305. \subsection TutorialSparse_Products Matrix products
  306. %Eigen supports various kind of sparse matrix products which are summarize below:
  307. - \b sparse-dense:
  308. \code
  309. dv2 = sm1 * dv1;
  310. dm2 = dm1 * sm1.adjoint();
  311. dm2 = 2. * sm1 * dm1;
  312. \endcode
  313. - \b symmetric \b sparse-dense. The product of a sparse symmetric matrix with a dense matrix (or vector) can also be optimized by specifying the symmetry with selfadjointView():
  314. \code
  315. dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored
  316. dm2 = A.selfadjointView<Upper>() * dm1; // if only the upper part of A is stored
  317. dm2 = A.selfadjointView<Lower>() * dm1; // if only the lower part of A is stored
  318. \endcode
  319. - \b sparse-sparse. For sparse-sparse products, two different algorithms are available. The default one is conservative and preserve the explicit zeros that might appear:
  320. \code
  321. sm3 = sm1 * sm2;
  322. sm3 = 4 * sm1.adjoint() * sm2;
  323. \endcode
  324. The second algorithm prunes on the fly the explicit zeros, or the values smaller than a given threshold. It is enabled and controlled through the prune() functions:
  325. \code
  326. sm3 = (sm1 * sm2).prune(); // removes numerical zeros
  327. sm3 = (sm1 * sm2).prune(ref); // removes elements much smaller than ref
  328. sm3 = (sm1 * sm2).prune(ref,epsilon); // removes elements smaller than ref*epsilon
  329. \endcode
  330. - \b permutations. Finally, permutations can be applied to sparse matrices too:
  331. \code
  332. PermutationMatrix<Dynamic,Dynamic> P = ...;
  333. sm2 = P * sm1;
  334. sm2 = sm1 * P.inverse();
  335. sm2 = sm1.transpose() * P;
  336. \endcode
  337. \subsection TutorialSparse_TriangularSelfadjoint Triangular and selfadjoint views
  338. Just as with dense matrices, the triangularView() function can be used to address a triangular part of the matrix, and perform triangular solves with a dense right hand side:
  339. \code
  340. dm2 = sm1.triangularView<Lower>(dm1);
  341. dv2 = sm1.transpose().triangularView<Upper>(dv1);
  342. \endcode
  343. The selfadjointView() function permits various operations:
  344. - optimized sparse-dense matrix products:
  345. \code
  346. dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored
  347. dm2 = A.selfadjointView<Upper>() * dm1; // if only the upper part of A is stored
  348. dm2 = A.selfadjointView<Lower>() * dm1; // if only the lower part of A is stored
  349. \endcode
  350. - copy of triangular parts:
  351. \code
  352. sm2 = sm1.selfadjointView<Upper>(); // makes a full selfadjoint matrix from the upper triangular part
  353. sm2.selfadjointView<Lower>() = sm1.selfadjointView<Upper>(); // copies the upper triangular part to the lower triangular part
  354. \endcode
  355. - application of symmetric permutations:
  356. \code
  357. PermutationMatrix<Dynamic,Dynamic> P = ...;
  358. sm2 = A.selfadjointView<Upper>().twistedBy(P); // compute P S P' from the upper triangular part of A, and make it a full matrix
  359. sm2.selfadjointView<Lower>() = A.selfadjointView<Lower>().twistedBy(P); // compute P S P' from the lower triangular part of A, and then only compute the lower part
  360. \endcode
  361. \subsection TutorialSparse_Submat Sub-matrices
  362. %Sparse matrices does not support yet the addressing of arbitrary sub matrices. Currently, one can only reference a set of contiguous \em inner vectors, i.e., a set of contiguous rows for a row-major matrix, or a set of contiguous columns for a column major matrix:
  363. \code
  364. sm1.innerVector(j); // returns an expression of the j-th column (resp. row) of the matrix if sm1 is col-major (resp. row-major)
  365. sm1.innerVectors(j, nb); // returns an expression of the nb columns (resp. row) starting from the j-th column (resp. row)
  366. // of the matrix if sm1 is col-major (resp. row-major)
  367. sm1.middleRows(j, nb); // for row major matrices only, get a range of nb rows
  368. sm1.middleCols(j, nb); // for column major matrices only, get a range of nb columns
  369. \endcode
  370. \li \b Next: \ref TutorialMapClass
  371. */
  372. }