Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions
Eigen::SparseLU< _MatrixType, _OrderingType > Class Template Reference

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

#include <SparseLU.h>

Inheritance diagram for Eigen::SparseLU< _MatrixType, _OrderingType >:
Inheritance graph
[legend]

List of all members.

Public Types

typedef internal::SparseLUImpl
< Scalar, Index
Base
typedef MatrixType::Index Index
typedef Matrix< Index, Dynamic, 1 > IndexVector
typedef _MatrixType MatrixType
typedef SparseMatrix< Scalar,
ColMajor, Index
NCMatrix
typedef _OrderingType OrderingType
typedef PermutationMatrix
< Dynamic, Dynamic, Index
PermutationType
typedef MatrixType::RealScalar RealScalar
typedef MatrixType::Scalar Scalar
typedef Matrix< Scalar,
Dynamic, 1 > 
ScalarVector
typedef
internal::MappedSuperNodalMatrix
< Scalar, Index
SCMatrix

Public Member Functions

template<typename Rhs , typename Dest >
bool _solve (const MatrixBase< Rhs > &B, MatrixBase< Dest > &_X) const
Scalar absDeterminant ()
void analyzePattern (const MatrixType &matrix)
Index cols () const
const PermutationTypecolsPermutation () const
void compute (const MatrixType &matrix)
void factorize (const MatrixType &matrix)
ComputationInfo info () const
 Reports whether previous computation was successful.
void isSymmetric (bool sym)
std::string lastErrorMessage () const
Scalar logAbsDeterminant () const
SparseLUMatrixLReturnType
< SCMatrix
matrixL () const
SparseLUMatrixUReturnType
< SCMatrix, MappedSparseMatrix
< Scalar, ColMajor, Index > > 
matrixU () const
Index rows () const
const PermutationTyperowsPermutation () const
void setPivotThreshold (const RealScalar &thresh)
Scalar signDeterminant ()
void simplicialfactorize (const MatrixType &matrix)
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
 SparseLU ()
 SparseLU (const MatrixType &matrix)
 ~SparseLU ()

Protected Member Functions

void initperfvalues ()

Protected Attributes

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

Private Member Functions

 SparseLU (const SparseLU &)

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 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:
TutorialSparseDirectSolvers
OrderingMethods_Module

Definition at line 73 of file SparseLU.h.


Member Typedef Documentation

template<typename _MatrixType , typename _OrderingType >
typedef internal::SparseLUImpl<Scalar, Index> Eigen::SparseLU< _MatrixType, _OrderingType >::Base

Definition at line 86 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
typedef MatrixType::Index Eigen::SparseLU< _MatrixType, _OrderingType >::Index

Definition at line 80 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
typedef Matrix<Index,Dynamic,1> Eigen::SparseLU< _MatrixType, _OrderingType >::IndexVector
template<typename _MatrixType , typename _OrderingType >
typedef _MatrixType Eigen::SparseLU< _MatrixType, _OrderingType >::MatrixType
template<typename _MatrixType , typename _OrderingType >
typedef SparseMatrix<Scalar,ColMajor,Index> Eigen::SparseLU< _MatrixType, _OrderingType >::NCMatrix

Definition at line 81 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
typedef _OrderingType Eigen::SparseLU< _MatrixType, _OrderingType >::OrderingType

Definition at line 77 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
typedef PermutationMatrix<Dynamic, Dynamic, Index> Eigen::SparseLU< _MatrixType, _OrderingType >::PermutationType

Definition at line 85 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
typedef MatrixType::RealScalar Eigen::SparseLU< _MatrixType, _OrderingType >::RealScalar
template<typename _MatrixType , typename _OrderingType >
typedef MatrixType::Scalar Eigen::SparseLU< _MatrixType, _OrderingType >::Scalar

Definition at line 78 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
typedef Matrix<Scalar,Dynamic,1> Eigen::SparseLU< _MatrixType, _OrderingType >::ScalarVector
template<typename _MatrixType , typename _OrderingType >
typedef internal::MappedSuperNodalMatrix<Scalar, Index> Eigen::SparseLU< _MatrixType, _OrderingType >::SCMatrix

Definition at line 82 of file SparseLU.h.


Constructor & Destructor Documentation

template<typename _MatrixType , typename _OrderingType >
Eigen::SparseLU< _MatrixType, _OrderingType >::SparseLU ( ) [inline]

Definition at line 89 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
Eigen::SparseLU< _MatrixType, _OrderingType >::SparseLU ( const MatrixType matrix) [inline]

Definition at line 93 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
Eigen::SparseLU< _MatrixType, _OrderingType >::~SparseLU ( ) [inline]

Definition at line 99 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
Eigen::SparseLU< _MatrixType, _OrderingType >::SparseLU ( const SparseLU< _MatrixType, _OrderingType > &  ) [private]

