Eigen::SuperLU< MatrixType_ > Class Template Reference

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

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

Public Types

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

Public Member Functions

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

Protected Member Functions

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

Protected Attributes

int m_analysisIsOk
 
bool m_extractedDataAreDirty
 
int m_factorizationIsOk
 
ComputationInfo m_info
 
bool m_isInitialized
 
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< int > m_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
 
- Protected Attributes inherited from Eigen::SuperLUBase< MatrixType_, SuperLU< MatrixType_ > >
int m_analysisIsOk
 
bool m_extractedDataAreDirty
 
int m_factorizationIsOk
 
ComputationInfo m_info
 
bool m_isInitialized
 
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< int > m_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
 
- 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.

This class follows the sparse solver concept .

See also
Sparse solver concept, class SparseLU

Definition at line 490 of file SuperLUSupport.h.

Member Typedef Documentation

◆ Base

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

Definition at line 493 of file SuperLUSupport.h.

◆ IntColVectorType

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

Definition at line 499 of file SuperLUSupport.h.

◆ IntRowVectorType

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

Definition at line 498 of file SuperLUSupport.h.

◆ LMatrixType

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

Definition at line 502 of file SuperLUSupport.h.

◆ LUMatrixType

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

Definition at line 501 of file SuperLUSupport.h.

◆ MatrixType

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

Definition at line 494 of file SuperLUSupport.h.

◆ PermutationMap

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

Definition at line 500 of file SuperLUSupport.h.

◆ RealScalar

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

Definition at line 496 of file SuperLUSupport.h.

◆ Scalar

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

Definition at line 495 of file SuperLUSupport.h.

◆ StorageIndex

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

Definition at line 497 of file SuperLUSupport.h.

◆ UMatrixType

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

Definition at line 503 of file SuperLUSupport.h.

Constructor & Destructor Documentation

◆ SuperLU() [1/3]

template<typename MatrixType_ >
Eigen::SuperLU< MatrixType_ >::SuperLU ( )
inline

Definition at line 508 of file SuperLUSupport.h.

508 : Base() { init(); }
SuperLUBase< MatrixType_, SuperLU > Base

◆ SuperLU() [2/3]

template<typename MatrixType_ >
Eigen::SuperLU< MatrixType_ >::SuperLU ( const MatrixType matrix)
inlineexplicit

Definition at line 510 of file SuperLUSupport.h.

510  : Base()
511  {
512  init();
513  Base::compute(matrix);
514  }
void compute(const MatrixType &matrix)

◆ ~SuperLU()

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

Definition at line 516 of file SuperLUSupport.h.

517  {
518  }

◆ SuperLU() [3/3]

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

Definition at line 611 of file SuperLUSupport.h.

611 { }

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

Definition at line 651 of file SuperLUSupport.h.

