MOAB
4.9.3pre
|
00001 // This file is part of Eigen, a lightweight C++ template library 00002 // for linear algebra. 00003 // 00004 // Copyright (C) 2008-2015 Gael Guennebaud <[email protected]> 00005 // 00006 // This Source Code Form is subject to the terms of the Mozilla 00007 // Public License v. 2.0. If a copy of the MPL was not distributed 00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 00009 00010 #ifndef EIGEN_SUPERLUSUPPORT_H 00011 #define EIGEN_SUPERLUSUPPORT_H 00012 00013 namespace Eigen { 00014 00015 #define DECL_GSSVX(PREFIX,FLOATTYPE,KEYTYPE) \ 00016 extern "C" { \ 00017 extern void PREFIX##gssvx(superlu_options_t *, SuperMatrix *, int *, int *, int *, \ 00018 char *, FLOATTYPE *, FLOATTYPE *, SuperMatrix *, SuperMatrix *, \ 00019 void *, int, SuperMatrix *, SuperMatrix *, \ 00020 FLOATTYPE *, FLOATTYPE *, FLOATTYPE *, FLOATTYPE *, \ 00021 mem_usage_t *, SuperLUStat_t *, int *); \ 00022 } \ 00023 inline float SuperLU_gssvx(superlu_options_t *options, SuperMatrix *A, \ 00024 int *perm_c, int *perm_r, int *etree, char *equed, \ 00025 FLOATTYPE *R, FLOATTYPE *C, SuperMatrix *L, \ 00026 SuperMatrix *U, void *work, int lwork, \ 00027 SuperMatrix *B, SuperMatrix *X, \ 00028 FLOATTYPE *recip_pivot_growth, \ 00029 FLOATTYPE *rcond, FLOATTYPE *ferr, FLOATTYPE *berr, \ 00030 SuperLUStat_t *stats, int *info, KEYTYPE) { \ 00031 mem_usage_t mem_usage; \ 00032 PREFIX##gssvx(options, A, perm_c, perm_r, etree, equed, R, C, L, \ 00033 U, work, lwork, B, X, recip_pivot_growth, rcond, \ 00034 ferr, berr, &mem_usage, stats, info); \ 00035 return mem_usage.for_lu; /* bytes used by the factor storage */ \ 00036 } 00037 00038 DECL_GSSVX(s,float,float) 00039 DECL_GSSVX(c,float,std::complex<float>) 00040 DECL_GSSVX(d,double,double) 00041 DECL_GSSVX(z,double,std::complex<double>) 00042 00043 #ifdef MILU_ALPHA 00044 #define EIGEN_SUPERLU_HAS_ILU 00045 #endif 00046 00047 #ifdef EIGEN_SUPERLU_HAS_ILU 00048 00049 // similarly for the incomplete factorization using gsisx 00050 #define DECL_GSISX(PREFIX,FLOATTYPE,KEYTYPE) \ 00051 extern "C" { \ 00052 extern void PREFIX##gsisx(superlu_options_t *, SuperMatrix *, int *, int *, int *, \ 00053 char *, FLOATTYPE *, FLOATTYPE *, SuperMatrix *, SuperMatrix *, \ 00054 void *, int, SuperMatrix *, SuperMatrix *, FLOATTYPE *, FLOATTYPE *, \ 00055 mem_usage_t *, SuperLUStat_t *, int *); \ 00056 } \ 00057 inline float SuperLU_gsisx(superlu_options_t *options, SuperMatrix *A, \ 00058 int *perm_c, int *perm_r, int *etree, char *equed, \ 00059 FLOATTYPE *R, FLOATTYPE *C, SuperMatrix *L, \ 00060 SuperMatrix *U, void *work, int lwork, \ 00061 SuperMatrix *B, SuperMatrix *X, \ 00062 FLOATTYPE *recip_pivot_growth, \ 00063 FLOATTYPE *rcond, \ 00064 SuperLUStat_t *stats, int *info, KEYTYPE) { \ 00065 mem_usage_t mem_usage; \ 00066 PREFIX##gsisx(options, A, perm_c, perm_r, etree, equed, R, C, L, \ 00067 U, work, lwork, B, X, recip_pivot_growth, rcond, \ 00068 &mem_usage, stats, info); \ 00069 return mem_usage.for_lu; /* bytes used by the factor storage */ \ 00070 } 00071 00072 DECL_GSISX(s,float,float) 00073 DECL_GSISX(c,float,std::complex<float>) 00074 DECL_GSISX(d,double,double) 00075 DECL_GSISX(z,double,std::complex<double>) 00076 00077 #endif 00078 00079 template<typename MatrixType> 00080 struct SluMatrixMapHelper; 00081 00089 struct SluMatrix : SuperMatrix 00090 { 00091 SluMatrix() 00092 { 00093 Store = &storage; 00094 } 00095 00096 SluMatrix(const SluMatrix& other) 00097 : SuperMatrix(other) 00098 { 00099 Store = &storage; 00100 storage = other.storage; 00101 } 00102 00103 SluMatrix& operator=(const SluMatrix& other) 00104 { 00105 SuperMatrix::operator=(static_cast<const SuperMatrix&>(other)); 00106 Store = &storage; 00107 storage = other.storage; 00108 return *this; 00109 } 00110 00111 struct 00112 { 00113 union {int nnz;int lda;}; 00114 void *values; 00115 int *innerInd; 00116 int *outerInd; 00117 } storage; 00118 00119 void setStorageType(Stype_t t) 00120 { 00121 Stype = t; 00122 if (t==SLU_NC || t==SLU_NR || t==SLU_DN) 00123 Store = &storage; 00124 else 00125 { 00126 eigen_assert(false && "storage type not supported"); 00127 Store = 0; 00128 } 00129 } 00130 00131 template<typename Scalar> 00132 void setScalarType() 00133 { 00134 if (internal::is_same<Scalar,float>::value) 00135 Dtype = SLU_S; 00136 else if (internal::is_same<Scalar,double>::value) 00137 Dtype = SLU_D; 00138 else if (internal::is_same<Scalar,std::complex<float> >::value) 00139 Dtype = SLU_C; 00140 else if (internal::is_same<Scalar,std::complex<double> >::value) 00141 Dtype = SLU_Z; 00142 else 00143 { 00144 eigen_assert(false && "Scalar type not supported by SuperLU"); 00145 } 00146 } 00147 00148 template<typename MatrixType> 00149 static SluMatrix Map(MatrixBase<MatrixType>& _mat) 00150 { 00151 MatrixType& mat(_mat.derived()); 00152 eigen_assert( ((MatrixType::Flags&RowMajorBit)!=RowMajorBit) && "row-major dense matrices are not supported by SuperLU"); 00153 SluMatrix res; 00154 res.setStorageType(SLU_DN); 00155 res.setScalarType<typename MatrixType::Scalar>(); 00156 res.Mtype = SLU_GE; 00157 00158 res.nrow = internal::convert_index<int>(mat.rows()); 00159 res.ncol = internal::convert_index<int>(mat.cols()); 00160 00161 res.storage.lda = internal::convert_index<int>(MatrixType::IsVectorAtCompileTime ? mat.size() : mat.outerStride()); 00162 res.storage.values = (void*)(mat.data()); 00163 return res; 00164 } 00165 00166 template<typename MatrixType> 00167 static SluMatrix Map(SparseMatrixBase<MatrixType>& a_mat) 00168 { 00169 MatrixType &mat(a_mat.derived()); 00170 SluMatrix res; 00171 if ((MatrixType::Flags&RowMajorBit)==RowMajorBit) 00172 { 00173 res.setStorageType(SLU_NR); 00174 res.nrow = internal::convert_index<int>(mat.cols()); 00175 res.ncol = internal::convert_index<int>(mat.rows()); 00176 } 00177 else 00178 { 00179 res.setStorageType(SLU_NC); 00180 res.nrow = internal::convert_index<int>(mat.rows()); 00181 res.ncol = internal::convert_index<int>(mat.cols()); 00182 } 00183 00184 res.Mtype = SLU_GE; 00185 00186 res.storage.nnz = internal::convert_index<int>(mat.nonZeros()); 00187 res.storage.values = mat.valuePtr(); 00188 res.storage.innerInd = mat.innerIndexPtr(); 00189 res.storage.outerInd = mat.outerIndexPtr(); 00190 00191 res.setScalarType<typename MatrixType::Scalar>(); 00192 00193 // FIXME the following is not very accurate 00194 if (MatrixType::Flags & Upper) 00195 res.Mtype = SLU_TRU; 00196 if (MatrixType::Flags & Lower) 00197 res.Mtype = SLU_TRL; 00198 00199 eigen_assert(((MatrixType::Flags & SelfAdjoint)==0) && "SelfAdjoint matrix shape not supported by SuperLU"); 00200 00201 return res; 00202 } 00203 }; 00204 00205 template<typename Scalar, int Rows, int Cols, int Options, int MRows, int MCols> 00206 struct SluMatrixMapHelper<Matrix<Scalar,Rows,Cols,Options,MRows,MCols> > 00207 { 00208 typedef Matrix<Scalar,Rows,Cols,Options,MRows,MCols> MatrixType; 00209 static void run(MatrixType& mat, SluMatrix& res) 00210 { 00211 eigen_assert( ((Options&RowMajor)!=RowMajor) && "row-major dense matrices is not supported by SuperLU"); 00212 res.setStorageType(SLU_DN); 00213 res.setScalarType<Scalar>(); 00214 res.Mtype = SLU_GE; 00215 00216 res.nrow = mat.rows(); 00217 res.ncol = mat.cols(); 00218 00219 res.storage.lda = mat.outerStride(); 00220 res.storage.values = mat.data(); 00221 } 00222 }; 00223 00224 template<typename Derived> 00225 struct SluMatrixMapHelper<SparseMatrixBase<Derived> > 00226 { 00227 typedef Derived MatrixType; 00228 static void run(MatrixType& mat, SluMatrix& res) 00229 { 00230 if ((MatrixType::Flags&RowMajorBit)==RowMajorBit) 00231 { 00232 res.setStorageType(SLU_NR); 00233 res.nrow = mat.cols(); 00234 res.ncol = mat.rows(); 00235 } 00236 else 00237 { 00238 res.setStorageType(SLU_NC); 00239 res.nrow = mat.rows(); 00240 res.ncol = mat.cols(); 00241 } 00242 00243 res.Mtype = SLU_GE; 00244 00245 res.storage.nnz = mat.nonZeros(); 00246 res.storage.values = mat.valuePtr(); 00247 res.storage.innerInd = mat.innerIndexPtr(); 00248 res.storage.outerInd = mat.outerIndexPtr(); 00249 00250 res.setScalarType<typename MatrixType::Scalar>(); 00251 00252 // FIXME the following is not very accurate 00253 if (MatrixType::Flags & Upper) 00254 res.Mtype = SLU_TRU; 00255 if (MatrixType::Flags & Lower) 00256 res.Mtype = SLU_TRL; 00257 00258 eigen_assert(((MatrixType::Flags & SelfAdjoint)==0) && "SelfAdjoint matrix shape not supported by SuperLU"); 00259 } 00260 }; 00261 00262 namespace internal { 00263 00264 template<typename MatrixType> 00265 SluMatrix asSluMatrix(MatrixType& mat) 00266 { 00267 return SluMatrix::Map(mat); 00268 } 00269 00271 template<typename Scalar, int Flags, typename Index> 00272 MappedSparseMatrix<Scalar,Flags,Index> map_superlu(SluMatrix& sluMat) 00273 { 00274 eigen_assert((Flags&RowMajor)==RowMajor && sluMat.Stype == SLU_NR 00275 || (Flags&ColMajor)==ColMajor && sluMat.Stype == SLU_NC); 00276 00277 Index outerSize = (Flags&RowMajor)==RowMajor ? sluMat.ncol : sluMat.nrow; 00278 00279 return MappedSparseMatrix<Scalar,Flags,Index>( 00280 sluMat.nrow, sluMat.ncol, sluMat.storage.outerInd[outerSize], 00281 sluMat.storage.outerInd, sluMat.storage.innerInd, reinterpret_cast<Scalar*>(sluMat.storage.values) ); 00282 } 00283 00284 } // end namespace internal 00285 00290 template<typename _MatrixType, typename Derived> 00291 class SuperLUBase : public SparseSolverBase<Derived> 00292 { 00293 protected: 00294 typedef SparseSolverBase<Derived> Base; 00295 using Base::derived; 00296 using Base::m_isInitialized; 00297 public: 00298 typedef _MatrixType MatrixType; 00299 typedef typename MatrixType::Scalar Scalar; 00300 typedef typename MatrixType::RealScalar RealScalar; 00301 typedef typename MatrixType::StorageIndex StorageIndex; 00302 typedef Matrix<Scalar,Dynamic,1> Vector; 00303 typedef Matrix<int, 1, MatrixType::ColsAtCompileTime> IntRowVectorType; 00304 typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> IntColVectorType; 00305 typedef Map<PermutationMatrix<Dynamic,Dynamic,int> > PermutationMap; 00306 typedef SparseMatrix<Scalar> LUMatrixType; 00307 enum { 00308 ColsAtCompileTime = MatrixType::ColsAtCompileTime, 00309 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime 00310 }; 00311 00312 public: 00313 00314 SuperLUBase() {} 00315 00316 ~SuperLUBase() 00317 { 00318 clearFactors(); 00319 } 00320 00321 inline Index rows() const { return m_matrix.rows(); } 00322 inline Index cols() const { return m_matrix.cols(); } 00323 00325 inline superlu_options_t& options() { return m_sluOptions; } 00326 00332 ComputationInfo info() const 00333 { 00334 eigen_assert(m_isInitialized && "Decomposition is not initialized."); 00335 return m_info; 00336 } 00337 00339 void compute(const MatrixType& matrix) 00340 { 00341 derived().analyzePattern(matrix); 00342 derived().factorize(matrix); 00343 } 00344 00351 void analyzePattern(const MatrixType& /*matrix*/) 00352 { 00353 m_isInitialized = true; 00354 m_info = Success; 00355 m_analysisIsOk = true; 00356 m_factorizationIsOk = false; 00357 } 00358 00359 template<typename Stream> 00360 void dumpMemory(Stream& /*s*/) 00361 {} 00362 00363 protected: 00364 00365 void initFactorization(const MatrixType& a) 00366 { 00367 set_default_options(&this->m_sluOptions); 00368 00369 const Index size = a.rows(); 00370 m_matrix = a; 00371 00372 m_sluA = internal::asSluMatrix(m_matrix); 00373 clearFactors(); 00374 00375 m_p.resize(size); 00376 m_q.resize(size); 00377 m_sluRscale.resize(size); 00378 m_sluCscale.resize(size); 00379 m_sluEtree.resize(size); 00380 00381 // set empty B and X 00382 m_sluB.setStorageType(SLU_DN); 00383 m_sluB.setScalarType<Scalar>(); 00384 m_sluB.Mtype = SLU_GE; 00385 m_sluB.storage.values = 0; 00386 m_sluB.nrow = 0; 00387 m_sluB.ncol = 0; 00388 m_sluB.storage.lda = internal::convert_index<int>(size); 00389 m_sluX = m_sluB; 00390 00391 m_extractedDataAreDirty = true; 00392 } 00393 00394 void init() 00395 { 00396 m_info = InvalidInput; 00397 m_isInitialized = false; 00398 m_sluL.Store = 0; 00399 m_sluU.Store = 0; 00400 } 00401 00402 void extractData() const; 00403 00404 void clearFactors() 00405 { 00406 if(m_sluL.Store) 00407 Destroy_SuperNode_Matrix(&m_sluL); 00408 if(m_sluU.Store) 00409 Destroy_CompCol_Matrix(&m_sluU); 00410 00411 m_sluL.Store = 0; 00412 m_sluU.Store = 0; 00413 00414 memset(&m_sluL,0,sizeof m_sluL); 00415 memset(&m_sluU,0,sizeof m_sluU); 00416 } 00417 00418 // cached data to reduce reallocation, etc. 00419 mutable LUMatrixType m_l; 00420 mutable LUMatrixType m_u; 00421 mutable IntColVectorType m_p; 00422 mutable IntRowVectorType m_q; 00423 00424 mutable LUMatrixType m_matrix; // copy of the factorized matrix 00425 mutable SluMatrix m_sluA; 00426 mutable SuperMatrix m_sluL, m_sluU; 00427 mutable SluMatrix m_sluB, m_sluX; 00428 mutable SuperLUStat_t m_sluStat; 00429 mutable superlu_options_t m_sluOptions; 00430 mutable std::vector<int> m_sluEtree; 00431 mutable Matrix<RealScalar,Dynamic,1> m_sluRscale, m_sluCscale; 00432 mutable Matrix<RealScalar,Dynamic,1> m_sluFerr, m_sluBerr; 00433 mutable char m_sluEqued; 00434 00435 mutable ComputationInfo m_info; 00436 int m_factorizationIsOk; 00437 int m_analysisIsOk; 00438 mutable bool m_extractedDataAreDirty; 00439 00440 private: 00441 SuperLUBase(SuperLUBase& ) { } 00442 }; 00443 00444 00461 template<typename _MatrixType> 00462 class SuperLU : public SuperLUBase<_MatrixType,SuperLU<_MatrixType> > 00463 { 00464 public: 00465 typedef SuperLUBase<_MatrixType,SuperLU> Base; 00466 typedef _MatrixType MatrixType; 00467 typedef typename Base::Scalar Scalar; 00468 typedef typename Base::RealScalar RealScalar; 00469 typedef typename Base::StorageIndex StorageIndex; 00470 typedef typename Base::IntRowVectorType IntRowVectorType; 00471 typedef typename Base::IntColVectorType IntColVectorType; 00472 typedef typename Base::PermutationMap PermutationMap; 00473 typedef typename Base::LUMatrixType LUMatrixType; 00474 typedef TriangularView<LUMatrixType, Lower|UnitDiag> LMatrixType; 00475 typedef TriangularView<LUMatrixType, Upper> UMatrixType; 00476 00477 public: 00478 using Base::_solve_impl; 00479 00480 SuperLU() : Base() { init(); } 00481 00482 explicit SuperLU(const MatrixType& matrix) : Base() 00483 { 00484 init(); 00485 Base::compute(matrix); 00486 } 00487 00488 ~SuperLU() 00489 { 00490 } 00491 00498 void analyzePattern(const MatrixType& matrix) 00499 { 00500 m_info = InvalidInput; 00501 m_isInitialized = false; 00502 Base::analyzePattern(matrix); 00503 } 00504 00511 void factorize(const MatrixType& matrix); 00512 00514 template<typename Rhs,typename Dest> 00515 void _solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest> &dest) const; 00516 00517 inline const LMatrixType& matrixL() const 00518 { 00519 if (m_extractedDataAreDirty) this->extractData(); 00520 return m_l; 00521 } 00522 00523 inline const UMatrixType& matrixU() const 00524 { 00525 if (m_extractedDataAreDirty) this->extractData(); 00526 return m_u; 00527 } 00528 00529 inline const IntColVectorType& permutationP() const 00530 { 00531 if (m_extractedDataAreDirty) this->extractData(); 00532 return m_p; 00533 } 00534 00535 inline const IntRowVectorType& permutationQ() const 00536 { 00537 if (m_extractedDataAreDirty) this->extractData(); 00538 return m_q; 00539 } 00540 00541 Scalar determinant() const; 00542 00543 protected: 00544 00545 using Base::m_matrix; 00546 using Base::m_sluOptions; 00547 using Base::m_sluA; 00548 using Base::m_sluB; 00549 using Base::m_sluX; 00550 using Base::m_p; 00551 using Base::m_q; 00552 using Base::m_sluEtree; 00553 using Base::m_sluEqued; 00554 using Base::m_sluRscale; 00555 using Base::m_sluCscale; 00556 using Base::m_sluL; 00557 using Base::m_sluU; 00558 using Base::m_sluStat; 00559 using Base::m_sluFerr; 00560 using Base::m_sluBerr; 00561 using Base::m_l; 00562 using Base::m_u; 00563 00564 using Base::m_analysisIsOk; 00565 using Base::m_factorizationIsOk; 00566 using Base::m_extractedDataAreDirty; 00567 using Base::m_isInitialized; 00568 using Base::m_info; 00569 00570 void init() 00571 { 00572 Base::init(); 00573 00574 set_default_options(&this->m_sluOptions); 00575 m_sluOptions.PrintStat = NO; 00576 m_sluOptions.ConditionNumber = NO; 00577 m_sluOptions.Trans = NOTRANS; 00578 m_sluOptions.ColPerm = COLAMD; 00579 } 00580 00581 00582 private: 00583 SuperLU(SuperLU& ) { } 00584 }; 00585 00586 template<typename MatrixType> 00587 void SuperLU<MatrixType>::factorize(const MatrixType& a) 00588 { 00589 eigen_assert(m_analysisIsOk && "You must first call analyzePattern()"); 00590 if(!m_analysisIsOk) 00591 { 00592 m_info = InvalidInput; 00593 return; 00594 } 00595 00596 this->initFactorization(a); 00597 00598 m_sluOptions.ColPerm = COLAMD; 00599 int info = 0; 00600 RealScalar recip_pivot_growth, rcond; 00601 RealScalar ferr, berr; 00602 00603 StatInit(&m_sluStat); 00604 SuperLU_gssvx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0], 00605 &m_sluEqued, &m_sluRscale[0], &m_sluCscale[0], 00606 &m_sluL, &m_sluU, 00607 NULL, 0, 00608 &m_sluB, &m_sluX, 00609 &recip_pivot_growth, &rcond, 00610 &ferr, &berr, 00611 &m_sluStat, &info, Scalar()); 00612 StatFree(&m_sluStat); 00613 00614 m_extractedDataAreDirty = true; 00615 00616 // FIXME how to better check for errors ??? 00617 m_info = info == 0 ? Success : NumericalIssue; 00618 m_factorizationIsOk = true; 00619 } 00620 00621 template<typename MatrixType> 00622 template<typename Rhs,typename Dest> 00623 void SuperLU<MatrixType>::_solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest>& x) const 00624 { 00625 eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or analyzePattern()/factorize()"); 00626 00627 const Index size = m_matrix.rows(); 00628 const Index rhsCols = b.cols(); 00629 eigen_assert(size==b.rows()); 00630 00631 m_sluOptions.Trans = NOTRANS; 00632 m_sluOptions.Fact = FACTORED; 00633 m_sluOptions.IterRefine = NOREFINE; 00634 00635 00636 m_sluFerr.resize(rhsCols); 00637 m_sluBerr.resize(rhsCols); 00638 00639 Ref<const Matrix<typename Rhs::Scalar,Dynamic,Dynamic,ColMajor> > b_ref(b); 00640 Ref<const Matrix<typename Dest::Scalar,Dynamic,Dynamic,ColMajor> > x_ref(x); 00641 00642 m_sluB = SluMatrix::Map(b_ref.const_cast_derived()); 00643 m_sluX = SluMatrix::Map(x_ref.const_cast_derived()); 00644 00645 typename Rhs::PlainObject b_cpy; 00646 if(m_sluEqued!='N') 00647 { 00648 b_cpy = b; 00649 m_sluB = SluMatrix::Map(b_cpy.const_cast_derived()); 00650 } 00651 00652 StatInit(&m_sluStat); 00653 int info = 0; 00654 RealScalar recip_pivot_growth, rcond; 00655 SuperLU_gssvx(&m_sluOptions, &m_sluA, 00656 m_q.data(), m_p.data(), 00657 &m_sluEtree[0], &m_sluEqued, 00658 &m_sluRscale[0], &m_sluCscale[0], 00659 &m_sluL, &m_sluU, 00660 NULL, 0, 00661 &m_sluB, &m_sluX, 00662 &recip_pivot_growth, &rcond, 00663 &m_sluFerr[0], &m_sluBerr[0], 00664 &m_sluStat, &info, Scalar()); 00665 StatFree(&m_sluStat); 00666 00667 if(x.derived().data() != x_ref.data()) 00668 x = x_ref; 00669 00670 m_info = info==0 ? Success : NumericalIssue; 00671 } 00672 00673 // the code of this extractData() function has been adapted from the SuperLU's Matlab support code, 00674 // 00675 // Copyright (c) 1994 by Xerox Corporation. All rights reserved. 00676 // 00677 // THIS MATERIAL IS PROVIDED AS IS, WITH ABSOLUTELY NO WARRANTY 00678 // EXPRESSED OR IMPLIED. ANY USE IS AT YOUR OWN RISK. 00679 // 00680 template<typename MatrixType, typename Derived> 00681 void SuperLUBase<MatrixType,Derived>::extractData() const 00682 { 00683 eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for extracting factors, you must first call either compute() or analyzePattern()/factorize()"); 00684 if (m_extractedDataAreDirty) 00685 { 00686 int upper; 00687 int fsupc, istart, nsupr; 00688 int lastl = 0, lastu = 0; 00689 SCformat *Lstore = static_cast<SCformat*>(m_sluL.Store); 00690 NCformat *Ustore = static_cast<NCformat*>(m_sluU.Store); 00691 Scalar *SNptr; 00692 00693 const Index size = m_matrix.rows(); 00694 m_l.resize(size,size); 00695 m_l.resizeNonZeros(Lstore->nnz); 00696 m_u.resize(size,size); 00697 m_u.resizeNonZeros(Ustore->nnz); 00698 00699 int* Lcol = m_l.outerIndexPtr(); 00700 int* Lrow = m_l.innerIndexPtr(); 00701 Scalar* Lval = m_l.valuePtr(); 00702 00703 int* Ucol = m_u.outerIndexPtr(); 00704 int* Urow = m_u.innerIndexPtr(); 00705 Scalar* Uval = m_u.valuePtr(); 00706 00707 Ucol[0] = 0; 00708 Ucol[0] = 0; 00709 00710 /* for each supernode */ 00711 for (int k = 0; k <= Lstore->nsuper; ++k) 00712 { 00713 fsupc = L_FST_SUPC(k); 00714 istart = L_SUB_START(fsupc); 00715 nsupr = L_SUB_START(fsupc+1) - istart; 00716 upper = 1; 00717 00718 /* for each column in the supernode */ 00719 for (int j = fsupc; j < L_FST_SUPC(k+1); ++j) 00720 { 00721 SNptr = &((Scalar*)Lstore->nzval)[L_NZ_START(j)]; 00722 00723 /* Extract U */ 00724 for (int i = U_NZ_START(j); i < U_NZ_START(j+1); ++i) 00725 { 00726 Uval[lastu] = ((Scalar*)Ustore->nzval)[i]; 00727 /* Matlab doesn't like explicit zero. */ 00728 if (Uval[lastu] != 0.0) 00729 Urow[lastu++] = U_SUB(i); 00730 } 00731 for (int i = 0; i < upper; ++i) 00732 { 00733 /* upper triangle in the supernode */ 00734 Uval[lastu] = SNptr[i]; 00735 /* Matlab doesn't like explicit zero. */ 00736 if (Uval[lastu] != 0.0) 00737 Urow[lastu++] = L_SUB(istart+i); 00738 } 00739 Ucol[j+1] = lastu; 00740 00741 /* Extract L */ 00742 Lval[lastl] = 1.0; /* unit diagonal */ 00743 Lrow[lastl++] = L_SUB(istart + upper - 1); 00744 for (int i = upper; i < nsupr; ++i) 00745 { 00746 Lval[lastl] = SNptr[i]; 00747 /* Matlab doesn't like explicit zero. */ 00748 if (Lval[lastl] != 0.0) 00749 Lrow[lastl++] = L_SUB(istart+i); 00750 } 00751 Lcol[j+1] = lastl; 00752 00753 ++upper; 00754 } /* for j ... */ 00755 00756 } /* for k ... */ 00757 00758 // squeeze the matrices : 00759 m_l.resizeNonZeros(lastl); 00760 m_u.resizeNonZeros(lastu); 00761 00762 m_extractedDataAreDirty = false; 00763 } 00764 } 00765 00766 template<typename MatrixType> 00767 typename SuperLU<MatrixType>::Scalar SuperLU<MatrixType>::determinant() const 00768 { 00769 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()"); 00770 00771 if (m_extractedDataAreDirty) 00772 this->extractData(); 00773 00774 Scalar det = Scalar(1); 00775 for (int j=0; j<m_u.cols(); ++j) 00776 { 00777 if (m_u.outerIndexPtr()[j+1]-m_u.outerIndexPtr()[j] > 0) 00778 { 00779 int lastId = m_u.outerIndexPtr()[j+1]-1; 00780 eigen_assert(m_u.innerIndexPtr()[lastId]<=j); 00781 if (m_u.innerIndexPtr()[lastId]==j) 00782 det *= m_u.valuePtr()[lastId]; 00783 } 00784 } 00785 if(PermutationMap(m_p.data(),m_p.size()).determinant()*PermutationMap(m_q.data(),m_q.size()).determinant()<0) 00786 det = -det; 00787 if(m_sluEqued!='N') 00788 return det/m_sluRscale.prod()/m_sluCscale.prod(); 00789 else 00790 return det; 00791 } 00792 00793 #ifdef EIGEN_PARSED_BY_DOXYGEN 00794 #define EIGEN_SUPERLU_HAS_ILU 00795 #endif 00796 00797 #ifdef EIGEN_SUPERLU_HAS_ILU 00798 00815 template<typename _MatrixType> 00816 class SuperILU : public SuperLUBase<_MatrixType,SuperILU<_MatrixType> > 00817 { 00818 public: 00819 typedef SuperLUBase<_MatrixType,SuperILU> Base; 00820 typedef _MatrixType MatrixType; 00821 typedef typename Base::Scalar Scalar; 00822 typedef typename Base::RealScalar RealScalar; 00823 00824 public: 00825 using Base::_solve_impl; 00826 00827 SuperILU() : Base() { init(); } 00828 00829 SuperILU(const MatrixType& matrix) : Base() 00830 { 00831 init(); 00832 Base::compute(matrix); 00833 } 00834 00835 ~SuperILU() 00836 { 00837 } 00838 00845 void analyzePattern(const MatrixType& matrix) 00846 { 00847 Base::analyzePattern(matrix); 00848 } 00849 00856 void factorize(const MatrixType& matrix); 00857 00858 #ifndef EIGEN_PARSED_BY_DOXYGEN 00859 00860 template<typename Rhs,typename Dest> 00861 void _solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest> &dest) const; 00862 #endif // EIGEN_PARSED_BY_DOXYGEN 00863 00864 protected: 00865 00866 using Base::m_matrix; 00867 using Base::m_sluOptions; 00868 using Base::m_sluA; 00869 using Base::m_sluB; 00870 using Base::m_sluX; 00871 using Base::m_p; 00872 using Base::m_q; 00873 using Base::m_sluEtree; 00874 using Base::m_sluEqued; 00875 using Base::m_sluRscale; 00876 using Base::m_sluCscale; 00877 using Base::m_sluL; 00878 using Base::m_sluU; 00879 using Base::m_sluStat; 00880 using Base::m_sluFerr; 00881 using Base::m_sluBerr; 00882 using Base::m_l; 00883 using Base::m_u; 00884 00885 using Base::m_analysisIsOk; 00886 using Base::m_factorizationIsOk; 00887 using Base::m_extractedDataAreDirty; 00888 using Base::m_isInitialized; 00889 using Base::m_info; 00890 00891 void init() 00892 { 00893 Base::init(); 00894 00895 ilu_set_default_options(&m_sluOptions); 00896 m_sluOptions.PrintStat = NO; 00897 m_sluOptions.ConditionNumber = NO; 00898 m_sluOptions.Trans = NOTRANS; 00899 m_sluOptions.ColPerm = MMD_AT_PLUS_A; 00900 00901 // no attempt to preserve column sum 00902 m_sluOptions.ILU_MILU = SILU; 00903 // only basic ILU(k) support -- no direct control over memory consumption 00904 // better to use ILU_DropRule = DROP_BASIC | DROP_AREA 00905 // and set ILU_FillFactor to max memory growth 00906 m_sluOptions.ILU_DropRule = DROP_BASIC; 00907 m_sluOptions.ILU_DropTol = NumTraits<Scalar>::dummy_precision()*10; 00908 } 00909 00910 private: 00911 SuperILU(SuperILU& ) { } 00912 }; 00913 00914 template<typename MatrixType> 00915 void SuperILU<MatrixType>::factorize(const MatrixType& a) 00916 { 00917 eigen_assert(m_analysisIsOk && "You must first call analyzePattern()"); 00918 if(!m_analysisIsOk) 00919 { 00920 m_info = InvalidInput; 00921 return; 00922 } 00923 00924 this->initFactorization(a); 00925 00926 int info = 0; 00927 RealScalar recip_pivot_growth, rcond; 00928 00929 StatInit(&m_sluStat); 00930 SuperLU_gsisx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0], 00931 &m_sluEqued, &m_sluRscale[0], &m_sluCscale[0], 00932 &m_sluL, &m_sluU, 00933 NULL, 0, 00934 &m_sluB, &m_sluX, 00935 &recip_pivot_growth, &rcond, 00936 &m_sluStat, &info, Scalar()); 00937 StatFree(&m_sluStat); 00938 00939 // FIXME how to better check for errors ??? 00940 m_info = info == 0 ? Success : NumericalIssue; 00941 m_factorizationIsOk = true; 00942 } 00943 00944 template<typename MatrixType> 00945 template<typename Rhs,typename Dest> 00946 void SuperILU<MatrixType>::_solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest>& x) const 00947 { 00948 eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or analyzePattern()/factorize()"); 00949 00950 const int size = m_matrix.rows(); 00951 const int rhsCols = b.cols(); 00952 eigen_assert(size==b.rows()); 00953 00954 m_sluOptions.Trans = NOTRANS; 00955 m_sluOptions.Fact = FACTORED; 00956 m_sluOptions.IterRefine = NOREFINE; 00957 00958 m_sluFerr.resize(rhsCols); 00959 m_sluBerr.resize(rhsCols); 00960 00961 Ref<const Matrix<typename Rhs::Scalar,Dynamic,Dynamic,ColMajor> > b_ref(b); 00962 Ref<const Matrix<typename Dest::Scalar,Dynamic,Dynamic,ColMajor> > x_ref(x); 00963 00964 m_sluB = SluMatrix::Map(b_ref.const_cast_derived()); 00965 m_sluX = SluMatrix::Map(x_ref.const_cast_derived()); 00966 00967 typename Rhs::PlainObject b_cpy; 00968 if(m_sluEqued!='N') 00969 { 00970 b_cpy = b; 00971 m_sluB = SluMatrix::Map(b_cpy.const_cast_derived()); 00972 } 00973 00974 int info = 0; 00975 RealScalar recip_pivot_growth, rcond; 00976 00977 StatInit(&m_sluStat); 00978 SuperLU_gsisx(&m_sluOptions, &m_sluA, 00979 m_q.data(), m_p.data(), 00980 &m_sluEtree[0], &m_sluEqued, 00981 &m_sluRscale[0], &m_sluCscale[0], 00982 &m_sluL, &m_sluU, 00983 NULL, 0, 00984 &m_sluB, &m_sluX, 00985 &recip_pivot_growth, &rcond, 00986 &m_sluStat, &info, Scalar()); 00987 StatFree(&m_sluStat); 00988 00989 if(&x.coeffRef(0) != x_ref.data()) 00990 x = x_ref; 00991 00992 m_info = info==0 ? Success : NumericalIssue; 00993 } 00994 #endif 00995 00996 } // end namespace Eigen 00997 00998 #endif // EIGEN_SUPERLUSUPPORT_H