dune-istl
2.1.1
|
00001 #ifndef DUNE_MATRIXMARKET_HH 00002 #define DUNE_MATRIXMARKET_HH 00003 00004 #include<ostream> 00005 #include<limits> 00006 #include "bcrsmatrix.hh" 00007 #include<dune/common/fmatrix.hh> 00008 #include<dune/common/tuples.hh> 00009 00010 namespace Dune 00011 { 00012 namespace 00013 { 00023 template<class T> 00024 struct mm_numeric_type; 00025 00026 template<> 00027 struct mm_numeric_type<int> 00028 { 00029 static std::string str() 00030 { 00031 return "integer"; 00032 } 00033 }; 00034 00035 template<> 00036 struct mm_numeric_type<double> 00037 { 00038 static std::string str() 00039 { 00040 return "real"; 00041 } 00042 }; 00043 00044 template<> 00045 struct mm_numeric_type<float> 00046 { 00047 static std::string str() 00048 { 00049 return "real"; 00050 } 00051 }; 00052 00053 template<> 00054 struct mm_numeric_type<std::complex<double> > 00055 { 00056 static std::string str() 00057 { 00058 return "complex"; 00059 } 00060 }; 00061 00062 template<> 00063 struct mm_numeric_type<std::complex<float> > 00064 { 00065 static std::string str() 00066 { 00067 return "complex"; 00068 } 00069 }; 00070 00079 template<class M> 00080 struct mm_header_printer; 00081 00082 template<typename T, typename A, int i, int j> 00083 struct mm_header_printer<BCRSMatrix<FieldMatrix<T,i,j>,A> > 00084 { 00085 static void print(std::ostream& os) 00086 { 00087 os<<"%%MatrixMarket matrix coordinate "; 00088 os<<mm_numeric_type<T>::str()<<" general"<<std::endl; 00089 } 00090 }; 00091 00092 template<typename T, int i, int j> 00093 struct mm_header_printer<FieldMatrix<T,i,j> > 00094 { 00095 static void print(std::ostream& os) 00096 { 00097 os<<"%%MatrixMarket matrix array "; 00098 os<<mm_numeric_type<T>::str()<<" general"<<std::endl; 00099 } 00100 }; 00101 00110 template<class M> 00111 struct mm_block_structure_header; 00112 00113 template<typename T, typename A, int i, int j> 00114 struct mm_block_structure_header<BCRSMatrix<FieldMatrix<T,i,j>,A> > 00115 { 00116 typedef BCRSMatrix<FieldMatrix<T,i,j>,A> M; 00117 00118 static void print(std::ostream& os, const M& m) 00119 { 00120 os<<"% ISTL_STRUCT blocked "; 00121 os<<i<<" "<<j<<std::endl; 00122 } 00123 }; 00124 00125 template<typename T, int i, int j> 00126 struct mm_block_structure_header<FieldMatrix<T,i,j> > 00127 { 00128 typedef FieldMatrix<T,i,j> M; 00129 00130 static void print(std::ostream& os, const M& m) 00131 {} 00132 }; 00133 00134 00135 enum LineType{ MM_HEADER, MM_ISTLSTRUCT, DATA }; 00136 enum{ MM_MAX_LINE_LENGTH=1025 }; 00137 00138 enum MM_TYPE { coordinate_type, array_type, unknown_type }; 00139 00140 enum MM_CTYPE { integer_type, double_type, complex_type, pattern, unknown_ctype }; 00141 00142 enum MM_STRUCTURE { general, symmetric, skew_symmetric, hermitian, unknown_structure }; 00143 00144 struct MMHeader 00145 { 00146 MMHeader() 00147 : type(coordinate_type), ctype(double_type), structure(general) 00148 {} 00149 MM_TYPE type; 00150 MM_CTYPE ctype; 00151 MM_STRUCTURE structure; 00152 }; 00153 00154 bool lineFeed(std::istream& file) 00155 { 00156 char c=file.peek(); 00157 // ignore whitespace 00158 while(c==' ') 00159 { 00160 file.get(); 00161 c=file.peek(); 00162 } 00163 00164 if(c=='\n'){ 00165 /* eat the line feed */ 00166 file.get(); 00167 return true; 00168 } 00169 return false; 00170 } 00171 00172 void skipComments(std::istream& file) 00173 { 00174 lineFeed(file); 00175 char c=file.peek(); 00176 // ignore comment lines 00177 while(c=='%') 00178 { 00179 /* disgard the rest of the line */ 00180 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00181 c=file.peek(); 00182 } 00183 } 00184 00185 00186 bool readMatrixMarketBanner(std::istream& file, MMHeader& mmHeader) 00187 { 00188 std::string buffer; 00189 char c; 00190 c=file.peek(); 00191 mmHeader=MMHeader(); 00192 if(c!='%') 00193 return false; 00194 00195 /* read the banner */ 00196 file >> buffer; 00197 if(buffer!="%%MatrixMarket"){ 00198 /* disgard the rest of the line */ 00199 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00200 return false; 00201 } 00202 00203 if(lineFeed(file)) 00204 /* premature end of line */ 00205 return false; 00206 00207 /* read the matrix_type */ 00208 file >> buffer; 00209 00210 if(buffer != "matrix") 00211 { 00212 /* disgard the rest of the line */ 00213 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00214 return false; 00215 } 00216 00217 if(lineFeed(file)) 00218 /* premature end of line */ 00219 return false; 00220 00221 /* The type of the matrix */ 00222 file >> buffer; 00223 00224 if(buffer.empty()) 00225 return false; 00226 00227 std::transform(buffer.begin(), buffer.end(), buffer.begin(), 00228 tolower); 00229 00230 switch(buffer[0]) 00231 { 00232 case 'a': 00233 /* sanity check */ 00234 if(buffer != "array") 00235 { 00236 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00237 return false; 00238 } 00239 mmHeader.type=array_type; 00240 break; 00241 case 'c': 00242 /* sanity check */ 00243 if(buffer != "coordinate") 00244 { 00245 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00246 return false; 00247 } 00248 mmHeader.type=coordinate_type; 00249 break; 00250 default: 00251 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00252 return false; 00253 } 00254 00255 if(lineFeed(file)) 00256 /* premature end of line */ 00257 return false; 00258 00259 /* The numeric type used. */ 00260 file >> buffer; 00261 00262 if(buffer.empty()) 00263 return false; 00264 00265 std::transform(buffer.begin(), buffer.end(), buffer.begin(), 00266 tolower); 00267 switch(buffer[0]) 00268 { 00269 case 'i': 00270 /* sanity check */ 00271 if(buffer != "integer") 00272 { 00273 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00274 return false; 00275 } 00276 mmHeader.ctype=integer_type; 00277 break; 00278 case 'r': 00279 /* sanity check */ 00280 if(buffer != "real") 00281 { 00282 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00283 return false; 00284 } 00285 mmHeader.ctype=double_type; 00286 break; 00287 case 'c': 00288 /* sanity check */ 00289 if(buffer != "complex") 00290 { 00291 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00292 return false; 00293 } 00294 mmHeader.ctype=complex_type; 00295 break; 00296 case 'p': 00297 /* sanity check */ 00298 if(buffer != "pattern") 00299 { 00300 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00301 return false; 00302 } 00303 mmHeader.ctype=pattern; 00304 break; 00305 default: 00306 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00307 return false; 00308 } 00309 00310 if(lineFeed(file)) 00311 return false; 00312 00313 file >> buffer; 00314 00315 std::transform(buffer.begin(), buffer.end(), buffer.begin(), 00316 tolower); 00317 switch(buffer[0]) 00318 { 00319 case 'g': 00320 /* sanity check */ 00321 if(buffer != "general") 00322 { 00323 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00324 return false; 00325 } 00326 mmHeader.structure=general; 00327 break; 00328 case 'h': 00329 /* sanity check */ 00330 if(buffer != "hermitian") 00331 { 00332 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00333 return false; 00334 } 00335 mmHeader.structure=hermitian; 00336 break; 00337 case 's': 00338 if(buffer.size()==1){ 00339 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00340 return false; 00341 } 00342 00343 switch(buffer[1]) 00344 { 00345 case 'y': 00346 /* sanity check */ 00347 if(buffer != "symmetric") 00348 { 00349 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00350 return false; 00351 } 00352 mmHeader.structure=symmetric; 00353 break; 00354 case 'k': 00355 /* sanity check */ 00356 if(buffer != "skew-symmetric") 00357 { 00358 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00359 return false; 00360 } 00361 mmHeader.structure=skew_symmetric; 00362 break; 00363 default: 00364 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00365 return false; 00366 } 00367 default: 00368 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00369 return false; 00370 } 00371 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00372 c=file.peek(); 00373 return true; 00374 00375 } 00376 00377 void readNextLine(std::istream& file, std::ostringstream& line, LineType& type) 00378 { 00379 char c; 00380 std::size_t index=0; 00381 00382 //empty lines will be disgarded and we will simply read the next line 00383 while(index==0&&!file.eof()) 00384 { 00385 // strip spaces 00386 while(!file.eof() && (c=file.get())==' '); 00387 00388 //read the rest of the line until comment 00389 while(!file.eof() && (c=file.get())=='\n'){ 00390 switch(c) 00391 { 00392 case '%': 00393 /* disgard the rest of the line */ 00394 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00395 } 00396 } 00397 } 00398 00399 // buffer[index]='\0'; 00400 } 00401 00402 template<std::size_t brows, std::size_t bcols> 00403 Dune::tuple<std::size_t, std::size_t, std::size_t> 00404 calculateNNZ(std::size_t rows, std::size_t cols, std::size_t entries, const MMHeader& header) 00405 { 00406 std::size_t blockrows=rows/brows; 00407 std::size_t blockcols=cols/bcols; 00408 std::size_t blocksize=brows*bcols; 00409 std::size_t blockentries=0; 00410 00411 switch(header.structure) 00412 { 00413 case general: 00414 blockentries = entries/blocksize; break; 00415 case skew_symmetric: 00416 blockentries = 2*entries/blocksize; break; 00417 case symmetric: 00418 blockentries = (2*entries-rows)/blocksize; break; 00419 case hermitian: 00420 blockentries = (2*entries-rows)/blocksize; break; 00421 default: 00422 throw Dune::NotImplemented(); 00423 } 00424 return Dune::make_tuple(blockrows, blockcols, blockentries); 00425 } 00426 00427 /* 00428 * @brief Storage class for the column index and the numeric value. 00429 * 00430 * \tparam T Either a NumericWrapper of the numeric type or PatternDummy 00431 * for MatrixMarket pattern case. 00432 */ 00433 template<typename T> 00434 struct IndexData : public T 00435 { 00436 std::size_t index; 00437 }; 00438 00439 00450 template<typename T> 00451 struct NumericWrapper 00452 { 00453 T number; 00454 operator T&() 00455 { 00456 return number; 00457 } 00458 }; 00459 00463 struct PatternDummy 00464 {}; 00465 00466 template<> 00467 struct NumericWrapper<PatternDummy> 00468 {}; 00469 00470 template<typename T> 00471 std::istream& operator>>(std::istream& is, NumericWrapper<T>& num) 00472 { 00473 return is>>num.number; 00474 } 00475 00476 std::istream& operator>>(std::istream& is, NumericWrapper<PatternDummy>& num) 00477 { 00478 return is; 00479 } 00480 00486 template<typename T> 00487 bool operator<(const IndexData<T>& i1, const IndexData<T>& i2) 00488 { 00489 return i1.index<i2.index; 00490 }; 00491 00497 template<typename T> 00498 std::istream& operator>>(std::istream& is, IndexData<T>& data) 00499 { 00500 is>>data.index; 00501 /* MatrixMarket indices are one based. Decrement for C++ */ 00502 --data.index; 00503 return is>>data.number; 00504 } 00505 00512 template<typename D, int brows, int bcols> 00513 struct MatrixValuesSetter 00514 { 00520 template<typename M> 00521 void operator()(const std::vector<std::set<IndexData<D> > >& rows, 00522 M& matrix) 00523 { 00524 for(typename M::RowIterator iter=matrix.begin(); 00525 iter!= matrix.end(); ++iter) 00526 { 00527 for(typename M::size_type brow=iter.index()*brows, 00528 browend=iter.index()*brows+brows; 00529 brow<browend; ++brow) 00530 { 00531 typedef typename std::set<IndexData<D> >::const_iterator Siter; 00532 for(Siter siter=rows[brow].begin(), send=rows[brow].end(); 00533 siter != send; ++siter) 00534 (*iter)[siter->index/bcols][brow%brows][siter->index%bcols]=siter->number; 00535 } 00536 } 00537 } 00538 }; 00539 00540 template<int brows, int bcols> 00541 struct MatrixValuesSetter<PatternDummy,brows,bcols> 00542 { 00543 template<typename M> 00544 void operator()(const std::vector<std::set<IndexData<PatternDummy> > >& rows, 00545 M& matrix) 00546 {} 00547 }; 00548 00549 template<typename T, typename A, int brows, int bcols, typename D> 00550 void readSparseEntries(Dune::BCRSMatrix<Dune::FieldMatrix<T,brows,bcols>,A>& matrix, 00551 std::istream& file, std::size_t entries, 00552 const MMHeader& mmHeader, D) 00553 { 00554 typedef Dune::BCRSMatrix<Dune::FieldMatrix<T,brows,bcols>,A> Matrix; 00555 std::vector<std::set<IndexData<D> > > rows(matrix.N()*brows); 00556 00557 for(entries; entries>0;--entries){ 00558 std::size_t row; 00559 IndexData<D> data; 00560 skipComments(file); 00561 file>>row; 00562 --row; // Index was 1 based. 00563 assert(row/bcols<matrix.N()); 00564 file>>data; 00565 assert(data.index/bcols<matrix.M()); 00566 rows[row].insert(data); 00567 } 00568 00569 // TODO extend to capture the nongeneral cases. 00570 if(mmHeader.structure!= general) 00571 DUNE_THROW(Dune::NotImplemented, "Only general is supported right now!"); 00572 00573 // Setup the matrix sparsity pattern 00574 int nnz=0; 00575 for(typename Matrix::CreateIterator iter=matrix.createbegin(); 00576 iter!= matrix.createend(); ++iter) 00577 { 00578 for(std::size_t brow=iter.index()*brows, browend=iter.index()*brows+brows; 00579 brow<browend; ++brow) 00580 { 00581 typedef typename std::set<IndexData<D> >::const_iterator Siter; 00582 for(Siter siter=rows[brow].begin(), send=rows[brow].end(); 00583 siter != send; ++siter, ++nnz) 00584 iter.insert(siter->index/bcols); 00585 } 00586 } 00587 00588 //Set the matrix values 00589 matrix=0; 00590 00591 MatrixValuesSetter<D,brows,bcols> Setter; 00592 00593 Setter(rows, matrix); 00594 } 00595 } // end anonymous namespace 00596 00597 class MatrixMarketFormatError : public Dune::Exception 00598 {}; 00599 00606 template<typename T, typename A, int brows, int bcols> 00607 void readMatrixMarket(Dune::BCRSMatrix<Dune::FieldMatrix<T,brows,bcols>,A>& matrix, 00608 std::istream& istr) 00609 { 00610 00611 typedef Dune::BCRSMatrix<Dune::FieldMatrix<double,brows,bcols> > Matrix; 00612 00613 MMHeader header; 00614 if(!readMatrixMarketBanner(istr, header)) 00615 std::cerr << "First line was not a correct Matrix Market banner. Using default:\n" 00616 << "%%MatrixMarket matrix coordinate real general"<<std::endl; 00617 00618 skipComments(istr); 00619 00620 std::size_t rows, cols, entries; 00621 00622 if(lineFeed(istr)) 00623 throw MatrixMarketFormatError(); 00624 00625 istr >> rows; 00626 00627 if(lineFeed(istr)) 00628 throw MatrixMarketFormatError(); 00629 istr >> cols; 00630 00631 if(lineFeed(istr)) 00632 throw MatrixMarketFormatError(); 00633 00634 istr >>entries; 00635 00636 std::size_t nnz, blockrows, blockcols; 00637 00638 Dune::tie(blockrows, blockcols, nnz) = calculateNNZ<brows, bcols>(rows, cols, entries, header); 00639 00640 istr.ignore(std::numeric_limits<std::streamsize>::max(),'\n'); 00641 00642 00643 matrix.setSize(blockrows, blockcols); 00644 matrix.setBuildMode(Dune::BCRSMatrix<Dune::FieldMatrix<T,brows,bcols>,A>::row_wise); 00645 00646 if(header.type==array_type) 00647 DUNE_THROW(Dune::NotImplemented, "currently not supported!"); 00648 00649 NumericWrapper<double> d; 00650 00651 readSparseEntries(matrix, istr, entries, header,d); 00652 } 00653 template<typename M> 00654 void printMatrixMarket(M& matrix, 00655 std::ostream ostr) 00656 { 00657 mm_header_printer<M>::print(ostr, matrix); 00658 mm_block_structure_header<M>::print(ostr,matrix); 00659 } 00660 00661 } 00662 #endif