Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
SparseLU< _MatrixType, _OrderingType > Class Template Reference

Detailed Description

template<typename _MatrixType, typename _OrderingType = COLAMDOrdering<typename _MatrixType::Index>>
class Eigen::SparseLU< _MatrixType, _OrderingType >

Sparse supernodal LU factorization for general matrices.

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 arithmetics 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);
*
Warning
The input matrix A should be in a compressed and column-major form. Otherwise an expensive copy will be made. You can call the inexpensive makeCompressed() to get a compressed matrix.
Note
Unlike the initial SuperLU implementation, there is no step to equilibrate the matrix. For badly scaled matrices, this step can be useful to reduce the pivoting during factorization. If this is the case for your matrices, you can try the basic scaling method at "unsupported/Eigen/src/IterativeSolvers/Scaling.h"
Template Parameters
_MatrixTypeThe type of the sparse matrix. It must be a column-major SparseMatrix<>
_OrderingTypeThe ordering method to use, either AMD, COLAMD or METIS. Default is COLMAD
See Also
Sparse solvers
OrderingMethods module

Inherits SparseLUImpl< _MatrixType::Scalar, _MatrixType::Index >.

Public Member Functions

Scalar absDeterminant ()
 
void analyzePattern (const MatrixType &matrix)
 
const PermutationTypecolsPermutation () const
 
void compute (const MatrixType &matrix)
 
Scalar determinant ()
 
void factorize (const MatrixType &matrix)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
void isSymmetric (bool sym)
 
std::string lastErrorMessage () const
 
Scalar logAbsDeterminant () const
 
SparseLUMatrixLReturnType
< SCMatrix > 
matrixL () const
 
SparseLUMatrixUReturnType
< SCMatrix, MappedSparseMatrix
< Scalar, ColMajor, Index > > 
matrixU () const
 
const PermutationTyperowsPermutation () const
 
void setPivotThreshold (const RealScalar &thresh)
 
Scalar signDeterminant ()
 
template<typename Rhs >
const internal::solve_retval
< SparseLU, Rhs > 
solve (const MatrixBase< Rhs > &B) const
 
template<typename Rhs >
const
internal::sparse_solve_retval
< SparseLU, Rhs > 
solve (const SparseMatrixBase< Rhs > &B) const
 

Member Function Documentation

Scalar absDeterminant ( )
inline
Returns
the absolute value of the determinant of the matrix of which *this is the QR decomposition.
Warning
a determinant can be very big or small, so for matrices of large enough dimension, there is a risk of overflow/underflow. One way to work around that is to use logAbsDeterminant() instead.
See Also
logAbsDeterminant(), signDeterminant()
void analyzePattern ( const MatrixType &  mat)

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

References Eigen::internal::coletree(), PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, IndexType >::indices(), PlainObjectBase< Derived >::resize(), and Eigen::internal::treePostorder().

Referenced by SparseLU< _MatrixType, _OrderingType >::compute().

const PermutationType& colsPermutation ( ) const
inline
Returns
a reference to the column matrix permutation $ P_c^T $ such that $P_r A P_c^T = L U$
See Also
rowsPermutation()
void compute ( const MatrixType &  matrix)
inline

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

References SparseLU< _MatrixType, _OrderingType >::analyzePattern(), and SparseLU< _MatrixType, _OrderingType >::factorize().

Scalar determinant ( )
inline
Returns
The determinant of the matrix.
See Also
absDeterminant(), logAbsDeterminant()
void 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.  

References Eigen::NumericalIssue, PlainObjectBase< Derived >::setConstant(), PlainObjectBase< Derived >::setZero(), and Eigen::Success.

Referenced by SparseLU< _MatrixType, _OrderingType >::compute().

ComputationInfo info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was succesful, NumericalIssue if the LU factorization reports a problem, zero diagonal for instance InvalidInput if the input matrix is invalid
See Also
iparm()
void isSymmetric ( bool  sym)
inline

Indicate that the pattern of the input matrix is symmetric

std::string lastErrorMessage ( ) const
inline
Returns
A string describing the type of error
Scalar logAbsDeterminant ( ) const
inline
Returns
the natural log of the absolute value of the determinant of the matrix of which **this is the QR decomposition
Note
This method is useful to work around the risk of overflow/underflow that's inherent to the determinant computation.
See Also
absDeterminant(), signDeterminant()
SparseLUMatrixLReturnType<SCMatrix> matrixL ( ) const
inline
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);
*
SparseLUMatrixUReturnType<SCMatrix,MappedSparseMatrix<Scalar,ColMajor,Index> > matrixU ( ) const
inline
Returns
an expression of the matrix U, The only operation available with this expression is the triangular solve
* y = b; matrixU().solveInPlace(y);
*
const PermutationType& rowsPermutation ( ) const
inline
Returns
a reference to the row matrix permutation $ P_r $ such that $P_r A P_c^T = L U$
See Also
colsPermutation()
void setPivotThreshold ( const RealScalar &  thresh)
inline

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

Scalar signDeterminant ( )
inline
Returns
A number representing the sign of the determinant
See Also
absDeterminant(), logAbsDeterminant()
const internal::solve_retval<SparseLU, Rhs> solve ( const MatrixBase< Rhs > &  B) const
inline
Returns
the solution X of $ A X = B $ using the current decomposition of A.
Warning
the destination matrix X in X = this->solve(B) must be colmun-major.
See Also
compute()
const internal::sparse_solve_retval<SparseLU, Rhs> solve ( const SparseMatrixBase< Rhs > &  B) const
inline
Returns
the solution X of $ A X = B $ using the current decomposition of A.
See Also
compute()

References EigenBase< Derived >::derived(), and SparseMatrixBase< Derived >::rows().


The documentation for this class was generated from the following file: