supermatrix.hh

Go to the documentation of this file.
00001 // -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 2 -*-
00002 // vi: set et ts=8 sw=2 sts=2:
00003 #ifndef DUNE_SUPERLUMATRIX_HH
00004 #define DUNE_SUPERLUMATRIX_HH
00005 
00006 #ifdef HAVE_SUPERLU
00007 #ifdef SUPERLU_POST_2005_VERSION
00008 #include "slu_ddefs.h"
00009 #else
00010 #include "dsp_defs.h"
00011 #endif
00012 #include"bcrsmatrix.hh"
00013 #include"bvector.hh"
00014 #include<dune/common/fmatrix.hh>
00015 #include<dune/common/fvector.hh>
00016 #include<limits>
00017 
00018 namespace Dune
00019 {
00020 
00026   template<class M>
00027   class MatrixRowSet
00028   {
00029   public:
00030     // @brief The type of the matrix.
00031     typedef M Matrix;
00032     // @brief The matrix row iterator type.
00033     typedef typename Matrix::ConstRowIterator const_iterator;
00034     
00039     MatrixRowSet(const Matrix& m)
00040       : m_(m)
00041     {}
00042     
00043     // @brief Get the row iterator at the first row.
00044     const_iterator begin() const
00045     {
00046       return m_.begin();
00047     }
00048     //@brief Get the row iterator at the end of all rows.
00049     const_iterator end() const
00050     {
00051       return m_.end();
00052     }
00053   private:
00054     const Matrix& m_;
00055   };
00056   
00064   template<class M, class S>
00065   class MatrixRowSubset
00066   {
00067   public:
00068     /* @brief the type of the matrix class. */
00069     typedef M Matrix;
00070     /* @brief the type of the set of valid row indices. */
00071     typedef S RowIndexSet;
00072     
00078     MatrixRowSubset(const Matrix& m, const RowIndexSet& s)
00079       : m_(m), s_(s)
00080     {}
00081 
00082     const Matrix& matrix() const
00083     {
00084       return m_;
00085     }
00086     
00087     const RowIndexSet& rowIndexSet() const
00088     {
00089       return s_;
00090     }
00091     
00092     // @brief The matrix row iterator type.
00093     class const_iterator
00094       : public ForwardIteratorFacade<const_iterator, const typename Matrix::row_type>
00095     {
00096     public:
00097       const_iterator(typename Matrix::const_iterator firstRow,
00098                      typename RowIndexSet::const_iterator pos)
00099         : firstRow_(firstRow), pos_(pos)
00100       {}
00101       
00102       
00103       const typename Matrix::row_type& dereference() const
00104       {
00105         return *(firstRow_+ *pos_);        
00106       }
00107       bool equals(const const_iterator& o) const
00108       {
00109         return pos_==o.pos_;
00110       }
00111       void increment()
00112       {
00113         ++pos_;
00114       }
00115       typename RowIndexSet::value_type index() const
00116       {
00117         return *pos_;
00118       }
00119       
00120     private:
00121       // @brief Iterator pointing to the first row of the matrix.
00122       typename Matrix::const_iterator firstRow_;
00123       // @brief Iterator pointing to the current row index.
00124       typename RowIndexSet::const_iterator pos_; 
00125     };
00126 
00127     // @brief Get the row iterator at the first row.
00128     const_iterator begin() const
00129     {
00130       return const_iterator(m_.begin(), s_.begin());
00131     }
00132     //@brief Get the row iterator at the end of all rows.
00133     const_iterator end() const
00134     {
00135       return const_iterator(m_.begin(), s_.end());
00136     }
00137     
00138   private:
00139     // @brief The matrix for which we manage the row subset.
00140     const Matrix& m_;
00141     // @brief The set of row indices we manage.
00142     const RowIndexSet& s_;
00143     
00144   };
00145   
00146   template<class T>
00147   struct GetSuperLUType
00148   {
00149   };
00150 
00151   template<>
00152   struct GetSuperLUType<double>
00153   {
00154     enum{ type = SLU_D};
00155   
00156   };
00157 
00158   template<>
00159   struct GetSuperLUType<float>
00160   {
00161     enum{ type = SLU_S};
00162   
00163   };
00164 
00165   template<>
00166   struct GetSuperLUType<std::complex<double> >
00167   {
00168     enum{ type = SLU_Z};
00169   
00170   };
00171 
00172   template<>
00173   struct GetSuperLUType<std::complex<float> >
00174   {
00175     enum{ type = SLU_C};
00176   
00177   };
00178 
00179 
00184   template<class M>
00185   struct SuperLUMatrix
00186   {
00187   
00188   };
00189 
00190   template<class M>
00191   struct SuperMatrixInitializer
00192   {};
00193 
00194   template<class M, class X, class TM, bool b, class T1>
00195   class SeqOverlappingSchwarz;
00196 
00200   template<class B, class TA, int n, int m>
00201   class SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
00202   {
00203     template<class M, class X, class TM, bool b, class T1>
00204     friend class SeqOverlappingSchwarz;
00205     friend class SuperMatrixInitializer<BCRSMatrix<FieldMatrix<B,n,m>,TA> >;
00206     
00207   public:
00209     typedef BCRSMatrix<FieldMatrix<B,n,m>,TA> Matrix;
00210   
00211     typedef typename Matrix::size_type size_type;
00212     
00217     SuperLUMatrix(const Matrix& mat);
00218 
00219     SuperLUMatrix();
00220     
00222     ~SuperLUMatrix();
00223     
00225     operator SuperMatrix&()
00226     {
00227       return A;
00228     }
00229     
00231     operator const SuperMatrix&()const
00232     {
00233       return A;
00234     }
00235     bool operator==(const Matrix& mat) const;
00236     
00241     size_type N() const
00242     {
00243       return N_;
00244     }
00245     
00250     size_type M() const
00251     {
00252       return M_;
00253     }
00254     
00255     SuperLUMatrix& operator=(const Matrix& mat);
00256     
00257     SuperLUMatrix& operator=(const SuperLUMatrix& mat);
00258 
00265     template<class S>
00266     void setMatrix(const Matrix& mat, const S& mrs);
00268     void free();
00269   private: 
00271     void setMatrix(const Matrix& mat);
00272 
00273     int N_, M_, Nnz_;
00274     B* values;
00275     int* rowindex;
00276     int* colstart;
00277     SuperMatrix A;
00278   };
00279 
00280   template<class T, class A, int n, int m>
00281   void writeCompColMatrixToMatlab(const SuperLUMatrix<BCRSMatrix<FieldMatrix<T,n,m>,A> >& mat,
00282                            std::ostream& os)
00283   {
00284     const SuperMatrix a=static_cast<const SuperMatrix&>(mat);
00285     NCformat *astore = (NCformat *) a.Store;
00286     double* dp = (double*)astore->nzval;
00287 
00288     // remember old flags
00289     std::ios_base::fmtflags oldflags = os.flags();
00290     // set the output format
00291     //os.setf(std::ios_base::scientific, std::ios_base::floatfield);
00292     int oldprec = os.precision();
00293     //os.precision(10);
00294     //dPrint_CompCol_Matrix("A", const_cast<SuperMatrix*>(&a));
00295     
00296     os <<"[";
00297     for(int row=0; row<a.nrow; ++row){
00298       for(int col= 0; col < a.ncol; ++col){
00299         // linear search for col
00300         int i;
00301         for(i=astore->colptr[col]; i < astore->colptr[col+1]; ++i)
00302           if(astore->rowind[i]==row){
00303             os<<dp[i]<<" ";
00304             break;
00305           }
00306         if(i==astore->colptr[col+1])
00307           // entry not present
00308           os<<0<<" ";
00309       }
00310       if(row==a.nrow-1)
00311         os<<"]";
00312       os<<std::endl;
00313     }
00314     // reset the output format
00315     os.flags(oldflags);
00316     os.precision(oldprec);
00317   }
00318   
00319   
00320   template<class T, class A, int n, int m>
00321   class SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >
00322   {
00323     template<class I, class S>
00324     friend class OverlappingSchwarzInitializer;
00325   public:
00326     typedef Dune::BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
00327     typedef Dune::SuperLUMatrix<Matrix> SuperLUMatrix;
00328     typedef typename Matrix::row_type::const_iterator CIter;
00329     typedef typename Matrix::size_type size_type;
00330     
00331     SuperMatrixInitializer(SuperLUMatrix& lum);
00332 
00333     SuperMatrixInitializer();
00334     
00335     ~SuperMatrixInitializer();
00336     
00337     template<typename Iter>
00338     void addRowNnz(const Iter& row)const;
00339     
00340     void allocate();
00341         
00342     template<typename Iter>
00343     void countEntries(const Iter& row, const CIter& col) const;
00344     
00345     void countEntries(size_type colidx)const;
00346 
00347     void calcColstart()const;
00348 
00349     template<typename Iter>
00350     void copyValue(const Iter& row, const CIter& col) const;
00351 
00352     void copyValue(const CIter& col, size_type rowindex, size_type colidx)const;
00353 
00354     void createMatrix() const;
00355     
00356   private:
00357     
00358     void allocateMatrixStorage()const;
00359 
00360     void allocateMarker();
00361 
00362     SuperLUMatrix* mat;
00363     int cols;
00364     typename Matrix::size_type *marker;
00365   };
00366   
00367   template<class T, class A, int n, int m>
00368   SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::SuperMatrixInitializer(SuperLUMatrix& mat_)
00369     : mat(&mat_), cols(mat_.N()), marker(0)
00370   {
00371     mat->Nnz_=0;
00372   }
00373 
00374   template<class T, class A, int n, int m>
00375   SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::SuperMatrixInitializer()
00376     : mat(0), cols(0), marker(0)
00377   {}
00378   
00379   template<class T, class A, int n, int m>
00380   SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::~SuperMatrixInitializer()
00381   {
00382     if(marker)
00383       delete[] marker;
00384   }
00385   
00386   template<class T, class A, int n, int m>
00387   template<typename Iter>
00388   void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::addRowNnz(const Iter& row)const
00389   {
00390     mat->Nnz_+=row->getsize();
00391   }
00392   
00393   template<class T, class A, int n, int m>
00394   void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::allocate()
00395   {
00396     allocateMatrixStorage();
00397     allocateMarker();
00398   }
00399   
00400   template<class T, class A, int n, int m>
00401   void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::allocateMatrixStorage()const
00402   {
00403     mat->Nnz_*=n*m;
00404     // initialize data
00405     mat->values=new T[mat->Nnz_];
00406     mat->rowindex=new int[mat->Nnz_];
00407     mat->colstart=new int[cols+1];
00408   }
00409 
00410   template<class T, class A, int n, int m>
00411   void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::allocateMarker()
00412   {          
00413     marker = new typename Matrix::size_type[cols];
00414 
00415     for(int i=0; i < cols; ++i)
00416       marker[i]=0;
00417   }
00418   
00419   template<class T, class A, int n, int m>
00420   template<typename Iter>
00421   void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::countEntries(const Iter& row, const CIter& col)const
00422   {
00423     countEntries(col.index());
00424     
00425   }
00426   
00427   template<class T, class A, int n, int m>
00428   void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::countEntries(size_type colindex)const
00429   {
00430     for(int i=0; i < m; ++i)
00431       marker[colindex*m+i]+=n;
00432   }
00433 
00434   template<class T, class A, int n, int m>
00435   void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::calcColstart()const
00436   {
00437     mat->colstart[0]=0;
00438     for(int i=0; i < cols; ++i){
00439         mat->colstart[i+1]=mat->colstart[i]+marker[i];
00440         marker[i]=mat->colstart[i];
00441     }
00442   }
00443   
00444   template<class T, class A, int n, int m>
00445   template<typename Iter>
00446   void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::copyValue(const Iter& row, const CIter& col)const
00447   {
00448     copyValue(col, row.index(), col.index());
00449   }
00450 
00451   template<class T, class A, int n, int m>
00452   void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::copyValue(const CIter& col, size_type rowindex, size_type colindex)const
00453   {
00454     for(int i=0; i<n;i++){
00455       for(int j=0; j<m; j++){
00456         mat->rowindex[marker[colindex*m+j]]=rowindex*n+i;
00457         mat->values[marker[colindex*m+j]]=(*col)[i][j];
00458         ++marker[colindex*m+j];
00459       }
00460     }
00461   }
00462   
00463   template<class T, class A, int n, int m>
00464   void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::createMatrix() const
00465   {
00466     dCreate_CompCol_Matrix(&mat->A, mat->N_, mat->M_, mat->colstart[cols], 
00467                            mat->values, mat->rowindex, mat->colstart, SLU_NC, static_cast<Dtype_t>(GetSuperLUType<T>::type), SLU_GE);
00468   }
00469 
00470   template<class F, class MRS>
00471   void copyToSuperMatrix(F& initializer, const MRS& mrs)
00472   {
00473     typedef typename MRS::const_iterator Iter;
00474     typedef typename  std::iterator_traits<Iter>::value_type::const_iterator CIter;
00475     for(Iter row=mrs.begin(); row!= mrs.end(); ++row)
00476       initializer.addRowNnz(row);
00477     
00478     initializer.allocate();
00479     
00480     for(Iter row=mrs.begin(); row!= mrs.end(); ++row){
00481       
00482       for(CIter col=row->begin(); col != row->end(); ++col)
00483         initializer.countEntries(row, col);
00484     }
00485     
00486     initializer.calcColstart();
00487     
00488     for(Iter row=mrs.begin(); row!= mrs.end(); ++row){
00489       for(CIter col=row->begin(); col != row->end(); ++col)
00490         initializer.copyValue(row, col);
00491     }
00492     initializer.createMatrix();
00493   }
00494 
00495   template<class F, class M,class S>
00496   void copyToSuperMatrix(F& initializer, const MatrixRowSubset<M,S>& mrs)
00497   {
00498     typedef MatrixRowSubset<M,S> MRS;
00499     typedef typename MRS::RowIndexSet SIS;
00500     typedef typename SIS::const_iterator SIter;
00501     typedef typename MRS::const_iterator Iter;
00502     typedef typename std::iterator_traits<Iter>::value_type row_type;
00503     typedef typename row_type::const_iterator CIter;
00504 
00505     // Calculate upper Bound for nonzeros
00506     for(Iter row=mrs.begin(); row!= mrs.end(); ++row)
00507       initializer.addRowNnz(row);
00508 
00509     initializer.allocate();
00510 
00511     typedef typename MRS::Matrix::size_type size_type;
00512 
00513     // A vector containing the corresponding indices in
00514     // the to create submatrix.
00515     // If an entry is the maximum of size_type then this index will not appear in 
00516     // the submatrix.
00517     std::vector<size_type> subMatrixIndex(mrs.matrix().N(), 
00518                                          std::numeric_limits<size_type>::max());
00519     size_type s=0;
00520     for(SIter index = mrs.rowIndexSet().begin(); index!=mrs.rowIndexSet().end(); ++index)
00521       subMatrixIndex[*index]=s++;
00522     
00523     for(Iter row=mrs.begin(); row!= mrs.end(); ++row)
00524       for(CIter col=row->begin(); col != row->end(); ++col){
00525         if(subMatrixIndex[col.index()]!=std::numeric_limits<size_type>::max())
00526           // This column is in our subset (use submatrix column index)
00527           initializer.countEntries(subMatrixIndex[col.index()]);
00528     }
00529     
00530     initializer.calcColstart();
00531     
00532     for(Iter row=mrs.begin(); row!= mrs.end(); ++row)
00533       for(CIter col=row->begin(); col != row->end(); ++col){
00534         if(subMatrixIndex[col.index()]!=std::numeric_limits<size_type>::max())
00535           // This value is in our submatrix -> copy (use submatrix indices
00536           initializer.copyValue(col, subMatrixIndex[row.index()], subMatrixIndex[col.index()]);
00537     }
00538     initializer.createMatrix();
00539   }
00540 
00541 #ifndef DOXYGEN
00542   
00543   template<class B, class TA, int n, int m>
00544   bool  SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator==(const BCRSMatrix<FieldMatrix<B,n,m>,TA>& mat) const
00545   {
00546     const NCformat* S=static_cast<const NCformat *>(A.Store);
00547     for(size_type col=0; col < M(); ++col){
00548       for(int j=S->colptr[col]; j < S->colptr[col+1]; ++j){
00549         int row=S->rowind[j];
00550         if((mat[row/n][col/m])[row%n][col%m]!=reinterpret_cast<B*>(S->nzval)[j]){
00551           std::cerr<<" bcrs["<<row/n<<"]["<<col/m<<"]["<<row%n<<"]["<<row%m
00552                    <<"]="<<(mat[row/n][col/m])[row%n][col%m]<<" super["<<row<<"]["<<col<<"]="<<reinterpret_cast<B*>(S->nzval)[j];
00553           
00554           return false;
00555         }
00556       }
00557     }
00558     
00559     return true;
00560   }
00561   
00562 #endif // DOYXGEN
00563 
00564   template<class B, class TA, int n, int m>
00565   bool operator==(SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >& sla, BCRSMatrix<FieldMatrix<B,n,m>,TA>& a)
00566   {
00567     return a==sla;
00568   }
00569 
00570 #ifndef DOXYGEN
00571   
00572   template<class B, class TA, int n, int m>
00573   SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::SuperLUMatrix()
00574     : N_(0), M_(0), Nnz_(0), values(0), rowindex(0), colstart(0)
00575   {}
00576 
00577   template<class B, class TA, int n, int m>
00578   SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
00579   ::SuperLUMatrix(const Matrix& mat)
00580     : N_(n*mat.N()), M_(m*mat.M()), Nnz_(n*m*mat.Nnz())
00581   {
00582     
00583   }
00584 
00585   template<class B, class TA, int n, int m>
00586   SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >&
00587   SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(const Matrix& mat)
00588   {
00589     if(N_+M_+Nnz_!=0)
00590       free();
00591     setMatrix(mat);
00592     return *this;
00593   }
00594 
00595   template<class B, class TA, int n, int m>
00596   SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >&
00597   SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(const SuperLUMatrix& mat)
00598   {
00599     if(N_+M_+Nnz_!=0)
00600       free();
00601     N_=mat.N_;
00602     M_=mat.M_;
00603     Nnz_= mat.Nnz_;
00604     if(M_>0){
00605       colstart=new int[M_+1];
00606       for(int i=0; i<=M_; ++i)
00607         colstart[i]=mat.colstart[i];
00608     }
00609     
00610     if(Nnz_>0){
00611       values = new B[Nnz_];
00612       rowindex = new int[Nnz_];
00613       
00614       for(int i=0; i<Nnz_; ++i)
00615         values[i]=mat.values[i];
00616       
00617       for(int i=0; i<Nnz_; ++i)
00618         rowindex[i]=mat.rowindex[i];
00619     }
00620     if(M_+Nnz_>0)
00621       dCreate_CompCol_Matrix(&A, N_, M_, Nnz_, 
00622                            values, rowindex, colstart, SLU_NC, static_cast<Dtype_t>(GetSuperLUType<B>::type), SLU_GE);
00623     return *this;
00624   }
00625 
00626   template<class B, class TA, int n, int m>
00627   void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
00628   ::setMatrix(const Matrix& mat)
00629   {
00630     N_=n*mat.N();
00631     M_=m*mat.M();
00632     SuperMatrixInitializer<Matrix> initializer(*this);
00633     
00634     copyToSuperMatrix(initializer, MatrixRowSet<Matrix>(mat));
00635     
00636 #ifdef DUNE_ISTL_WITH_CHECKING
00637     char name[] = {'A',0};
00638     if(N_<0)
00639       dPrint_CompCol_Matrix(name,&A);
00640     assert(*this==mat);
00641 #endif
00642   }
00643 
00644    template<class B, class TA, int n, int m>
00645    template<class S>
00646   void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
00647    ::setMatrix(const Matrix& mat, const S& mrs)
00648   {
00649     if(N_+M_+Nnz_!=0)
00650       free();
00651     N_=mrs.size()*n;
00652     M_=mrs.size()*m;
00653     SuperMatrixInitializer<Matrix> initializer(*this);
00654     
00655     copyToSuperMatrix(initializer, MatrixRowSubset<Matrix,S>(mat,mrs));
00656     
00657 #ifdef DUNE_ISTL_WITH_CHECKING
00658     char name[] = {'A',0};
00659     if(N_<0)
00660       dPrint_CompCol_Matrix(name,&A);
00661 #endif
00662   }
00663 
00664   template<class B, class TA, int n, int m>
00665   SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::~SuperLUMatrix()
00666   {
00667     if(N_+M_+Nnz_!=0)
00668       free();
00669   }
00670 
00671   template<class B, class TA, int n, int m>
00672   void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::free()
00673   {
00674     delete[] values;
00675     delete[] rowindex;
00676     delete[] colstart;
00677     SUPERLU_FREE(A.Store);
00678     N_=M_=Nnz_=0;
00679   }
00680 
00681 #endif // DOXYGEN
00682 
00683 }
00684 #endif  
00685 #endif
Generated on Sat Apr 24 11:13:47 2010 for dune-istl by  doxygen 1.6.3