io.hh

Go to the documentation of this file.
00001 #ifndef DUNE_ISTLIO_HH
00002 #define DUNE_ISTLIO_HH
00003 
00004 #include<cmath>
00005 #include<complex>
00006 #include <limits>
00007 #include<ios>
00008 #include<iomanip>
00009 #include<fstream>
00010 #include<string>
00011 
00012 #include "istlexception.hh"
00013 #include <dune/common/fvector.hh>
00014 #include <dune/common/fmatrix.hh>
00015 #include "bcrsmatrix.hh"
00016 
00017 
00018 namespace Dune {
00019    
00030   //==== pretty printing of vectors
00031 
00032   // recursively print all the blocks
00033   template<class V>
00034   void recursive_printvector (std::ostream& s, const V& v, std::string rowtext, int& counter, 
00035                                                           int columns, int width, int precision)
00036   {
00037         for (typename V::ConstIterator i=v.begin(); i!=v.end(); ++i)
00038           recursive_printvector(s,*i,rowtext,counter,columns,width,precision);
00039   }
00040 
00041   // specialization for FieldVector
00042   template<class K, int n>
00043   void recursive_printvector (std::ostream& s, const FieldVector<K,n>& v, std::string rowtext, int& counter, 
00044                                                           int columns, int width, int precision)
00045   {
00046         // we now can print n numbers
00047         for (int i=0; i<n; i++)
00048           {
00049                 if (counter%columns==0)
00050                   { 
00051                         s << rowtext; // start a new row
00052                         s << " ";     // space in front of each entry
00053                         s.width(4);   // set width for counter
00054                         s << counter; // number of first entry in a line
00055                   }
00056                 s << " ";         // space in front of each entry
00057                 s.width(width);   // set width for each entry anew
00058                 s << v[i];        // yeah, the number !
00059                 counter++;        // increment the counter
00060                 if (counter%columns==0)
00061                   s << std::endl; // start a new line
00062           }
00063   }
00064 
00065 
00066   template<class V>
00067   void printvector (std::ostream& s, const V& v, std::string title, std::string rowtext, 
00068                                         int columns=1, int width=10, int precision=2)
00069   {
00070         // count the numbers printed to make columns
00071         int counter=0;
00072 
00073         // remember old flags
00074         std::ios_base::fmtflags oldflags = s.flags();
00075         
00076         // set the output format
00077         s.setf(std::ios_base::scientific, std::ios_base::floatfield);
00078         int oldprec = s.precision();
00079         s.precision(precision);
00080 
00081         // print title
00082         s << title << " [blocks=" << v.N() << ",dimension=" << v.dim() << "]" << std::endl;
00083 
00084         // print data from all blocks
00085         recursive_printvector(s,v,rowtext,counter,columns,width,precision);
00086 
00087         // check if new line is required
00088         if (counter%columns!=0)
00089           s << std::endl;
00090 
00091         // reset the output format
00092         s.flags(oldflags);
00093         s.precision(oldprec);
00094   }
00095 
00096 
00097   //==== pretty printing of matrices
00098 
00099 
00101   inline void fill_row (std::ostream& s, int m, int width, int precision)
00102   {
00103         for (int j=0; j<m; j++)
00104           {
00105                 s << " ";         // space in front of each entry
00106                 s.width(width);   // set width for each entry anew
00107                 s << ".";         // yeah, the number !
00108           }
00109   }
00110 
00111   template<typename B, typename A>
00112   class BCRSMatrix;
00113   
00114   template<typename K, int n, int m>
00115   class FieldMatrix;
00116   
00117   template<typename M>
00118   struct MatrixDimension
00119   {
00120   };
00121     
00122     
00123   template<typename B, typename TA>
00124   struct MatrixDimension<BCRSMatrix<B,TA> >
00125   {
00126     typedef BCRSMatrix<B,TA> Matrix;
00127     typedef typename Matrix::block_type block_type;
00128     typedef typename Matrix::size_type size_type;
00129       
00130     static size_type rowdim (const Matrix& A, size_type i)
00131     {
00132       const B* row = A.r[i].getptr();
00133       if(row)
00134         return MatrixDimension<block_type>::rowdim(*row);
00135       else
00136         return 0;
00137     }
00138       
00139     static size_type coldim (const Matrix& A, size_type c)
00140     {
00141       // find an entry in column j
00142       if (A.nnz>0)
00143         {
00144           for (size_type k=0; k<A.nnz; k++) {
00145             if (A.j[k]==c) {
00146               return MatrixDimension<block_type>::coldim(A.a[k]);
00147             }
00148           }
00149         }
00150       else
00151         {
00152           for (size_type i=0; i<A.N(); i++)
00153             {
00154               size_type* j = A.r[i].getindexptr();
00155               B*   a = A.r[i].getptr();
00156               for (size_type k=0; k<A.r[i].getsize(); k++)
00157                 if (j[k]==c) {
00158                   return MatrixDimension<block_type>::coldim(a[k]);
00159                 }
00160             }
00161         }
00162 
00163       // not found
00164       return 0;
00165     }
00166 
00167     static size_type rowdim (const Matrix& A){
00168       size_type nn=0;
00169       for (size_type i=0; i<A.N(); i++)
00170         nn += rowdim(A,i);
00171       return nn;
00172     }
00173       
00174     static size_type coldim (const Matrix& A){
00175       typedef typename Matrix::ConstRowIterator ConstRowIterator;
00176       typedef typename Matrix::ConstColIterator ConstColIterator;
00177           
00178       // The following code has a complexity of nnz, and
00179       // typically a very small constant.
00180       //
00181       std::vector<size_type> coldims(A.M(),std::numeric_limits<size_type>::max());
00182           
00183       for (ConstRowIterator row=A.begin(); row!=A.end(); ++row) 
00184         for (ConstColIterator col=row->begin(); col!=row->end(); ++col)
00185           // only compute blocksizes we don't already have
00186           if (coldims[col.index()]==std::numeric_limits<size_type>::max()) 
00187             coldims[col.index()] = MatrixDimension<block_type>::coldim(*col);
00188 
00189       size_type sum = 0;
00190       for (typename std::vector<size_type>::iterator it=coldims.begin(); it!=coldims.end(); ++it)
00191         // skip rows for which no coldim could be determined
00192         if ((*it)>=0)
00193           sum += *it;
00194 
00195       return sum;
00196     }
00197   };
00198     
00199 
00200     template<typename B, int n, int m, typename TA>
00201   struct MatrixDimension<BCRSMatrix<FieldMatrix<B,n,m> ,TA> >
00202   {
00203     typedef BCRSMatrix<FieldMatrix<B,n,m> ,TA> Matrix;
00204     typedef typename Matrix::size_type size_type;
00205       
00206     static size_type rowdim (const Matrix& A, size_type i)
00207     {
00208         return n;
00209     }
00210       
00211     static size_type coldim (const Matrix& A, size_type c)
00212     {
00213         return m;
00214     }
00215 
00216     static size_type rowdim (const Matrix& A){
00217         return A.N()*n;
00218     }
00219       
00220     static size_type coldim (const Matrix& A){
00221         return A.M()*m;
00222     }
00223   };
00224     
00225   template<typename K, int n, int m>
00226   struct MatrixDimension<FieldMatrix<K,n,m> >
00227   {
00228     typedef FieldMatrix<K,n,m> Matrix;
00229     typedef typename Matrix::size_type size_type;
00230       
00231     static size_type rowdim(const Matrix& A, size_type r)
00232     {
00233       return 1;
00234     }
00235 
00236     static size_type coldim(const Matrix& A, size_type r)
00237     {
00238       return 1;
00239     }
00240 
00241     static size_type rowdim(const Matrix& A)
00242     {
00243       return n;
00244     }
00245 
00246     static size_type coldim(const Matrix& A)
00247     {
00248       return m;
00249     }
00250   };
00251   
00253   template<class M>
00254   void print_row (std::ostream& s, const M& A, typename M::size_type I,
00255                   typename M::size_type J, typename M::size_type therow,
00256                   int width, int precision)
00257   {
00258         typename M::size_type i0=I;
00259         for (typename M::size_type i=0; i<A.N(); i++)
00260           {
00261                 if (therow>=i0 && therow<i0+MatrixDimension<M>::rowdim(A,i))
00262                   {
00263                         // the row is in this block row !
00264                         typename M::size_type j0=J;
00265                         for (typename M::size_type j=0; j<A.M(); j++)
00266                           {
00267                                 // find this block
00268                                 typename M::ConstColIterator it = A[i].find(j);
00269                                 
00270                                 // print row or filler
00271                                 if (it!=A[i].end())
00272                                   print_row(s,*it,i0,j0,therow,width,precision);
00273                                 else
00274                                   fill_row(s,MatrixDimension<M>::coldim(A,j),width,precision);
00275                                 
00276                                 // advance columns
00277                                 j0 += MatrixDimension<M>::coldim(A,j);
00278                           }
00279                   }
00280                 // advance rows
00281                 i0 += MatrixDimension<M>::rowdim(A,i);
00282           }
00283   }
00284 
00286   template<class K, int n, int m>
00287   void print_row (std::ostream& s, const FieldMatrix<K,n,m>& A, 
00288                   typename FieldMatrix<K,n,m>::size_type I, typename FieldMatrix<K,n,m>::size_type J,
00289                   typename FieldMatrix<K,n,m>::size_type therow, int width, int precision)
00290   {
00291     typedef typename FieldMatrix<K,n,m>::size_type size_type;
00292     
00293         for (size_type i=0; i<n; i++)
00294           if (I+i==therow)
00295                   for (int j=0; j<m; j++)
00296                         {
00297                           s << " ";         // space in front of each entry
00298                           s.width(width);   // set width for each entry anew
00299                           s << A[i][j];     // yeah, the number !
00300                         }
00301   }
00302 
00304   template<class K>
00305   void print_row (std::ostream& s, const FieldMatrix<K,1,1>& A, typename FieldMatrix<K,1,1>::size_type I, 
00306                   typename FieldMatrix<K,1,1>::size_type J, typename FieldMatrix<K,1,1>::size_type therow, 
00307                   int width, int precision)
00308   {
00309         if (I==therow)
00310           {
00311                 s << " ";         // space in front of each entry
00312                 s.width(width);   // set width for each entry anew
00313                 s << static_cast<K>(A);         // yeah, the number !
00314           }
00315   }
00316 
00320   template<class M>
00321   void printmatrix (std::ostream& s, const M& A, std::string title, std::string rowtext, 
00322                                         int width=10, int precision=2)
00323   {
00324     
00325         // remember old flags
00326         std::ios_base::fmtflags oldflags = s.flags();
00327 
00328         // set the output format
00329         s.setf(std::ios_base::scientific, std::ios_base::floatfield);
00330         int oldprec = s.precision();
00331         s.precision(precision);
00332 
00333         // print title
00334         s << title 
00335           << " [n=" << A.N() 
00336           << ",m=" << A.M() 
00337           << ",rowdim=" << MatrixDimension<M>::rowdim(A) 
00338           << ",coldim=" << MatrixDimension<M>::coldim(A) 
00339           << "]" << std::endl;
00340 
00341         // print all rows
00342         for (typename M::size_type i=0; i<MatrixDimension<M>::rowdim(A); i++)
00343           {
00344                 s << rowtext;  // start a new row
00345                 s << " ";      // space in front of each entry
00346                 s.width(4);    // set width for counter
00347                 s << i;        // number of first entry in a line
00348                 print_row(s,A,0,0,i,width,precision); // generic print
00349                 s << std::endl;// start a new line
00350           }
00351 
00352         // reset the output format
00353         s.flags(oldflags);
00354         s.precision(oldprec);
00355   }
00356 
00374   template<class B, int n, int m, class A>
00375   void printSparseMatrix(std::ostream& s, 
00376                          const BCRSMatrix<FieldMatrix<B,n,m>,A>& mat, 
00377                          std::string title, std::string rowtext, 
00378                          int width=3, int precision=2)
00379   {
00380     typedef BCRSMatrix<FieldMatrix<B,n,m>,A> Matrix;
00381     // remember old flags
00382     std::ios_base::fmtflags oldflags = s.flags();
00383     // set the output format
00384     s.setf(std::ios_base::scientific, std::ios_base::floatfield);
00385     int oldprec = s.precision();
00386     s.precision(precision);
00387     // print title
00388     s << title 
00389       << " [n=" << mat.N() 
00390       << ",m=" << mat.M() 
00391       << ",rowdim=" << MatrixDimension<Matrix>::rowdim(mat) 
00392       << ",coldim=" << MatrixDimension<Matrix>::coldim(mat) 
00393       << "]" << std::endl;
00394     
00395     typedef typename Matrix::ConstRowIterator Row;
00396     
00397     for(Row row=mat.begin(); row != mat.end();++row){
00398       int skipcols=0;
00399       bool reachedEnd=false;
00400       
00401       while(!reachedEnd){
00402         for(int innerrow=0; innerrow<n; ++innerrow){
00403           int count=0;
00404           typedef typename Matrix::ConstColIterator Col;
00405           Col col=row->begin();   
00406           for(; col != row->end(); ++col,++count){
00407             if(count<skipcols)
00408               continue;
00409             if(count>=skipcols+width)
00410               break;
00411             if(innerrow==0){
00412               if(count==skipcols){
00413                 s << rowtext;  // start a new row
00414                 s << " ";      // space in front of each entry
00415                 s.width(4);    // set width for counter
00416                 s << row.index()<<": ";        // number of first entry in a line
00417               }
00418               s.width(4);
00419               s<<col.index()<<": |";
00420             }else{
00421               if(count==skipcols){
00422                 for(int i=0; i < rowtext.length(); i++)
00423                   s<<" ";
00424                 s<<"       ";
00425               }
00426               s<<"      |";
00427             }
00428             for(int innercol=0; innercol < m; ++innercol){
00429               s.width(9);
00430               s<<(*col)[innerrow][innercol]<<" ";
00431             }
00432             
00433             s<<"|";
00434           }
00435           if(innerrow==n-1 && col==row->end())
00436             reachedEnd=true;
00437           else
00438             s<<std::endl;
00439         }
00440         skipcols+=width;
00441         s<<std::endl;
00442       }
00443       s<<std::endl;
00444     }
00445     
00446     // reset the output format
00447     s.flags(oldflags);
00448     s.precision(oldprec);
00449   }
00450   
00451   namespace
00452   {
00453     template<typename T>
00454     struct MatlabPODWriter
00455     {
00456       static std::ostream& write(const T& t,  std::ostream& s)
00457       {
00458         s<<t;
00459         return s;
00460       }
00461     };
00462     template<typename T>
00463     struct MatlabPODWriter<std::complex<T> >
00464     {
00465       static std::ostream& write(const std::complex<T>& t,  std::ostream& s)
00466       {
00467         s<<t.real()<<" "<<t.imag();
00468         return s;
00469       }
00470     }; 
00471   }  
00476     template <class FieldType, int rows, int cols>
00477     void writeMatrixToMatlabHelper(const FieldMatrix<FieldType,rows,cols>& matrix,
00478                                    int rowOffset, int colOffset,
00479                                    std::ostream& s)
00480     {
00481         for (int i=0; i<rows; i++)
00482             for (int j=0; j<cols; j++){
00483                 //+1 for Matlab numbering
00484                 s << rowOffset + i + 1 << " " << colOffset + j + 1 << " ";
00485                 MatlabPODWriter<FieldType>::write(matrix[i][j], s)<< std::endl;
00486             }
00487     }
00488     
00489     template <class MatrixType>
00490     void writeMatrixToMatlabHelper(const MatrixType& matrix, 
00491                                    int externalRowOffset, int externalColOffset,
00492                                    std::ostream& s)
00493 {
00494     // Precompute the accumulated sizes of the columns
00495     std::vector<typename MatrixType::size_type> colOffset(matrix.M());
00496     if (colOffset.size() > 0)
00497         colOffset[0] = 0;
00498 
00499     for (typename MatrixType::size_type i=0; i<matrix.M()-1; i++)
00500         colOffset[i+1] = colOffset[i] + MatrixDimension<MatrixType>::coldim(matrix,i);
00501     
00502     typename MatrixType::size_type rowOffset = 0;
00503 
00504     // Loop over all matrix rows
00505     for (typename MatrixType::size_type rowIdx=0; rowIdx<matrix.N(); rowIdx++) {
00506         
00507         const typename MatrixType::row_type& row = matrix[rowIdx];
00508 
00509         typename MatrixType::row_type::ConstIterator cIt   = row.begin();
00510         typename MatrixType::row_type::ConstIterator cEndIt = row.end();
00511         
00512         // Loop over all columns in this row
00513         for (; cIt!=cEndIt; ++cIt) 
00514             writeMatrixToMatlabHelper(*cIt, 
00515                                       externalRowOffset+rowOffset, 
00516                                       externalColOffset + colOffset[cIt.index()], 
00517                                       s);
00518 
00519         rowOffset += MatrixDimension<MatrixType>::rowdim(matrix, rowIdx);
00520     }
00521     
00522 }
00523 
00535 template <class MatrixType>
00536 void writeMatrixToMatlab(const MatrixType& matrix, 
00537                          const std::string& filename)
00538 {
00539     std::ofstream outStream(filename.c_str());
00540 
00541     writeMatrixToMatlabHelper(matrix, 0, 0, outStream);
00542 }
00543 
00546 } // end namespace
00547 
00548 #endif

Generated on Sun Nov 15 22:29:36 2009 for dune-istl by  doxygen 1.5.6