Public Types | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
Eigen::SparseLU Class Reference

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

#include <SparseLU.h>

Public Types

enum  { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef internal::SparseLUImpl< Scalar, StorageIndexBase
 
typedef Matrix< StorageIndex, Dynamic, 1 > IndexVector
 
typedef _MatrixType MatrixType
 
typedef SparseMatrix< Scalar, ColMajor, StorageIndexNCMatrix
 
typedef _OrderingType OrderingType
 
typedef PermutationMatrix< Dynamic, Dynamic, StorageIndexPermutationType
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
typedef Matrix< Scalar, Dynamic, 1 > ScalarVector
 
typedef internal::MappedSuperNodalMatrix< Scalar, StorageIndexSCMatrix
 
typedef MatrixType::StorageIndex StorageIndex
 

Public Member Functions

template<typename Rhs , typename Dest >
bool _solve_impl (const MatrixBase< Rhs > &B, MatrixBase< Dest > &X_base) const
 
template<typename Rhs , typename Dest >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 
Scalar absDeterminant ()
 
const SparseLUTransposeView< true, SparseLU< _MatrixType, _OrderingType > > adjoint ()
 
void analyzePattern (const MatrixType &matrix)
 
Index cols () const
 
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< SCMatrixmatrixL () const
 
SparseLUMatrixUReturnType< SCMatrix, MappedSparseMatrix< Scalar, ColMajor, StorageIndex > > matrixU () const
 
Index nnzL () const
 
Index nnzU () const
 
Index rows () const
 
const PermutationTyperowsPermutation () const
 
void setPivotThreshold (const RealScalar &thresh)
 
Scalar signDeterminant ()
 
void simplicialfactorize (const MatrixType &matrix)
 
 SparseLU ()
 
 SparseLU (const MatrixType &matrix)
 
const SparseLUTransposeView< false, SparseLU< _MatrixType, _OrderingType > > transpose ()
 
 ~SparseLU ()
 

Protected Types

typedef SparseSolverBase< SparseLU< _MatrixType, _OrderingType > > APIBase
 

Protected Member Functions

void initperfvalues ()
 

Protected Attributes

bool m_analysisIsOk
 
Index m_detPermC
 
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 m_perfv
 
PermutationType m_perm_c
 
PermutationType m_perm_r
 
bool m_symmetricmode
 
MappedSparseMatrix< Scalar, ColMajor, StorageIndexm_Ustore
 

Private Member Functions

