Eigen::SuperLU< MatrixType_ > Class Template Reference

A sparse direct LU factorization and solver based on the SuperLU library. More...

#include <SuperLUSupport.h>

+ Inheritance diagram for Eigen::SuperLU< MatrixType_ >:

Public Types

typedef SuperLUBase< MatrixType_, SuperLUBase
 
typedef MatrixType_ MatrixType
 
typedef Base::Scalar Scalar
 
typedef Base::RealScalar RealScalar
 
typedef Base::StorageIndex StorageIndex
 
typedef Base::IntRowVectorType IntRowVectorType
 
typedef Base::IntColVectorType IntColVectorType
 
typedef Base::PermutationMap PermutationMap
 
typedef Base::LUMatrixType LUMatrixType
 
typedef TriangularView< LUMatrixType, Lower|UnitDiagLMatrixType
 
typedef TriangularView< LUMatrixType, UpperUMatrixType
 
- Public Types inherited from Eigen::SuperLUBase< MatrixType_, SuperLU< MatrixType_ > >
enum  
 
typedef MatrixType_ MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef Matrix< Scalar, Dynamic, 1 > Vector
 
typedef Matrix< int, 1, MatrixType::ColsAtCompileTime > IntRowVectorType
 
typedef Matrix< int, MatrixType::RowsAtCompileTime, 1 > IntColVectorType
 
typedef Map< PermutationMatrix< Dynamic, Dynamic, int > > PermutationMap
 
typedef SparseMatrix< ScalarLUMatrixType
 

Public Member Functions

 SuperLU ()
 
 SuperLU (const MatrixType &matrix)
 
 ~SuperLU ()
 
void analyzePattern (const MatrixType &matrix)
 
void factorize (const MatrixType &matrix)
 
template<typename Rhs , typename Dest >
void _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
 
const LMatrixTypematrixL () const
 
const UMatrixTypematrixU () const
 
const IntColVectorTypepermutationP () const
 
const IntRowVectorTypepermutationQ () const
 
Scalar determinant () const
 
- Public Member Functions inherited from Eigen::SuperLUBase< MatrixType_, SuperLU< MatrixType_ > >
 SuperLUBase ()
 
 ~SuperLUBase ()
 
Index rows () const
 
Index cols () const
 
superlu_options_toptions ()
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
void compute (const MatrixType &matrix)
 
void analyzePattern (const MatrixType &)
 
void dumpMemory (Stream &)
 
- Public Member Functions inherited from Eigen::SparseSolverBase< Derived >
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
template<typename Rhs , typename Dest >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 

Protected Member Functions

void init ()
 
- Protected Member Functions inherited from Eigen::SuperLUBase< MatrixType_, SuperLU< MatrixType_ > >
void initFactorization (const MatrixType &a)
 
void init ()
 
void extractData () const
 
void clearFactors ()
 
SuperLU< MatrixType_ > & derived ()
 
const SuperLU< MatrixType_ > & derived () const
 

Protected Attributes

LUMatrixType m_l
 
LUMatrixType m_matrix
 
IntColVectorType m_p
 
IntRowVectorType m_q
 
SluMatrix m_sluA
 
SluMatrix m_sluB
 
Matrix< RealScalar, Dynamic, 1 > m_sluBerr
 
Matrix< RealScalar, Dynamic, 1 > m_sluCscale
 
char m_sluEqued
 
std::vector< intm_sluEtree
 
Matrix< RealScalar, Dynamic, 1 > m_sluFerr
 
SuperMatrix m_sluL
 
superlu_options_t m_sluOptions
 
Matrix< RealScalar, Dynamic, 1 > m_sluRscale
 
SuperLUStat_t m_sluStat
 
SuperMatrix m_sluU
 
SluMatrix m_sluX
 
LUMatrixType m_u
 
int m_analysisIsOk
 
bool m_extractedDataAreDirty
 
int m_factorizationIsOk
 
ComputationInfo m_info
 
