3 #ifndef DUNE_REPARTITION_HH
4 #define DUNE_REPARTITION_HH
14 #include <dune/common/timer.hh>
15 #include <dune/common/unused.hh>
16 #include <dune/common/enumset.hh>
17 #include <dune/common/stdstreams.hh>
18 #include <dune/common/parallel/mpitraits.hh>
19 #include <dune/common/parallel/communicator.hh>
20 #include <dune/common/parallel/indexset.hh>
21 #include <dune/common/parallel/indicessyncer.hh>
22 #include <dune/common/parallel/remoteindices.hh>
51 template<
class G,
class T1,
class T2>
55 typedef typename IndexSet::LocalIndex::Attribute Attribute;
57 IndexSet& indexSet = oocomm.
indexSet();
61 typedef typename G::ConstVertexIterator VertexIterator;
64 std::size_t sum=0, needed = graph.noVertices()-indexSet.size();
65 std::vector<std::size_t> neededall(oocomm.
communicator().size(), 0);
67 MPI_Allgather(&needed, 1, MPITraits<std::size_t>::getType() , &(neededall[0]), 1, MPITraits<std::size_t>::getType(), oocomm.
communicator());
77 typedef typename IndexSet::const_iterator Iterator;
80 for(Iterator it = indexSet.begin(); it != end; ++it)
81 maxgi=std::max(maxgi,it->global());
90 maxgi=maxgi+neededall[i];
93 std::map<int,SLList<std::pair<T1,Attribute> > > globalIndices;
94 storeGlobalIndicesOfRemoteIndices(globalIndices, oocomm.
remoteIndices());
95 indexSet.beginResize();
97 for(VertexIterator vertex = graph.begin(), vend=graph.end(); vertex != vend; ++vertex) {
98 const typename IndexSet::IndexPair* pair=lookup.pair(*vertex);
106 indexSet.endResize();
108 repairLocalIndexPointers(globalIndices, oocomm.
remoteIndices(), indexSet);
113 std::cout<<
"Holes are filled!"<<std::endl;
121 class ParmetisDuneIndexMap
124 template<
class Graph,
class OOComm>
125 ParmetisDuneIndexMap(
const Graph& graph,
const OOComm& com);
126 int toParmetis(
int i)
const
130 int toLocalParmetis(
int i)
const
134 int operator[](
int i)
const
138 int toDune(
int i)
const
142 std::vector<int>::size_type numOfOwnVtx()
const
159 template<
class G,
class OOComm>
160 ParmetisDuneIndexMap::ParmetisDuneIndexMap(
const G& graph,
const OOComm& oocomm)
163 int npes=oocomm.communicator().size(), mype=oocomm.communicator().rank();
165 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
166 typedef typename OOComm::OwnerSet OwnerSet;
169 Iterator end = oocomm.indexSet().end();
170 for(Iterator
index = oocomm.indexSet().begin();
index != end; ++
index) {
171 if (OwnerSet::contains(
index->local().attribute())) {
176 std::vector<int> globalNumOfVtx(npes);
178 MPI_Allgather(&numOfOwnVtx, 1, MPI_INT, &(globalNumOfVtx[0]), 1, MPI_INT, oocomm.communicator());
182 for(
int i=0; i<npes; i++) {
184 base += globalNumOfVtx[i];
192 typedef typename G::ConstVertexIterator VertexIterator;
194 std::cout << oocomm.communicator().rank()<<
" vtxDist: ";
195 for(
int i=0; i<= npes; ++i)
197 std::cout<<std::endl;
204 VertexIterator vend = graph.end();
205 for(VertexIterator vertex = graph.begin(); vertex != vend; ++vertex) {
206 const typename OOComm::ParallelIndexSet::IndexPair*
index=oocomm.globalLookup().pair(*vertex);
208 if (OwnerSet::contains(index->local().attribute())) {
221 std::cout <<oocomm.communicator().rank()<<
": before ";
224 std::cout<<std::endl;
228 std::cout <<oocomm.communicator().rank()<<
": after ";
231 std::cout<<std::endl;
239 void setCommunicator(MPI_Comm comm)
243 template<
class Flags,
class IS>
244 void buildSendInterface(
const std::vector<int>& toPart,
const IS& idxset)
246 std::map<int,int> sizes;
248 typedef typename IS::const_iterator IIter;
249 for(IIter i=idxset.begin(), end=idxset.end(); i!=end; ++i)
250 if(Flags::contains(i->local().attribute()))
251 ++sizes[toPart[i->local()]];
254 typedef std::map<int,int>::const_iterator MIter;
255 for(MIter i=sizes.begin(), end=sizes.end(); i!=end; ++i)
256 interfaces()[i->first].first.reserve(i->second);
259 typedef typename IS::const_iterator IIter;
260 for(IIter i=idxset.begin(), end=idxset.end(); i!=end; ++i)
261 if(Flags::contains(i->local().attribute()))
262 interfaces()[toPart[i->local()]].first.add(i->local());
265 void reserveSpaceForReceiveInterface(
int proc,
int size)
267 interfaces()[proc].second.reserve(size);
269 void addReceiveIndex(
int proc, std::size_t idx)
271 interfaces()[proc].second.add(idx);
273 template<
typename TG>
274 void buildReceiveInterface(std::vector<std::pair<TG,int> >& indices)
276 typedef typename std::vector<std::pair<TG,int> >::const_iterator VIter;
278 for(VIter idx=indices.begin(); idx!= indices.end(); ++idx) {
279 interfaces()[idx->second].second.add(i++);
300 void createSendBuf(std::vector<GI>& ownerVec, std::set<GI>& overlapVec, std::set<int>& neighbors,
char *sendBuf,
int buffersize, MPI_Comm comm) {
302 std::size_t s=ownerVec.size();
306 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
307 MPI_Pack(&(ownerVec[0]), s, MPITraits<GI>::getType(), sendBuf, buffersize, &pos, comm);
308 s = overlapVec.size();
309 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
310 typedef typename std::set<GI>::iterator Iter;
311 for(Iter i=overlapVec.begin(), end= overlapVec.end(); i != end; ++i)
312 MPI_Pack(const_cast<GI*>(&(*i)), 1, MPITraits<GI>::getType(), sendBuf, buffersize, &pos, comm);
315 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
316 typedef typename std::set<int>::iterator IIter;
318 for(IIter i=neighbors.begin(), end= neighbors.end(); i != end; ++i)
319 MPI_Pack(const_cast<int*>(&(*i)), 1, MPI_INT, sendBuf, buffersize, &pos, comm);
330 void saveRecvBuf(
char *recvBuf,
int bufferSize, std::vector<std::pair<GI,int> >& ownerVec,
331 std::set<GI>& overlapVec, std::set<int>& neighbors, RedistributeInterface& inf,
int from, MPI_Comm comm) {
335 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
336 inf.reserveSpaceForReceiveInterface(from, size);
337 ownerVec.reserve(ownerVec.size()+size);
338 for(; size!=0; --size) {
340 MPI_Unpack(recvBuf, bufferSize, &pos, &gi, 1, MPITraits<GI>::getType(), comm);
341 ownerVec.push_back(std::make_pair(gi,from));
344 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
345 typename std::set<GI>::iterator ipos = overlapVec.begin();
346 Dune::dverb <<
"unpacking "<<size<<
" overlap"<<std::endl;
347 for(; size!=0; --size) {
349 MPI_Unpack(recvBuf, bufferSize, &pos, &gi, 1, MPITraits<GI>::getType(), comm);
350 ipos=overlapVec.insert(ipos, gi);
353 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
354 Dune::dverb <<
"unpacking "<<size<<
" neighbors"<<std::endl;
355 typename std::set<int>::iterator npos = neighbors.begin();
356 for(; size!=0; --size) {
358 MPI_Unpack(recvBuf, bufferSize, &pos, &n, 1, MPI_INT, comm);
359 npos=neighbors.insert(npos, n);
377 void getDomain(
const MPI_Comm& comm, T *part,
int numOfOwnVtx,
int nparts,
int *myDomain, std::vector<int> &domainMapping) {
379 MPI_Comm_size(comm, &npes);
380 MPI_Comm_rank(comm, &mype);
387 std::vector<int> domain(nparts);
388 std::vector<int> assigned(npes);
390 for (i=0; i<nparts; i++) {
391 domainMapping[i] = -1;
394 for (i=0; i<npes; i++) {
398 for (i=0; i<numOfOwnVtx; i++) {
402 int *domainMatrix =
new int[npes * nparts];
404 for(i=0; i<npes*nparts; i++) {
409 int *buf =
new int[nparts];
410 for (i=0; i<nparts; i++) {
412 domainMatrix[mype*nparts+i] = domain[i];
415 int src = (mype-1+npes)%npes;
416 int dest = (mype+1)%npes;
418 for (i=0; i<npes-1; i++) {
419 MPI_Sendrecv_replace(buf, nparts, MPI_INT, dest, 0, src, 0, comm, &status);
421 pe = ((mype-1-i)+npes)%npes;
422 for(j=0; j<nparts; j++) {
424 domainMatrix[pe*nparts+j] = buf[j];
432 int maxOccurance = 0;
434 for(i=0; i<nparts; i++) {
435 for(j=0; j<npes; j++) {
437 if (assigned[j]==0) {
438 if (maxOccurance < domainMatrix[j*nparts+i]) {
439 maxOccurance = domainMatrix[j*nparts+i];
447 domainMapping[i] = pe;
458 delete[] domainMatrix;
465 bool operator()(
const T& t1,
const T& t2)
const
483 void mergeVec(std::vector<std::pair<GI, int> >& ownerVec, std::set<GI>& overlapSet) {
485 typedef typename std::vector<std::pair<GI,int> >::const_iterator VIter;
488 if(ownerVec.size()>0)
490 VIter old=ownerVec.begin();
491 for(VIter i=old+1, end=ownerVec.end(); i != end; old=i++)
493 if(i->first==old->first)
495 std::cerr<<
"Value at indes"<<old-ownerVec.begin()<<
" is the same as at index "
496 <<i-ownerVec.begin()<<
" ["<<old->first<<
","<<old->second<<
"]==["
497 <<i->first<<
","<<i->second<<
"]"<<std::endl;
505 typedef typename std::set<GI>::iterator SIter;
506 VIter v=ownerVec.begin(), vend=ownerVec.end();
507 for(SIter s=overlapSet.begin(), send=overlapSet.end(); s!=send;)
509 while(v!=vend && v->first<*s) ++v;
510 if(v!=vend && v->first==*s) {
515 overlapSet.erase(tmp);
535 template<
class OwnerSet,
class Graph,
class IS,
class GI>
536 void getNeighbor(
const Graph& g, std::vector<int>& part,
537 typename Graph::VertexDescriptor vtx,
const IS& indexSet,
538 int toPe, std::set<GI>& neighbor, std::set<int>& neighborProcs) {
539 typedef typename Graph::ConstEdgeIterator Iter;
540 for(Iter edge=g.beginEdges(vtx), end=g.endEdges(vtx); edge!=end; ++edge)
542 const typename IS::IndexPair* pindex = indexSet.pair(edge.target());
544 if(part[pindex->local()]!=toPe || !OwnerSet::contains(pindex->local().attribute()))
547 neighbor.insert(pindex->global());
548 neighborProcs.insert(part[pindex->local()]);
553 template<
class T,
class I>
554 void my_push_back(std::vector<T>& ownerVec,
const I& index,
int proc)
556 DUNE_UNUSED_PARAMETER(proc);
557 ownerVec.push_back(index);
560 template<
class T,
class I>
561 void my_push_back(std::vector<std::pair<T,int> >& ownerVec,
const I& index,
int proc)
563 ownerVec.push_back(std::make_pair(index,proc));
566 void reserve(std::vector<T>&, RedistributeInterface&,
int)
569 void reserve(std::vector<std::pair<T,int> >& ownerVec, RedistributeInterface& redist,
int proc)
571 redist.reserveSpaceForReceiveInterface(proc, ownerVec.size());
592 template<
class OwnerSet,
class G,
class IS,
class T,
class GI>
593 void getOwnerOverlapVec(
const G& graph, std::vector<int>& part, IS& indexSet,
594 int myPe,
int toPe, std::vector<T>& ownerVec, std::set<GI>& overlapSet,
595 RedistributeInterface& redist, std::set<int>& neighborProcs) {
596 DUNE_UNUSED_PARAMETER(myPe);
598 typedef typename IS::const_iterator Iterator;
599 for(Iterator index = indexSet.begin(); index != indexSet.end(); ++
index) {
601 if(OwnerSet::contains(index->local().attribute()))
603 if(part[index->local()]==toPe)
605 getNeighbor<OwnerSet>(graph, part, index->local(), indexSet,
606 toPe, overlapSet, neighborProcs);
607 my_push_back(ownerVec, index->global(), toPe);
611 reserve(ownerVec, redist, toPe);
622 template<
class F,
class IS>
623 inline bool isOwner(IS& indexSet,
int index) {
625 const typename IS::IndexPair* pindex=indexSet.pair(index);
628 return F::contains(pindex->local().attribute());
632 class BaseEdgeFunctor
635 BaseEdgeFunctor(
int* adj,
const ParmetisDuneIndexMap& data)
640 void operator()(
const T& edge)
660 :
public BaseEdgeFunctor
662 EdgeFunctor(
int* adj,
const ParmetisDuneIndexMap& data, std::size_t s)
663 : BaseEdgeFunctor(adj, data)
673 template<
class G,
class V,
class E,
class VM,
class EM>
674 class EdgeFunctor<Dune::Amg::PropertiesGraph<G,V,E,VM,EM> >
675 :
public BaseEdgeFunctor
678 EdgeFunctor(
int* adj,
const ParmetisDuneIndexMap& data, std::size_t s)
679 : BaseEdgeFunctor(adj, data)
685 void operator()(
const T& edge)
687 weight_[
index()]=edge.properties().depends() ? 3 : 1;
688 BaseEdgeFunctor::operator()(edge);
719 template<
class F,
class G,
class IS,
class EW>
720 void getAdjArrays(G& graph, IS& indexSet,
int *xadj,
726 typedef typename G::ConstVertexIterator VertexIterator;
728 typedef typename IS::const_iterator Iterator;
730 VertexIterator vend = graph.end();
733 for(VertexIterator vertex = graph.begin(); vertex != vend; ++vertex) {
734 if (isOwner<F>(indexSet,*vertex)) {
736 typedef typename G::ConstEdgeIterator EdgeIterator;
737 EdgeIterator eend = vertex.end();
738 xadj[j] = ew.index();
740 for(EdgeIterator edge = vertex.begin(); edge != eend; ++edge) {
745 xadj[j] = ew.index();
749 template<
class G,
class T1,
class T2>
753 RedistributeInterface& redistInf,
758 #if PARMETIS_MAJOR_VERSION > 3
761 void METIS_PartGraphKway(
int *nvtxs, idxtype *xadj, idxtype *adjncy, idxtype *vwgt,
762 idxtype *adjwgt,
int *wgtflag,
int *numflag,
int *nparts,
763 int *options,
int *edgecut, idxtype *part);
765 void METIS_PartGraphRecursive(
int *nvtxs, idxtype *xadj, idxtype *adjncy, idxtype *vwgt,
766 idxtype *adjwgt,
int *wgtflag,
int *numflag,
int *nparts,
767 int *options,
int *edgecut, idxtype *part);
773 template<
class S,
class T>
776 for(T *cur=array, *end=array+l; cur!=end; ++cur)
781 idxtype* adjncy,
bool checkSymmetry)
786 if(xadj[vtx]>noEdges||xadj[vtx]<0) {
787 std::cerr <<
"Check graph: xadj["<<vtx<<
"]="<<xadj[vtx]<<
" (>"
788 <<noEdges<<
") out of range!"<<std::endl;
791 if(xadj[vtx+1]>noEdges||xadj[vtx+1]<0) {
792 std::cerr <<
"Check graph: xadj["<<vtx+1<<
"]="<<xadj[vtx+1]<<
" (>"
793 <<noEdges<<
") out of range!"<<std::endl;
797 for(
idxtype i=xadj[vtx]; i< xadj[vtx+1]; ++i) {
798 if(adjncy[i]<0||((std::size_t)adjncy[i])>gnoVtx) {
799 std::cerr<<
" Edge "<<adjncy[i]<<
" out of range ["<<0<<
","<<noVtx<<
")"
805 for(
idxtype i=xadj[vtx]; i< xadj[vtx+1]; ++i) {
809 for(
idxtype j=xadj[target]; j< xadj[target+1]; ++j)
813 std::cerr<<
"Edge ("<<target<<
","<<vtx<<
") "<<i<<
" time"<<std::endl;
822 template<
class M,
class T1,
class T2>
829 std::cout<<
"Repartitioning from "<<oocomm.
communicator().size()
830 <<
" to "<<nparts<<
" parts"<<std::endl;
834 int* part =
new int[1];
849 typedef typename RemoteIndices::const_iterator
878 xadj[1]=noNeighbours;
908 typedef typename RemoteIndices::const_iterator NeighbourIterator;
916 adjwgt =
new idxtype[noNeighbours];
922 if(n->first != rank) {
932 noNeighbours, xadj, adjncy,
false));
934 int wgtflag=0, numflag=0, edgecut;
938 float *tpwgts =
new float[nparts];
939 for(
int i=0; i<nparts; ++i)
940 tpwgts[i]=1.0/nparts;
941 int options[5] ={ 0,1,15,0,0};
944 Dune::dinfo<<rank<<
" vtxdist: ";
946 Dune::dinfo<<std::endl<<rank<<
" xadj: ";
948 Dune::dinfo<<std::endl<<rank<<
" adjncy: ";
952 Dune::dinfo<<std::endl<<rank<<
" vwgt: ";
954 Dune::dinfo<<std::endl<<rank<<
" adwgt: ";
957 Dune::dinfo<<std::endl;
960 std::cout<<
"Creating comm graph took "<<time.elapsed()<<std::endl;
963 #ifdef PARALLEL_PARTITION
970 ParMETIS_V3_PartKway(vtxdist, xadj, adjncy,
971 vwgt, adjwgt, &wgtflag,
972 &numflag, &ncon, &nparts, tpwgts, &ubvec, options, &edgecut, part,
975 std::cout<<
"ParMETIS took "<<time.elapsed()<<std::endl;
979 std::size_t gnoedges=0;
982 Dune::dverb<<
"noNeighbours: "<<noNeighbours<<std::endl;
984 MPI_Allgather(&noNeighbours,1,MPI_INT,noedges,1, MPI_INT,oocomm.
communicator());
987 std::cout<<
"Gathering noedges took "<<time1.elapsed()<<std::endl;
1002 std::size_t localNoVtx=vtxdist[rank+1]-vtxdist[rank];
1007 Dune::dinfo<<
"noedges: ";
1009 Dune::dinfo<<std::endl;
1017 noxs[i]=vtxdist[i+1]-vtxdist[i]+1;
1018 novs[i]=vtxdist[i+1]-vtxdist[i];
1023 for(
int *xcurr = xdispl, *vcurr = vdispl, *end=vdispl+oocomm.
communicator().size();
1024 vcurr!=end; ++vcurr, ++xcurr, ++so, ++offset) {
1026 *xcurr = offset + *so;
1032 for(
int *curr=noedges, *end=noedges+oocomm.
communicator().size()-1;
1033 curr!=end; ++curr) {
1038 Dune::dinfo<<
"displ: ";
1040 Dune::dinfo<<std::endl;
1044 for(
int *curr=noedges, *end=noedges+oocomm.
communicator().size();
1049 Dune::dinfo<<
"gxadjlen: "<<gxadjlen<<
" noVertices: "<<noVertices
1050 <<
" gnoedges: "<<gnoedges<<std::endl;
1051 gxadj =
new idxtype[gxadjlen];
1052 gpart =
new idxtype[noVertices];
1054 gvwgt =
new idxtype[noVertices];
1055 gadjwgt =
new idxtype[gnoedges];
1057 gadjncy =
new idxtype[gnoedges];
1061 std::cout<<
"Preparing global graph took "<<time1.elapsed()<<std::endl;
1065 MPI_Allgatherv(xadj,2,MPITraits<idxtype>::getType(),
1066 gxadj,noxs,xdispl,MPITraits<idxtype>::getType(),
1068 MPI_Allgatherv(adjncy,noNeighbours,MPITraits<idxtype>::getType(),
1069 gadjncy,noedges,displ,MPITraits<idxtype>::getType(),
1072 MPI_Allgatherv(adjwgt,noNeighbours,MPITraits<idxtype>::getType(),
1073 gadjwgt,noedges,displ,MPITraits<idxtype>::getType(),
1075 MPI_Allgatherv(vwgt,localNoVtx,MPITraits<idxtype>::getType(),
1076 gvwgt,novs,vdispl,MPITraits<idxtype>::getType(),
1080 std::cout<<
"Gathering global graph data took "<<time1.elapsed()<<std::endl;
1090 idxtype increment = vtxdist[1];
1094 int lprev = vtxdist[i]-vtxdist[i-1];
1095 int l = vtxdist[i+1]-vtxdist[i];
1097 assert((start+l+offset)-gxadj<=static_cast<idxtype>(gxadjlen));
1098 increment = *(start-1);
1099 std::transform(start+offset, start+l+offset, start, std::bind2nd(std::plus<idxtype>(), increment));
1101 Dune::dinfo<<std::endl<<
"shifted xadj:";
1103 Dune::dinfo<<std::endl<<
" gadjncy: ";
1106 Dune::dinfo<<std::endl<<
" gvwgt: ";
1108 Dune::dinfo<<std::endl<<
"adjwgt: ";
1110 Dune::dinfo<<std::endl;
1114 std::cout<<
"Postprocesing global graph data took "<<time1.elapsed()<<std::endl;
1118 gxadj, gadjncy,
true));
1122 std::cout<<
"Creating grah one 1 process took "<<time.elapsed()<<std::endl;
1124 options[0]=0; options[1]=1; options[2]=1; options[3]=3; options[4]=3;
1126 METIS_PartGraphRecursive(&noVertices, gxadj, gadjncy, gvwgt, gadjwgt, &wgtflag,
1127 &numflag, &nparts, options, &edgecut, gpart);
1130 std::cout<<
"METIS took "<<time.elapsed()<<std::endl;
1133 Dune::dinfo<<std::endl<<
"part:";
1144 MPI_Scatter(gpart, 1, MPITraits<idxtype>::getType(), part, 1,
1145 MPITraits<idxtype>::getType(), 0, comm);
1169 Dune::dinfo<<
" repart "<<rank <<
" -> "<< part[0]<<std::endl;
1171 std::vector<int> realpart(mat.N(), part[0]);
1177 std::cout<<
"Scattering repartitioning took "<<time.elapsed()<<std::endl;
1185 std::cout<<
"Filling index set took "<<time.elapsed()<<std::endl;
1193 std::cout<<
"Average no neighbours was "<<noNeighbours<<std::endl;
1198 std::cout<<
"Building index sets took "<<time.elapsed()<<std::endl;
1220 template<
class G,
class T1,
class T2>
1233 std::cout<<
"Filling holes took "<<time.elapsed()<<std::endl;
1240 double t1=0.0, t2=0.0, t3=0.0, t4=0.0, tSum=0.0;
1268 typedef typename OOComm::OwnerSet OwnerSet;
1273 ParmetisDuneIndexMap indexMap(graph,oocomm);
1277 std::size_t *part =
new std::size_t[indexMap.numOfOwnVtx()];
1279 for(std::size_t i=0; i < indexMap.numOfOwnVtx(); ++i)
1284 std::cerr<<
"ParMETIS not activated. Will repartition to 1 domain instead of requested "
1285 <<nparts<<
" domains."<<std::endl;
1294 EdgeFunctor<G> ef(adjncy, indexMap, graph.noEdges());
1295 getAdjArrays<OwnerSet>(graph, oocomm.
globalLookup(), xadj, ef);
1301 int numflag=0, wgtflag=0, options[3], edgecut=0, ncon=1;
1303 float *tpwgts =
new float[nparts];
1304 for(
int i=0; i<nparts; ++i)
1305 tpwgts[i]=1.0/nparts;
1314 wgtflag = (ef.getWeights()!=NULL) ? 1 : 0;
1322 std::cout<<std::endl;
1323 std::cout<<
"Testing ParMETIS_V3_PartKway with options[1-2] = {"
1324 <<options[1]<<
" "<<options[2]<<
"}, Ncon: "
1325 <<ncon<<
", Nparts: "<<nparts<<std::endl;
1338 std::cout<<
"Preparing for parmetis took "<<time.elapsed()<<std::endl;
1345 ParMETIS_V3_PartKway(indexMap.vtxDist(), xadj, adjncy,
1346 NULL, ef.getWeights(), &wgtflag,
1347 &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgecut, part, &
const_cast<MPI_Comm&
>(comm));
1358 std::cout<<std::endl;
1359 std::cout<<
"ParMETIS_V3_PartKway reported a cut of "<<edgecut<<std::endl;
1360 std::cout<<std::endl;
1362 std::cout<<mype<<
": PARMETIS-Result: ";
1363 for(
int i=0; i < indexMap.vtxDist()[mype+1]-indexMap.vtxDist()[mype]; ++i) {
1364 std::cout<<part[i]<<
" ";
1366 std::cout<<std::endl;
1367 std::cout<<
"Testing ParMETIS_V3_PartKway with options[1-2] = {"
1368 <<options[1]<<
" "<<options[2]<<
"}, Ncon: "
1369 <<ncon<<
", Nparts: "<<nparts<<std::endl;
1382 std::cout<<
"Parmetis took "<<time.elapsed()<<std::endl;
1389 for(std::size_t i=0; i<indexMap.numOfOwnVtx(); ++i)
1399 std::vector<int> domainMapping(nparts);
1401 getDomain(comm, part, indexMap.numOfOwnVtx(), nparts, &myDomain, domainMapping);
1406 std::cout<<mype<<
": myDomain: "<<myDomain<<std::endl;
1407 std::cout<<mype<<
": DomainMapping: ";
1408 for(
int j=0; j<nparts; j++) {
1409 std::cout<<
" do: "<<j<<
" pe: "<<domainMapping[j]<<
" ";
1411 std::cout<<std::endl;
1418 std::vector<int> setPartition(oocomm.
indexSet().size(), -1);
1420 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
1423 if(OwnerSet::contains(index->local().attribute())) {
1424 setPartition[index->local()]=domainMapping[part[i++]];
1432 static_cast<int>(SolverCategory::nonoverlapping))
1439 std::cout<<
"Creating indexsets took "<<time.elapsed()<<std::endl;
1446 template<
class G,
class T1,
class T2>
1454 typedef typename OOComm::OwnerSet OwnerSet;
1488 std::set<int> recvFrom;
1494 typedef typename std::vector<int>::const_iterator VIter;
1498 std::set<int> tsendTo;
1499 for(VIter i=setPartition.begin(), iend = setPartition.end(); i!=iend; ++i)
1502 noSendTo = tsendTo.size();
1503 sendTo =
new int[noSendTo];
1504 typedef std::set<int>::const_iterator iterator;
1506 for(iterator i=tsendTo.begin(); i != tsendTo.end(); ++i, ++idx)
1512 int* gsendToDispl =
new int[oocomm.
communicator().size()+1];
1514 MPI_Allgather(&noSendTo, 1, MPI_INT, gnoSend, 1,
1518 int totalNoRecv = 0;
1519 for(
int i=0; i<npes; ++i)
1520 totalNoRecv += gnoSend[i];
1522 int *gsendTo =
new int[totalNoRecv];
1526 for(
int i=0; i<npes; ++i)
1527 gsendToDispl[i+1]=gsendToDispl[i]+gnoSend[i];
1530 MPI_Allgatherv(sendTo, noSendTo, MPI_INT, gsendTo, gnoSend, gsendToDispl,
1534 for(
int proc=0; proc < npes; ++proc)
1535 for(
int i=gsendToDispl[proc]; i < gsendToDispl[proc+1]; ++i)
1536 if(gsendTo[i]==mype)
1537 recvFrom.insert(proc);
1539 bool existentOnNextLevel = recvFrom.size()>0;
1543 delete[] gsendToDispl;
1548 if(recvFrom.size()) {
1549 std::cout<<mype<<
": recvFrom: ";
1550 typedef typename std::set<int>::const_iterator siter;
1551 for(siter i=recvFrom.begin(); i!= recvFrom.end(); ++i) {
1556 std::cout<<std::endl<<std::endl;
1557 std::cout<<mype<<
": sendTo: ";
1558 for(
int i=0; i<noSendTo; i++) {
1559 std::cout<<sendTo[i]<<
" ";
1561 std::cout<<std::endl<<std::endl;
1566 std::cout<<
" Communicating the receive information took "<<
1567 time.elapsed()<<std::endl;
1581 typedef typename OOComm::ParallelIndexSet::GlobalIndex GI;
1582 typedef std::vector<GI> GlobalVector;
1583 std::vector<std::pair<GI,int> > myOwnerVec;
1584 std::set<GI> myOverlapSet;
1585 GlobalVector sendOwnerVec;
1586 std::set<GI> sendOverlapSet;
1587 std::set<int> myNeighbors;
1592 char **sendBuffers=
new char*[noSendTo];
1593 MPI_Request *requests =
new MPI_Request[noSendTo];
1596 for(
int i=0; i < noSendTo; ++i) {
1598 sendOwnerVec.clear();
1599 sendOverlapSet.clear();
1602 std::set<int> neighbors;
1603 getOwnerOverlapVec<OwnerSet>(graph, setPartition, oocomm.
globalLookup(),
1604 mype, sendTo[i], sendOwnerVec, sendOverlapSet, redistInf,
1610 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.
communicator(), &buffersize);
1611 MPI_Pack_size(sendOwnerVec.size(), MPITraits<GI>::getType(), oocomm.
communicator(), &tsize);
1613 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.
communicator(), &tsize);
1615 MPI_Pack_size(sendOverlapSet.size(), MPITraits<GI>::getType(), oocomm.
communicator(), &tsize);
1616 buffersize += tsize;
1617 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.
communicator(), &tsize);
1618 buffersize += tsize;
1619 MPI_Pack_size(neighbors.size(), MPI_INT, oocomm.
communicator(), &tsize);
1620 buffersize += tsize;
1622 sendBuffers[i] =
new char[buffersize];
1625 std::cout<<mype<<
" sending "<<sendOwnerVec.size()<<
" owner and "<<
1626 sendOverlapSet.size()<<
" overlap to "<<sendTo[i]<<
" buffersize="<<buffersize<<std::endl;
1628 createSendBuf(sendOwnerVec, sendOverlapSet, neighbors, sendBuffers[i], buffersize, oocomm.
communicator());
1629 MPI_Issend(sendBuffers[i], buffersize, MPI_PACKED, sendTo[i], 99, oocomm.
communicator(), requests+i);
1635 std::cout<<
" Creating sends took "<<
1636 time.elapsed()<<std::endl;
1641 int noRecv = recvFrom.size();
1642 int oldbuffersize=0;
1647 MPI_Probe(MPI_ANY_SOURCE, 99, oocomm.
communicator(), &stat);
1649 MPI_Get_count(&stat, MPI_PACKED, &buffersize);
1651 if(oldbuffersize<buffersize) {
1654 recvBuf =
new char[buffersize];
1655 oldbuffersize = buffersize;
1657 MPI_Recv(recvBuf, buffersize, MPI_PACKED, stat.MPI_SOURCE, 99, oocomm.
communicator(), &stat);
1658 saveRecvBuf(recvBuf, buffersize, myOwnerVec, myOverlapSet, myNeighbors, redistInf,
1668 MPI_Status *statuses =
new MPI_Status[noSendTo];
1669 int send = MPI_Waitall(noSendTo, requests, statuses);
1672 if(send==MPI_ERR_IN_STATUS) {
1673 std::cerr<<mype<<
": Error in sending :"<<std::endl;
1675 for(
int i=0; i< noSendTo; i++)
1676 if(statuses[i].MPI_ERROR!=MPI_SUCCESS) {
1679 MPI_Error_string(statuses[i].MPI_ERROR, message, &messageLength);
1680 std::cerr<<
" source="<<statuses[i].MPI_SOURCE<<
" message: ";
1681 for(
int i=0; i< messageLength; i++)
1682 std::cout<<message[i];
1684 std::cerr<<std::endl;
1690 std::cout<<
" Receiving and saving took "<<
1691 time.elapsed()<<std::endl;
1695 for(
int i=0; i < noSendTo; ++i)
1696 delete[] sendBuffers[i];
1698 delete[] sendBuffers;
1713 if (!existentOnNextLevel) {
1715 color= MPI_UNDEFINED;
1717 MPI_Comm outputComm;
1725 std::vector<int> tneighbors;
1726 tneighbors.reserve(myNeighbors.size());
1728 typename OOComm::ParallelIndexSet& outputIndexSet = outcomm->indexSet();
1730 MPI_Allgather(&newrank, 1, MPI_INT, newranks, 1,
1732 typedef typename std::set<int>::const_iterator IIter;
1736 for(IIter i=myNeighbors.begin(), end=myNeighbors.end();
1738 assert(newranks[*i]>=0);
1739 std::cout<<*i<<
"->"<<newranks[*i]<<
" ";
1740 tneighbors.push_back(newranks[*i]);
1742 std::cout<<std::endl;
1744 for(IIter i=myNeighbors.begin(), end=myNeighbors.end();
1746 tneighbors.push_back(newranks[*i]);
1750 myNeighbors.clear();
1755 std::cout<<
" Calculating new neighbours ("<<tneighbors.size()<<
") took "<<
1756 time.elapsed()<<std::endl;
1761 outputIndexSet.beginResize();
1764 std::sort(myOwnerVec.begin(), myOwnerVec.end(), SortFirst());
1768 typedef typename OOComm::ParallelIndexSet::LocalIndex LocalIndex;
1769 typedef typename std::vector<std::pair<GI,int> >::const_iterator VPIter;
1771 for(VPIter g=myOwnerVec.begin(), end =myOwnerVec.end(); g!=end; ++g, ++i ) {
1772 outputIndexSet.add(g->first,LocalIndex(i, OwnerOverlapCopyAttributeSet::owner,
true));
1779 std::cout<<
" Adding owner indices took "<<
1780 time.elapsed()<<std::endl;
1789 mergeVec(myOwnerVec, myOverlapSet);
1793 myOwnerVec.swap(myOwnerVec);
1798 std::cout<<
" Merging indices took "<<
1799 time.elapsed()<<std::endl;
1805 typedef typename std::set<GI>::const_iterator SIter;
1806 for(SIter g=myOverlapSet.begin(), end=myOverlapSet.end(); g!=end; ++g, i++) {
1807 outputIndexSet.add(*g,LocalIndex(i, OwnerOverlapCopyAttributeSet::copy,
true));
1809 myOverlapSet.clear();
1810 outputIndexSet.endResize();
1812 #ifdef DUNE_ISTL_WITH_CHECKING
1814 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
1815 Iterator end = outputIndexSet.end();
1816 for(Iterator index = outputIndexSet.begin(); index != end; ++
index) {
1817 if (OwnerSet::contains(index->local().attribute())) {
1828 Iterator index=outputIndexSet.begin();
1831 for(Iterator old = outputIndexSet.begin(); index != end; old=index++) {
1832 if(old->global()>index->global())
1833 DUNE_THROW(
ISTLError,
"Index set's globalindex not sorted correctly");
1840 std::cout<<
" Adding overlap indices took "<<
1841 time.elapsed()<<std::endl;
1846 if(color != MPI_UNDEFINED) {
1847 outcomm->remoteIndices().setNeighbours(tneighbors);
1848 outcomm->remoteIndices().template rebuild<true>();
1858 std::cout<<
" Storing indexsets took "<<
1859 time.elapsed()<<std::endl;
1865 tSum = t1 + t2 + t3 + t4;
1866 std::cout<<std::endl
1867 <<mype<<
": WTime for step 1): "<<t1
1875 return color!=MPI_UNDEFINED;
1879 template<
class G,
class P,
class T1,
class T2,
class R>
1885 if(nparts!=oocomm.size())
1886 DUNE_THROW(NotImplemented,
"only available for MPI programs");
1890 template<
class G,
class P,
class T1,
class T2,
class R>
1896 if(nparts!=oocomm.size())
1897 DUNE_THROW(NotImplemented,
"only available for MPI programs");