Dune Core Modules (2.3.1)

matrixmarket.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_MATRIXMARKET_HH
4#define DUNE_MATRIXMARKET_HH
5
6#include <ostream>
7#include <istream>
8#include <fstream>
9#include <sstream>
10#include <limits>
11#include <ios>
12#include "matrixutils.hh"
13#include "bcrsmatrix.hh"
14#include "owneroverlapcopy.hh"
16#include <dune/common/tuples.hh>
17#include <dune/common/unused.hh>
18
19namespace Dune
20{
21
47 namespace
48 {
58 template<class T>
59 struct mm_numeric_type {
60 enum {
64 is_numeric=false
65 };
66 };
67
68 template<>
69 struct mm_numeric_type<int>
70 {
71 enum {
75 is_numeric=true
76 };
77
78 static std::string str()
79 {
80 return "integer";
81 }
82 };
83
84 template<>
85 struct mm_numeric_type<double>
86 {
87 enum {
91 is_numeric=true
92 };
93
94 static std::string str()
95 {
96 return "real";
97 }
98 };
99
100 template<>
101 struct mm_numeric_type<float>
102 {
103 enum {
107 is_numeric=true
108 };
109
110 static std::string str()
111 {
112 return "real";
113 }
114 };
115
116 template<>
117 struct mm_numeric_type<std::complex<double> >
118 {
119 enum {
123 is_numeric=true
124 };
125
126 static std::string str()
127 {
128 return "complex";
129 }
130 };
131
132 template<>
133 struct mm_numeric_type<std::complex<float> >
134 {
135 enum {
139 is_numeric=true
140 };
141
142 static std::string str()
143 {
144 return "complex";
145 }
146 };
147
156 template<class M>
157 struct mm_header_printer;
158
159 template<typename T, typename A, int i, int j>
160 struct mm_header_printer<BCRSMatrix<FieldMatrix<T,i,j>,A> >
161 {
162 static void print(std::ostream& os)
163 {
164 os<<"%%MatrixMarket matrix coordinate ";
165 os<<mm_numeric_type<T>::str()<<" general"<<std::endl;
166 }
167 };
168
169 template<typename B, typename A>
170 struct mm_header_printer<BlockVector<B,A> >
171 {
172 static void print(std::ostream& os)
173 {
174 os<<"%%MatrixMarket matrix array ";
175 os<<mm_numeric_type<typename B::field_type>::str()<<" general"<<std::endl;
176 }
177 };
178
179 template<typename T, int j>
180 struct mm_header_printer<FieldVector<T,j> >
181 {
182 static void print(std::ostream& os)
183 {
184 os<<"%%MatrixMarket matrix array ";
185 os<<mm_numeric_type<T>::str()<<" general"<<std::endl;
186 }
187 };
188
189 template<typename T, int i, int j>
190 struct mm_header_printer<FieldMatrix<T,i,j> >
191 {
192 static void print(std::ostream& os)
193 {
194 os<<"%%MatrixMarket matrix array ";
195 os<<mm_numeric_type<T>::str()<<" general"<<std::endl;
196 }
197 };
198
207 template<class M>
208 struct mm_block_structure_header;
209
210 template<typename T, typename A, int i>
211 struct mm_block_structure_header<BlockVector<FieldVector<T,i>,A> >
212 {
213 typedef BlockVector<FieldVector<T,i>,A> M;
214
215 static void print(std::ostream& os, const M& m)
216 {
217 os<<"% ISTL_STRUCT blocked ";
218 os<<i<<" "<<1<<std::endl;
219 }
220 };
221
222 template<typename T, typename A, int i, int j>
223 struct mm_block_structure_header<BCRSMatrix<FieldMatrix<T,i,j>,A> >
224 {
225 typedef BCRSMatrix<FieldMatrix<T,i,j>,A> M;
226
227 static void print(std::ostream& os, const M& m)
228 {
229 os<<"% ISTL_STRUCT blocked ";
230 os<<i<<" "<<j<<std::endl;
231 }
232 };
233
234
235 template<typename T, int i, int j>
236 struct mm_block_structure_header<FieldMatrix<T,i,j> >
237 {
238 typedef FieldMatrix<T,i,j> M;
239
240 static void print(std::ostream& os, const M& m)
241 {}
242 };
243
244 template<typename T, int i>
245 struct mm_block_structure_header<FieldVector<T,i> >
246 {
247 typedef FieldVector<T,i> M;
248
249 static void print(std::ostream& os, const M& m)
250 {}
251 };
252
253 enum LineType { MM_HEADER, MM_ISTLSTRUCT, DATA };
254 enum { MM_MAX_LINE_LENGTH=1025 };
255
256 enum MM_TYPE { coordinate_type, array_type, unknown_type };
257
258 enum MM_CTYPE { integer_type, double_type, complex_type, pattern, unknown_ctype };
259
260 enum MM_STRUCTURE { general, symmetric, skew_symmetric, hermitian, unknown_structure };
261
262 struct MMHeader
263 {
264 MMHeader()
265 : type(coordinate_type), ctype(double_type), structure(general)
266 {}
267 MM_TYPE type;
268 MM_CTYPE ctype;
269 MM_STRUCTURE structure;
270 };
271
272 bool lineFeed(std::istream& file)
273 {
274 char c;
275 if(!file.eof())
276 c=file.peek();
277 else
278 return false;
279 // ignore whitespace
280 while(c==' ')
281 {
282 file.get();
283 if(file.eof())
284 return false;
285 c=file.peek();
286 }
287
288 if(c=='\n') {
289 /* eat the line feed */
290 file.get();
291 return true;
292 }
293 return false;
294 }
295
296 void skipComments(std::istream& file)
297 {
298 lineFeed(file);
299 char c=file.peek();
300 // ignore comment lines
301 while(c=='%')
302 {
303 /* disgard the rest of the line */
304 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
305 c=file.peek();
306 }
307 }
308
309
310 bool readMatrixMarketBanner(std::istream& file, MMHeader& mmHeader)
311 {
312 std::string buffer;
313 char c;
314 file >> buffer;
315 c=buffer[0];
316 mmHeader=MMHeader();
317 if(c!='%')
318 return false;
319 std::cout<<buffer<<std::endl;
320 /* read the banner */
321 if(buffer!="%%MatrixMarket") {
322 /* disgard the rest of the line */
323 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
324 return false;
325 }
326
327 if(lineFeed(file))
328 /* premature end of line */
329 return false;
330
331 /* read the matrix_type */
332 file >> buffer;
333
334 if(buffer != "matrix")
335 {
336 /* disgard the rest of the line */
337 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
338 return false;
339 }
340
341 if(lineFeed(file))
342 /* premature end of line */
343 return false;
344
345 /* The type of the matrix */
346 file >> buffer;
347
348 if(buffer.empty())
349 return false;
350
351 std::transform(buffer.begin(), buffer.end(), buffer.begin(),
352 tolower);
353
354 switch(buffer[0])
355 {
356 case 'a' :
357 /* sanity check */
358 if(buffer != "array")
359 {
360 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
361 return false;
362 }
363 mmHeader.type=array_type;
364 break;
365 case 'c' :
366 /* sanity check */
367 if(buffer != "coordinate")
368 {
369 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
370 return false;
371 }
372 mmHeader.type=coordinate_type;
373 break;
374 default :
375 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
376 return false;
377 }
378
379 if(lineFeed(file))
380 /* premature end of line */
381 return false;
382
383 /* The numeric type used. */
384 file >> buffer;
385
386 if(buffer.empty())
387 return false;
388
389 std::transform(buffer.begin(), buffer.end(), buffer.begin(),
390 tolower);
391 switch(buffer[0])
392 {
393 case 'i' :
394 /* sanity check */
395 if(buffer != "integer")
396 {
397 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
398 return false;
399 }
400 mmHeader.ctype=integer_type;
401 break;
402 case 'r' :
403 /* sanity check */
404 if(buffer != "real")
405 {
406 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
407 return false;
408 }
409 mmHeader.ctype=double_type;
410 break;
411 case 'c' :
412 /* sanity check */
413 if(buffer != "complex")
414 {
415 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
416 return false;
417 }
418 mmHeader.ctype=complex_type;
419 break;
420 case 'p' :
421 /* sanity check */
422 if(buffer != "pattern")
423 {
424 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
425 return false;
426 }
427 mmHeader.ctype=pattern;
428 break;
429 default :
430 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
431 return false;
432 }
433
434 if(lineFeed(file))
435 return false;
436
437 file >> buffer;
438
439 std::transform(buffer.begin(), buffer.end(), buffer.begin(),
440 tolower);
441 switch(buffer[0])
442 {
443 case 'g' :
444 /* sanity check */
445 if(buffer != "general")
446 {
447 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
448 return false;
449 }
450 mmHeader.structure=general;
451 break;
452 case 'h' :
453 /* sanity check */
454 if(buffer != "hermitian")
455 {
456 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
457 return false;
458 }
459 mmHeader.structure=hermitian;
460 break;
461 case 's' :
462 if(buffer.size()==1) {
463 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
464 return false;
465 }
466
467 switch(buffer[1])
468 {
469 case 'y' :
470 /* sanity check */
471 if(buffer != "symmetric")
472 {
473 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
474 return false;
475 }
476 mmHeader.structure=symmetric;
477 break;
478 case 'k' :
479 /* sanity check */
480 if(buffer != "skew-symmetric")
481 {
482 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
483 return false;
484 }
485 mmHeader.structure=skew_symmetric;
486 break;
487 default :
488 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
489 return false;
490 }
491 default :
492 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
493 return false;
494 }
495 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
496 c=file.peek();
497 return true;
498
499 }
500
501 void readNextLine(std::istream& file, std::ostringstream&, LineType& type)
502 {
504 char c;
505 std::size_t index=0;
506
507 //empty lines will be disgarded and we will simply read the next line
508 while(index==0&&!file.eof())
509 {
510 // strip spaces
511 while(!file.eof() && (c=file.get())==' ') ;
512
513 //read the rest of the line until comment
514 while(!file.eof() && (c=file.get())=='\n') {
515 switch(c)
516 {
517 case '%' :
518 /* disgard the rest of the line */
519 file.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
520 }
521 }
522 }
523
524 // buffer[index]='\0';
525 }
526
527 template<std::size_t brows, std::size_t bcols>
529 calculateNNZ(std::size_t rows, std::size_t cols, std::size_t entries, const MMHeader& header)
530 {
531 std::size_t blockrows=rows/brows;
532 std::size_t blockcols=cols/bcols;
533 std::size_t blocksize=brows*bcols;
534 std::size_t blockentries=0;
535
536 switch(header.structure)
537 {
538 case general :
539 blockentries = entries/blocksize; break;
540 case skew_symmetric :
541 blockentries = 2*entries/blocksize; break;
542 case symmetric :
543 blockentries = (2*entries-rows)/blocksize; break;
544 case hermitian :
545 blockentries = (2*entries-rows)/blocksize; break;
546 default :
547 throw Dune::NotImplemented();
548 }
549 return Dune::make_tuple(blockrows, blockcols, blockentries);
550 }
551
552 /*
553 * @brief Storage class for the column index and the numeric value.
554 *
555 * \tparam T Either a NumericWrapper of the numeric type or PatternDummy
556 * for MatrixMarket pattern case.
557 */
558 template<typename T>
559 struct IndexData : public T
560 {
561 std::size_t index;
562 };
563
564
575 template<typename T>
576 struct NumericWrapper
577 {
578 T number;
579 operator T&()
580 {
581 return number;
582 }
583 };
584
588 struct PatternDummy
589 {};
590
591 template<>
592 struct NumericWrapper<PatternDummy>
593 {};
594
595 template<typename T>
596 std::istream& operator>>(std::istream& is, NumericWrapper<T>& num)
597 {
598 return is>>num.number;
599 }
600
601 std::istream& operator>>(std::istream& is, NumericWrapper<PatternDummy>& num)
602 {
604 return is;
605 }
606
612 template<typename T>
613 bool operator<(const IndexData<T>& i1, const IndexData<T>& i2)
614 {
615 return i1.index<i2.index;
616 }
617
623 template<typename T>
624 std::istream& operator>>(std::istream& is, IndexData<T>& data)
625 {
626 is>>data.index;
627 /* MatrixMarket indices are one based. Decrement for C++ */
628 --data.index;
629 return is>>data.number;
630 }
631
638 template<typename D, int brows, int bcols>
639 struct MatrixValuesSetter
640 {
646 template<typename M>
647 void operator()(const std::vector<std::set<IndexData<D> > >& rows,
648 M& matrix)
649 {
650 for(typename M::RowIterator iter=matrix.begin();
651 iter!= matrix.end(); ++iter)
652 {
653 for(typename M::size_type brow=iter.index()*brows,
654 browend=iter.index()*brows+brows;
655 brow<browend; ++brow)
656 {
657 typedef typename std::set<IndexData<D> >::const_iterator Siter;
658 for(Siter siter=rows[brow].begin(), send=rows[brow].end();
659 siter != send; ++siter)
660 (*iter)[siter->index/bcols][brow%brows][siter->index%bcols]=siter->number;
661 }
662 }
663 }
664 };
665
666 template<int brows, int bcols>
667 struct MatrixValuesSetter<PatternDummy,brows,bcols>
668 {
669 template<typename M>
670 void operator()(const std::vector<std::set<IndexData<PatternDummy> > >& rows,
671 M& matrix)
672 {}
673 };
674
675 template<typename T, typename A, int brows, int bcols, typename D>
676 void readSparseEntries(Dune::BCRSMatrix<Dune::FieldMatrix<T,brows,bcols>,A>& matrix,
677 std::istream& file, std::size_t entries,
678 const MMHeader& mmHeader, const D&)
679 {
681 // First path
682 // store entries together with column index in a speparate
683 // data structure
684 std::vector<std::set<IndexData<D> > > rows(matrix.N()*brows);
685
686 for(; entries>0; --entries) {
687 std::size_t row;
688 IndexData<D> data;
689 skipComments(file);
690 file>>row;
691 --row; // Index was 1 based.
692 assert(row/bcols<matrix.N());
693 file>>data;
694 assert(data.index/bcols<matrix.M());
695 rows[row].insert(data);
696 }
697
698 // TODO extend to capture the nongeneral cases.
699 if(mmHeader.structure!= general)
700 DUNE_THROW(Dune::NotImplemented, "Only general is supported right now!");
701
702 // Setup the matrix sparsity pattern
703 int nnz=0;
704 for(typename Matrix::CreateIterator iter=matrix.createbegin();
705 iter!= matrix.createend(); ++iter)
706 {
707 for(std::size_t brow=iter.index()*brows, browend=iter.index()*brows+brows;
708 brow<browend; ++brow)
709 {
710 typedef typename std::set<IndexData<D> >::const_iterator Siter;
711 for(Siter siter=rows[brow].begin(), send=rows[brow].end();
712 siter != send; ++siter, ++nnz)
713 iter.insert(siter->index/bcols);
714 }
715 }
716
717 //Set the matrix values
718 matrix=0;
719
720 MatrixValuesSetter<D,brows,bcols> Setter;
721
722 Setter(rows, matrix);
723 }
724 } // end anonymous namespace
725
726 class MatrixMarketFormatError : public Dune::Exception
727 {};
728
729
730 void mm_read_header(std::size_t& rows, std::size_t& cols, MMHeader& header, std::istream& istr,
731 bool isVector)
732 {
733 if(!readMatrixMarketBanner(istr, header)) {
734 std::cerr << "First line was not a correct Matrix Market banner. Using default:\n"
735 << "%%MatrixMarket matrix coordinate real general"<<std::endl;
736 // Go to the beginning of the file
737 istr.clear() ;
738 istr.seekg(0, std::ios::beg);
739 if(isVector)
740 header.type=array_type;
741 }
742
743 skipComments(istr);
744
745 if(lineFeed(istr))
746 throw MatrixMarketFormatError();
747
748 istr >> rows;
749
750 if(lineFeed(istr))
751 throw MatrixMarketFormatError();
752 istr >> cols;
753 }
754
755 template<typename T, typename A, int entries>
756 void mm_read_vector_entries(Dune::BlockVector<Dune::FieldVector<T,entries>,A>& vector,
757 std::size_t size,
758 std::istream& istr)
759 {
760 for(int i=0; size>0; ++i, --size) {
761 T val;
762 istr>>val;
763 vector[i/entries][i%entries]=val;
764 }
765 }
766
767
774 template<typename T, typename A, int entries>
776 std::istream& istr)
777 {
778 MMHeader header;
779 std::size_t rows, cols;
780 mm_read_header(rows,cols,header,istr, true);
781 if(cols!=1)
782 DUNE_THROW(MatrixMarketFormatError, "cols!=1, therefore this is no vector!");
783
784 if(header.type!=array_type)
785 DUNE_THROW(MatrixMarketFormatError, "Vectors have to be stored in array format!");
786
787 std::size_t size=rows/entries;
788 if(size*entries!=rows)
789 DUNE_THROW(MatrixMarketFormatError, "Block size of vector is not correct1");
790
791 vector.resize(size);
792
793 istr.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
794
795 mm_read_vector_entries(vector, rows, istr);
796 }
797
798
805 template<typename T, typename A, int brows, int bcols>
807 std::istream& istr)
808 {
809 MMHeader header;
810 if(!readMatrixMarketBanner(istr, header)) {
811 std::cerr << "First line was not a correct Matrix Market banner. Using default:\n"
812 << "%%MatrixMarket matrix coordinate real general"<<std::endl;
813 // Go to the beginning of the file
814 istr.clear() ;
815 istr.seekg(0, std::ios::beg);
816 }
817 skipComments(istr);
818
819 std::size_t rows, cols, entries;
820
821 if(lineFeed(istr))
822 throw MatrixMarketFormatError();
823
824 istr >> rows;
825
826 if(lineFeed(istr))
827 throw MatrixMarketFormatError();
828 istr >> cols;
829
830 if(lineFeed(istr))
831 throw MatrixMarketFormatError();
832
833 istr >>entries;
834
835 std::size_t nnz, blockrows, blockcols;
836
837 Dune::tie(blockrows, blockcols, nnz) = calculateNNZ<brows, bcols>(rows, cols, entries, header);
838
839 istr.ignore(std::numeric_limits<std::streamsize>::max(),'\n');
840
841
842 matrix.setSize(blockrows, blockcols);
843 matrix.setBuildMode(Dune::BCRSMatrix<Dune::FieldMatrix<T,brows,bcols>,A>::row_wise);
844
845 if(header.type==array_type)
846 DUNE_THROW(Dune::NotImplemented, "Array format currently not supported for matrices!");
847
848 readSparseEntries(matrix, istr, entries, header, NumericWrapper<T>());
849 }
850
851 template<typename M>
852 struct mm_multipliers
853 {};
854
855 template<typename B, int i, int j, typename A>
856 struct mm_multipliers<BCRSMatrix<FieldMatrix<B,i,j>,A> >
857 {
858 enum {
859 rows = i,
860 cols = j
861 };
862 };
863
864 template<typename B, int i, int j>
865 void mm_print_entry(const FieldMatrix<B,i,j>& entry,
866 typename FieldMatrix<B,i,j>::size_type rowidx,
867 typename FieldMatrix<B,i,j>::size_type colidx,
868 std::ostream& ostr)
869 {
870 typedef typename FieldMatrix<B,i,j>::const_iterator riterator;
871 typedef typename FieldMatrix<B,i,j>::ConstColIterator citerator;
872 for(riterator row=entry.begin(); row != entry.end(); ++row, ++rowidx) {
873 int coli=colidx;
874 for(citerator col = row->begin(); col != row->end(); ++col, ++coli)
875 ostr<< rowidx<<" "<<coli<<" "<<*col<<std::endl;
876 }
877 }
878
879 // Write a vector entry
880 template<typename V>
881 void mm_print_vector_entry(const V& entry, std::ostream& ostr,
882 const integral_constant<int,1>&)
883 {
884 ostr<<entry<<std::endl;
885 }
886
887 // Write a vector
888 template<typename V>
889 void mm_print_vector_entry(const V& vector, std::ostream& ostr,
890 const integral_constant<int,0>&)
891 {
892 // Is the entry a supported numeric type?
893 const int isnumeric = mm_numeric_type<typename V::block_type>::is_numeric;
894 typedef typename V::const_iterator VIter;
895
896 for(VIter i=vector.begin(); i != vector.end(); ++i)
897
898 mm_print_vector_entry(*i, ostr,
899 integral_constant<int,isnumeric>());
900 }
901
902 template<typename T, typename A, int i>
903 std::size_t countEntries(const BlockVector<FieldVector<T,i>,A>& vector)
904 {
905 return vector.size()*i;
906 }
907
908 // Version for writing vectors.
909 template<typename V>
910 void writeMatrixMarket(const V& vector, std::ostream& ostr,
911 const integral_constant<int,0>&)
912 {
913 ostr<<countEntries(vector)<<" "<<1<<std::endl;
914 const int isnumeric = mm_numeric_type<typename V::block_type>::is_numeric;
915 mm_print_vector_entry(vector,ostr, integral_constant<int,isnumeric>());
916 }
917
918 // Versions for writing matrices
919 template<typename M>
920 void writeMatrixMarket(const M& matrix,
921 std::ostream& ostr,
922 const integral_constant<int,1>&)
923 {
924 ostr<<matrix.M()*mm_multipliers<M>::rows<<" "
925 <<matrix.N()*mm_multipliers<M>::cols<<" "
926 <<countNonZeros(matrix)<<std::endl;
927
928 typedef typename M::const_iterator riterator;
929 typedef typename M::ConstColIterator citerator;
930 for(riterator row=matrix.begin(); row != matrix.end(); ++row)
931 for(citerator col = row->begin(); col != row->end(); ++col)
932 // Matrix Market indexing start with 1!
933 mm_print_entry(*col, row.index()*mm_multipliers<M>::rows+1,
934 col.index()*mm_multipliers<M>::cols+1, ostr);
935 }
936
937
941 template<typename M>
942 void writeMatrixMarket(const M& matrix,
943 std::ostream& ostr)
944 {
945 // Write header information
946 mm_header_printer<M>::print(ostr);
947 mm_block_structure_header<M>::print(ostr,matrix);
948 // Choose the correct function for matrix and vector
949 writeMatrixMarket(matrix,ostr,integral_constant<int,IsMatrix<M>::value>());
950 }
951
952
963 template<typename M>
964 void storeMatrixMarket(const M& matrix,
965 std::string filename)
966 {
967 std::ofstream file(filename.c_str());
968 file.setf(std::ios::scientific,std::ios::floatfield);
969 writeMatrixMarket(matrix, file);
970 file.close();
971 }
972
973#if HAVE_MPI
987 template<typename M, typename G, typename L>
988 void storeMatrixMarket(const M& matrix,
989 std::string filename,
991 bool storeIndices=true)
992 {
993 // Get our rank
994 int rank = comm.communicator().rank();
995 // Write the local matrix
996 std::ostringstream rfilename;
997 rfilename<<filename <<"_"<<rank<<".mm";
998 std::cout<< rfilename.str()<<std::endl;
999 std::ofstream file(rfilename.str().c_str());
1000 file.setf(std::ios::scientific,std::ios::floatfield);
1001 writeMatrixMarket(matrix, file);
1002 file.close();
1003
1004 if(!storeIndices)
1005 return;
1006
1007 // Write the global to local index mapping
1008 rfilename.str("");
1009 rfilename<<filename<<"_"<<rank<<".idx";
1010 file.open(rfilename.str().c_str());
1011 file.setf(std::ios::scientific,std::ios::floatfield);
1013 typedef typename IndexSet::const_iterator Iterator;
1014 for(Iterator iter = comm.indexSet().begin();
1015 iter != comm.indexSet().end(); ++iter) {
1016 file << iter->global()<<" "<<(std::size_t)iter->local()<<" "
1017 <<(int)iter->local().attribute()<<" "<<(int)iter->local().isPublic()<<std::endl;
1018 }
1019 // Store neighbour information for efficient remote indices setup.
1020 file<<"neighbours:";
1021 const std::set<int>& neighbours=comm.remoteIndices().getNeighbours();
1022 typedef std::set<int>::const_iterator SIter;
1023 for(SIter neighbour=neighbours.begin(); neighbour != neighbours.end(); ++neighbour) {
1024 file<<" "<< *neighbour;
1025 }
1026 file.close();
1027 }
1028
1043 template<typename M, typename G, typename L>
1044 void loadMatrixMarket(M& matrix,
1045 const std::string& filename,
1047 bool readIndices=true)
1048 {
1050 typedef typename LocalIndex::Attribute Attribute;
1051 // Get our rank
1052 int rank = comm.communicator().rank();
1053 // load local matrix
1054 std::ostringstream rfilename;
1055 rfilename<<filename <<"_"<<rank<<".mm";
1056 std::ifstream file;
1057 file.open(rfilename.str().c_str(), std::ios::in);
1058 if(!file)
1059 DUNE_THROW(IOError, "Could not open file" << rfilename.str().c_str());
1060 //if(!file.is_open())
1061 //
1062 readMatrixMarket(matrix,file);
1063 file.close();
1064
1065 if(!readIndices)
1066 return;
1067
1068 // read indices
1070 IndexSet& pis=comm.pis;
1071 rfilename.str("");
1072 rfilename<<filename<<"_"<<rank<<".idx";
1073 file.open(rfilename.str().c_str());
1074 if(pis.size()!=0)
1075 DUNE_THROW(InvalidIndexSetState, "Index set is not empty!");
1076
1077 pis.beginResize();
1078 while(!file.eof() && file.peek()!='n') {
1079 G g;
1080 file >>g;
1081 std::size_t l;
1082 file >>l;
1083 int c;
1084 file >>c;
1085 bool b;
1086 file >> b;
1087 pis.add(g,LocalIndex(l,Attribute(c),b));
1088 lineFeed(file);
1089 }
1090 pis.endResize();
1091 if(!file.eof()) {
1092 // read neighbours
1093 std::string s;
1094 file>>s;
1095 if(s!="neighbours:")
1096 DUNE_THROW(MatrixMarketFormatError, "was exspecting the string: \"neighbours:\"");
1097 std::set<int> nb;
1098 while(!file.eof()) {
1099 int i;
1100 file >> i;
1101 nb.insert(i);
1102 }
1103 file.close();
1104 comm.ri.setNeighbours(nb);
1105 }
1106 comm.ri.template rebuild<false>();
1107 }
1108
1109 #endif
1110
1121 template<typename M>
1122 void loadMatrixMarket(M& matrix,
1123 const std::string& filename)
1124 {
1125 std::ifstream file;
1126 file.open(filename.c_str(), std::ios::in);
1127 if(!file)
1128 DUNE_THROW(IOError, "Could not open file" << filename);
1129 readMatrixMarket(matrix,file);
1130 file.close();
1131 }
1132
1134}
1135#endif
Implementation of the BCRSMatrix class.
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:414
A vector of blocks with memory management.
Definition: bvector.hh:254
int rank() const
Return rank, is between 0 and size()-1.
Definition: mpicollectivecommunication.hh:162
ConstIterator const_iterator
typedef for stl compliant access
Definition: densematrix.hh:290
row_type::ConstIterator ConstColIterator
rename the iterators for easier access
Definition: densematrix.hh:294
Base class for Dune-Exceptions.
Definition: exceptions.hh:92
A dense n x m matrix.
Definition: fmatrix.hh:67
vector space out of a tensor product of fields.
Definition: fvector.hh:92
Default exception class for I/O errors.
Definition: exceptions.hh:257
Index Set Interface base class.
Definition: indexidset.hh:78
IndexType size(GeometryType type) const
Return total number of entities of given geometry type in entity set .
Definition: indexidset.hh:204
Exception indicating that the index set is not in the expected state.
Definition: indexset.hh:205
An index present on the local process.
Definition: localindex.hh:34
Default exception for dummy implementations.
Definition: exceptions.hh:289
A class setting up standard communication for a two-valued attribute set with owner/overlap/copy sema...
Definition: owneroverlapcopy.hh:173
const ParallelIndexSet & indexSet() const
Get the underlying parallel index set.
Definition: owneroverlapcopy.hh:466
const RemoteIndices & remoteIndices() const
Get the underlying remote indices.
Definition: owneroverlapcopy.hh:475
A Tuple of objects.
Definition: tuples.hh:294
Implements a matrix constructed from a given type representing a field and compile-time given number ...
iterator begin()
Get an iterator over the indices positioned at the first index.
iterator end()
Get an iterator over the indices positioned after the last index.
LI LocalIndex
The type of the local index, e.g. ParallelLocalIndex.
Definition: indexset.hh:239
std::istream & operator>>(std::istream &is, Pair< T1, T2 > &pair)
Read a pair or tuple.
Definition: tuples.hh:949
#define DUNE_THROW(E, m)
Definition: exceptions.hh:244
void readMatrixMarket(Dune::BlockVector< Dune::FieldVector< T, entries >, A > &vector, std::istream &istr)
Reads a BlockVector from a matrix market file.
Definition: matrixmarket.hh:775
void loadMatrixMarket(M &matrix, const std::string &filename, OwnerOverlapCopyCommunication< G, L > &comm, bool readIndices=true)
Load a parallel matrix/vector stored in matrix market format.
Definition: matrixmarket.hh:1044
void storeMatrixMarket(const M &matrix, std::string filename)
Stores a parallel matrix/vector in matrix market format in a file.
Definition: matrixmarket.hh:964
int countNonZeros(const M &matrix)
Get the number of nonzero fields in the matrix.
Definition: matrixutils.hh:155
Some handy generic functions for ISTL matrices.
Dune namespace.
Definition: alignment.hh:14
STL namespace.
Classes providing communication interfaces for overlapping Schwarz methods.
Test whether a type is an ISTL Matrix.
Definition: matrixutils.hh:453
Generate a type for a given integral constant.
Definition: typetraits.hh:457
Fallback implementation of the std::tuple class.
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentional unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)