- Home
- About DUNE
- Download
- Documentation
- Community
- Development
00001 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*- 00002 // vi: set et ts=4 sw=2 sts=2: 00003 #ifndef DUNE_ISTL_ILUSUBDOMAIN_HH 00004 #define DUNE_ISTL_ILUSUBDOMAIN_HH 00005 00006 #include<map> 00007 #include<dune/common/typetraits.hh> 00008 #include"matrix.hh" 00009 #include<cmath> 00010 #include<cstdlib> 00011 00012 namespace Dune{ 00013 00032 template<class M, class X, class Y> 00033 class ILUSubdomainSolver { 00034 public: 00036 typedef typename Dune::remove_const<M>::type matrix_type; 00038 typedef X domain_type; 00040 typedef Y range_type; 00041 00048 virtual void apply (X& v, const Y& d) =0; 00049 00050 virtual ~ILUSubdomainSolver() 00051 {} 00052 00053 protected: 00059 template<class S> 00060 std::size_t copyToLocalMatrix(const M& A, S& rowset); 00061 00063 // for ILUN 00064 matrix_type ILU; 00065 }; 00066 00073 template<class M, class X, class Y> 00074 class ILU0SubdomainSolver 00075 : public ILUSubdomainSolver<M,X,Y>{ 00076 public: 00078 typedef typename Dune::remove_const<M>::type matrix_type; 00079 typedef typename Dune::remove_const<M>::type rilu_type; 00081 typedef X domain_type; 00083 typedef Y range_type; 00084 00085 00090 void apply (X& v, const Y& d) 00091 { 00092 bilu_backsolve(this->ILU,v,d); 00093 } 00101 template<class S> 00102 void setSubMatrix(const M& A, S& rowset); 00103 00104 }; 00105 00106 template<class M, class X, class Y> 00107 class ILUNSubdomainSolver 00108 : public ILUSubdomainSolver<M,X,Y>{ 00109 public: 00111 typedef typename Dune::remove_const<M>::type matrix_type; 00112 typedef typename Dune::remove_const<M>::type rilu_type; 00114 typedef X domain_type; 00116 typedef Y range_type; 00117 00122 void apply (X& v, const Y& d) 00123 { 00124 bilu_backsolve(RILU,v,d); 00125 } 00126 00134 template<class S> 00135 void setSubMatrix(const M& A, S& rowset); 00136 00137 private: 00141 rilu_type RILU; 00142 }; 00143 00144 00145 00146 template<class M, class X, class Y> 00147 template<class S> 00148 std::size_t ILUSubdomainSolver<M,X,Y>::copyToLocalMatrix(const M& A, S& rowSet) 00149 { 00150 // Calculate consecutive indices for local problem 00151 // while perserving the ordering 00152 typedef typename M::size_type size_type; 00153 typedef std::map<typename S::value_type,size_type> IndexMap; 00154 typedef typename IndexMap::iterator IMIter; 00155 IndexMap indexMap; 00156 IMIter guess = indexMap.begin(); 00157 size_type localIndex=0; 00158 00159 typedef typename S::const_iterator SIter; 00160 for(SIter rowIdx = rowSet.begin(), rowEnd=rowSet.end(); 00161 rowIdx!= rowEnd; ++rowIdx, ++localIndex) 00162 guess = indexMap.insert(guess, 00163 std::make_pair(*rowIdx,localIndex)); 00164 00165 00166 // Build Matrix for local subproblem 00167 ILU.setSize(rowSet.size(),rowSet.size()); 00168 ILU.setBuildMode(matrix_type::row_wise); 00169 00170 // Create sparsity pattern 00171 typedef typename matrix_type::CreateIterator CIter; 00172 CIter rowCreator = ILU.createbegin(); 00173 typedef typename S::const_iterator RIter; 00174 std::size_t offset=0; 00175 for(SIter rowIdx = rowSet.begin(), rowEnd=rowSet.end(); 00176 rowIdx!= rowEnd; ++rowIdx, ++rowCreator){ 00177 // See wich row entries are in our subset and add them to 00178 // the sparsity pattern 00179 guess = indexMap.begin(); 00180 00181 for(typename matrix_type::ConstColIterator col=A[*rowIdx].begin(), 00182 endcol=A[*rowIdx].end(); col != endcol; ++col){ 00183 // search for the entry in the row set 00184 guess = indexMap.find(col.index()); 00185 if(guess!=indexMap.end()){ 00186 // add local index to row 00187 rowCreator.insert(guess->second); 00188 offset=std::max(offset,(std::size_t)std::abs((int)(guess->second-rowCreator.index()))); 00189 } 00190 } 00191 00192 } 00193 00194 // Insert the matrix values for the local problem 00195 typename matrix_type::iterator iluRow=ILU.begin(); 00196 00197 for(SIter rowIdx = rowSet.begin(), rowEnd=rowSet.end(); 00198 rowIdx!= rowEnd; ++rowIdx, ++iluRow){ 00199 // See wich row entries are in our subset and add them to 00200 // the sparsity pattern 00201 typename matrix_type::ColIterator localCol=iluRow->begin(); 00202 for(typename matrix_type::ConstColIterator col=A[*rowIdx].begin(), 00203 endcol=A[*rowIdx].end(); col != endcol; ++col){ 00204 // search for the entry in the row set 00205 guess = indexMap.find(col.index()); 00206 if(guess!=indexMap.end()){ 00207 // set local value 00208 (*localCol)=(*col); 00209 ++localCol; 00210 } 00211 } 00212 } 00213 return offset; 00214 } 00215 00216 00217 template<class M, class X, class Y> 00218 template<class S> 00219 void ILU0SubdomainSolver<M,X,Y>::setSubMatrix(const M& A, S& rowSet) 00220 { 00221 this->copyToLocalMatrix(A,rowSet); 00222 bilu0_decomposition(this->ILU); 00223 } 00224 00225 template<class M, class X, class Y> 00226 template<class S> 00227 void ILUNSubdomainSolver<M,X,Y>::setSubMatrix(const M& A, S& rowSet) 00228 { 00229 std::size_t offset=copyToLocalMatrix(A,rowSet); 00230 RILU.setSize(rowSet.size(),rowSet.size(), (1+2*offset)*rowSet.size()); 00231 RILU.setBuildMode(matrix_type::row_wise); 00232 bilu_decomposition(this->ILU, (offset+1)/2, RILU); 00233 } 00234 00236 } // end name space DUNE 00237 00238 00239 #endif
Generated on Fri Apr 29 2011 with Doxygen (ver 1.7.1) [doxygen-log,error-log].