blob: 9a212ae02c8f268d8621831d536765faaa15486b [file] [log] [blame]
Gael Guennebaud86a19262009-01-16 17:20:12 +00001namespace Eigen {
2
Gael Guennebaud41ea92d2010-07-04 10:14:47 +02003/** \page TutorialSparse Tutorial page 9 - Sparse Matrix
Gael Guennebaud86a19262009-01-16 17:20:12 +00004 \ingroup Tutorial
5
Gael Guennebaud41ea92d2010-07-04 10:14:47 +02006\li \b Previous: \ref TutorialGeometry
Tim Holy44b19b42012-02-08 22:11:12 +01007\li \b Next: \ref TutorialMapClass
Gael Guennebaud86a19262009-01-16 17:20:12 +00008
Gael Guennebaud2ea1e492012-12-28 18:58:07 +01009\tableofcontents
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020010
Gael Guennebaud86a19262009-01-16 17:20:12 +000011<hr>
12
Gael Guennebaud70206ab2011-11-24 17:32:30 +010013Manipulating and solving sparse problems involves various modules which are summarized below:
Gael Guennebaud86a19262009-01-16 17:20:12 +000014
Gael Guennebaudf66fe262010-10-19 11:40:49 +020015<table class="manual">
Gael Guennebaud70206ab2011-11-24 17:32:30 +010016<tr><th>Module</th><th>Header file</th><th>Contents</th></tr>
Gael Guennebaudf41d96d2012-12-24 13:33:22 +010017<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>
Gael Guennebaud70206ab2011-11-24 17:32:30 +010018<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>
19<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>
Gael Guennebaudf41d96d2012-12-24 13:33:22 +010020<tr><td>\link Sparse_modules Sparse \endlink</td><td>\code#include <Eigen/Sparse>\endcode</td><td>Includes all the above modules</td></tr>
Gael Guennebaud86a19262009-01-16 17:20:12 +000021</table>
22
Gael Guennebaud3f56de22011-12-03 10:26:00 +010023\section TutorialSparseIntro Sparse matrix representation
Gael Guennebaud70206ab2011-11-24 17:32:30 +010024
25In 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.
26
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020027\b The \b %SparseMatrix \b class
28
Gael Guennebaud70206ab2011-11-24 17:32:30 +010029The class SparseMatrix is the main sparse matrix representation of Eigen's sparse module; it offers high performance and low memory usage.
30It implements a more versatile variant of the widely-used Compressed Column (or Row) Storage scheme.
31It 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.
Jitse Niesen8c71d732012-06-20 09:52:45 +010034 - \c OuterStarts: stores for each column (resp. row) the index of the first non-zero in the previous two arrays.
Gael Guennebaud3f56de22011-12-03 10:26:00 +010035 - \c InnerNNZs: stores the number of non-zeros of each column (resp. row).
Gael Guennebaud70206ab2011-11-24 17:32:30 +010036The 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.
37The word \c outer refers to the other direction.
38
39This storage scheme is better explained on an example. The following matrix
40<table class="manual">
41<tr><td> 0</td><td>3</td><td> 0</td><td>0</td><td> 0</td></tr>
42<tr><td>22</td><td>0</td><td> 0</td><td>0</td><td>17</td></tr>
43<tr><td> 7</td><td>5</td><td> 0</td><td>1</td><td> 0</td></tr>
44<tr><td> 0</td><td>0</td><td> 0</td><td>0</td><td> 0</td></tr>
45<tr><td> 0</td><td>0</td><td>14</td><td>0</td><td> 8</td></tr>
46</table>
47
Gael Guennebaud3f56de22011-12-03 10:26:00 +010048and one of its possible sparse, \b column \b major representation:
Gael Guennebaud70206ab2011-11-24 17:32:30 +010049<table class="manual">
50<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>
51<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>
52</table>
53<table class="manual">
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020054<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 Guennebaud3f56de22011-12-03 10:26:00 +010055<tr><td>InnerNNZs:</td> <td>2</td><td>2</td><td>1</td><td>1</td><td> 2</td><td></td></tr>
Gael Guennebaud70206ab2011-11-24 17:32:30 +010056</table>
57
Gael Guennebaud3f56de22011-12-03 10:26:00 +010058Currently the elements of a given inner vector are guaranteed to be always sorted by increasing inner indices.
59The \c "_" indicates available free space to quickly insert new elements.
60Assuming 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.
61On 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 Guennebaud70206ab2011-11-24 17:32:30 +010062
Gael Guennebaud3f56de22011-12-03 10:26:00 +010063The case where no empty space is available is a special case, and is refered as the \em compressed mode.
64It corresponds to the widely used Compressed Column (or Row) Storage schemes (CCS or CRS).
65Any SparseMatrix can be turned to this form by calling the SparseMatrix::makeCompressed() function.
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020066In 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 Guennebaud3f56de22011-12-03 10:26:00 +010067Therefore, in practice a call to SparseMatrix::makeCompressed() frees this buffer.
68
69It is worth noting that most of our wrappers to external libraries requires compressed matrices as inputs.
70
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020071The results of %Eigen's operations always produces \b compressed sparse matrices.
Gael Guennebaud3f56de22011-12-03 10:26:00 +010072On the other hand, the insertion of a new element into a SparseMatrix converts this later to the \b uncompressed mode.
Gael Guennebaud70206ab2011-11-24 17:32:30 +010073
74Here is the previous matrix represented in compressed mode:
Gael Guennebaudf66fe262010-10-19 11:40:49 +020075<table class="manual">
Gael Guennebaud86a19262009-01-16 17:20:12 +000076<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 Guennebaud70206ab2011-11-24 17:32:30 +010077<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 Guennebaud86a19262009-01-16 17:20:12 +000078</table>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020079<table class="manual">
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020080<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 Guennebaud25f16582009-01-21 17:44:58 +000081</table>
Gael Guennebaud70206ab2011-11-24 17:32:30 +010082
Gael Guennebaud3f56de22011-12-03 10:26:00 +010083A SparseVector is a special case of a SparseMatrix where only the \c Values and \c InnerIndices arrays are stored.
84There is no notion of compressed/uncompressed mode for a SparseVector.
Gael Guennebaud25f16582009-01-21 17:44:58 +000085
Gael Guennebaud86a19262009-01-16 17:20:12 +000086
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020087\section TutorialSparseExample First example
Gael Guennebaud86a19262009-01-16 17:20:12 +000088
Jitse Niesen8c71d732012-06-20 09:52:45 +010089Before 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.
90Such 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.
Gael Guennebaud25f16582009-01-21 17:44:58 +000091
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020092<table class="manual">
93<tr><td>
94\include Tutorial_sparse_example.cpp
95</td>
96<td>
97\image html Tutorial_sparse_example.jpeg
98</td></tr></table>
99
100In 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.
101
102In 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.
103The raw and flat list of non-zero entries is then converted to a true SparseMatrix object \c A.
104Note that the elements of the list do not have to be sorted, and possible duplicate entries will be summed up.
105
Jitse Niesen8c71d732012-06-20 09:52:45 +0100106The last step consists of effectively solving the assembled problem.
107Since 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.
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200108
109The 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.
110
Jitse Niesen8c71d732012-06-20 09:52:45 +0100111Describing 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.
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200112
113
114
115
116\section TutorialSparseSparseMatrix The SparseMatrix class
117
118\b %Matrix \b and \b vector \b properties \n
119
120The SparseMatrix and SparseVector classes take three template arguments:
121 * the scalar type (e.g., double)
122 * the storage order (ColMajor or RowMajor, the default is RowMajor)
123 * the inner index type (default is \c int).
124
125As for dense Matrix objects, constructors takes the size of the object.
126Here are some examples:
127
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100128\code
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200129SparseMatrix<std::complex<float> > mat(1000,2000); // declares a 1000x2000 column-major compressed sparse matrix of complex<float>
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100130SparseMatrix<double,RowMajor> mat(1000,2000); // declares a 1000x2000 row-major compressed sparse matrix of double
131SparseVector<std::complex<float> > vec(1000); // declares a column sparse vector of complex<float> of size 1000
132SparseVector<double,RowMajor> vec(1000); // declares a row sparse vector of double of size 1000
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100133\endcode
134
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200135In the rest of the tutorial, \c mat and \c vec represent any sparse-matrix and sparse-vector objects, respectively.
136
137The dimensions of a matrix can be queried using the following functions:
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200138<table class="manual">
Gael Guennebaud86a19262009-01-16 17:20:12 +0000139<tr><td>Standard \n dimensions</td><td>\code
140mat.rows()
141mat.cols()\endcode</td>
142<td>\code
143vec.size() \endcode</td>
144</tr>
145<tr><td>Sizes along the \n inner/outer dimensions</td><td>\code
146mat.innerSize()
147mat.outerSize()\endcode</td>
148<td></td>
149</tr>
Tim Holy16a2d892011-06-20 22:47:58 -0500150<tr><td>Number of non \n zero coefficients</td><td>\code
Gael Guennebaud86a19262009-01-16 17:20:12 +0000151mat.nonZeros() \endcode</td>
152<td>\code
153vec.nonZeros() \endcode</td></tr>
154</table>
155
Gael Guennebaud25f16582009-01-21 17:44:58 +0000156
Gael Guennebaud86a19262009-01-16 17:20:12 +0000157\b Iterating \b over \b the \b nonzero \b coefficients \n
158
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200159Random access to the elements of a sparse object can be done through the \c coeffRef(i,j) function.
160However, this function involves a quite expensive binary search.
161In 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.
162Here is an example:
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200163<table class="manual">
Gael Guennebaud86a19262009-01-16 17:20:12 +0000164<tr><td>
165\code
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200166SparseMatrix<double> mat(rows,cols);
Jitse Niesen59b83c12011-09-10 09:18:18 +0100167for (int k=0; k<mat.outerSize(); ++k)
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200168 for (SparseMatrix<double>::InnerIterator it(mat,k); it; ++it)
Gael Guennebaud86a19262009-01-16 17:20:12 +0000169 {
170 it.value();
171 it.row(); // row index
172 it.col(); // col index (here it is equal to k)
173 it.index(); // inner index, here it is equal to it.row()
174 }
175\endcode
176</td><td>
177\code
Gael Guennebaud25f16582009-01-21 17:44:58 +0000178SparseVector<double> vec(size);
179for (SparseVector<double>::InnerIterator it(vec); it; ++it)
Gael Guennebaud86a19262009-01-16 17:20:12 +0000180{
Gael Guennebaud25f16582009-01-21 17:44:58 +0000181 it.value(); // == vec[ it.index() ]
Gael Guennebaudc532f422009-09-25 16:30:44 +0200182 it.index();
Gael Guennebaud86a19262009-01-16 17:20:12 +0000183}
184\endcode
185</td></tr>
186</table>
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200187For a writable expression, the referenced value can be modified using the valueRef() function.
Jitse Niesen59b83c12011-09-10 09:18:18 +0100188If the type of the sparse matrix or vector depends on a template parameter, then the \c typename keyword is
189required to indicate that \c InnerIterator denotes a type; see \ref TopicTemplateKeyword for details.
190
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100191
Gael Guennebaud86a19262009-01-16 17:20:12 +0000192\section TutorialSparseFilling Filling a sparse matrix
193
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100194Because of the special storage scheme of a SparseMatrix, special care has to be taken when adding new nonzero entries.
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200195For 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 Guennebaudc532f422009-09-25 16:30:44 +0200196
Jitse Niesen8c71d732012-06-20 09:52:45 +0100197The 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.
Gael Guennebaud1763f862012-02-04 10:44:07 +0100198
199Here is a typical usage example:
200\code
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200201typedef Eigen::Triplet<double> T;
Gael Guennebaud1763f862012-02-04 10:44:07 +0100202std::vector<T> tripletList;
Gael Guennebaud1edb3962012-09-26 19:24:41 +0200203tripletList.reserve(estimation_of_entries);
Gael Guennebaud1763f862012-02-04 10:44:07 +0100204for(...)
205{
206 // ...
207 tripletList.push_back(T(i,j,v_ij));
208}
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200209SparseMatrixType mat(rows,cols);
210mat.setFromTriplets(tripletList.begin(), tripletList.end());
211// mat is ready to go!
Gael Guennebaud1763f862012-02-04 10:44:07 +0100212\endcode
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200213The \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 Guennebaud1763f862012-02-04 10:44:07 +0100214See the SparseMatrix::setFromTriplets() function and class Triplet for more details.
215
216
Jitse Niesen8c71d732012-06-20 09:52:45 +0100217In some cases, however, slightly higher performance, and lower memory consumption can be reached by directly inserting the non-zeros into the destination matrix.
Gael Guennebaud1763f862012-02-04 10:44:07 +0100218A typical scenario of this approach is illustrated bellow:
Gael Guennebaud25f16582009-01-21 17:44:58 +0000219\code
Gael Guennebaudc861e052011-12-03 18:14:51 +01002201: SparseMatrix<double> mat(rows,cols); // default is column major
Gael Guennebaud1aa6c7f2011-12-11 17:18:14 +01002212: mat.reserve(VectorXi::Constant(cols,6));
Gael Guennebaud70206ab2011-11-24 17:32:30 +01002223: for each i,j such that v_ij != 0
Gael Guennebaud3f56de22011-12-03 10:26:00 +01002234: mat.insert(i,j) = v_ij; // alternative: mat.coeffRef(i,j) += v_ij;
2245: mat.makeCompressed(); // optional
Gael Guennebaud25f16582009-01-21 17:44:58 +0000225\endcode
226
Jitse Niesen8c71d732012-06-20 09:52:45 +0100227- 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.
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100228- 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.
229- 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.
230- The line 5 suppresses the remaining empty space and transforms the matrix into a compressed column storage.
Gael Guennebaud25f16582009-01-21 17:44:58 +0000231
Gael Guennebaud25f16582009-01-21 17:44:58 +0000232
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100233\section TutorialSparseDirectSolvers Solving linear problems
Gael Guennebaud86a19262009-01-16 17:20:12 +0000234
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200235%Eigen currently provides a limited set of built-in solvers, as well as wrappers to external solver libraries.
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100236They are summarized in the following table:
Gael Guennebaudab41c182010-08-18 15:33:58 +0200237
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100238<table class="manual">
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200239<tr><th>Class</th><th>Module</th><th>Solver kind</th><th>Matrix kind</th><th>Features related to performance</th>
240 <th>Dependencies,License</th><th class="width20em"><p>Notes</p></th></tr>
241<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 Guennebaud70206ab2011-11-24 17:32:30 +0100242 <td>built-in, LGPL</td>
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200243 <td>SimplicialLDLT is often preferable</td></tr>
244<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 Guennebaud70206ab2011-11-24 17:32:30 +0100245 <td>built-in, LGPL</td>
246 <td>Recommended for very sparse and not too large problems (e.g., 2D Poisson eq.)</td></tr>
247<tr><td>ConjugateGradient</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>Classic iterative CG</td><td>SPD</td><td>Preconditionning</td>
248 <td>built-in, LGPL</td>
249 <td>Recommended for large symmetric problems (e.g., 3D Poisson eq.)</td></tr>
250<tr><td>BiCGSTAB</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>Iterative stabilized bi-conjugate gradient</td><td>Square</td><td>Preconditionning</td>
251 <td>built-in, LGPL</td>
252 <td>Might not always converge</td></tr>
253
254
Gael Guennebaud6f3057f2012-06-21 10:51:22 +0200255<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>
256 <td>Requires the <a href="http://pastix.gforge.inria.fr">PaStiX</a> package, \b CeCILL-C </td>
257 <td>optimized for tough problems and symmetric patterns</td></tr>
258<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 Guennebaud70206ab2011-11-24 17:32:30 +0100259 <td>Requires the <a href="http://www.cise.ufl.edu/research/sparse/SuiteSparse/">SuiteSparse</a> package, \b GPL </td>
260 <td></td></tr>
261<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>
262 <td>Requires the <a href="http://www.cise.ufl.edu/research/sparse/SuiteSparse/">SuiteSparse</a> package, \b GPL </td>
263 <td></td></tr>
264<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 Guennebaud3f56de22011-12-03 10:26:00 +0100265 <td>Requires the <a href="http://crd-legacy.lbl.gov/~xiaoye/SuperLU/">SuperLU</a> library, (BSD-like)</td>
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100266 <td></td></tr>
267</table>
268
269Here \c SPD means symmetric positive definite.
270
271All these solvers follow the same general concept.
272Here is a typical and general example:
Gael Guennebaudab41c182010-08-18 15:33:58 +0200273\code
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100274#include <Eigen/RequiredModuleName>
Gael Guennebaudab41c182010-08-18 15:33:58 +0200275// ...
276SparseMatrix<double> A;
277// fill A
278VectorXd b, x;
279// fill b
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100280// solve Ax = b
281SolverClassName<SparseMatrix<double> > solver;
282solver.compute(A);
Gael Guennebaudf23dd7c2012-07-28 13:07:45 +0200283if(solver.info()!=Success) {
Tim Holy16a2d892011-06-20 22:47:58 -0500284 // decomposition failed
Gael Guennebaudab41c182010-08-18 15:33:58 +0200285 return;
286}
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100287x = solver.solve(b);
Gael Guennebaudf23dd7c2012-07-28 13:07:45 +0200288if(solver.info()!=Success) {
Gael Guennebaudab41c182010-08-18 15:33:58 +0200289 // solving failed
290 return;
291}
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100292// solve for another right hand side:
293x1 = solver.solve(b1);
Gael Guennebaudab41c182010-08-18 15:33:58 +0200294\endcode
295
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100296For \c SPD solvers, a second optional template argument allows to specify which triangular part have to be used, e.g.:
297
298\code
299#include <Eigen/IterativeLinearSolvers>
300
301ConjugateGradient<SparseMatrix<double>, Eigen::Upper> solver;
302x = solver.compute(A).solve(b);
303\endcode
304In 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.
305
306In the case where multiple problems with the same sparcity pattern have to be solved, then the "compute" step can be decomposed as follow:
307\code
308SolverClassName<SparseMatrix<double> > solver;
309solver.analyzePattern(A); // for this step the numerical values of A are not used
310solver.factorize(A);
311x1 = solver.solve(b1);
312x2 = solver.solve(b2);
313...
314A = ...; // modify the values of the nonzeros of A, the nonzeros pattern must stay unchanged
315solver.factorize(A);
316x1 = solver.solve(b1);
317x2 = solver.solve(b2);
318...
319\endcode
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100320The compute() method is equivalent to calling both analyzePattern() and factorize().
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100321
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100322Finally, each solver provides some specific features, such as determinant, access to the factors, controls of the iterations, and so on.
323More details are availble in the documentations of the respective classes.
Gael Guennebaud86a19262009-01-16 17:20:12 +0000324
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200325
326\section TutorialSparseFeatureSet Supported operators and functions
327
328Because of their special storage format, sparse matrices cannot offer the same level of flexbility than dense matrices.
329In Eigen's sparse module we chose to expose only the subset of the dense matrix API which can be efficiently implemented.
330In the following \em sm denotes a sparse matrix, \em sv a sparse vector, \em dm a dense matrix, and \em dv a dense vector.
331
332\subsection TutorialSparse_BasicOps Basic operations
333
334%Sparse expressions support most of the unary and binary coefficient wise operations:
335\code
336sm1.real() sm1.imag() -sm1 0.5*sm1
337sm1+sm2 sm1-sm2 sm1.cwiseProduct(sm2)
338\endcode
339However, a strong restriction is that the storage orders must match. For instance, in the following example:
340\code
341sm4 = sm1 + sm2 + sm3;
342\endcode
343sm1, sm2, and sm3 must all be row-major or all column major.
344On the other hand, there is no restriction on the target matrix sm4.
345For 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:
346\code
347SparseMatrix<double> A, B;
348B = SparseMatrix<double>(A.transpose()) + A;
349\endcode
350
351Binary coefficient wise operators can also mix sparse and dense expressions:
352\code
353sm2 = sm1.cwiseProduct(dm1);
354dm2 = sm1 + dm1;
355\endcode
356
357
358%Sparse expressions also support transposition:
359\code
360sm1 = sm2.transpose();
361sm1 = sm2.adjoint();
362\endcode
363However, there is no transposeInPlace() method.
364
365
366\subsection TutorialSparse_Products Matrix products
367
368%Eigen supports various kind of sparse matrix products which are summarize below:
369 - \b sparse-dense:
370 \code
371dv2 = sm1 * dv1;
372dm2 = dm1 * sm1.adjoint();
373dm2 = 2. * sm1 * dm1;
374 \endcode
375 - \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():
376 \code
377dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored
378dm2 = A.selfadjointView<Upper>() * dm1; // if only the upper part of A is stored
379dm2 = A.selfadjointView<Lower>() * dm1; // if only the lower part of A is stored
380 \endcode
381 - \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:
382 \code
383sm3 = sm1 * sm2;
384sm3 = 4 * sm1.adjoint() * sm2;
385 \endcode
Jitse Niesen8c71d732012-06-20 09:52:45 +0100386 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:
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200387 \code
388sm3 = (sm1 * sm2).prune(); // removes numerical zeros
389sm3 = (sm1 * sm2).prune(ref); // removes elements much smaller than ref
Jitse Niesen8c71d732012-06-20 09:52:45 +0100390sm3 = (sm1 * sm2).prune(ref,epsilon); // removes elements smaller than ref*epsilon
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200391 \endcode
392
393 - \b permutations. Finally, permutations can be applied to sparse matrices too:
394 \code
395PermutationMatrix<Dynamic,Dynamic> P = ...;
396sm2 = P * sm1;
397sm2 = sm1 * P.inverse();
398sm2 = sm1.transpose() * P;
399 \endcode
400
401
402\subsection TutorialSparse_TriangularSelfadjoint Triangular and selfadjoint views
403
Jitse Niesen8c71d732012-06-20 09:52:45 +0100404Just 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:
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200405\code
406dm2 = sm1.triangularView<Lower>(dm1);
407dv2 = sm1.transpose().triangularView<Upper>(dv1);
408\endcode
409
410The selfadjointView() function permits various operations:
411 - optimized sparse-dense matrix products:
412 \code
413dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored
414dm2 = A.selfadjointView<Upper>() * dm1; // if only the upper part of A is stored
415dm2 = A.selfadjointView<Lower>() * dm1; // if only the lower part of A is stored
416 \endcode
417 - copy of triangular parts:
418 \code
419sm2 = sm1.selfadjointView<Upper>(); // makes a full selfadjoint matrix from the upper triangular part
Jitse Niesen8c71d732012-06-20 09:52:45 +0100420sm2.selfadjointView<Lower>() = sm1.selfadjointView<Upper>(); // copies the upper triangular part to the lower triangular part
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200421 \endcode
422 - application of symmetric permutations:
423 \code
424PermutationMatrix<Dynamic,Dynamic> P = ...;
425sm2 = A.selfadjointView<Upper>().twistedBy(P); // compute P S P' from the upper triangular part of A, and make it a full matrix
426sm2.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
427 \endcode
428
429\subsection TutorialSparse_Submat Sub-matrices
430
431%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:
432\code
433 sm1.innerVector(j); // returns an expression of the j-th column (resp. row) of the matrix if sm1 is col-major (resp. row-major)
434 sm1.innerVectors(j, nb); // returns an expression of the nb columns (resp. row) starting from the j-th column (resp. row)
435 // of the matrix if sm1 is col-major (resp. row-major)
436 sm1.middleRows(j, nb); // for row major matrices only, get a range of nb rows
437 sm1.middleCols(j, nb); // for column major matrices only, get a range of nb columns
438\endcode
439
Tim Holy44b19b42012-02-08 22:11:12 +0100440\li \b Next: \ref TutorialMapClass
Gael Guennebaud86a19262009-01-16 17:20:12 +0000441
442*/
443
444}