00001 #ifndef DUNE_MATRIXREDIST_HH
00002 #define DUNE_MATRIXREDIST_HH
00003 #include"repartition.hh"
00004 #include<dune/common/exceptions.hh>
00005 #include<dune/istl/indexset.hh>
00006 #include<dune/istl/owneroverlapcopy.hh>
00012 namespace Dune
00013 {
00014 template<typename T>
00015 struct RedistributeInformation
00016 {
00017 bool isSetup() const
00018 {
00019 return false;
00020 }
00021 template<class D>
00022 void redistribute(const D& from, D& to) const
00023 {}
00024
00025 template<class D>
00026 void redistributeBackward(D& from, const D& to) const
00027 {}
00028 };
00029
00030 #if HAVE_MPI
00031 template<typename T, typename T1>
00032 class RedistributeInformation<OwnerOverlapCopyCommunication<T,T1> >
00033 {
00034 public:
00035 typedef OwnerOverlapCopyCommunication<T,T1> Comm;
00036
00037 RedistributeInformation()
00038 : interface(), setup_(false)
00039 {}
00040
00041 RedistributeInterface& getInterface()
00042 {
00043 return interface;
00044 }
00045 template<typename IS>
00046 void checkInterface(const IS& source,
00047 const IS& target, MPI_Comm comm)
00048 {
00049 RemoteIndices<IS> *ri=new RemoteIndices<IS>(source, target, comm);
00050 ri->template rebuild<true>();
00051 Interface inf;
00052 typename OwnerOverlapCopyCommunication<int>::OwnerSet flags;
00053 int rank;
00054 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
00055 inf.free();
00056 inf.build(*ri, flags, flags);
00057
00058
00059 #ifdef DEBUG_REPART
00060 if(inf!=interface){
00061
00062 int rank;
00063 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
00064 if(rank==0)
00065 std::cout<<"Interfaces do not match!"<<std::endl;
00066 std::cout<<rank<<": redist interface new :"<<inf<<std::endl;
00067 std::cout<<rank<<": redist interface :"<<interface<<std::endl;
00068
00069 throw "autsch!";
00070 }
00071
00072 #endif
00073 }
00074 void setSetup()
00075 {
00076 setup_=true;
00077 interface.strip();
00078 }
00079
00080 template<class GatherScatter, class D>
00081 void redistribute(const D& from, D& to) const
00082 {
00083 BufferedCommunicator communicator;
00084 communicator.template build<D>(from,to, interface);
00085 communicator.template forward<GatherScatter>(from, to);
00086 communicator.free();
00087 }
00088 template<class GatherScatter, class D>
00089 void redistributeBackward(D& from, const D& to) const
00090 {
00091
00092 BufferedCommunicator communicator;
00093 communicator.template build<D>(from,to, interface);
00094 communicator.template backward<GatherScatter>(from, to);
00095 communicator.free();
00096 }
00097
00098 template<class D>
00099 void redistribute(const D& from, D& to) const
00100 {
00101 redistribute<CopyGatherScatter<D> >(from,to);
00102 }
00103 template<class D>
00104 void redistributeBackward(D& from, const D& to) const
00105 {
00106 redistributeBackward<CopyGatherScatter<D> >(from,to);
00107 }
00108 bool isSetup() const
00109 {
00110 return setup_;
00111 }
00112
00113 private:
00114 RedistributeInterface interface;
00115 bool setup_;
00116 };
00117
00126 template<class M>
00127 struct CommMatrixRowSize
00128 {
00129
00130 typedef typename M::size_type value_type;
00131 typedef typename M::size_type size_type;
00132
00138 CommMatrixRowSize(const M& m_, std::vector<size_type>& rowsize_)
00139 : matrix(m_), rowsize(rowsize_)
00140 {}
00141 const M& matrix;
00142 std::vector<size_type>& rowsize;
00143
00144 };
00145
00146
00155 template<class M, class I>
00156 struct CommMatrixSparsityPattern
00157 {
00158 typedef typename M::size_type size_type;
00159
00166 CommMatrixSparsityPattern(M& m_, const Dune::GlobalLookupIndexSet<I>& idxset_, const I& aggidxset_)
00167 : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), rowsize()
00168 {}
00169
00177 CommMatrixSparsityPattern(M& m_, const Dune::GlobalLookupIndexSet<I>& idxset_, const I& aggidxset_,
00178 const std::vector<typename M::size_type>& rowsize_)
00179 : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), sparsity(aggidxset_.size()), rowsize(&rowsize_)
00180 {}
00181
00188 void storeSparsityPattern(M& m)
00189 {
00190
00191 typedef typename Dune::GlobalLookupIndexSet<I>::const_iterator IIter;
00192 typedef typename Dune::OwnerOverlapCopyCommunication<int>::OwnerSet OwnerSet;
00193 std::size_t nnz=0;
00194 #ifdef DEBUG_REPART
00195 int rank;
00196
00197 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
00198 #endif
00199 for(IIter i= aggidxset.begin(), end=aggidxset.end(); i!=end; ++i){
00200 if(!OwnerSet::contains(i->local().attribute())){
00201 #ifdef DEBUG_REPART
00202 std::cout<<rank<<" Inserting diagonal for"<<i->local()<<std::endl;
00203 #endif
00204 sparsity[i->local()].insert(i->local());
00205 }
00206
00207 nnz+=sparsity[i->local()].size();
00208 }
00209 assert( aggidxset.size()==sparsity.size());
00210
00211 if(nnz>0){
00212 m.setSize(aggidxset.size(), aggidxset.size(), nnz);
00213 m.setBuildMode(M::row_wise);
00214 typename M::CreateIterator citer=m.createbegin();
00215 #ifdef DEBUG_REPART
00216 std::size_t idx=0;
00217 bool correct=true;
00218 Dune::GlobalLookupIndexSet<I> global(aggidxset);
00219 #endif
00220 typedef typename std::vector<std::set<size_type> >::const_iterator Iter;
00221 for(Iter i=sparsity.begin(), end=sparsity.end(); i!=end; ++i, ++citer)
00222 {
00223 typedef typename std::set<size_type>::const_iterator SIter;
00224 for(SIter si=i->begin(), send=i->end(); si!=send; ++si)
00225 citer.insert(*si);
00226 #ifdef DEBUG_REPART
00227 if(i->find(idx)==i->end()){
00228 const typename I::IndexPair* gi=global.pair(idx);
00229 assert(gi);
00230 std::cout<<rank<<": row "<<idx<<" is missing a diagonal entry! global="<<gi->global()<<" attr="<<gi->local().attribute()<<" "<<
00231 OwnerSet::contains(gi->local().attribute())<<
00232 " row size="<<i->size()<<std::endl;
00233 correct=false;
00234 }
00235 ++idx;
00236 #endif
00237 }
00238 #ifdef DEBUG_REPART
00239 if(!correct)
00240 throw "bla";
00241 #endif
00242 }
00243 }
00244
00245 M& matrix;
00246 typedef Dune::GlobalLookupIndexSet<I> LookupIndexSet;
00247 const Dune::GlobalLookupIndexSet<I>& idxset;
00248 const I& aggidxset;
00249 std::vector<std::set<size_type> > sparsity;
00250 const std::vector<size_type>* rowsize;
00251 };
00252
00253 template<class M, class I>
00254 struct CommPolicy<CommMatrixSparsityPattern<M,I> >
00255 {
00256 typedef CommMatrixSparsityPattern<M,I> Type;
00257
00262 typedef typename I::GlobalIndex IndexedType;
00263
00265 typedef VariableSize IndexedTypeFlag;
00266
00267 static typename M::size_type getSize(const Type& t, std::size_t i)
00268 {
00269 if(!t.rowsize)
00270 return t.matrix[i].size();
00271 else
00272 {
00273 assert((*t.rowsize)[i]>0);
00274 return (*t.rowsize)[i];
00275 }
00276 }
00277 };
00278
00285 template<class M, class I>
00286 struct CommMatrixRow
00287 {
00296 CommMatrixRow(M& m_, const Dune::GlobalLookupIndexSet<I>& idxset_, const I& aggidxset_)
00297 : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), rowsize()
00298 {}
00299
00303 CommMatrixRow(M& m_, const Dune::GlobalLookupIndexSet<I>& idxset_, const I& aggidxset_,
00304 std::vector<size_t>& rowsize_)
00305 : matrix(m_), idxset(idxset_), aggidxset(aggidxset_), rowsize(&rowsize_)
00306 {}
00312 void setOverlapRowsToDirichlet()
00313 {
00314 typedef typename Dune::GlobalLookupIndexSet<I>::const_iterator Iter;
00315 typedef typename Dune::OwnerOverlapCopyCommunication<int>::OwnerSet OwnerSet;
00316
00317 for(Iter i= aggidxset.begin(), end=aggidxset.end(); i!=end; ++i)
00318 if(!OwnerSet::contains(i->local().attribute())){
00319
00320 typedef typename M::ColIterator CIter;
00321 for(CIter c=matrix[i->local()].begin(), cend= matrix[i->local()].end();
00322 c!= cend; ++c)
00323 {
00324 *c=0;
00325 if(c.index()==i->local()){
00326 typedef typename M::block_type::RowIterator RIter;
00327 for(RIter r=c->begin(), rend=c->end();
00328 r != rend; ++r)
00329 (*r)[r.index()]=1;
00330 }
00331 }
00332 }
00333 }
00335 M& matrix;
00337 const Dune::GlobalLookupIndexSet<I>& idxset;
00339 const I& aggidxset;
00341 std::vector<size_t>* rowsize;
00342 };
00343
00344 template<class M, class I>
00345 struct CommPolicy<CommMatrixRow<M,I> >
00346 {
00347 typedef CommMatrixRow<M,I> Type;
00348
00353 typedef std::pair<typename I::GlobalIndex,typename M::block_type> IndexedType;
00354
00356 typedef VariableSize IndexedTypeFlag;
00357
00358 static std::size_t getSize(const Type& t, std::size_t i)
00359 {
00360 if(!t.rowsize)
00361 return t.matrix[i].size();
00362 else
00363 {
00364 assert((*t.rowsize)[i]>0);
00365 return (*t.rowsize)[i];
00366 }
00367 }
00368 };
00369
00370 template<class M, class I>
00371 struct MatrixRowSizeGatherScatter
00372 {
00373 typedef CommMatrixRowSize<M> Container;
00374
00375 static const typename M::size_type gather(const Container& cont, std::size_t i)
00376 {
00377 return cont.matrix[i].size();
00378 }
00379 static void scatter(Container& cont, const typename M::size_type& rowsize, std::size_t i)
00380 {
00381 assert(rowsize);
00382 cont.rowsize[i]=rowsize;
00383 }
00384
00385 };
00386 template<class M, class I>
00387 struct MatrixSparsityPatternGatherScatter
00388 {
00389 typedef typename I::GlobalIndex GlobalIndex;
00390 typedef CommMatrixSparsityPattern<M,I> Container;
00391 typedef typename M::ConstColIterator ColIter;
00392
00393 static ColIter col;
00394
00395 static const GlobalIndex& gather(const Container& cont, std::size_t i, std::size_t j)
00396 {
00397 if(j==0)
00398 col=cont.matrix[i].begin();
00399 else
00400 ++col;
00401 assert(col!=cont.matrix[i].end());
00402 const typename I::IndexPair* index=cont.idxset.pair(col.index());
00403 assert(index);
00404 return index->global();
00405 }
00406 static void scatter(Container& cont, const GlobalIndex& gi, std::size_t i, std::size_t j)
00407 {
00408 try{
00409 const typename I::IndexPair& ip=cont.aggidxset.at(gi);
00410 assert(ip.global()==gi);
00411 std::size_t col = ip.local();
00412 cont.sparsity[i].insert(col);
00413
00414 typedef typename Dune::OwnerOverlapCopyCommunication<int>::OwnerSet OwnerSet;
00415 if(!OwnerSet::contains(ip.local().attribute()))
00416
00417 cont.sparsity[col].insert(i);
00418 }
00419 catch(Dune::RangeError er){
00420
00421 #ifdef DEBUG_REPART
00422 typedef typename Container::LookupIndexSet GlobalLookup;
00423 typedef typename GlobalLookup::IndexPair IndexPair;
00424 typedef typename Dune::OwnerOverlapCopyCommunication<int>::OwnerSet OwnerSet;
00425
00426 GlobalLookup lookup(cont.aggidxset);
00427 const IndexPair* pi=lookup.pair(i);
00428 assert(pi);
00429 if(OwnerSet::contains(pi->local().attribute())){
00430 int rank;
00431 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
00432 std::cout<<rank<<cont.aggidxset<<std::endl;
00433 std::cout<<rank<<": row "<<i<<" (global="<<gi <<") not in index set for owner index "<<pi->global()<<std::endl;
00434 throw er;
00435 }
00436 #endif
00437 }
00438 }
00439
00440 };
00441 template<class M, class I>
00442 typename MatrixSparsityPatternGatherScatter<M,I>::ColIter MatrixSparsityPatternGatherScatter<M,I>::col;
00443
00444 template<class M, class I>
00445 struct MatrixRowGatherScatter
00446 {
00447 typedef typename I::GlobalIndex GlobalIndex;
00448 typedef CommMatrixRow<M,I> Container;
00449 typedef typename M::ConstColIterator ColIter;
00450 typedef typename std::pair<GlobalIndex,typename M::block_type> Data;
00451 static ColIter col;
00452 static Data datastore;
00453
00454 static const Data& gather(const Container& cont, std::size_t i, std::size_t j)
00455 {
00456 if(j==0)
00457 col=cont.matrix[i].begin();
00458 else
00459 ++col;
00460
00461 const typename I::IndexPair* index=cont.idxset.pair(col.index());
00462 assert(index);
00463
00464 datastore = std::make_pair(index->global(),*col);
00465 return datastore;
00466 }
00467 static void scatter(Container& cont, const Data& data, std::size_t i, std::size_t j)
00468 {
00469 try{
00470 typename M::size_type column=cont.aggidxset.at(data.first).local();
00471 cont.matrix[i][column]=data.second;
00472 }
00473 catch(Dune::RangeError er){
00474
00475 }
00476
00477 }
00478 };
00479
00480 template<class M, class I>
00481 typename MatrixRowGatherScatter<M,I>::ColIter MatrixRowGatherScatter<M,I>::col;
00482
00483 template<class M, class I>
00484 typename MatrixRowGatherScatter<M,I>::Data MatrixRowGatherScatter<M,I>::datastore;
00485
00502 template<typename M, typename C>
00503 void redistributeMatrix(M& origMatrix, M& newMatrix, C& origComm, C& newComm,
00504 RedistributeInformation<C>& ri)
00505 {
00506 std::vector<typename M::size_type> rowsize(newComm.indexSet().size(), 0);
00507
00508 typedef typename C::ParallelIndexSet IndexSet;
00509 CommMatrixRowSize<M> commRowSize(origMatrix, rowsize);
00510 ri.template redistribute<MatrixRowSizeGatherScatter<M,IndexSet> >(commRowSize,commRowSize);
00511
00512 origComm.buildGlobalLookup();
00513
00514 CommMatrixSparsityPattern<M,IndexSet> origsp(origMatrix, origComm.globalLookup(), newComm.indexSet());
00515 CommMatrixSparsityPattern<M,IndexSet> newsp(origMatrix, origComm.globalLookup(), newComm.indexSet(), rowsize);
00516
00517 ri.template redistribute<MatrixSparsityPatternGatherScatter<M,IndexSet> >(origsp,newsp);
00518
00519 newsp.storeSparsityPattern(newMatrix);
00520
00521 #ifdef DUNE_ISTL_WITH_CHECKING
00522
00523 int ret=0;
00524 typedef typename M::ConstRowIterator RIter;
00525 for(RIter row=newMatrix.begin(), rend=newMatrix.end(); row != rend; ++row){
00526 typedef typename M::ConstColIterator CIter;
00527 for(CIter col=row->begin(), cend=row->end(); col!=cend; ++col)
00528 {
00529 try{
00530 newMatrix[col.index()][row.index()];
00531 }catch(Dune::ISTLError e){
00532 std::cerr<<newComm.communicator().rank()<<": entry ("
00533 <<col.index()<<","<<row.index()<<") missing! for symmetry!"<<std::endl;
00534 ret=1;
00535
00536 }
00537
00538 }
00539 }
00540
00541 if(ret)
00542 DUNE_THROW(ISTLError, "Matrix not symmetric!");
00543 #endif
00544
00545 CommMatrixRow<M,IndexSet> origrow(origMatrix, origComm.globalLookup(), newComm.indexSet());
00546 CommMatrixRow<M,IndexSet> newrow(newMatrix, origComm.globalLookup(), newComm.indexSet(),rowsize);
00547
00548 ri.template redistribute<MatrixRowGatherScatter<M,IndexSet> >(origrow,newrow);
00549 newrow.setOverlapRowsToDirichlet();
00550 if(newMatrix.N()>0&&newMatrix.N()<20)
00551 printmatrix(std::cout, newMatrix, "redist", "row");
00552 }
00553
00554 #endif
00555 }
00556 #endif