1 #ifndef DUNE_REPARTITION_HH
2 #define DUNE_REPARTITION_HH
12 #include <dune/common/timer.hh>
13 #include <dune/common/enumset.hh>
14 #include <dune/common/mpitraits.hh>
15 #include <dune/common/stdstreams.hh>
16 #include <dune/common/parallel/communicator.hh>
17 #include <dune/common/parallel/indexset.hh>
18 #include <dune/common/parallel/indicessyncer.hh>
19 #include <dune/common/parallel/remoteindices.hh>
48 template<
class G,
class T1,
class T2>
52 typedef typename IndexSet::LocalIndex::Attribute Attribute;
54 IndexSet& indexSet = oocomm.
indexSet();
58 typedef typename G::ConstVertexIterator VertexIterator;
61 std::size_t sum=0, needed = graph.noVertices()-indexSet.size();
62 std::vector<std::size_t> neededall(oocomm.
communicator().size(), 0);
64 MPI_Allgather(&needed, 1, MPITraits<std::size_t>::getType() , &(neededall[0]), 1, MPITraits<std::size_t>::getType(), oocomm.
communicator());
74 typedef typename IndexSet::const_iterator Iterator;
77 for(Iterator it = indexSet.begin(); it != end; ++it)
78 maxgi=std::max(maxgi,it->global());
87 maxgi=maxgi+neededall[i];
90 std::map<int,SLList<std::pair<T1,Attribute> > > globalIndices;
91 storeGlobalIndicesOfRemoteIndices(globalIndices, oocomm.
remoteIndices(), indexSet);
92 indexSet.beginResize();
94 for(VertexIterator vertex = graph.begin(), vend=graph.end(); vertex != vend; ++vertex){
95 const typename IndexSet::IndexPair* pair=lookup.pair(*vertex);
103 indexSet.endResize();
105 repairLocalIndexPointers(globalIndices, oocomm.
remoteIndices(), indexSet);
110 std::cout<<
"Holes are filled!"<<std::endl;
118 class ParmetisDuneIndexMap
121 template<
class Graph,
class OOComm>
122 ParmetisDuneIndexMap(
const Graph& graph,
const OOComm& com);
123 int toParmetis(
int i)
const
127 int toLocalParmetis(
int i)
const
131 int operator[](
int i)
const
135 int toDune(
int i)
const
139 std::vector<int>::size_type numOfOwnVtx()
const
156 template<
class G,
class OOComm>
157 ParmetisDuneIndexMap::ParmetisDuneIndexMap(
const G& graph,
const OOComm& oocomm)
160 int npes=oocomm.communicator().size(), mype=oocomm.communicator().rank();
162 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
163 typedef typename OOComm::OwnerSet OwnerSet;
166 Iterator end = oocomm.indexSet().end();
167 for(Iterator
index = oocomm.indexSet().begin();
index != end; ++
index) {
168 if (OwnerSet::contains(
index->local().attribute())) {
173 std::vector<int> globalNumOfVtx(npes);
175 MPI_Allgather(&numOfOwnVtx, 1, MPI_INT, &(globalNumOfVtx[0]), 1, MPI_INT, oocomm.communicator());
179 for(
int i=0; i<npes; i++) {
181 base += globalNumOfVtx[i];
189 typedef typename G::ConstVertexIterator VertexIterator;
191 std::cout << oocomm.communicator().rank()<<
" vtxDist: ";
192 for(
int i=0; i<= npes;++i)
194 std::cout<<std::endl;
201 VertexIterator vend = graph.end();
202 for(VertexIterator vertex = graph.begin(); vertex != vend; ++vertex) {
203 const typename OOComm::ParallelIndexSet::IndexPair*
index=oocomm.globalLookup().pair(*vertex);
205 if (OwnerSet::contains(index->local().attribute())) {
218 std::cout <<oocomm.communicator().rank()<<
": before ";
221 std::cout<<std::endl;
225 std::cout <<oocomm.communicator().rank()<<
": after ";
228 std::cout<<std::endl;
236 void setCommunicator(MPI_Comm comm)
240 template<
class Flags,
class IS>
241 void buildSendInterface(
const std::vector<int>& toPart,
const IS& idxset)
243 std::map<int,int> sizes;
245 typedef typename IS::const_iterator IIter;
246 for(IIter i=idxset.begin(), end=idxset.end();i!=end; ++i)
247 if(Flags::contains(i->local().attribute()))
248 ++sizes[toPart[i->local()]];
251 typedef std::map<int,int>::const_iterator MIter;
252 for(MIter i=sizes.begin(), end=sizes.end(); i!=end; ++i)
253 interfaces()[i->first].first.reserve(i->second);
256 typedef typename IS::const_iterator IIter;
257 for(IIter i=idxset.begin(), end=idxset.end();i!=end; ++i)
258 if(Flags::contains(i->local().attribute()))
259 interfaces()[toPart[i->local()]].first.add(i->local());
262 void reserveSpaceForReceiveInterface(
int proc,
int size)
264 interfaces()[proc].second.reserve(size);
266 void addReceiveIndex(
int proc, std::size_t idx)
268 interfaces()[proc].second.add(idx);
270 template<
typename TG>
271 void buildReceiveInterface(std::vector<std::pair<TG,int> >& indices)
273 typedef typename std::vector<std::pair<TG,int> >::const_iterator VIter;
275 for(VIter idx=indices.begin(); idx!= indices.end(); ++idx){
276 interfaces()[idx->second].second.add(i++);
298 void createSendBuf(std::vector<GI>& ownerVec, std::set<GI>& overlapVec, std::set<int>& neighbors,
char *sendBuf,
int buffersize, MPI_Comm comm) {
300 std::size_t s=ownerVec.size();
304 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
305 MPI_Pack(&(ownerVec[0]), s, MPITraits<GI>::getType(), sendBuf, buffersize, &pos, comm);
306 s = overlapVec.size();
307 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
308 typedef typename std::set<GI>::iterator Iter;
309 for(Iter i=overlapVec.begin(), end= overlapVec.end(); i != end; ++i)
310 MPI_Pack(const_cast<GI*>(&(*i)), 1, MPITraits<GI>::getType(), sendBuf, buffersize, &pos, comm);
313 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
314 typedef typename std::set<int>::iterator IIter;
316 for(IIter i=neighbors.begin(), end= neighbors.end(); i != end; ++i)
317 MPI_Pack(const_cast<int*>(&(*i)), 1, MPI_INT, sendBuf, buffersize, &pos, comm);
328 void saveRecvBuf(
char *recvBuf,
int bufferSize, std::vector<std::pair<GI,int> >& ownerVec,
329 std::set<GI>& overlapVec, std::set<int>& neighbors, RedistributeInterface& inf,
int from, MPI_Comm comm) {
333 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
334 inf.reserveSpaceForReceiveInterface(from, size);
335 ownerVec.reserve(ownerVec.size()+size);
336 for(;size!=0;--size){
338 MPI_Unpack(recvBuf, bufferSize, &pos, &gi, 1, MPITraits<GI>::getType(), comm);
339 ownerVec.push_back(std::make_pair(gi,from));
342 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
343 typename std::set<GI>::iterator ipos = overlapVec.begin();
344 Dune::dverb <<
"unpacking "<<size<<
" overlap"<<std::endl;
345 for(;size!=0;--size){
347 MPI_Unpack(recvBuf, bufferSize, &pos, &gi, 1, MPITraits<GI>::getType(), comm);
348 ipos=overlapVec.insert(ipos, gi);
351 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
352 Dune::dverb <<
"unpacking "<<size<<
" neighbors"<<std::endl;
353 typename std::set<int>::iterator npos = neighbors.begin();
354 for(;size!=0;--size){
356 MPI_Unpack(recvBuf, bufferSize, &pos, &n, 1, MPI_INT, comm);
357 npos=neighbors.insert(npos, n);
375 void getDomain(
const MPI_Comm& comm, T *part,
int numOfOwnVtx,
int nparts,
int *myDomain, std::vector<int> &domainMapping) {
377 MPI_Comm_size(comm, &npes);
378 MPI_Comm_rank(comm, &mype);
385 std::vector<int> domain(nparts);
386 std::vector<int> assigned(npes);
388 for (i=0; i<nparts; i++) {
389 domainMapping[i] = -1;
392 for (i=0; i<npes; i++) {
396 for (i=0; i<numOfOwnVtx; i++) {
400 int *domainMatrix =
new int[npes * nparts];
402 for(i=0; i<npes*nparts; i++) {
407 int *buf =
new int[nparts];
408 for (i=0; i<nparts; i++) {
410 domainMatrix[mype*nparts+i] = domain[i];
413 int src = (mype-1+npes)%npes;
414 int dest = (mype+1)%npes;
416 for (i=0; i<npes-1; i++) {
417 MPI_Sendrecv_replace(buf, nparts, MPI_INT, dest, 0, src, 0, comm, &status);
419 pe = ((mype-1-i)+npes)%npes;
420 for(j=0; j<nparts; j++) {
422 domainMatrix[pe*nparts+j] = buf[j];
430 int maxOccurance = 0;
432 for(i=0; i<nparts; i++) {
433 for(j=0; j<npes; j++) {
435 if (assigned[j]==0) {
436 if (maxOccurance < domainMatrix[j*nparts+i]) {
437 maxOccurance = domainMatrix[j*nparts+i];
445 domainMapping[i] = pe;
456 delete[] domainMatrix;
463 bool operator()(
const T& t1,
const T& t2)
const
481 void mergeVec(std::vector<std::pair<GI, int> >& ownerVec, std::set<GI>& overlapSet) {
483 typedef typename std::vector<std::pair<GI,int> >::const_iterator VIter;
486 if(ownerVec.size()>0)
488 VIter old=ownerVec.begin();
489 for(VIter i=old+1, end=ownerVec.end(); i != end; old=i++)
491 if(i->first==old->first)
493 std::cerr<<
"Value at indes"<<old-ownerVec.begin()<<
" is the same as at index "
494 <<i-ownerVec.begin()<<
" ["<<old->first<<
","<<old->second<<
"]==["
495 <<i->first<<
","<<i->second<<
"]"<<std::endl;
503 typedef typename std::set<GI>::iterator SIter;
504 VIter v=ownerVec.begin(), vend=ownerVec.end();
505 for(SIter s=overlapSet.begin(), send=overlapSet.end(); s!=send;)
507 while(v!=vend && v->first<*s) ++v;
508 if(v!=vend && v->first==*s){
513 overlapSet.erase(tmp);
533 template<
class OwnerSet,
class Graph,
class IS,
class GI>
534 void getNeighbor(
const Graph& g, std::vector<int>& part,
535 typename Graph::VertexDescriptor vtx,
const IS& indexSet,
536 int toPe, std::set<GI>& neighbor, std::set<int>& neighborProcs) {
537 typedef typename Graph::ConstEdgeIterator Iter;
538 for(Iter edge=g.beginEdges(vtx), end=g.endEdges(vtx); edge!=end; ++edge)
540 const typename IS::IndexPair* pindex = indexSet.pair(edge.target());
542 if(part[pindex->local()]!=toPe || !OwnerSet::contains(pindex->local().attribute()))
545 neighbor.insert(pindex->global());
546 neighborProcs.insert(part[pindex->local()]);
551 template<
class T,
class I>
552 void my_push_back(std::vector<T>& ownerVec,
const I& index,
int proc)
554 ownerVec.push_back(index);
557 template<
class T,
class I>
558 void my_push_back(std::vector<std::pair<T,int> >& ownerVec,
const I& index,
int proc)
560 ownerVec.push_back(std::make_pair(index,proc));
563 void reserve(std::vector<T>&, RedistributeInterface&,
int)
567 void reserve(std::vector<std::pair<T,int> >& ownerVec, RedistributeInterface& redist,
int proc)
569 redist.reserveSpaceForReceiveInterface(proc, ownerVec.size());
590 template<
class OwnerSet,
class G,
class IS,
class T,
class GI>
591 void getOwnerOverlapVec(
const G& graph, std::vector<int>& part, IS& indexSet,
592 int myPe,
int toPe, std::vector<T>& ownerVec, std::set<GI>& overlapSet,
593 RedistributeInterface& redist, std::set<int>& neighborProcs) {
596 typedef typename IS::const_iterator Iterator;
597 for(Iterator index = indexSet.begin(); index != indexSet.end(); ++
index) {
599 if(OwnerSet::contains(index->local().attribute()))
601 if(part[index->local()]==toPe)
603 getNeighbor<OwnerSet>(graph, part, index->local(), indexSet,
604 toPe, overlapSet, neighborProcs);
605 my_push_back(ownerVec, index->global(), toPe);
609 reserve(ownerVec, redist, toPe);
620 template<
class F,
class IS>
621 inline bool isOwner(IS& indexSet,
int index) {
623 const typename IS::IndexPair* pindex=indexSet.pair(index);
626 return F::contains(pindex->local().attribute());
630 class BaseEdgeFunctor
633 BaseEdgeFunctor(
int* adj,
const ParmetisDuneIndexMap& data)
638 void operator()(
const T& edge)
658 :
public BaseEdgeFunctor
660 EdgeFunctor(
int* adj,
const ParmetisDuneIndexMap& data, std::size_t s)
661 : BaseEdgeFunctor(adj, data)
671 template<
class G,
class V,
class E,
class VM,
class EM>
672 class EdgeFunctor<Dune::Amg::PropertiesGraph<G,V,E,VM,EM> >
673 :
public BaseEdgeFunctor
676 EdgeFunctor(
int* adj,
const ParmetisDuneIndexMap& data, std::size_t s)
677 :BaseEdgeFunctor(adj, data)
683 void operator()(
const T& edge)
685 weight_[
index()]=edge.properties().depends()?3:1;
686 BaseEdgeFunctor::operator()(edge);
717 template<
class F,
class G,
class IS,
class EW>
718 void getAdjArrays(G& graph, IS& indexSet,
int *xadj,
724 typedef typename G::ConstVertexIterator VertexIterator;
726 typedef typename IS::const_iterator Iterator;
728 VertexIterator vend = graph.end();
731 for(VertexIterator vertex = graph.begin(); vertex != vend; ++vertex){
732 if (isOwner<F>(indexSet,*vertex)) {
734 typedef typename G::ConstEdgeIterator EdgeIterator;
735 EdgeIterator eend = vertex.end();
736 xadj[j] = ew.index();
738 for(EdgeIterator edge = vertex.begin(); edge != eend; ++edge){
743 xadj[j] = ew.index();
747 template<
class G,
class T1,
class T2>
751 RedistributeInterface& redistInf,
756 idxtype *adjwgt,
int *wgtflag,
int *numflag,
int *nparts,
757 int *options,
int *edgecut,
idxtype *part);
760 idxtype *adjwgt,
int *wgtflag,
int *numflag,
int *nparts,
761 int *options,
int *edgecut,
idxtype *part);
765 template<
class S,
class T>
768 for(T *cur=array, *end=array+l; cur!=end; ++cur)
776 idxtype* adjncy,
bool checkSymmetry)
781 if(xadj[vtx]>noEdges||xadj[vtx]<0){
782 std::cerr <<
"Check graph: xadj["<<vtx<<
"]="<<xadj[vtx]<<
" (>"
783 <<noEdges<<
") out of range!"<<std::endl;
786 if(xadj[vtx+1]>noEdges||xadj[vtx+1]<0){
787 std::cerr <<
"Check graph: xadj["<<vtx+1<<
"]="<<xadj[vtx+1]<<
" (>"
788 <<noEdges<<
") out of range!"<<std::endl;
792 for(
idxtype i=xadj[vtx]; i< xadj[vtx+1];++i){
793 if(adjncy[i]<0||((std::size_t)adjncy[i])>gnoVtx){
794 std::cerr<<
" Edge "<<adjncy[i]<<
" out of range ["<<0<<
","<<noVtx<<
")"
800 for(
idxtype i=xadj[vtx]; i< xadj[vtx+1];++i){
804 for(
idxtype j=xadj[target]; j< xadj[target+1];++j)
808 std::cerr<<
"Edge ("<<target<<
","<<vtx<<
") "<<i<<
" time"<<std::endl;
817 template<
class M,
class T1,
class T2>
824 std::cout<<
"Repartitioning from "<<oocomm.
communicator().size()
825 <<
" to "<<nparts<<
" parts"<<std::endl;
829 int* part =
new int[1];
844 typedef typename RemoteIndices::const_iterator
869 xadj[1]=noNeighbours;
899 typedef typename RemoteIndices::const_iterator NeighbourIterator;
901 typedef typename IndexSet::LocalIndex LocalIndex;
909 adjwgt =
new idxtype[noNeighbours];
915 if(n->first != rank){
925 noNeighbours, xadj, adjncy,
false));
927 int wgtflag=0, numflag=0, edgecut;
931 float *tpwgts =
new float[nparts];
932 for(
int i=0; i<nparts; ++i)
933 tpwgts[i]=1.0/nparts;
934 int options[5] ={ 0,1,15,0,0};
937 Dune::dinfo<<rank<<
" vtxdist: ";
939 Dune::dinfo<<std::endl<<rank<<
" xadj: ";
941 Dune::dinfo<<std::endl<<rank<<
" adjncy: ";
945 Dune::dinfo<<std::endl<<rank<<
" vwgt: ";
947 Dune::dinfo<<std::endl<<rank<<
" adwgt: ";
950 Dune::dinfo<<std::endl;
953 std::cout<<
"Creating comm graph took "<<time.elapsed()<<std::endl;
956 #ifdef PARALLEL_PARTITION
963 ParMETIS_V3_PartKway(vtxdist, xadj, adjncy,
964 vwgt, adjwgt, &wgtflag,
965 &numflag, &ncon, &nparts, tpwgts, &ubvec, options, &edgecut, part,
968 std::cout<<
"ParMETIS took "<<time.elapsed()<<std::endl;
972 std::size_t gnoedges=0;
975 Dune::dverb<<
"noNeighbours: "<<noNeighbours<<std::endl;
977 MPI_Allgather(&noNeighbours,1,MPI_INT,noedges,1, MPI_INT,oocomm.
communicator());
980 std::cout<<
"Gathering noedges took "<<time1.elapsed()<<std::endl;
994 std::size_t localNoVtx=vtxdist[rank+1]-vtxdist[rank];
998 Dune::dinfo<<
"noedges: ";
1000 Dune::dinfo<<std::endl;
1008 noxs[i]=vtxdist[i+1]-vtxdist[i]+1;
1009 novs[i]=vtxdist[i+1]-vtxdist[i];
1014 for(
int *xcurr = xdispl, *vcurr = vdispl, *end=vdispl+oocomm.
communicator().size();
1015 vcurr!=end; ++vcurr, ++xcurr, ++so, ++offset){
1017 *xcurr = offset + *so;
1023 for(
int *curr=noedges, *end=noedges+oocomm.
communicator().size()-1;
1029 Dune::dinfo<<
"displ: ";
1031 Dune::dinfo<<std::endl;
1035 for(
int *curr=noedges, *end=noedges+oocomm.
communicator().size();
1040 Dune::dinfo<<
"gxadjlen: "<<gxadjlen<<
" noVertices: "<<noVertices
1041 <<
" gnoedges: "<<gnoedges<<std::endl;
1042 gxadj =
new idxtype[gxadjlen];
1043 gpart =
new idxtype[noVertices];
1045 gvwgt =
new idxtype[noVertices];
1046 gadjwgt =
new idxtype[gnoedges];
1048 gadjncy =
new idxtype[gnoedges];
1052 std::cout<<
"Preparing global graph took "<<time1.elapsed()<<std::endl;
1056 MPI_Allgatherv(xadj,2,MPITraits<idxtype>::getType(),
1057 gxadj,noxs,xdispl,MPITraits<idxtype>::getType(),
1059 MPI_Allgatherv(adjncy,noNeighbours,MPITraits<idxtype>::getType(),
1060 gadjncy,noedges,displ,MPITraits<idxtype>::getType(),
1063 MPI_Allgatherv(adjwgt,noNeighbours,MPITraits<idxtype>::getType(),
1064 gadjwgt,noedges,displ,MPITraits<idxtype>::getType(),
1066 MPI_Allgatherv(vwgt,localNoVtx,MPITraits<idxtype>::getType(),
1067 gvwgt,novs,vdispl,MPITraits<idxtype>::getType(),
1071 std::cout<<
"Gathering global graph data took "<<time1.elapsed()<<std::endl;
1081 idxtype increment = vtxdist[1];
1085 int lprev = vtxdist[i]-vtxdist[i-1];
1086 int l = vtxdist[i+1]-vtxdist[i];
1088 assert((start+l+offset)-gxadj<=static_cast<idxtype>(gxadjlen));
1089 increment = *(start-1);
1090 std::transform(start+offset, start+l+offset, start, std::bind2nd(std::plus<idxtype>(), increment));
1092 Dune::dinfo<<std::endl<<
"shifted xadj:";
1094 Dune::dinfo<<std::endl<<
" gadjncy: ";
1097 Dune::dinfo<<std::endl<<
" gvwgt: ";
1099 Dune::dinfo<<std::endl<<
"adjwgt: ";
1101 Dune::dinfo<<std::endl;
1105 std::cout<<
"Postprocesing global graph data took "<<time1.elapsed()<<std::endl;
1109 gxadj, gadjncy,
true));
1113 std::cout<<
"Creating grah one 1 process took "<<time.elapsed()<<std::endl;
1115 options[0]=0; options[1]=1; options[2]=1; options[3]=3; options[4]=3;
1117 METIS_PartGraphRecursive(&noVertices, gxadj, gadjncy, gvwgt, gadjwgt, &wgtflag,
1118 &numflag, &nparts, options, &edgecut, gpart);
1121 std::cout<<
"METIS took "<<time.elapsed()<<std::endl;
1124 Dune::dinfo<<std::endl<<
"part:";
1135 MPI_Scatter(gpart, 1, MPITraits<idxtype>::getType(), part, 1,
1136 MPITraits<idxtype>::getType(), 0, comm);
1160 Dune::dinfo<<
" repart "<<rank <<
" -> "<< part[0]<<std::endl;
1162 std::vector<int> realpart(mat.N(), part[0]);
1168 std::cout<<
"Scattering repartitioning took "<<time.elapsed()<<std::endl;
1176 std::cout<<
"Filling index set took "<<time.elapsed()<<std::endl;
1184 std::cout<<
"Average no neighbours was "<<noNeighbours<<std::endl;
1189 std::cout<<
"Building index sets took "<<time.elapsed()<<std::endl;
1211 template<
class G,
class T1,
class T2>
1224 std::cout<<
"Filling holes took "<<time.elapsed()<<std::endl;
1231 double t1=0.0, t2=0.0, t3=0.0, t4=0.0, tSum=0.0;
1259 typedef typename OOComm::OwnerSet OwnerSet;
1264 ParmetisDuneIndexMap indexMap(graph,oocomm);
1268 std::size_t *part =
new std::size_t[indexMap.numOfOwnVtx()];
1270 for(std::size_t i=0; i < indexMap.numOfOwnVtx(); ++i)
1275 std::cerr<<
"ParMETIS not activated. Will repartition to 1 domain instead of requested "
1276 <<nparts<<
" domains."<<std::endl;
1285 EdgeFunctor<G> ef(adjncy, indexMap, graph.noEdges());
1286 getAdjArrays<OwnerSet>(graph, oocomm.
globalLookup(), xadj, ef);
1292 int numflag=0, wgtflag=0, options[3], edgecut=0, ncon=1;
1294 float *tpwgts =
new float[nparts];
1295 for(
int i=0; i<nparts; ++i)
1296 tpwgts[i]=1.0/nparts;
1305 wgtflag = (ef.getWeights()!=NULL)?1:0;
1313 std::cout<<std::endl;
1314 std::cout<<
"Testing ParMETIS_V3_PartKway with options[1-2] = {"
1315 <<options[1]<<
" "<<options[2]<<
"}, Ncon: "
1316 <<ncon<<
", Nparts: "<<nparts<<std::endl;
1329 std::cout<<
"Preparing for parmetis took "<<time.elapsed()<<std::endl;
1336 ParMETIS_V3_PartKway(indexMap.vtxDist(), xadj, adjncy,
1337 NULL, ef.getWeights(), &wgtflag,
1338 &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgecut, part, &
const_cast<MPI_Comm&
>(comm));
1349 std::cout<<std::endl;
1350 std::cout<<
"ParMETIS_V3_PartKway reported a cut of "<<edgecut<<std::endl;
1351 std::cout<<std::endl;
1353 std::cout<<mype<<
": PARMETIS-Result: ";
1354 for(
int i=0; i < indexMap.vtxDist()[mype+1]-indexMap.vtxDist()[mype]; ++i) {
1355 std::cout<<part[i]<<
" ";
1357 std::cout<<std::endl;
1358 std::cout<<
"Testing ParMETIS_V3_PartKway with options[1-2] = {"
1359 <<options[1]<<
" "<<options[2]<<
"}, Ncon: "
1360 <<ncon<<
", Nparts: "<<nparts<<std::endl;
1373 std::cout<<
"Parmetis took "<<time.elapsed()<<std::endl;
1380 for(std::size_t i=0; i<indexMap.numOfOwnVtx();++i)
1390 std::vector<int> domainMapping(nparts);
1392 getDomain(comm, part, indexMap.numOfOwnVtx(), nparts, &myDomain, domainMapping);
1397 std::cout<<mype<<
": myDomain: "<<myDomain<<std::endl;
1398 std::cout<<mype<<
": DomainMapping: ";
1399 for(
int j=0; j<nparts; j++) {
1400 std::cout<<
" do: "<<j<<
" pe: "<<domainMapping[j]<<
" ";
1402 std::cout<<std::endl;
1409 std::vector<int> setPartition(oocomm.
indexSet().size(), -1);
1411 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
1414 if(OwnerSet::contains(index->local().attribute())){
1415 setPartition[index->local()]=domainMapping[part[i++]];
1423 static_cast<int>(SolverCategory::nonoverlapping))
1430 std::cout<<
"Creating indexsets took "<<time.elapsed()<<std::endl;
1437 template<
class G,
class T1,
class T2>
1445 typedef typename OOComm::OwnerSet OwnerSet;
1479 std::set<int> recvFrom;
1485 typedef typename std::vector<int>::const_iterator VIter;
1489 std::set<int> tsendTo;
1490 for(VIter i=setPartition.begin(), iend = setPartition.end(); i!=iend; ++i)
1493 noSendTo = tsendTo.size();
1494 sendTo =
new int[noSendTo];
1495 typedef std::set<int>::const_iterator iterator;
1497 for(iterator i=tsendTo.begin(); i != tsendTo.end(); ++i, ++idx)
1503 int* gsendToDispl =
new int[oocomm.
communicator().size()+1];
1505 MPI_Allgather(&noSendTo, 1, MPI_INT, gnoSend, 1,
1509 int totalNoRecv = 0;
1510 for(
int i=0; i<npes; ++i)
1511 totalNoRecv += gnoSend[i];
1513 int *gsendTo =
new int[totalNoRecv];
1517 for(
int i=0; i<npes; ++i)
1518 gsendToDispl[i+1]=gsendToDispl[i]+gnoSend[i];
1521 MPI_Allgatherv(sendTo, noSendTo, MPI_INT, gsendTo, gnoSend, gsendToDispl,
1525 for(
int proc=0; proc < npes; ++proc)
1526 for(
int i=gsendToDispl[proc]; i < gsendToDispl[proc+1]; ++i)
1527 if(gsendTo[i]==mype)
1528 recvFrom.insert(proc);
1530 bool existentOnNextLevel = recvFrom.size()>0;
1534 delete[] gsendToDispl;
1539 if(recvFrom.size()){
1540 std::cout<<mype<<
": recvFrom: ";
1541 typedef typename std::set<int>::const_iterator siter;
1542 for(siter i=recvFrom.begin(); i!= recvFrom.end(); ++i) {
1547 std::cout<<std::endl<<std::endl;
1548 std::cout<<mype<<
": sendTo: ";
1549 for(
int i=0; i<noSendTo; i++) {
1550 std::cout<<sendTo[i]<<
" ";
1552 std::cout<<std::endl<<std::endl;
1557 std::cout<<
" Communicating the receive information took "<<
1558 time.elapsed()<<std::endl;
1572 typedef typename OOComm::ParallelIndexSet::GlobalIndex GI;
1573 typedef std::vector<GI> GlobalVector;
1574 std::vector<std::pair<GI,int> > myOwnerVec;
1575 std::set<GI> myOverlapSet;
1576 GlobalVector sendOwnerVec;
1577 std::set<GI> sendOverlapSet;
1578 std::set<int> myNeighbors;
1583 char **sendBuffers=
new char*[noSendTo];
1584 MPI_Request *requests =
new MPI_Request[noSendTo];
1587 for(
int i=0; i < noSendTo; ++i){
1589 sendOwnerVec.clear();
1590 sendOverlapSet.clear();
1593 std::set<int> neighbors;
1594 getOwnerOverlapVec<OwnerSet>(graph, setPartition, oocomm.
globalLookup(),
1595 mype, sendTo[i], sendOwnerVec, sendOverlapSet, redistInf,
1601 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.
communicator(), &buffersize);
1602 MPI_Pack_size(sendOwnerVec.size(), MPITraits<GI>::getType(), oocomm.
communicator(), &tsize);
1604 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.
communicator(), &tsize);
1606 MPI_Pack_size(sendOverlapSet.size(), MPITraits<GI>::getType(), oocomm.
communicator(), &tsize);
1607 buffersize += tsize;
1608 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.
communicator(), &tsize);
1609 buffersize += tsize;
1610 MPI_Pack_size(neighbors.size(), MPI_INT, oocomm.
communicator(), &tsize);
1611 buffersize += tsize;
1613 sendBuffers[i] =
new char[buffersize];
1616 std::cout<<mype<<
" sending "<<sendOwnerVec.size()<<
" owner and "<<
1617 sendOverlapSet.size()<<
" overlap to "<<sendTo[i]<<
" buffersize="<<buffersize<<std::endl;
1619 createSendBuf(sendOwnerVec, sendOverlapSet, neighbors, sendBuffers[i], buffersize, oocomm.
communicator());
1620 MPI_Issend(sendBuffers[i], buffersize, MPI_PACKED, sendTo[i], 99, oocomm.
communicator(), requests+i);
1626 std::cout<<
" Creating sends took "<<
1627 time.elapsed()<<std::endl;
1632 int noRecv = recvFrom.size();
1633 int oldbuffersize=0;
1638 MPI_Probe(MPI_ANY_SOURCE, 99, oocomm.
communicator(), &stat);
1640 MPI_Get_count(&stat, MPI_PACKED, &buffersize);
1642 if(oldbuffersize<buffersize){
1645 recvBuf =
new char[buffersize];
1646 oldbuffersize = buffersize;
1648 MPI_Recv(recvBuf, buffersize, MPI_PACKED, stat.MPI_SOURCE, 99, oocomm.
communicator(), &stat);
1649 saveRecvBuf(recvBuf, buffersize, myOwnerVec, myOverlapSet, myNeighbors, redistInf,
1659 MPI_Status *statuses =
new MPI_Status[noSendTo];
1660 int send = MPI_Waitall(noSendTo, requests, statuses);
1663 if(send==MPI_ERR_IN_STATUS){
1664 std::cerr<<mype<<
": Error in sending :"<<std::endl;
1666 for(
int i=0; i< noSendTo; i++)
1667 if(statuses[i].MPI_ERROR!=MPI_SUCCESS){
1670 MPI_Error_string(statuses[i].MPI_ERROR, message, &messageLength);
1671 std::cerr<<
" source="<<statuses[i].MPI_SOURCE<<
" message: ";
1672 for(
int i=0; i< messageLength; i++)
1673 std::cout<<message[i];
1675 std::cerr<<std::endl;
1681 std::cout<<
" Receiving and saving took "<<
1682 time.elapsed()<<std::endl;
1686 for(
int i=0; i < noSendTo; ++i)
1687 delete[] sendBuffers[i];
1689 delete[] sendBuffers;
1704 if (!existentOnNextLevel) {
1706 color= MPI_UNDEFINED;
1708 MPI_Comm outputComm;
1716 std::vector<int> tneighbors;
1717 tneighbors.reserve(myNeighbors.size());
1719 typename OOComm::ParallelIndexSet& outputIndexSet = outcomm->indexSet();
1721 MPI_Allgather(&newrank, 1, MPI_INT, newranks, 1,
1723 typedef typename std::set<int>::const_iterator IIter;
1727 for(IIter i=myNeighbors.begin(), end=myNeighbors.end();
1729 assert(newranks[*i]>=0);
1730 std::cout<<*i<<
"->"<<newranks[*i]<<
" ";
1731 tneighbors.push_back(newranks[*i]);
1733 std::cout<<std::endl;
1735 for(IIter i=myNeighbors.begin(), end=myNeighbors.end();
1737 tneighbors.push_back(newranks[*i]);
1741 myNeighbors.clear();
1746 std::cout<<
" Calculating new neighbours ("<<tneighbors.size()<<
") took "<<
1747 time.elapsed()<<std::endl;
1752 outputIndexSet.beginResize();
1755 std::sort(myOwnerVec.begin(), myOwnerVec.end(), SortFirst());
1759 typedef typename OOComm::ParallelIndexSet::LocalIndex LocalIndex;
1760 typedef typename std::vector<std::pair<GI,int> >::const_iterator VPIter;
1762 for(VPIter g=myOwnerVec.begin(), end =myOwnerVec.end(); g!=end; ++g, ++i ) {
1763 outputIndexSet.add(g->first,LocalIndex(i, OwnerOverlapCopyAttributeSet::owner,
true));
1770 std::cout<<
" Adding owner indices took "<<
1771 time.elapsed()<<std::endl;
1780 mergeVec(myOwnerVec, myOverlapSet);
1784 myOwnerVec.swap(myOwnerVec);
1789 std::cout<<
" Merging indices took "<<
1790 time.elapsed()<<std::endl;
1796 typedef typename std::set<GI>::const_iterator SIter;
1797 for(SIter g=myOverlapSet.begin(), end=myOverlapSet.end(); g!=end; ++g, i++) {
1798 outputIndexSet.add(*g,LocalIndex(i, OwnerOverlapCopyAttributeSet::copy,
true));
1800 myOverlapSet.clear();
1801 outputIndexSet.endResize();
1803 #ifdef DUNE_ISTL_WITH_CHECKING
1805 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
1806 Iterator end = outputIndexSet.end();
1807 for(Iterator index = outputIndexSet.begin(); index != end; ++
index) {
1808 if (OwnerSet::contains(index->local().attribute())) {
1819 Iterator index=outputIndexSet.begin();
1822 for(Iterator old = outputIndexSet.begin(); index != end; old=index++) {
1823 if(old->global()>index->global())
1824 DUNE_THROW(
ISTLError,
"Index set's globalindex not sorted correctly");
1831 std::cout<<
" Adding overlap indices took "<<
1832 time.elapsed()<<std::endl;
1837 if(color != MPI_UNDEFINED){
1838 outcomm->remoteIndices().setNeighbours(tneighbors);
1839 outcomm->remoteIndices().template rebuild<true>();
1849 std::cout<<
" Storing indexsets took "<<
1850 time.elapsed()<<std::endl;
1856 tSum = t1 + t2 + t3 + t4;
1857 std::cout<<std::endl
1858 <<mype<<
": WTime for step 1): "<<t1
1866 return color!=MPI_UNDEFINED;
1870 template<
class G,
class P,
class T1,
class T2,
class R>
1876 if(nparts!=oocomm.size())
1877 DUNE_THROW(NotImplemented,
"only available for MPI programs");
1881 template<
class G,
class P,
class T1,
class T2,
class R>
1887 if(nparts!=oocomm.size())
1888 DUNE_THROW(NotImplemented,
"only available for MPI programs");