Member Function Documentation

template<typename _MatrixType , typename _OrderingType >
template<typename Rhs , typename Dest >
bool Eigen::SparseLU< _MatrixType, _OrderingType >::_solve ( const MatrixBase< Rhs > &  B,
MatrixBase< Dest > &  _X 
) const [inline]

Definition at line 222 of file SparseLU.h.

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()

Definition at line 256 of file SparseLU.h.

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

Definition at line 369 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseLU< _MatrixType, _OrderingType >::cols ( void  ) const [inline]

Definition at line 121 of file SparseLU.h.

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()

Definition at line 161 of file SparseLU.h.

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.

Definition at line 112 of file SparseLU.h.

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.

Definition at line 452 of file SparseLU.h.

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

Definition at line 207 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
void Eigen::SparseLU< _MatrixType, _OrderingType >::initperfvalues ( ) [inline, protected]

Definition at line 317 of file SparseLU.h.

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

Indicate that the pattern of the input matrix is symmetric

Definition at line 123 of file SparseLU.h.

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

Definition at line 216 of file SparseLU.h.

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()

Definition at line 286 of file SparseLU.h.

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);

Definition at line 134 of file SparseLU.h.

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

Definition at line 144 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseLU< _MatrixType, _OrderingType >::rows ( void  ) const [inline]

Definition at line 120 of file SparseLU.h.

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()

Definition at line 153 of file SparseLU.h.

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.

Definition at line 166 of file SparseLU.h.

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()

Definition at line 309 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
void Eigen::SparseLU< _MatrixType, _OrderingType >::simplicialfactorize ( const MatrixType matrix)
template<typename _MatrixType , typename _OrderingType >
template<typename Rhs >
const internal::solve_retval<SparseLU, Rhs> Eigen::SparseLU< _MatrixType, _OrderingType >::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()

Definition at line 178 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
template<typename Rhs >
const internal::sparse_solve_retval<SparseLU, Rhs> Eigen::SparseLU< _MatrixType, _OrderingType >::solve ( const SparseMatrixBase< Rhs > &  B) const [inline]
Returns:
the solution X of $ A X = B $ using the current decomposition of A.
See also:
compute()

Definition at line 191 of file SparseLU.h.


Member Data Documentation

template<typename _MatrixType , typename _OrderingType >
bool Eigen::SparseLU< _MatrixType, _OrderingType >::m_analysisIsOk [protected]

Definition at line 331 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseLU< _MatrixType, _OrderingType >::m_detPermR [protected]

Definition at line 348 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
RealScalar Eigen::SparseLU< _MatrixType, _OrderingType >::m_diagpivotthresh [protected]

Definition at line 346 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
IndexVector Eigen::SparseLU< _MatrixType, _OrderingType >::m_etree [protected]

Definition at line 338 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
bool Eigen::SparseLU< _MatrixType, _OrderingType >::m_factorizationIsOk [protected]

Definition at line 330 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
Base::GlobalLU_t Eigen::SparseLU< _MatrixType, _OrderingType >::m_glu [protected]

Definition at line 340 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
ComputationInfo Eigen::SparseLU< _MatrixType, _OrderingType >::m_info [mutable, protected]

Definition at line 328 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
bool Eigen::SparseLU< _MatrixType, _OrderingType >::m_isInitialized [protected]

Definition at line 329 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
std::string Eigen::SparseLU< _MatrixType, _OrderingType >::m_lastError [protected]

Definition at line 332 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
SCMatrix Eigen::SparseLU< _MatrixType, _OrderingType >::m_Lstore [protected]

Definition at line 334 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
NCMatrix Eigen::SparseLU< _MatrixType, _OrderingType >::m_mat [protected]

Definition at line 333 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseLU< _MatrixType, _OrderingType >::m_nnzL [protected]

Definition at line 347 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
Index Eigen::SparseLU< _MatrixType, _OrderingType >::m_nnzU [protected]

Definition at line 347 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
internal::perfvalues<Index> Eigen::SparseLU< _MatrixType, _OrderingType >::m_perfv [protected]

Definition at line 345 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
PermutationType Eigen::SparseLU< _MatrixType, _OrderingType >::m_perm_c [protected]

Definition at line 336 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
PermutationType Eigen::SparseLU< _MatrixType, _OrderingType >::m_perm_r [protected]

Definition at line 337 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
bool Eigen::SparseLU< _MatrixType, _OrderingType >::m_symmetricmode [protected]

Definition at line 343 of file SparseLU.h.

template<typename _MatrixType , typename _OrderingType >
MappedSparseMatrix<Scalar,ColMajor,Index> Eigen::SparseLU< _MatrixType, _OrderingType >::m_Ustore [protected]

Definition at line 335 of file SparseLU.h.


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


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Thu Aug 27 2015 12:02:03