652 {
653  eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or analyzePattern()/factorize()");
654 
655  const Index rhsCols = b.cols();
656  eigen_assert(m_matrix.rows()==b.rows());
657 
658  m_sluOptions.Trans = NOTRANS;
659  m_sluOptions.Fact = FACTORED;
660  m_sluOptions.IterRefine = NOREFINE;
661 
662 
663  m_sluFerr.resize(rhsCols);
664  m_sluBerr.resize(rhsCols);
665 
666  Ref<const Matrix<typename Rhs::Scalar,Dynamic,Dynamic,ColMajor> > b_ref(b);
667  Ref<const Matrix<typename Dest::Scalar,Dynamic,Dynamic,ColMajor> > x_ref(x);
668 
669  m_sluB = SluMatrix::Map(b_ref.const_cast_derived());
670  m_sluX = SluMatrix::Map(x_ref.const_cast_derived());
671 
672  typename Rhs::PlainObject b_cpy;
673  if(m_sluEqued!='N')
674  {
675  b_cpy = b;
676  m_sluB = SluMatrix::Map(b_cpy.const_cast_derived());
677  }
678 
679  StatInit(&m_sluStat);
680  int info = 0;
681  RealScalar recip_pivot_growth, rcond;
682  SuperLU_gssvx(&m_sluOptions, &m_sluA,
683  m_q.data(), m_p.data(),
684  &m_sluEtree[0], &m_sluEqued,
685  &m_sluRscale[0], &m_sluCscale[0],
686  &m_sluL, &m_sluU,
687  NULL, 0,
688  &m_sluB, &m_sluX,
689  &recip_pivot_growth, &rcond,
690  &m_sluFerr[0], &m_sluBerr[0],
691  &m_sluStat, &info, Scalar());
692  StatFree(&m_sluStat);
693 
694  if(x.derived().data() != x_ref.data())
695  x = x_ref;
696 
698 }
Array< int, 3, 1 > b
#define eigen_assert(x)
Definition: Macros.h:902
const Scalar * data() const
constexpr void resize(Index rows, Index cols)
Index rows() const
Definition: SparseMatrix.h:165
ComputationInfo info() const
Reports whether previous computation was successful.
std::vector< int > m_sluEtree
Matrix< RealScalar, Dynamic, 1 > m_sluRscale
SuperMatrix m_sluL
SuperLUStat_t m_sluStat
Base::Scalar Scalar
SluMatrix m_sluB
Base::RealScalar RealScalar
SuperMatrix m_sluU
IntRowVectorType m_q
Matrix< RealScalar, Dynamic, 1 > m_sluCscale
SluMatrix m_sluX
IntColVectorType m_p
Matrix< RealScalar, Dynamic, 1 > m_sluFerr
Matrix< RealScalar, Dynamic, 1 > m_sluBerr
superlu_options_t m_sluOptions
ComputationInfo m_info
SluMatrix m_sluA
LUMatrixType m_matrix
@ NumericalIssue
Definition: Constants.h:448
@ Success
Definition: Constants.h:446
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
static SluMatrix Map(MatrixBase< MatrixType > &_mat)

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

Definition at line 526 of file SuperLUSupport.h.

527  {
529  m_isInitialized = false;
530  Base::analyzePattern(matrix);
531  }
void analyzePattern(const MatrixType &)
@ InvalidInput
Definition: Constants.h:453

◆ determinant()

template<typename MatrixType >
SuperLU< MatrixType >::Scalar Eigen::SuperLU< MatrixType >::determinant

Definition at line 794 of file SuperLUSupport.h.

795 {
796  eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for computing the determinant, you must first call either compute() or analyzePattern()/factorize()");
797 
799  this->extractData();
800 
801  Scalar det = Scalar(1);
802  for (int j=0; j<m_u.cols(); ++j)
803  {
804  if (m_u.outerIndexPtr()[j+1]-m_u.outerIndexPtr()[j] > 0)
805  {
806  int lastId = m_u.outerIndexPtr()[j+1]-1;
807  eigen_assert(m_u.innerIndexPtr()[lastId]<=j);
808  if (m_u.innerIndexPtr()[lastId]==j)
809  det *= m_u.valuePtr()[lastId];
810  }
811  }
812  if(PermutationMap(m_p.data(),m_p.size()).determinant()*PermutationMap(m_q.data(),m_q.size()).determinant()<0)
813  det = -det;
814  if(m_sluEqued!='N')
815  return det/m_sluRscale.prod()/m_sluCscale.prod();
816  else
817  return det;
818 }
Index cols() const
Definition: SparseMatrix.h:167
const Scalar * valuePtr() const
Definition: SparseMatrix.h:177
const StorageIndex * outerIndexPtr() const
Definition: SparseMatrix.h:195
const StorageIndex * innerIndexPtr() const
Definition: SparseMatrix.h:186
bool m_extractedDataAreDirty
LUMatrixType m_u
Base::PermutationMap PermutationMap
std::ptrdiff_t 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()

Definition at line 615 of file SuperLUSupport.h.