bool m_isInitialized
 
- Protected Attributes inherited from Eigen::SuperLUBase< MatrixType_, SuperLU< MatrixType_ > >
LUMatrixType m_l
 
LUMatrixType m_u
 
IntColVectorType m_p
 
IntRowVectorType m_q
 
LUMatrixType m_matrix
 
SluMatrix m_sluA
 
SuperMatrix m_sluL
 
SuperMatrix m_sluU
 
SluMatrix m_sluB
 
SluMatrix m_sluX
 
SuperLUStat_t m_sluStat
 
superlu_options_t m_sluOptions
 
std::vector< intm_sluEtree
 
Matrix< RealScalar, Dynamic, 1 > m_sluRscale
 
Matrix< RealScalar, Dynamic, 1 > m_sluCscale
 
Matrix< RealScalar, Dynamic, 1 > m_sluFerr
 
Matrix< RealScalar, Dynamic, 1 > m_sluBerr
 
char m_sluEqued
 
ComputationInfo m_info
 
int m_factorizationIsOk
 
int m_analysisIsOk
 
bool m_extractedDataAreDirty
 
bool m_isInitialized
 
- Protected Attributes inherited from Eigen::SparseSolverBase< Derived >
bool m_isInitialized
 

Private Member Functions

 SuperLU (SuperLU &)
 

Additional Inherited Members

- Protected Types inherited from Eigen::SuperLUBase< MatrixType_, SuperLU< MatrixType_ > >
typedef SparseSolverBase< SuperLU< MatrixType_ > > Base
 

Detailed Description

template<typename MatrixType_>
class Eigen::SuperLU< MatrixType_ >

A sparse direct LU factorization and solver based on the SuperLU library.

This class allows to solve for A.X = B sparse linear problems via a direct LU factorization using the SuperLU library. The sparse matrix A must be squared and invertible. The vectors or matrices X and B can be either dense or sparse.

Template Parameters
MatrixType_the type of the sparse matrix A, it must be a SparseMatrix<>
Warning
This class is only for the 4.x versions of SuperLU. The 3.x and 5.x versions are not supported.

\implsparsesolverconcept

See also
Sparse solver concept, class SparseLU

Member Typedef Documentation

◆ Base

template<typename MatrixType_ >
typedef SuperLUBase<MatrixType_, SuperLU> Eigen::SuperLU< MatrixType_ >::Base

◆ IntColVectorType

template<typename MatrixType_ >
typedef Base::IntColVectorType Eigen::SuperLU< MatrixType_ >::IntColVectorType

◆ IntRowVectorType

template<typename MatrixType_ >
typedef Base::IntRowVectorType Eigen::SuperLU< MatrixType_ >::IntRowVectorType

◆ LMatrixType

template<typename MatrixType_ >
typedef TriangularView<LUMatrixType, Lower | UnitDiag> Eigen::SuperLU< MatrixType_ >::LMatrixType

◆ LUMatrixType

template<typename MatrixType_ >
typedef Base::LUMatrixType Eigen::SuperLU< MatrixType_ >::LUMatrixType

◆ MatrixType

template<typename MatrixType_ >
typedef MatrixType_ Eigen::SuperLU< MatrixType_ >::MatrixType

◆ PermutationMap

template<typename MatrixType_ >
typedef Base::PermutationMap Eigen::SuperLU< MatrixType_ >::PermutationMap

◆ RealScalar

template<typename MatrixType_ >
typedef Base::RealScalar Eigen::SuperLU< MatrixType_ >::RealScalar

◆ Scalar

template<typename MatrixType_ >
typedef Base::Scalar Eigen::SuperLU< MatrixType_ >::Scalar

◆ StorageIndex

template<typename MatrixType_ >
typedef Base::StorageIndex Eigen::SuperLU< MatrixType_ >::StorageIndex

◆ UMatrixType

