blob: 9f06005fa5140efb7fecf4b447623057f4facfee [file] [log] [blame]
Gael Guennebaud86a19262009-01-16 17:20:12 +00001namespace Eigen {
2
Gael Guennebaud93ee82b2013-01-05 16:37:11 +01003/** \eigenManualPage TutorialSparse Sparse matrix manipulations
Gael Guennebaud86a19262009-01-16 17:20:12 +00004
Gael Guennebaud93ee82b2013-01-05 16:37:11 +01005\eigeneigenAutoToc
Gael Guennebaud52dce0c2012-06-20 09:28:32 +02006
Gael Guennebaud70206ab2011-11-24 17:32:30 +01007Manipulating and solving sparse problems involves various modules which are summarized below:
Gael Guennebaud86a19262009-01-16 17:20:12 +00008
Gael Guennebaudf66fe262010-10-19 11:40:49 +02009<table class="manual">
Gael Guennebaud70206ab2011-11-24 17:32:30 +010010<tr><th>Module</th><th>Header file</th><th>Contents</th></tr>
Gael Guennebaudf41d96d2012-12-24 13:33:22 +010011<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 +010012<tr><td>\link SparseCholesky_Module SparseCholesky \endlink</td><td>\code#include <Eigen/SparseCholesky>\endcode</td><td>Direct sparse LLT and LDLT Cholesky factorization to solve sparse self-adjoint positive definite problems</td></tr>
13<tr><td>\link 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 +010014<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 +000015</table>
16
Gael Guennebaud3f56de22011-12-03 10:26:00 +010017\section TutorialSparseIntro Sparse matrix representation
Gael Guennebaud70206ab2011-11-24 17:32:30 +010018
19In 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.
20
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020021\b The \b %SparseMatrix \b class
22
Gael Guennebaud70206ab2011-11-24 17:32:30 +010023The class SparseMatrix is the main sparse matrix representation of Eigen's sparse module; it offers high performance and low memory usage.
24It implements a more versatile variant of the widely-used Compressed Column (or Row) Storage scheme.
25It consists of four compact arrays:
26 - \c Values: stores the coefficient values of the non-zeros.
27 - \c InnerIndices: stores the row (resp. column) indices of the non-zeros.
Jitse Niesen8c71d732012-06-20 09:52:45 +010028 - \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 +010029 - \c InnerNNZs: stores the number of non-zeros of each column (resp. row).
Gael Guennebaud70206ab2011-11-24 17:32:30 +010030The 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.
31The word \c outer refers to the other direction.
32
33This storage scheme is better explained on an example. The following matrix
34<table class="manual">
35<tr><td> 0</td><td>3</td><td> 0</td><td>0</td><td> 0</td></tr>
36<tr><td>22</td><td>0</td><td> 0</td><td>0</td><td>17</td></tr>
37<tr><td> 7</td><td>5</td><td> 0</td><td>1</td><td> 0</td></tr>
38<tr><td> 0</td><td>0</td><td> 0</td><td>0</td><td> 0</td></tr>
39<tr><td> 0</td><td>0</td><td>14</td><td>0</td><td> 8</td></tr>
40</table>
41
Gael Guennebaud3f56de22011-12-03 10:26:00 +010042and one of its possible sparse, \b column \b major representation:
Gael Guennebaud70206ab2011-11-24 17:32:30 +010043<table class="manual">
44<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>
45<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>
46</table>
47<table class="manual">
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020048<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 +010049<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 +010050</table>
51
Gael Guennebaud3f56de22011-12-03 10:26:00 +010052Currently the elements of a given inner vector are guaranteed to be always sorted by increasing inner indices.
53The \c "_" indicates available free space to quickly insert new elements.
54Assuming 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.
55On 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 +010056
Gael Guennebaud3f56de22011-12-03 10:26:00 +010057The case where no empty space is available is a special case, and is refered as the \em compressed mode.
58It corresponds to the widely used Compressed Column (or Row) Storage schemes (CCS or CRS).
59Any SparseMatrix can be turned to this form by calling the SparseMatrix::makeCompressed() function.
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020060In 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 +010061Therefore, in practice a call to SparseMatrix::makeCompressed() frees this buffer.
62
63It is worth noting that most of our wrappers to external libraries requires compressed matrices as inputs.
64
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020065The results of %Eigen's operations always produces \b compressed sparse matrices.
Gael Guennebaud3f56de22011-12-03 10:26:00 +010066On 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 +010067
68Here is the previous matrix represented in compressed mode:
Gael Guennebaudf66fe262010-10-19 11:40:49 +020069<table class="manual">
Gael Guennebaud86a19262009-01-16 17:20:12 +000070<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 +010071<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 +000072</table>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020073<table class="manual">
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020074<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 +000075</table>
Gael Guennebaud70206ab2011-11-24 17:32:30 +010076
Gael Guennebaud3f56de22011-12-03 10:26:00 +010077A SparseVector is a special case of a SparseMatrix where only the \c Values and \c InnerIndices arrays are stored.
78There is no notion of compressed/uncompressed mode for a SparseVector.
Gael Guennebaud25f16582009-01-21 17:44:58 +000079
Gael Guennebaud86a19262009-01-16 17:20:12 +000080
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020081\section TutorialSparseExample First example
Gael Guennebaud86a19262009-01-16 17:20:12 +000082
Jitse Niesen8c71d732012-06-20 09:52:45 +010083Before 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.
84Such 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 +000085
Gael Guennebaud52dce0c2012-06-20 09:28:32 +020086<table class="manual">
87<tr><td>
88\include Tutorial_sparse_example.cpp
89</td>
90<td>
91\image html Tutorial_sparse_example.jpeg
92</td></tr></table>
93
94In 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.
95
96In 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.
97The raw and flat list of non-zero entries is then converted to a true SparseMatrix object \c A.
98Note that the elements of the list do not have to be sorted, and possible duplicate entries will be summed up.
99
Jitse Niesen8c71d732012-06-20 09:52:45 +0100100The last step consists of effectively solving the assembled problem.
101Since 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 +0200102
103The 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.
104
Jitse Niesen8c71d732012-06-20 09:52:45 +0100105Describing 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 +0200106
107
108
109
110\section TutorialSparseSparseMatrix The SparseMatrix class
111
112\b %Matrix \b and \b vector \b properties \n
113
114The SparseMatrix and SparseVector classes take three template arguments:
115 * the scalar type (e.g., double)
116 * the storage order (ColMajor or RowMajor, the default is RowMajor)
117 * the inner index type (default is \c int).
118
119As for dense Matrix objects, constructors takes the size of the object.
120Here are some examples:
121
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100122\code
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200123SparseMatrix<std::complex<float> > mat(1000,2000); // declares a 1000x2000 column-major compressed sparse matrix of complex<float>
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100124SparseMatrix<double,RowMajor> mat(1000,2000); // declares a 1000x2000 row-major compressed sparse matrix of double
125SparseVector<std::complex<float> > vec(1000); // declares a column sparse vector of complex<float> of size 1000
126SparseVector<double,RowMajor> vec(1000); // declares a row sparse vector of double of size 1000
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100127\endcode
128
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200129In the rest of the tutorial, \c mat and \c vec represent any sparse-matrix and sparse-vector objects, respectively.
130
131The dimensions of a matrix can be queried using the following functions:
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200132<table class="manual">
Gael Guennebaud86a19262009-01-16 17:20:12 +0000133<tr><td>Standard \n dimensions</td><td>\code
134mat.rows()
135mat.cols()\endcode</td>
136<td>\code
137vec.size() \endcode</td>
138</tr>
139<tr><td>Sizes along the \n inner/outer dimensions</td><td>\code
140mat.innerSize()
141mat.outerSize()\endcode</td>
142<td></td>
143</tr>
Tim Holy16a2d892011-06-20 22:47:58 -0500144<tr><td>Number of non \n zero coefficients</td><td>\code
Gael Guennebaud86a19262009-01-16 17:20:12 +0000145mat.nonZeros() \endcode</td>
146<td>\code
147vec.nonZeros() \endcode</td></tr>
148</table>
149
Gael Guennebaud25f16582009-01-21 17:44:58 +0000150
Gael Guennebaud86a19262009-01-16 17:20:12 +0000151\b Iterating \b over \b the \b nonzero \b coefficients \n
152
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200153Random access to the elements of a sparse object can be done through the \c coeffRef(i,j) function.
154However, this function involves a quite expensive binary search.
155In 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.
156Here is an example:
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200157<table class="manual">
Gael Guennebaud86a19262009-01-16 17:20:12 +0000158<tr><td>
159\code
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200160SparseMatrix<double> mat(rows,cols);
Jitse Niesen59b83c12011-09-10 09:18:18 +0100161for (int k=0; k<mat.outerSize(); ++k)
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200162 for (SparseMatrix<double>::InnerIterator it(mat,k); it; ++it)
Gael Guennebaud86a19262009-01-16 17:20:12 +0000163 {
164 it.value();
165 it.row(); // row index
166 it.col(); // col index (here it is equal to k)
167 it.index(); // inner index, here it is equal to it.row()
168 }
169\endcode
170</td><td>
171\code
Gael Guennebaud25f16582009-01-21 17:44:58 +0000172SparseVector<double> vec(size);
173for (SparseVector<double>::InnerIterator it(vec); it; ++it)
Gael Guennebaud86a19262009-01-16 17:20:12 +0000174{
Gael Guennebaud25f16582009-01-21 17:44:58 +0000175 it.value(); // == vec[ it.index() ]
Gael Guennebaudc532f422009-09-25 16:30:44 +0200176 it.index();
Gael Guennebaud86a19262009-01-16 17:20:12 +0000177}
178\endcode
179</td></tr>
180</table>
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200181For a writable expression, the referenced value can be modified using the valueRef() function.
Jitse Niesen59b83c12011-09-10 09:18:18 +0100182If the type of the sparse matrix or vector depends on a template parameter, then the \c typename keyword is
183required to indicate that \c InnerIterator denotes a type; see \ref TopicTemplateKeyword for details.
184
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100185
Gael Guennebaud86a19262009-01-16 17:20:12 +0000186\section TutorialSparseFilling Filling a sparse matrix
187
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100188Because 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 +0200189For 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 +0200190
Jitse Niesen8c71d732012-06-20 09:52:45 +0100191The 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 +0100192
193Here is a typical usage example:
194\code
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200195typedef Eigen::Triplet<double> T;
Gael Guennebaud1763f862012-02-04 10:44:07 +0100196std::vector<T> tripletList;
Gael Guennebaud1edb3962012-09-26 19:24:41 +0200197tripletList.reserve(estimation_of_entries);
Gael Guennebaud1763f862012-02-04 10:44:07 +0100198for(...)
199{
200 // ...
201 tripletList.push_back(T(i,j,v_ij));
202}
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200203SparseMatrixType mat(rows,cols);
204mat.setFromTriplets(tripletList.begin(), tripletList.end());
205// mat is ready to go!
Gael Guennebaud1763f862012-02-04 10:44:07 +0100206\endcode
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200207The \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 +0100208See the SparseMatrix::setFromTriplets() function and class Triplet for more details.
209
210
Jitse Niesen8c71d732012-06-20 09:52:45 +0100211In 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 +0100212A typical scenario of this approach is illustrated bellow:
Gael Guennebaud25f16582009-01-21 17:44:58 +0000213\code
Gael Guennebaudc861e052011-12-03 18:14:51 +01002141: SparseMatrix<double> mat(rows,cols); // default is column major
Gael Guennebaud1aa6c7f2011-12-11 17:18:14 +01002152: mat.reserve(VectorXi::Constant(cols,6));
Gael Guennebaud70206ab2011-11-24 17:32:30 +01002163: for each i,j such that v_ij != 0
Gael Guennebaud3f56de22011-12-03 10:26:00 +01002174: mat.insert(i,j) = v_ij; // alternative: mat.coeffRef(i,j) += v_ij;
2185: mat.makeCompressed(); // optional
Gael Guennebaud25f16582009-01-21 17:44:58 +0000219\endcode
220
Jitse Niesen8c71d732012-06-20 09:52:45 +0100221- 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 +0100222- 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.
223- 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.
224- The line 5 suppresses the remaining empty space and transforms the matrix into a compressed column storage.
Gael Guennebaud25f16582009-01-21 17:44:58 +0000225
Gael Guennebaud25f16582009-01-21 17:44:58 +0000226
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100227\section TutorialSparseDirectSolvers Solving linear problems
Gael Guennebaud86a19262009-01-16 17:20:12 +0000228
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200229%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 +0100230They are summarized in the following table:
Gael Guennebaudab41c182010-08-18 15:33:58 +0200231
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100232<table class="manual">
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200233<tr><th>Class</th><th>Module</th><th>Solver kind</th><th>Matrix kind</th><th>Features related to performance</th>
234 <th>Dependencies,License</th><th class="width20em"><p>Notes</p></th></tr>
235<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 +0100236 <td>built-in, LGPL</td>
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200237 <td>SimplicialLDLT is often preferable</td></tr>
238<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 +0100239 <td>built-in, LGPL</td>
240 <td>Recommended for very sparse and not too large problems (e.g., 2D Poisson eq.)</td></tr>
241<tr><td>ConjugateGradient</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>Classic iterative CG</td><td>SPD</td><td>Preconditionning</td>
242 <td>built-in, LGPL</td>
243 <td>Recommended for large symmetric problems (e.g., 3D Poisson eq.)</td></tr>
244<tr><td>BiCGSTAB</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>Iterative stabilized bi-conjugate gradient</td><td>Square</td><td>Preconditionning</td>
245 <td>built-in, LGPL</td>
246 <td>Might not always converge</td></tr>
247
248
Gael Guennebaud6f3057f2012-06-21 10:51:22 +0200249<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>
250 <td>Requires the <a href="http://pastix.gforge.inria.fr">PaStiX</a> package, \b CeCILL-C </td>
251 <td>optimized for tough problems and symmetric patterns</td></tr>
252<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 +0100253 <td>Requires the <a href="http://www.cise.ufl.edu/research/sparse/SuiteSparse/">SuiteSparse</a> package, \b GPL </td>
254 <td></td></tr>
255<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>
256 <td>Requires the <a href="http://www.cise.ufl.edu/research/sparse/SuiteSparse/">SuiteSparse</a> package, \b GPL </td>
257 <td></td></tr>
258<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 +0100259 <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 +0100260 <td></td></tr>
261</table>
262
263Here \c SPD means symmetric positive definite.
264
265All these solvers follow the same general concept.
266Here is a typical and general example:
Gael Guennebaudab41c182010-08-18 15:33:58 +0200267\code
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100268#include <Eigen/RequiredModuleName>
Gael Guennebaudab41c182010-08-18 15:33:58 +0200269// ...
270SparseMatrix<double> A;
271// fill A
272VectorXd b, x;
273// fill b
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100274// solve Ax = b
275SolverClassName<SparseMatrix<double> > solver;
276solver.compute(A);
Gael Guennebaudf23dd7c2012-07-28 13:07:45 +0200277if(solver.info()!=Success) {
Tim Holy16a2d892011-06-20 22:47:58 -0500278 // decomposition failed
Gael Guennebaudab41c182010-08-18 15:33:58 +0200279 return;
280}
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100281x = solver.solve(b);
Gael Guennebaudf23dd7c2012-07-28 13:07:45 +0200282if(solver.info()!=Success) {
Gael Guennebaudab41c182010-08-18 15:33:58 +0200283 // solving failed
284 return;
285}
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100286// solve for another right hand side:
287x1 = solver.solve(b1);
Gael Guennebaudab41c182010-08-18 15:33:58 +0200288\endcode
289
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100290For \c SPD solvers, a second optional template argument allows to specify which triangular part have to be used, e.g.:
291
292\code
293#include <Eigen/IterativeLinearSolvers>
294
295ConjugateGradient<SparseMatrix<double>, Eigen::Upper> solver;
296x = solver.compute(A).solve(b);
297\endcode
298In 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.
299
300In the case where multiple problems with the same sparcity pattern have to be solved, then the "compute" step can be decomposed as follow:
301\code
302SolverClassName<SparseMatrix<double> > solver;
303solver.analyzePattern(A); // for this step the numerical values of A are not used
304solver.factorize(A);
305x1 = solver.solve(b1);
306x2 = solver.solve(b2);
307...
308A = ...; // modify the values of the nonzeros of A, the nonzeros pattern must stay unchanged
309solver.factorize(A);
310x1 = solver.solve(b1);
311x2 = solver.solve(b2);
312...
313\endcode
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100314The compute() method is equivalent to calling both analyzePattern() and factorize().
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100315
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100316Finally, each solver provides some specific features, such as determinant, access to the factors, controls of the iterations, and so on.
317More details are availble in the documentations of the respective classes.
Gael Guennebaud86a19262009-01-16 17:20:12 +0000318
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200319
320\section TutorialSparseFeatureSet Supported operators and functions
321
322Because of their special storage format, sparse matrices cannot offer the same level of flexbility than dense matrices.
323In Eigen's sparse module we chose to expose only the subset of the dense matrix API which can be efficiently implemented.
324In the following \em sm denotes a sparse matrix, \em sv a sparse vector, \em dm a dense matrix, and \em dv a dense vector.
325
326\subsection TutorialSparse_BasicOps Basic operations
327
328%Sparse expressions support most of the unary and binary coefficient wise operations:
329\code
330sm1.real() sm1.imag() -sm1 0.5*sm1
331sm1+sm2 sm1-sm2 sm1.cwiseProduct(sm2)
332\endcode
333However, a strong restriction is that the storage orders must match. For instance, in the following example:
334\code
335sm4 = sm1 + sm2 + sm3;
336\endcode
337sm1, sm2, and sm3 must all be row-major or all column major.
338On the other hand, there is no restriction on the target matrix sm4.
339For 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:
340\code
341SparseMatrix<double> A, B;
342B = SparseMatrix<double>(A.transpose()) + A;
343\endcode
344
345Binary coefficient wise operators can also mix sparse and dense expressions:
346\code
347sm2 = sm1.cwiseProduct(dm1);
348dm2 = sm1 + dm1;
349\endcode
350
351
352%Sparse expressions also support transposition:
353\code
354sm1 = sm2.transpose();
355sm1 = sm2.adjoint();
356\endcode
357However, there is no transposeInPlace() method.
358
359
360\subsection TutorialSparse_Products Matrix products
361
362%Eigen supports various kind of sparse matrix products which are summarize below:
363 - \b sparse-dense:
364 \code
365dv2 = sm1 * dv1;
366dm2 = dm1 * sm1.adjoint();
367dm2 = 2. * sm1 * dm1;
368 \endcode
369 - \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():
370 \code
371dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored
372dm2 = A.selfadjointView<Upper>() * dm1; // if only the upper part of A is stored
373dm2 = A.selfadjointView<Lower>() * dm1; // if only the lower part of A is stored
374 \endcode
375 - \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:
376 \code
377sm3 = sm1 * sm2;
378sm3 = 4 * sm1.adjoint() * sm2;
379 \endcode
Jitse Niesen8c71d732012-06-20 09:52:45 +0100380 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 +0200381 \code
382sm3 = (sm1 * sm2).prune(); // removes numerical zeros
383sm3 = (sm1 * sm2).prune(ref); // removes elements much smaller than ref
Jitse Niesen8c71d732012-06-20 09:52:45 +0100384sm3 = (sm1 * sm2).prune(ref,epsilon); // removes elements smaller than ref*epsilon
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200385 \endcode
386
387 - \b permutations. Finally, permutations can be applied to sparse matrices too:
388 \code
389PermutationMatrix<Dynamic,Dynamic> P = ...;
390sm2 = P * sm1;
391sm2 = sm1 * P.inverse();
392sm2 = sm1.transpose() * P;
393 \endcode
394
395
396\subsection TutorialSparse_TriangularSelfadjoint Triangular and selfadjoint views
397
Jitse Niesen8c71d732012-06-20 09:52:45 +0100398Just 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 +0200399\code
400dm2 = sm1.triangularView<Lower>(dm1);
401dv2 = sm1.transpose().triangularView<Upper>(dv1);
402\endcode
403
404The selfadjointView() function permits various operations:
405 - optimized sparse-dense matrix products:
406 \code
407dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored
408dm2 = A.selfadjointView<Upper>() * dm1; // if only the upper part of A is stored
409dm2 = A.selfadjointView<Lower>() * dm1; // if only the lower part of A is stored
410 \endcode
411 - copy of triangular parts:
412 \code
413sm2 = sm1.selfadjointView<Upper>(); // makes a full selfadjoint matrix from the upper triangular part
Jitse Niesen8c71d732012-06-20 09:52:45 +0100414sm2.selfadjointView<Lower>() = sm1.selfadjointView<Upper>(); // copies the upper triangular part to the lower triangular part
Gael Guennebaud52dce0c2012-06-20 09:28:32 +0200415 \endcode
416 - application of symmetric permutations:
417 \code
418PermutationMatrix<Dynamic,Dynamic> P = ...;
419sm2 = A.selfadjointView<Upper>().twistedBy(P); // compute P S P' from the upper triangular part of A, and make it a full matrix
420sm2.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
421 \endcode
422
423\subsection TutorialSparse_Submat Sub-matrices
424
425%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:
426\code
427 sm1.innerVector(j); // returns an expression of the j-th column (resp. row) of the matrix if sm1 is col-major (resp. row-major)
428 sm1.innerVectors(j, nb); // returns an expression of the nb columns (resp. row) starting from the j-th column (resp. row)
429 // of the matrix if sm1 is col-major (resp. row-major)
430 sm1.middleRows(j, nb); // for row major matrices only, get a range of nb rows
431 sm1.middleCols(j, nb); // for column major matrices only, get a range of nb columns
432\endcode
433
Gael Guennebaud86a19262009-01-16 17:20:12 +0000434*/
435
436}