![]() |
|
A sparse direct LU factorization and solver based on the SuperLU library. More...
#include <SuperLUSupport.h>
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 LMatrixType & | matrixL () const |
const UMatrixType & | matrixU () const |
const IntColVectorType & | permutationP () const |
const IntRowVectorType & | permutationQ () const |
Scalar | determinant () const |
![]() | |
SuperLUBase () | |
~SuperLUBase () | |
Index | rows () const |
Index | cols () const |
superlu_options_t & | options () |
ComputationInfo | info () const |
Reports whether previous computation was successful. More... | |
void | compute (const MatrixType &matrix) |
void | analyzePattern (const MatrixType &) |
void | dumpMemory (Stream &) |
![]() | |
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 () |
![]() | |
void | initFactorization (const MatrixType &a) |
void | init () |
void | extractData () const |
void | clearFactors () |
SuperLU< MatrixType_ > & | derived () |
const SuperLU< MatrixType_ > & | derived () const |
Private Member Functions | |
SuperLU (SuperLU &) | |
Additional Inherited Members | |
![]() | |
typedef SparseSolverBase< SuperLU< MatrixType_ > > | Base |
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.
MatrixType_ | the type of the sparse matrix A, it must be a SparseMatrix<> |
\implsparsesolverconcept
typedef SuperLUBase<MatrixType_, SuperLU> Eigen::SuperLU< MatrixType_ >::Base |
typedef Base::IntColVectorType Eigen::SuperLU< MatrixType_ >::IntColVectorType |
typedef Base::IntRowVectorType Eigen::SuperLU< MatrixType_ >::IntRowVectorType |
typedef TriangularView<LUMatrixType, Lower | UnitDiag> Eigen::SuperLU< MatrixType_ >::LMatrixType |
typedef Base::LUMatrixType Eigen::SuperLU< MatrixType_ >::LUMatrixType |
typedef MatrixType_ Eigen::SuperLU< MatrixType_ >::MatrixType |
typedef Base::PermutationMap Eigen::SuperLU< MatrixType_ >::PermutationMap |
typedef Base::RealScalar Eigen::SuperLU< MatrixType_ >::RealScalar |
typedef Base::Scalar Eigen::SuperLU< MatrixType_ >::Scalar |
typedef Base::StorageIndex Eigen::SuperLU< MatrixType_ >::StorageIndex |
typedef TriangularView<LUMatrixType, Upper> Eigen::SuperLU< MatrixType_ >::UMatrixType |
|
inline |
References Eigen::SuperLU< MatrixType_ >::init().
|
inlineexplicit |
References Eigen::SuperLUBase< MatrixType_, Derived >::compute(), Eigen::SuperLU< MatrixType_ >::init(), and matrix().
|
inline |
|
inlineprivate |
void Eigen::SuperLU< MatrixType >::_solve_impl | ( | const MatrixBase< Rhs > & | b, |
MatrixBase< Dest > & | dest | ||
) | const |
References b, eigen_assert, FACTORED, info, Eigen::SluMatrix::Map(), NOREFINE, NOTRANS, Eigen::NumericalIssue, StatFree(), StatInit(), Eigen::Success, and plotDoE::x.
|
inline |
Performs a symbolic decomposition on the sparcity of matrix.
This function is particularly useful when solving for several problems having the same structure.
References Eigen::SuperLUBase< MatrixType_, Derived >::analyzePattern(), Eigen::InvalidInput, Eigen::SuperLU< MatrixType_ >::m_info, Eigen::SuperLU< MatrixType_ >::m_isInitialized, and matrix().
SuperLU< MatrixType >::Scalar Eigen::SuperLU< MatrixType >::determinant |
References eigen_assert, and j.
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.
References a, COLAMD, eigen_assert, info, Eigen::InvalidInput, Eigen::NumericalIssue, StatFree(), StatInit(), and Eigen::Success.
|
inlineprotected |
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().
|
inline |
|
inline |
|
inline |
|
inline |
|
protected |
|
mutableprotected |
|
protected |
|
mutableprotected |
Referenced by Eigen::SuperLU< MatrixType_ >::analyzePattern().
|
mutableprotected |
Referenced by Eigen::SuperLU< MatrixType_ >::analyzePattern().
|
mutableprotected |
Referenced by Eigen::SuperLU< MatrixType_ >::matrixL().
|
mutableprotected |
|
mutableprotected |
Referenced by Eigen::SuperLU< MatrixType_ >::permutationP().
|
mutableprotected |
Referenced by Eigen::SuperLU< MatrixType_ >::permutationQ().
|
mutableprotected |
|
mutableprotected |
|
protected |
|
protected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
Referenced by Eigen::SuperLU< MatrixType_ >::init().
|
mutableprotected |
|
mutableprotected |
|
protected |
|
protected |
|
mutableprotected |
Referenced by Eigen::SuperLU< MatrixType_ >::matrixU().