template<typename MatrixType_ >
typedef TriangularView<LUMatrixType, Upper> Eigen::SuperLU< MatrixType_ >::UMatrixType

Constructor & Destructor Documentation

◆ SuperLU() [1/3]

template<typename MatrixType_ >
Eigen::SuperLU< MatrixType_ >::SuperLU ( )
inline
448 : Base() { init(); }
void init()
Definition: SuperLUSupport.h:529
SuperLUBase< MatrixType_, SuperLU > Base
Definition: SuperLUSupport.h:433

References Eigen::SuperLU< MatrixType_ >::init().

◆ SuperLU() [2/3]

template<typename MatrixType_ >
Eigen::SuperLU< MatrixType_ >::SuperLU ( const MatrixType matrix)
inlineexplicit
450  : Base() {
451  init();
453  }
void compute(const MatrixType &matrix)
Definition: SuperLUSupport.h:318
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85

References Eigen::SuperLUBase< MatrixType_, Derived >::compute(), Eigen::SuperLU< MatrixType_ >::init(), and matrix().

◆ ~SuperLU()

template<typename MatrixType_ >
Eigen::SuperLU< MatrixType_ >::~SuperLU ( )
inline
455 {}

◆ SuperLU() [3/3]

template<typename MatrixType_ >
Eigen::SuperLU< MatrixType_ >::SuperLU ( SuperLU< MatrixType_ > &  )
inlineprivate
540 {}

Member Function Documentation

◆ _solve_impl()

template<typename MatrixType >
template<typename Rhs , typename Dest >
void Eigen::SuperLU< MatrixType >::_solve_impl ( const MatrixBase< Rhs > &  b,
MatrixBase< Dest > &  dest 
) const
573  {
575  "The decomposition is not in a valid state for solving, you must first call either compute() or "
576  "analyzePattern()/factorize()");
577 
578  const Index rhsCols = b.cols();
579  eigen_assert(m_matrix.rows() == b.rows());
580 
584 
585  m_sluFerr.resize(rhsCols);
586  m_sluBerr.resize(rhsCols);
587 
588  Ref<const Matrix<typename Rhs::Scalar, Dynamic, Dynamic, ColMajor> > b_ref(b);
589  Ref<const Matrix<typename Dest::Scalar, Dynamic, Dynamic, ColMajor> > x_ref(x);
590 
591  m_sluB = SluMatrix::Map(b_ref.const_cast_derived());
592  m_sluX = SluMatrix::Map(x_ref.const_cast_derived());
593 
594  typename Rhs::PlainObject b_cpy;
595  if (m_sluEqued != 'N') {
596  b_cpy = b;
597  m_sluB = SluMatrix::Map(b_cpy.const_cast_derived());
598  }
599 
601  int info = 0;
602  RealScalar recip_pivot_growth, rcond;
603  SuperLU_gssvx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0], &m_sluEqued, &m_sluRscale[0],
604  &m_sluCscale[0], &m_sluL, &m_sluU, NULL, 0, &m_sluB, &m_sluX, &recip_pivot_growth, &rcond,
605  &m_sluFerr[0], &m_sluBerr[0], &m_sluStat, &info, Scalar());
607 
608  if (x.derived().data() != x_ref.data()) x = x_ref;
609 
610  m_info = info == 0 ? Success : NumericalIssue;
611 }
#define eigen_assert(x)
Definition: Macros.h:910
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
constexpr EIGEN_DEVICE_FUNC const Scalar * data() const
Definition: PlainObjectBase.h:273
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:294
Index rows() const
Definition: SparseMatrix.h:159
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SuperLUSupport.h:312
Matrix< RealScalar, Dynamic, 1 > m_sluRscale
Definition: SuperLUSupport.h:401
std::vector< int > m_sluEtree
Definition: SuperLUSupport.h:400
SuperMatrix m_sluL
Definition: SuperLUSupport.h:396
SuperLUStat_t m_sluStat
Definition: SuperLUSupport.h:398
Matrix< RealScalar, Dynamic, 1 > m_sluBerr
Definition: SuperLUSupport.h:402
SluMatrix m_sluB
Definition: SuperLUSupport.h:397
SuperMatrix m_sluU
Definition: SuperLUSupport.h:396
IntRowVectorType m_q
Definition: SuperLUSupport.h:392
SluMatrix m_sluX
Definition: SuperLUSupport.h:397
Matrix< RealScalar, Dynamic, 1 > m_sluFerr
Definition: SuperLUSupport.h:402
Matrix< RealScalar, Dynamic, 1 > m_sluCscale
Definition: SuperLUSupport.h:401
int m_factorizationIsOk
Definition: SuperLUSupport.h:406
IntColVectorType m_p
Definition: SuperLUSupport.h:391
superlu_options_t m_sluOptions
Definition: SuperLUSupport.h:399
char m_sluEqued
Definition: SuperLUSupport.h:403
ComputationInfo m_info
Definition: SuperLUSupport.h:405
SluMatrix m_sluA
Definition: SuperLUSupport.h:395
LUMatrixType m_matrix
Definition: SuperLUSupport.h:394
@ NumericalIssue
Definition: Constants.h:442
@ Success
Definition: Constants.h:440
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
list x
Definition: plotDoE.py:28
@ NOTRANS
Definition: oomph_superlu_4.3/superlu_enum_consts.h:21
@ NOREFINE
Definition: oomph_superlu_4.3/superlu_enum_consts.h:23
@ FACTORED
Definition: oomph_superlu_4.3/superlu_enum_consts.h:17
void StatInit(SuperLUStat_t *)
void StatFree(SuperLUStat_t *)
static SluMatrix Map(MatrixBase< MatrixType > &_mat)
Definition: SuperLUSupport.h:152
trans_t Trans
Definition: slu_util.h:250
IterRefine_t IterRefine
Definition: slu_util.h:251
fact_t Fact
Definition: slu_util.h:247

