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-2010 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_CHOLMODSUPPORT_H 00011 #define EIGEN_CHOLMODSUPPORT_H 00012 00013 namespace Eigen { 00014 00015 namespace internal { 00016 00017 template<typename Scalar, typename CholmodType> 00018 void cholmod_configure_matrix(CholmodType& mat) 00019 { 00020 if (internal::is_same<Scalar,float>::value) 00021 { 00022 mat.xtype = CHOLMOD_REAL; 00023 mat.dtype = CHOLMOD_SINGLE; 00024 } 00025 else if (internal::is_same<Scalar,double>::value) 00026 { 00027 mat.xtype = CHOLMOD_REAL; 00028 mat.dtype = CHOLMOD_DOUBLE; 00029 } 00030 else if (internal::is_same<Scalar,std::complex<float> >::value) 00031 { 00032 mat.xtype = CHOLMOD_COMPLEX; 00033 mat.dtype = CHOLMOD_SINGLE; 00034 } 00035 else if (internal::is_same<Scalar,std::complex<double> >::value) 00036 { 00037 mat.xtype = CHOLMOD_COMPLEX; 00038 mat.dtype = CHOLMOD_DOUBLE; 00039 } 00040 else 00041 { 00042 eigen_assert(false && "Scalar type not supported by CHOLMOD"); 00043 } 00044 } 00045 00046 } // namespace internal 00047 00051 template<typename _Scalar, int _Options, typename _StorageIndex> 00052 cholmod_sparse viewAsCholmod(SparseMatrix<_Scalar,_Options,_StorageIndex>& mat) 00053 { 00054 cholmod_sparse res; 00055 res.nzmax = mat.nonZeros(); 00056 res.nrow = mat.rows();; 00057 res.ncol = mat.cols(); 00058 res.p = mat.outerIndexPtr(); 00059 res.i = mat.innerIndexPtr(); 00060 res.x = mat.valuePtr(); 00061 res.z = 0; 00062 res.sorted = 1; 00063 if(mat.isCompressed()) 00064 { 00065 res.packed = 1; 00066 res.nz = 0; 00067 } 00068 else 00069 { 00070 res.packed = 0; 00071 res.nz = mat.innerNonZeroPtr(); 00072 } 00073 00074 res.dtype = 0; 00075 res.stype = -1; 00076 00077 if (internal::is_same<_StorageIndex,int>::value) 00078 { 00079 res.itype = CHOLMOD_INT; 00080 } 00081 else if (internal::is_same<_StorageIndex,long>::value) 00082 { 00083 res.itype = CHOLMOD_LONG; 00084 } 00085 else 00086 { 00087 eigen_assert(false && "Index type not supported yet"); 00088 } 00089 00090 // setup res.xtype 00091 internal::cholmod_configure_matrix<_Scalar>(res); 00092 00093 res.stype = 0; 00094 00095 return res; 00096 } 00097 00098 template<typename _Scalar, int _Options, typename _Index> 00099 const cholmod_sparse viewAsCholmod(const SparseMatrix<_Scalar,_Options,_Index>& mat) 00100 { 00101 cholmod_sparse res = viewAsCholmod(mat.const_cast_derived()); 00102 return res; 00103 } 00104 00107 template<typename _Scalar, int _Options, typename _Index, unsigned int UpLo> 00108 cholmod_sparse viewAsCholmod(const SparseSelfAdjointView<const SparseMatrix<_Scalar,_Options,_Index>, UpLo>& mat) 00109 { 00110 cholmod_sparse res = viewAsCholmod(mat.matrix().const_cast_derived()); 00111 00112 if(UpLo==Upper) res.stype = 1; 00113 if(UpLo==Lower) res.stype = -1; 00114 00115 return res; 00116 } 00117 00120 template<typename Derived> 00121 cholmod_dense viewAsCholmod(MatrixBase<Derived>& mat) 00122 { 00123 EIGEN_STATIC_ASSERT((internal::traits<Derived>::Flags&RowMajorBit)==0,THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES); 00124 typedef typename Derived::Scalar Scalar; 00125 00126 cholmod_dense res; 00127 res.nrow = mat.rows(); 00128 res.ncol = mat.cols(); 00129 res.nzmax = res.nrow * res.ncol; 00130 res.d = Derived::IsVectorAtCompileTime ? mat.derived().size() : mat.derived().outerStride(); 00131 res.x = (void*)(mat.derived().data()); 00132 res.z = 0; 00133 00134 internal::cholmod_configure_matrix<Scalar>(res); 00135 00136 return res; 00137 } 00138 00141 template<typename Scalar, int Flags, typename StorageIndex> 00142 MappedSparseMatrix<Scalar,Flags,StorageIndex> viewAsEigen(cholmod_sparse& cm) 00143 { 00144 return MappedSparseMatrix<Scalar,Flags,StorageIndex> 00145 (cm.nrow, cm.ncol, static_cast<StorageIndex*>(cm.p)[cm.ncol], 00146 static_cast<StorageIndex*>(cm.p), static_cast<StorageIndex*>(cm.i),static_cast<Scalar*>(cm.x) ); 00147 } 00148 00149 enum CholmodMode { 00150 CholmodAuto, CholmodSimplicialLLt, CholmodSupernodalLLt, CholmodLDLt 00151 }; 00152 00153 00159 template<typename _MatrixType, int _UpLo, typename Derived> 00160 class CholmodBase : public SparseSolverBase<Derived> 00161 { 00162 protected: 00163 typedef SparseSolverBase<Derived> Base; 00164 using Base::derived; 00165 using Base::m_isInitialized; 00166 public: 00167 typedef _MatrixType MatrixType; 00168 enum { UpLo = _UpLo }; 00169 typedef typename MatrixType::Scalar Scalar; 00170 typedef typename MatrixType::RealScalar RealScalar; 00171 typedef MatrixType CholMatrixType; 00172 typedef typename MatrixType::StorageIndex StorageIndex; 00173 enum { 00174 ColsAtCompileTime = MatrixType::ColsAtCompileTime, 00175 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime 00176 }; 00177 00178 public: 00179 00180 CholmodBase() 00181 : m_cholmodFactor(0), m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false) 00182 { 00183 m_shiftOffset[0] = m_shiftOffset[1] = RealScalar(0.0); 00184 cholmod_start(&m_cholmod); 00185 } 00186 00187 explicit CholmodBase(const MatrixType& matrix) 00188 : m_cholmodFactor(0), m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false) 00189 { 00190 m_shiftOffset[0] = m_shiftOffset[1] = RealScalar(0.0); 00191 cholmod_start(&m_cholmod); 00192 compute(matrix); 00193 } 00194 00195 ~CholmodBase() 00196 { 00197 if(m_cholmodFactor) 00198 cholmod_free_factor(&m_cholmodFactor, &m_cholmod); 00199 cholmod_finish(&m_cholmod); 00200 } 00201 00202 inline StorageIndex cols() const { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); } 00203 inline StorageIndex rows() const { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); } 00204 00210 ComputationInfo info() const 00211 { 00212 eigen_assert(m_isInitialized && "Decomposition is not initialized."); 00213 return m_info; 00214 } 00215 00217 Derived& compute(const MatrixType& matrix) 00218 { 00219 analyzePattern(matrix); 00220 factorize(matrix); 00221 return derived(); 00222 } 00223 00230 void analyzePattern(const MatrixType& matrix) 00231 { 00232 if(m_cholmodFactor) 00233 { 00234 cholmod_free_factor(&m_cholmodFactor, &m_cholmod); 00235 m_cholmodFactor = 0; 00236 } 00237 cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>()); 00238 m_cholmodFactor = cholmod_analyze(&A, &m_cholmod); 00239 00240 this->m_isInitialized = true; 00241 this->m_info = Success; 00242 m_analysisIsOk = true; 00243 m_factorizationIsOk = false; 00244 } 00245 00252 void factorize(const MatrixType& matrix) 00253 { 00254 eigen_assert(m_analysisIsOk && "You must first call analyzePattern()"); 00255 cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>()); 00256 cholmod_factorize_p(&A, m_shiftOffset, 0, 0, m_cholmodFactor, &m_cholmod); 00257 00258 // If the factorization failed, minor is the column at which it did. On success minor == n. 00259 this->m_info = (m_cholmodFactor->minor == m_cholmodFactor->n ? Success : NumericalIssue); 00260 m_factorizationIsOk = true; 00261 } 00262 00265 cholmod_common& cholmod() { return m_cholmod; } 00266 00267 #ifndef EIGEN_PARSED_BY_DOXYGEN 00268 00269 template<typename Rhs,typename Dest> 00270 void _solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest> &dest) const 00271 { 00272 eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()"); 00273 const Index size = m_cholmodFactor->n; 00274 EIGEN_UNUSED_VARIABLE(size); 00275 eigen_assert(size==b.rows()); 00276 00277 // Cholmod needs column-major stoarge without inner-stride, which corresponds to the default behavior of Ref. 00278 Ref<const Matrix<typename Rhs::Scalar,Dynamic,Dynamic,ColMajor> > b_ref(b.derived()); 00279 00280 cholmod_dense b_cd = viewAsCholmod(b_ref); 00281 cholmod_dense* x_cd = cholmod_solve(CHOLMOD_A, m_cholmodFactor, &b_cd, &m_cholmod); 00282 if(!x_cd) 00283 { 00284 this->m_info = NumericalIssue; 00285 return; 00286 } 00287 // TODO optimize this copy by swapping when possible (be careful with alignment, etc.) 00288 dest = Matrix<Scalar,Dest::RowsAtCompileTime,Dest::ColsAtCompileTime>::Map(reinterpret_cast<Scalar*>(x_cd->x),b.rows(),b.cols()); 00289 cholmod_free_dense(&x_cd, &m_cholmod); 00290 } 00291 00293 template<typename RhsScalar, int RhsOptions, typename RhsIndex, typename DestScalar, int DestOptions, typename DestIndex> 00294 void _solve_impl(const SparseMatrix<RhsScalar,RhsOptions,RhsIndex> &b, SparseMatrix<DestScalar,DestOptions,DestIndex> &dest) const 00295 { 00296 eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()"); 00297 const Index size = m_cholmodFactor->n; 00298 EIGEN_UNUSED_VARIABLE(size); 00299 eigen_assert(size==b.rows()); 00300 00301 // note: cs stands for Cholmod Sparse 00302 cholmod_sparse b_cs = viewAsCholmod(b); 00303 cholmod_sparse* x_cs = cholmod_spsolve(CHOLMOD_A, m_cholmodFactor, &b_cs, &m_cholmod); 00304 if(!x_cs) 00305 { 00306 this->m_info = NumericalIssue; 00307 return; 00308 } 00309 // TODO optimize this copy by swapping when possible (be careful with alignment, etc.) 00310 dest = viewAsEigen<DestScalar,DestOptions,DestIndex>(*x_cs); 00311 cholmod_free_sparse(&x_cs, &m_cholmod); 00312 } 00313 #endif // EIGEN_PARSED_BY_DOXYGEN 00314 00315 00325 Derived& setShift(const RealScalar& offset) 00326 { 00327 m_shiftOffset[0] = offset; 00328 return derived(); 00329 } 00330 00332 Scalar determinant() const 00333 { 00334 using std::exp; 00335 return exp(logDeterminant()); 00336 } 00337 00339 Scalar logDeterminant() const 00340 { 00341 using std::log; 00342 using numext::real; 00343 eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()"); 00344 00345 RealScalar logDet = 0; 00346 Scalar *x = static_cast<Scalar*>(m_cholmodFactor->x); 00347 if (m_cholmodFactor->is_super) 00348 { 00349 // Supernodal factorization stored as a packed list of dense column-major blocs, 00350 // as described by the following structure: 00351 00352 // super[k] == index of the first column of the j-th super node 00353 StorageIndex *super = static_cast<StorageIndex*>(m_cholmodFactor->super); 00354 // pi[k] == offset to the description of row indices 00355 StorageIndex *pi = static_cast<StorageIndex*>(m_cholmodFactor->pi); 00356 // px[k] == offset to the respective dense block 00357 StorageIndex *px = static_cast<StorageIndex*>(m_cholmodFactor->px); 00358 00359 Index nb_super_nodes = m_cholmodFactor->nsuper; 00360 for (Index k=0; k < nb_super_nodes; ++k) 00361 { 00362 StorageIndex ncols = super[k + 1] - super[k]; 00363 StorageIndex nrows = pi[k + 1] - pi[k]; 00364 00365 Map<const Array<Scalar,1,Dynamic>, 0, InnerStride<> > sk(x + px[k], ncols, InnerStride<>(nrows+1)); 00366 logDet += sk.real().log().sum(); 00367 } 00368 } 00369 else 00370 { 00371 // Simplicial factorization stored as standard CSC matrix. 00372 StorageIndex *p = static_cast<StorageIndex*>(m_cholmodFactor->p); 00373 Index size = m_cholmodFactor->n; 00374 for (Index k=0; k<size; ++k) 00375 logDet += log(real( x[p[k]] )); 00376 } 00377 if (m_cholmodFactor->is_ll) 00378 logDet *= 2.0; 00379 return logDet; 00380 }; 00381 00382 template<typename Stream> 00383 void dumpMemory(Stream& /*s*/) 00384 {} 00385 00386 protected: 00387 mutable cholmod_common m_cholmod; 00388 cholmod_factor* m_cholmodFactor; 00389 RealScalar m_shiftOffset[2]; 00390 mutable ComputationInfo m_info; 00391 int m_factorizationIsOk; 00392 int m_analysisIsOk; 00393 }; 00394 00415 template<typename _MatrixType, int _UpLo = Lower> 00416 class CholmodSimplicialLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLLT<_MatrixType, _UpLo> > 00417 { 00418 typedef CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLLT> Base; 00419 using Base::m_cholmod; 00420 00421 public: 00422 00423 typedef _MatrixType MatrixType; 00424 00425 CholmodSimplicialLLT() : Base() { init(); } 00426 00427 CholmodSimplicialLLT(const MatrixType& matrix) : Base() 00428 { 00429 init(); 00430 this->compute(matrix); 00431 } 00432 00433 ~CholmodSimplicialLLT() {} 00434 protected: 00435 void init() 00436 { 00437 m_cholmod.final_asis = 0; 00438 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL; 00439 m_cholmod.final_ll = 1; 00440 } 00441 }; 00442 00443 00464 template<typename _MatrixType, int _UpLo = Lower> 00465 class CholmodSimplicialLDLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLDLT<_MatrixType, _UpLo> > 00466 { 00467 typedef CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLDLT> Base; 00468 using Base::m_cholmod; 00469 00470 public: 00471 00472 typedef _MatrixType MatrixType; 00473 00474 CholmodSimplicialLDLT() : Base() { init(); } 00475 00476 CholmodSimplicialLDLT(const MatrixType& matrix) : Base() 00477 { 00478 init(); 00479 this->compute(matrix); 00480 } 00481 00482 ~CholmodSimplicialLDLT() {} 00483 protected: 00484 void init() 00485 { 00486 m_cholmod.final_asis = 1; 00487 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL; 00488 } 00489 }; 00490 00511 template<typename _MatrixType, int _UpLo = Lower> 00512 class CholmodSupernodalLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSupernodalLLT<_MatrixType, _UpLo> > 00513 { 00514 typedef CholmodBase<_MatrixType, _UpLo, CholmodSupernodalLLT> Base; 00515 using Base::m_cholmod; 00516 00517 public: 00518 00519 typedef _MatrixType MatrixType; 00520 00521 CholmodSupernodalLLT() : Base() { init(); } 00522 00523 CholmodSupernodalLLT(const MatrixType& matrix) : Base() 00524 { 00525 init(); 00526 this->compute(matrix); 00527 } 00528 00529 ~CholmodSupernodalLLT() {} 00530 protected: 00531 void init() 00532 { 00533 m_cholmod.final_asis = 1; 00534 m_cholmod.supernodal = CHOLMOD_SUPERNODAL; 00535 } 00536 }; 00537 00560 template<typename _MatrixType, int _UpLo = Lower> 00561 class CholmodDecomposition : public CholmodBase<_MatrixType, _UpLo, CholmodDecomposition<_MatrixType, _UpLo> > 00562 { 00563 typedef CholmodBase<_MatrixType, _UpLo, CholmodDecomposition> Base; 00564 using Base::m_cholmod; 00565 00566 public: 00567 00568 typedef _MatrixType MatrixType; 00569 00570 CholmodDecomposition() : Base() { init(); } 00571 00572 CholmodDecomposition(const MatrixType& matrix) : Base() 00573 { 00574 init(); 00575 this->compute(matrix); 00576 } 00577 00578 ~CholmodDecomposition() {} 00579 00580 void setMode(CholmodMode mode) 00581 { 00582 switch(mode) 00583 { 00584 case CholmodAuto: 00585 m_cholmod.final_asis = 1; 00586 m_cholmod.supernodal = CHOLMOD_AUTO; 00587 break; 00588 case CholmodSimplicialLLt: 00589 m_cholmod.final_asis = 0; 00590 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL; 00591 m_cholmod.final_ll = 1; 00592 break; 00593 case CholmodSupernodalLLt: 00594 m_cholmod.final_asis = 1; 00595 m_cholmod.supernodal = CHOLMOD_SUPERNODAL; 00596 break; 00597 case CholmodLDLt: 00598 m_cholmod.final_asis = 1; 00599 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL; 00600 break; 00601 default: 00602 break; 00603 } 00604 } 00605 protected: 00606 void init() 00607 { 00608 m_cholmod.final_asis = 1; 00609 m_cholmod.supernodal = CHOLMOD_AUTO; 00610 } 00611 }; 00612 00613 } // end namespace Eigen 00614 00615 #endif // EIGEN_CHOLMODSUPPORT_H