10#ifndef EIGEN_SUPERLUSUPPORT_H
11#define EIGEN_SUPERLUSUPPORT_H
14#include "./InternalHeaderCheck.h"
18#if defined(SUPERLU_MAJOR_VERSION) && (SUPERLU_MAJOR_VERSION >= 5)
19#define DECL_GSSVX(PREFIX, FLOATTYPE, KEYTYPE) \
21 extern void PREFIX##gssvx(superlu_options_t *, SuperMatrix *, int *, int *, int *, char *, FLOATTYPE *, FLOATTYPE *, \
22 SuperMatrix *, SuperMatrix *, void *, int, SuperMatrix *, SuperMatrix *, FLOATTYPE *, \
23 FLOATTYPE *, FLOATTYPE *, FLOATTYPE *, GlobalLU_t *, mem_usage_t *, SuperLUStat_t *, \
26 inline float SuperLU_gssvx(superlu_options_t *options, SuperMatrix *A, int *perm_c, int *perm_r, int *etree, \
27 char *equed, FLOATTYPE *R, FLOATTYPE *C, SuperMatrix *L, SuperMatrix *U, void *work, \
28 int lwork, SuperMatrix *B, SuperMatrix *X, FLOATTYPE *recip_pivot_growth, \
29 FLOATTYPE *rcond, FLOATTYPE *ferr, FLOATTYPE *berr, SuperLUStat_t *stats, int *info, \
31 mem_usage_t mem_usage; \
33 PREFIX##gssvx(options, A, perm_c, perm_r, etree, equed, R, C, L, U, work, lwork, B, X, recip_pivot_growth, rcond, \
34 ferr, berr, &gLU, &mem_usage, stats, info); \
35 return mem_usage.for_lu; \
38#define DECL_GSSVX(PREFIX, FLOATTYPE, KEYTYPE) \
40 extern void PREFIX##gssvx(superlu_options_t *, SuperMatrix *, int *, int *, int *, char *, FLOATTYPE *, FLOATTYPE *, \
41 SuperMatrix *, SuperMatrix *, void *, int, SuperMatrix *, SuperMatrix *, FLOATTYPE *, \
42 FLOATTYPE *, FLOATTYPE *, FLOATTYPE *, mem_usage_t *, SuperLUStat_t *, int *); \
44 inline float SuperLU_gssvx(superlu_options_t *options, SuperMatrix *A, int *perm_c, int *perm_r, int *etree, \
45 char *equed, FLOATTYPE *R, FLOATTYPE *C, SuperMatrix *L, SuperMatrix *U, void *work, \
46 int lwork, SuperMatrix *B, SuperMatrix *X, FLOATTYPE *recip_pivot_growth, \
47 FLOATTYPE *rcond, FLOATTYPE *ferr, FLOATTYPE *berr, SuperLUStat_t *stats, int *info, \
49 mem_usage_t mem_usage; \
50 PREFIX##gssvx(options, A, perm_c, perm_r, etree, equed, R, C, L, U, work, lwork, B, X, recip_pivot_growth, rcond, \
51 ferr, berr, &mem_usage, stats, info); \
52 return mem_usage.for_lu; \
56DECL_GSSVX(s,
float,
float)
57DECL_GSSVX(c,
float, std::complex<float>)
58DECL_GSSVX(d,
double,
double)
59DECL_GSSVX(z,
double, std::complex<double>)
62#define EIGEN_SUPERLU_HAS_ILU
65#ifdef EIGEN_SUPERLU_HAS_ILU
68#if defined(SUPERLU_MAJOR_VERSION) && (SUPERLU_MAJOR_VERSION >= 5)
69#define DECL_GSISX(PREFIX, FLOATTYPE, KEYTYPE) \
71 extern void PREFIX##gsisx(superlu_options_t *, SuperMatrix *, int *, int *, int *, char *, FLOATTYPE *, FLOATTYPE *, \
72 SuperMatrix *, SuperMatrix *, void *, int, SuperMatrix *, SuperMatrix *, FLOATTYPE *, \
73 FLOATTYPE *, GlobalLU_t *, mem_usage_t *, SuperLUStat_t *, int *); \
75 inline float SuperLU_gsisx(superlu_options_t *options, SuperMatrix *A, int *perm_c, int *perm_r, int *etree, \
76 char *equed, FLOATTYPE *R, FLOATTYPE *C, SuperMatrix *L, SuperMatrix *U, void *work, \
77 int lwork, SuperMatrix *B, SuperMatrix *X, FLOATTYPE *recip_pivot_growth, \
78 FLOATTYPE *rcond, SuperLUStat_t *stats, int *info, KEYTYPE) { \
79 mem_usage_t mem_usage; \
81 PREFIX##gsisx(options, A, perm_c, perm_r, etree, equed, R, C, L, U, work, lwork, B, X, recip_pivot_growth, rcond, \
82 &gLU, &mem_usage, stats, info); \
83 return mem_usage.for_lu; \
86#define DECL_GSISX(PREFIX, FLOATTYPE, KEYTYPE) \
88 extern void PREFIX##gsisx(superlu_options_t *, SuperMatrix *, int *, int *, int *, char *, FLOATTYPE *, FLOATTYPE *, \
89 SuperMatrix *, SuperMatrix *, void *, int, SuperMatrix *, SuperMatrix *, FLOATTYPE *, \
90 FLOATTYPE *, mem_usage_t *, SuperLUStat_t *, int *); \
92 inline float SuperLU_gsisx(superlu_options_t *options, SuperMatrix *A, int *perm_c, int *perm_r, int *etree, \
93 char *equed, FLOATTYPE *R, FLOATTYPE *C, SuperMatrix *L, SuperMatrix *U, void *work, \
94 int lwork, SuperMatrix *B, SuperMatrix *X, FLOATTYPE *recip_pivot_growth, \
95 FLOATTYPE *rcond, SuperLUStat_t *stats, int *info, KEYTYPE) { \
96 mem_usage_t mem_usage; \
97 PREFIX##gsisx(options, A, perm_c, perm_r, etree, equed, R, C, L, U, work, lwork, B, X, recip_pivot_growth, rcond, \
98 &mem_usage, stats, info); \
99 return mem_usage.for_lu; \
103DECL_GSISX(s,
float,
float)
104DECL_GSISX(c,
float, std::complex<float>)
105DECL_GSISX(d,
double,
double)
106DECL_GSISX(z,
double, std::complex<double>)
110template <
typename MatrixType>
111struct SluMatrixMapHelper;
120struct SluMatrix : SuperMatrix {
121 SluMatrix() { Store = &storage; }
123 SluMatrix(
const SluMatrix &other) : SuperMatrix(other) {
125 storage = other.storage;
128 SluMatrix &operator=(
const SluMatrix &other) {
129 SuperMatrix::operator=(
static_cast<const SuperMatrix &
>(other));
131 storage = other.storage;
145 void setStorageType(Stype_t t) {
147 if (t == SLU_NC || t == SLU_NR || t == SLU_DN)
150 eigen_assert(
false &&
"storage type not supported");
155 template <
typename Scalar>
156 void setScalarType() {
157 if (internal::is_same<Scalar, float>::value)
159 else if (internal::is_same<Scalar, double>::value)
161 else if (internal::is_same<Scalar, std::complex<float> >::value)
163 else if (internal::is_same<Scalar, std::complex<double> >::value)
166 eigen_assert(
false &&
"Scalar type not supported by SuperLU");
170 template <
typename MatrixType>
171 static SluMatrix Map(MatrixBase<MatrixType> &_mat) {
172 MatrixType &mat(_mat.derived());
174 "row-major dense matrices are not supported by SuperLU");
176 res.setStorageType(SLU_DN);
177 res.setScalarType<
typename MatrixType::Scalar>();
180 res.nrow = internal::convert_index<int>(mat.rows());
181 res.ncol = internal::convert_index<int>(mat.cols());
183 res.storage.lda = internal::convert_index<int>(MatrixType::IsVectorAtCompileTime ? mat.size() : mat.outerStride());
184 res.storage.values = (
void *)(mat.data());
188 template <
typename MatrixType>
189 static SluMatrix Map(SparseMatrixBase<MatrixType> &a_mat) {
190 MatrixType &mat(a_mat.derived());
193 res.setStorageType(SLU_NR);
194 res.nrow = internal::convert_index<int>(mat.cols());
195 res.ncol = internal::convert_index<int>(mat.rows());
197 res.setStorageType(SLU_NC);
198 res.nrow = internal::convert_index<int>(mat.rows());
199 res.ncol = internal::convert_index<int>(mat.cols());
204 res.storage.nnz = internal::convert_index<int>(mat.nonZeros());
205 res.storage.values = mat.valuePtr();
206 res.storage.innerInd = mat.innerIndexPtr();
207 res.storage.outerInd = mat.outerIndexPtr();
209 res.setScalarType<
typename MatrixType::Scalar>();
212 if (
int(MatrixType::Flags) &
int(
Upper)) res.Mtype = SLU_TRU;
213 if (
int(MatrixType::Flags) &
int(
Lower)) res.Mtype = SLU_TRL;
215 eigen_assert(((
int(MatrixType::Flags) &
int(
SelfAdjoint)) == 0) &&
216 "SelfAdjoint matrix shape not supported by SuperLU");
222template <
typename Scalar,
int Rows,
int Cols,
int Options,
int MRows,
int MCols>
223struct SluMatrixMapHelper<
Matrix<
Scalar, Rows, Cols, Options, MRows, MCols> > {
224 typedef Matrix<Scalar, Rows, Cols, Options, MRows, MCols> MatrixType;
225 static void run(MatrixType &mat, SluMatrix &res) {
226 eigen_assert(((Options &
RowMajor) !=
RowMajor) &&
"row-major dense matrices is not supported by SuperLU");
227 res.setStorageType(SLU_DN);
228 res.setScalarType<Scalar>();
231 res.nrow = mat.rows();
232 res.ncol = mat.cols();
234 res.storage.lda = mat.outerStride();
235 res.storage.values = mat.data();
239template <
typename Derived>
241 typedef Derived MatrixType;
242 static void run(MatrixType &mat, SluMatrix &res) {
244 res.setStorageType(SLU_NR);
245 res.nrow = mat.cols();
246 res.ncol = mat.rows();
248 res.setStorageType(SLU_NC);
249 res.nrow = mat.rows();
250 res.ncol = mat.cols();
255 res.storage.nnz = mat.nonZeros();
256 res.storage.values = mat.valuePtr();
257 res.storage.innerInd = mat.innerIndexPtr();
258 res.storage.outerInd = mat.outerIndexPtr();
260 res.setScalarType<
typename MatrixType::Scalar>();
263 if (MatrixType::Flags &
Upper) res.Mtype = SLU_TRU;
264 if (MatrixType::Flags &
Lower) res.Mtype = SLU_TRL;
266 eigen_assert(((MatrixType::Flags &
SelfAdjoint) == 0) &&
"SelfAdjoint matrix shape not supported by SuperLU");
272template <
typename MatrixType>
273SluMatrix asSluMatrix(MatrixType &mat) {
274 return SluMatrix::Map(mat);
278template <
typename Scalar,
int Flags,
typename Index>
279Map<SparseMatrix<Scalar, Flags, Index> > map_superlu(SluMatrix &sluMat) {
285 return Map<SparseMatrix<Scalar, Flags, Index> >(sluMat.nrow, sluMat.ncol, sluMat.storage.outerInd[outerSize],
286 sluMat.storage.outerInd, sluMat.storage.innerInd,
287 reinterpret_cast<Scalar *
>(sluMat.storage.values));
296template <
typename MatrixType_,
typename Derived>
301 using Base::m_isInitialized;
304 typedef MatrixType_ MatrixType;
305 typedef typename MatrixType::Scalar Scalar;
306 typedef typename MatrixType::RealScalar RealScalar;
307 typedef typename MatrixType::StorageIndex StorageIndex;
313 enum { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime };
318 ~SuperLUBase() { clearFactors(); }
320 inline Index rows()
const {
return m_matrix.rows(); }
321 inline Index cols()
const {
return m_matrix.cols(); }
324 inline superlu_options_t &
options() {
return m_sluOptions; }
332 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
338 derived().analyzePattern(matrix);
339 derived().factorize(matrix);
349 m_isInitialized =
true;
351 m_analysisIsOk =
true;
352 m_factorizationIsOk =
false;
355 template <
typename Stream>
356 void dumpMemory(Stream & ) {}
359 void initFactorization(
const MatrixType &a) {
360 set_default_options(&this->m_sluOptions);
362 const Index size = a.rows();
365 m_sluA = internal::asSluMatrix(m_matrix);
372 m_sluEtree.resize(size);
375 m_sluB.setStorageType(SLU_DN);
376 m_sluB.setScalarType<Scalar>();
377 m_sluB.Mtype = SLU_GE;
378 m_sluB.storage.values = 0;
381 m_sluB.storage.lda = internal::convert_index<int>(size);
384 m_extractedDataAreDirty =
true;
389 m_isInitialized =
false;
394 void extractData()
const;
396 void clearFactors() {
397 if (m_sluL.Store) Destroy_SuperNode_Matrix(&m_sluL);
398 if (m_sluU.Store) Destroy_CompCol_Matrix(&m_sluU);
403 memset(&m_sluL, 0,
sizeof m_sluL);
404 memset(&m_sluU, 0,
sizeof m_sluU);
408 mutable LUMatrixType m_l;
409 mutable LUMatrixType m_u;
410 mutable IntColVectorType m_p;
411 mutable IntRowVectorType m_q;
413 mutable LUMatrixType m_matrix;
414 mutable SluMatrix m_sluA;
415 mutable SuperMatrix m_sluL, m_sluU;
416 mutable SluMatrix m_sluB, m_sluX;
417 mutable SuperLUStat_t m_sluStat;
418 mutable superlu_options_t m_sluOptions;
419 mutable std::vector<int> m_sluEtree;
420 mutable Matrix<RealScalar, Dynamic, 1> m_sluRscale, m_sluCscale;
421 mutable Matrix<RealScalar, Dynamic, 1> m_sluFerr, m_sluBerr;
422 mutable char m_sluEqued;
425 int m_factorizationIsOk;
427 mutable bool m_extractedDataAreDirty;
430 SuperLUBase(SuperLUBase &) {}
449template <
typename MatrixType_>
450class SuperLU :
public SuperLUBase<MatrixType_, SuperLU<MatrixType_> > {
452 typedef SuperLUBase<MatrixType_, SuperLU> Base;
453 typedef MatrixType_ MatrixType;
454 typedef typename Base::Scalar Scalar;
455 typedef typename Base::RealScalar RealScalar;
456 typedef typename Base::StorageIndex StorageIndex;
457 typedef typename Base::IntRowVectorType IntRowVectorType;
458 typedef typename Base::IntColVectorType IntColVectorType;
459 typedef typename Base::PermutationMap PermutationMap;
460 typedef typename Base::LUMatrixType LUMatrixType;
465 using Base::_solve_impl;
467 SuperLU() : Base() { init(); }
469 explicit SuperLU(
const MatrixType &matrix) : Base() {
484 m_isInitialized =
false;
498 template <
typename Rhs,
typename Dest>
501 inline const LMatrixType &matrixL()
const {
502 if (m_extractedDataAreDirty) this->extractData();
506 inline const UMatrixType &matrixU()
const {
507 if (m_extractedDataAreDirty) this->extractData();
511 inline const IntColVectorType &permutationP()
const {
512 if (m_extractedDataAreDirty) this->extractData();
516 inline const IntRowVectorType &permutationQ()
const {
517 if (m_extractedDataAreDirty) this->extractData();
521 Scalar determinant()
const;
525 using Base::m_matrix;
530 using Base::m_sluBerr;
531 using Base::m_sluCscale;
532 using Base::m_sluEqued;
533 using Base::m_sluEtree;
534 using Base::m_sluFerr;
536 using Base::m_sluOptions;
537 using Base::m_sluRscale;
538 using Base::m_sluStat;
543 using Base::m_analysisIsOk;
544 using Base::m_extractedDataAreDirty;
545 using Base::m_factorizationIsOk;
547 using Base::m_isInitialized;
552 set_default_options(&this->m_sluOptions);
553 m_sluOptions.PrintStat = NO;
554 m_sluOptions.ConditionNumber = NO;
555 m_sluOptions.Trans = NOTRANS;
556 m_sluOptions.ColPerm = COLAMD;
560 SuperLU(SuperLU &) {}
563template <
typename MatrixType>
565 eigen_assert(m_analysisIsOk &&
"You must first call analyzePattern()");
566 if (!m_analysisIsOk) {
571 this->initFactorization(a);
573 m_sluOptions.ColPerm = COLAMD;
575 RealScalar recip_pivot_growth, rcond;
576 RealScalar ferr, berr;
578 StatInit(&m_sluStat);
579 SuperLU_gssvx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0], &m_sluEqued, &m_sluRscale[0],
580 &m_sluCscale[0], &m_sluL, &m_sluU, NULL, 0, &m_sluB, &m_sluX, &recip_pivot_growth, &rcond, &ferr, &berr,
581 &m_sluStat, &
info, Scalar());
582 StatFree(&m_sluStat);
584 m_extractedDataAreDirty =
true;
588 m_factorizationIsOk =
true;
591template <
typename MatrixType>
592template <
typename Rhs,
typename Dest>
594 eigen_assert(m_factorizationIsOk &&
595 "The decomposition is not in a valid state for solving, you must first call either compute() or "
596 "analyzePattern()/factorize()");
598 const Index rhsCols = b.cols();
599 eigen_assert(m_matrix.rows() == b.rows());
601 m_sluOptions.Trans = NOTRANS;
602 m_sluOptions.Fact = FACTORED;
603 m_sluOptions.IterRefine = NOREFINE;
605 m_sluFerr.resize(rhsCols);
606 m_sluBerr.resize(rhsCols);
611 m_sluB = SluMatrix::Map(b_ref.const_cast_derived());
612 m_sluX = SluMatrix::Map(x_ref.const_cast_derived());
614 typename Rhs::PlainObject b_cpy;
615 if (m_sluEqued !=
'N') {
617 m_sluB = SluMatrix::Map(b_cpy.const_cast_derived());
620 StatInit(&m_sluStat);
622 RealScalar recip_pivot_growth, rcond;
623 SuperLU_gssvx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0], &m_sluEqued, &m_sluRscale[0],
624 &m_sluCscale[0], &m_sluL, &m_sluU, NULL, 0, &m_sluB, &m_sluX, &recip_pivot_growth, &rcond,
625 &m_sluFerr[0], &m_sluBerr[0], &m_sluStat, &info, Scalar());
626 StatFree(&m_sluStat);
628 if (x.derived().data() != x_ref.data()) x = x_ref;
640template <
typename MatrixType,
typename Derived>
641void SuperLUBase<MatrixType, Derived>::extractData()
const {
642 eigen_assert(m_factorizationIsOk &&
643 "The decomposition is not in a valid state for extracting factors, you must first call either compute() "
644 "or analyzePattern()/factorize()");
645 if (m_extractedDataAreDirty) {
647 int fsupc, istart, nsupr;
648 int lastl = 0, lastu = 0;
649 SCformat *Lstore =
static_cast<SCformat *
>(m_sluL.Store);
650 NCformat *Ustore =
static_cast<NCformat *
>(m_sluU.Store);
653 const Index size = m_matrix.rows();
654 m_l.resize(size, size);
655 m_l.resizeNonZeros(Lstore->nnz);
656 m_u.resize(size, size);
657 m_u.resizeNonZeros(Ustore->nnz);
659 int *Lcol = m_l.outerIndexPtr();
660 int *Lrow = m_l.innerIndexPtr();
661 Scalar *Lval = m_l.valuePtr();
663 int *Ucol = m_u.outerIndexPtr();
664 int *Urow = m_u.innerIndexPtr();
665 Scalar *Uval = m_u.valuePtr();
671 for (
int k = 0; k <= Lstore->nsuper; ++k) {
672 fsupc = L_FST_SUPC(k);
673 istart = L_SUB_START(fsupc);
674 nsupr = L_SUB_START(fsupc + 1) - istart;
678 for (
int j = fsupc; j < L_FST_SUPC(k + 1); ++j) {
679 SNptr = &((
Scalar *)Lstore->nzval)[L_NZ_START(j)];
682 for (
int i = U_NZ_START(j); i < U_NZ_START(j + 1); ++i) {
683 Uval[lastu] = ((
Scalar *)Ustore->nzval)[i];
685 if (Uval[lastu] != 0.0) Urow[lastu++] = U_SUB(i);
687 for (
int i = 0; i < upper; ++i) {
689 Uval[lastu] = SNptr[i];
691 if (Uval[lastu] != 0.0) Urow[lastu++] = L_SUB(istart + i);
697 Lrow[lastl++] = L_SUB(istart + upper - 1);
698 for (
int i = upper; i < nsupr; ++i) {
699 Lval[lastl] = SNptr[i];
701 if (Lval[lastl] != 0.0) Lrow[lastl++] = L_SUB(istart + i);
711 m_l.resizeNonZeros(lastl);
712 m_u.resizeNonZeros(lastu);
714 m_extractedDataAreDirty =
false;
718template <
typename MatrixType>
719typename SuperLU<MatrixType>::Scalar SuperLU<MatrixType>::determinant()
const {
720 eigen_assert(m_factorizationIsOk &&
721 "The decomposition is not in a valid state for computing the determinant, you must first call either "
722 "compute() or analyzePattern()/factorize()");
724 if (m_extractedDataAreDirty) this->extractData();
727 for (
int j = 0; j < m_u.cols(); ++j) {
728 if (m_u.outerIndexPtr()[j + 1] - m_u.outerIndexPtr()[j] > 0) {
729 int lastId = m_u.outerIndexPtr()[j + 1] - 1;
730 eigen_assert(m_u.innerIndexPtr()[lastId] <= j);
731 if (m_u.innerIndexPtr()[lastId] == j) det *= m_u.valuePtr()[lastId];
734 if (PermutationMap(m_p.data(), m_p.size()).determinant() * PermutationMap(m_q.data(), m_q.size()).determinant() < 0)
736 if (m_sluEqued !=
'N')
737 return det / m_sluRscale.prod() / m_sluCscale.prod();
742#ifdef EIGEN_PARSED_BY_DOXYGEN
743#define EIGEN_SUPERLU_HAS_ILU
746#ifdef EIGEN_SUPERLU_HAS_ILU
765template <
typename MatrixType_>
766class SuperILU :
public SuperLUBase<MatrixType_, SuperILU<MatrixType_> > {
768 typedef SuperLUBase<MatrixType_, SuperILU> Base;
769 typedef MatrixType_ MatrixType;
770 typedef typename Base::Scalar Scalar;
771 typedef typename Base::RealScalar RealScalar;
774 using Base::_solve_impl;
776 SuperILU() : Base() { init(); }
778 SuperILU(
const MatrixType &matrix) : Base() {
802#ifndef EIGEN_PARSED_BY_DOXYGEN
804 template <
typename Rhs,
typename Dest>
810 using Base::m_matrix;
815 using Base::m_sluBerr;
816 using Base::m_sluCscale;
817 using Base::m_sluEqued;
818 using Base::m_sluEtree;
819 using Base::m_sluFerr;
821 using Base::m_sluOptions;
822 using Base::m_sluRscale;
823 using Base::m_sluStat;
828 using Base::m_analysisIsOk;
829 using Base::m_extractedDataAreDirty;
830 using Base::m_factorizationIsOk;
832 using Base::m_isInitialized;
837 ilu_set_default_options(&m_sluOptions);
838 m_sluOptions.PrintStat = NO;
839 m_sluOptions.ConditionNumber = NO;
840 m_sluOptions.Trans = NOTRANS;
841 m_sluOptions.ColPerm = MMD_AT_PLUS_A;
844 m_sluOptions.ILU_MILU = SILU;
848 m_sluOptions.ILU_DropRule = DROP_BASIC;
849 m_sluOptions.ILU_DropTol = NumTraits<Scalar>::dummy_precision() * 10;
853 SuperILU(SuperILU &) {}
856template <
typename MatrixType>
858 eigen_assert(m_analysisIsOk &&
"You must first call analyzePattern()");
859 if (!m_analysisIsOk) {
864 this->initFactorization(a);
867 RealScalar recip_pivot_growth, rcond;
869 StatInit(&m_sluStat);
870 SuperLU_gsisx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0], &m_sluEqued, &m_sluRscale[0],
871 &m_sluCscale[0], &m_sluL, &m_sluU, NULL, 0, &m_sluB, &m_sluX, &recip_pivot_growth, &rcond, &m_sluStat,
873 StatFree(&m_sluStat);
877 m_factorizationIsOk =
true;
880#ifndef EIGEN_PARSED_BY_DOXYGEN
881template <
typename MatrixType>
882template <
typename Rhs,
typename Dest>
884 eigen_assert(m_factorizationIsOk &&
885 "The decomposition is not in a valid state for solving, you must first call either compute() or "
886 "analyzePattern()/factorize()");
888 const int rhsCols = b.cols();
889 eigen_assert(m_matrix.rows() == b.rows());
891 m_sluOptions.Trans = NOTRANS;
892 m_sluOptions.Fact = FACTORED;
893 m_sluOptions.IterRefine = NOREFINE;
895 m_sluFerr.resize(rhsCols);
896 m_sluBerr.resize(rhsCols);
901 m_sluB = SluMatrix::Map(b_ref.const_cast_derived());
902 m_sluX = SluMatrix::Map(x_ref.const_cast_derived());
904 typename Rhs::PlainObject b_cpy;
905 if (m_sluEqued !=
'N') {
907 m_sluB = SluMatrix::Map(b_cpy.const_cast_derived());
911 RealScalar recip_pivot_growth, rcond;
913 StatInit(&m_sluStat);
914 SuperLU_gsisx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0], &m_sluEqued, &m_sluRscale[0],
915 &m_sluCscale[0], &m_sluL, &m_sluU, NULL, 0, &m_sluB, &m_sluX, &recip_pivot_growth, &rcond, &m_sluStat,
917 StatFree(&m_sluStat);
919 if (x.derived().data() != x_ref.data()) x = x_ref;
A matrix or vector expression mapping an existing array of data.
Definition Map.h:96
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
constexpr void resize(Index rows, Index cols)
Definition PlainObjectBase.h:268
A matrix or vector expression mapping an existing expression.
Definition Ref.h:264
Base class of any sparse matrices or sparse expressions.
Definition SparseMatrixBase.h:30
A versatible sparse matrix representation.
Definition SparseMatrix.h:121
SparseSolverBase()
Definition SparseSolverBase.h:70
A sparse direct incomplete LU factorization and solver based on the SuperLU library.
Definition SuperLUSupport.h:766
void analyzePattern(const MatrixType &matrix)
Definition SuperLUSupport.h:791
void factorize(const MatrixType &matrix)
Definition SuperLUSupport.h:857
void compute(const MatrixType &matrix)
Definition SuperLUSupport.h:337
ComputationInfo info() const
Reports whether previous computation was successful.
Definition SuperLUSupport.h:331
void analyzePattern(const MatrixType &)
Definition SuperLUSupport.h:348
superlu_options_t & options()
Definition SuperLUSupport.h:324
void factorize(const MatrixType &matrix)
Definition SuperLUSupport.h:564
void analyzePattern(const MatrixType &matrix)
Definition SuperLUSupport.h:482
Expression of a triangular part in a matrix.
Definition TriangularMatrix.h:167
ComputationInfo
Definition Constants.h:438
@ SelfAdjoint
Definition Constants.h:227
@ Lower
Definition Constants.h:211
@ Upper
Definition Constants.h:213
@ NumericalIssue
Definition Constants.h:442
@ InvalidInput
Definition Constants.h:447
@ Success
Definition Constants.h:440
@ ColMajor
Definition Constants.h:318
@ RowMajor
Definition Constants.h:320
const unsigned int RowMajorBit
Definition Constants.h:70
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82