namespace Eigen { /** \page TutorialSparse Tutorial page 9 - Sparse Matrix \ingroup Tutorial \li \b Previous: \ref TutorialGeometry \li \b Next: TODO \b Table \b of \b contents \n - \ref TutorialSparseIntro - \ref TutorialSparseFilling - \ref TutorialSparseFeatureSet - \ref TutorialSparseDirectSolvers
Manipulating and solving sparse problems involves various modules which are summarized below:
ModuleHeader fileContents
\link Sparse_Module SparseCore \endlink\code#include \endcodeSparseMatrix and SparseVector classes, matrix assembly, basic sparse linear algebra (including sparse triangular solvers)
\link SparseCholesky_Module SparseCholesky \endlink\code#include \endcodeDirect sparse LLT and LDLT Cholesky factorization to solve sparse self-adjoint positive definite problems
\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlink\code#include \endcodeIterative solvers to solve large general linear square problems (including self-adjoint positive definite problems)
\code#include \endcodeIncludes all the above modules
\section TutorialSparseIntro Representing sparse matrices In many applications (e.g., finite element methods) it is common to deal with very large matrices where only a few coefficients are different from zero. In such cases, memory consumption can be reduced and performance increased by using a specialized representation storing only the nonzero coefficients. Such a matrix is called a sparse matrix. \b The \b SparseMatrix \b class The class SparseMatrix is the main sparse matrix representation of Eigen's sparse module; it offers high performance and low memory usage. It implements a more versatile variant of the widely-used Compressed Column (or Row) Storage scheme. It consists of four compact arrays: - \c Values: stores the coefficient values of the non-zeros. - \c InnerIndices: stores the row (resp. column) indices of the non-zeros. - \c OuterIndexPtrs: stores for each colmun (resp. row) the index of the first non zero in the previous arrays. - \c InnerSizes: stores the number of non-zeros of each column (resp. row). The word \c inner refers to an \em inner \em vector that is a column for a column-major matrix, or a row for a row-major matrix. The word \c outer refers to the other direction. This storage scheme is better explained on an example. The following matrix
03 00 0
220 0017
75 01 0
00 00 0
00140 8
and its sparse, \b column \b major representation:
Values: 227_3514__1_178
InnerIndices: 12_02 4__2_ 14
OuterIndexPtrs:035810\em 12
InnerSizes: 2211 2
By default the elements of a given inner vector are always sorted by increasing inner indices. The \c _ indicates available free space to quickly insert new elements. Assuming no reallocation is needed, the insertion of a random element of coordinate is therefore in O(nnz) where nnz is the number of nonzeros of the respective inner vector. On the other hand, inserting elements with increasing inner indices in a given inner vector is much more efficient since this only requires to incresae the respective \c InnerSizes entry that is a O(1) operation. The case where no empty space is available is a special case, and is refered as the \em compressed mode and it corresponds to the widely used Compressed Column (or Row) Storage scheme. Any SparseMatrix can be turned in this form by calling the SparseMatrix::makeCompressed() function. In this case one can remark that the \c InnerSizes array is superfluous because \c InnerSizes[j] is simply equals to \c OuterIndexPtrs[j+1]-\c OuterIndexPtrs[j]. Therefore, a call to SparseMatrix::makeCompressed() frees this buffer. Here is the previous matrix represented in compressed mode:
Values: 22735141178
InnerIndices: 1202 42 14
OuterIndexPtrs:02456\em 8
In this mode any insertion of new non zero elements would be extremely costly. \b Matrix \b and \b vector \b properties \n Here mat and vec represent any sparse-matrix and sparse-vector type, respectively. Declarations: \code SparseMatrix > mat(1000,2000); // declare a 1000x2000 col-major compressed sparse matrix of complex SparseMatrix mat(1000,2000); // declare a 1000x2000 row-major compressed sparse matrix of double SparseVector > vec(1000); // declare a column sparse vector of complex of size 1000 SparseVector vec(1000); // declare a row sparse vector of double of size 1000 \endcode Although a sparse matrix could also be used to represent a sparse vector, for that purpose it is better to use the specialized SparseVector class: \code \endcode
Standard \n dimensions\code mat.rows() mat.cols()\endcode \code vec.size() \endcode
Sizes along the \n inner/outer dimensions\code mat.innerSize() mat.outerSize()\endcode
Number of non \n zero coefficients\code mat.nonZeros() \endcode \code vec.nonZeros() \endcode
\b Iterating \b over \b the \b nonzero \b coefficients \n Iterating over the coefficients of a sparse matrix can be done only in the same order as the storage order. Here is an example:
\code SparseMatrixType mat(rows,cols); for (int k=0; k \code SparseVector vec(size); for (SparseVector::InnerIterator it(vec); it; ++it) { it.value(); // == vec[ it.index() ] it.index(); } \endcode
If the type of the sparse matrix or vector depends on a template parameter, then the \c typename keyword is required to indicate that \c InnerIterator denotes a type; see \ref TopicTemplateKeyword for details. \section TutorialSparseFilling Filling a sparse matrix Because of the special storage scheme of a SparseMatrix, special care has to be taken when adding new nonzero entries. For 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. A typical scenario to insert nonzeros is illustrated bellow: \code 1: SparseMatrix m(rows,cols); // default is column major 2: aux.reserve(VectorXi::Constant(rows,6)); 3: for each i,j such that v_ij != 0 4: mat.insert(i,j) = v_ij; 5: mat.makeCompressed(); // optional \endcode - 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 provifing a VectorXi or std::vector compatible object. 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. - The line 4 performs a sorted insertion. In this example, the ideal case is when the j-th column is not full and contains non-zeros whose inner-indices are smaller than i. In this case, this operation boils down to trivial O(1) operation. - The line 5 suppresses the left empty room and transform the matrix to a compressed column storage. \section TutorialSparseFeatureSet Supported operators and functions In the following \em sm denotes a sparse matrix, \em sv a sparse vector, \em dm a dense matrix, and \em dv a dense vector. In 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: \code SparseMatrixType sm1, sm2, sm3; sm3 = sm1.transpose() + sm2; // invalid, because transpose() changes the storage order sm3 = SparseMatrixType(sm1.transpose()) + sm2; // correct, because evaluation reformats as column-major \endcode Here are some examples of supported operations: \code sm1 *= 0.5; sm4 = sm1 + sm2 + sm3; // only if sm1, sm2 and sm3 have the same storage order sm3 = sm1 * sm2; dv3 = sm1 * dv2; dm3 = sm1 * dm2; dm3 = dm2 * sm1; sm3 = sm1.cwiseProduct(sm2); // only if sm1 and sm2 have the same storage order dv2 = sm1.triangularView().solve(dv2); \endcode The 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: \code res = A.selfadjointView<>() * d; // if all coefficients of A are stored res = A.selfadjointView() * d; // if only the upper part of A is stored res = A.selfadjointView() * d; // if only the lower part of A is stored \endcode \section TutorialSparseDirectSolvers Solving linear problems Eigen currently provides a limited set of built-in solvers as well as wrappers to external solver libraries. They are summarized in the following table:
Class ModuleSolver kindMatrix kindFeatures related to performance Dependencies,License Notes
SimplicialLLt \link SparseCholesky_Module SparseCholesky \endlinkDirect LLt factorizationSPDFill-in reducing built-in, LGPL SimplicialLDLt is often preferable
SimplicialLDLt \link SparseCholesky_Module SparseCholesky \endlinkDirect LDLt factorizationSPDFill-in reducing built-in, LGPL Recommended for very sparse and not too large problems (e.g., 2D Poisson eq.)
ConjugateGradient\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlinkClassic iterative CGSPDPreconditionning built-in, LGPL Recommended for large symmetric problems (e.g., 3D Poisson eq.)
BiCGSTAB\link IterativeLinearSolvers_Module IterativeLinearSolvers \endlinkIterative stabilized bi-conjugate gradientSquarePreconditionning built-in, LGPL Might not always converge
CholmodDecomposition\link CholmodSupport_Module CholmodSupport \endlinkDirect LLT factorizationSPDFill-in reducing, Leverage fast dense algebra Requires the SuiteSparse package, \b GPL
UmfPackLU\link UmfPackSupport_Module UmfPackSupport \endlinkDirect LU factorizationSquareFill-in reducing, Leverage fast dense algebra Requires the SuiteSparse package, \b GPL
SuperLU\link SuperLUSupport_Module SuperLUSupport \endlinkDirect LU factorizationSquareFill-in reducing, Leverage fast dense algebra Requires the SuperLU library, custom (BSD-like)
Here \c SPD means symmetric positive definite. All these solvers follow the same general concept. Here is a typical and general example: \code #include // ... SparseMatrix A; // fill A VectorXd b, x; // fill b // solve Ax = b SolverClassName > solver; solver.compute(A); if(solver.info()!=Succeeded) { // decomposition failed return; } x = solver.solve(b); if(solver.info()!=Succeeded) { // solving failed return; } // solve for another right hand side: x1 = solver.solve(b1); \endcode For \c SPD solvers, a second optional template argument allows to specify which triangular part have to be used, e.g.: \code #include ConjugateGradient, Eigen::Upper> solver; x = solver.compute(A).solve(b); \endcode In the above example, only the upper triangular part of the input matrix A is considered for solving. The opposite triangle might either be empty or contain arbitrary values. In the case where multiple problems with the same sparcity pattern have to be solved, then the "compute" step can be decomposed as follow: \code SolverClassName > solver; solver.analyzePattern(A); // for this step the numerical values of A are not used solver.factorize(A); x1 = solver.solve(b1); x2 = solver.solve(b2); ... A = ...; // modify the values of the nonzeros of A, the nonzeros pattern must stay unchanged solver.factorize(A); x1 = solver.solve(b1); x2 = solver.solve(b2); ... \endcode The compute() methode is equivalent to calling both analyzePattern() and factorize(). Finally, each solver provides some specific features, such as determinant, access to the factors, controls of the iterations, etc. More details are availble on the documentation of the respective classes. \li \b Next: TODO */ }