10 #ifndef EIGEN_SUPERLUSUPPORT_H 11 #define EIGEN_SUPERLUSUPPORT_H 15 #define DECL_GSSVX(PREFIX,FLOATTYPE,KEYTYPE) \ 17 typedef struct { FLOATTYPE for_lu; FLOATTYPE total_needed; int expansions; } PREFIX##mem_usage_t; \ 18 extern void PREFIX##gssvx(superlu_options_t *, SuperMatrix *, int *, int *, int *, \ 19 char *, FLOATTYPE *, FLOATTYPE *, SuperMatrix *, SuperMatrix *, \ 20 void *, int, SuperMatrix *, SuperMatrix *, \ 21 FLOATTYPE *, FLOATTYPE *, FLOATTYPE *, FLOATTYPE *, \ 22 PREFIX##mem_usage_t *, SuperLUStat_t *, int *); \ 24 inline float SuperLU_gssvx(superlu_options_t *options, SuperMatrix *A, \ 25 int *perm_c, int *perm_r, int *etree, char *equed, \ 26 FLOATTYPE *R, FLOATTYPE *C, SuperMatrix *L, \ 27 SuperMatrix *U, void *work, int lwork, \ 28 SuperMatrix *B, SuperMatrix *X, \ 29 FLOATTYPE *recip_pivot_growth, \ 30 FLOATTYPE *rcond, FLOATTYPE *ferr, FLOATTYPE *berr, \ 31 SuperLUStat_t *stats, int *info, KEYTYPE) { \ 32 PREFIX##mem_usage_t mem_usage; \ 33 PREFIX##gssvx(options, A, perm_c, perm_r, etree, equed, R, C, L, \ 34 U, work, lwork, B, X, recip_pivot_growth, rcond, \ 35 ferr, berr, &mem_usage, stats, info); \ 36 return mem_usage.for_lu; \ 45 #define EIGEN_SUPERLU_HAS_ILU 48 #ifdef EIGEN_SUPERLU_HAS_ILU 51 #define DECL_GSISX(PREFIX,FLOATTYPE,KEYTYPE) \ 53 extern void PREFIX##gsisx(superlu_options_t *, SuperMatrix *, int *, int *, int *, \ 54 char *, FLOATTYPE *, FLOATTYPE *, SuperMatrix *, SuperMatrix *, \ 55 void *, int, SuperMatrix *, SuperMatrix *, FLOATTYPE *, FLOATTYPE *, \ 56 PREFIX##mem_usage_t *, SuperLUStat_t *, int *); \ 58 inline float SuperLU_gsisx(superlu_options_t *options, SuperMatrix *A, \ 59 int *perm_c, int *perm_r, int *etree, char *equed, \ 60 FLOATTYPE *R, FLOATTYPE *C, SuperMatrix *L, \ 61 SuperMatrix *U, void *work, int lwork, \ 62 SuperMatrix *B, SuperMatrix *X, \ 63 FLOATTYPE *recip_pivot_growth, \ 65 SuperLUStat_t *stats, int *info, KEYTYPE) { \ 66 PREFIX##mem_usage_t mem_usage; \ 67 PREFIX##gsisx(options, A, perm_c, perm_r, etree, equed, R, C, L, \ 68 U, work, lwork, B, X, recip_pivot_growth, rcond, \ 69 &mem_usage, stats, info); \ 70 return mem_usage.for_lu; \ 73 DECL_GSISX(s,
float,
float)
74 DECL_GSISX(c,
float,
std::complex<
float>)
75 DECL_GSISX(d,
double,
double)
76 DECL_GSISX(z,
double,
std::complex<
double>)
80 template<
typename MatrixType>
106 SuperMatrix::operator=(static_cast<const SuperMatrix&>(other));
123 if (t==SLU_NC || t==SLU_NR || t==SLU_DN)
132 template<
typename Scalar>
145 eigen_assert(
false &&
"Scalar type not supported by SuperLU");
149 template<
typename MatrixType>
152 MatrixType& mat(_mat.derived());
153 eigen_assert( ((MatrixType::Flags&
RowMajorBit)!=RowMajorBit) &&
"row-major dense matrices are not supported by SuperLU");
159 res.nrow = mat.rows();
160 res.ncol = mat.cols();
162 res.
storage.
lda = MatrixType::IsVectorAtCompileTime ? mat.size() : mat.outerStride();
167 template<
typename MatrixType>
174 res.nrow = mat.
cols();
175 res.ncol = mat.
rows();
180 res.nrow = mat.
rows();
181 res.ncol = mat.
cols();
194 if (MatrixType::Flags &
Upper)
196 if (MatrixType::Flags &
Lower)
205 template<
typename Scalar,
int Rows,
int Cols,
int Options,
int MRows,
int MCols>
211 eigen_assert( ((Options&
RowMajor)!=RowMajor) &&
"row-major dense matrices is not supported by SuperLU");
216 res.nrow = mat.
rows();
217 res.ncol = mat.
cols();
224 template<
typename Derived>
233 res.nrow = mat.cols();
234 res.ncol = mat.rows();
239 res.nrow = mat.rows();
240 res.ncol = mat.cols();
253 if (MatrixType::Flags &
Upper)
255 if (MatrixType::Flags &
Lower)
264 template<
typename MatrixType>
271 template<
typename Scalar,
int Flags,
typename Index>
275 || (Flags&
ColMajor)==ColMajor && sluMat.Stype == SLU_NC);
277 Index outerSize = (Flags&
RowMajor)==RowMajor ? sluMat.ncol : sluMat.nrow;
290 template<
typename _MatrixType,
typename Derived>
295 typedef typename MatrixType::Scalar
Scalar;
297 typedef typename MatrixType::Index
Index;
312 Derived&
derived() {
return *
static_cast<Derived*
>(
this); }
313 const Derived&
derived()
const {
return *
static_cast<const Derived*
>(
this); }
315 inline Index
rows()
const {
return m_matrix.rows(); }
316 inline Index
cols()
const {
return m_matrix.cols(); }
319 inline superlu_options_t&
options() {
return m_sluOptions; }
328 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
335 derived().analyzePattern(matrix);
336 derived().factorize(matrix);
343 template<
typename Rhs>
346 eigen_assert(m_isInitialized &&
"SuperLU is not initialized.");
348 &&
"SuperLU::solve(): invalid number of rows of the right hand side matrix b");
356 template<
typename Rhs>
359 eigen_assert(m_isInitialized &&
"SuperLU is not initialized.");
361 &&
"SuperLU::solve(): invalid number of rows of the right hand side matrix b");
373 m_isInitialized =
true;
375 m_analysisIsOk =
true;
376 m_factorizationIsOk =
false;
379 template<
typename Stream>
387 set_default_options(&this->m_sluOptions);
389 const int size = a.rows();
397 m_sluRscale.resize(size);
398 m_sluCscale.resize(size);
399 m_sluEtree.resize(size);
402 m_sluB.setStorageType(SLU_DN);
403 m_sluB.setScalarType<Scalar>();
404 m_sluB.Mtype = SLU_GE;
405 m_sluB.storage.values = 0;
408 m_sluB.storage.lda = size;
411 m_extractedDataAreDirty =
true;
417 m_isInitialized =
false;
422 void extractData()
const;
427 Destroy_SuperNode_Matrix(&m_sluL);
429 Destroy_CompCol_Matrix(&m_sluU);
434 memset(&m_sluL,0,
sizeof m_sluL);
435 memset(&m_sluU,0,
sizeof m_sluU);
441 mutable IntColVectorType
m_p;
442 mutable IntRowVectorType
m_q;
478 template<
typename _MatrixType>
500 Base::compute(matrix);
516 m_isInitialized =
false;
517 Base::analyzePattern(matrix);
526 void factorize(
const MatrixType& matrix);
528 #ifndef EIGEN_PARSED_BY_DOXYGEN 530 template<
typename Rhs,
typename Dest>
532 #endif // EIGEN_PARSED_BY_DOXYGEN 536 if (m_extractedDataAreDirty) this->extractData();
542 if (m_extractedDataAreDirty) this->extractData();
548 if (m_extractedDataAreDirty) this->extractData();
554 if (m_extractedDataAreDirty) this->extractData();
558 Scalar determinant()
const;
562 using Base::m_matrix;
563 using Base::m_sluOptions;
569 using Base::m_sluEtree;
570 using Base::m_sluEqued;
571 using Base::m_sluRscale;
572 using Base::m_sluCscale;
575 using Base::m_sluStat;
576 using Base::m_sluFerr;
577 using Base::m_sluBerr;
581 using Base::m_analysisIsOk;
582 using Base::m_factorizationIsOk;
583 using Base::m_extractedDataAreDirty;
584 using Base::m_isInitialized;
591 set_default_options(&this->m_sluOptions);
592 m_sluOptions.PrintStat = NO;
593 m_sluOptions.ConditionNumber = NO;
594 m_sluOptions.Trans = NOTRANS;
595 m_sluOptions.ColPerm = COLAMD;
603 template<
typename MatrixType>
606 eigen_assert(m_analysisIsOk &&
"You must first call analyzePattern()");
613 this->initFactorization(a);
615 m_sluOptions.ColPerm = COLAMD;
620 StatInit(&m_sluStat);
621 SuperLU_gssvx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0],
622 &m_sluEqued, &m_sluRscale[0], &m_sluCscale[0],
626 &recip_pivot_growth, &rcond,
628 &m_sluStat, &info,
Scalar());
629 StatFree(&m_sluStat);
631 m_extractedDataAreDirty =
true;
635 m_factorizationIsOk =
true;
638 template<
typename MatrixType>
639 template<
typename Rhs,
typename Dest>
642 eigen_assert(m_factorizationIsOk &&
"The decomposition is not in a valid state for solving, you must first call either compute() or analyzePattern()/factorize()");
644 const int size = m_matrix.rows();
645 const int rhsCols = b.cols();
648 m_sluOptions.Trans = NOTRANS;
649 m_sluOptions.Fact = FACTORED;
650 m_sluOptions.IterRefine = NOREFINE;
653 m_sluFerr.resize(rhsCols);
654 m_sluBerr.resize(rhsCols);
658 typename Rhs::PlainObject b_cpy;
665 StatInit(&m_sluStat);
668 SuperLU_gssvx(&m_sluOptions, &m_sluA,
669 m_q.data(), m_p.data(),
670 &m_sluEtree[0], &m_sluEqued,
671 &m_sluRscale[0], &m_sluCscale[0],
675 &recip_pivot_growth, &rcond,
676 &m_sluFerr[0], &m_sluBerr[0],
677 &m_sluStat, &info,
Scalar());
678 StatFree(&m_sluStat);
689 template<
typename MatrixType,
typename Derived>
692 eigen_assert(m_factorizationIsOk &&
"The decomposition is not in a valid state for extracting factors, you must first call either compute() or analyzePattern()/factorize()");
693 if (m_extractedDataAreDirty)
696 int fsupc, istart, nsupr;
697 int lastl = 0, lastu = 0;
698 SCformat *Lstore =
static_cast<SCformat*
>(m_sluL.Store);
699 NCformat *Ustore =
static_cast<NCformat*
>(m_sluU.Store);
702 const int size = m_matrix.rows();
703 m_l.resize(size,size);
704 m_l.resizeNonZeros(Lstore->nnz);
705 m_u.resize(size,size);
706 m_u.resizeNonZeros(Ustore->nnz);
708 int* Lcol = m_l.outerIndexPtr();
709 int* Lrow = m_l.innerIndexPtr();
710 Scalar* Lval = m_l.valuePtr();
712 int* Ucol = m_u.outerIndexPtr();
713 int* Urow = m_u.innerIndexPtr();
714 Scalar* Uval = m_u.valuePtr();
720 for (
int k = 0; k <= Lstore->nsuper; ++k)
722 fsupc = L_FST_SUPC(k);
723 istart = L_SUB_START(fsupc);
724 nsupr = L_SUB_START(fsupc+1) - istart;
728 for (
int j = fsupc; j < L_FST_SUPC(k+1); ++j)
730 SNptr = &((
Scalar*)Lstore->nzval)[L_NZ_START(j)];
733 for (
int i = U_NZ_START(j); i < U_NZ_START(j+1); ++i)
735 Uval[lastu] = ((
Scalar*)Ustore->nzval)[i];
737 if (Uval[lastu] != 0.0)
738 Urow[lastu++] = U_SUB(i);
740 for (
int i = 0; i < upper; ++i)
743 Uval[lastu] = SNptr[i];
745 if (Uval[lastu] != 0.0)
746 Urow[lastu++] = L_SUB(istart+i);
752 Lrow[lastl++] = L_SUB(istart + upper - 1);
753 for (
int i = upper; i < nsupr; ++i)
755 Lval[lastl] = SNptr[i];
757 if (Lval[lastl] != 0.0)
758 Lrow[lastl++] = L_SUB(istart+i);
768 m_l.resizeNonZeros(lastl);
769 m_u.resizeNonZeros(lastu);
771 m_extractedDataAreDirty =
false;
775 template<
typename MatrixType>
778 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()");
780 if (m_extractedDataAreDirty)
784 for (
int j=0; j<m_u.cols(); ++j)
786 if (m_u.outerIndexPtr()[j+1]-m_u.outerIndexPtr()[j] > 0)
788 int lastId = m_u.outerIndexPtr()[j+1]-1;
790 if (m_u.innerIndexPtr()[lastId]==j)
791 det *= m_u.valuePtr()[lastId];
795 return det/m_sluRscale.prod()/m_sluCscale.prod();
800 #ifdef EIGEN_PARSED_BY_DOXYGEN 801 #define EIGEN_SUPERLU_HAS_ILU 804 #ifdef EIGEN_SUPERLU_HAS_ILU 820 template<
typename _MatrixType>
821 class SuperILU :
public SuperLUBase<_MatrixType,SuperILU<_MatrixType> >
832 SuperILU() : Base() { init(); }
834 SuperILU(
const MatrixType& matrix) : Base()
837 Base::compute(matrix);
850 void analyzePattern(
const MatrixType& matrix)
852 Base::analyzePattern(matrix);
861 void factorize(
const MatrixType& matrix);
863 #ifndef EIGEN_PARSED_BY_DOXYGEN 865 template<
typename Rhs,
typename Dest>
867 #endif // EIGEN_PARSED_BY_DOXYGEN 871 using Base::m_matrix;
872 using Base::m_sluOptions;
878 using Base::m_sluEtree;
879 using Base::m_sluEqued;
880 using Base::m_sluRscale;
881 using Base::m_sluCscale;
884 using Base::m_sluStat;
885 using Base::m_sluFerr;
886 using Base::m_sluBerr;
890 using Base::m_analysisIsOk;
891 using Base::m_factorizationIsOk;
892 using Base::m_extractedDataAreDirty;
893 using Base::m_isInitialized;
900 ilu_set_default_options(&m_sluOptions);
901 m_sluOptions.PrintStat = NO;
902 m_sluOptions.ConditionNumber = NO;
903 m_sluOptions.Trans = NOTRANS;
904 m_sluOptions.ColPerm = MMD_AT_PLUS_A;
907 m_sluOptions.ILU_MILU = SILU;
911 m_sluOptions.ILU_DropRule = DROP_BASIC;
916 SuperILU(SuperILU& ) { }
919 template<
typename MatrixType>
920 void SuperILU<MatrixType>::factorize(
const MatrixType& a)
922 eigen_assert(m_analysisIsOk &&
"You must first call analyzePattern()");
929 this->initFactorization(a);
934 StatInit(&m_sluStat);
935 SuperLU_gsisx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0],
936 &m_sluEqued, &m_sluRscale[0], &m_sluCscale[0],
940 &recip_pivot_growth, &rcond,
941 &m_sluStat, &info,
Scalar());
942 StatFree(&m_sluStat);
946 m_factorizationIsOk =
true;
949 template<
typename MatrixType>
950 template<
typename Rhs,
typename Dest>
953 eigen_assert(m_factorizationIsOk &&
"The decomposition is not in a valid state for solving, you must first call either compute() or analyzePattern()/factorize()");
955 const int size = m_matrix.rows();
956 const int rhsCols = b.cols();
959 m_sluOptions.Trans = NOTRANS;
960 m_sluOptions.Fact = FACTORED;
961 m_sluOptions.IterRefine = NOREFINE;
963 m_sluFerr.resize(rhsCols);
964 m_sluBerr.resize(rhsCols);
968 typename Rhs::PlainObject b_cpy;
978 StatInit(&m_sluStat);
979 SuperLU_gsisx(&m_sluOptions, &m_sluA,
980 m_q.data(), m_p.data(),
981 &m_sluEtree[0], &m_sluEqued,
982 &m_sluRscale[0], &m_sluCscale[0],
986 &recip_pivot_growth, &rcond,
987 &m_sluStat, &info,
Scalar());
988 StatFree(&m_sluStat);
996 template<
typename _MatrixType,
typename Derived,
typename Rhs>
1003 template<typename Dest>
void evalTo(Dest& dst)
const 1005 dec().derived()._solve(rhs(),dst);
1009 template<
typename _MatrixType,
typename Derived,
typename Rhs>
1016 template<typename Dest>
void evalTo(Dest& dst)
const 1018 this->defaultEvalTo(dst);
1026 #endif // EIGEN_SUPERLUSUPPORT_H
const IntColVectorType & permutationP() const
MappedSparseMatrix< Scalar, Flags, Index > map_superlu(SluMatrix &sluMat)
SluMatrix asSluMatrix(MatrixType &mat)
void compute(const MatrixType &matrix)
Matrix< int, 1, MatrixType::ColsAtCompileTime > IntRowVectorType
Index outerStride() const
const LMatrixType & matrixL() const
const internal::sparse_solve_retval< SuperLUBase, Rhs > solve(const SparseMatrixBase< Rhs > &b) const
A sparse direct LU factorization and solver based on the SuperLU library.
MatrixType::Scalar Scalar
bool m_extractedDataAreDirty
#define DECL_GSSVX(PREFIX, FLOATTYPE, KEYTYPE)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
ComputationInfo info() const
Reports whether previous computation was successful.
const Derived & derived() const
const IntRowVectorType & permutationQ() const
const internal::solve_retval< SuperLUBase, Rhs > solve(const MatrixBase< Rhs > &b) const
const unsigned int RowMajorBit
void _solve(const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
Matrix< Scalar, Rows, Cols, Options, MRows, MCols > MatrixType
void analyzePattern(const MatrixType &)
Base::IntColVectorType IntColVectorType
EIGEN_STRONG_INLINE Index rows() const
SuperLU(const MatrixType &matrix)
Base class of any sparse matrices or sparse expressions.
Base::RealScalar RealScalar
SluMatrix & operator=(const SluMatrix &other)
struct Eigen::SluMatrix::@438 storage
superlu_options_t m_sluOptions
void dumpMemory(Stream &)
The base class for the direct and incomplete LU factorization of SuperLU.
SuperLUBase(SuperLUBase &)
Matrix< RealScalar, Dynamic, 1 > m_sluRscale
TFSIMD_FORCE_INLINE const tfScalar & x() const
std::vector< int > m_sluEtree
static void run(MatrixType &mat, SluMatrix &res)
TriangularView< LUMatrixType, Lower|UnitDiag > LMatrixType
#define EIGEN_MAKE_SPARSE_SOLVE_HELPERS(DecompositionType, Rhs)
void initFactorization(const MatrixType &a)
Base::LUMatrixType LUMatrixType
EIGEN_STRONG_INLINE const Scalar * data() const
TriangularView< LUMatrixType, Upper > UMatrixType
Matrix< int, MatrixType::RowsAtCompileTime, 1 > IntColVectorType
const Derived & derived() const
static SluMatrix Map(SparseMatrixBase< MatrixType > &mat)
Matrix< RealScalar, Dynamic, 1 > m_sluFerr
SuperLUBase< _MatrixType, SuperLU > Base
Base class for triangular part in a matrix.
SparseMatrix< Scalar > LUMatrixType
MatrixType::RealScalar RealScalar
void setStorageType(Stype_t t)
SuperLUBase< _MatrixType, Derived > Dec
#define EIGEN_MAKE_SOLVE_HELPERS(DecompositionType, Rhs)
const UMatrixType & matrixU() const
static SluMatrix Map(MatrixBase< MatrixType > &_mat)
void factorize(const MatrixType &matrix)
superlu_options_t & options()
The matrix class, also used for vectors and row-vectors.
EIGEN_STRONG_INLINE Index cols() const
SluMatrix(const SluMatrix &other)
Matrix< Scalar, Dynamic, 1 > Vector
Scalar determinant() const
Base class for all dense matrices, vectors, and expressions.
Base::IntRowVectorType IntRowVectorType
void analyzePattern(const MatrixType &matrix)
SuperLUBase< _MatrixType, Derived > Dec
static void run(MatrixType &mat, SluMatrix &res)