References b, eigen_assert, FACTORED, info, Eigen::SluMatrix::Map(), NOREFINE, NOTRANS, Eigen::NumericalIssue, StatFree(), StatInit(), Eigen::Success, and plotDoE::x.

◆ analyzePattern()

template<typename MatrixType_ >
void Eigen::SuperLU< MatrixType_ >::analyzePattern ( const MatrixType matrix)
inline

Performs a symbolic decomposition on the sparcity of matrix.

This function is particularly useful when solving for several problems having the same structure.

See also
factorize()
463  {
465  m_isInitialized = false;
467  }
void analyzePattern(const MatrixType &)
Definition: SuperLUSupport.h:329
bool m_isInitialized
Definition: SparseSolverBase.h:110
@ InvalidInput
Definition: Constants.h:447

References Eigen::SuperLUBase< MatrixType_, Derived >::analyzePattern(), Eigen::InvalidInput, Eigen::SuperLU< MatrixType_ >::m_info, Eigen::SuperLU< MatrixType_ >::m_isInitialized, and matrix().

◆ determinant()

template<typename MatrixType >
SuperLU< MatrixType >::Scalar Eigen::SuperLU< MatrixType >::determinant
699  {
701  "The decomposition is not in a valid state for computing the determinant, you must first call either "
702  "compute() or analyzePattern()/factorize()");
703 
705 
706  Scalar det = Scalar(1);
707  for (int j = 0; j < m_u.cols(); ++j) {
708  if (m_u.outerIndexPtr()[j + 1] - m_u.outerIndexPtr()[j] > 0) {
709  int lastId = m_u.outerIndexPtr()[j + 1] - 1;
710  eigen_assert(m_u.innerIndexPtr()[lastId] <= j);
711  if (m_u.innerIndexPtr()[lastId] == j) det *= m_u.valuePtr()[lastId];
712  }
713  }
714  if (PermutationMap(m_p.data(), m_p.size()).determinant() * PermutationMap(m_q.data(), m_q.size()).determinant() < 0)
715  det = -det;
716  if (m_sluEqued != 'N')
717  return det / m_sluRscale.prod() / m_sluCscale.prod();
718  else
719  return det;
720 }
Index cols() const
Definition: SparseMatrix.h:161
const Scalar * valuePtr() const
Definition: SparseMatrix.h:171
const StorageIndex * outerIndexPtr() const
Definition: SparseMatrix.h:189
const StorageIndex * innerIndexPtr() const
Definition: SparseMatrix.h:180
void extractData() const
Definition: SuperLUSupport.h:621
bool m_extractedDataAreDirty
Definition: SuperLUSupport.h:408
Base::Scalar Scalar
Definition: SuperLUSupport.h:435
LUMatrixType m_u
Definition: SuperLUSupport.h:390
Base::PermutationMap PermutationMap
Definition: SuperLUSupport.h:440
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References eigen_assert, and j.

