Path Tracer
Public Types | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
Eigen::SparseLU< _MatrixType, _OrderingType > Class Template Reference

Sparse supernodal LU factorization for general matrices. More...

#include <SparseLU.h>

Public Types

enum  { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef _MatrixType MatrixType
 
typedef _OrderingType OrderingType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef SparseMatrix< Scalar, ColMajor, StorageIndex > NCMatrix
 
typedef internal::MappedSuperNodalMatrix< Scalar, StorageIndex > SCMatrix
 
typedef Matrix< Scalar, Dynamic, 1 > ScalarVector
 
typedef Matrix< StorageIndex, Dynamic, 1 > IndexVector
 
typedef PermutationMatrix< Dynamic, Dynamic, StorageIndex > PermutationType
 
typedef internal::SparseLUImpl< Scalar, StorageIndex > Base
 

Public Member Functions

 SparseLU (const MatrixType &matrix)
 
void analyzePattern (const MatrixType &matrix)
 
void factorize (const MatrixType &matrix)
 
void simplicialfactorize (const MatrixType &matrix)
 
void compute (const MatrixType &matrix)
 
Index rows () const
 
Index cols () const
 
void isSymmetric (bool sym)
 
SparseLUMatrixLReturnType< SCMatrixmatrixL () const
 
SparseLUMatrixUReturnType< SCMatrix, MappedSparseMatrix< Scalar, ColMajor, StorageIndex > > matrixU () const
 
const PermutationTyperowsPermutation () const
 
const PermutationTypecolsPermutation () const
 
void setPivotThreshold (const RealScalar &thresh)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
std::string lastErrorMessage () const
 
template<typename Rhs , typename Dest >
bool _solve_impl (const MatrixBase< Rhs > &B, MatrixBase< Dest > &X_base) const
 
Scalar absDeterminant ()
 
Scalar logAbsDeterminant () const
 
Scalar signDeterminant ()
 
Scalar determinant ()
 
Index nnzL () const
 
Index nnzU () const
 
template<typename Rhs , typename Dest >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 

Protected Types

typedef SparseSolverBase< SparseLU< _MatrixType, _OrderingType > > APIBase
 

Protected Member Functions

void initperfvalues ()
 

Protected Attributes

ComputationInfo m_info
 
bool m_factorizationIsOk
 
bool m_analysisIsOk
 
std::string m_lastError
 
NCMatrix m_mat
 
SCMatrix m_Lstore
 
MappedSparseMatrix< Scalar, ColMajor, StorageIndex > m_Ustore
 
PermutationType m_perm_c
 
PermutationType m_perm_r
 
IndexVector m_etree
 
Base::GlobalLU_t m_glu
 
bool m_symmetricmode
 
internal::perfvalues m_perfv
 
RealScalar m_diagpivotthresh
 
Index m_nnzL
 
Index m_nnzU
 
Index m_detPermR
 
Index m_detPermC
 
bool m_isInitialized
 

Detailed Description

template<typename _MatrixType, typename _OrderingType>
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 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> A;
SparseLU<SparseMatrix<double>, COLAMDOrdering<int> > 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

\implsparsesolverconcept

See also
TutorialSparseSolverConcept
OrderingMethods_Module

Member Function Documentation

◆ absDeterminant()

template<typename _MatrixType , typename _OrderingType >
Scalar Eigen::SparseLU< _MatrixType, _OrderingType >::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()

◆ analyzePattern()

template<typename MatrixType , typename OrderingType >
void Eigen::SparseLU< MatrixType, OrderingType >::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

◆ colsPermutation()

template<typename _MatrixType , typename _OrderingType >
const PermutationType& Eigen::SparseLU< _MatrixType, _OrderingType >::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()

◆ compute()

template<typename _MatrixType , typename _OrderingType >
void Eigen::SparseLU< _MatrixType, _OrderingType >::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.

◆ determinant()

template<typename _MatrixType , typename _OrderingType >
Scalar Eigen::SparseLU< _MatrixType, _OrderingType >::determinant ( )
inline
Returns
The determinant of the matrix.
See also
absDeterminant(), logAbsDeterminant()

◆ factorize()

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.  

◆ info()

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

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
See also
iparm()

◆ isSymmetric()

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

Indicate that the pattern of the input matrix is symmetric

◆ lastErrorMessage()

template<typename _MatrixType , typename _OrderingType >
std::string Eigen::SparseLU< _MatrixType, _OrderingType >::lastErrorMessage ( ) const
inline
Returns
A string describing the type of error

◆ logAbsDeterminant()

template<typename _MatrixType , typename _OrderingType >
Scalar Eigen::SparseLU< _MatrixType, _OrderingType >::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()

◆ matrixL()

template<typename _MatrixType , typename _OrderingType >
SparseLUMatrixLReturnType<SCMatrix> Eigen::SparseLU< _MatrixType, _OrderingType >::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);

◆ matrixU()

template<typename _MatrixType , typename _OrderingType >
SparseLUMatrixUReturnType<SCMatrix,MappedSparseMatrix<Scalar,ColMajor,StorageIndex> > Eigen::SparseLU< _MatrixType, _OrderingType >::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);

◆ rowsPermutation()

template<typename _MatrixType , typename _OrderingType >
const PermutationType& Eigen::SparseLU< _MatrixType, _OrderingType >::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()

◆ setPivotThreshold()

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

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

◆ signDeterminant()

template<typename _MatrixType , typename _OrderingType >
Scalar Eigen::SparseLU< _MatrixType, _OrderingType >::signDeterminant ( )
inline
Returns
A number representing the sign of the determinant
See also
absDeterminant(), logAbsDeterminant()

The documentation for this class was generated from the following file:
Eigen::SparseLU::matrixL
SparseLUMatrixLReturnType< SCMatrix > matrixL() const
Definition: SparseLU.h:146
Eigen::SparseLU::matrixU
SparseLUMatrixUReturnType< SCMatrix, MappedSparseMatrix< Scalar, ColMajor, StorageIndex > > matrixU() const
Definition: SparseLU.h:156