namespace Eigen { /** \page TutorialSparse Tutorial page 9 - Sparse Matrix \ingroup Tutorial \li \b Previous: \ref TutorialGeometry \li \b Next: \ref TutorialMapClass \b Table \b of \b contents \n - \ref TutorialSparseIntro - \ref TutorialSparseExample "Example" - \ref TutorialSparseSparseMatrix - \ref TutorialSparseFilling - \ref TutorialSparseDirectSolvers - \ref TutorialSparseFeatureSet - \ref TutorialSparse_BasicOps - \ref TutorialSparse_Products - \ref TutorialSparse_TriangularSelfadjoint - \ref TutorialSparse_Submat
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 Sparse matrix representation 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 OuterStarts: stores for each column (resp. row) the index of the first non-zero in the previous two arrays. - \c InnerNNZs: 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 one of its possible sparse, \b column \b major representation:
Values: 227_3514__1_178
InnerIndices: 12_02 4__2_ 14
OuterStarts:035810\em 12
InnerNNZs: 2211 2
Currently the elements of a given inner vector are guaranteed to be 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 is therefore in O(nnz_j) where nnz_j 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 increase the respective \c InnerNNZs 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. It corresponds to the widely used Compressed Column (or Row) Storage schemes (CCS or CRS). Any SparseMatrix can be turned to this form by calling the SparseMatrix::makeCompressed() function. In 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]. Therefore, in practice a call to SparseMatrix::makeCompressed() frees this buffer. It is worth noting that most of our wrappers to external libraries requires compressed matrices as inputs. The results of %Eigen's operations always produces \b compressed sparse matrices. On the other hand, the insertion of a new element into a SparseMatrix converts this later to the \b uncompressed mode. Here is the previous matrix represented in compressed mode:
Values: 22735141178
InnerIndices: 1202 42 14
OuterStarts:02456\em 8
A SparseVector is a special case of a SparseMatrix where only the \c Values and \c InnerIndices arrays are stored. There is no notion of compressed/uncompressed mode for a SparseVector. \section TutorialSparseExample First example Before 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. Such 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.
\include Tutorial_sparse_example.cpp \image html Tutorial_sparse_example.jpeg
In this example, we start by defining a column-major sparse matrix type of double \c SparseMatrix, and a triplet list of the same scalar type \c Triplet. A triplet is a simple object representing a non-zero entry as the triplet: \c row index, \c column index, \c value. In 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. The raw and flat list of non-zero entries is then converted to a true SparseMatrix object \c A. Note that the elements of the list do not have to be sorted, and possible duplicate entries will be summed up. The last step consists of effectively solving the assembled problem. Since 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. The 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. Describing 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. \section TutorialSparseSparseMatrix The SparseMatrix class \b %Matrix \b and \b vector \b properties \n The SparseMatrix and SparseVector classes take three template arguments: * the scalar type (e.g., double) * the storage order (ColMajor or RowMajor, the default is RowMajor) * the inner index type (default is \c int). As for dense Matrix objects, constructors takes the size of the object. Here are some examples: \code SparseMatrix > mat(1000,2000); // declares a 1000x2000 column-major compressed sparse matrix of complex SparseMatrix mat(1000,2000); // declares a 1000x2000 row-major compressed sparse matrix of double SparseVector > vec(1000); // declares a column sparse vector of complex of size 1000 SparseVector vec(1000); // declares a row sparse vector of double of size 1000 \endcode In the rest of the tutorial, \c mat and \c vec represent any sparse-matrix and sparse-vector objects, respectively. The dimensions of a matrix can be queried using the following functions:
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 Random access to the elements of a sparse object can be done through the \c coeffRef(i,j) function. However, this function involves a quite expensive binary search. In 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. Here is an example:
\code SparseMatrix mat(rows,cols); for (int k=0; k::InnerIterator it(mat,k); it; ++it) { it.value(); it.row(); // row index it.col(); // col index (here it is equal to k) it.index(); // inner index, here it is equal to it.row() } \endcode \code SparseVector vec(size); for (SparseVector::InnerIterator it(vec); it; ++it) { it.value(); // == vec[ it.index() ] it.index(); } \endcode
For a writable expression, the referenced value can be modified using the valueRef() function. 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 a single purely random insertion into a SparseMatrix is \c O(nnz), where \c nnz is the current number of non-zero coefficients. The 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. Here is a typical usage example: \code typedef Eigen::Triplet T; std::vector tripletList; triplets.reserve(estimation_of_entries); for(...) { // ... tripletList.push_back(T(i,j,v_ij)); } SparseMatrixType mat(rows,cols); mat.setFromTriplets(tripletList.begin(), tripletList.end()); // mat is ready to go! \endcode The \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(). See the SparseMatrix::setFromTriplets() function and class Triplet for more details. In some cases, however, slightly higher performance, and lower memory consumption can be reached by directly inserting the non-zeros into the destination matrix. A typical scenario of this approach is illustrated bellow: \code 1: SparseMatrix mat(rows,cols); // default is column major 2: mat.reserve(VectorXi::Constant(cols,6)); 3: for each i,j such that v_ij != 0 4: mat.insert(i,j) = v_ij; // alternative: mat.coeffRef(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-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). 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. - 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. - 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. - The line 5 suppresses the remaining empty space and transforms the matrix into a compressed column storage. \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:
ClassModuleSolver 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
PastixLLT \n PastixLDLT \n PastixLU\link PaStiXSupport_Module PaStiXSupport \endlinkDirect LLt, LDLt, LU factorizationsSPD \n SPD \n SquareFill-in reducing, Leverage fast dense algebra, Multithreading Requires the PaStiX package, \b CeCILL-C optimized for tough problems and symmetric patterns
CholmodSupernodalLLT\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, (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() method 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, and so on. More details are availble in the documentations of the respective classes. \section TutorialSparseFeatureSet Supported operators and functions Because of their special storage format, sparse matrices cannot offer the same level of flexbility than dense matrices. In Eigen's sparse module we chose to expose only the subset of the dense matrix API which can be efficiently implemented. 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. \subsection TutorialSparse_BasicOps Basic operations %Sparse expressions support most of the unary and binary coefficient wise operations: \code sm1.real() sm1.imag() -sm1 0.5*sm1 sm1+sm2 sm1-sm2 sm1.cwiseProduct(sm2) \endcode However, a strong restriction is that the storage orders must match. For instance, in the following example: \code sm4 = sm1 + sm2 + sm3; \endcode sm1, sm2, and sm3 must all be row-major or all column major. On the other hand, there is no restriction on the target matrix sm4. For 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: \code SparseMatrix A, B; B = SparseMatrix(A.transpose()) + A; \endcode Binary coefficient wise operators can also mix sparse and dense expressions: \code sm2 = sm1.cwiseProduct(dm1); dm2 = sm1 + dm1; \endcode %Sparse expressions also support transposition: \code sm1 = sm2.transpose(); sm1 = sm2.adjoint(); \endcode However, there is no transposeInPlace() method. \subsection TutorialSparse_Products Matrix products %Eigen supports various kind of sparse matrix products which are summarize below: - \b sparse-dense: \code dv2 = sm1 * dv1; dm2 = dm1 * sm1.adjoint(); dm2 = 2. * sm1 * dm1; \endcode - \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(): \code dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored dm2 = A.selfadjointView() * dm1; // if only the upper part of A is stored dm2 = A.selfadjointView() * dm1; // if only the lower part of A is stored \endcode - \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: \code sm3 = sm1 * sm2; sm3 = 4 * sm1.adjoint() * sm2; \endcode 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: \code sm3 = (sm1 * sm2).prune(); // removes numerical zeros sm3 = (sm1 * sm2).prune(ref); // removes elements much smaller than ref sm3 = (sm1 * sm2).prune(ref,epsilon); // removes elements smaller than ref*epsilon \endcode - \b permutations. Finally, permutations can be applied to sparse matrices too: \code PermutationMatrix P = ...; sm2 = P * sm1; sm2 = sm1 * P.inverse(); sm2 = sm1.transpose() * P; \endcode \subsection TutorialSparse_TriangularSelfadjoint Triangular and selfadjoint views Just 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: \code dm2 = sm1.triangularView(dm1); dv2 = sm1.transpose().triangularView(dv1); \endcode The selfadjointView() function permits various operations: - optimized sparse-dense matrix products: \code dm2 = sm1.selfadjointView<>() * dm1; // if all coefficients of A are stored dm2 = A.selfadjointView() * dm1; // if only the upper part of A is stored dm2 = A.selfadjointView() * dm1; // if only the lower part of A is stored \endcode - copy of triangular parts: \code sm2 = sm1.selfadjointView(); // makes a full selfadjoint matrix from the upper triangular part sm2.selfadjointView() = sm1.selfadjointView(); // copies the upper triangular part to the lower triangular part \endcode - application of symmetric permutations: \code PermutationMatrix P = ...; sm2 = A.selfadjointView().twistedBy(P); // compute P S P' from the upper triangular part of A, and make it a full matrix sm2.selfadjointView() = A.selfadjointView().twistedBy(P); // compute P S P' from the lower triangular part of A, and then only compute the lower part \endcode \subsection TutorialSparse_Submat Sub-matrices %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: \code sm1.innerVector(j); // returns an expression of the j-th column (resp. row) of the matrix if sm1 is col-major (resp. row-major) sm1.innerVectors(j, nb); // returns an expression of the nb columns (resp. row) starting from the j-th column (resp. row) // of the matrix if sm1 is col-major (resp. row-major) sm1.middleRows(j, nb); // for row major matrices only, get a range of nb rows sm1.middleCols(j, nb); // for column major matrices only, get a range of nb columns \endcode \li \b Next: \ref TutorialMapClass */ }