blob: 792b72c459f265ddb1197f24f23a559a7d2a69b2 [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
9\b Table \b of \b contents \n
10 - \ref TutorialSparseIntro
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020011 - \ref TutorialSparseExample "Example"
12 - \ref TutorialSparseSparseMatrix
Gael Guennebaud86a19262009-01-16 17:20:12 +000013 - \ref TutorialSparseFilling
Gael Guennebaud86a19262009-01-16 17:20:12 +000014 - \ref TutorialSparseDirectSolvers
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020015 - \ref TutorialSparseFeatureSet
16 - \ref TutorialSparse_BasicOps
17 - \ref TutorialSparse_Products
18 - \ref TutorialSparse_TriangularSelfadjoint
19 - \ref TutorialSparse_Submat
20
21
Gael Guennebaud86a19262009-01-16 17:20:12 +000022<hr>
23
Gael Guennebaud70206ab2011-11-24 17:32:30 +010024Manipulating and solving sparse problems involves various modules which are summarized below:
Gael Guennebaud86a19262009-01-16 17:20:12 +000025
Gael Guennebaudf66fe262010-10-19 11:40:49 +020026<table class="manual">
Gael Guennebaud70206ab2011-11-24 17:32:30 +010027<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 Guennebaud86a19262009-01-16 17:20:12 +000032</table>
33
Gael Guennebaud3f56de22011-12-03 10:26:00 +010034\section TutorialSparseIntro Sparse matrix representation
Gael Guennebaud70206ab2011-11-24 17:32:30 +010035
36In 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 Guennebaud52dce0c2012-06-20 09:28:32 +020038\b The \b %SparseMatrix \b class
39
Gael Guennebaud70206ab2011-11-24 17:32:30 +010040The class SparseMatrix is the main sparse matrix representation of Eigen's sparse module; it offers high performance and low memory usage.
41It implements a more versatile variant of the widely-used Compressed Column (or Row) Storage scheme.
42It 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.
Jitse Niesen8c71d732012-06-20 09:52:45 +010045 - \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 +010046 - \c InnerNNZs: stores the number of non-zeros of each column (resp. row).
Gael Guennebaud70206ab2011-11-24 17:32:30 +010047The 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.
48The word \c outer refers to the other direction.
49
50This 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 Guennebaud3f56de22011-12-03 10:26:00 +010059and one of its possible sparse, \b column \b major representation:
Gael Guennebaud70206ab2011-11-24 17:32:30 +010060<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 Guennebaud52dce0c2012-06-20 09:28:32 +020065<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 +010066<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 +010067</table>
68
Gael Guennebaud3f56de22011-12-03 10:26:00 +010069Currently the elements of a given inner vector are guaranteed to be always sorted by increasing inner indices.
70The \c "_" indicates available free space to quickly insert new elements.
71Assuming 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.
72On 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 +010073
Gael Guennebaud3f56de22011-12-03 10:26:00 +010074The case where no empty space is available is a special case, and is refered as the \em compressed mode.
75It corresponds to the widely used Compressed Column (or Row) Storage schemes (CCS or CRS).
76Any SparseMatrix can be turned to this form by calling the SparseMatrix::makeCompressed() function.
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020077In 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 +010078Therefore, in practice a call to SparseMatrix::makeCompressed() frees this buffer.
79
80It is worth noting that most of our wrappers to external libraries requires compressed matrices as inputs.
81
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020082The results of %Eigen's operations always produces \b compressed sparse matrices.
Gael Guennebaud3f56de22011-12-03 10:26:00 +010083On 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 +010084
85Here is the previous matrix represented in compressed mode:
Gael Guennebaudf66fe262010-10-19 11:40:49 +020086<table class="manual">
Gael Guennebaud86a19262009-01-16 17:20:12 +000087<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 +010088<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 +000089</table>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020090<table class="manual">
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020091<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 +000092</table>
Gael Guennebaud70206ab2011-11-24 17:32:30 +010093
Gael Guennebaud3f56de22011-12-03 10:26:00 +010094A SparseVector is a special case of a SparseMatrix where only the \c Values and \c InnerIndices arrays are stored.
95There is no notion of compressed/uncompressed mode for a SparseVector.
Gael Guennebaud25f16582009-01-21 17:44:58 +000096
Gael Guennebaud86a19262009-01-16 17:20:12 +000097
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020098\section TutorialSparseExample First example
Gael Guennebaud86a19262009-01-16 17:20:12 +000099
Jitse Niesen8c71d732012-06-20 09:52:45 +0100100Before 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.
101Such 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 +0000102
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200103<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
111In 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
113In 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.
114The raw and flat list of non-zero entries is then converted to a true SparseMatrix object \c A.
115Note that the elements of the list do not have to be sorted, and possible duplicate entries will be summed up.
116
Jitse Niesen8c71d732012-06-20 09:52:45 +0100117The last step consists of effectively solving the assembled problem.
118Since 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 +0200119
120The 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
Jitse Niesen8c71d732012-06-20 09:52:45 +0100122Describing 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 +0200123
124
125
126
127\section TutorialSparseSparseMatrix The SparseMatrix class
128
129\b %Matrix \b and \b vector \b properties \n
130
131The 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
136As for dense Matrix objects, constructors takes the size of the object.
137Here are some examples:
138
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100139\code
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200140SparseMatrix<std::complex<float> > mat(1000,2000); // declares a 1000x2000 column-major compressed sparse matrix of complex<float>
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100141SparseMatrix<double,RowMajor> mat(1000,2000); // declares a 1000x2000 row-major compressed sparse matrix of double
142SparseVector<std::complex<float> > vec(1000); // declares a column sparse vector of complex<float> of size 1000
143SparseVector<double,RowMajor> vec(1000); // declares a row sparse vector of double of size 1000
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100144\endcode
145
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200146In the rest of the tutorial, \c mat and \c vec represent any sparse-matrix and sparse-vector objects, respectively.
147
148The dimensions of a matrix can be queried using the following functions:
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200149<table class="manual">
Gael Guennebaud86a19262009-01-16 17:20:12 +0000150<tr><td>Standard \n dimensions</td><td>\code
151mat.rows()
152mat.cols()\endcode</td>
153<td>\code
154vec.size() \endcode</td>
155</tr>
156<tr><td>Sizes along the \n inner/outer dimensions</td><td>\code
157mat.innerSize()
158mat.outerSize()\endcode</td>
159<td></td>
160</tr>
Tim Holy16a2d892011-06-20 22:47:58 -0500161<tr><td>Number of non \n zero coefficients</td><td>\code
Gael Guennebaud86a19262009-01-16 17:20:12 +0000162mat.nonZeros() \endcode</td>
163<td>\code
164vec.nonZeros() \endcode</td></tr>
165</table>
166
Gael Guennebaud25f16582009-01-21 17:44:58 +0000167
Gael Guennebaud86a19262009-01-16 17:20:12 +0000168\b Iterating \b over \b the \b nonzero \b coefficients \n
169
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200170Random access to the elements of a sparse object can be done through the \c coeffRef(i,j) function.
171However, this function involves a quite expensive binary search.
172In 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.
173Here is an example:
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200174<table class="manual">
Gael Guennebaud86a19262009-01-16 17:20:12 +0000175<tr><td>
176\code
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200177SparseMatrix<double> mat(rows,cols);
Jitse Niesen59b83c12011-09-10 09:18:18 +0100178for (int k=0; k<mat.outerSize(); ++k)
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200179 for (SparseMatrix<double>::InnerIterator it(mat,k); it; ++it)
Gael Guennebaud86a19262009-01-16 17:20:12 +0000180 {
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 Guennebaud25f16582009-01-21 17:44:58 +0000189SparseVector<double> vec(size);
190for (SparseVector<double>::InnerIterator it(vec); it; ++it)
Gael Guennebaud86a19262009-01-16 17:20:12 +0000191{
Gael Guennebaud25f16582009-01-21 17:44:58 +0000192 it.value(); // == vec[ it.index() ]
Gael Guennebaudc532f422009-09-25 16:30:44 +0200193 it.index();
Gael Guennebaud86a19262009-01-16 17:20:12 +0000194}
195\endcode
196</td></tr>
197</table>
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200198For a writable expression, the referenced value can be modified using the valueRef() function.
Jitse Niesen59b83c12011-09-10 09:18:18 +0100199If the type of the sparse matrix or vector depends on a template parameter, then the \c typename keyword is
200required to indicate that \c InnerIterator denotes a type; see \ref TopicTemplateKeyword for details.
201
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100202
Gael Guennebaud86a19262009-01-16 17:20:12 +0000203\section TutorialSparseFilling Filling a sparse matrix
204
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100205Because 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 +0200206For 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 +0200207
Jitse Niesen8c71d732012-06-20 09:52:45 +0100208The 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 +0100209
210Here is a typical usage example:
211\code
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200212typedef Eigen::Triplet<double> T;
Gael Guennebaud1763f862012-02-04 10:44:07 +0100213std::vector<T> tripletList;
214triplets.reserve(estimation_of_entries);
215for(...)
216{
217 // ...
218 tripletList.push_back(T(i,j,v_ij));
219}
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200220SparseMatrixType mat(rows,cols);
221mat.setFromTriplets(tripletList.begin(), tripletList.end());
222// mat is ready to go!
Gael Guennebaud1763f862012-02-04 10:44:07 +0100223\endcode
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200224The \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 +0100225See the SparseMatrix::setFromTriplets() function and class Triplet for more details.
226
227
Jitse Niesen8c71d732012-06-20 09:52:45 +0100228In 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 +0100229A typical scenario of this approach is illustrated bellow:
Gael Guennebaud25f16582009-01-21 17:44:58 +0000230\code
Gael Guennebaudc861e052011-12-03 18:14:51 +01002311: SparseMatrix<double> mat(rows,cols); // default is column major
Gael Guennebaud1aa6c7f2011-12-11 17:18:14 +01002322: mat.reserve(VectorXi::Constant(cols,6));
Gael Guennebaud70206ab2011-11-24 17:32:30 +01002333: for each i,j such that v_ij != 0
Gael Guennebaud3f56de22011-12-03 10:26:00 +01002344: mat.insert(i,j) = v_ij; // alternative: mat.coeffRef(i,j) += v_ij;
2355: mat.makeCompressed(); // optional
Gael Guennebaud25f16582009-01-21 17:44:58 +0000236\endcode
237
Jitse Niesen8c71d732012-06-20 09:52:45 +0100238- 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 +0100239- 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 Guennebaud25f16582009-01-21 17:44:58 +0000242
Gael Guennebaud25f16582009-01-21 17:44:58 +0000243
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100244\section TutorialSparseDirectSolvers Solving linear problems
Gael Guennebaud86a19262009-01-16 17:20:12 +0000245
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200246%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 +0100247They are summarized in the following table:
Gael Guennebaudab41c182010-08-18 15:33:58 +0200248
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100249<table class="manual">
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200250<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 Guennebaud70206ab2011-11-24 17:32:30 +0100253 <td>built-in, LGPL</td>
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200254 <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 Guennebaud70206ab2011-11-24 17:32:30 +0100256 <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 Guennebaud52dce0c2012-06-20 09:28:32 +0200267<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 +0100268 <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 Guennebaud3f56de22011-12-03 10:26:00 +0100274 <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 +0100275 <td></td></tr>
276</table>
277
278Here \c SPD means symmetric positive definite.
279
280All these solvers follow the same general concept.
281Here is a typical and general example:
Gael Guennebaudab41c182010-08-18 15:33:58 +0200282\code
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100283#include <Eigen/RequiredModuleName>
Gael Guennebaudab41c182010-08-18 15:33:58 +0200284// ...
285SparseMatrix<double> A;
286// fill A
287VectorXd b, x;
288// fill b
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100289// solve Ax = b
290SolverClassName<SparseMatrix<double> > solver;
291solver.compute(A);
292if(solver.info()!=Succeeded) {
Tim Holy16a2d892011-06-20 22:47:58 -0500293 // decomposition failed
Gael Guennebaudab41c182010-08-18 15:33:58 +0200294 return;
295}
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100296x = solver.solve(b);
297if(solver.info()!=Succeeded) {
Gael Guennebaudab41c182010-08-18 15:33:58 +0200298 // solving failed
299 return;
300}
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100301// solve for another right hand side:
302x1 = solver.solve(b1);
Gael Guennebaudab41c182010-08-18 15:33:58 +0200303\endcode
304
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100305For \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
310ConjugateGradient<SparseMatrix<double>, Eigen::Upper> solver;
311x = solver.compute(A).solve(b);
312\endcode
313In 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
315In 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
317SolverClassName<SparseMatrix<double> > solver;
318solver.analyzePattern(A); // for this step the numerical values of A are not used
319solver.factorize(A);
320x1 = solver.solve(b1);
321x2 = solver.solve(b2);
322...
323A = ...; // modify the values of the nonzeros of A, the nonzeros pattern must stay unchanged
324solver.factorize(A);
325x1 = solver.solve(b1);
326x2 = solver.solve(b2);
327...
328\endcode
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100329The compute() method is equivalent to calling both analyzePattern() and factorize().
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100330
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100331Finally, each solver provides some specific features, such as determinant, access to the factors, controls of the iterations, and so on.
332More details are availble in the documentations of the respective classes.
Gael Guennebaud86a19262009-01-16 17:20:12 +0000333
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200334
335\section TutorialSparseFeatureSet Supported operators and functions
336
337Because of their special storage format, sparse matrices cannot offer the same level of flexbility than dense matrices.
338In Eigen's sparse module we chose to expose only the subset of the dense matrix API which can be efficiently implemented.
339In 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
345sm1.real() sm1.imag() -sm1 0.5*sm1
346sm1+sm2 sm1-sm2 sm1.cwiseProduct(sm2)
347\endcode
348However, a strong restriction is that the storage orders must match. For instance, in the following example:
349\code
350sm4 = sm1 + sm2 + sm3;
351\endcode
352sm1, sm2, and sm3 must all be row-major or all column major.
353On the other hand, there is no restriction on the target matrix sm4.
354For 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
356SparseMatrix<double> A, B;
357B = SparseMatrix<double>(A.transpose()) + A;
358\endcode
359
360Binary coefficient wise operators can also mix sparse and dense expressions:
361\code
362sm2 = sm1.cwiseProduct(dm1);
363dm2 = sm1 + dm1;
364\endcode
365
366
367%Sparse expressions also support transposition:
368\code
369sm1 = sm2.transpose();
370sm1 = sm2.adjoint();
371\endcode
372However, 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
380dv2 = sm1 * dv1;
381dm2 = dm1 * sm1.adjoint();
382dm2 = 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
386dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored
387dm2 = A.selfadjointView<Upper>() * dm1; // if only the upper part of A is stored
388dm2 = 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
392sm3 = sm1 * sm2;
393sm3 = 4 * sm1.adjoint() * sm2;
394 \endcode
Jitse Niesen8c71d732012-06-20 09:52:45 +0100395 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 +0200396 \code
397sm3 = (sm1 * sm2).prune(); // removes numerical zeros
398sm3 = (sm1 * sm2).prune(ref); // removes elements much smaller than ref
Jitse Niesen8c71d732012-06-20 09:52:45 +0100399sm3 = (sm1 * sm2).prune(ref,epsilon); // removes elements smaller than ref*epsilon
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200400 \endcode
401
402 - \b permutations. Finally, permutations can be applied to sparse matrices too:
403 \code
404PermutationMatrix<Dynamic,Dynamic> P = ...;
405sm2 = P * sm1;
406sm2 = sm1 * P.inverse();
407sm2 = sm1.transpose() * P;
408 \endcode
409
410
411\subsection TutorialSparse_TriangularSelfadjoint Triangular and selfadjoint views
412
Jitse Niesen8c71d732012-06-20 09:52:45 +0100413Just 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 +0200414\code
415dm2 = sm1.triangularView<Lower>(dm1);
416dv2 = sm1.transpose().triangularView<Upper>(dv1);
417\endcode
418
419The selfadjointView() function permits various operations:
420 - optimized sparse-dense matrix products:
421 \code
422dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored
423dm2 = A.selfadjointView<Upper>() * dm1; // if only the upper part of A is stored
424dm2 = A.selfadjointView<Lower>() * dm1; // if only the lower part of A is stored
425 \endcode
426 - copy of triangular parts:
427 \code
428sm2 = sm1.selfadjointView<Upper>(); // makes a full selfadjoint matrix from the upper triangular part
Jitse Niesen8c71d732012-06-20 09:52:45 +0100429sm2.selfadjointView<Lower>() = sm1.selfadjointView<Upper>(); // copies the upper triangular part to the lower triangular part
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200430 \endcode
431 - application of symmetric permutations:
432 \code
433PermutationMatrix<Dynamic,Dynamic> P = ...;
434sm2 = A.selfadjointView<Upper>().twistedBy(P); // compute P S P' from the upper triangular part of A, and make it a full matrix
435sm2.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 Holy44b19b42012-02-08 22:11:12 +0100449\li \b Next: \ref TutorialMapClass
Gael Guennebaud86a19262009-01-16 17:20:12 +0000450
451*/
452
453}