616 {
617  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
618  if(!m_analysisIsOk)
619  {
621  return;
622  }
623 
624  this->initFactorization(a);
625 
626  m_sluOptions.ColPerm = COLAMD;
627  int info = 0;
628  RealScalar recip_pivot_growth, rcond;
629  RealScalar ferr, berr;
630 
631  StatInit(&m_sluStat);
632  SuperLU_gssvx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0],
634  &m_sluL, &m_sluU,
635  NULL, 0,
636  &m_sluB, &m_sluX,
637  &recip_pivot_growth, &rcond,
638  &ferr, &berr,
639  &m_sluStat, &info, Scalar());
640  StatFree(&m_sluStat);
641 
643 
644  // FIXME how to better check for errors ???
645  m_info = info == 0 ? Success : NumericalIssue;
646  m_factorizationIsOk = true;
647 }

◆ init()

template<typename MatrixType_ >
void Eigen::SuperLU< MatrixType_ >::init ( )
inlineprotected

Definition at line 598 of file SuperLUSupport.h.

599  {
600  Base::init();
601 
602  set_default_options(&this->m_sluOptions);
603  m_sluOptions.PrintStat = NO;
604  m_sluOptions.ConditionNumber = NO;
605  m_sluOptions.Trans = NOTRANS;
606  m_sluOptions.ColPerm = COLAMD;
607  }

◆ matrixL()

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

Definition at line 545 of file SuperLUSupport.h.

546  {
548  return m_l;
549  }
LUMatrixType m_l

◆ matrixU()

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

Definition at line 551 of file SuperLUSupport.h.

552  {
554  return m_u;
555  }

◆ permutationP()

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

Definition at line 557 of file SuperLUSupport.h.

558  {
560  return m_p;
561  }

◆ permutationQ()

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

Definition at line 563 of file SuperLUSupport.h.

564  {
566  return m_q;
567  }

Member Data Documentation

◆ m_analysisIsOk

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

Definition at line 465 of file SuperLUSupport.h.

◆ m_extractedDataAreDirty

template<typename MatrixType_ >
bool Eigen::SuperLUBase< MatrixType_, Derived >::m_extractedDataAreDirty
mutableprotected

Definition at line 466 of file SuperLUSupport.h.

◆ m_factorizationIsOk

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

Definition at line 464 of file SuperLUSupport.h.

◆ m_info

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

Definition at line 463 of file SuperLUSupport.h.

◆ m_isInitialized

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

Definition at line 324 of file SparseSolverBase.h.

◆ m_l

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

Definition at line 447 of file SuperLUSupport.h.

◆ m_matrix

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

Definition at line 452 of file SuperLUSupport.h.

◆ m_p

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

Definition at line 449 of file SuperLUSupport.h.

◆ m_q

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

Definition at line 450 of file SuperLUSupport.h.

◆ m_sluA

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

Definition at line 453 of file SuperLUSupport.h.

◆ m_sluB

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

Definition at line 455 of file SuperLUSupport.h.

◆ m_sluBerr

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

Definition at line 460 of file SuperLUSupport.h.

◆ m_sluCscale

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

Definition at line 459 of file SuperLUSupport.h.

◆ m_sluEqued

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

Definition at line 461 of file SuperLUSupport.h.

◆ m_sluEtree

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

Definition at line 458 of file SuperLUSupport.h.

◆ m_sluFerr

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

Definition at line 460 of file SuperLUSupport.h.

◆ m_sluL

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

Definition at line 454 of file SuperLUSupport.h.

◆ m_sluOptions

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

Definition at line 457 of file SuperLUSupport.h.

◆ m_sluRscale

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

Definition at line 459 of file SuperLUSupport.h.

◆ m_sluStat

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

Definition at line 456 of file SuperLUSupport.h.

◆ m_sluU

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

Definition at line 454 of file SuperLUSupport.h.

◆ m_sluX

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

Definition at line 455 of file SuperLUSupport.h.

◆ m_u

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

Definition at line 448 of file SuperLUSupport.h.


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