◆ factorize()

template<typename MatrixType >
void Eigen::SuperLU< MatrixType >::factorize ( const MatrixType matrix)

Performs a numeric decomposition of matrix

The given matrix must has the same sparcity than the matrix on which the symbolic decomposition has been performed.

See also
analyzePattern()
544  {
545  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
546  if (!m_analysisIsOk) {
548  return;
549  }
550 
551  this->initFactorization(a);
552 
554  int info = 0;
555  RealScalar recip_pivot_growth, rcond;
556  RealScalar ferr, berr;
557 
559  SuperLU_gssvx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0], &m_sluEqued, &m_sluRscale[0],
560  &m_sluCscale[0], &m_sluL, &m_sluU, NULL, 0, &m_sluB, &m_sluX, &recip_pivot_growth, &rcond, &ferr, &berr,
561  &m_sluStat, &info, Scalar());
563 
565 
566  // FIXME how to better check for errors ???
567  m_info = info == 0 ? Success : NumericalIssue;
568  m_factorizationIsOk = true;
569 }
void initFactorization(const MatrixType &a)
Definition: SuperLUSupport.h:340
int m_analysisIsOk
Definition: SuperLUSupport.h:407
const Scalar * a
Definition: level2_cplx_impl.h:32
@ COLAMD
Definition: oomph_superlu_4.3/superlu_enum_consts.h:19
colperm_t ColPerm
Definition: slu_util.h:249

References a, COLAMD, eigen_assert, info, Eigen::InvalidInput, Eigen::NumericalIssue, StatFree(), StatInit(), and Eigen::Success.

◆ init()

template<typename MatrixType_ >
void Eigen::SuperLU< MatrixType_ >::init ( )
inlineprotected
529  {
530  Base::init();
531 
537  }
void init()
Definition: SuperLUSupport.h:368
@ NO
Definition: oomph_superlu_4.3/superlu_enum_consts.h:16
void set_default_options(superlu_options_t *options)
yes_no_t ConditionNumber
Definition: slu_util.h:255
yes_no_t PrintStat
Definition: slu_util.h:268

References COLAMD, superlu_options_t::ColPerm, superlu_options_t::ConditionNumber, Eigen::SuperLUBase< MatrixType_, Derived >::init(), Eigen::SuperLU< MatrixType_ >::m_sluOptions, NO, NOTRANS, superlu_options_t::PrintStat, set_default_options(), and superlu_options_t::Trans.

Referenced by Eigen::SuperLU< MatrixType_ >::SuperLU().

◆ matrixL()

template<typename MatrixType_ >
const LMatrixType& Eigen::SuperLU< MatrixType_ >::matrixL ( ) const
inline

◆ matrixU()

template<typename MatrixType_ >
const UMatrixType& Eigen::SuperLU< MatrixType_ >::matrixU ( ) const
inline

◆ permutationP()

template<typename MatrixType_ >
const IntColVectorType& Eigen::SuperLU< MatrixType_ >::permutationP ( ) const
inline

◆ permutationQ()

