template<typename _MatrixType, typename _OrderingType>
Eigen::SparseLU class

Sparse supernodal LU factorization for general matrices.

Template parameters
_MatrixType The type of the sparse matrix. It must be a column-major SparseMatrix<>
_OrderingType The ordering method to use, either AMD, COLAMD or METIS. Default is COLMAD

Contents

This class implements the supernodal LU factorization for general matrices. It uses the main techniques from the sequential SuperLU package (http://crd-legacy.lbl.gov/~xiaoye/SuperLU/). It handles transparently real and complex arithmetic with single and double precision, depending on the scalar type of your input matrix. The code has been optimized to provide BLAS-3 operations during supernode-panel updates. It benefits directly from the built-in high-performant Eigen BLAS routines. Moreover, when the size of a supernode is very small, the BLAS calls are avoided to enable a better optimization from the compiler. For best performance, you should compile it with NDEBUG flag to avoid the numerous bounds checking on vectors.

An important parameter of this class is the ordering method. It is used to reorder the columns (and eventually the rows) of the matrix to reduce the number of new elements that are created during numerical factorization. The cheapest method available is COLAMD. See the OrderingMethods module for the list of built-in and external ordering methods.

Simple example with key steps

VectorXd x(n), b(n);
SparseMatrix<double, ColMajor> A;
SparseLU<SparseMatrix<scalar, ColMajor>, COLAMDOrdering<Index> >   solver;
// fill A and b;
// Compute the ordering permutation vector from the structural pattern of A
solver.analyzePattern(A); 
// Compute the numerical factorization 
solver.factorize(A); 
//Use the factors to solve the linear system 
x = solver.solve(b);

This class follows the sparse solver concept.

Base classes

template<typename Derived>
class SparseSolverBase
A base class for sparse solvers.

Public functions

auto absDeterminant() -> Scalar
void analyzePattern(const MatrixType& matrix)
auto colsPermutation() const -> const PermutationType&
void compute(const MatrixType& matrix)
auto determinant() -> Scalar
void factorize(const MatrixType& matrix)
auto info() const -> ComputationInfo
Reports whether previous computation was successful.
void isSymmetric(bool sym)
auto lastErrorMessage() const -> std::string
auto logAbsDeterminant() const -> Scalar
auto matrixL() const -> SparseLUMatrixLReturnType<SCMatrix>
auto matrixU() const -> SparseLUMatrixUReturnType<SCMatrix, MappedSparseMatrix<Scalar, ColMajor, StorageIndex>>
auto rowsPermutation() const -> const PermutationType&
void setPivotThreshold(const RealScalar& thresh)
auto signDeterminant() -> Scalar
template<typename Rhs>
auto solve(const MatrixBase<Rhs>& B) const -> const Solve<SparseLU, Rhs>

Function documentation

template<typename _MatrixType, typename _OrderingType>
Scalar Eigen::SparseLU<_MatrixType, _OrderingType>::absDeterminant()

Returns the absolute value of the determinant of the matrix of which *this is the QR decomposition.

template<typename _MatrixType, typename _OrderingType>
void Eigen::SparseLU<_MatrixType, _OrderingType>::analyzePattern(const MatrixType& matrix)

Compute the column permutation to minimize the fill-in

  • Apply this permutation to the input matrix -
  • Compute the column elimination tree on the permuted matrix
  • Postorder the elimination tree and the column permutation

template<typename _MatrixType, typename _OrderingType>
const PermutationType& Eigen::SparseLU<_MatrixType, _OrderingType>::colsPermutation() const

Returns a reference to the column matrix permutation $ P_c^T $ such that $P_r A P_c^T = L U$

template<typename _MatrixType, typename _OrderingType>
void Eigen::SparseLU<_MatrixType, _OrderingType>::compute(const MatrixType& matrix)

Compute the symbolic and numeric factorization of the input sparse matrix. The input matrix should be in column-major storage.

template<typename _MatrixType, typename _OrderingType>
Scalar Eigen::SparseLU<_MatrixType, _OrderingType>::determinant()

Returns The determinant of the matrix.

template<typename _MatrixType, typename _OrderingType>
void Eigen::SparseLU<_MatrixType, _OrderingType>::factorize(const MatrixType& matrix)

  • Numerical factorization
  • Interleaved with the symbolic factorization On exit, info is

    = 0: successful factorization

0: if info = i, and i is

  <= A->ncol: U(i,i) is exactly zero. The factorization has
     been completed, but the factor U is exactly singular,
     and division by zero will occur if it is used to solve a
     system of equations.

  > A->ncol: number of bytes allocated when memory allocation
    failure occurred, plus A->ncol. If lwork = -1, it is
    the estimated amount of space needed, plus A->ncol.  

template<typename _MatrixType, typename _OrderingType>
ComputationInfo Eigen::SparseLU<_MatrixType, _OrderingType>::info() const

Reports whether previous computation was successful.

Returns Success if computation was successful, NumericalIssue if the LU factorization reports a problem, zero diagonal for instance InvalidInput if the input matrix is invalid

template<typename _MatrixType, typename _OrderingType>
void Eigen::SparseLU<_MatrixType, _OrderingType>::isSymmetric(bool sym)

Indicate that the pattern of the input matrix is symmetric

template<typename _MatrixType, typename _OrderingType>
std::string Eigen::SparseLU<_MatrixType, _OrderingType>::lastErrorMessage() const

Returns A string describing the type of error

template<typename _MatrixType, typename _OrderingType>
Scalar Eigen::SparseLU<_MatrixType, _OrderingType>::logAbsDeterminant() const

Returns the natural log of the absolute value of the determinant of the matrix of which **this is the QR decomposition

template<typename _MatrixType, typename _OrderingType>
SparseLUMatrixLReturnType<SCMatrix> Eigen::SparseLU<_MatrixType, _OrderingType>::matrixL() const

Returns an expression of the matrix L, internally stored as supernodes The only operation available with this expression is the triangular solve y = b; matrixL().solveInPlace(y);

template<typename _MatrixType, typename _OrderingType>
SparseLUMatrixUReturnType<SCMatrix, MappedSparseMatrix<Scalar, ColMajor, StorageIndex>> Eigen::SparseLU<_MatrixType, _OrderingType>::matrixU() const

Returns an expression of the matrix U, The only operation available with this expression is the triangular solve y = b; matrixU().solveInPlace(y);

template<typename _MatrixType, typename _OrderingType>
const PermutationType& Eigen::SparseLU<_MatrixType, _OrderingType>::rowsPermutation() const

Returns a reference to the row matrix permutation $ P_r $ such that $P_r A P_c^T = L U$

template<typename _MatrixType, typename _OrderingType>
void Eigen::SparseLU<_MatrixType, _OrderingType>::setPivotThreshold(const RealScalar& thresh)

Set the threshold used for a diagonal entry to be an acceptable pivot.

template<typename _MatrixType, typename _OrderingType>
Scalar Eigen::SparseLU<_MatrixType, _OrderingType>::signDeterminant()

Returns A number representing the sign of the determinant

template<typename _MatrixType, typename _OrderingType> template<typename Rhs>
const Solve<SparseLU, Rhs> Eigen::SparseLU<_MatrixType, _OrderingType>::solve(const MatrixBase<Rhs>& B) const

Returns the solution X of $ A X = B $ using the current decomposition of A.