3 #ifndef DUNE_ISTL_REPARTITION_HH 4 #define DUNE_ISTL_REPARTITION_HH 20 #include <dune/common/timer.hh> 21 #include <dune/common/unused.hh> 22 #include <dune/common/enumset.hh> 23 #include <dune/common/stdstreams.hh> 24 #include <dune/common/parallel/mpitraits.hh> 25 #include <dune/common/parallel/communicator.hh> 26 #include <dune/common/parallel/indexset.hh> 27 #include <dune/common/parallel/indicessyncer.hh> 28 #include <dune/common/parallel/remoteindices.hh> 29 #include <dune/common/rangeutilities.hh> 48 #if HAVE_PARMETIS && defined(REALTYPEWIDTH) 54 #if HAVE_PARMETIS && defined(IDXTYPEWIDTH) 56 #elif HAVE_PARMETIS && defined(HAVE_SCOTCH_NUM_TYPE) 57 using idx_t = SCOTCH_Num;
80 template<
class G,
class T1,
class T2>
84 typedef typename IndexSet::LocalIndex::Attribute Attribute;
86 IndexSet& indexSet = oocomm.
indexSet();
90 typedef typename G::ConstVertexIterator VertexIterator;
93 std::size_t sum=0, needed = graph.noVertices()-indexSet.size();
94 std::vector<std::size_t> neededall(oocomm.
communicator().size(), 0);
96 MPI_Allgather(&needed, 1, MPITraits<std::size_t>::getType() , &(neededall[0]), 1, MPITraits<std::size_t>::getType(), oocomm.
communicator());
106 typedef typename IndexSet::const_iterator Iterator;
108 end = indexSet.end();
109 for(Iterator it = indexSet.begin(); it != end; ++it)
110 maxgi=std::max(maxgi,it->global());
119 maxgi=maxgi+neededall[i];
122 std::map<int,SLList<std::pair<T1,Attribute> > > globalIndices;
123 storeGlobalIndicesOfRemoteIndices(globalIndices, oocomm.
remoteIndices());
124 indexSet.beginResize();
126 for(VertexIterator vertex = graph.begin(), vend=graph.end(); vertex != vend; ++vertex) {
127 const typename IndexSet::IndexPair* pair=lookup.pair(*vertex);
135 indexSet.endResize();
137 repairLocalIndexPointers(globalIndices, oocomm.
remoteIndices(), indexSet);
142 std::cout<<
"Holes are filled!"<<std::endl;
150 class ParmetisDuneIndexMap
153 template<
class Graph,
class OOComm>
154 ParmetisDuneIndexMap(
const Graph& graph,
const OOComm& com);
155 int toParmetis(
int i)
const 157 return duneToParmetis[i];
159 int toLocalParmetis(
int i)
const 161 return duneToParmetis[i]-base_;
163 int operator[](
int i)
const 165 return duneToParmetis[i];
167 int toDune(
int i)
const 169 return parmetisToDune[i];
171 std::vector<int>::size_type numOfOwnVtx()
const 173 return parmetisToDune.size();
182 std::vector<int> duneToParmetis;
183 std::vector<int> parmetisToDune;
185 std::vector<Metis::idx_t> vtxDist_;
188 template<
class G,
class OOComm>
189 ParmetisDuneIndexMap::ParmetisDuneIndexMap(
const G& graph,
const OOComm& oocomm)
190 : duneToParmetis(graph.noVertices(), -1), vtxDist_(oocomm.communicator().size()+1)
192 int npes=oocomm.communicator().size(), mype=oocomm.communicator().rank();
194 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
195 typedef typename OOComm::OwnerSet OwnerSet;
198 Iterator end = oocomm.indexSet().end();
199 for(Iterator index = oocomm.indexSet().begin(); index != end; ++index) {
200 if (OwnerSet::contains(index->local().attribute())) {
204 parmetisToDune.resize(numOfOwnVtx);
205 std::vector<int> globalNumOfVtx(npes);
207 MPI_Allgather(&numOfOwnVtx, 1, MPI_INT, &(globalNumOfVtx[0]), 1, MPI_INT, oocomm.communicator());
211 for(
int i=0; i<npes; i++) {
213 base += globalNumOfVtx[i];
215 vtxDist_[i+1] = vtxDist_[i] + globalNumOfVtx[i];
221 typedef typename G::ConstVertexIterator VertexIterator;
223 std::cout << oocomm.communicator().rank()<<
" vtxDist: ";
224 for(
int i=0; i<= npes; ++i)
225 std::cout << vtxDist_[i]<<
" ";
226 std::cout<<std::endl;
233 VertexIterator vend = graph.end();
234 for(VertexIterator vertex = graph.begin(); vertex != vend; ++vertex) {
235 const typename OOComm::ParallelIndexSet::IndexPair* index=oocomm.globalLookup().pair(*vertex);
237 if (OwnerSet::contains(index->local().attribute())) {
239 parmetisToDune[base-base_]=index->local();
240 duneToParmetis[index->local()] = base++;
250 std::cout <<oocomm.communicator().rank()<<
": before ";
251 for(std::size_t i=0; i<duneToParmetis.size(); ++i)
252 std::cout<<duneToParmetis[i]<<
" ";
253 std::cout<<std::endl;
255 oocomm.copyOwnerToAll(duneToParmetis,duneToParmetis);
257 std::cout <<oocomm.communicator().rank()<<
": after ";
258 for(std::size_t i=0; i<duneToParmetis.size(); ++i)
259 std::cout<<duneToParmetis[i]<<
" ";
260 std::cout<<std::endl;
272 template<
class Flags,
class IS>
275 std::map<int,int> sizes;
277 typedef typename IS::const_iterator IIter;
278 for(IIter i=idxset.begin(), end=idxset.end(); i!=end; ++i)
279 if(Flags::contains(i->local().attribute()))
280 ++sizes[toPart[i->local()]];
283 typedef std::map<int,int>::const_iterator MIter;
284 for(MIter i=sizes.begin(), end=sizes.end(); i!=end; ++i)
285 interfaces()[i->first].first.reserve(i->second);
288 typedef typename IS::const_iterator IIter;
289 for(IIter i=idxset.begin(), end=idxset.end(); i!=end; ++i)
290 if(Flags::contains(i->local().attribute()))
291 interfaces()[toPart[i->local()]].first.add(i->local());
296 interfaces()[proc].second.reserve(size);
300 interfaces()[proc].second.add(idx);
302 template<
typename TG>
305 typedef typename std::vector<std::pair<TG,int> >::const_iterator VIter;
307 for(VIter idx=indices.begin(); idx!= indices.end(); ++idx) {
308 interfaces()[idx->second].second.add(i++);
329 void createSendBuf(std::vector<GI>& ownerVec, std::set<GI>& overlapVec, std::set<int>& neighbors,
char *sendBuf,
int buffersize, MPI_Comm comm) {
331 std::size_t s=ownerVec.size();
335 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
336 MPI_Pack(&(ownerVec[0]), s, MPITraits<GI>::getType(), sendBuf, buffersize, &pos, comm);
337 s = overlapVec.size();
338 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
339 typedef typename std::set<GI>::iterator Iter;
340 for(Iter i=overlapVec.begin(), end= overlapVec.end(); i != end; ++i)
341 MPI_Pack(const_cast<GI*>(&(*i)), 1, MPITraits<GI>::getType(), sendBuf, buffersize, &pos, comm);
344 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
345 typedef typename std::set<int>::iterator IIter;
347 for(IIter i=neighbors.begin(), end= neighbors.end(); i != end; ++i)
348 MPI_Pack(const_cast<int*>(&(*i)), 1, MPI_INT, sendBuf, buffersize, &pos, comm);
359 void saveRecvBuf(
char *recvBuf,
int bufferSize, std::vector<std::pair<GI,int> >& ownerVec,
360 std::set<GI>& overlapVec, std::set<int>& neighbors,
RedistributeInterface& inf,
int from, MPI_Comm comm) {
364 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
366 ownerVec.reserve(ownerVec.size()+size);
367 for(; size!=0; --size) {
369 MPI_Unpack(recvBuf, bufferSize, &pos, &gi, 1, MPITraits<GI>::getType(), comm);
370 ownerVec.push_back(std::make_pair(gi,from));
373 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
374 typename std::set<GI>::iterator ipos = overlapVec.begin();
375 Dune::dverb <<
"unpacking "<<size<<
" overlap"<<std::endl;
376 for(; size!=0; --size) {
378 MPI_Unpack(recvBuf, bufferSize, &pos, &gi, 1, MPITraits<GI>::getType(), comm);
379 ipos=overlapVec.insert(ipos, gi);
382 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
383 Dune::dverb <<
"unpacking "<<size<<
" neighbors"<<std::endl;
384 typename std::set<int>::iterator npos = neighbors.begin();
385 for(; size!=0; --size) {
387 MPI_Unpack(recvBuf, bufferSize, &pos, &n, 1, MPI_INT, comm);
388 npos=neighbors.insert(npos, n);
406 void getDomain(
const MPI_Comm& comm, T *part,
int numOfOwnVtx,
int nparts,
int *myDomain, std::vector<int> &domainMapping) {
408 MPI_Comm_size(comm, &npes);
409 MPI_Comm_rank(comm, &mype);
416 std::vector<int> domain(nparts, 0);
417 std::vector<int> assigned(npes, 0);
419 domainMapping.assign(domainMapping.size(), -1);
422 for (i=0; i<numOfOwnVtx; i++) {
426 std::vector<int> domainMatrix(npes * nparts, -1);
429 int *buf =
new int[nparts];
430 for (i=0; i<nparts; i++) {
432 domainMatrix[mype*nparts+i] = domain[i];
435 int src = (mype-1+npes)%npes;
436 int dest = (mype+1)%npes;
438 for (i=0; i<npes-1; i++) {
439 MPI_Sendrecv_replace(buf, nparts, MPI_INT, dest, 0, src, 0, comm, &status);
441 pe = ((mype-1-i)+npes)%npes;
442 for(j=0; j<nparts; j++) {
444 domainMatrix[pe*nparts+j] = buf[j];
452 int maxOccurance = 0;
454 std::set<std::size_t> unassigned;
456 for(i=0; i<nparts; i++) {
457 for(j=0; j<npes; j++) {
459 if (assigned[j]==0) {
460 if (maxOccurance < domainMatrix[j*nparts+i]) {
461 maxOccurance = domainMatrix[j*nparts+i];
469 domainMapping[i] = pe;
479 unassigned.insert(i);
484 typename std::vector<int>::iterator next_free = assigned.begin();
486 for(
typename std::set<std::size_t>::iterator domain = unassigned.begin(),
487 end = unassigned.end(); domain != end; ++domain)
489 next_free = std::find_if(next_free, assigned.end(), std::bind(std::less<int>(), std::placeholders::_1, 1));
490 assert(next_free != assigned.end());
491 domainMapping[*domain] = next_free-assigned.begin();
499 bool operator()(
const T& t1,
const T& t2)
const 517 void mergeVec(std::vector<std::pair<GI, int> >& ownerVec, std::set<GI>& overlapSet) {
519 typedef typename std::vector<std::pair<GI,int> >::const_iterator VIter;
522 if(ownerVec.size()>0)
524 VIter old=ownerVec.begin();
525 for(VIter i=old+1, end=ownerVec.end(); i != end; old=i++)
527 if(i->first==old->first)
529 std::cerr<<
"Value at indes"<<old-ownerVec.begin()<<
" is the same as at index " 530 <<i-ownerVec.begin()<<
" ["<<old->first<<
","<<old->second<<
"]==[" 531 <<i->first<<
","<<i->second<<
"]"<<std::endl;
539 typedef typename std::set<GI>::iterator SIter;
540 VIter v=ownerVec.begin(), vend=ownerVec.end();
541 for(SIter s=overlapSet.begin(), send=overlapSet.end(); s!=send;)
543 while(v!=vend && v->first<*s) ++v;
544 if(v!=vend && v->first==*s) {
549 overlapSet.erase(tmp);
569 template<
class OwnerSet,
class Graph,
class IS,
class GI>
570 void getNeighbor(
const Graph& g, std::vector<int>& part,
571 typename Graph::VertexDescriptor vtx,
const IS& indexSet,
572 int toPe, std::set<GI>& neighbor, std::set<int>& neighborProcs) {
573 typedef typename Graph::ConstEdgeIterator Iter;
574 for(Iter edge=g.beginEdges(vtx), end=g.endEdges(vtx); edge!=end; ++edge)
576 const typename IS::IndexPair* pindex = indexSet.pair(edge.target());
578 if(part[pindex->local()]!=toPe || !OwnerSet::contains(pindex->local().attribute()))
581 neighbor.insert(pindex->global());
582 neighborProcs.insert(part[pindex->local()]);
587 template<
class T,
class I>
588 void my_push_back(std::vector<T>& ownerVec,
const I& index,
int proc)
590 DUNE_UNUSED_PARAMETER(proc);
591 ownerVec.push_back(index);
594 template<
class T,
class I>
595 void my_push_back(std::vector<std::pair<T,int> >& ownerVec,
const I& index,
int proc)
597 ownerVec.push_back(std::make_pair(index,proc));
626 template<
class OwnerSet,
class G,
class IS,
class T,
class GI>
627 void getOwnerOverlapVec(
const G& graph, std::vector<int>& part, IS& indexSet,
628 int myPe,
int toPe, std::vector<T>& ownerVec, std::set<GI>& overlapSet,
630 DUNE_UNUSED_PARAMETER(myPe);
632 typedef typename IS::const_iterator Iterator;
633 for(Iterator index = indexSet.begin(); index != indexSet.end(); ++index) {
635 if(OwnerSet::contains(index->local().attribute()))
637 if(part[index->local()]==toPe)
639 getNeighbor<OwnerSet>(graph, part, index->local(), indexSet,
640 toPe, overlapSet, neighborProcs);
641 my_push_back(ownerVec, index->global(), toPe);
645 reserve(ownerVec, redist, toPe);
656 template<
class F,
class IS>
657 inline bool isOwner(IS& indexSet,
int index) {
659 const typename IS::IndexPair* pindex=indexSet.pair(index);
662 return F::contains(pindex->local().attribute());
666 class BaseEdgeFunctor
669 BaseEdgeFunctor(
Metis::idx_t* adj,
const ParmetisDuneIndexMap& data)
670 : i_(), adj_(adj), data_(data)
674 void operator()(
const T& edge)
678 adj_[i_] = data_.toParmetis(edge.target());
689 const ParmetisDuneIndexMap& data_;
694 :
public BaseEdgeFunctor
696 EdgeFunctor(
Metis::idx_t* adj,
const ParmetisDuneIndexMap& data, std::size_t)
697 : BaseEdgeFunctor(adj, data)
707 template<
class G,
class V,
class E,
class VM,
class EM>
709 :
public BaseEdgeFunctor
712 EdgeFunctor(
Metis::idx_t* adj,
const ParmetisDuneIndexMap& data, std::size_t s)
713 : BaseEdgeFunctor(adj, data)
719 void operator()(
const T& edge)
721 weight_[index()]=edge.properties().depends() ? 3 : 1;
722 BaseEdgeFunctor::operator()(edge);
753 template<
class F,
class G,
class IS,
class EW>
754 void getAdjArrays(G& graph, IS& indexSet,
Metis::idx_t *xadj,
760 typedef typename G::ConstVertexIterator VertexIterator;
762 typedef typename IS::const_iterator Iterator;
764 VertexIterator vend = graph.end();
767 for(VertexIterator vertex = graph.begin(); vertex != vend; ++vertex) {
768 if (isOwner<F>(indexSet,*vertex)) {
770 typedef typename G::ConstEdgeIterator EdgeIterator;
771 EdgeIterator eend = vertex.end();
772 xadj[j] = ew.index();
774 for(EdgeIterator edge = vertex.begin(); edge != eend; ++edge) {
779 xadj[j] = ew.index();
783 template<
class G,
class T1,
class T2>
790 #ifndef METIS_VER_MAJOR 795 Metis::idx_t *adjwgt,
int *wgtflag,
int *numflag,
int *nparts,
799 Metis::idx_t *adjwgt,
int *wgtflag,
int *numflag,
int *nparts,
803 #endif // HAVE_PARMETIS 805 template<
class S,
class T>
808 for(T *cur=array, *end=array+l; cur!=end; ++cur)
812 template<
class S,
class T>
813 inline bool isValidGraph(std::size_t noVtx, std::size_t gnoVtx, S noEdges, T* xadj,
814 T* adjncy,
bool checkSymmetry)
819 if(xadj[vtx]>noEdges||xadj[vtx]<0) {
820 std::cerr <<
"Check graph: xadj["<<vtx<<
"]="<<xadj[vtx]<<
" (>" 821 <<noEdges<<
") out of range!"<<std::endl;
824 if(xadj[vtx+1]>noEdges||xadj[vtx+1]<0) {
825 std::cerr <<
"Check graph: xadj["<<vtx+1<<
"]="<<xadj[vtx+1]<<
" (>" 826 <<noEdges<<
") out of range!"<<std::endl;
831 if(adjncy[i]<0||((std::size_t)adjncy[i])>gnoVtx) {
832 std::cerr<<
" Edge "<<adjncy[i]<<
" out of range ["<<0<<
","<<noVtx<<
")" 842 for(
Metis::idx_t j=xadj[target]; j< xadj[target+1]; ++j)
846 std::cerr<<
"Edge ("<<target<<
","<<vtx<<
") "<<i<<
" time"<<std::endl;
855 template<
class M,
class T1,
class T2>
863 std::cout<<
"Repartitioning from "<<oocomm.
communicator().size()
864 <<
" to "<<nparts<<
" parts"<<std::endl;
868 int* part =
new int[1];
883 typedef typename RemoteIndices::const_iterator
912 xadj[1]=noNeighbours;
942 typedef typename RemoteIndices::const_iterator NeighbourIterator;
956 if(n->first != rank) {
966 noNeighbours, xadj, adjncy,
false));
974 for(
int i=0; i<nparts; ++i)
975 tpwgts[i]=1.0/nparts;
976 int options[5] ={ 0,1,15,0,0};
979 Dune::dinfo<<rank<<
" vtxdist: ";
981 Dune::dinfo<<std::endl<<rank<<
" xadj: ";
983 Dune::dinfo<<std::endl<<rank<<
" adjncy: ";
987 Dune::dinfo<<std::endl<<rank<<
" vwgt: ";
989 Dune::dinfo<<std::endl<<rank<<
" adwgt: ";
992 Dune::dinfo<<std::endl;
995 std::cout<<
"Creating comm graph took "<<time.elapsed()<<std::endl;
998 #ifdef PARALLEL_PARTITION 1005 ParMETIS_V3_PartKway(vtxdist, xadj, adjncy,
1006 vwgt, adjwgt, &wgtflag,
1007 &numflag, &ncon, &nparts, tpwgts, &ubvec, options, &edgecut, part,
1010 std::cout<<
"ParMETIS took "<<time.elapsed()<<std::endl;
1014 std::size_t gnoedges=0;
1017 Dune::dverb<<
"noNeighbours: "<<noNeighbours<<std::endl;
1019 MPI_Allgather(&noNeighbours,1,MPI_INT,noedges,1, MPI_INT,oocomm.
communicator());
1022 std::cout<<
"Gathering noedges took "<<time1.elapsed()<<std::endl;
1037 std::size_t localNoVtx=vtxdist[rank+1]-vtxdist[rank];
1042 Dune::dinfo<<
"noedges: ";
1044 Dune::dinfo<<std::endl;
1052 noxs[i]=vtxdist[i+1]-vtxdist[i]+1;
1053 novs[i]=vtxdist[i+1]-vtxdist[i];
1058 for(
int *xcurr = xdispl, *vcurr = vdispl, *end=vdispl+oocomm.
communicator().size();
1059 vcurr!=end; ++vcurr, ++xcurr, ++so, ++offset) {
1061 *xcurr = offset + *so;
1067 for(
int *curr=noedges, *end=noedges+oocomm.
communicator().size()-1;
1068 curr!=end; ++curr) {
1073 Dune::dinfo<<
"displ: ";
1075 Dune::dinfo<<std::endl;
1079 for(
int *curr=noedges, *end=noedges+oocomm.
communicator().size();
1084 Dune::dinfo<<
"gxadjlen: "<<gxadjlen<<
" noVertices: "<<noVertices
1085 <<
" gnoedges: "<<gnoedges<<std::endl;
1096 std::cout<<
"Preparing global graph took "<<time1.elapsed()<<std::endl;
1100 MPI_Allgatherv(xadj,2,MPITraits<Metis::idx_t>::getType(),
1101 gxadj,noxs,xdispl,MPITraits<Metis::idx_t>::getType(),
1103 MPI_Allgatherv(adjncy,noNeighbours,MPITraits<Metis::idx_t>::getType(),
1104 gadjncy,noedges,displ,MPITraits<Metis::idx_t>::getType(),
1107 MPI_Allgatherv(adjwgt,noNeighbours,MPITraits<Metis::idx_t>::getType(),
1108 gadjwgt,noedges,displ,MPITraits<Metis::idx_t>::getType(),
1110 MPI_Allgatherv(vwgt,localNoVtx,MPITraits<Metis::idx_t>::getType(),
1111 gvwgt,novs,vdispl,MPITraits<Metis::idx_t>::getType(),
1115 std::cout<<
"Gathering global graph data took "<<time1.elapsed()<<std::endl;
1129 int lprev = vtxdist[i]-vtxdist[i-1];
1130 int l = vtxdist[i+1]-vtxdist[i];
1132 assert((start+l+offset)-gxadj<=static_cast<Metis::idx_t>(gxadjlen));
1133 increment = *(start-1);
1134 std::transform(start+offset, start+l+offset, start, std::bind(std::plus<Metis::idx_t>(), std::placeholders::_1, increment));
1136 Dune::dinfo<<std::endl<<
"shifted xadj:";
1138 Dune::dinfo<<std::endl<<
" gadjncy: ";
1141 Dune::dinfo<<std::endl<<
" gvwgt: ";
1143 Dune::dinfo<<std::endl<<
"adjwgt: ";
1145 Dune::dinfo<<std::endl;
1149 std::cout<<
"Postprocesing global graph data took "<<time1.elapsed()<<std::endl;
1153 gxadj, gadjncy,
true));
1157 std::cout<<
"Creating grah one 1 process took "<<time.elapsed()<<std::endl;
1159 options[0]=0; options[1]=1; options[2]=1; options[3]=3; options[4]=3;
1160 #if METIS_VER_MAJOR >= 5 1163 METIS_SetDefaultOptions(moptions);
1164 moptions[METIS_OPTION_NUMBERING] = numflag;
1165 METIS_PartGraphRecursive(&noVertices, &ncon, gxadj, gadjncy, gvwgt, NULL, gadjwgt,
1166 &nparts, NULL, NULL, moptions, &edgecut, gpart);
1169 METIS_PartGraphRecursive(&noVertices, gxadj, gadjncy, gvwgt, gadjwgt, &wgtflag,
1170 &numflag, &nparts, options, &edgecut, gpart);
1174 std::cout<<
"METIS took "<<time.elapsed()<<std::endl;
1177 Dune::dinfo<<std::endl<<
"part:";
1188 MPI_Scatter(gpart, 1, MPITraits<Metis::idx_t>::getType(), part, 1,
1189 MPITraits<Metis::idx_t>::getType(), 0, comm);
1213 Dune::dinfo<<
" repart "<<rank <<
" -> "<< part[0]<<std::endl;
1215 std::vector<int> realpart(mat.N(), part[0]);
1221 std::cout<<
"Scattering repartitioning took "<<time.elapsed()<<std::endl;
1229 std::cout<<
"Filling index set took "<<time.elapsed()<<std::endl;
1237 std::cout<<
"Average no neighbours was "<<noNeighbours<<std::endl;
1242 std::cout<<
"Building index sets took "<<time.elapsed()<<std::endl;
1264 template<
class G,
class T1,
class T2>
1277 std::cout<<
"Filling holes took "<<time.elapsed()<<std::endl;
1284 double t1=0.0, t2=0.0, t3=0.0, t4=0.0, tSum=0.0;
1291 assert(nparts<=static_cast<Metis::idx_t>(oocomm.
communicator().size()));
1312 typedef typename OOComm::OwnerSet OwnerSet;
1317 ParmetisDuneIndexMap indexMap(graph,oocomm);
1319 for(std::size_t i=0; i < indexMap.numOfOwnVtx(); ++i)
1324 std::cerr<<
"ParMETIS not activated. Will repartition to 1 domain instead of requested " 1325 <<nparts<<
" domains."<<std::endl;
1334 EdgeFunctor<G> ef(adjncy, indexMap, graph.noEdges());
1335 getAdjArrays<OwnerSet>(graph, oocomm.
globalLookup(), xadj, ef);
1341 Metis::idx_t numflag=0, wgtflag=0, options[3], edgecut=0, ncon=1;
1344 for(
int i=0; i<nparts; ++i)
1345 tpwgts[i]=1.0/nparts;
1354 wgtflag = (ef.getWeights()!=NULL) ? 1 : 0;
1362 std::cout<<std::endl;
1363 std::cout<<
"Testing ParMETIS_V3_PartKway with options[1-2] = {" 1364 <<options[1]<<
" "<<options[2]<<
"}, Ncon: " 1365 <<ncon<<
", Nparts: "<<nparts<<std::endl;
1378 std::cout<<
"Preparing for parmetis took "<<time.elapsed()<<std::endl;
1385 ParMETIS_V3_PartKway(indexMap.vtxDist(), xadj, adjncy,
1386 NULL, ef.getWeights(), &wgtflag,
1387 &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgecut, part, &
const_cast<MPI_Comm&
>(comm));
1398 std::cout<<std::endl;
1399 std::cout<<
"ParMETIS_V3_PartKway reported a cut of "<<edgecut<<std::endl;
1400 std::cout<<std::endl;
1402 std::cout<<mype<<
": PARMETIS-Result: ";
1403 for(
int i=0; i < indexMap.vtxDist()[mype+1]-indexMap.vtxDist()[mype]; ++i) {
1404 std::cout<<part[i]<<
" ";
1406 std::cout<<std::endl;
1407 std::cout<<
"Testing ParMETIS_V3_PartKway with options[1-2] = {" 1408 <<options[1]<<
" "<<options[2]<<
"}, Ncon: " 1409 <<ncon<<
", Nparts: "<<nparts<<std::endl;
1422 std::cout<<
"Parmetis took "<<time.elapsed()<<std::endl;
1429 for(std::size_t i=0; i<indexMap.numOfOwnVtx(); ++i)
1439 std::vector<int> domainMapping(nparts);
1441 getDomain(comm, part, indexMap.numOfOwnVtx(), nparts, &myDomain, domainMapping);
1446 std::cout<<mype<<
": myDomain: "<<myDomain<<std::endl;
1447 std::cout<<mype<<
": DomainMapping: ";
1448 for(
auto j : range(nparts)) {
1449 std::cout<<
" do: "<<j<<
" pe: "<<domainMapping[j]<<
" ";
1451 std::cout<<std::endl;
1458 std::vector<int> setPartition(oocomm.
indexSet().size(), -1);
1460 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
1462 for(Iterator index = oocomm.
indexSet().begin(); index != oocomm.
indexSet().end(); ++index)
1463 if(OwnerSet::contains(index->local().attribute())) {
1464 setPartition[index->local()]=domainMapping[part[i++]];
1472 static_cast<int>(SolverCategory::nonoverlapping))
1479 std::cout<<
"Creating indexsets took "<<time.elapsed()<<std::endl;
1486 template<
class G,
class T1,
class T2>
1494 typedef typename OOComm::OwnerSet OwnerSet;
1528 std::set<int> recvFrom;
1534 typedef typename std::vector<int>::const_iterator VIter;
1538 std::set<int> tsendTo;
1539 for(VIter i=setPartition.begin(), iend = setPartition.end(); i!=iend; ++i)
1542 noSendTo = tsendTo.size();
1543 sendTo =
new int[noSendTo];
1544 typedef std::set<int>::const_iterator iterator;
1546 for(iterator i=tsendTo.begin(); i != tsendTo.end(); ++i, ++idx)
1552 int* gsendToDispl =
new int[oocomm.
communicator().size()+1];
1554 MPI_Allgather(&noSendTo, 1, MPI_INT, gnoSend, 1,
1558 int totalNoRecv = 0;
1559 for(
int i=0; i<npes; ++i)
1560 totalNoRecv += gnoSend[i];
1562 int *gsendTo =
new int[totalNoRecv];
1566 for(
int i=0; i<npes; ++i)
1567 gsendToDispl[i+1]=gsendToDispl[i]+gnoSend[i];
1570 MPI_Allgatherv(sendTo, noSendTo, MPI_INT, gsendTo, gnoSend, gsendToDispl,
1574 for(
int proc=0; proc < npes; ++proc)
1575 for(
int i=gsendToDispl[proc]; i < gsendToDispl[proc+1]; ++i)
1576 if(gsendTo[i]==mype)
1577 recvFrom.insert(proc);
1579 bool existentOnNextLevel = recvFrom.size()>0;
1583 delete[] gsendToDispl;
1588 if(recvFrom.size()) {
1589 std::cout<<mype<<
": recvFrom: ";
1590 typedef typename std::set<int>::const_iterator siter;
1591 for(siter i=recvFrom.begin(); i!= recvFrom.end(); ++i) {
1596 std::cout<<std::endl<<std::endl;
1597 std::cout<<mype<<
": sendTo: ";
1598 for(
int i=0; i<noSendTo; i++) {
1599 std::cout<<sendTo[i]<<
" ";
1601 std::cout<<std::endl<<std::endl;
1606 std::cout<<
" Communicating the receive information took "<<
1607 time.elapsed()<<std::endl;
1621 typedef typename OOComm::ParallelIndexSet::GlobalIndex GI;
1622 typedef std::vector<GI> GlobalVector;
1623 std::vector<std::pair<GI,int> > myOwnerVec;
1624 std::set<GI> myOverlapSet;
1625 GlobalVector sendOwnerVec;
1626 std::set<GI> sendOverlapSet;
1627 std::set<int> myNeighbors;
1632 char **sendBuffers=
new char*[noSendTo];
1633 MPI_Request *requests =
new MPI_Request[noSendTo];
1636 for(
int i=0; i < noSendTo; ++i) {
1638 sendOwnerVec.clear();
1639 sendOverlapSet.clear();
1642 std::set<int> neighbors;
1643 getOwnerOverlapVec<OwnerSet>(graph, setPartition, oocomm.
globalLookup(),
1644 mype, sendTo[i], sendOwnerVec, sendOverlapSet, redistInf,
1650 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.
communicator(), &buffersize);
1651 MPI_Pack_size(sendOwnerVec.size(), MPITraits<GI>::getType(), oocomm.
communicator(), &tsize);
1653 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.
communicator(), &tsize);
1655 MPI_Pack_size(sendOverlapSet.size(), MPITraits<GI>::getType(), oocomm.
communicator(), &tsize);
1656 buffersize += tsize;
1657 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.
communicator(), &tsize);
1658 buffersize += tsize;
1659 MPI_Pack_size(neighbors.size(), MPI_INT, oocomm.
communicator(), &tsize);
1660 buffersize += tsize;
1662 sendBuffers[i] =
new char[buffersize];
1665 std::cout<<mype<<
" sending "<<sendOwnerVec.size()<<
" owner and "<<
1666 sendOverlapSet.size()<<
" overlap to "<<sendTo[i]<<
" buffersize="<<buffersize<<std::endl;
1668 createSendBuf(sendOwnerVec, sendOverlapSet, neighbors, sendBuffers[i], buffersize, oocomm.
communicator());
1669 MPI_Issend(sendBuffers[i], buffersize, MPI_PACKED, sendTo[i], 99, oocomm.
communicator(), requests+i);
1675 std::cout<<
" Creating sends took "<<
1676 time.elapsed()<<std::endl;
1681 int noRecv = recvFrom.size();
1682 int oldbuffersize=0;
1687 MPI_Probe(MPI_ANY_SOURCE, 99, oocomm.
communicator(), &stat);
1689 MPI_Get_count(&stat, MPI_PACKED, &buffersize);
1691 if(oldbuffersize<buffersize) {
1694 recvBuf =
new char[buffersize];
1695 oldbuffersize = buffersize;
1697 MPI_Recv(recvBuf, buffersize, MPI_PACKED, stat.MPI_SOURCE, 99, oocomm.
communicator(), &stat);
1698 saveRecvBuf(recvBuf, buffersize, myOwnerVec, myOverlapSet, myNeighbors, redistInf,
1708 MPI_Status *statuses =
new MPI_Status[noSendTo];
1709 int send = MPI_Waitall(noSendTo, requests, statuses);
1712 if(send==MPI_ERR_IN_STATUS) {
1713 std::cerr<<mype<<
": Error in sending :"<<std::endl;
1715 for(
int i=0; i< noSendTo; i++)
1716 if(statuses[i].MPI_ERROR!=MPI_SUCCESS) {
1719 MPI_Error_string(statuses[i].MPI_ERROR, message, &messageLength);
1720 std::cerr<<
" source="<<statuses[i].MPI_SOURCE<<
" message: ";
1721 for(
int j = 0; j < messageLength; j++)
1722 std::cout<<message[j];
1724 std::cerr<<std::endl;
1730 std::cout<<
" Receiving and saving took "<<
1731 time.elapsed()<<std::endl;
1735 for(
int i=0; i < noSendTo; ++i)
1736 delete[] sendBuffers[i];
1738 delete[] sendBuffers;
1753 if (!existentOnNextLevel) {
1755 color= MPI_UNDEFINED;
1757 MPI_Comm outputComm;
1765 std::vector<int> tneighbors;
1766 tneighbors.reserve(myNeighbors.size());
1768 typename OOComm::ParallelIndexSet& outputIndexSet = outcomm->indexSet();
1770 MPI_Allgather(&newrank, 1, MPI_INT, newranks, 1,
1772 typedef typename std::set<int>::const_iterator IIter;
1776 for(IIter i=myNeighbors.begin(), end=myNeighbors.end();
1778 assert(newranks[*i]>=0);
1779 std::cout<<*i<<
"->"<<newranks[*i]<<
" ";
1780 tneighbors.push_back(newranks[*i]);
1782 std::cout<<std::endl;
1784 for(IIter i=myNeighbors.begin(), end=myNeighbors.end();
1786 tneighbors.push_back(newranks[*i]);
1790 myNeighbors.clear();
1795 std::cout<<
" Calculating new neighbours ("<<tneighbors.size()<<
") took "<<
1796 time.elapsed()<<std::endl;
1801 outputIndexSet.beginResize();
1804 std::sort(myOwnerVec.begin(), myOwnerVec.end(), SortFirst());
1808 typedef typename OOComm::ParallelIndexSet::LocalIndex LocalIndex;
1809 typedef typename std::vector<std::pair<GI,int> >::const_iterator VPIter;
1811 for(VPIter g=myOwnerVec.begin(), end =myOwnerVec.end(); g!=end; ++g, ++i ) {
1812 outputIndexSet.add(g->first,LocalIndex(i, OwnerOverlapCopyAttributeSet::owner,
true));
1819 std::cout<<
" Adding owner indices took "<<
1820 time.elapsed()<<std::endl;
1829 mergeVec(myOwnerVec, myOverlapSet);
1833 myOwnerVec.swap(myOwnerVec);
1838 std::cout<<
" Merging indices took "<<
1839 time.elapsed()<<std::endl;
1845 typedef typename std::set<GI>::const_iterator SIter;
1846 for(SIter g=myOverlapSet.begin(), end=myOverlapSet.end(); g!=end; ++g, i++) {
1847 outputIndexSet.add(*g,LocalIndex(i, OwnerOverlapCopyAttributeSet::copy,
true));
1849 myOverlapSet.clear();
1850 outputIndexSet.endResize();
1852 #ifdef DUNE_ISTL_WITH_CHECKING 1854 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
1855 Iterator end = outputIndexSet.end();
1856 for(Iterator index = outputIndexSet.begin(); index != end; ++index) {
1857 if (OwnerSet::contains(index->local().attribute())) {
1868 Iterator index=outputIndexSet.begin();
1871 for(Iterator old = outputIndexSet.begin(); index != end; old=index++) {
1872 if(old->global()>index->global())
1873 DUNE_THROW(
ISTLError,
"Index set's globalindex not sorted correctly");
1880 std::cout<<
" Adding overlap indices took "<<
1881 time.elapsed()<<std::endl;
1886 if(color != MPI_UNDEFINED) {
1887 outcomm->remoteIndices().setNeighbours(tneighbors);
1888 outcomm->remoteIndices().template rebuild<true>();
1898 std::cout<<
" Storing indexsets took "<<
1899 time.elapsed()<<std::endl;
1905 tSum = t1 + t2 + t3 + t4;
1906 std::cout<<std::endl
1907 <<mype<<
": WTime for step 1): "<<t1
1915 return color!=MPI_UNDEFINED;
1919 template<
class G,
class P,
class T1,
class T2,
class R>
1925 if(nparts!=oocomm.size())
1926 DUNE_THROW(NotImplemented,
"only available for MPI programs");
1930 template<
class G,
class P,
class T1,
class T2,
class R>
1936 if(nparts!=oocomm.size())
1937 DUNE_THROW(NotImplemented,
"only available for MPI programs");
void fillIndexSetHoles(const G &graph, Dune::OwnerOverlapCopyCommunication< T1, T2 > &oocomm)
Fills the holes in an index set.
Definition: repartition.hh:81
Definition: allocator.hh:7
void print_carray(S &os, T *array, std::size_t l)
Definition: repartition.hh:806
void buildSendInterface(const std::vector< int > &toPart, const IS &idxset)
Definition: repartition.hh:273
void buildReceiveInterface(std::vector< std::pair< TG, int > > &indices)
Definition: repartition.hh:303
void copyOwnerToAll(const T &source, T &dest) const
Communicate values from owner data points to all other data points.
Definition: owneroverlapcopy.hh:315
bool isValidGraph(std::size_t noVtx, std::size_t gnoVtx, S noEdges, T *xadj, T *adjncy, bool checkSymmetry)
Definition: repartition.hh:813
const ParallelIndexSet & indexSet() const
Get the underlying parallel index set.
Definition: owneroverlapcopy.hh:463
Dune::ParallelIndexSet< GlobalIdType, LI, 512 > ParallelIndexSet
The type of the parallel index set.
Definition: owneroverlapcopy.hh:450
const GlobalLookupIndexSet & globalLookup() const
Definition: owneroverlapcopy.hh:527
const RemoteIndices & remoteIndices() const
Get the underlying remote indices.
Definition: owneroverlapcopy.hh:472
void freeGlobalLookup()
Definition: owneroverlapcopy.hh:521
Dune::GlobalLookupIndexSet< ParallelIndexSet > GlobalLookupIndexSet
The type of the reverse lookup of indices.
Definition: owneroverlapcopy.hh:457
Matrix & mat
Definition: matrixmatrix.hh:345
The (undirected) graph of a matrix.
Definition: graph.hh:48
Attaches properties to the edges and vertices of a graph.
Definition: graph.hh:975
bool buildCommunication(const G &graph, std::vector< int > &realparts, Dune::OwnerOverlapCopyCommunication< T1, T2 > &oocomm, Dune::OwnerOverlapCopyCommunication< T1, T2 > *&outcomm, RedistributeInterface &redistInf, bool verbose=false)
Definition: repartition.hh:1487
bool graphRepartition(const G &graph, Dune::OwnerOverlapCopyCommunication< T1, T2 > &oocomm, Metis::idx_t nparts, Dune::OwnerOverlapCopyCommunication< T1, T2 > *&outcomm, RedistributeInterface &redistInf, bool verbose=false)
execute a graph repartition for a giving graph and indexset.
Definition: repartition.hh:1265
derive error class from the base class in common
Definition: istlexception.hh:16
void reserveSpaceForReceiveInterface(int proc, int size)
Definition: repartition.hh:294
void buildGlobalLookup()
Definition: owneroverlapcopy.hh:496
void addReceiveIndex(int proc, std::size_t idx)
Definition: repartition.hh:298
float real_t
Definition: repartition.hh:51
bool commGraphRepartition(const M &mat, Dune::OwnerOverlapCopyCommunication< T1, T2 > &oocomm, Metis::idx_t nparts, Dune::OwnerOverlapCopyCommunication< T1, T2 > *&outcomm, RedistributeInterface &redistInf, bool verbose=false)
Definition: repartition.hh:856
Definition: repartition.hh:265
void copyCopyToAll(const T &source, T &dest) const
Communicate values from copy data points to all other data points.
Definition: owneroverlapcopy.hh:332
A class setting up standard communication for a two-valued attribute set with owner/overlap/copy sema...
Definition: owneroverlapcopy.hh:171
std::size_t idx_t
Definition: repartition.hh:61
Classes providing communication interfaces for overlapping Schwarz methods.
Definition: owneroverlapcopy.hh:59
void setCommunicator(MPI_Comm comm)
Definition: repartition.hh:268
SolverCategory::Category getSolverCategory() const
Get Solver Category.
Definition: owneroverlapcopy.hh:299
const CollectiveCommunication< MPI_Comm > & communicator() const
Definition: owneroverlapcopy.hh:303
Dune::RemoteIndices< PIS > RemoteIndices
The type of the remote indices.
Definition: owneroverlapcopy.hh:453
~RedistributeInterface()
Definition: repartition.hh:312
int globalOwnerVertices
Definition: repartition.hh:179
Provides classes for building the matrix graph.