Merged in f5soh/librepilot/update_credits (pull request #529)
[librepilot.git] / ground / gcs / src / libs / eigen / doc / SparseQuickReference.dox
bloba25622e800e0df278456d6d99badc47d0807a993
1 namespace Eigen {
2 /** \eigenManualPage SparseQuickRefPage Quick reference guide for sparse matrices
3 \eigenAutoToc
5 <hr>
7 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 : 
8 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. 
10 \section SparseMatrixInit Sparse Matrix Initialization
11 <table class="manual">
12 <tr><th> Category </th> <th> Operations</th> <th>Notes</th></tr>
13 <tr><td>Constructor</td>
14 <td>
15 \code
16   SparseMatrix<double> sm1(1000,1000); 
17   SparseMatrix<std::complex<double>,RowMajor> sm2;
18 \endcode
19 </td> <td> Default is ColMajor</td> </tr>
20 <tr class="alt">
21 <td> Resize/Reserve</td>
22 <td> 
23  \code
24     sm1.resize(m,n);      // Change sm1 to a m x n matrix.
25     sm1.reserve(nnz);     // Allocate room for nnz nonzeros elements.   
26   \endcode 
27 </td>
28 <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>
29 </tr>
30 <tr> 
31 <td> Assignment </td>
32 <td> 
33 \code 
34   SparseMatrix<double,Colmajor> sm1;
35  // Initialize sm2 with sm1.
36   SparseMatrix<double,Rowmajor> sm2(sm1), sm3;        
37   // Assignment and evaluations modify the storage order.
38   sm3 = sm1; 
39  \endcode
40 </td>
41 <td> The copy constructor can be used to convert from a storage order to another</td>
42 </tr>
43 <tr class="alt">
44 <td> Element-wise Insertion</td>
45 <td>
46 \code 
47 // Insert a new element; 
48  sm1.insert(i, j) = v_ij;  
50 // Update the value v_ij
51  sm1.coeffRef(i,j) = v_ij;
52  sm1.coeffRef(i,j) += v_ij;
53  sm1.coeffRef(i,j) -= v_ij;
54 \endcode
55 </td>
56 <td> insert() assumes that the element does not already exist; otherwise, use coeffRef()</td>
57 </tr>
58 <tr> 
59 <td> Batch insertion</td>
60 <td>
61 \code
62   std::vector< Eigen::Triplet<double> > tripletList;
63   tripletList.reserve(estimation_of_entries);
64   // -- Fill tripletList with nonzero elements...
65   sm1.setFromTriplets(TripletList.begin(), TripletList.end());
66 \endcode
67 </td>
68 <td>A complete example is available at \link TutorialSparseFilling Triplet Insertion \endlink.</td>
69 </tr>
70 <tr class="alt"> 
71 <td> Constant or Random Insertion</td>
72 <td>
73 \code
74 sm1.setZero();
75 \endcode
76 </td>
77 <td>Remove all non-zero coefficients</td>
78 </tr>
79 </table>
82 \section SparseBasicInfos Matrix properties
83 Beyond the basic functions rows() and cols(), there are some useful functions that are available to easily get some informations from the matrix. 
84 <table class="manual">
85 <tr>
86   <td> \code
87   sm1.rows();         // Number of rows
88   sm1.cols();         // Number of columns 
89   sm1.nonZeros();     // Number of non zero values   
90   sm1.outerSize();    // Number of columns (resp. rows) for a column major (resp. row major )
91   sm1.innerSize();    // Number of rows (resp. columns) for a row major (resp. column major)
92   sm1.norm();         // Euclidian norm of the matrix
93   sm1.squaredNorm();  // Squared norm of the matrix
94   sm1.blueNorm();
95   sm1.isVector();     // Check if sm1 is a sparse vector or a sparse matrix
96   sm1.isCompressed(); // Check if sm1 is in compressed form
97   ...
98   \endcode </td>
99 </tr>
100 </table>
102 \section SparseBasicOps Arithmetic operations
103 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.
104 <table class="manual">
105 <tr><th> Operations </th> <th> Code </th> <th> Notes </th></tr>
107 <tr>
108   <td> add subtract </td> 
109   <td> \code
110   sm3 = sm1 + sm2; 
111   sm3 = sm1 - sm2;
112   sm2 += sm1; 
113   sm2 -= sm1; \endcode
114   </td>
115   <td> 
116   sm1 and sm2 should have the same storage order
117   </td> 
118 </tr>
120 <tr class="alt"><td>
121   scalar product</td><td>\code
122   sm3 = sm1 * s1;   sm3 *= s1; 
123   sm3 = s1 * sm1 + s2 * sm2; sm3 /= s1;\endcode
124   </td>
125   <td>
126     Many combinations are possible if the dimensions and the storage order agree.
127 </tr>
129 <tr>
130   <td> %Sparse %Product </td>
131   <td> \code
132   sm3 = sm1 * sm2;
133   dm2 = sm1 * dm1;
134   dv2 = sm1 * dv1;
135   \endcode </td>
136   <td>
137   </td>
138 </tr> 
140 <tr class='alt'>
141   <td> transposition, adjoint</td>
142   <td> \code
143   sm2 = sm1.transpose();
144   sm2 = sm1.adjoint();
145   \endcode </td>
146   <td>
147   Note that the transposition change the storage order. There is no support for transposeInPlace().
148   </td>
149 </tr> 
150 <tr>
151 <td> Permutation </td>
152 <td> 
153 \code 
154 perm.indices();      // Reference to the vector of indices
155 sm1.twistedBy(perm); // Permute rows and columns
156 sm2 = sm1 * perm;    // Permute the columns
157 sm2 = perm * sm1;    // Permute the columns
158 \endcode 
159 </td>
160 <td> 
162 </td>
163 </tr>
164 <tr>
165   <td>
166   Component-wise ops
167   </td>
168   <td>\code 
169   sm1.cwiseProduct(sm2);
170   sm1.cwiseQuotient(sm2);
171   sm1.cwiseMin(sm2);
172   sm1.cwiseMax(sm2);
173   sm1.cwiseAbs();
174   sm1.cwiseSqrt();
175   \endcode</td>
176   <td>
177   sm1 and sm2 should have the same storage order
178   </td>
179 </tr>
180 </table>
182 \section sparseotherops Other supported operations
183 <table class="manual">
184 <tr><th style="min-width:initial"> Code </th> <th> Notes</th> </tr>
185 <tr><td colspan="2">Sub-matrices</td></tr>
186 <tr>
187 <td> 
188 \code 
189   sm1.block(startRow, startCol, rows, cols); 
190   sm1.block(startRow, startCol); 
191   sm1.topLeftCorner(rows, cols); 
192   sm1.topRightCorner(rows, cols);
193   sm1.bottomLeftCorner( rows, cols);
194   sm1.bottomRightCorner( rows, cols);
195   \endcode
196 </td><td>
197 Contrary to dense matrices, here <strong>all these methods are read-only</strong>.\n
198 See \ref TutorialSparse_SubMatrices and below for read-write sub-matrices.
199 </td>
200 </tr>
201 <tr class="alt"><td colspan="2"> Range </td></tr>
202 <tr class="alt">
203 <td> 
204 \code 
205   sm1.innerVector(outer);           // RW
206   sm1.innerVectors(start, size);    // RW
207   sm1.leftCols(size);               // RW
208   sm2.rightCols(size);              // RO because sm2 is row-major
209   sm1.middleRows(start, numRows);   // RO because sm1 is column-major
210   sm1.middleCols(start, numCols);   // RW
211   sm1.col(j);                       // RW
212 \endcode
213 </td>
214 <td>
215 A inner vector is either a row (for row-major) or a column (for column-major).\n
216 As stated earlier, for a read-write sub-matrix (RW), the evaluation can be done in a matrix with different storage order.
217 </td>
218 </tr>
219 <tr><td colspan="2"> Triangular and selfadjoint views</td></tr>
220 <tr>
221 <td> 
222 \code
223   sm2 = sm1.triangularview<Lower>();
224   sm2 = sm1.selfadjointview<Lower>();
225 \endcode
226 </td>
227 <td> Several combination between triangular views and blocks views are possible
228 \code 
229   \endcode </td>
230 </tr>
231 <tr class="alt"><td colspan="2">Triangular solve </td></tr>
232 <tr class="alt">
233 <td> 
234 \code 
235  dv2 = sm1.triangularView<Upper>().solve(dv1);
236  dv2 = sm1.topLeftCorner(size, size)
237           .triangularView<Lower>().solve(dv1);
238 \endcode 
239 </td>
240 <td> For general sparse solve, Use any suitable module described at \ref TopicSparseSystems </td>
241 </tr>
242 <tr><td colspan="2"> Low-level API</td></tr>
243 <tr>
244 <td>
245 \code
246 sm1.valuePtr();      // Pointer to the values
247 sm1.innerIndextr();  // Pointer to the indices.
248 sm1.outerIndexPtr(); // Pointer to the beginning of each inner vector
249 \endcode
250 </td>
251 <td>
252 If the matrix is not in compressed form, makeCompressed() should be called before.\n
253 Note that these functions are mostly provided for interoperability purposes with external libraries.\n
254 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>
255 </tr>
256 <tr class="alt"><td colspan="2">Mapping external buffers</td></tr>
257 <tr class="alt">
258 <td>
259 \code
260 int outerIndexPtr[cols+1];
261 int innerIndices[nnz];
262 double values[nnz];
263 Map<SparseMatrix<double> > sm1(rows,cols,nnz,outerIndexPtr, // read-write
264                                innerIndices,values);
265 Map<const SparseMatrix<double> > sm2(...);                  // read-only
266 \endcode
267 </td>
268 <td>As for dense matrices, class Map<SparseMatrixType> can be used to see external buffers as an %Eigen's SparseMatrix object. </td>
269 </tr>
270 </table>