blob: 86a6d1430e4fe6ffcbff6bb34d553a06a08c3b01 [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
11 - \ref TutorialSparseFilling
12 - \ref TutorialSparseFeatureSet
13 - \ref TutorialSparseDirectSolvers
14<hr>
15
Gael Guennebaud70206ab2011-11-24 17:32:30 +010016Manipulating and solving sparse problems involves various modules which are summarized below:
Gael Guennebaud86a19262009-01-16 17:20:12 +000017
Gael Guennebaudf66fe262010-10-19 11:40:49 +020018<table class="manual">
Gael Guennebaud70206ab2011-11-24 17:32:30 +010019<tr><th>Module</th><th>Header file</th><th>Contents</th></tr>
20<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>
21<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>
22<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>
23<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 +000024</table>
25
Gael Guennebaud3f56de22011-12-03 10:26:00 +010026\section TutorialSparseIntro Sparse matrix representation
Gael Guennebaud70206ab2011-11-24 17:32:30 +010027
28In 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.
29
30\b The \b SparseMatrix \b class
31The class SparseMatrix is the main sparse matrix representation of Eigen's sparse module; it offers high performance and low memory usage.
32It implements a more versatile variant of the widely-used Compressed Column (or Row) Storage scheme.
33It consists of four compact arrays:
34 - \c Values: stores the coefficient values of the non-zeros.
35 - \c InnerIndices: stores the row (resp. column) indices of the non-zeros.
36 - \c OuterIndexPtrs: stores for each colmun (resp. row) the index of the first non zero in the previous arrays.
Gael Guennebaud3f56de22011-12-03 10:26:00 +010037 - \c InnerNNZs: stores the number of non-zeros of each column (resp. row).
Gael Guennebaud70206ab2011-11-24 17:32:30 +010038The 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.
39The word \c outer refers to the other direction.
40
41This storage scheme is better explained on an example. The following matrix
42<table class="manual">
43<tr><td> 0</td><td>3</td><td> 0</td><td>0</td><td> 0</td></tr>
44<tr><td>22</td><td>0</td><td> 0</td><td>0</td><td>17</td></tr>
45<tr><td> 7</td><td>5</td><td> 0</td><td>1</td><td> 0</td></tr>
46<tr><td> 0</td><td>0</td><td> 0</td><td>0</td><td> 0</td></tr>
47<tr><td> 0</td><td>0</td><td>14</td><td>0</td><td> 8</td></tr>
48</table>
49
Gael Guennebaud3f56de22011-12-03 10:26:00 +010050and one of its possible sparse, \b column \b major representation:
Gael Guennebaud70206ab2011-11-24 17:32:30 +010051<table class="manual">
52<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>
53<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>
54</table>
55<table class="manual">
56<tr><td>OuterIndexPtrs:</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 +010057<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 +010058</table>
59
Gael Guennebaud3f56de22011-12-03 10:26:00 +010060Currently the elements of a given inner vector are guaranteed to be always sorted by increasing inner indices.
61The \c "_" indicates available free space to quickly insert new elements.
62Assuming 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.
63On 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 +010064
Gael Guennebaud3f56de22011-12-03 10:26:00 +010065The case where no empty space is available is a special case, and is refered as the \em compressed mode.
66It corresponds to the widely used Compressed Column (or Row) Storage schemes (CCS or CRS).
67Any SparseMatrix can be turned to this form by calling the SparseMatrix::makeCompressed() function.
68In this case, one can remark that the \c InnerNNZs array is redundant with \c OuterIndexPtrs because we the equality: \c InnerNNZs[j] = \c OuterIndexPtrs[j+1]-\c OuterIndexPtrs[j].
69Therefore, in practice a call to SparseMatrix::makeCompressed() frees this buffer.
70
71It is worth noting that most of our wrappers to external libraries requires compressed matrices as inputs.
72
73The results of Eigen's operations always produces \b compressed sparse matrices.
74On 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 +010075
76Here is the previous matrix represented in compressed mode:
Gael Guennebaudf66fe262010-10-19 11:40:49 +020077<table class="manual">
Gael Guennebaud86a19262009-01-16 17:20:12 +000078<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 +010079<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 +000080</table>
Gael Guennebaudf66fe262010-10-19 11:40:49 +020081<table class="manual">
Gael Guennebaud70206ab2011-11-24 17:32:30 +010082<tr><td>OuterIndexPtrs:</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 +000083</table>
Gael Guennebaud70206ab2011-11-24 17:32:30 +010084
Gael Guennebaud3f56de22011-12-03 10:26:00 +010085A SparseVector is a special case of a SparseMatrix where only the \c Values and \c InnerIndices arrays are stored.
86There is no notion of compressed/uncompressed mode for a SparseVector.
Gael Guennebaud25f16582009-01-21 17:44:58 +000087
Gael Guennebaud86a19262009-01-16 17:20:12 +000088
89\b Matrix \b and \b vector \b properties \n
90
Tim Holy16a2d892011-06-20 22:47:58 -050091Here mat and vec represent any sparse-matrix and sparse-vector type, respectively.
Gael Guennebaud25f16582009-01-21 17:44:58 +000092
Gael Guennebaud70206ab2011-11-24 17:32:30 +010093Declarations:
94\code
Gael Guennebaud3f56de22011-12-03 10:26:00 +010095SparseMatrix<std::complex<float> > mat(1000,2000); // declares a 1000x2000 col-major compressed sparse matrix of complex<float>
96SparseMatrix<double,RowMajor> mat(1000,2000); // declares a 1000x2000 row-major compressed sparse matrix of double
97SparseVector<std::complex<float> > vec(1000); // declares a column sparse vector of complex<float> of size 1000
98SparseVector<double,RowMajor> vec(1000); // declares a row sparse vector of double of size 1000
Gael Guennebaud70206ab2011-11-24 17:32:30 +010099\endcode
100
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200101<table class="manual">
Gael Guennebaud86a19262009-01-16 17:20:12 +0000102<tr><td>Standard \n dimensions</td><td>\code
103mat.rows()
104mat.cols()\endcode</td>
105<td>\code
106vec.size() \endcode</td>
107</tr>
108<tr><td>Sizes along the \n inner/outer dimensions</td><td>\code
109mat.innerSize()
110mat.outerSize()\endcode</td>
111<td></td>
112</tr>
Tim Holy16a2d892011-06-20 22:47:58 -0500113<tr><td>Number of non \n zero coefficients</td><td>\code
Gael Guennebaud86a19262009-01-16 17:20:12 +0000114mat.nonZeros() \endcode</td>
115<td>\code
116vec.nonZeros() \endcode</td></tr>
117</table>
118
Gael Guennebaud25f16582009-01-21 17:44:58 +0000119
Gael Guennebaud86a19262009-01-16 17:20:12 +0000120\b Iterating \b over \b the \b nonzero \b coefficients \n
121
Tim Holy16a2d892011-06-20 22:47:58 -0500122Iterating over the coefficients of a sparse matrix can be done only in the same order as the storage order. Here is an example:
Gael Guennebaudf66fe262010-10-19 11:40:49 +0200123<table class="manual">
Gael Guennebaud86a19262009-01-16 17:20:12 +0000124<tr><td>
125\code
Gael Guennebaud25f16582009-01-21 17:44:58 +0000126SparseMatrixType mat(rows,cols);
Jitse Niesen59b83c12011-09-10 09:18:18 +0100127for (int k=0; k<mat.outerSize(); ++k)
Gael Guennebaud25f16582009-01-21 17:44:58 +0000128 for (SparseMatrixType::InnerIterator it(mat,k); it; ++it)
Gael Guennebaud86a19262009-01-16 17:20:12 +0000129 {
130 it.value();
131 it.row(); // row index
132 it.col(); // col index (here it is equal to k)
133 it.index(); // inner index, here it is equal to it.row()
134 }
135\endcode
136</td><td>
137\code
Gael Guennebaud25f16582009-01-21 17:44:58 +0000138SparseVector<double> vec(size);
139for (SparseVector<double>::InnerIterator it(vec); it; ++it)
Gael Guennebaud86a19262009-01-16 17:20:12 +0000140{
Gael Guennebaud25f16582009-01-21 17:44:58 +0000141 it.value(); // == vec[ it.index() ]
Gael Guennebaudc532f422009-09-25 16:30:44 +0200142 it.index();
Gael Guennebaud86a19262009-01-16 17:20:12 +0000143}
144\endcode
145</td></tr>
146</table>
147
Jitse Niesen59b83c12011-09-10 09:18:18 +0100148If the type of the sparse matrix or vector depends on a template parameter, then the \c typename keyword is
149required to indicate that \c InnerIterator denotes a type; see \ref TopicTemplateKeyword for details.
150
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100151
Gael Guennebaud86a19262009-01-16 17:20:12 +0000152\section TutorialSparseFilling Filling a sparse matrix
153
Gael Guennebaud1763f862012-02-04 10:44:07 +0100154
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100155Because of the special storage scheme of a SparseMatrix, special care has to be taken when adding new nonzero entries.
156For instance, the cost of inserting nnz non zeros in a a single purely random insertion into a SparseMatrix is O(nnz), where nnz is the current number of nonzero coefficients.
Gael Guennebaudc532f422009-09-25 16:30:44 +0200157
Gael Guennebaud1763f862012-02-04 10:44:07 +0100158The simplest way to create a sparse matrix while guarantying good performance is to first build a list of so called \em triplets, and then convert it to a SparseMatrix.
159
160Here is a typical usage example:
161\code
162typedef Triplet<double> T;
163std::vector<T> tripletList;
164triplets.reserve(estimation_of_entries);
165for(...)
166{
167 // ...
168 tripletList.push_back(T(i,j,v_ij));
169}
170SparseMatrixType m(rows,cols);
171m.setFromTriplets(tripletList.begin(), tripletList.end());
172// m is ready to go!
173\endcode
174The std::vector triplets might contain the elements in arbitrary order, and might even contain duplicated elements that will be summed up by setFromTriplets().
175See the SparseMatrix::setFromTriplets() function and class Triplet for more details.
176
177
178In some cases, however, slightly higher performance, and lower memory consumption can be reached by directly inserting the non zeros into the destination matrix.
179A typical scenario of this approach is illustrated bellow:
Gael Guennebaud25f16582009-01-21 17:44:58 +0000180\code
Gael Guennebaudc861e052011-12-03 18:14:51 +01001811: SparseMatrix<double> mat(rows,cols); // default is column major
Gael Guennebaud1aa6c7f2011-12-11 17:18:14 +01001822: mat.reserve(VectorXi::Constant(cols,6));
Gael Guennebaud70206ab2011-11-24 17:32:30 +01001833: for each i,j such that v_ij != 0
Gael Guennebaud3f56de22011-12-03 10:26:00 +01001844: mat.insert(i,j) = v_ij; // alternative: mat.coeffRef(i,j) += v_ij;
1855: mat.makeCompressed(); // optional
Gael Guennebaud25f16582009-01-21 17:44:58 +0000186\endcode
187
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100188- The key ingredient here is the line 2 where we reserve room for 6 non zeros per column. In many cases, the number of non zero per column or row can be easily 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.
189- 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.
190- 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.
191- The line 5 suppresses the remaining empty space and transforms the matrix into a compressed column storage.
Gael Guennebaud25f16582009-01-21 17:44:58 +0000192
Gael Guennebaud25f16582009-01-21 17:44:58 +0000193
Gael Guennebaud86a19262009-01-16 17:20:12 +0000194\section TutorialSparseFeatureSet Supported operators and functions
195
Tim Holy16a2d892011-06-20 22:47:58 -0500196In the following \em sm denotes a sparse matrix, \em sv a sparse vector, \em dm a dense matrix, and \em dv a dense vector.
197In Eigen's sparse module we chose to expose only the subset of the dense matrix API which can be efficiently implemented. Moreover, not every combination is allowed; for instance, it is not possible to add two sparse matrices having two different storage orders. On the other hand, it is perfectly fine to evaluate a sparse matrix or expression to a matrix having a different storage order:
Gael Guennebaud86a19262009-01-16 17:20:12 +0000198\code
199SparseMatrixType sm1, sm2, sm3;
Tim Holy16a2d892011-06-20 22:47:58 -0500200sm3 = sm1.transpose() + sm2; // invalid, because transpose() changes the storage order
201sm3 = SparseMatrixType(sm1.transpose()) + sm2; // correct, because evaluation reformats as column-major
Gael Guennebaud86a19262009-01-16 17:20:12 +0000202\endcode
203
Tim Holy16a2d892011-06-20 22:47:58 -0500204Here are some examples of supported operations:
Gael Guennebaud86a19262009-01-16 17:20:12 +0000205\code
Tim Holy16a2d892011-06-20 22:47:58 -0500206sm1 *= 0.5;
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100207sm1 = sm2 * 0.5;
208sm1 = sm2.transpose();
209sm1 = sm2.adjoint();
210sm4 = sm1 + sm2 + sm3; // only if sm1, sm2 and sm3 have the same storage order
211sm3 = sm1 * sm2; // conservative sparse * sparse product preserving numerical zeros
212sm3 = (sm1 * sm2).pruned(); // sparse * sparse product that removes numerical zeros (triggers a different algorithm)
213sm3 = (sm1 * sm2).pruned(ref); // sparse * sparse product that removes elements much smaller than ref
214sm3 = (sm1 * sm2).pruned(ref,epsilon); // sparse * sparse product that removes elements smaller than ref*epsilon
Gael Guennebaud86a19262009-01-16 17:20:12 +0000215dv3 = sm1 * dv2;
216dm3 = sm1 * dm2;
217dm3 = dm2 * sm1;
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100218sm3 = sm1.cwiseProduct(sm2); // only if sm1 and sm2 have the same storage order
Gael Guennebaud41ea92d2010-07-04 10:14:47 +0200219dv2 = sm1.triangularView<Upper>().solve(dv2);
Gael Guennebaud86a19262009-01-16 17:20:12 +0000220\endcode
221
Tim Holy16a2d892011-06-20 22:47:58 -0500222The product of a sparse \em symmetric matrix A with a dense matrix (or vector) d can be optimized by specifying the symmetry of A using selfadjointView:
Gael Guennebaud86a19262009-01-16 17:20:12 +0000223\code
Tim Holy16a2d892011-06-20 22:47:58 -0500224res = A.selfadjointView<>() * d; // if all coefficients of A are stored
225res = A.selfadjointView<Upper>() * d; // if only the upper part of A is stored
226res = A.selfadjointView<Lower>() * d; // if only the lower part of A is stored
Gael Guennebaud86a19262009-01-16 17:20:12 +0000227\endcode
228
229
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100230
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100231\section TutorialSparseDirectSolvers Solving linear problems
Gael Guennebaud86a19262009-01-16 17:20:12 +0000232
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100233Eigen currently provides a limited set of built-in solvers as well as wrappers to external solver libraries.
234They are summarized in the following table:
Gael Guennebaudab41c182010-08-18 15:33:58 +0200235
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100236<table class="manual">
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100237<tr><td>Class</td><td>Module</td><td>Solver kind</td><td>Matrix kind</td><td>Features related to performance</td>
238 <td>Dependencies,License</td><td class="width20em"><p>Notes</p></td></tr>
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100239<tr><td>SimplicialLLt </td><td>\link SparseCholesky_Module SparseCholesky \endlink</td><td>Direct LLt factorization</td><td>SPD</td><td>Fill-in reducing</td>
240 <td>built-in, LGPL</td>
241 <td>SimplicialLDLt is often preferable</td></tr>
242<tr><td>SimplicialLDLt </td><td>\link SparseCholesky_Module SparseCholesky \endlink</td><td>Direct LDLt factorization</td><td>SPD</td><td>Fill-in reducing</td>
243 <td>built-in, LGPL</td>
244 <td>Recommended for very sparse and not too large problems (e.g., 2D Poisson eq.)</td></tr>
245<tr><td>ConjugateGradient</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>Classic iterative CG</td><td>SPD</td><td>Preconditionning</td>
246 <td>built-in, LGPL</td>
247 <td>Recommended for large symmetric problems (e.g., 3D Poisson eq.)</td></tr>
248<tr><td>BiCGSTAB</td><td>\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink</td><td>Iterative stabilized bi-conjugate gradient</td><td>Square</td><td>Preconditionning</td>
249 <td>built-in, LGPL</td>
250 <td>Might not always converge</td></tr>
251
252
253
254<tr><td>CholmodDecomposition</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>
255 <td>Requires the <a href="http://www.cise.ufl.edu/research/sparse/SuiteSparse/">SuiteSparse</a> package, \b GPL </td>
256 <td></td></tr>
257<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>
258 <td>Requires the <a href="http://www.cise.ufl.edu/research/sparse/SuiteSparse/">SuiteSparse</a> package, \b GPL </td>
259 <td></td></tr>
260<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 +0100261 <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 +0100262 <td></td></tr>
263</table>
264
265Here \c SPD means symmetric positive definite.
266
267All these solvers follow the same general concept.
268Here is a typical and general example:
Gael Guennebaudab41c182010-08-18 15:33:58 +0200269\code
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100270#include <Eigen/RequiredModuleName>
Gael Guennebaudab41c182010-08-18 15:33:58 +0200271// ...
272SparseMatrix<double> A;
273// fill A
274VectorXd b, x;
275// fill b
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100276// solve Ax = b
277SolverClassName<SparseMatrix<double> > solver;
278solver.compute(A);
279if(solver.info()!=Succeeded) {
Tim Holy16a2d892011-06-20 22:47:58 -0500280 // decomposition failed
Gael Guennebaudab41c182010-08-18 15:33:58 +0200281 return;
282}
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100283x = solver.solve(b);
284if(solver.info()!=Succeeded) {
Gael Guennebaudab41c182010-08-18 15:33:58 +0200285 // solving failed
286 return;
287}
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100288// solve for another right hand side:
289x1 = solver.solve(b1);
Gael Guennebaudab41c182010-08-18 15:33:58 +0200290\endcode
291
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100292For \c SPD solvers, a second optional template argument allows to specify which triangular part have to be used, e.g.:
293
294\code
295#include <Eigen/IterativeLinearSolvers>
296
297ConjugateGradient<SparseMatrix<double>, Eigen::Upper> solver;
298x = solver.compute(A).solve(b);
299\endcode
300In 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.
301
302In the case where multiple problems with the same sparcity pattern have to be solved, then the "compute" step can be decomposed as follow:
303\code
304SolverClassName<SparseMatrix<double> > solver;
305solver.analyzePattern(A); // for this step the numerical values of A are not used
306solver.factorize(A);
307x1 = solver.solve(b1);
308x2 = solver.solve(b2);
309...
310A = ...; // modify the values of the nonzeros of A, the nonzeros pattern must stay unchanged
311solver.factorize(A);
312x1 = solver.solve(b1);
313x2 = solver.solve(b2);
314...
315\endcode
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100316The compute() method is equivalent to calling both analyzePattern() and factorize().
Gael Guennebaud70206ab2011-11-24 17:32:30 +0100317
Gael Guennebaud3f56de22011-12-03 10:26:00 +0100318Finally, each solver provides some specific features, such as determinant, access to the factors, controls of the iterations, and so on.
319More details are availble in the documentations of the respective classes.
Gael Guennebaud86a19262009-01-16 17:20:12 +0000320
Tim Holy44b19b42012-02-08 22:11:12 +0100321\li \b Next: \ref TutorialMapClass
Gael Guennebaud86a19262009-01-16 17:20:12 +0000322
323*/
324
325}