10 #ifndef EIGEN_SUPERLUSUPPORT_H
11 #define EIGEN_SUPERLUSUPPORT_H
15 #if defined(SUPERLU_MAJOR_VERSION) && (SUPERLU_MAJOR_VERSION >= 5)
16 #define DECL_GSSVX(PREFIX,FLOATTYPE,KEYTYPE) \
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 GlobalLU_t *, 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 mem_usage_t mem_usage; \
34 PREFIX##gssvx(options, A, perm_c, perm_r, etree, equed, R, C, L, \
35 U, work, lwork, B, X, recip_pivot_growth, rcond, \
36 ferr, berr, &gLU, &mem_usage, stats, info); \
37 return mem_usage.for_lu; \
39 #else // version < 5.0
40 #define DECL_GSSVX(PREFIX,FLOATTYPE,KEYTYPE) \
42 extern void PREFIX##gssvx(superlu_options_t *, SuperMatrix *, int *, int *, int *, \
43 char *, FLOATTYPE *, FLOATTYPE *, SuperMatrix *, SuperMatrix *, \
44 void *, int, SuperMatrix *, SuperMatrix *, \
45 FLOATTYPE *, FLOATTYPE *, FLOATTYPE *, FLOATTYPE *, \
46 mem_usage_t *, SuperLUStat_t *, int *); \
48 inline float SuperLU_gssvx(superlu_options_t *options, SuperMatrix *A, \
49 int *perm_c, int *perm_r, int *etree, char *equed, \
50 FLOATTYPE *R, FLOATTYPE *C, SuperMatrix *L, \
51 SuperMatrix *U, void *work, int lwork, \
52 SuperMatrix *B, SuperMatrix *X, \
53 FLOATTYPE *recip_pivot_growth, \
54 FLOATTYPE *rcond, FLOATTYPE *ferr, FLOATTYPE *berr, \
55 SuperLUStat_t *stats, int *info, KEYTYPE) { \
56 mem_usage_t mem_usage; \
57 PREFIX##gssvx(options, A, perm_c, perm_r, etree, equed, R, C, L, \
58 U, work, lwork, B, X, recip_pivot_growth, rcond, \
59 ferr, berr, &mem_usage, stats, info); \
60 return mem_usage.for_lu; \
70 #define EIGEN_SUPERLU_HAS_ILU
73 #ifdef EIGEN_SUPERLU_HAS_ILU
76 #define DECL_GSISX(PREFIX,FLOATTYPE,KEYTYPE) \
78 extern void PREFIX##gsisx(superlu_options_t *, SuperMatrix *, int *, int *, int *, \
79 char *, FLOATTYPE *, FLOATTYPE *, SuperMatrix *, SuperMatrix *, \
80 void *, int, SuperMatrix *, SuperMatrix *, FLOATTYPE *, FLOATTYPE *, \
81 mem_usage_t *, SuperLUStat_t *, int *); \
83 inline float SuperLU_gsisx(superlu_options_t *options, SuperMatrix *A, \
84 int *perm_c, int *perm_r, int *etree, char *equed, \
85 FLOATTYPE *R, FLOATTYPE *C, SuperMatrix *L, \
86 SuperMatrix *U, void *work, int lwork, \
87 SuperMatrix *B, SuperMatrix *X, \
88 FLOATTYPE *recip_pivot_growth, \
90 SuperLUStat_t *stats, int *info, KEYTYPE) { \
91 mem_usage_t mem_usage; \
92 PREFIX##gsisx(options, A, perm_c, perm_r, etree, equed, R, C, L, \
93 U, work, lwork, B, X, recip_pivot_growth, rcond, \
94 &mem_usage, stats, info); \
95 return mem_usage.for_lu; \
98 DECL_GSISX(
s,
float,
float)
99 DECL_GSISX(
c,
float,std::complex<float>)
100 DECL_GSISX(d,
double,
double)
101 DECL_GSISX(z,
double,std::complex<double>)
105 template<
typename MatrixType>
131 SuperMatrix::operator=(
static_cast<const SuperMatrix&
>(other));
148 if (t==SLU_NC || t==SLU_NR || t==SLU_DN)
157 template<
typename Scalar>
170 eigen_assert(
false &&
"Scalar type not supported by SuperLU");
174 template<
typename MatrixType>
184 res.nrow = internal::convert_index<int>(
mat.rows());
185 res.ncol = internal::convert_index<int>(
mat.cols());
187 res.
storage.
lda = internal::convert_index<int>(MatrixType::IsVectorAtCompileTime ?
mat.size() :
mat.outerStride());
192 template<
typename MatrixType>
200 res.nrow = internal::convert_index<int>(
mat.cols());
201 res.ncol = internal::convert_index<int>(
mat.rows());
206 res.nrow = internal::convert_index<int>(
mat.rows());
207 res.ncol = internal::convert_index<int>(
mat.cols());
212 res.
storage.
nnz = internal::convert_index<int>(
mat.nonZeros());
220 if (MatrixType::Flags &
Upper)
222 if (MatrixType::Flags &
Lower)
231 template<
typename Scalar,
int Rows,
int Cols,
int Options,
int MRows,
int MCols>
242 res.nrow =
mat.rows();
243 res.ncol =
mat.cols();
250 template<
typename Derived>
259 res.nrow =
mat.cols();
260 res.ncol =
mat.rows();
265 res.nrow =
mat.rows();
266 res.ncol =
mat.cols();
279 if (MatrixType::Flags &
Upper)
281 if (MatrixType::Flags &
Lower)
290 template<
typename MatrixType>
297 template<
typename Scalar,
int Flags,
typename Index>
316 template<
typename _MatrixType,
typename Derived>
385 template<
typename Stream>
433 Destroy_SuperNode_Matrix(&
m_sluL);
435 Destroy_CompCol_Matrix(&
m_sluU);
487 template<
typename _MatrixType>
540 template<
typename Rhs,
typename Dest>
612 template<
typename MatrixType>
615 eigen_assert(m_analysisIsOk &&
"You must first call analyzePattern()");
622 this->initFactorization(
a);
624 m_sluOptions.ColPerm = COLAMD;
629 StatInit(&m_sluStat);
630 SuperLU_gssvx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0],
631 &m_sluEqued, &m_sluRscale[0], &m_sluCscale[0],
635 &recip_pivot_growth, &rcond,
638 StatFree(&m_sluStat);
640 m_extractedDataAreDirty =
true;
644 m_factorizationIsOk =
true;
647 template<
typename MatrixType>
648 template<
typename Rhs,
typename Dest>
651 eigen_assert(m_factorizationIsOk &&
"The decomposition is not in a valid state for solving, you must first call either compute() or analyzePattern()/factorize()");
654 const Index rhsCols =
b.cols();
657 m_sluOptions.Trans = NOTRANS;
658 m_sluOptions.Fact = FACTORED;
659 m_sluOptions.IterRefine = NOREFINE;
662 m_sluFerr.resize(rhsCols);
663 m_sluBerr.resize(rhsCols);
671 typename Rhs::PlainObject b_cpy;
678 StatInit(&m_sluStat);
681 SuperLU_gssvx(&m_sluOptions, &m_sluA,
682 m_q.data(), m_p.data(),
683 &m_sluEtree[0], &m_sluEqued,
684 &m_sluRscale[0], &m_sluCscale[0],
688 &recip_pivot_growth, &rcond,
689 &m_sluFerr[0], &m_sluBerr[0],
691 StatFree(&m_sluStat);
693 if(
x.derived().data() != x_ref.data())
706 template<
typename MatrixType,
typename Derived>
709 eigen_assert(m_factorizationIsOk &&
"The decomposition is not in a valid state for extracting factors, you must first call either compute() or analyzePattern()/factorize()");
710 if (m_extractedDataAreDirty)
713 int fsupc, istart, nsupr;
714 int lastl = 0, lastu = 0;
715 SCformat *Lstore =
static_cast<SCformat*
>(m_sluL.Store);
716 NCformat *Ustore =
static_cast<NCformat*
>(m_sluU.Store);
721 m_l.resizeNonZeros(Lstore->nnz);
723 m_u.resizeNonZeros(Ustore->nnz);
725 int* Lcol = m_l.outerIndexPtr();
726 int* Lrow = m_l.innerIndexPtr();
727 Scalar* Lval = m_l.valuePtr();
729 int* Ucol = m_u.outerIndexPtr();
730 int* Urow = m_u.innerIndexPtr();
731 Scalar* Uval = m_u.valuePtr();
737 for (
int k = 0; k <= Lstore->nsuper; ++k)
739 fsupc = L_FST_SUPC(k);
740 istart = L_SUB_START(fsupc);
741 nsupr = L_SUB_START(fsupc+1) - istart;
745 for (
int j = fsupc; j < L_FST_SUPC(k+1); ++j)
747 SNptr = &((
Scalar*)Lstore->nzval)[L_NZ_START(j)];
750 for (
int i = U_NZ_START(j); i < U_NZ_START(j+1); ++i)
752 Uval[lastu] = ((
Scalar*)Ustore->nzval)[i];
754 if (Uval[lastu] != 0.0)
755 Urow[lastu++] = U_SUB(i);
757 for (
int i = 0; i < upper; ++i)
760 Uval[lastu] = SNptr[i];
762 if (Uval[lastu] != 0.0)
763 Urow[lastu++] = L_SUB(istart+i);
769 Lrow[lastl++] = L_SUB(istart + upper - 1);
770 for (
int i = upper; i < nsupr; ++i)
772 Lval[lastl] = SNptr[i];
774 if (Lval[lastl] != 0.0)
775 Lrow[lastl++] = L_SUB(istart+i);
785 m_l.resizeNonZeros(lastl);
786 m_u.resizeNonZeros(lastu);
788 m_extractedDataAreDirty =
false;
792 template<
typename MatrixType>
795 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 if (m_extractedDataAreDirty)
801 for (
int j=0; j<m_u.cols(); ++j)
803 if (m_u.outerIndexPtr()[j+1]-m_u.outerIndexPtr()[j] > 0)
805 int lastId = m_u.outerIndexPtr()[j+1]-1;
807 if (m_u.innerIndexPtr()[lastId]==j)
808 det *= m_u.valuePtr()[lastId];
814 return det/m_sluRscale.prod()/m_sluCscale.prod();
819 #ifdef EIGEN_PARSED_BY_DOXYGEN
820 #define EIGEN_SUPERLU_HAS_ILU
823 #ifdef EIGEN_SUPERLU_HAS_ILU
841 template<
typename _MatrixType>
842 class SuperILU :
public SuperLUBase<_MatrixType,SuperILU<_MatrixType> >
851 using Base::_solve_impl;
853 SuperILU() : Base() { init(); }
873 Base::analyzePattern(
matrix);
884 #ifndef EIGEN_PARSED_BY_DOXYGEN
886 template<
typename Rhs,
typename Dest>
887 void _solve_impl(
const MatrixBase<Rhs> &
b, MatrixBase<Dest> &dest)
const;
888 #endif // EIGEN_PARSED_BY_DOXYGEN
892 using Base::m_matrix;
893 using Base::m_sluOptions;
899 using Base::m_sluEtree;
900 using Base::m_sluEqued;
901 using Base::m_sluRscale;
902 using Base::m_sluCscale;
905 using Base::m_sluStat;
906 using Base::m_sluFerr;
907 using Base::m_sluBerr;
911 using Base::m_analysisIsOk;
912 using Base::m_factorizationIsOk;
913 using Base::m_extractedDataAreDirty;
914 using Base::m_isInitialized;
921 ilu_set_default_options(&m_sluOptions);
922 m_sluOptions.PrintStat = NO;
923 m_sluOptions.ConditionNumber = NO;
924 m_sluOptions.Trans = NOTRANS;
925 m_sluOptions.ColPerm = MMD_AT_PLUS_A;
928 m_sluOptions.ILU_MILU = SILU;
932 m_sluOptions.ILU_DropRule = DROP_BASIC;
933 m_sluOptions.ILU_DropTol = NumTraits<Scalar>::dummy_precision()*10;
937 SuperILU(SuperILU& ) { }
940 template<
typename MatrixType>
941 void SuperILU<MatrixType>::factorize(
const MatrixType&
a)
943 eigen_assert(m_analysisIsOk &&
"You must first call analyzePattern()");
950 this->initFactorization(
a);
955 StatInit(&m_sluStat);
956 SuperLU_gsisx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0],
957 &m_sluEqued, &m_sluRscale[0], &m_sluCscale[0],
961 &recip_pivot_growth, &rcond,
963 StatFree(&m_sluStat);
967 m_factorizationIsOk =
true;
970 #ifndef EIGEN_PARSED_BY_DOXYGEN
971 template<
typename MatrixType>
972 template<
typename Rhs,
typename Dest>
973 void SuperILU<MatrixType>::_solve_impl(
const MatrixBase<Rhs> &
b, MatrixBase<Dest>&
x)
const
975 eigen_assert(m_factorizationIsOk &&
"The decomposition is not in a valid state for solving, you must first call either compute() or analyzePattern()/factorize()");
977 const int size = m_matrix.rows();
978 const int rhsCols =
b.cols();
981 m_sluOptions.Trans = NOTRANS;
982 m_sluOptions.Fact = FACTORED;
983 m_sluOptions.IterRefine = NOREFINE;
985 m_sluFerr.resize(rhsCols);
986 m_sluBerr.resize(rhsCols);
988 Ref<const Matrix<typename Rhs::Scalar,Dynamic,Dynamic,ColMajor> > b_ref(
b);
989 Ref<const Matrix<typename Dest::Scalar,Dynamic,Dynamic,ColMajor> > x_ref(
x);
994 typename Rhs::PlainObject b_cpy;
1004 StatInit(&m_sluStat);
1005 SuperLU_gsisx(&m_sluOptions, &m_sluA,
1006 m_q.data(), m_p.data(),
1007 &m_sluEtree[0], &m_sluEqued,
1008 &m_sluRscale[0], &m_sluCscale[0],
1012 &recip_pivot_growth, &rcond,
1014 StatFree(&m_sluStat);
1016 if(
x.derived().data() != x_ref.data())
1027 #endif // EIGEN_SUPERLUSUPPORT_H