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 typedef std::vector<int>::const_iterator Iter;
244 std::map<int,int> sizes;
246 typedef typename IS::const_iterator IIter;
247 for(IIter i=idxset.begin(), end=idxset.end();i!=end; ++i)
248 if(Flags::contains(i->local().attribute()))
249 ++sizes[toPart[i->local()]];
252 typedef std::map<int,int>::const_iterator MIter;
253 for(MIter i=sizes.begin(), end=sizes.end(); i!=end; ++i)
254 interfaces()[i->first].first.reserve(i->second);
257 typedef typename IS::const_iterator IIter;
258 for(IIter i=idxset.begin(), end=idxset.end();i!=end; ++i)
259 if(Flags::contains(i->local().attribute()))
260 interfaces()[toPart[i->local()]].first.add(i->local());
263 void reserveSpaceForReceiveInterface(
int proc,
int size)
265 interfaces()[proc].second.reserve(size);
267 void addReceiveIndex(
int proc, std::size_t idx)
269 interfaces()[proc].second.add(idx);
271 template<
typename TG>
272 void buildReceiveInterface(std::vector<std::pair<TG,int> >& indices)
274 typedef typename std::vector<std::pair<TG,int> >::const_iterator VIter;
276 for(VIter idx=indices.begin(); idx!= indices.end(); ++idx){
277 interfaces()[idx->second].second.add(i++);
299 void createSendBuf(std::vector<GI>& ownerVec, std::set<GI>& overlapVec, std::set<int>& neighbors,
char *sendBuf,
int buffersize, MPI_Comm comm) {
301 std::size_t s=ownerVec.size();
305 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
306 MPI_Pack(&(ownerVec[0]), s, MPITraits<GI>::getType(), sendBuf, buffersize, &pos, comm);
307 s = overlapVec.size();
308 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
309 typedef typename std::set<GI>::iterator Iter;
310 for(Iter i=overlapVec.begin(), end= overlapVec.end(); i != end; ++i)
311 MPI_Pack(const_cast<GI*>(&(*i)), 1, MPITraits<GI>::getType(), sendBuf, buffersize, &pos, comm);
314 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
315 typedef typename std::set<int>::iterator IIter;
317 for(IIter i=neighbors.begin(), end= neighbors.end(); i != end; ++i)
318 MPI_Pack(const_cast<int*>(&(*i)), 1, MPI_INT, sendBuf, buffersize, &pos, comm);
329 void saveRecvBuf(
char *recvBuf,
int bufferSize, std::vector<std::pair<GI,int> >& ownerVec,
330 std::set<GI>& overlapVec, std::set<int>& neighbors, RedistributeInterface& inf,
int from, MPI_Comm comm) {
334 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
335 inf.reserveSpaceForReceiveInterface(from, size);
336 ownerVec.reserve(ownerVec.size()+size);
337 for(;size!=0;--size){
339 MPI_Unpack(recvBuf, bufferSize, &pos, &gi, 1, MPITraits<GI>::getType(), comm);
340 ownerVec.push_back(std::make_pair(gi,from));
343 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
344 typename std::set<GI>::iterator ipos = overlapVec.begin();
345 Dune::dverb <<
"unpacking "<<size<<
" overlap"<<std::endl;
346 for(;size!=0;--size){
348 MPI_Unpack(recvBuf, bufferSize, &pos, &gi, 1, MPITraits<GI>::getType(), comm);
349 ipos=overlapVec.insert(ipos, gi);
352 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
353 Dune::dverb <<
"unpacking "<<size<<
" neighbors"<<std::endl;
354 typename std::set<int>::iterator npos = neighbors.begin();
355 for(;size!=0;--size){
357 MPI_Unpack(recvBuf, bufferSize, &pos, &n, 1, MPI_INT, comm);
358 npos=neighbors.insert(npos, n);
376 void getDomain(
const MPI_Comm& comm, T *part,
int numOfOwnVtx,
int nparts,
int *myDomain,
int domainMapping[]) {
378 MPI_Comm_size(comm, &npes);
379 MPI_Comm_rank(comm, &mype);
389 for (i=0; i<nparts; i++) {
390 domainMapping[i] = -1;
393 for (i=0; i<npes; i++) {
397 for (i=0; i<numOfOwnVtx; i++) {
401 int *domainMatrix =
new int[npes * nparts];
403 for(i=0; i<npes*nparts; i++) {
408 int *buf =
new int[nparts];
409 for (i=0; i<nparts; i++) {
411 domainMatrix[mype*nparts+i] = domain[i];
414 int src = (mype-1+npes)%npes;
415 int dest = (mype+1)%npes;
417 for (i=0; i<npes-1; i++) {
418 MPI_Sendrecv_replace(buf, nparts, MPI_INT, dest, 0, src, 0, comm, &status);
420 pe = ((mype-1-i)+npes)%npes;
421 for(j=0; j<nparts; j++) {
423 domainMatrix[pe*nparts+j] = buf[j];
431 int maxOccurance = 0;
433 for(i=0; i<nparts; i++) {
434 for(j=0; j<npes; j++) {
436 if (assigned[j]==0) {
437 if (maxOccurance < domainMatrix[j*nparts+i]) {
438 maxOccurance = domainMatrix[j*nparts+i];
446 domainMapping[i] = pe;
457 delete[] domainMatrix;
464 bool operator()(
const T& t1,
const T& t2)
const
482 void mergeVec(std::vector<std::pair<GI, int> >& ownerVec, std::set<GI>& overlapSet) {
484 typedef typename std::vector<std::pair<GI,int> >::const_iterator VIter;
487 if(ownerVec.size()>0)
489 VIter old=ownerVec.begin();
490 for(VIter i=old+1, end=ownerVec.end(); i != end; old=i++)
492 if(i->first==old->first)
494 std::cerr<<
"Value at indes"<<old-ownerVec.begin()<<
" is the same as at index "
495 <<i-ownerVec.begin()<<
" ["<<old->first<<
","<<old->second<<
"]==["
496 <<i->first<<
","<<i->second<<
"]"<<std::endl;
504 typedef typename std::set<GI>::iterator SIter;
505 VIter v=ownerVec.begin(), vend=ownerVec.end();
506 for(SIter s=overlapSet.begin(), send=overlapSet.end(); s!=send;)
508 while(v!=vend && v->first<*s) ++v;
509 if(v!=vend && v->first==*s){
514 overlapSet.erase(tmp);
534 template<
class OwnerSet,
class Graph,
class IS,
class GI>
535 void getNeighbor(
const Graph& g, std::vector<int>& part,
536 typename Graph::VertexDescriptor vtx,
const IS& indexSet,
537 int toPe, std::set<GI>& neighbor, std::set<int>& neighborProcs) {
538 typedef typename Graph::ConstEdgeIterator Iter;
539 for(Iter edge=g.beginEdges(vtx), end=g.endEdges(vtx); edge!=end; ++edge)
541 const typename IS::IndexPair* pindex = indexSet.pair(edge.target());
543 if(part[pindex->local()]!=toPe || !OwnerSet::contains(pindex->local().attribute()))
546 neighbor.insert(pindex->global());
547 neighborProcs.insert(part[pindex->local()]);
552 template<
class T,
class I>
553 void my_push_back(std::vector<T>& ownerVec,
const I& index,
int proc)
555 ownerVec.push_back(index);
558 template<
class T,
class I>
559 void my_push_back(std::vector<std::pair<T,int> >& ownerVec,
const I& index,
int proc)
561 ownerVec.push_back(std::make_pair(index,proc));
564 void reserve(std::vector<T>&, RedistributeInterface&,
int)
568 void reserve(std::vector<std::pair<T,int> >& ownerVec, RedistributeInterface& redist,
int proc)
570 redist.reserveSpaceForReceiveInterface(proc, ownerVec.size());
591 template<
class OwnerSet,
class G,
class IS,
class T,
class GI>
592 void getOwnerOverlapVec(
const G& graph, std::vector<int>& part, IS& indexSet,
593 int myPe,
int toPe, std::vector<T>& ownerVec, std::set<GI>& overlapSet,
594 RedistributeInterface& redist, std::set<int>& neighborProcs) {
597 typedef typename IS::const_iterator Iterator;
598 for(Iterator index = indexSet.begin(); index != indexSet.end(); ++
index) {
600 if(OwnerSet::contains(index->local().attribute()))
602 if(part[index->local()]==toPe)
604 getNeighbor<OwnerSet>(graph, part, index->local(), indexSet,
605 toPe, overlapSet, neighborProcs);
606 my_push_back(ownerVec, index->global(), toPe);
610 reserve(ownerVec, redist, toPe);
621 template<
class F,
class IS>
622 inline bool isOwner(IS& indexSet,
int index) {
624 const typename IS::IndexPair* pindex=indexSet.pair(index);
627 return F::contains(pindex->local().attribute());
631 class BaseEdgeFunctor
634 BaseEdgeFunctor(
int* adj,
const ParmetisDuneIndexMap& data)
639 void operator()(
const T& edge)
659 :
public BaseEdgeFunctor
661 EdgeFunctor(
int* adj,
const ParmetisDuneIndexMap& data, std::size_t s)
662 : BaseEdgeFunctor(adj, data)
672 template<
class G,
class V,
class E,
class VM,
class EM>
673 class EdgeFunctor<Dune::Amg::PropertiesGraph<G,V,E,VM,EM> >
674 :
public BaseEdgeFunctor
677 EdgeFunctor(
int* adj,
const ParmetisDuneIndexMap& data, std::size_t s)
678 :BaseEdgeFunctor(adj, data)
684 void operator()(
const T& edge)
686 weight_[
index()]=edge.properties().depends()?3:1;
687 BaseEdgeFunctor::operator()(edge);
718 template<
class F,
class G,
class IS,
class EW>
719 void getAdjArrays(G& graph, IS& indexSet,
int *xadj,
725 typedef typename G::ConstVertexIterator VertexIterator;
727 typedef typename IS::const_iterator Iterator;
729 VertexIterator vend = graph.end();
732 for(VertexIterator vertex = graph.begin(); vertex != vend; ++vertex){
733 if (isOwner<F>(indexSet,*vertex)) {
735 typedef typename G::ConstEdgeIterator EdgeIterator;
736 EdgeIterator eend = vertex.end();
737 xadj[j] = ew.index();
739 for(EdgeIterator edge = vertex.begin(); edge != eend; ++edge){
744 xadj[j] = ew.index();
748 template<
class G,
class T1,
class T2>
752 RedistributeInterface& redistInf,
757 idxtype *adjwgt,
int *wgtflag,
int *numflag,
int *nparts,
758 int *options,
int *edgecut,
idxtype *part);
761 idxtype *adjwgt,
int *wgtflag,
int *numflag,
int *nparts,
762 int *options,
int *edgecut,
idxtype *part);
766 template<
class S,
class T>
769 for(T *cur=array, *end=array+l; cur!=end; ++cur)
777 idxtype* adjncy,
bool checkSymmetry)
782 if(xadj[vtx]>noEdges||xadj[vtx]<0){
783 std::cerr <<
"Check graph: xadj["<<vtx<<
"]="<<xadj[vtx]<<
" (>"
784 <<noEdges<<
") out of range!"<<std::endl;
787 if(xadj[vtx+1]>noEdges||xadj[vtx+1]<0){
788 std::cerr <<
"Check graph: xadj["<<vtx+1<<
"]="<<xadj[vtx+1]<<
" (>"
789 <<noEdges<<
") out of range!"<<std::endl;
793 for(
idxtype i=xadj[vtx]; i< xadj[vtx+1];++i){
794 if(adjncy[i]<0||((std::size_t)adjncy[i])>gnoVtx){
795 std::cerr<<
" Edge "<<adjncy[i]<<
" out of range ["<<0<<
","<<noVtx<<
")"
801 for(
idxtype i=xadj[vtx]; i< xadj[vtx+1];++i){
805 for(
idxtype j=xadj[target]; j< xadj[target+1];++j)
809 std::cerr<<
"Edge ("<<target<<
","<<vtx<<
") "<<i<<
" time"<<std::endl;
818 template<
class M,
class T1,
class T2>
825 std::cout<<
"Repartitioning from "<<oocomm.
communicator().size()
826 <<
" to "<<nparts<<
" parts"<<std::endl;
830 int* part =
new int[1];
845 typedef typename RemoteIndices::const_iterator
870 xadj[1]=noNeighbours;
900 typedef typename RemoteIndices::const_iterator NeighbourIterator;
902 typedef typename IndexSet::LocalIndex LocalIndex;
910 adjwgt =
new idxtype[noNeighbours];
916 if(n->first != rank){
926 noNeighbours, xadj, adjncy,
false));
928 int wgtflag=0, numflag=0, edgecut;
932 float *tpwgts =
new float[nparts];
933 for(
int i=0; i<nparts; ++i)
934 tpwgts[i]=1.0/nparts;
935 int options[5] ={ 0,1,15,0,0};
938 Dune::dinfo<<rank<<
" vtxdist: ";
940 Dune::dinfo<<std::endl<<rank<<
" xadj: ";
942 Dune::dinfo<<std::endl<<rank<<
" adjncy: ";
946 Dune::dinfo<<std::endl<<rank<<
" vwgt: ";
948 Dune::dinfo<<std::endl<<rank<<
" adwgt: ";
951 Dune::dinfo<<std::endl;
954 std::cout<<
"Creating comm graph took "<<time.elapsed()<<std::endl;
957 #ifdef PARALLEL_PARTITION
964 ParMETIS_V3_PartKway(vtxdist, xadj, adjncy,
965 vwgt, adjwgt, &wgtflag,
966 &numflag, &ncon, &nparts, tpwgts, &ubvec, options, &edgecut, part,
969 std::cout<<
"ParMETIS took "<<time.elapsed()<<std::endl;
973 std::size_t gnoedges=0;
976 Dune::dverb<<
"noNeighbours: "<<noNeighbours<<std::endl;
978 MPI_Allgather(&noNeighbours,1,MPI_INT,noedges,1, MPI_INT,oocomm.
communicator());
981 std::cout<<
"Gathering noedges took "<<time1.elapsed()<<std::endl;
995 std::size_t localNoVtx=vtxdist[rank+1]-vtxdist[rank];
999 Dune::dinfo<<
"noedges: ";
1001 Dune::dinfo<<std::endl;
1009 noxs[i]=vtxdist[i+1]-vtxdist[i]+1;
1010 novs[i]=vtxdist[i+1]-vtxdist[i];
1015 for(
int *xcurr = xdispl, *vcurr = vdispl, *end=vdispl+oocomm.
communicator().size();
1016 vcurr!=end; ++vcurr, ++xcurr, ++so, ++offset){
1018 *xcurr = offset + *so;
1024 for(
int *curr=noedges, *end=noedges+oocomm.
communicator().size()-1;
1030 Dune::dinfo<<
"displ: ";
1032 Dune::dinfo<<std::endl;
1036 for(
int *curr=noedges, *end=noedges+oocomm.
communicator().size();
1041 Dune::dinfo<<
"gxadjlen: "<<gxadjlen<<
" noVertices: "<<noVertices
1042 <<
" gnoedges: "<<gnoedges<<std::endl;
1043 gxadj =
new idxtype[gxadjlen];
1044 gpart =
new idxtype[noVertices];
1046 gvwgt =
new idxtype[noVertices];
1047 gadjwgt =
new idxtype[gnoedges];
1049 gadjncy =
new idxtype[gnoedges];
1053 std::cout<<
"Preparing global graph took "<<time1.elapsed()<<std::endl;
1057 MPI_Allgatherv(xadj,2,MPITraits<idxtype>::getType(),
1058 gxadj,noxs,xdispl,MPITraits<idxtype>::getType(),
1060 MPI_Allgatherv(adjncy,noNeighbours,MPITraits<idxtype>::getType(),
1061 gadjncy,noedges,displ,MPITraits<idxtype>::getType(),
1064 MPI_Allgatherv(adjwgt,noNeighbours,MPITraits<idxtype>::getType(),
1065 gadjwgt,noedges,displ,MPITraits<idxtype>::getType(),
1067 MPI_Allgatherv(vwgt,localNoVtx,MPITraits<idxtype>::getType(),
1068 gvwgt,novs,vdispl,MPITraits<idxtype>::getType(),
1072 std::cout<<
"Gathering global graph data took "<<time1.elapsed()<<std::endl;
1082 idxtype increment = vtxdist[1];
1086 int lprev = vtxdist[i]-vtxdist[i-1];
1087 int l = vtxdist[i+1]-vtxdist[i];
1089 assert((start+l+offset)-gxadj<=static_cast<idxtype>(gxadjlen));
1090 increment = *(start-1);
1091 std::transform(start+offset, start+l+offset, start, std::bind2nd(std::plus<idxtype>(), increment));
1093 Dune::dinfo<<std::endl<<
"shifted xadj:";
1095 Dune::dinfo<<std::endl<<
" gadjncy: ";
1098 Dune::dinfo<<std::endl<<
" gvwgt: ";
1100 Dune::dinfo<<std::endl<<
"adjwgt: ";
1102 Dune::dinfo<<std::endl;
1106 std::cout<<
"Postprocesing global graph data took "<<time1.elapsed()<<std::endl;
1110 gxadj, gadjncy,
true));
1114 std::cout<<
"Creating grah one 1 process took "<<time.elapsed()<<std::endl;
1116 options[0]=0; options[1]=1; options[2]=1; options[3]=3; options[4]=3;
1118 METIS_PartGraphRecursive(&noVertices, gxadj, gadjncy, gvwgt, gadjwgt, &wgtflag,
1119 &numflag, &nparts, options, &edgecut, gpart);
1122 std::cout<<
"METIS took "<<time.elapsed()<<std::endl;
1125 Dune::dinfo<<std::endl<<
"part:";
1136 MPI_Scatter(gpart, 1, MPITraits<idxtype>::getType(), part, 1,
1137 MPITraits<idxtype>::getType(), 0, comm);
1161 Dune::dinfo<<
" repart "<<rank <<
" -> "<< part[0]<<std::endl;
1163 std::vector<int> realpart(mat.N(), part[0]);
1169 std::cout<<
"Scattering repartitioning took "<<time.elapsed()<<std::endl;
1177 std::cout<<
"Filling index set took "<<time.elapsed()<<std::endl;
1185 std::cout<<
"Average no neighbours was "<<noNeighbours<<std::endl;
1190 std::cout<<
"Building index sets took "<<time.elapsed()<<std::endl;
1212 template<
class G,
class T1,
class T2>
1225 std::cout<<
"Filling holes took "<<time.elapsed()<<std::endl;
1232 double t1=0.0, t2=0.0, t3=0.0, t4=0.0, tSum=0.0;
1242 typedef std::vector<GI> GlobalVector;
1262 typedef typename OOComm::OwnerSet OwnerSet;
1267 ParmetisDuneIndexMap indexMap(graph,oocomm);
1271 std::size_t *part =
new std::size_t[indexMap.numOfOwnVtx()];
1273 for(std::size_t i=0; i < indexMap.numOfOwnVtx(); ++i)
1278 std::cerr<<
"ParMETIS not activated. Will repartition to 1 domain instead of requested "
1279 <<nparts<<
" domains."<<std::endl;
1287 idxtype *xadj =
new idxtype[indexMap.numOfOwnVtx()+1];
1288 idxtype *adjncy =
new idxtype[graph.noEdges()];
1289 EdgeFunctor<G> ef(adjncy, indexMap, graph.noEdges());
1290 getAdjArrays<OwnerSet>(graph, oocomm.
globalLookup(), xadj, ef);
1296 int numflag=0, wgtflag=0, options[3], edgecut=0, ncon=1;
1297 float *tpwgts = NULL;
1306 wgtflag = (ef.getWeights()!=NULL)?1:0;
1314 std::cout<<std::endl;
1315 std::cout<<
"Testing ParMETIS_V3_PartKway with options[1-2] = {"
1316 <<options[1]<<
" "<<options[2]<<
"}, Ncon: "
1317 <<ncon<<
", Nparts: "<<nparts<<std::endl;
1330 std::cout<<
"Preparing for parmetis took "<<time.elapsed()<<std::endl;
1337 ParMETIS_V3_PartKway(indexMap.vtxDist(), xadj, adjncy,
1338 NULL, ef.getWeights(), &wgtflag,
1339 &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgecut, part, &
const_cast<MPI_Comm&
>(comm));
1348 std::cout<<std::endl;
1349 std::cout<<
"ParMETIS_V3_PartKway reported a cut of "<<edgecut<<std::endl;
1350 std::cout<<std::endl;
1352 std::cout<<mype<<
": PARMETIS-Result: ";
1353 for(
int i=0; i < indexMap.vtxDist()[mype+1]-indexMap.vtxDist()[mype]; ++i) {
1354 std::cout<<part[i]<<
" ";
1356 std::cout<<std::endl;
1357 std::cout<<
"Testing ParMETIS_V3_PartKway with options[1-2] = {"
1358 <<options[1]<<
" "<<options[2]<<
"}, Ncon: "
1359 <<ncon<<
", Nparts: "<<nparts<<std::endl;
1372 std::cout<<
"Parmetis took "<<time.elapsed()<<std::endl;
1379 for(std::size_t i=0; i<indexMap.numOfOwnVtx();++i)
1389 int domainMapping[nparts];
1391 getDomain(comm, part, indexMap.numOfOwnVtx(), nparts, &myDomain, domainMapping);
1396 std::cout<<mype<<
": myDomain: "<<myDomain<<std::endl;
1397 std::cout<<mype<<
": DomainMapping: ";
1398 for(
int j=0; j<nparts; j++) {
1399 std::cout<<
" do: "<<j<<
" pe: "<<domainMapping[j]<<
" ";
1401 std::cout<<std::endl;
1408 std::vector<int> setPartition(oocomm.
indexSet().size(), -1);
1410 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
1413 if(OwnerSet::contains(index->local().attribute())){
1414 setPartition[index->local()]=domainMapping[part[i++]];
1422 static_cast<int>(SolverCategory::nonoverlapping))
1429 std::cout<<
"Creating indexsets took "<<time.elapsed()<<std::endl;
1436 template<
class G,
class T1,
class T2>
1444 typedef typename OOComm::OwnerSet OwnerSet;
1478 std::set<int> recvFrom;
1484 typedef typename std::vector<int>::const_iterator VIter;
1488 std::set<int> tsendTo;
1489 for(VIter i=setPartition.begin(), iend = setPartition.end(); i!=iend; ++i)
1492 noSendTo = tsendTo.size();
1493 sendTo =
new int[noSendTo];
1494 typedef std::set<int>::const_iterator iterator;
1496 for(iterator i=tsendTo.begin(); i != tsendTo.end(); ++i, ++idx)
1502 int* gsendToDispl =
new int[oocomm.
communicator().size()+1];
1504 MPI_Allgather(&noSendTo, 1, MPI_INT, gnoSend, 1,
1508 int totalNoRecv = 0;
1509 for(
int i=0; i<npes; ++i)
1510 totalNoRecv += gnoSend[i];
1512 int *gsendTo =
new int[totalNoRecv];
1516 for(
int i=0; i<npes; ++i)
1517 gsendToDispl[i+1]=gsendToDispl[i]+gnoSend[i];
1520 MPI_Allgatherv(sendTo, noSendTo, MPI_INT, gsendTo, gnoSend, gsendToDispl,
1524 for(
int proc=0; proc < npes; ++proc)
1525 for(
int i=gsendToDispl[proc]; i < gsendToDispl[proc+1]; ++i)
1526 if(gsendTo[i]==mype)
1527 recvFrom.insert(proc);
1529 bool existentOnNextLevel = recvFrom.size()>0;
1533 delete[] gsendToDispl;
1538 if(recvFrom.size()){
1539 std::cout<<mype<<
": recvFrom: ";
1540 typedef typename std::set<int>::const_iterator siter;
1541 for(siter i=recvFrom.begin(); i!= recvFrom.end(); ++i) {
1546 std::cout<<std::endl<<std::endl;
1547 std::cout<<mype<<
": sendTo: ";
1548 for(
int i=0; i<noSendTo; i++) {
1549 std::cout<<sendTo[i]<<
" ";
1551 std::cout<<std::endl<<std::endl;
1556 std::cout<<
" Communicating the receive information took "<<
1557 time.elapsed()<<std::endl;
1571 typedef typename OOComm::ParallelIndexSet::GlobalIndex GI;
1572 typedef std::vector<GI> GlobalVector;
1573 std::vector<std::pair<GI,int> > myOwnerVec;
1574 std::set<GI> myOverlapSet;
1575 GlobalVector sendOwnerVec;
1576 std::set<GI> sendOverlapSet;
1577 std::set<int> myNeighbors;
1582 char **sendBuffers=
new char*[noSendTo];
1583 MPI_Request *requests =
new MPI_Request[noSendTo];
1586 for(
int i=0; i < noSendTo; ++i){
1588 sendOwnerVec.clear();
1589 sendOverlapSet.clear();
1592 std::set<int> neighbors;
1593 getOwnerOverlapVec<OwnerSet>(graph, setPartition, oocomm.
globalLookup(),
1594 mype, sendTo[i], sendOwnerVec, sendOverlapSet, redistInf,
1600 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.
communicator(), &buffersize);
1601 MPI_Pack_size(sendOwnerVec.size(), MPITraits<GI>::getType(), oocomm.
communicator(), &tsize);
1603 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.
communicator(), &tsize);
1605 MPI_Pack_size(sendOverlapSet.size(), MPITraits<GI>::getType(), oocomm.
communicator(), &tsize);
1606 buffersize += tsize;
1607 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.
communicator(), &tsize);
1608 buffersize += tsize;
1609 MPI_Pack_size(neighbors.size(), MPI_INT, oocomm.
communicator(), &tsize);
1610 buffersize += tsize;
1612 sendBuffers[i] =
new char[buffersize];
1615 std::cout<<mype<<
" sending "<<sendOwnerVec.size()<<
" owner and "<<
1616 sendOverlapSet.size()<<
" overlap to "<<sendTo[i]<<
" buffersize="<<buffersize<<std::endl;
1618 createSendBuf(sendOwnerVec, sendOverlapSet, neighbors, sendBuffers[i], buffersize, oocomm.
communicator());
1619 MPI_Issend(sendBuffers[i], buffersize, MPI_PACKED, sendTo[i], 99, oocomm.
communicator(), requests+i);
1625 std::cout<<
" Creating sends took "<<
1626 time.elapsed()<<std::endl;
1631 int noRecv = recvFrom.size();
1632 int oldbuffersize=0;
1637 MPI_Probe(MPI_ANY_SOURCE, 99, oocomm.
communicator(), &stat);
1639 MPI_Get_count(&stat, MPI_PACKED, &buffersize);
1641 if(oldbuffersize<buffersize){
1644 recvBuf =
new char[buffersize];
1645 oldbuffersize = buffersize;
1647 MPI_Recv(recvBuf, buffersize, MPI_PACKED, stat.MPI_SOURCE, 99, oocomm.
communicator(), &stat);
1648 saveRecvBuf(recvBuf, buffersize, myOwnerVec, myOverlapSet, myNeighbors, redistInf,
1658 MPI_Status *statuses =
new MPI_Status[noSendTo];
1659 int send = MPI_Waitall(noSendTo, requests, statuses);
1662 if(send==MPI_ERR_IN_STATUS){
1663 std::cerr<<mype<<
": Error in sending :"<<std::endl;
1665 for(
int i=0; i< noSendTo; i++)
1666 if(statuses[i].MPI_ERROR!=MPI_SUCCESS){
1669 MPI_Error_string(statuses[i].MPI_ERROR, message, &messageLength);
1670 std::cerr<<
" source="<<statuses[i].MPI_SOURCE<<
" message: ";
1671 for(
int i=0; i< messageLength; i++)
1672 std::cout<<message[i];
1674 std::cerr<<std::endl;
1680 std::cout<<
" Receiving and saving took "<<
1681 time.elapsed()<<std::endl;
1685 for(
int i=0; i < noSendTo; ++i)
1686 delete[] sendBuffers[i];
1688 delete[] sendBuffers;
1703 if (!existentOnNextLevel) {
1705 color= MPI_UNDEFINED;
1707 MPI_Comm outputComm;
1715 std::vector<int> tneighbors;
1716 tneighbors.reserve(myNeighbors.size());
1718 typename OOComm::ParallelIndexSet& outputIndexSet = outcomm->indexSet();
1720 MPI_Allgather(&newrank, 1, MPI_INT, newranks, 1,
1722 typedef typename std::set<int>::const_iterator IIter;
1726 for(IIter i=myNeighbors.begin(), end=myNeighbors.end();
1728 assert(newranks[*i]>=0);
1729 std::cout<<*i<<
"->"<<newranks[*i]<<
" ";
1730 tneighbors.push_back(newranks[*i]);
1732 std::cout<<std::endl;
1734 for(IIter i=myNeighbors.begin(), end=myNeighbors.end();
1736 tneighbors.push_back(newranks[*i]);
1740 myNeighbors.clear();
1745 std::cout<<
" Calculating new neighbours ("<<tneighbors.size()<<
") took "<<
1746 time.elapsed()<<std::endl;
1751 outputIndexSet.beginResize();
1754 std::sort(myOwnerVec.begin(), myOwnerVec.end(), SortFirst());
1758 typedef typename OOComm::ParallelIndexSet::LocalIndex LocalIndex;
1759 typedef typename std::vector<std::pair<GI,int> >::const_iterator VPIter;
1761 for(VPIter g=myOwnerVec.begin(), end =myOwnerVec.end(); g!=end; ++g, ++i ) {
1762 outputIndexSet.add(g->first,LocalIndex(i, OwnerOverlapCopyAttributeSet::owner,
true));
1769 std::cout<<
" Adding owner indices took "<<
1770 time.elapsed()<<std::endl;
1779 mergeVec(myOwnerVec, myOverlapSet);
1783 myOwnerVec.swap(myOwnerVec);
1788 std::cout<<
" Merging indices took "<<
1789 time.elapsed()<<std::endl;
1795 typedef typename std::set<GI>::const_iterator SIter;
1796 for(SIter g=myOverlapSet.begin(), end=myOverlapSet.end(); g!=end; ++g, i++) {
1797 outputIndexSet.add(*g,LocalIndex(i, OwnerOverlapCopyAttributeSet::copy,
true));
1799 myOverlapSet.clear();
1800 outputIndexSet.endResize();
1802 #ifdef DUNE_ISTL_WITH_CHECKING
1804 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
1805 Iterator end = outputIndexSet.end();
1806 for(Iterator index = outputIndexSet.begin(); index != end; ++
index) {
1807 if (OwnerSet::contains(index->local().attribute())) {
1818 Iterator index=outputIndexSet.begin();
1821 for(Iterator old = outputIndexSet.begin(); index != end; old=index++) {
1822 if(old->global()>index->global())
1823 DUNE_THROW(
ISTLError,
"Index set's globalindex not sorted correctly");
1830 std::cout<<
" Adding overlap indices took "<<
1831 time.elapsed()<<std::endl;
1836 if(color != MPI_UNDEFINED){
1837 outcomm->remoteIndices().setNeighbours(tneighbors);
1838 outcomm->remoteIndices().template rebuild<true>();
1848 std::cout<<
" Storing indexsets took "<<
1849 time.elapsed()<<std::endl;
1855 tSum = t1 + t2 + t3 + t4;
1856 std::cout<<std::endl
1857 <<mype<<
": WTime for step 1): "<<t1
1865 return color!=MPI_UNDEFINED;
1869 template<
class G,
class P,
class T1,
class T2,
class R>
1875 if(nparts!=oocomm.size())
1876 DUNE_THROW(NotImplemented,
"only available for MPI programs");
1880 template<
class G,
class P,
class T1,
class T2,
class R>
1886 if(nparts!=oocomm.size())
1887 DUNE_THROW(NotImplemented,
"only available for MPI programs");