Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 1 | namespace Eigen { |
| 2 | |
Gael Guennebaud | 41ea92d | 2010-07-04 10:14:47 +0200 | [diff] [blame] | 3 | /** \page TutorialSparse Tutorial page 9 - Sparse Matrix |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 4 | \ingroup Tutorial |
| 5 | |
Gael Guennebaud | 41ea92d | 2010-07-04 10:14:47 +0200 | [diff] [blame] | 6 | \li \b Previous: \ref TutorialGeometry |
Tim Holy | 44b19b4 | 2012-02-08 22:11:12 +0100 | [diff] [blame] | 7 | \li \b Next: \ref TutorialMapClass |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 8 | |
| 9 | \b Table \b of \b contents \n |
| 10 | - \ref TutorialSparseIntro |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 11 | - \ref TutorialSparseExample "Example" |
| 12 | - \ref TutorialSparseSparseMatrix |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 13 | - \ref TutorialSparseFilling |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 14 | - \ref TutorialSparseDirectSolvers |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 15 | - \ref TutorialSparseFeatureSet |
| 16 | - \ref TutorialSparse_BasicOps |
| 17 | - \ref TutorialSparse_Products |
| 18 | - \ref TutorialSparse_TriangularSelfadjoint |
| 19 | - \ref TutorialSparse_Submat |
| 20 | |
| 21 | |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 22 | <hr> |
| 23 | |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 24 | Manipulating and solving sparse problems involves various modules which are summarized below: |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 25 | |
Gael Guennebaud | f66fe26 | 2010-10-19 11:40:49 +0200 | [diff] [blame] | 26 | <table class="manual"> |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 27 | <tr><th>Module</th><th>Header file</th><th>Contents</th></tr> |
| 28 | <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> |
| 29 | <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> |
| 30 | <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> |
| 31 | <tr><td></td><td>\code#include <Eigen/Sparse>\endcode</td><td>Includes all the above modules</td></tr> |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 32 | </table> |
| 33 | |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 34 | \section TutorialSparseIntro Sparse matrix representation |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 35 | |
| 36 | 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. |
| 37 | |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 38 | \b The \b %SparseMatrix \b class |
| 39 | |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 40 | The class SparseMatrix is the main sparse matrix representation of Eigen's sparse module; it offers high performance and low memory usage. |
| 41 | It implements a more versatile variant of the widely-used Compressed Column (or Row) Storage scheme. |
| 42 | It consists of four compact arrays: |
| 43 | - \c Values: stores the coefficient values of the non-zeros. |
| 44 | - \c InnerIndices: stores the row (resp. column) indices of the non-zeros. |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 45 | - \c OuterStarts: stores for each colmun (resp. row) the index of the first non zero in the previous two arrays. |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 46 | - \c InnerNNZs: stores the number of non-zeros of each column (resp. row). |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 47 | 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. |
| 48 | The word \c outer refers to the other direction. |
| 49 | |
| 50 | This storage scheme is better explained on an example. The following matrix |
| 51 | <table class="manual"> |
| 52 | <tr><td> 0</td><td>3</td><td> 0</td><td>0</td><td> 0</td></tr> |
| 53 | <tr><td>22</td><td>0</td><td> 0</td><td>0</td><td>17</td></tr> |
| 54 | <tr><td> 7</td><td>5</td><td> 0</td><td>1</td><td> 0</td></tr> |
| 55 | <tr><td> 0</td><td>0</td><td> 0</td><td>0</td><td> 0</td></tr> |
| 56 | <tr><td> 0</td><td>0</td><td>14</td><td>0</td><td> 8</td></tr> |
| 57 | </table> |
| 58 | |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 59 | and one of its possible sparse, \b column \b major representation: |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 60 | <table class="manual"> |
| 61 | <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> |
| 62 | <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> |
| 63 | </table> |
| 64 | <table class="manual"> |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 65 | <tr><td>OuterStarts:</td><td>0</td><td>3</td><td>5</td><td>8</td><td>10</td><td>\em 12 </td></tr> |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 66 | <tr><td>InnerNNZs:</td> <td>2</td><td>2</td><td>1</td><td>1</td><td> 2</td><td></td></tr> |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 67 | </table> |
| 68 | |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 69 | Currently the elements of a given inner vector are guaranteed to be always sorted by increasing inner indices. |
| 70 | The \c "_" indicates available free space to quickly insert new elements. |
| 71 | 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. |
| 72 | 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. |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 73 | |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 74 | The case where no empty space is available is a special case, and is refered as the \em compressed mode. |
| 75 | It corresponds to the widely used Compressed Column (or Row) Storage schemes (CCS or CRS). |
| 76 | Any SparseMatrix can be turned to this form by calling the SparseMatrix::makeCompressed() function. |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 77 | 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]. |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 78 | Therefore, in practice a call to SparseMatrix::makeCompressed() frees this buffer. |
| 79 | |
| 80 | It is worth noting that most of our wrappers to external libraries requires compressed matrices as inputs. |
| 81 | |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 82 | The results of %Eigen's operations always produces \b compressed sparse matrices. |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 83 | On the other hand, the insertion of a new element into a SparseMatrix converts this later to the \b uncompressed mode. |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 84 | |
| 85 | Here is the previous matrix represented in compressed mode: |
Gael Guennebaud | f66fe26 | 2010-10-19 11:40:49 +0200 | [diff] [blame] | 86 | <table class="manual"> |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 87 | <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> |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 88 | <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> |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 89 | </table> |
Gael Guennebaud | f66fe26 | 2010-10-19 11:40:49 +0200 | [diff] [blame] | 90 | <table class="manual"> |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 91 | <tr><td>OuterStarts:</td><td>0</td><td>2</td><td>4</td><td>5</td><td>6</td><td>\em 8 </td></tr> |
Gael Guennebaud | 25f1658 | 2009-01-21 17:44:58 +0000 | [diff] [blame] | 92 | </table> |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 93 | |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 94 | A SparseVector is a special case of a SparseMatrix where only the \c Values and \c InnerIndices arrays are stored. |
| 95 | There is no notion of compressed/uncompressed mode for a SparseVector. |
Gael Guennebaud | 25f1658 | 2009-01-21 17:44:58 +0000 | [diff] [blame] | 96 | |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 97 | |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 98 | \section TutorialSparseExample First example |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 99 | |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 100 | Before describing each individual classes, lets start with the following typical example solving for the Lapace equation \f$ \nabla u = 0 \f$ onto a regular 2D grid using a finite difference scheme and Dirichlet boundary conditions. |
| 101 | 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 unknows (in our case, the values of the pixels), \f$ b \f$ is the right hand side vector resuting from the boundary conditions, and \f$ A \f$ is a \f$ m \times m \f$ matrix containing only a few non-zeros elements resulting from the discretization of the Laplacian operator. |
Gael Guennebaud | 25f1658 | 2009-01-21 17:44:58 +0000 | [diff] [blame] | 102 | |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 103 | <table class="manual"> |
| 104 | <tr><td> |
| 105 | \include Tutorial_sparse_example.cpp |
| 106 | </td> |
| 107 | <td> |
| 108 | \image html Tutorial_sparse_example.jpeg |
| 109 | </td></tr></table> |
| 110 | |
| 111 | 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. |
| 112 | |
| 113 | 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. |
| 114 | The raw and flat list of non-zero entries is then converted to a true SparseMatrix object \c A. |
| 115 | Note that the elements of the list do not have to be sorted, and possible duplicate entries will be summed up. |
| 116 | |
| 117 | The last step consists in effectively solving the such assembled problem. |
| 118 | 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 counter part for dense objects. |
| 119 | |
| 120 | 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. |
| 121 | |
| 122 | Commenting 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. |
| 123 | |
| 124 | |
| 125 | |
| 126 | |
| 127 | \section TutorialSparseSparseMatrix The SparseMatrix class |
| 128 | |
| 129 | \b %Matrix \b and \b vector \b properties \n |
| 130 | |
| 131 | The SparseMatrix and SparseVector classes take three template arguments: |
| 132 | * the scalar type (e.g., double) |
| 133 | * the storage order (ColMajor or RowMajor, the default is RowMajor) |
| 134 | * the inner index type (default is \c int). |
| 135 | |
| 136 | As for dense Matrix objects, constructors takes the size of the object. |
| 137 | Here are some examples: |
| 138 | |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 139 | \code |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 140 | SparseMatrix<std::complex<float> > mat(1000,2000); // declares a 1000x2000 column-major compressed sparse matrix of complex<float> |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 141 | SparseMatrix<double,RowMajor> mat(1000,2000); // declares a 1000x2000 row-major compressed sparse matrix of double |
| 142 | SparseVector<std::complex<float> > vec(1000); // declares a column sparse vector of complex<float> of size 1000 |
| 143 | SparseVector<double,RowMajor> vec(1000); // declares a row sparse vector of double of size 1000 |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 144 | \endcode |
| 145 | |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 146 | In the rest of the tutorial, \c mat and \c vec represent any sparse-matrix and sparse-vector objects, respectively. |
| 147 | |
| 148 | The dimensions of a matrix can be queried using the following functions: |
Gael Guennebaud | f66fe26 | 2010-10-19 11:40:49 +0200 | [diff] [blame] | 149 | <table class="manual"> |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 150 | <tr><td>Standard \n dimensions</td><td>\code |
| 151 | mat.rows() |
| 152 | mat.cols()\endcode</td> |
| 153 | <td>\code |
| 154 | vec.size() \endcode</td> |
| 155 | </tr> |
| 156 | <tr><td>Sizes along the \n inner/outer dimensions</td><td>\code |
| 157 | mat.innerSize() |
| 158 | mat.outerSize()\endcode</td> |
| 159 | <td></td> |
| 160 | </tr> |
Tim Holy | 16a2d89 | 2011-06-20 22:47:58 -0500 | [diff] [blame] | 161 | <tr><td>Number of non \n zero coefficients</td><td>\code |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 162 | mat.nonZeros() \endcode</td> |
| 163 | <td>\code |
| 164 | vec.nonZeros() \endcode</td></tr> |
| 165 | </table> |
| 166 | |
Gael Guennebaud | 25f1658 | 2009-01-21 17:44:58 +0000 | [diff] [blame] | 167 | |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 168 | \b Iterating \b over \b the \b nonzero \b coefficients \n |
| 169 | |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 170 | Random access to the elements of a sparse object can be done through the \c coeffRef(i,j) function. |
| 171 | However, this function involves a quite expensive binary search. |
| 172 | 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. |
| 173 | Here is an example: |
Gael Guennebaud | f66fe26 | 2010-10-19 11:40:49 +0200 | [diff] [blame] | 174 | <table class="manual"> |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 175 | <tr><td> |
| 176 | \code |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 177 | SparseMatrix<double> mat(rows,cols); |
Jitse Niesen | 59b83c1 | 2011-09-10 09:18:18 +0100 | [diff] [blame] | 178 | for (int k=0; k<mat.outerSize(); ++k) |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 179 | for (SparseMatrix<double>::InnerIterator it(mat,k); it; ++it) |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 180 | { |
| 181 | it.value(); |
| 182 | it.row(); // row index |
| 183 | it.col(); // col index (here it is equal to k) |
| 184 | it.index(); // inner index, here it is equal to it.row() |
| 185 | } |
| 186 | \endcode |
| 187 | </td><td> |
| 188 | \code |
Gael Guennebaud | 25f1658 | 2009-01-21 17:44:58 +0000 | [diff] [blame] | 189 | SparseVector<double> vec(size); |
| 190 | for (SparseVector<double>::InnerIterator it(vec); it; ++it) |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 191 | { |
Gael Guennebaud | 25f1658 | 2009-01-21 17:44:58 +0000 | [diff] [blame] | 192 | it.value(); // == vec[ it.index() ] |
Gael Guennebaud | c532f42 | 2009-09-25 16:30:44 +0200 | [diff] [blame] | 193 | it.index(); |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 194 | } |
| 195 | \endcode |
| 196 | </td></tr> |
| 197 | </table> |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 198 | For a writable expression, the referenced value can be modified using the valueRef() function. |
Jitse Niesen | 59b83c1 | 2011-09-10 09:18:18 +0100 | [diff] [blame] | 199 | If the type of the sparse matrix or vector depends on a template parameter, then the \c typename keyword is |
| 200 | required to indicate that \c InnerIterator denotes a type; see \ref TopicTemplateKeyword for details. |
| 201 | |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 202 | |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 203 | \section TutorialSparseFilling Filling a sparse matrix |
| 204 | |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 205 | Because of the special storage scheme of a SparseMatrix, special care has to be taken when adding new nonzero entries. |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 206 | 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. |
Gael Guennebaud | c532f42 | 2009-09-25 16:30:44 +0200 | [diff] [blame] | 207 | |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 208 | The simplest way to create a sparse matrix while guarantying good performance is thus to first build a list of so called \em triplets, and then convert it to a SparseMatrix. |
Gael Guennebaud | 1763f86 | 2012-02-04 10:44:07 +0100 | [diff] [blame] | 209 | |
| 210 | Here is a typical usage example: |
| 211 | \code |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 212 | typedef Eigen::Triplet<double> T; |
Gael Guennebaud | 1763f86 | 2012-02-04 10:44:07 +0100 | [diff] [blame] | 213 | std::vector<T> tripletList; |
| 214 | triplets.reserve(estimation_of_entries); |
| 215 | for(...) |
| 216 | { |
| 217 | // ... |
| 218 | tripletList.push_back(T(i,j,v_ij)); |
| 219 | } |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 220 | SparseMatrixType mat(rows,cols); |
| 221 | mat.setFromTriplets(tripletList.begin(), tripletList.end()); |
| 222 | // mat is ready to go! |
Gael Guennebaud | 1763f86 | 2012-02-04 10:44:07 +0100 | [diff] [blame] | 223 | \endcode |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 224 | 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(). |
Gael Guennebaud | 1763f86 | 2012-02-04 10:44:07 +0100 | [diff] [blame] | 225 | See the SparseMatrix::setFromTriplets() function and class Triplet for more details. |
| 226 | |
| 227 | |
| 228 | In some cases, however, slightly higher performance, and lower memory consumption can be reached by directly inserting the non zeros into the destination matrix. |
| 229 | A typical scenario of this approach is illustrated bellow: |
Gael Guennebaud | 25f1658 | 2009-01-21 17:44:58 +0000 | [diff] [blame] | 230 | \code |
Gael Guennebaud | c861e05 | 2011-12-03 18:14:51 +0100 | [diff] [blame] | 231 | 1: SparseMatrix<double> mat(rows,cols); // default is column major |
Gael Guennebaud | 1aa6c7f | 2011-12-11 17:18:14 +0100 | [diff] [blame] | 232 | 2: mat.reserve(VectorXi::Constant(cols,6)); |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 233 | 3: for each i,j such that v_ij != 0 |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 234 | 4: mat.insert(i,j) = v_ij; // alternative: mat.coeffRef(i,j) += v_ij; |
| 235 | 5: mat.makeCompressed(); // optional |
Gael Guennebaud | 25f1658 | 2009-01-21 17:44:58 +0000 | [diff] [blame] | 236 | \endcode |
| 237 | |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 238 | - 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 zero 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. |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 239 | - 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. |
| 240 | - 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. |
| 241 | - The line 5 suppresses the remaining empty space and transforms the matrix into a compressed column storage. |
Gael Guennebaud | 25f1658 | 2009-01-21 17:44:58 +0000 | [diff] [blame] | 242 | |
Gael Guennebaud | 25f1658 | 2009-01-21 17:44:58 +0000 | [diff] [blame] | 243 | |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 244 | \section TutorialSparseDirectSolvers Solving linear problems |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 245 | |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 246 | %Eigen currently provides a limited set of built-in solvers, as well as wrappers to external solver libraries. |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 247 | They are summarized in the following table: |
Gael Guennebaud | ab41c18 | 2010-08-18 15:33:58 +0200 | [diff] [blame] | 248 | |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 249 | <table class="manual"> |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 250 | <tr><th>Class</th><th>Module</th><th>Solver kind</th><th>Matrix kind</th><th>Features related to performance</th> |
| 251 | <th>Dependencies,License</th><th class="width20em"><p>Notes</p></th></tr> |
| 252 | <tr><td>SimplicialLLT </td><td>\link SparseCholesky_Module SparseCholesky \endlink</td><td>Direct LLt factorization</td><td>SPD</td><td>Fill-in reducing</td> |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 253 | <td>built-in, LGPL</td> |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 254 | <td>SimplicialLDLT is often preferable</td></tr> |
| 255 | <tr><td>SimplicialLDLT </td><td>\link SparseCholesky_Module SparseCholesky \endlink</td><td>Direct LDLt factorization</td><td>SPD</td><td>Fill-in reducing</td> |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 256 | <td>built-in, LGPL</td> |
| 257 | <td>Recommended for very sparse and not too large problems (e.g., 2D Poisson eq.)</td></tr> |
| 258 | <tr><td>ConjugateGradient</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>Classic iterative CG</td><td>SPD</td><td>Preconditionning</td> |
| 259 | <td>built-in, LGPL</td> |
| 260 | <td>Recommended for large symmetric problems (e.g., 3D Poisson eq.)</td></tr> |
| 261 | <tr><td>BiCGSTAB</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>Iterative stabilized bi-conjugate gradient</td><td>Square</td><td>Preconditionning</td> |
| 262 | <td>built-in, LGPL</td> |
| 263 | <td>Might not always converge</td></tr> |
| 264 | |
| 265 | |
| 266 | |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 267 | <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> |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 268 | <td>Requires the <a href="http://www.cise.ufl.edu/research/sparse/SuiteSparse/">SuiteSparse</a> package, \b GPL </td> |
| 269 | <td></td></tr> |
| 270 | <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> |
| 271 | <td>Requires the <a href="http://www.cise.ufl.edu/research/sparse/SuiteSparse/">SuiteSparse</a> package, \b GPL </td> |
| 272 | <td></td></tr> |
| 273 | <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> |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 274 | <td>Requires the <a href="http://crd-legacy.lbl.gov/~xiaoye/SuperLU/">SuperLU</a> library, (BSD-like)</td> |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 275 | <td></td></tr> |
| 276 | </table> |
| 277 | |
| 278 | Here \c SPD means symmetric positive definite. |
| 279 | |
| 280 | All these solvers follow the same general concept. |
| 281 | Here is a typical and general example: |
Gael Guennebaud | ab41c18 | 2010-08-18 15:33:58 +0200 | [diff] [blame] | 282 | \code |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 283 | #include <Eigen/RequiredModuleName> |
Gael Guennebaud | ab41c18 | 2010-08-18 15:33:58 +0200 | [diff] [blame] | 284 | // ... |
| 285 | SparseMatrix<double> A; |
| 286 | // fill A |
| 287 | VectorXd b, x; |
| 288 | // fill b |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 289 | // solve Ax = b |
| 290 | SolverClassName<SparseMatrix<double> > solver; |
| 291 | solver.compute(A); |
| 292 | if(solver.info()!=Succeeded) { |
Tim Holy | 16a2d89 | 2011-06-20 22:47:58 -0500 | [diff] [blame] | 293 | // decomposition failed |
Gael Guennebaud | ab41c18 | 2010-08-18 15:33:58 +0200 | [diff] [blame] | 294 | return; |
| 295 | } |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 296 | x = solver.solve(b); |
| 297 | if(solver.info()!=Succeeded) { |
Gael Guennebaud | ab41c18 | 2010-08-18 15:33:58 +0200 | [diff] [blame] | 298 | // solving failed |
| 299 | return; |
| 300 | } |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 301 | // solve for another right hand side: |
| 302 | x1 = solver.solve(b1); |
Gael Guennebaud | ab41c18 | 2010-08-18 15:33:58 +0200 | [diff] [blame] | 303 | \endcode |
| 304 | |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 305 | For \c SPD solvers, a second optional template argument allows to specify which triangular part have to be used, e.g.: |
| 306 | |
| 307 | \code |
| 308 | #include <Eigen/IterativeLinearSolvers> |
| 309 | |
| 310 | ConjugateGradient<SparseMatrix<double>, Eigen::Upper> solver; |
| 311 | x = solver.compute(A).solve(b); |
| 312 | \endcode |
| 313 | 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. |
| 314 | |
| 315 | In the case where multiple problems with the same sparcity pattern have to be solved, then the "compute" step can be decomposed as follow: |
| 316 | \code |
| 317 | SolverClassName<SparseMatrix<double> > solver; |
| 318 | solver.analyzePattern(A); // for this step the numerical values of A are not used |
| 319 | solver.factorize(A); |
| 320 | x1 = solver.solve(b1); |
| 321 | x2 = solver.solve(b2); |
| 322 | ... |
| 323 | A = ...; // modify the values of the nonzeros of A, the nonzeros pattern must stay unchanged |
| 324 | solver.factorize(A); |
| 325 | x1 = solver.solve(b1); |
| 326 | x2 = solver.solve(b2); |
| 327 | ... |
| 328 | \endcode |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 329 | The compute() method is equivalent to calling both analyzePattern() and factorize(). |
Gael Guennebaud | 70206ab | 2011-11-24 17:32:30 +0100 | [diff] [blame] | 330 | |
Gael Guennebaud | 3f56de2 | 2011-12-03 10:26:00 +0100 | [diff] [blame] | 331 | Finally, each solver provides some specific features, such as determinant, access to the factors, controls of the iterations, and so on. |
| 332 | More details are availble in the documentations of the respective classes. |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 333 | |
Gael Guennebaud | 52dce0c | 2012-06-20 09:28:32 +0200 | [diff] [blame] | 334 | |
| 335 | \section TutorialSparseFeatureSet Supported operators and functions |
| 336 | |
| 337 | Because of their special storage format, sparse matrices cannot offer the same level of flexbility than dense matrices. |
| 338 | In Eigen's sparse module we chose to expose only the subset of the dense matrix API which can be efficiently implemented. |
| 339 | 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. |
| 340 | |
| 341 | \subsection TutorialSparse_BasicOps Basic operations |
| 342 | |
| 343 | %Sparse expressions support most of the unary and binary coefficient wise operations: |
| 344 | \code |
| 345 | sm1.real() sm1.imag() -sm1 0.5*sm1 |
| 346 | sm1+sm2 sm1-sm2 sm1.cwiseProduct(sm2) |
| 347 | \endcode |
| 348 | However, a strong restriction is that the storage orders must match. For instance, in the following example: |
| 349 | \code |
| 350 | sm4 = sm1 + sm2 + sm3; |
| 351 | \endcode |
| 352 | sm1, sm2, and sm3 must all be row-major or all column major. |
| 353 | On the other hand, there is no restriction on the target matrix sm4. |
| 354 | 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: |
| 355 | \code |
| 356 | SparseMatrix<double> A, B; |
| 357 | B = SparseMatrix<double>(A.transpose()) + A; |
| 358 | \endcode |
| 359 | |
| 360 | Binary coefficient wise operators can also mix sparse and dense expressions: |
| 361 | \code |
| 362 | sm2 = sm1.cwiseProduct(dm1); |
| 363 | dm2 = sm1 + dm1; |
| 364 | \endcode |
| 365 | |
| 366 | |
| 367 | %Sparse expressions also support transposition: |
| 368 | \code |
| 369 | sm1 = sm2.transpose(); |
| 370 | sm1 = sm2.adjoint(); |
| 371 | \endcode |
| 372 | However, there is no transposeInPlace() method. |
| 373 | |
| 374 | |
| 375 | \subsection TutorialSparse_Products Matrix products |
| 376 | |
| 377 | %Eigen supports various kind of sparse matrix products which are summarize below: |
| 378 | - \b sparse-dense: |
| 379 | \code |
| 380 | dv2 = sm1 * dv1; |
| 381 | dm2 = dm1 * sm1.adjoint(); |
| 382 | dm2 = 2. * sm1 * dm1; |
| 383 | \endcode |
| 384 | - \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(): |
| 385 | \code |
| 386 | dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored |
| 387 | dm2 = A.selfadjointView<Upper>() * dm1; // if only the upper part of A is stored |
| 388 | dm2 = A.selfadjointView<Lower>() * dm1; // if only the lower part of A is stored |
| 389 | \endcode |
| 390 | - \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: |
| 391 | \code |
| 392 | sm3 = sm1 * sm2; |
| 393 | sm3 = 4 * sm1.adjoint() * sm2; |
| 394 | \endcode |
| 395 | The second algorithm punes on the fly the explicit zeros, or the values smaller than a given threshold. It is enabled and controlled through the prune() functions: |
| 396 | \code |
| 397 | sm3 = (sm1 * sm2).prune(); // removes numerical zeros |
| 398 | sm3 = (sm1 * sm2).prune(ref); // removes elements much smaller than ref |
| 399 | sm3 = (sm1 * sm2).prune(ref,epsilon); // removes elements much smaller than ref*epsilon |
| 400 | \endcode |
| 401 | |
| 402 | - \b permutations. Finally, permutations can be applied to sparse matrices too: |
| 403 | \code |
| 404 | PermutationMatrix<Dynamic,Dynamic> P = ...; |
| 405 | sm2 = P * sm1; |
| 406 | sm2 = sm1 * P.inverse(); |
| 407 | sm2 = sm1.transpose() * P; |
| 408 | \endcode |
| 409 | |
| 410 | |
| 411 | \subsection TutorialSparse_TriangularSelfadjoint Triangular and selfadjoint views |
| 412 | |
| 413 | Just as dense matrices, the triangularView() function can be used to address a triangular part of the matrix, and perform triangular solves with a dense right and side: |
| 414 | \code |
| 415 | dm2 = sm1.triangularView<Lower>(dm1); |
| 416 | dv2 = sm1.transpose().triangularView<Upper>(dv1); |
| 417 | \endcode |
| 418 | |
| 419 | The selfadjointView() function permits various operations: |
| 420 | - optimized sparse-dense matrix products: |
| 421 | \code |
| 422 | dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored |
| 423 | dm2 = A.selfadjointView<Upper>() * dm1; // if only the upper part of A is stored |
| 424 | dm2 = A.selfadjointView<Lower>() * dm1; // if only the lower part of A is stored |
| 425 | \endcode |
| 426 | - copy of triangular parts: |
| 427 | \code |
| 428 | sm2 = sm1.selfadjointView<Upper>(); // makes a full selfadjoint matrix from the upper triangular part |
| 429 | sm2.selfadjointView<Lower>() = sm1.selfadjointView<Upper>(); // copie the upper triangular part to the lower triangular part |
| 430 | \endcode |
| 431 | - application of symmetric permutations: |
| 432 | \code |
| 433 | PermutationMatrix<Dynamic,Dynamic> P = ...; |
| 434 | sm2 = A.selfadjointView<Upper>().twistedBy(P); // compute P S P' from the upper triangular part of A, and make it a full matrix |
| 435 | 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 |
| 436 | \endcode |
| 437 | |
| 438 | \subsection TutorialSparse_Submat Sub-matrices |
| 439 | |
| 440 | %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: |
| 441 | \code |
| 442 | sm1.innerVector(j); // returns an expression of the j-th column (resp. row) of the matrix if sm1 is col-major (resp. row-major) |
| 443 | sm1.innerVectors(j, nb); // returns an expression of the nb columns (resp. row) starting from the j-th column (resp. row) |
| 444 | // of the matrix if sm1 is col-major (resp. row-major) |
| 445 | sm1.middleRows(j, nb); // for row major matrices only, get a range of nb rows |
| 446 | sm1.middleCols(j, nb); // for column major matrices only, get a range of nb columns |
| 447 | \endcode |
| 448 | |
Tim Holy | 44b19b4 | 2012-02-08 22:11:12 +0100 | [diff] [blame] | 449 | \li \b Next: \ref TutorialMapClass |
Gael Guennebaud | 86a1926 | 2009-01-16 17:20:12 +0000 | [diff] [blame] | 450 | |
| 451 | */ |
| 452 | |
| 453 | } |