00001 #ifndef DUNE_SUPERLUMATRIX_HH
00002 #define DUNE_SUPERLUMATRIX_HH
00003
00004 #ifdef HAVE_SUPERLU
00005 #ifdef SUPERLU_POST_2005_VERSION
00006 #include "slu_ddefs.h"
00007 #else
00008 #include "dsp_defs.h"
00009 #endif
00010 #include"bcrsmatrix.hh"
00011 #include"bvector.hh"
00012 #include<dune/common/fmatrix.hh>
00013 #include<dune/common/fvector.hh>
00014 namespace Dune
00015 {
00016
00017
00018 template<class T>
00019 struct GetSuperLUType
00020 {
00021 };
00022
00023 template<>
00024 struct GetSuperLUType<double>
00025 {
00026 enum{ type = SLU_D};
00027
00028 };
00029
00030 template<>
00031 struct GetSuperLUType<float>
00032 {
00033 enum{ type = SLU_S};
00034
00035 };
00036
00037 template<>
00038 struct GetSuperLUType<std::complex<double> >
00039 {
00040 enum{ type = SLU_Z};
00041
00042 };
00043
00044 template<>
00045 struct GetSuperLUType<std::complex<float> >
00046 {
00047 enum{ type = SLU_C};
00048
00049 };
00050
00051
00056 template<class M>
00057 struct SuperLUMatrix
00058 {
00059
00060 };
00061
00062 template<class M>
00063 struct SuperMatrixInitializer
00064 {};
00065
00066 template<class M, class X, class TM, class T1>
00067 class SeqOverlappingSchwarz;
00068
00072 template<class B, class TA, int n, int m>
00073 class SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
00074 {
00075 template<class M, class X, class TM, class T1>
00076 friend class SeqOverlappingSchwarz;
00077 friend class SuperMatrixInitializer<BCRSMatrix<FieldMatrix<B,n,m>,TA> >;
00078
00079 public:
00081 typedef BCRSMatrix<FieldMatrix<B,n,m>,TA> Matrix;
00082
00083 typedef typename Matrix::size_type size_type;
00084
00089 SuperLUMatrix(const Matrix& mat);
00090
00091 SuperLUMatrix();
00092
00094 ~SuperLUMatrix();
00095
00097 operator SuperMatrix&()
00098 {
00099 return A;
00100 }
00101
00103 operator const SuperMatrix&()const
00104 {
00105 return A;
00106 }
00107 bool operator==(const Matrix& mat) const;
00108
00113 size_type N() const
00114 {
00115 return N_;
00116 }
00117
00122 size_type M() const
00123 {
00124 return M_;
00125 }
00126
00127 SuperLUMatrix& operator=(const Matrix& mat);
00128
00129 SuperLUMatrix& operator=(const SuperLUMatrix& mat);
00130
00131 private:
00133 void setMatrix(const Matrix& mat);
00135 void free();
00136
00137 int N_, M_, Nnz_;
00138 B* values;
00139 int* rowindex;
00140 int* colstart;
00141 SuperMatrix A;
00142 };
00143
00144 template<class T, class A, int n, int m>
00145 void writeCompColMatrixToMatlab(const SuperLUMatrix<BCRSMatrix<FieldMatrix<T,n,m>,A> >& mat,
00146 std::ostream& os)
00147 {
00148 const SuperMatrix a=static_cast<const SuperMatrix&>(mat);
00149 NCformat *astore = (NCformat *) a.Store;
00150 double* dp = (double*)astore->nzval;
00151
00152
00153 std::ios_base::fmtflags oldflags = os.flags();
00154
00155
00156 int oldprec = os.precision();
00157
00158
00159
00160 os <<"[";
00161 for(int row=0; row<a.nrow; ++row){
00162 for(int col= 0; col < a.ncol; ++col){
00163
00164 int i;
00165 for(i=astore->colptr[col]; i < astore->colptr[col+1]; ++i)
00166 if(astore->rowind[i]==row){
00167 os<<dp[i]<<" ";
00168 break;
00169 }
00170 if(i==astore->colptr[col+1])
00171
00172 os<<0<<" ";
00173 }
00174 if(row==a.nrow-1)
00175 os<<"]";
00176 os<<std::endl;
00177 }
00178
00179 os.flags(oldflags);
00180 os.precision(oldprec);
00181 }
00182
00183
00184 template<class T, class A, int n, int m>
00185 class SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >
00186 {
00187 template<class I, class S>
00188 friend class OverlappingSchwarzInitializer;
00189 public:
00190 typedef Dune::BCRSMatrix<FieldMatrix<T,n,m>,A> Matrix;
00191 typedef Dune::SuperLUMatrix<Matrix> SuperLUMatrix;
00192 typedef typename Matrix::const_iterator Iter;
00193 typedef typename Matrix::row_type::const_iterator CIter;
00194 typedef typename Matrix::size_type size_type;
00195
00196 SuperMatrixInitializer(SuperLUMatrix& lum);
00197
00198 SuperMatrixInitializer();
00199
00200 ~SuperMatrixInitializer();
00201
00202 void addRowNnz(const Iter& row)const;
00203
00204 void allocate();
00205
00206 void countEntries(const Iter& row, const CIter& col) const;
00207
00208 void countEntries(size_type colidx)const;
00209
00210 void calcColstart()const;
00211
00212 void copyValue(const Iter& row, const CIter& col) const;
00213
00214 void copyValue(const CIter& col, size_type rowindex, size_type colidx)const;
00215
00216 void createMatrix() const;
00217
00218 private:
00219
00220 void allocateMatrixStorage()const;
00221
00222 void allocateMarker();
00223
00224 SuperLUMatrix* mat;
00225 int cols;
00226 typename Matrix::size_type *marker;
00227 };
00228
00229 template<class T, class A, int n, int m>
00230 SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::SuperMatrixInitializer(SuperLUMatrix& mat_)
00231 : mat(&mat_), cols(mat_.M()), marker(0)
00232 {
00233 mat->Nnz_=0;
00234 }
00235
00236 template<class T, class A, int n, int m>
00237 SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::SuperMatrixInitializer()
00238 : mat(0), cols(0), marker(0)
00239 {}
00240
00241 template<class T, class A, int n, int m>
00242 SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::~SuperMatrixInitializer()
00243 {
00244 if(marker)
00245 delete[] marker;
00246 }
00247
00248 template<class T, class A, int n, int m>
00249 void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::addRowNnz(const Iter& row)const
00250 {
00251 mat->Nnz_+=row->getsize();
00252 }
00253
00254 template<class T, class A, int n, int m>
00255 void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::allocate()
00256 {
00257 allocateMatrixStorage();
00258 allocateMarker();
00259 }
00260
00261 template<class T, class A, int n, int m>
00262 void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::allocateMatrixStorage()const
00263 {
00264 mat->Nnz_*=n*m;
00265
00266 mat->values=new T[mat->Nnz_];
00267 mat->rowindex=new int[mat->Nnz_];
00268 mat->colstart=new int[cols+1];
00269 }
00270
00271 template<class T, class A, int n, int m>
00272 void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::allocateMarker()
00273 {
00274 marker = new typename Matrix::size_type[cols];
00275
00276 for(int i=0; i < cols; ++i)
00277 marker[i]=0;
00278 }
00279
00280 template<class T, class A, int n, int m>
00281 void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::countEntries(const Iter& row, const CIter& col)const
00282 {
00283 countEntries(col.index());
00284
00285 }
00286
00287 template<class T, class A, int n, int m>
00288 void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::countEntries(size_type colindex)const
00289 {
00290 for(int i=0; i < m; ++i)
00291 marker[colindex*m+i]+=n;
00292 }
00293
00294 template<class T, class A, int n, int m>
00295 void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::calcColstart()const
00296 {
00297 mat->colstart[0]=0;
00298 for(int i=0; i < cols; ++i){
00299 mat->colstart[i+1]=mat->colstart[i]+marker[i];
00300 marker[i]=mat->colstart[i];
00301 }
00302 }
00303
00304 template<class T, class A, int n, int m>
00305 void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::copyValue(const Iter& row, const CIter& col)const
00306 {
00307 copyValue(col, row.index(), col.index());
00308 }
00309
00310 template<class T, class A, int n, int m>
00311 void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::copyValue(const CIter& col, size_type rowindex, size_type colindex)const
00312 {
00313 for(int i=0; i<n;i++){
00314 for(int j=0; j<m; j++){
00315 mat->rowindex[marker[colindex*m+j]]=rowindex*n+i;
00316 mat->values[marker[colindex*m+j]]=(*col)[i][j];
00317 ++marker[colindex*m+j];
00318 }
00319 }
00320 }
00321
00322 template<class T, class A, int n, int m>
00323 void SuperMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A> >::createMatrix() const
00324 {
00325 dCreate_CompCol_Matrix(&mat->A, mat->N_, mat->M_, mat->Nnz_,
00326 mat->values, mat->rowindex, mat->colstart, SLU_NC, static_cast<Dtype_t>(GetSuperLUType<T>::type), SLU_GE);
00327 }
00328
00329 template<class F, class Matrix>
00330 void copyToSuperMatrix(F& initializer, const Matrix& mat)
00331 {
00332 typedef typename Matrix::const_iterator Iter;
00333
00334 for(Iter row=mat.begin(); row!= mat.end(); ++row)
00335 initializer.addRowNnz(row);
00336
00337 initializer.allocate();
00338
00339 for(Iter row=mat.begin(); row!= mat.end(); ++row){
00340 typedef typename Matrix::row_type::const_iterator CIter;
00341 for(CIter col=row->begin(); col != row->end(); ++col)
00342 initializer.countEntries(row, col);
00343 }
00344
00345 initializer.calcColstart();
00346
00347 for(Iter row=mat.begin(); row!= mat.end(); ++row){
00348 typedef typename Matrix::row_type::const_iterator CIter;
00349 for(CIter col=row->begin(); col != row->end(); ++col)
00350 initializer.copyValue(row, col);
00351 }
00352 initializer.createMatrix();
00353 }
00354
00355
00356 template<class B, class TA, int n, int m>
00357 bool SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator==(const BCRSMatrix<FieldMatrix<B,n,m>,TA>& mat) const
00358 {
00359 const NCformat* S=static_cast<const NCformat *>(A.Store);
00360 for(size_type col=0; col < M(); ++col){
00361 for(int j=S->colptr[col]; j < S->colptr[col+1]; ++j){
00362 int row=S->rowind[j];
00363 if((mat[row/n][col/m])[row%n][col%m]!=reinterpret_cast<B*>(S->nzval)[j]){
00364 std::cerr<<" bcrs["<<row/n<<"]["<<col/m<<"]["<<row%n<<"]["<<row%m
00365 <<"]="<<(mat[row/n][col/m])[row%n][col%m]<<" super["<<row<<"]["<<col<<"]="<<reinterpret_cast<B*>(S->nzval)[j];
00366
00367 return false;
00368 }
00369 }
00370 }
00371
00372 return true;
00373 }
00374
00375 template<class B, class TA, int n, int m>
00376 bool operator==(SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >& sla, BCRSMatrix<FieldMatrix<B,n,m>,TA>& a)
00377 {
00378 return a==sla;
00379 }
00380
00381 template<class B, class TA, int n, int m>
00382 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::SuperLUMatrix()
00383 : N_(0), M_(0), Nnz_(0), values(0), rowindex(0), colstart(0)
00384 {}
00385
00386 template<class B, class TA, int n, int m>
00387 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
00388 ::SuperLUMatrix(const Matrix& mat)
00389 : N_(n*mat.N()), M_(m*mat.M()), Nnz_(n*m*mat.Nnz())
00390 {
00391
00392 }
00393
00394 template<class B, class TA, int n, int m>
00395 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >&
00396 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(const Matrix& mat)
00397 {
00398 if(N_+M_+Nnz_!=0)
00399 free();
00400 setMatrix(mat);
00401 return *this;
00402 }
00403
00404 template<class B, class TA, int n, int m>
00405 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >&
00406 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::operator=(const SuperLUMatrix& mat)
00407 {
00408 if(N_+M_+Nnz_!=0)
00409 free();
00410 N_=mat.N_;
00411 M_=mat.M_;
00412 Nnz_= mat.Nnz_;
00413 if(M_>0){
00414 colstart=new int[M_+1];
00415 for(int i=0; i<=M_; ++i)
00416 colstart[i]=mat.colstart[i];
00417 }
00418
00419 if(Nnz_>0){
00420 values = new B[Nnz_];
00421 rowindex = new int[Nnz_];
00422
00423 for(int i=0; i<Nnz_; ++i)
00424 values[i]=mat.values[i];
00425
00426 for(int i=0; i<Nnz_; ++i)
00427 rowindex[i]=mat.rowindex[i];
00428 }
00429 if(M_+Nnz_>0)
00430 dCreate_CompCol_Matrix(&A, N_, M_, Nnz_,
00431 values, rowindex, colstart, SLU_NC, static_cast<Dtype_t>(GetSuperLUType<B>::type), SLU_GE);
00432 return *this;
00433 }
00434
00435 template<class B, class TA, int n, int m>
00436 void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >
00437 ::setMatrix(const Matrix& mat)
00438 {
00439 N_=n*mat.N();
00440 M_=m*mat.M();
00441 SuperMatrixInitializer<Matrix> initializer(*this);
00442
00443 copyToSuperMatrix(initializer,mat);
00444
00445 #ifdef DUNE_ISTL_WITH_CHECKING
00446 char name[] = {'A',0};
00447 if(N_<0)
00448 dPrint_CompCol_Matrix(name,&A);
00449 assert(*this==mat);
00450 #endif
00451 }
00452
00453 template<class B, class TA, int n, int m>
00454 SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::~SuperLUMatrix()
00455 {
00456 if(N_+M_+Nnz_!=0)
00457 free();
00458 }
00459
00460 template<class B, class TA, int n, int m>
00461 void SuperLUMatrix<BCRSMatrix<FieldMatrix<B,n,m>,TA> >::free()
00462 {
00463 delete[] values;
00464 delete[] rowindex;
00465 delete[] colstart;
00466 SUPERLU_FREE(A.Store);
00467 }
00468
00469 }
00470 #endif
00471 #endif