template<typename MatrixType_ >
const IntRowVectorType& Eigen::SuperLU< MatrixType_ >::permutationQ ( ) const
inline

Member Data Documentation

◆ m_analysisIsOk

template<typename MatrixType_ >
int Eigen::SuperLUBase< MatrixType_, Derived >::m_analysisIsOk
protected

◆ m_extractedDataAreDirty

◆ m_factorizationIsOk

template<typename MatrixType_ >
int Eigen::SuperLUBase< MatrixType_, Derived >::m_factorizationIsOk
protected

◆ m_info

template<typename MatrixType_ >
ComputationInfo Eigen::SuperLUBase< MatrixType_, Derived >::m_info
mutableprotected

◆ m_isInitialized

template<typename MatrixType_ >
bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprotected

◆ m_l

template<typename MatrixType_ >
LUMatrixType Eigen::SuperLUBase< MatrixType_, Derived >::m_l
mutableprotected

◆ m_matrix

template<typename MatrixType_ >
LUMatrixType Eigen::SuperLUBase< MatrixType_, Derived >::m_matrix
mutableprotected

◆ m_p

template<typename MatrixType_ >
IntColVectorType Eigen::SuperLUBase< MatrixType_, Derived >::m_p
mutableprotected

◆ m_q

template<typename MatrixType_ >
IntRowVectorType Eigen::SuperLUBase< MatrixType_, Derived >::m_q
mutableprotected

◆ m_sluA

template<typename MatrixType_ >
SluMatrix Eigen::SuperLUBase< MatrixType_, Derived >::m_sluA
mutableprotected

◆ m_sluB

template<typename MatrixType_ >
SluMatrix Eigen::SuperLUBase< MatrixType_, Derived >::m_sluB
mutableprotected

◆ m_sluBerr

template<typename MatrixType_ >
Matrix<RealScalar, Dynamic, 1> Eigen::SuperLUBase< MatrixType_, Derived >::m_sluBerr
protected

◆ m_sluCscale

template<typename MatrixType_ >
Matrix<RealScalar, Dynamic, 1> Eigen::SuperLUBase< MatrixType_, Derived >::m_sluCscale
protected

◆ m_sluEqued

template<typename MatrixType_ >
char Eigen::SuperLUBase< MatrixType_, Derived >::m_sluEqued
mutableprotected

◆ m_sluEtree

template<typename MatrixType_ >
std::vector<int> Eigen::SuperLUBase< MatrixType_, Derived >::m_sluEtree
mutableprotected

◆ m_sluFerr

template<typename MatrixType_ >
Matrix<RealScalar, Dynamic, 1> Eigen::SuperLUBase< MatrixType_, Derived >::m_sluFerr
mutableprotected

◆ m_sluL

template<typename MatrixType_ >
SuperMatrix Eigen::SuperLUBase< MatrixType_, Derived >::m_sluL
mutableprotected

◆ m_sluOptions

template<typename MatrixType_ >
superlu_options_t Eigen::SuperLUBase< MatrixType_, Derived >::m_sluOptions
mutableprotected

◆ m_sluRscale

template<typename MatrixType_ >
Matrix<RealScalar, Dynamic, 1> Eigen::SuperLUBase< MatrixType_, Derived >::m_sluRscale
mutableprotected

◆ m_sluStat

template<typename MatrixType_ >
SuperLUStat_t Eigen::SuperLUBase< MatrixType_, Derived >::m_sluStat
mutableprotected

◆ m_sluU

template<typename MatrixType_ >
SuperMatrix Eigen::SuperLUBase< MatrixType_, Derived >::m_sluU
protected

◆ m_sluX

template<typename MatrixType_ >
SluMatrix Eigen::SuperLUBase< MatrixType_, Derived >::m_sluX
protected

◆ m_u

template<typename MatrixType_ >
LUMatrixType Eigen::SuperLUBase< MatrixType_, Derived >::m_u
mutableprotected

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