 SparseLU (const SparseLU &)
 

Detailed Description

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
Sparse solver concept
OrderingMethods_Module

Definition at line 17 of file SparseLU.h.

Member Typedef Documentation

◆ APIBase

typedef SparseSolverBase<SparseLU<_MatrixType,_OrderingType> > Eigen::SparseLU::APIBase
protected

Definition at line 134 of file SparseLU.h.

◆ Base

Definition at line 149 of file SparseLU.h.

◆ IndexVector

Definition at line 147 of file SparseLU.h.

◆ MatrixType

typedef _MatrixType Eigen::SparseLU::MatrixType

Definition at line 139 of file SparseLU.h.

◆ NCMatrix

Definition at line 144 of file SparseLU.h.

◆ OrderingType

typedef _OrderingType Eigen::SparseLU::OrderingType

Definition at line 140 of file SparseLU.h.

◆ PermutationType

Definition at line 148 of file SparseLU.h.

◆ RealScalar

typedef MatrixType::RealScalar Eigen::SparseLU::RealScalar

Definition at line 142 of file SparseLU.h.

◆ Scalar

typedef MatrixType::Scalar Eigen::SparseLU::Scalar

Definition at line 141 of file SparseLU.h.

◆ ScalarVector

Definition at line 146 of file SparseLU.h.

◆ SCMatrix

Definition at line 145 of file SparseLU.h.

◆ StorageIndex

typedef MatrixType::StorageIndex Eigen::SparseLU::StorageIndex

Definition at line 143 of file SparseLU.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
ColsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 151 of file SparseLU.h.

Constructor & Destructor Documentation

◆ SparseLU() [1/3]

Eigen::SparseLU::SparseLU ( )
inline

Definition at line 158 of file SparseLU.h.

◆ SparseLU() [2/3]

Eigen::SparseLU::SparseLU ( const MatrixType matrix)
inlineexplicit

Definition at line 162 of file SparseLU.h.

◆ ~SparseLU()

Eigen::SparseLU::~SparseLU ( )
inline

Definition at line 169 of file SparseLU.h.

◆ SparseLU() [3/3]

Eigen::SparseLU::SparseLU ( const SparseLU )
private

Member Function Documentation

◆ _solve_impl() [1/2]

template<typename Rhs , typename Dest >
bool Eigen::SparseLU::_solve_impl ( const MatrixBase< Rhs > &  B,
MatrixBase< Dest > &  X_base 
) const
inline

Definition at line 314 of file SparseLU.h.

◆ _solve_impl() [2/2]

template<typename Rhs , typename Dest >
void Eigen::SparseSolverBase< Derived >::_solve_impl ( typename Rhs  ,
typename Dest   
)
inline

Definition at line 111 of file SparseSolverBase.h.

◆ absDeterminant()

Scalar Eigen::SparseLU::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 350 of file SparseLU.h.

◆ adjoint()

const SparseLUTransposeView<true, SparseLU<_MatrixType,_OrderingType> > Eigen::SparseLU::adjoint ( )
inline
Returns
an expression of the adjoint of the factored matrix

A typical usage is to solve for the adjoint problem A' x = b:

solver.compute(A);
x = solver.adjoint().solve(b);

For real scalar types, this function is equivalent to transpose().

See also
transpose(), solve()

Definition at line 221 of file SparseLU.h.

◆ analyzePattern()

void Eigen::SparseLU::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 510 of file SparseLU.h.

◆ cols()

Index Eigen::SparseLU::cols ( ) const
inline

Definition at line 230 of file SparseLU.h.

◆ colsPermutation()

const PermutationType& Eigen::SparseLU::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 270 of file SparseLU.h.

◆ compute()

void Eigen::SparseLU::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 182 of file SparseLU.h.

◆ determinant()

Scalar Eigen::SparseLU::determinant ( )
inline
Returns
The determinant of the matrix.
See also
absDeterminant(), logAbsDeterminant()

Definition at line 434 of file SparseLU.h.

◆ factorize()

void Eigen::SparseLU::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 595 of file SparseLU.h.

◆ info()

ComputationInfo Eigen::SparseLU::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()

Definition at line 299 of file SparseLU.h.

◆ initperfvalues()

void Eigen::SparseLU::initperfvalues ( )
inlineprotected

Definition at line 460 of file SparseLU.h.

◆ isSymmetric()

void Eigen::SparseLU::isSymmetric ( bool  sym)
inline

Indicate that the pattern of the input matrix is symmetric

Definition at line 232 of file SparseLU.h.

◆ lastErrorMessage()

std::string Eigen::SparseLU::lastErrorMessage ( ) const
inline
Returns
A string describing the type of error

Definition at line 308 of file SparseLU.h.

◆ logAbsDeterminant()

Scalar Eigen::SparseLU::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 380 of file SparseLU.h.

◆ matrixL()

SparseLUMatrixLReturnType<SCMatrix> Eigen::SparseLU::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 243 of file SparseLU.h.

◆ matrixU()

SparseLUMatrixUReturnType<SCMatrix,MappedSparseMatrix<Scalar,ColMajor,StorageIndex> > Eigen::SparseLU::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 253 of file SparseLU.h.

◆ nnzL()

Index Eigen::SparseLU::nnzL ( ) const
inline

Definition at line 455 of file SparseLU.h.

◆ nnzU()

Index Eigen::SparseLU::nnzU ( ) const
inline

Definition at line 456 of file SparseLU.h.

◆ rows()

Index Eigen::SparseLU::rows ( ) const
inline

Definition at line 229 of file SparseLU.h.

◆ rowsPermutation()

const PermutationType& Eigen::SparseLU::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 262 of file SparseLU.h.

◆ setPivotThreshold()

void Eigen::SparseLU::setPivotThreshold ( const RealScalar thresh)
inline

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

Definition at line 275 of file SparseLU.h.

◆ signDeterminant()

Scalar Eigen::SparseLU::signDeterminant ( )
inline
Returns
A number representing the sign of the determinant
See also
absDeterminant(), logAbsDeterminant()

Definition at line 406 of file SparseLU.h.

◆ simplicialfactorize()

void Eigen::SparseLU::simplicialfactorize ( const MatrixType matrix)

◆ transpose()

const SparseLUTransposeView<false,SparseLU<_MatrixType,_OrderingType> > Eigen::SparseLU::transpose ( )
inline
Returns
an expression of the transposed of the factored matrix.

A typical usage is to solve for the transposed problem A^T x = b:

solver.compute(A);
x = solver.transpose().solve(b);
See also
adjoint(), solve()

Definition at line 200 of file SparseLU.h.

Member Data Documentation

◆ m_analysisIsOk

bool Eigen::SparseLU::m_analysisIsOk
protected

Definition at line 473 of file SparseLU.h.

◆ m_detPermC

Index Eigen::SparseLU::m_detPermC
protected

Definition at line 490 of file SparseLU.h.

◆ m_detPermR

Index Eigen::SparseLU::m_detPermR
protected

Definition at line 490 of file SparseLU.h.

◆ m_diagpivotthresh

RealScalar Eigen::SparseLU::m_diagpivotthresh
protected

Definition at line 488 of file SparseLU.h.

◆ m_etree

IndexVector Eigen::SparseLU::m_etree
protected

Definition at line 480 of file SparseLU.h.

◆ m_factorizationIsOk

bool Eigen::SparseLU::m_factorizationIsOk
protected

Definition at line 472 of file SparseLU.h.

◆ m_glu

Base::GlobalLU_t Eigen::SparseLU::m_glu
protected

Definition at line 482 of file SparseLU.h.

◆ m_info

ComputationInfo Eigen::SparseLU::m_info
mutableprotected

Definition at line 471 of file SparseLU.h.

◆ m_isInitialized

bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprotected

Definition at line 119 of file SparseSolverBase.h.

◆ m_lastError

std::string Eigen::SparseLU::m_lastError
protected

Definition at line 474 of file SparseLU.h.

◆ m_Lstore

SCMatrix Eigen::SparseLU::m_Lstore
protected

Definition at line 476 of file SparseLU.h.

◆ m_mat

NCMatrix Eigen::SparseLU::m_mat
protected

Definition at line 475 of file SparseLU.h.

◆ m_nnzL

Index Eigen::SparseLU::m_nnzL
protected

Definition at line 489 of file SparseLU.h.

◆ m_nnzU

Index Eigen::SparseLU::m_nnzU
protected

Definition at line 489 of file SparseLU.h.

◆ m_perfv

internal::perfvalues Eigen::SparseLU::m_perfv
protected

Definition at line 487 of file SparseLU.h.

◆ m_perm_c

PermutationType Eigen::SparseLU::m_perm_c
protected

Definition at line 478 of file SparseLU.h.

◆ m_perm_r

PermutationType Eigen::SparseLU::m_perm_r
protected

Definition at line 479 of file SparseLU.h.

◆ m_symmetricmode

bool Eigen::SparseLU::m_symmetricmode
protected

Definition at line 485 of file SparseLU.h.

◆ m_Ustore

MappedSparseMatrix<Scalar,ColMajor,StorageIndex> Eigen::SparseLU::m_Ustore
protected

Definition at line 477 of file SparseLU.h.


The documentation for this class was generated from the following file:
b
Scalar * b
Definition: benchVecAdd.cpp:17
x
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition: gnuplot_common_settings.hh:12
solver
BiCGSTAB< SparseMatrix< double > > solver
Definition: BiCGSTAB_simple.cpp:5
n
int n
Definition: BiCGSTAB_simple.cpp:1
Eigen::SparseLU::matrixU
SparseLUMatrixUReturnType< SCMatrix, MappedSparseMatrix< Scalar, ColMajor, StorageIndex > > matrixU() const
Definition: SparseLU.h:253
A
Definition: test_numpy_dtypes.cpp:298
y
Scalar * y
Definition: level1_cplx_impl.h:124
Eigen::SparseLU::matrixL
SparseLUMatrixLReturnType< SCMatrix > matrixL() const
Definition: SparseLU.h:243


gtsam
Author(s):
autogenerated on Fri Nov 1 2024 03:46:25