source: pacpussensors/trunk/Vislab/lib3dv-1.2.0/lib3dv/eigen/doc/TutorialSparse.dox

Last change on this file was 136, checked in by ldecherf, 8 years ago

Doc

File size: 19.8 KB
Line 
1namespace Eigen {
2
3/** \eigenManualPage TutorialSparse Sparse matrix manipulations
4
5\eigenAutoToc
6
7Manipulating and solving sparse problems involves various modules which are summarized below:
8
9<table class="manual">
10<tr><th>Module</th><th>Header file</th><th>Contents</th></tr>
11<tr><td>\link SparseCore_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>
12<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>
13<tr><td>\link SparseLU_Module SparseLU \endlink</td><td>\code #include<Eigen/SparseLU> \endcode</td>
14<td>%Sparse LU factorization to solve general square sparse systems</td></tr>
15<tr><td>\link SparseQR_Module SparseQR \endlink</td><td>\code #include<Eigen/SparseQR>\endcode </td><td>%Sparse QR factorization for solving sparse linear least-squares problems</td></tr>
16<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>
17<tr><td>\link Sparse_Module Sparse \endlink</td><td>\code#include <Eigen/Sparse>\endcode</td><td>Includes all the above modules</td></tr>
18</table>
19
20\section TutorialSparseIntro Sparse matrix format
21
22In 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.
23
24\b The \b %SparseMatrix \b class
25
26The class SparseMatrix is the main sparse matrix representation of Eigen's sparse module; it offers high performance and low memory usage.
27It implements a more versatile variant of the widely-used Compressed Column (or Row) Storage scheme.
28It consists of four compact arrays:
29 - \c Values: stores the coefficient values of the non-zeros.
30 - \c InnerIndices: stores the row (resp. column) indices of the non-zeros.
31 - \c OuterStarts: stores for each column (resp. row) the index of the first non-zero in the previous two arrays.
32 - \c InnerNNZs: stores the number of non-zeros of each column (resp. row).
33The 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.
34The word \c outer refers to the other direction.
35
36This storage scheme is better explained on an example. The following matrix
37<table class="manual">
38<tr><td> 0</td><td>3</td><td> 0</td><td>0</td><td> 0</td></tr>
39<tr><td>22</td><td>0</td><td> 0</td><td>0</td><td>17</td></tr>
40<tr><td> 7</td><td>5</td><td> 0</td><td>1</td><td> 0</td></tr>
41<tr><td> 0</td><td>0</td><td> 0</td><td>0</td><td> 0</td></tr>
42<tr><td> 0</td><td>0</td><td>14</td><td>0</td><td> 8</td></tr>
43</table>
44
45and one of its possible sparse, \b column \b major representation:
46<table class="manual">
47<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>
48<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>
49</table>
50<table class="manual">
51<tr><td>OuterStarts:</td><td>0</td><td>3</td><td>5</td><td>8</td><td>10</td><td>\em 12 </td></tr>
52<tr><td>InnerNNZs:</td> <td>2</td><td>2</td><td>1</td><td>1</td><td> 2</td><td></td></tr>
53</table>
54
55Currently the elements of a given inner vector are guaranteed to be always sorted by increasing inner indices.
56The \c "_" indicates available free space to quickly insert new elements.
57Assuming 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.
58On 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
60The case where no empty space is available is a special case, and is refered as the \em compressed mode.
61It corresponds to the widely used Compressed Column (or Row) Storage schemes (CCS or CRS).
62Any SparseMatrix can be turned to this form by calling the SparseMatrix::makeCompressed() function.
63In 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].
64Therefore, in practice a call to SparseMatrix::makeCompressed() frees this buffer.
65
66It is worth noting that most of our wrappers to external libraries requires compressed matrices as inputs.
67
68The results of %Eigen's operations always produces \b compressed sparse matrices.
69On the other hand, the insertion of a new element into a SparseMatrix converts this later to the \b uncompressed mode.
70
71Here is the previous matrix represented in compressed mode:
72<table class="manual">
73<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>
74<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>
75</table>
76<table class="manual">
77<tr><td>OuterStarts:</td><td>0</td><td>2</td><td>4</td><td>5</td><td>6</td><td>\em 8 </td></tr>
78</table>
79
80A SparseVector is a special case of a SparseMatrix where only the \c Values and \c InnerIndices arrays are stored.
81There is no notion of compressed/uncompressed mode for a SparseVector.
82
83
84\section TutorialSparseExample First example
85
86Before describing each individual class, let's start with the following typical example: solving the Laplace equation \f$ \nabla u = 0 \f$ on a regular 2D grid using a finite difference scheme and Dirichlet boundary conditions.
87Such 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.
88
89<table class="manual">
90<tr><td>
91\include Tutorial_sparse_example.cpp
92</td>
93<td>
94\image html Tutorial_sparse_example.jpeg
95</td></tr></table>
96
97In 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.
98
99In 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.
100The raw and flat list of non-zero entries is then converted to a true SparseMatrix object \c A.
101Note that the elements of the list do not have to be sorted, and possible duplicate entries will be summed up.
102
103The last step consists of effectively solving the assembled problem.
104Since 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.
105
106The 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.
107
108Describing 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.
109
110
111
112
113\section TutorialSparseSparseMatrix The SparseMatrix class
114
115\b %Matrix \b and \b vector \b properties \n
116
117The SparseMatrix and SparseVector classes take three template arguments:
118 * the scalar type (e.g., double)
119 * the storage order (ColMajor or RowMajor, the default is ColMajor)
120 * the inner index type (default is \c int).
121
122As for dense Matrix objects, constructors takes the size of the object.
123Here are some examples:
124
125\code
126SparseMatrix<std::complex<float> > mat(1000,2000); // declares a 1000x2000 column-major compressed sparse matrix of complex<float>
127SparseMatrix<double,RowMajor> mat(1000,2000); // declares a 1000x2000 row-major compressed sparse matrix of double
128SparseVector<std::complex<float> > vec(1000); // declares a column sparse vector of complex<float> of size 1000
129SparseVector<double,RowMajor> vec(1000); // declares a row sparse vector of double of size 1000
130\endcode
131
132In the rest of the tutorial, \c mat and \c vec represent any sparse-matrix and sparse-vector objects, respectively.
133
134The dimensions of a matrix can be queried using the following functions:
135<table class="manual">
136<tr><td>Standard \n dimensions</td><td>\code
137mat.rows()
138mat.cols()\endcode</td>
139<td>\code
140vec.size() \endcode</td>
141</tr>
142<tr><td>Sizes along the \n inner/outer dimensions</td><td>\code
143mat.innerSize()
144mat.outerSize()\endcode</td>
145<td></td>
146</tr>
147<tr><td>Number of non \n zero coefficients</td><td>\code
148mat.nonZeros() \endcode</td>
149<td>\code
150vec.nonZeros() \endcode</td></tr>
151</table>
152
153
154\b Iterating \b over \b the \b nonzero \b coefficients \n
155
156Random access to the elements of a sparse object can be done through the \c coeffRef(i,j) function.
157However, this function involves a quite expensive binary search.
158In 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.
159Here is an example:
160<table class="manual">
161<tr><td>
162\code
163SparseMatrix<double> mat(rows,cols);
164for (int k=0; k<mat.outerSize(); ++k)
165 for (SparseMatrix<double>::InnerIterator it(mat,k); it; ++it)
166 {
167 it.value();
168 it.row(); // row index
169 it.col(); // col index (here it is equal to k)
170 it.index(); // inner index, here it is equal to it.row()
171 }
172\endcode
173</td><td>
174\code
175SparseVector<double> vec(size);
176for (SparseVector<double>::InnerIterator it(vec); it; ++it)
177{
178 it.value(); // == vec[ it.index() ]
179 it.index();
180}
181\endcode
182</td></tr>
183</table>
184For a writable expression, the referenced value can be modified using the valueRef() function.
185If the type of the sparse matrix or vector depends on a template parameter, then the \c typename keyword is
186required to indicate that \c InnerIterator denotes a type; see \ref TopicTemplateKeyword for details.
187
188
189\section TutorialSparseFilling Filling a sparse matrix
190
191Because of the special storage scheme of a SparseMatrix, special care has to be taken when adding new nonzero entries.
192For 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.
193
194The 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.
195
196Here is a typical usage example:
197\code
198typedef Eigen::Triplet<double> T;
199std::vector<T> tripletList;
200tripletList.reserve(estimation_of_entries);
201for(...)
202{
203 // ...
204 tripletList.push_back(T(i,j,v_ij));
205}
206SparseMatrixType mat(rows,cols);
207mat.setFromTriplets(tripletList.begin(), tripletList.end());
208// mat is ready to go!
209\endcode
210The \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().
211See the SparseMatrix::setFromTriplets() function and class Triplet for more details.
212
213
214In some cases, however, slightly higher performance, and lower memory consumption can be reached by directly inserting the non-zeros into the destination matrix.
215A typical scenario of this approach is illustrated bellow:
216\code
2171: SparseMatrix<double> mat(rows,cols); // default is column major
2182: mat.reserve(VectorXi::Constant(cols,6));
2193: for each i,j such that v_ij != 0
2204: mat.insert(i,j) = v_ij; // alternative: mat.coeffRef(i,j) += v_ij;
2215: mat.makeCompressed(); // optional
222\endcode
223
224- 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.
225- 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.
226- 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.
227- The line 5 suppresses the remaining empty space and transforms the matrix into a compressed column storage.
228
229
230
231\section TutorialSparseFeatureSet Supported operators and functions
232
233Because of their special storage format, sparse matrices cannot offer the same level of flexibility than dense matrices.
234In Eigen's sparse module we chose to expose only the subset of the dense matrix API which can be efficiently implemented.
235In the following \em sm denotes a sparse matrix, \em sv a sparse vector, \em dm a dense matrix, and \em dv a dense vector.
236
237\subsection TutorialSparse_BasicOps Basic operations
238
239%Sparse expressions support most of the unary and binary coefficient wise operations:
240\code
241sm1.real() sm1.imag() -sm1 0.5*sm1
242sm1+sm2 sm1-sm2 sm1.cwiseProduct(sm2)
243\endcode
244However, <strong>a strong restriction is that the storage orders must match</strong>. For instance, in the following example:
245\code
246sm4 = sm1 + sm2 + sm3;
247\endcode
248sm1, sm2, and sm3 must all be row-major or all column-major.
249On the other hand, there is no restriction on the target matrix sm4.
250For 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:
251\code
252SparseMatrix<double> A, B;
253B = SparseMatrix<double>(A.transpose()) + A;
254\endcode
255
256Some binary coefficient-wise operators can also mix sparse and dense expressions:
257\code
258sm2 = sm1.cwiseProduct(dm1);
259dm1 += sm1;
260\endcode
261
262However, it is not yet possible to add a sparse and a dense matrix as in <tt>dm2 = sm1 + dm1</tt>.
263Please write this as the equivalent <tt>dm2 = dm1; dm2 += sm1</tt> (we plan to lift this restriction
264in the next release of %Eigen).
265
266%Sparse expressions also support transposition:
267\code
268sm1 = sm2.transpose();
269sm1 = sm2.adjoint();
270\endcode
271However, there is no transposeInPlace() method.
272
273
274\subsection TutorialSparse_Products Matrix products
275
276%Eigen supports various kind of sparse matrix products which are summarize below:
277 - \b sparse-dense:
278 \code
279dv2 = sm1 * dv1;
280dm2 = dm1 * sm1.adjoint();
281dm2 = 2. * sm1 * dm1;
282 \endcode
283 - \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():
284 \code
285dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored
286dm2 = A.selfadjointView<Upper>() * dm1; // if only the upper part of A is stored
287dm2 = A.selfadjointView<Lower>() * dm1; // if only the lower part of A is stored
288 \endcode
289 - \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:
290 \code
291sm3 = sm1 * sm2;
292sm3 = 4 * sm1.adjoint() * sm2;
293 \endcode
294 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:
295 \code
296sm3 = (sm1 * sm2).pruned(); // removes numerical zeros
297sm3 = (sm1 * sm2).pruned(ref); // removes elements much smaller than ref
298sm3 = (sm1 * sm2).pruned(ref,epsilon); // removes elements smaller than ref*epsilon
299 \endcode
300
301 - \b permutations. Finally, permutations can be applied to sparse matrices too:
302 \code
303PermutationMatrix<Dynamic,Dynamic> P = ...;
304sm2 = P * sm1;
305sm2 = sm1 * P.inverse();
306sm2 = sm1.transpose() * P;
307 \endcode
308
309
310\subsection TutorialSparse_SubMatrices Block operations
311
312Regarding read-access, sparse matrices expose the same API than for dense matrices to access to sub-matrices such as blocks, columns, and rows. See \ref TutorialBlockOperations for a detailed introduction.
313However, for performance reasons, writing to a sub-sparse-matrix is much more limited, and currently only contiguous sets of columns (resp. rows) of a column-major (resp. row-major) SparseMatrix are writable. Moreover, this information has to be known at compile-time, leaving out methods such as <tt>block(...)</tt> and <tt>corner*(...)</tt>. The available API for write-access to a SparseMatrix are summarized below:
314\code
315SparseMatrix<double,ColMajor> sm1;
316sm1.col(j) = ...;
317sm1.leftCols(ncols) = ...;
318sm1.middleCols(j,ncols) = ...;
319sm1.rightCols(ncols) = ...;
320
321SparseMatrix<double,RowMajor> sm2;
322sm2.row(i) = ...;
323sm2.topRows(nrows) = ...;
324sm2.middleRows(i,nrows) = ...;
325sm2.bottomRows(nrows) = ...;
326\endcode
327
328In addition, sparse matrices expose the SparseMatrixBase::innerVector() and SparseMatrixBase::innerVectors() methods, which are aliases to the col/middleCols methods for a column-major storage, and to the row/middleRows methods for a row-major storage.
329
330\subsection TutorialSparse_TriangularSelfadjoint Triangular and selfadjoint views
331
332Just 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:
333\code
334dm2 = sm1.triangularView<Lower>(dm1);
335dv2 = sm1.transpose().triangularView<Upper>(dv1);
336\endcode
337
338The selfadjointView() function permits various operations:
339 - optimized sparse-dense matrix products:
340 \code
341dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored
342dm2 = A.selfadjointView<Upper>() * dm1; // if only the upper part of A is stored
343dm2 = A.selfadjointView<Lower>() * dm1; // if only the lower part of A is stored
344 \endcode
345 - copy of triangular parts:
346 \code
347sm2 = sm1.selfadjointView<Upper>(); // makes a full selfadjoint matrix from the upper triangular part
348sm2.selfadjointView<Lower>() = sm1.selfadjointView<Upper>(); // copies the upper triangular part to the lower triangular part
349 \endcode
350 - application of symmetric permutations:
351 \code
352PermutationMatrix<Dynamic,Dynamic> P = ...;
353sm2 = A.selfadjointView<Upper>().twistedBy(P); // compute P S P' from the upper triangular part of A, and make it a full matrix
354sm2.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
355 \endcode
356
357Please, refer to the \link SparseQuickRefPage Quick Reference \endlink guide for the list of supported operations. The list of linear solvers available is \link TopicSparseSystems here. \endlink
358
359*/
360
361}
Note: See TracBrowser for help on using the repository browser.