Dune Core Modules (2.7.0)

repartition.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_ISTL_REPARTITION_HH
4#define DUNE_ISTL_REPARTITION_HH
5
6#include <cassert>
7#include <map>
8#include <utility>
9#include <cmath>
10
11#if HAVE_PARMETIS
12// Explicitly use C linkage as scotch does not extern "C" in its headers.
13// Works because ParMETIS/METIS checks whether compiler is C++ and otherwise
14// does not use extern "C". Therfore no nested extern "C" will be created
15extern "C"
16{
17#include <parmetis.h>
18}
19#endif
20
21#include <dune/common/timer.hh>
22#include <dune/common/unused.hh>
31
34
43namespace Dune
44{
45 namespace Metis
46 {
47 // Explicitly specify a real_t and idx_t for older (Par)METIS versions that do not
48 // provide these typedefs
49#if HAVE_PARMETIS && defined(REALTYPEWIDTH)
50 using real_t = ::real_t;
51#else
52 using real_t = float;
53#endif
54
55#if HAVE_PARMETIS && defined(IDXTYPEWIDTH)
56 using idx_t = ::idx_t;
57#elif HAVE_PARMETIS && defined(HAVE_SCOTCH_NUM_TYPE)
58 using idx_t = SCOTCH_Num;
59#elif HAVE_PARMETIS
60 using idx_t = int;
61#else
62 using idx_t = std::size_t;
63#endif
64 }
65
66
67#if HAVE_MPI
81 template<class G, class T1, class T2>
83 {
85 typedef typename IndexSet::LocalIndex::Attribute Attribute;
86
87 IndexSet& indexSet = oocomm.indexSet();
88 const typename Dune::OwnerOverlapCopyCommunication<T1,T2>::GlobalLookupIndexSet& lookup =oocomm.globalLookup();
89
90 // The type of the const vertex iterator.
91 typedef typename G::ConstVertexIterator VertexIterator;
92
93
94 std::size_t sum=0, needed = graph.noVertices()-indexSet.size();
95 std::vector<std::size_t> neededall(oocomm.communicator().size(), 0);
96
97 MPI_Allgather(&needed, 1, MPITraits<std::size_t>::getType() , &(neededall[0]), 1, MPITraits<std::size_t>::getType(), oocomm.communicator());
98 for(int i=0; i<oocomm.communicator().size(); ++i)
99 sum=sum+neededall[i]; // MAke this for generic
100
101 if(sum==0)
102 // Nothing to do
103 return;
104
105 //Compute Maximum Global Index
106 T1 maxgi=0;
107 typedef typename IndexSet::const_iterator Iterator;
108 Iterator end;
109 end = indexSet.end();
110 for(Iterator it = indexSet.begin(); it != end; ++it)
111 maxgi=std::max(maxgi,it->global());
112
113 //Process p creates global indices consecutively
114 //starting atmaxgi+\sum_{i=1}^p neededall[i]
115 // All created indices are owned by the process
116 maxgi=oocomm.communicator().max(maxgi);
117 ++maxgi; //Sart with the next free index.
118
119 for(int i=0; i<oocomm.communicator().rank(); ++i)
120 maxgi=maxgi+neededall[i]; // TODO: make this more generic
121
122 // Store the global index information for repairing the remote index information
123 std::map<int,SLList<std::pair<T1,Attribute> > > globalIndices;
124 storeGlobalIndicesOfRemoteIndices(globalIndices, oocomm.remoteIndices());
125 indexSet.beginResize();
126
127 for(VertexIterator vertex = graph.begin(), vend=graph.end(); vertex != vend; ++vertex) {
128 const typename IndexSet::IndexPair* pair=lookup.pair(*vertex);
129 if(pair==0) {
130 // No index yet, add new one
131 indexSet.add(maxgi, typename IndexSet::LocalIndex(*vertex, OwnerOverlapCopyAttributeSet::owner, false));
132 ++maxgi;
133 }
134 }
135
136 indexSet.endResize();
137
138 repairLocalIndexPointers(globalIndices, oocomm.remoteIndices(), indexSet);
139
140 oocomm.freeGlobalLookup();
141 oocomm.buildGlobalLookup();
142#ifdef DEBUG_REPART
143 std::cout<<"Holes are filled!"<<std::endl;
144 std::cout<<oocomm.communicator().rank()<<": "<<oocomm.indexSet()<<std::endl;
145#endif
146 }
147
148 namespace
149 {
150
151 class ParmetisDuneIndexMap
152 {
153 public:
154 template<class Graph, class OOComm>
155 ParmetisDuneIndexMap(const Graph& graph, const OOComm& com);
156 int toParmetis(int i) const
157 {
158 return duneToParmetis[i];
159 }
160 int toLocalParmetis(int i) const
161 {
162 return duneToParmetis[i]-base_;
163 }
164 int operator[](int i) const
165 {
166 return duneToParmetis[i];
167 }
168 int toDune(int i) const
169 {
170 return parmetisToDune[i];
171 }
172 std::vector<int>::size_type numOfOwnVtx() const
173 {
174 return parmetisToDune.size();
175 }
176 Metis::idx_t* vtxDist()
177 {
178 return &vtxDist_[0];
179 }
180 int globalOwnerVertices;
181 private:
182 int base_;
183 std::vector<int> duneToParmetis;
184 std::vector<int> parmetisToDune;
185 // range of vertices for processor i: vtxdist[i] to vtxdist[i+1] (parmetis global)
186 std::vector<Metis::idx_t> vtxDist_;
187 };
188
189 template<class G, class OOComm>
190 ParmetisDuneIndexMap::ParmetisDuneIndexMap(const G& graph, const OOComm& oocomm)
191 : duneToParmetis(graph.noVertices(), -1), vtxDist_(oocomm.communicator().size()+1)
192 {
193 int npes=oocomm.communicator().size(), mype=oocomm.communicator().rank();
194
195 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
196 typedef typename OOComm::OwnerSet OwnerSet;
197
198 int numOfOwnVtx=0;
199 Iterator end = oocomm.indexSet().end();
200 for(Iterator index = oocomm.indexSet().begin(); index != end; ++index) {
201 if (OwnerSet::contains(index->local().attribute())) {
202 numOfOwnVtx++;
203 }
204 }
205 parmetisToDune.resize(numOfOwnVtx);
206 std::vector<int> globalNumOfVtx(npes);
207 // make this number available to all processes
208 MPI_Allgather(&numOfOwnVtx, 1, MPI_INT, &(globalNumOfVtx[0]), 1, MPI_INT, oocomm.communicator());
209
210 int base=0;
211 vtxDist_[0] = 0;
212 for(int i=0; i<npes; i++) {
213 if (i<mype) {
214 base += globalNumOfVtx[i];
215 }
216 vtxDist_[i+1] = vtxDist_[i] + globalNumOfVtx[i];
217 }
218 globalOwnerVertices=vtxDist_[npes];
219 base_=base;
220
221 // The type of the const vertex iterator.
222 typedef typename G::ConstVertexIterator VertexIterator;
223#ifdef DEBUG_REPART
224 std::cout << oocomm.communicator().rank()<<" vtxDist: ";
225 for(int i=0; i<= npes; ++i)
226 std::cout << vtxDist_[i]<<" ";
227 std::cout<<std::endl;
228#endif
229
230 // Traverse the graph and assign a new consecutive number/index
231 // starting by "base" to all owner vertices.
232 // The new index is used as the ParMETIS global index and is
233 // stored in the vector "duneToParmetis"
234 VertexIterator vend = graph.end();
235 for(VertexIterator vertex = graph.begin(); vertex != vend; ++vertex) {
236 const typename OOComm::ParallelIndexSet::IndexPair* index=oocomm.globalLookup().pair(*vertex);
237 assert(index);
238 if (OwnerSet::contains(index->local().attribute())) {
239 // assign and count the index
240 parmetisToDune[base-base_]=index->local();
241 duneToParmetis[index->local()] = base++;
242 }
243 }
244
245 // At this point, every process knows the ParMETIS global index
246 // of it's owner vertices. The next step is to get the
247 // ParMETIS global index of the overlap vertices from the
248 // associated processes. To do this, the Dune::Interface class
249 // is used.
250#ifdef DEBUG_REPART
251 std::cout <<oocomm.communicator().rank()<<": before ";
252 for(std::size_t i=0; i<duneToParmetis.size(); ++i)
253 std::cout<<duneToParmetis[i]<<" ";
254 std::cout<<std::endl;
255#endif
256 oocomm.copyOwnerToAll(duneToParmetis,duneToParmetis);
257#ifdef DEBUG_REPART
258 std::cout <<oocomm.communicator().rank()<<": after ";
259 for(std::size_t i=0; i<duneToParmetis.size(); ++i)
260 std::cout<<duneToParmetis[i]<<" ";
261 std::cout<<std::endl;
262#endif
263 }
264 }
265
266 struct RedistributeInterface
267 : public Interface
268 {
269 void setCommunicator(MPI_Comm comm)
270 {
271 communicator_=comm;
272 }
273 template<class Flags,class IS>
274 void buildSendInterface(const std::vector<int>& toPart, const IS& idxset)
275 {
276 std::map<int,int> sizes;
277
278 typedef typename IS::const_iterator IIter;
279 for(IIter i=idxset.begin(), end=idxset.end(); i!=end; ++i)
280 if(Flags::contains(i->local().attribute()))
281 ++sizes[toPart[i->local()]];
282
283 // Allocate the necessary space
284 typedef std::map<int,int>::const_iterator MIter;
285 for(MIter i=sizes.begin(), end=sizes.end(); i!=end; ++i)
286 interfaces()[i->first].first.reserve(i->second);
287
288 //Insert the interface information
289 typedef typename IS::const_iterator IIter;
290 for(IIter i=idxset.begin(), end=idxset.end(); i!=end; ++i)
291 if(Flags::contains(i->local().attribute()))
292 interfaces()[toPart[i->local()]].first.add(i->local());
293 }
294
295 void reserveSpaceForReceiveInterface(int proc, int size)
296 {
297 interfaces()[proc].second.reserve(size);
298 }
299 void addReceiveIndex(int proc, std::size_t idx)
300 {
301 interfaces()[proc].second.add(idx);
302 }
303 template<typename TG>
304 void buildReceiveInterface(std::vector<std::pair<TG,int> >& indices)
305 {
306 typedef typename std::vector<std::pair<TG,int> >::const_iterator VIter;
307 std::size_t i=0;
308 for(VIter idx=indices.begin(); idx!= indices.end(); ++idx) {
309 interfaces()[idx->second].second.add(i++);
310 }
311 }
312
313 ~RedistributeInterface()
314 {}
315
316 };
317
318 namespace
319 {
329 template<class GI>
330 void createSendBuf(std::vector<GI>& ownerVec, std::set<GI>& overlapVec, std::set<int>& neighbors, char *sendBuf, int buffersize, MPI_Comm comm) {
331 // Pack owner vertices
332 std::size_t s=ownerVec.size();
333 int pos=0;
334 if(s==0)
335 ownerVec.resize(1); // otherwise would read beyond the memory bound
336 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
337 MPI_Pack(&(ownerVec[0]), s, MPITraits<GI>::getType(), sendBuf, buffersize, &pos, comm);
338 s = overlapVec.size();
339 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
340 typedef typename std::set<GI>::iterator Iter;
341 for(Iter i=overlapVec.begin(), end= overlapVec.end(); i != end; ++i)
342 MPI_Pack(const_cast<GI*>(&(*i)), 1, MPITraits<GI>::getType(), sendBuf, buffersize, &pos, comm);
343
344 s=neighbors.size();
345 MPI_Pack(&s, 1, MPITraits<std::size_t>::getType(), sendBuf, buffersize, &pos, comm);
346 typedef typename std::set<int>::iterator IIter;
347
348 for(IIter i=neighbors.begin(), end= neighbors.end(); i != end; ++i)
349 MPI_Pack(const_cast<int*>(&(*i)), 1, MPI_INT, sendBuf, buffersize, &pos, comm);
350 }
359 template<class GI>
360 void saveRecvBuf(char *recvBuf, int bufferSize, std::vector<std::pair<GI,int> >& ownerVec,
361 std::set<GI>& overlapVec, std::set<int>& neighbors, RedistributeInterface& inf, int from, MPI_Comm comm) {
362 std::size_t size;
363 int pos=0;
364 // unpack owner vertices
365 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
366 inf.reserveSpaceForReceiveInterface(from, size);
367 ownerVec.reserve(ownerVec.size()+size);
368 for(; size!=0; --size) {
369 GI gi;
370 MPI_Unpack(recvBuf, bufferSize, &pos, &gi, 1, MPITraits<GI>::getType(), comm);
371 ownerVec.push_back(std::make_pair(gi,from));
372 }
373 // unpack overlap vertices
374 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
375 typename std::set<GI>::iterator ipos = overlapVec.begin();
376 Dune::dverb << "unpacking "<<size<<" overlap"<<std::endl;
377 for(; size!=0; --size) {
378 GI gi;
379 MPI_Unpack(recvBuf, bufferSize, &pos, &gi, 1, MPITraits<GI>::getType(), comm);
380 ipos=overlapVec.insert(ipos, gi);
381 }
382 //unpack neighbors
383 MPI_Unpack(recvBuf, bufferSize, &pos, &size, 1, MPITraits<std::size_t>::getType(), comm);
384 Dune::dverb << "unpacking "<<size<<" neighbors"<<std::endl;
385 typename std::set<int>::iterator npos = neighbors.begin();
386 for(; size!=0; --size) {
387 int n;
388 MPI_Unpack(recvBuf, bufferSize, &pos, &n, 1, MPI_INT, comm);
389 npos=neighbors.insert(npos, n);
390 }
391 }
392
406 template<typename T>
407 void getDomain(const MPI_Comm& comm, T *part, int numOfOwnVtx, int nparts, int *myDomain, std::vector<int> &domainMapping) {
408 int npes, mype;
409 MPI_Comm_size(comm, &npes);
410 MPI_Comm_rank(comm, &mype);
411 MPI_Status status;
412
413 *myDomain = -1;
414 int i=0;
415 int j=0;
416
417 std::vector<int> domain(nparts, 0);
418 std::vector<int> assigned(npes, 0);
419 // init domain Mapping
420 domainMapping.assign(domainMapping.size(), -1);
421
422 // count the occurrence of domains
423 for (i=0; i<numOfOwnVtx; i++) {
424 domain[part[i]]++;
425 }
426
427 std::vector<int> domainMatrix(npes * nparts, -1);
428
429 // init buffer with the own domain
430 int *buf = new int[nparts];
431 for (i=0; i<nparts; i++) {
432 buf[i] = domain[i];
433 domainMatrix[mype*nparts+i] = domain[i];
434 }
435 int pe=0;
436 int src = (mype-1+npes)%npes;
437 int dest = (mype+1)%npes;
438 // ring communication, we need n-1 communications for n processors
439 for (i=0; i<npes-1; i++) {
440 MPI_Sendrecv_replace(buf, nparts, MPI_INT, dest, 0, src, 0, comm, &status);
441 // pe is the process of the actual received buffer
442 pe = ((mype-1-i)+npes)%npes;
443 for(j=0; j<nparts; j++) {
444 // save the values to the domain matrix
445 domainMatrix[pe*nparts+j] = buf[j];
446 }
447 }
448 delete[] buf;
449
450 // Start the domain calculation.
451 // The process which contains the maximum number of vertices of a
452 // particular domain is selected to choose it's favorate domain
453 int maxOccurance = 0;
454 pe = -1;
455 std::set<std::size_t> unassigned;
456
457 for(i=0; i<nparts; i++) {
458 for(j=0; j<npes; j++) {
459 // process has no domain assigned
460 if (assigned[j]==0) {
461 if (maxOccurance < domainMatrix[j*nparts+i]) {
462 maxOccurance = domainMatrix[j*nparts+i];
463 pe = j;
464 }
465 }
466
467 }
468 if (pe!=-1) {
469 // process got a domain, ...
470 domainMapping[i] = pe;
471 // ...mark as assigned
472 assigned[pe] = 1;
473 if (pe==mype) {
474 *myDomain = i;
475 }
476 pe = -1;
477 }
478 else
479 {
480 unassigned.insert(i);
481 }
482 maxOccurance = 0;
483 }
484
485 typename std::vector<int>::iterator next_free = assigned.begin();
486
487 for(typename std::set<std::size_t>::iterator domain = unassigned.begin(),
488 end = unassigned.end(); domain != end; ++domain)
489 {
490 next_free = std::find_if(next_free, assigned.end(), std::bind(std::less<int>(), std::placeholders::_1, 1));
491 assert(next_free != assigned.end());
492 domainMapping[*domain] = next_free-assigned.begin();
493 *next_free = 1;
494 }
495 }
496
497 struct SortFirst
498 {
499 template<class T>
500 bool operator()(const T& t1, const T& t2) const
501 {
502 return t1<t2;
503 }
504 };
505
506
517 template<class GI>
518 void mergeVec(std::vector<std::pair<GI, int> >& ownerVec, std::set<GI>& overlapSet) {
519
520 typedef typename std::vector<std::pair<GI,int> >::const_iterator VIter;
521#ifdef DEBUG_REPART
522 // Safety check for duplicates.
523 if(ownerVec.size()>0)
524 {
525 VIter old=ownerVec.begin();
526 for(VIter i=old+1, end=ownerVec.end(); i != end; old=i++)
527 {
528 if(i->first==old->first)
529 {
530 std::cerr<<"Value at indes"<<old-ownerVec.begin()<<" is the same as at index "
531 <<i-ownerVec.begin()<<" ["<<old->first<<","<<old->second<<"]==["
532 <<i->first<<","<<i->second<<"]"<<std::endl;
533 throw "Huch!";
534 }
535 }
536 }
537
538#endif
539
540 typedef typename std::set<GI>::iterator SIter;
541 VIter v=ownerVec.begin(), vend=ownerVec.end();
542 for(SIter s=overlapSet.begin(), send=overlapSet.end(); s!=send;)
543 {
544 while(v!=vend && v->first<*s) ++v;
545 if(v!=vend && v->first==*s) {
546 // Move to the next element before erasing
547 // thus s stays valid!
548 SIter tmp=s;
549 ++s;
550 overlapSet.erase(tmp);
551 }else
552 ++s;
553 }
554 }
555
556
570 template<class OwnerSet, class Graph, class IS, class GI>
571 void getNeighbor(const Graph& g, std::vector<int>& part,
572 typename Graph::VertexDescriptor vtx, const IS& indexSet,
573 int toPe, std::set<GI>& neighbor, std::set<int>& neighborProcs) {
574 typedef typename Graph::ConstEdgeIterator Iter;
575 for(Iter edge=g.beginEdges(vtx), end=g.endEdges(vtx); edge!=end; ++edge)
576 {
577 const typename IS::IndexPair* pindex = indexSet.pair(edge.target());
578 assert(pindex);
579 if(part[pindex->local()]!=toPe || !OwnerSet::contains(pindex->local().attribute()))
580 {
581 // is sent to another process and therefore becomes overlap
582 neighbor.insert(pindex->global());
583 neighborProcs.insert(part[pindex->local()]);
584 }
585 }
586 }
587
588 template<class T, class I>
589 void my_push_back(std::vector<T>& ownerVec, const I& index, int proc)
590 {
592 ownerVec.push_back(index);
593 }
594
595 template<class T, class I>
596 void my_push_back(std::vector<std::pair<T,int> >& ownerVec, const I& index, int proc)
597 {
598 ownerVec.push_back(std::make_pair(index,proc));
599 }
600 template<class T>
601 void reserve(std::vector<T>&, RedistributeInterface&, int)
602 {}
603 template<class T>
604 void reserve(std::vector<std::pair<T,int> >& ownerVec, RedistributeInterface& redist, int proc)
605 {
606 redist.reserveSpaceForReceiveInterface(proc, ownerVec.size());
607 }
608
609
627 template<class OwnerSet, class G, class IS, class T, class GI>
628 void getOwnerOverlapVec(const G& graph, std::vector<int>& part, IS& indexSet,
629 int myPe, int toPe, std::vector<T>& ownerVec, std::set<GI>& overlapSet,
630 RedistributeInterface& redist, std::set<int>& neighborProcs) {
632 //typedef typename IndexSet::const_iterator Iterator;
633 typedef typename IS::const_iterator Iterator;
634 for(Iterator index = indexSet.begin(); index != indexSet.end(); ++index) {
635 // Only Process owner vertices, the others are not in the parmetis graph.
636 if(OwnerSet::contains(index->local().attribute()))
637 {
638 if(part[index->local()]==toPe)
639 {
640 getNeighbor<OwnerSet>(graph, part, index->local(), indexSet,
641 toPe, overlapSet, neighborProcs);
642 my_push_back(ownerVec, index->global(), toPe);
643 }
644 }
645 }
646 reserve(ownerVec, redist, toPe);
647
648 }
649
650
657 template<class F, class IS>
658 inline bool isOwner(IS& indexSet, int index) {
659
660 const typename IS::IndexPair* pindex=indexSet.pair(index);
661
662 assert(pindex);
663 return F::contains(pindex->local().attribute());
664 }
665
666
667 class BaseEdgeFunctor
668 {
669 public:
670 BaseEdgeFunctor(Metis::idx_t* adj,const ParmetisDuneIndexMap& data)
671 : i_(), adj_(adj), data_(data)
672 {}
673
674 template<class T>
675 void operator()(const T& edge)
676 {
677 // Get the egde weight
678 // const Weight& weight=edge.weight();
679 adj_[i_] = data_.toParmetis(edge.target());
680 i_++;
681 }
682 std::size_t index()
683 {
684 return i_;
685 }
686
687 private:
688 std::size_t i_;
689 Metis::idx_t* adj_;
690 const ParmetisDuneIndexMap& data_;
691 };
692
693 template<typename G>
694 struct EdgeFunctor
695 : public BaseEdgeFunctor
696 {
697 EdgeFunctor(Metis::idx_t* adj, const ParmetisDuneIndexMap& data, std::size_t)
698 : BaseEdgeFunctor(adj, data)
699 {}
700
701 Metis::idx_t* getWeights()
702 {
703 return NULL;
704 }
705 void free(){}
706 };
707
708 template<class G, class V, class E, class VM, class EM>
709 class EdgeFunctor<Dune::Amg::PropertiesGraph<G,V,E,VM,EM> >
710 : public BaseEdgeFunctor
711 {
712 public:
713 EdgeFunctor(Metis::idx_t* adj, const ParmetisDuneIndexMap& data, std::size_t s)
714 : BaseEdgeFunctor(adj, data)
715 {
716 weight_=new Metis::idx_t[s];
717 }
718
719 template<class T>
720 void operator()(const T& edge)
721 {
722 weight_[index()]=edge.properties().depends() ? 3 : 1;
723 BaseEdgeFunctor::operator()(edge);
724 }
725 Metis::idx_t* getWeights()
726 {
727 return weight_;
728 }
729 void free(){
730 if(weight_!=0) {
731 delete weight_;
732 weight_=0;
733 }
734 }
735 private:
736 Metis::idx_t* weight_;
737 };
738
739
740
754 template<class F, class G, class IS, class EW>
755 void getAdjArrays(G& graph, IS& indexSet, Metis::idx_t *xadj,
756 EW& ew)
757 {
758 int j=0;
759
760 // The type of the const vertex iterator.
761 typedef typename G::ConstVertexIterator VertexIterator;
762 //typedef typename IndexSet::const_iterator Iterator;
763 typedef typename IS::const_iterator Iterator;
764
765 VertexIterator vend = graph.end();
766 Iterator end;
767
768 for(VertexIterator vertex = graph.begin(); vertex != vend; ++vertex) {
769 if (isOwner<F>(indexSet,*vertex)) {
770 // The type of const edge iterator.
771 typedef typename G::ConstEdgeIterator EdgeIterator;
772 EdgeIterator eend = vertex.end();
773 xadj[j] = ew.index();
774 j++;
775 for(EdgeIterator edge = vertex.begin(); edge != eend; ++edge) {
776 ew(edge);
777 }
778 }
779 }
780 xadj[j] = ew.index();
781 }
782 } // end anonymous namespace
783
784 template<class G, class T1, class T2>
785 bool buildCommunication(const G& graph, std::vector<int>& realparts,
787 std::shared_ptr<Dune::OwnerOverlapCopyCommunication<T1,T2>>& outcomm,
788 RedistributeInterface& redistInf,
789 bool verbose=false);
790#if HAVE_PARMETIS
791#ifndef METIS_VER_MAJOR
792 extern "C"
793 {
794 // backwards compatibility to parmetis < 4.0.0
795 void METIS_PartGraphKway(int *nvtxs, Metis::idx_t *xadj, Metis::idx_t *adjncy, Metis::idx_t *vwgt,
796 Metis::idx_t *adjwgt, int *wgtflag, int *numflag, int *nparts,
797 int *options, int *edgecut, Metis::idx_t *part);
798
799 void METIS_PartGraphRecursive(int *nvtxs, Metis::idx_t *xadj, Metis::idx_t *adjncy, Metis::idx_t *vwgt,
800 Metis::idx_t *adjwgt, int *wgtflag, int *numflag, int *nparts,
801 int *options, int *edgecut, Metis::idx_t *part);
802 }
803#endif
804#endif // HAVE_PARMETIS
805
806 template<class S, class T>
807 inline void print_carray(S& os, T* array, std::size_t l)
808 {
809 for(T *cur=array, *end=array+l; cur!=end; ++cur)
810 os<<*cur<<" ";
811 }
812
813 template<class S, class T>
814 inline bool isValidGraph(std::size_t noVtx, std::size_t gnoVtx, S noEdges, T* xadj,
815 T* adjncy, bool checkSymmetry)
816 {
817 bool correct=true;
818
819 using std::signbit;
820 for(Metis::idx_t vtx=0; vtx<(Metis::idx_t)noVtx; ++vtx) {
821 if(xadj[vtx]>noEdges || signbit(xadj[vtx])) {
822 std::cerr <<"Check graph: xadj["<<vtx<<"]="<<xadj[vtx]<<" (>"
823 <<noEdges<<") out of range!"<<std::endl;
824 correct=false;
825 }
826 if(xadj[vtx+1]>noEdges || signbit(xadj[vtx+1])) {
827 std::cerr <<"Check graph: xadj["<<vtx+1<<"]="<<xadj[vtx+1]<<" (>"
828 <<noEdges<<") out of range!"<<std::endl;
829 correct=false;
830 }
831 // Check numbers in adjncy
832 for(Metis::idx_t i=xadj[vtx]; i< xadj[vtx+1]; ++i) {
833 if(adjncy[i]<0||((std::size_t)adjncy[i])>gnoVtx) {
834 std::cerr<<" Edge "<<adjncy[i]<<" out of range ["<<0<<","<<noVtx<<")"
835 <<std::endl;
836 correct=false;
837 }
838 }
839 if(checkSymmetry) {
840 for(Metis::idx_t i=xadj[vtx]; i< xadj[vtx+1]; ++i) {
841 Metis::idx_t target=adjncy[i];
842 // search for symmetric edge
843 int found=0;
844 for(Metis::idx_t j=xadj[target]; j< xadj[target+1]; ++j)
845 if(adjncy[j]==vtx)
846 found++;
847 if(found!=1) {
848 std::cerr<<"Edge ("<<target<<","<<vtx<<") "<<i<<" time"<<std::endl;
849 correct=false;
850 }
851 }
852 }
853 }
854 return correct;
855 }
856
857 template<class M, class T1, class T2>
858 bool commGraphRepartition(const M& mat, Dune::OwnerOverlapCopyCommunication<T1,T2>& oocomm,
859 Metis::idx_t nparts,
860 std::shared_ptr<Dune::OwnerOverlapCopyCommunication<T1,T2>>& outcomm,
861 RedistributeInterface& redistInf,
862 bool verbose=false)
863 {
864 if(verbose && oocomm.communicator().rank()==0)
865 std::cout<<"Repartitioning from "<<oocomm.communicator().size()
866 <<" to "<<nparts<<" parts"<<std::endl;
867 Timer time;
868 int rank = oocomm.communicator().rank();
869#if !HAVE_PARMETIS
870 int* part = new int[1];
871 part[0]=0;
872#else
873 Metis::idx_t* part = new Metis::idx_t[1]; // where all our data moves to
874
875 if(nparts>1) {
876
877 part[0]=rank;
878
879 { // sublock for automatic memory deletion
880
881 // Build the graph of the communication scheme and create an appropriate indexset.
882 // calculate the neighbour vertices
883 int noNeighbours = oocomm.remoteIndices().neighbours();
885 typedef typename RemoteIndices::const_iterator
886 NeighbourIterator;
887
888 for(NeighbourIterator n= oocomm.remoteIndices().begin(); n != oocomm.remoteIndices().end();
889 ++n)
890 if(n->first==rank) {
891 //do not include ourselves.
892 --noNeighbours;
893 break;
894 }
895
896 // A parmetis graph representing the communication graph.
897 // The diagonal entries are the number of nodes on the process.
898 // The offdiagonal entries are the number of edges leading to other processes.
899
900 Metis::idx_t *xadj=new Metis::idx_t[2];
901 Metis::idx_t *vtxdist=new Metis::idx_t[oocomm.communicator().size()+1];
902 Metis::idx_t *adjncy=new Metis::idx_t[noNeighbours];
903#ifdef USE_WEIGHTS
904 Metis::idx_t *vwgt = 0;
905 Metis::idx_t *adjwgt = 0;
906#endif
907
908 // each process has exactly one vertex!
909 for(int i=0; i<oocomm.communicator().size(); ++i)
910 vtxdist[i]=i;
911 vtxdist[oocomm.communicator().size()]=oocomm.communicator().size();
912
913 xadj[0]=0;
914 xadj[1]=noNeighbours;
915
916 // count edges to other processor
917 // a vector mapping the index to the owner
918 // std::vector<int> owner(mat.N(), oocomm.communicator().rank());
919 // for(NeighbourIterator n= oocomm.remoteIndices().begin(); n != oocomm.remoteIndices().end();
920 // ++n)
921 // {
922 // if(n->first!=oocomm.communicator().rank()){
923 // typedef typename RemoteIndices::RemoteIndexList RIList;
924 // const RIList& rlist = *(n->second.first);
925 // typedef typename RIList::const_iterator LIter;
926 // for(LIter entry=rlist.begin(); entry!=rlist.end(); ++entry){
927 // if(entry->attribute()==OwnerOverlapCopyAttributeSet::owner)
928 // owner[entry->localIndexPair().local()] = n->first;
929 // }
930 // }
931 // }
932
933 // std::map<int,Metis::idx_t> edgecount; // edges to other processors
934 // typedef typename M::ConstRowIterator RIter;
935 // typedef typename M::ConstColIterator CIter;
936
937 // // calculate edge count
938 // for(RIter row=mat.begin(), endr=mat.end(); row != endr; ++row)
939 // if(owner[row.index()]==OwnerOverlapCopyAttributeSet::owner)
940 // for(CIter entry= row->begin(), ende = row->end(); entry != ende; ++entry)
941 // ++edgecount[owner[entry.index()]];
942
943 // setup edge and weight pattern
944 typedef typename RemoteIndices::const_iterator NeighbourIterator;
945
946 Metis::idx_t* adjp=adjncy;
947
948#ifdef USE_WEIGHTS
949 vwgt = new Metis::idx_t[1];
950 vwgt[0]= mat.N(); // weight is numer of rows TODO: Should actually be the nonzeros.
951
952 adjwgt = new Metis::idx_t[noNeighbours];
953 Metis::idx_t* adjwp=adjwgt;
954#endif
955
956 for(NeighbourIterator n= oocomm.remoteIndices().begin(); n != oocomm.remoteIndices().end();
957 ++n)
958 if(n->first != rank) {
959 *adjp=n->first;
960 ++adjp;
961#ifdef USE_WEIGHTS
962 *adjwp=1; //edgecount[n->first];
963 ++adjwp;
964#endif
965 }
966 assert(isValidGraph(vtxdist[rank+1]-vtxdist[rank],
967 vtxdist[oocomm.communicator().size()],
968 noNeighbours, xadj, adjncy, false));
969
970 DUNE_UNUSED Metis::idx_t wgtflag=0;
971 Metis::idx_t numflag=0;
972 Metis::idx_t edgecut;
973#ifdef USE_WEIGHTS
974 wgtflag=3;
975#endif
976 Metis::real_t *tpwgts = new Metis::real_t[nparts];
977 for(int i=0; i<nparts; ++i)
978 tpwgts[i]=1.0/nparts;
979 MPI_Comm comm=oocomm.communicator();
980
981 Dune::dinfo<<rank<<" vtxdist: ";
982 print_carray(Dune::dinfo, vtxdist, oocomm.communicator().size()+1);
983 Dune::dinfo<<std::endl<<rank<<" xadj: ";
984 print_carray(Dune::dinfo, xadj, 2);
985 Dune::dinfo<<std::endl<<rank<<" adjncy: ";
986 print_carray(Dune::dinfo, adjncy, noNeighbours);
987
988#ifdef USE_WEIGHTS
989 Dune::dinfo<<std::endl<<rank<<" vwgt: ";
990 print_carray(Dune::dinfo, vwgt, 1);
991 Dune::dinfo<<std::endl<<rank<<" adwgt: ";
992 print_carray(Dune::dinfo, adjwgt, noNeighbours);
993#endif
994 Dune::dinfo<<std::endl;
995 oocomm.communicator().barrier();
996 if(verbose && oocomm.communicator().rank()==0)
997 std::cout<<"Creating comm graph took "<<time.elapsed()<<std::endl;
998 time.reset();
999
1000#ifdef PARALLEL_PARTITION
1001 Metis::real_t ubvec = 1.15;
1002 int ncon=1;
1003 int options[5] ={ 0,1,15,0,0};
1004
1005 //=======================================================
1006 // ParMETIS_V3_PartKway
1007 //=======================================================
1008 ParMETIS_V3_PartKway(vtxdist, xadj, adjncy,
1009 vwgt, adjwgt, &wgtflag,
1010 &numflag, &ncon, &nparts, tpwgts, &ubvec, options, &edgecut, part,
1011 &comm);
1012 if(verbose && oocomm.communicator().rank()==0)
1013 std::cout<<"ParMETIS took "<<time.elapsed()<<std::endl;
1014 time.reset();
1015#else
1016 Timer time1;
1017 std::size_t gnoedges=0;
1018 int* noedges = 0;
1019 noedges = new int[oocomm.communicator().size()];
1020 Dune::dverb<<"noNeighbours: "<<noNeighbours<<std::endl;
1021 // gather number of edges for each vertex.
1022 MPI_Allgather(&noNeighbours,1,MPI_INT,noedges,1, MPI_INT,oocomm.communicator());
1023
1024 if(verbose && oocomm.communicator().rank()==0)
1025 std::cout<<"Gathering noedges took "<<time1.elapsed()<<std::endl;
1026 time1.reset();
1027
1028 Metis::idx_t noVertices = vtxdist[oocomm.communicator().size()];
1029 Metis::idx_t *gxadj = 0;
1030 Metis::idx_t *gvwgt = 0;
1031 Metis::idx_t *gadjncy = 0;
1032 Metis::idx_t *gadjwgt = 0;
1033 Metis::idx_t *gpart = 0;
1034 int* displ = 0;
1035 int* noxs = 0;
1036 int* xdispl = 0; // displacement for xadj
1037 int* novs = 0;
1038 int* vdispl=0; // real vertex displacement
1039#ifdef USE_WEIGHTS
1040 std::size_t localNoVtx=vtxdist[rank+1]-vtxdist[rank];
1041#endif
1042 std::size_t gxadjlen = vtxdist[oocomm.communicator().size()]-vtxdist[0]+oocomm.communicator().size();
1043
1044 {
1045 Dune::dinfo<<"noedges: ";
1046 print_carray(Dune::dinfo, noedges, oocomm.communicator().size());
1047 Dune::dinfo<<std::endl;
1048 displ = new int[oocomm.communicator().size()];
1049 xdispl = new int[oocomm.communicator().size()];
1050 noxs = new int[oocomm.communicator().size()];
1051 vdispl = new int[oocomm.communicator().size()];
1052 novs = new int[oocomm.communicator().size()];
1053
1054 for(int i=0; i < oocomm.communicator().size(); ++i) {
1055 noxs[i]=vtxdist[i+1]-vtxdist[i]+1;
1056 novs[i]=vtxdist[i+1]-vtxdist[i];
1057 }
1058
1059 Metis::idx_t *so= vtxdist;
1060 int offset = 0;
1061 for(int *xcurr = xdispl, *vcurr = vdispl, *end=vdispl+oocomm.communicator().size();
1062 vcurr!=end; ++vcurr, ++xcurr, ++so, ++offset) {
1063 *vcurr = *so;
1064 *xcurr = offset + *so;
1065 }
1066
1067 int *pdispl =displ;
1068 int cdispl = 0;
1069 *pdispl = 0;
1070 for(int *curr=noedges, *end=noedges+oocomm.communicator().size()-1;
1071 curr!=end; ++curr) {
1072 ++pdispl; // next displacement
1073 cdispl += *curr; // next value
1074 *pdispl = cdispl;
1075 }
1076 Dune::dinfo<<"displ: ";
1077 print_carray(Dune::dinfo, displ, oocomm.communicator().size());
1078 Dune::dinfo<<std::endl;
1079
1080 // calculate global number of edges
1081 // It is bigger than the actual one as we habe size-1 additional end entries
1082 for(int *curr=noedges, *end=noedges+oocomm.communicator().size();
1083 curr!=end; ++curr)
1084 gnoedges += *curr;
1085
1086 // alocate gobal graph
1087 Dune::dinfo<<"gxadjlen: "<<gxadjlen<<" noVertices: "<<noVertices
1088 <<" gnoedges: "<<gnoedges<<std::endl;
1089 gxadj = new Metis::idx_t[gxadjlen];
1090 gpart = new Metis::idx_t[noVertices];
1091#ifdef USE_WEIGHTS
1092 gvwgt = new Metis::idx_t[noVertices];
1093 gadjwgt = new Metis::idx_t[gnoedges];
1094#endif
1095 gadjncy = new Metis::idx_t[gnoedges];
1096 }
1097
1098 if(verbose && oocomm.communicator().rank()==0)
1099 std::cout<<"Preparing global graph took "<<time1.elapsed()<<std::endl;
1100 time1.reset();
1101 // Communicate data
1102
1103 MPI_Allgatherv(xadj,2,MPITraits<Metis::idx_t>::getType(),
1104 gxadj,noxs,xdispl,MPITraits<Metis::idx_t>::getType(),
1105 comm);
1106 MPI_Allgatherv(adjncy,noNeighbours,MPITraits<Metis::idx_t>::getType(),
1107 gadjncy,noedges,displ,MPITraits<Metis::idx_t>::getType(),
1108 comm);
1109#ifdef USE_WEIGHTS
1110 MPI_Allgatherv(adjwgt,noNeighbours,MPITraits<Metis::idx_t>::getType(),
1111 gadjwgt,noedges,displ,MPITraits<Metis::idx_t>::getType(),
1112 comm);
1113 MPI_Allgatherv(vwgt,localNoVtx,MPITraits<Metis::idx_t>::getType(),
1114 gvwgt,novs,vdispl,MPITraits<Metis::idx_t>::getType(),
1115 comm);
1116#endif
1117 if(verbose && oocomm.communicator().rank()==0)
1118 std::cout<<"Gathering global graph data took "<<time1.elapsed()<<std::endl;
1119 time1.reset();
1120
1121 {
1122 // create the real gxadj array
1123 // i.e. shift entries and add displacements.
1124
1125 print_carray(Dune::dinfo, gxadj, gxadjlen);
1126
1127 int offset = 0;
1128 Metis::idx_t increment = vtxdist[1];
1129 Metis::idx_t *start=gxadj+1;
1130 for(int i=1; i<oocomm.communicator().size(); ++i) {
1131 offset+=1;
1132 int lprev = vtxdist[i]-vtxdist[i-1];
1133 int l = vtxdist[i+1]-vtxdist[i];
1134 start+=lprev;
1135 assert((start+l+offset)-gxadj<=static_cast<Metis::idx_t>(gxadjlen));
1136 increment = *(start-1);
1137 std::transform(start+offset, start+l+offset, start, std::bind(std::plus<Metis::idx_t>(), std::placeholders::_1, increment));
1138 }
1139 Dune::dinfo<<std::endl<<"shifted xadj:";
1140 print_carray(Dune::dinfo, gxadj, noVertices+1);
1141 Dune::dinfo<<std::endl<<" gadjncy: ";
1142 print_carray(Dune::dinfo, gadjncy, gnoedges);
1143#ifdef USE_WEIGHTS
1144 Dune::dinfo<<std::endl<<" gvwgt: ";
1145 print_carray(Dune::dinfo, gvwgt, noVertices);
1146 Dune::dinfo<<std::endl<<"adjwgt: ";
1147 print_carray(Dune::dinfo, gadjwgt, gnoedges);
1148 Dune::dinfo<<std::endl;
1149#endif
1150 // everything should be fine now!!!
1151 if(verbose && oocomm.communicator().rank()==0)
1152 std::cout<<"Postprocesing global graph data took "<<time1.elapsed()<<std::endl;
1153 time1.reset();
1154#ifndef NDEBUG
1155 assert(isValidGraph(noVertices, noVertices, gnoedges,
1156 gxadj, gadjncy, true));
1157#endif
1158
1159 if(verbose && oocomm.communicator().rank()==0)
1160 std::cout<<"Creating grah one 1 process took "<<time.elapsed()<<std::endl;
1161 time.reset();
1162#if METIS_VER_MAJOR >= 5
1163 Metis::idx_t ncon = 1;
1164 Metis::idx_t moptions[METIS_NOPTIONS];
1165 METIS_SetDefaultOptions(moptions);
1166 moptions[METIS_OPTION_NUMBERING] = numflag;
1167 METIS_PartGraphRecursive(&noVertices, &ncon, gxadj, gadjncy, gvwgt, NULL, gadjwgt,
1168 &nparts, NULL, NULL, moptions, &edgecut, gpart);
1169#else
1170 int options[5] = {0, 1, 1, 3, 3};
1171 // Call metis
1172 METIS_PartGraphRecursive(&noVertices, gxadj, gadjncy, gvwgt, gadjwgt, &wgtflag,
1173 &numflag, &nparts, options, &edgecut, gpart);
1174#endif
1175
1176 if(verbose && oocomm.communicator().rank()==0)
1177 std::cout<<"METIS took "<<time.elapsed()<<std::endl;
1178 time.reset();
1179
1180 Dune::dinfo<<std::endl<<"part:";
1181 print_carray(Dune::dinfo, gpart, noVertices);
1182
1183 delete[] gxadj;
1184 delete[] gadjncy;
1185#ifdef USE_WEIGHTS
1186 delete[] gvwgt;
1187 delete[] gadjwgt;
1188#endif
1189 }
1190 // Scatter result
1191 MPI_Scatter(gpart, 1, MPITraits<Metis::idx_t>::getType(), part, 1,
1192 MPITraits<Metis::idx_t>::getType(), 0, comm);
1193
1194 {
1195 // release remaining memory
1196 delete[] gpart;
1197 delete[] noedges;
1198 delete[] displ;
1199 }
1200
1201
1202#endif
1203 delete[] xadj;
1204 delete[] vtxdist;
1205 delete[] adjncy;
1206#ifdef USE_WEIGHTS
1207 delete[] vwgt;
1208 delete[] adjwgt;
1209#endif
1210 delete[] tpwgts;
1211 }
1212 }else{
1213 part[0]=0;
1214 }
1215#endif
1216 Dune::dinfo<<" repart "<<rank <<" -> "<< part[0]<<std::endl;
1217
1218 std::vector<int> realpart(mat.N(), part[0]);
1219 delete[] part;
1220
1221 oocomm.copyOwnerToAll(realpart, realpart);
1222
1223 if(verbose && oocomm.communicator().rank()==0)
1224 std::cout<<"Scattering repartitioning took "<<time.elapsed()<<std::endl;
1225 time.reset();
1226
1227
1228 oocomm.buildGlobalLookup(mat.N());
1229 Dune::Amg::MatrixGraph<M> graph(const_cast<M&>(mat));
1230 fillIndexSetHoles(graph, oocomm);
1231 if(verbose && oocomm.communicator().rank()==0)
1232 std::cout<<"Filling index set took "<<time.elapsed()<<std::endl;
1233 time.reset();
1234
1235 if(verbose) {
1236 int noNeighbours=oocomm.remoteIndices().neighbours();
1237 noNeighbours = oocomm.communicator().sum(noNeighbours)
1238 / oocomm.communicator().size();
1239 if(oocomm.communicator().rank()==0)
1240 std::cout<<"Average no neighbours was "<<noNeighbours<<std::endl;
1241 }
1242 bool ret = buildCommunication(graph, realpart, oocomm, outcomm, redistInf,
1243 verbose);
1244 if(verbose && oocomm.communicator().rank()==0)
1245 std::cout<<"Building index sets took "<<time.elapsed()<<std::endl;
1246 time.reset();
1247
1248
1249 return ret;
1250
1251 }
1252
1267 template<class G, class T1, class T2>
1268 bool graphRepartition(const G& graph, Dune::OwnerOverlapCopyCommunication<T1,T2>& oocomm, Metis::idx_t nparts,
1269 std::shared_ptr<Dune::OwnerOverlapCopyCommunication<T1,T2>>& outcomm,
1270 RedistributeInterface& redistInf,
1271 bool verbose=false)
1272 {
1273 Timer time;
1274
1275 MPI_Comm comm=oocomm.communicator();
1276 oocomm.buildGlobalLookup(graph.noVertices());
1277 fillIndexSetHoles(graph, oocomm);
1278
1279 if(verbose && oocomm.communicator().rank()==0)
1280 std::cout<<"Filling holes took "<<time.elapsed()<<std::endl;
1281 time.reset();
1282
1283 // simple precondition checks
1284
1285#ifdef PERF_REPART
1286 // Profiling variables
1287 double t1=0.0, t2=0.0, t3=0.0, t4=0.0, tSum=0.0;
1288#endif
1289
1290
1291 // MPI variables
1292 int mype = oocomm.communicator().rank();
1293
1294 assert(nparts<=static_cast<Metis::idx_t>(oocomm.communicator().size()));
1295
1296 int myDomain = -1;
1297
1298 //
1299 // 1) Prepare the required parameters for using ParMETIS
1300 // Especially the arrays that represent the graph must be
1301 // generated by the DUNE Graph and IndexSet input variables.
1302 // These are the arrays:
1303 // - vtxdist
1304 // - xadj
1305 // - adjncy
1306 //
1307 //
1308#ifdef PERF_REPART
1309 // reset timer for step 1)
1310 t1=MPI_Wtime();
1311#endif
1312
1313
1314 typedef typename Dune::OwnerOverlapCopyCommunication<T1,T2> OOComm;
1315 typedef typename OOComm::OwnerSet OwnerSet;
1316
1317 // Create the vtxdist array and parmetisVtxMapping.
1318 // Global communications are necessary
1319 // The parmetis global identifiers for the owner vertices.
1320 ParmetisDuneIndexMap indexMap(graph,oocomm);
1321 Metis::idx_t *part = new Metis::idx_t[indexMap.numOfOwnVtx()];
1322 for(std::size_t i=0; i < indexMap.numOfOwnVtx(); ++i)
1323 part[i]=mype;
1324
1325#if !HAVE_PARMETIS
1326 if(oocomm.communicator().rank()==0 && nparts>1)
1327 std::cerr<<"ParMETIS not activated. Will repartition to 1 domain instead of requested "
1328 <<nparts<<" domains."<<std::endl;
1329 nparts=1; // No parmetis available, fallback to agglomerating to 1 process
1330
1331#else
1332
1333 if(nparts>1) {
1334 // Create the xadj and adjncy arrays
1335 Metis::idx_t *xadj = new Metis::idx_t[indexMap.numOfOwnVtx()+1];
1336 Metis::idx_t *adjncy = new Metis::idx_t[graph.noEdges()];
1337 EdgeFunctor<G> ef(adjncy, indexMap, graph.noEdges());
1338 getAdjArrays<OwnerSet>(graph, oocomm.globalLookup(), xadj, ef);
1339
1340 //
1341 // 2) Call ParMETIS
1342 //
1343 //
1344 Metis::idx_t numflag=0, wgtflag=0, options[3], edgecut=0, ncon=1;
1345 //float *tpwgts = NULL;
1346 Metis::real_t *tpwgts = new Metis::real_t[nparts];
1347 for(int i=0; i<nparts; ++i)
1348 tpwgts[i]=1.0/nparts;
1349 Metis::real_t ubvec[1];
1350 options[0] = 0; // 0=default, 1=options are defined in [1]+[2]
1351#ifdef DEBUG_REPART
1352 options[1] = 3; // show info: 0=no message
1353#else
1354 options[1] = 0; // show info: 0=no message
1355#endif
1356 options[2] = 1; // random number seed, default is 15
1357 wgtflag = (ef.getWeights()!=NULL) ? 1 : 0;
1358 numflag = 0;
1359 edgecut = 0;
1360 ncon=1;
1361 ubvec[0]=1.05; // recommended by ParMETIS
1362
1363#ifdef DEBUG_REPART
1364 if (mype == 0) {
1365 std::cout<<std::endl;
1366 std::cout<<"Testing ParMETIS_V3_PartKway with options[1-2] = {"
1367 <<options[1]<<" "<<options[2]<<"}, Ncon: "
1368 <<ncon<<", Nparts: "<<nparts<<std::endl;
1369 }
1370#endif
1371#ifdef PERF_REPART
1372 // stop the time for step 1)
1373 t1=MPI_Wtime()-t1;
1374 // reset timer for step 2)
1375 t2=MPI_Wtime();
1376#endif
1377
1378 if(verbose) {
1379 oocomm.communicator().barrier();
1380 if(oocomm.communicator().rank()==0)
1381 std::cout<<"Preparing for parmetis took "<<time.elapsed()<<std::endl;
1382 }
1383 time.reset();
1384
1385 //=======================================================
1386 // ParMETIS_V3_PartKway
1387 //=======================================================
1388 ParMETIS_V3_PartKway(indexMap.vtxDist(), xadj, adjncy,
1389 NULL, ef.getWeights(), &wgtflag,
1390 &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgecut, part, &const_cast<MPI_Comm&>(comm));
1391
1392
1393 delete[] xadj;
1394 delete[] adjncy;
1395 delete[] tpwgts;
1396
1397 ef.free();
1398
1399#ifdef DEBUG_REPART
1400 if (mype == 0) {
1401 std::cout<<std::endl;
1402 std::cout<<"ParMETIS_V3_PartKway reported a cut of "<<edgecut<<std::endl;
1403 std::cout<<std::endl;
1404 }
1405 std::cout<<mype<<": PARMETIS-Result: ";
1406 for(int i=0; i < indexMap.vtxDist()[mype+1]-indexMap.vtxDist()[mype]; ++i) {
1407 std::cout<<part[i]<<" ";
1408 }
1409 std::cout<<std::endl;
1410 std::cout<<"Testing ParMETIS_V3_PartKway with options[1-2] = {"
1411 <<options[1]<<" "<<options[2]<<"}, Ncon: "
1412 <<ncon<<", Nparts: "<<nparts<<std::endl;
1413#endif
1414#ifdef PERF_REPART
1415 // stop the time for step 2)
1416 t2=MPI_Wtime()-t2;
1417 // reset timer for step 3)
1418 t3=MPI_Wtime();
1419#endif
1420
1421
1422 if(verbose) {
1423 oocomm.communicator().barrier();
1424 if(oocomm.communicator().rank()==0)
1425 std::cout<<"Parmetis took "<<time.elapsed()<<std::endl;
1426 }
1427 time.reset();
1428 }else
1429#endif
1430 {
1431 // Everything goes to process 0!
1432 for(std::size_t i=0; i<indexMap.numOfOwnVtx(); ++i)
1433 part[i]=0;
1434 }
1435
1436
1437 //
1438 // 3) Find a optimal domain based on the ParMETIS repartitioning
1439 // result
1440 //
1441
1442 std::vector<int> domainMapping(nparts);
1443 if(nparts>1)
1444 getDomain(comm, part, indexMap.numOfOwnVtx(), nparts, &myDomain, domainMapping);
1445 else
1446 domainMapping[0]=0;
1447
1448#ifdef DEBUG_REPART
1449 std::cout<<mype<<": myDomain: "<<myDomain<<std::endl;
1450 std::cout<<mype<<": DomainMapping: ";
1451 for(auto j : range(nparts)) {
1452 std::cout<<" do: "<<j<<" pe: "<<domainMapping[j]<<" ";
1453 }
1454 std::cout<<std::endl;
1455#endif
1456
1457 // Make a domain mapping for the indexset and translate
1458 //domain number to real process number
1459 // domainMapping is the one of parmetis, that is without
1460 // the overlap/copy vertices
1461 std::vector<int> setPartition(oocomm.indexSet().size(), -1);
1462
1463 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
1464 std::size_t i=0; // parmetis index
1465 for(Iterator index = oocomm.indexSet().begin(); index != oocomm.indexSet().end(); ++index)
1466 if(OwnerSet::contains(index->local().attribute())) {
1467 setPartition[index->local()]=domainMapping[part[i++]];
1468 }
1469
1470 delete[] part;
1471 oocomm.copyOwnerToAll(setPartition, setPartition);
1472 // communication only needed for ALU
1473 // (ghosts with same global id as owners on the same process)
1474 if (SolverCategory::category(oocomm) ==
1475 static_cast<int>(SolverCategory::nonoverlapping))
1476 oocomm.copyCopyToAll(setPartition, setPartition);
1477 bool ret = buildCommunication(graph, setPartition, oocomm, outcomm, redistInf,
1478 verbose);
1479 if(verbose) {
1480 oocomm.communicator().barrier();
1481 if(oocomm.communicator().rank()==0)
1482 std::cout<<"Creating indexsets took "<<time.elapsed()<<std::endl;
1483 }
1484 return ret;
1485 }
1486
1487
1488
1489 template<class G, class T1, class T2>
1490 bool buildCommunication(const G& graph,
1491 std::vector<int>& setPartition, Dune::OwnerOverlapCopyCommunication<T1,T2>& oocomm,
1492 std::shared_ptr<Dune::OwnerOverlapCopyCommunication<T1,T2>>& outcomm,
1493 RedistributeInterface& redistInf,
1494 bool verbose)
1495 {
1496 typedef typename Dune::OwnerOverlapCopyCommunication<T1,T2> OOComm;
1497 typedef typename OOComm::OwnerSet OwnerSet;
1498
1499 Timer time;
1500
1501 // Build the send interface
1502 redistInf.buildSendInterface<OwnerSet>(setPartition, oocomm.indexSet());
1503
1504#ifdef PERF_REPART
1505 // stop the time for step 3)
1506 t3=MPI_Wtime()-t3;
1507 // reset timer for step 4)
1508 t4=MPI_Wtime();
1509#endif
1510
1511
1512 //
1513 // 4) Create the output IndexSet and RemoteIndices
1514 // 4.1) Determine the "send to" and "receive from" relation
1515 // according to the new partition using a MPI ring
1516 // communication.
1517 //
1518 // 4.2) Depends on the "send to" and "receive from" vector,
1519 // the processes will exchange the vertices each other
1520 //
1521 // 4.3) Create the IndexSet, RemoteIndices and the new MPI
1522 // communicator
1523 //
1524
1525 //
1526 // 4.1) Let's start...
1527 //
1528 int npes = oocomm.communicator().size();
1529 int *sendTo = 0;
1530 int noSendTo = 0;
1531 std::set<int> recvFrom;
1532
1533 // the max number of vertices is stored in the sendTo buffer,
1534 // not the number of vertices to send! Because the max number of Vtx
1535 // is used as the fixed buffer size by the MPI send/receive calls
1536
1537 typedef typename std::vector<int>::const_iterator VIter;
1538 int mype = oocomm.communicator().rank();
1539
1540 {
1541 std::set<int> tsendTo;
1542 for(VIter i=setPartition.begin(), iend = setPartition.end(); i!=iend; ++i)
1543 tsendTo.insert(*i);
1544
1545 noSendTo = tsendTo.size();
1546 sendTo = new int[noSendTo];
1547 typedef std::set<int>::const_iterator iterator;
1548 int idx=0;
1549 for(iterator i=tsendTo.begin(); i != tsendTo.end(); ++i, ++idx)
1550 sendTo[idx]=*i;
1551 }
1552
1553 //
1554 int* gnoSend= new int[oocomm.communicator().size()];
1555 int* gsendToDispl = new int[oocomm.communicator().size()+1];
1556
1557 MPI_Allgather(&noSendTo, 1, MPI_INT, gnoSend, 1,
1558 MPI_INT, oocomm.communicator());
1559
1560 // calculate total receive message size
1561 int totalNoRecv = 0;
1562 for(int i=0; i<npes; ++i)
1563 totalNoRecv += gnoSend[i];
1564
1565 int *gsendTo = new int[totalNoRecv];
1566
1567 // calculate displacement for allgatherv
1568 gsendToDispl[0]=0;
1569 for(int i=0; i<npes; ++i)
1570 gsendToDispl[i+1]=gsendToDispl[i]+gnoSend[i];
1571
1572 // gather the data
1573 MPI_Allgatherv(sendTo, noSendTo, MPI_INT, gsendTo, gnoSend, gsendToDispl,
1574 MPI_INT, oocomm.communicator());
1575
1576 // Extract from which processes we will receive data
1577 for(int proc=0; proc < npes; ++proc)
1578 for(int i=gsendToDispl[proc]; i < gsendToDispl[proc+1]; ++i)
1579 if(gsendTo[i]==mype)
1580 recvFrom.insert(proc);
1581
1582 bool existentOnNextLevel = recvFrom.size()>0;
1583
1584 // Delete memory
1585 delete[] gnoSend;
1586 delete[] gsendToDispl;
1587 delete[] gsendTo;
1588
1589
1590#ifdef DEBUG_REPART
1591 if(recvFrom.size()) {
1592 std::cout<<mype<<": recvFrom: ";
1593 typedef typename std::set<int>::const_iterator siter;
1594 for(siter i=recvFrom.begin(); i!= recvFrom.end(); ++i) {
1595 std::cout<<*i<<" ";
1596 }
1597 }
1598
1599 std::cout<<std::endl<<std::endl;
1600 std::cout<<mype<<": sendTo: ";
1601 for(int i=0; i<noSendTo; i++) {
1602 std::cout<<sendTo[i]<<" ";
1603 }
1604 std::cout<<std::endl<<std::endl;
1605#endif
1606
1607 if(verbose)
1608 if(oocomm.communicator().rank()==0)
1609 std::cout<<" Communicating the receive information took "<<
1610 time.elapsed()<<std::endl;
1611 time.reset();
1612
1613 //
1614 // 4.2) Start the communication
1615 //
1616
1617 // Get all the owner and overlap vertices for myself ans save
1618 // it in the vectors myOwnerVec and myOverlapVec.
1619 // The received vertices from the other processes are simple
1620 // added to these vector.
1621 //
1622
1623
1624 typedef typename OOComm::ParallelIndexSet::GlobalIndex GI;
1625 typedef std::vector<GI> GlobalVector;
1626 std::vector<std::pair<GI,int> > myOwnerVec;
1627 std::set<GI> myOverlapSet;
1628 GlobalVector sendOwnerVec;
1629 std::set<GI> sendOverlapSet;
1630 std::set<int> myNeighbors;
1631
1632 // getOwnerOverlapVec<OwnerSet>(graph, setPartition, oocomm.globalLookup(),
1633 // mype, mype, myOwnerVec, myOverlapSet, redistInf, myNeighbors);
1634
1635 char **sendBuffers=new char*[noSendTo];
1636 MPI_Request *requests = new MPI_Request[noSendTo];
1637
1638 // Create all messages to be sent
1639 for(int i=0; i < noSendTo; ++i) {
1640 // clear the vector for sending
1641 sendOwnerVec.clear();
1642 sendOverlapSet.clear();
1643 // get all owner and overlap vertices for process j and save these
1644 // in the vectors sendOwnerVec and sendOverlapSet
1645 std::set<int> neighbors;
1646 getOwnerOverlapVec<OwnerSet>(graph, setPartition, oocomm.globalLookup(),
1647 mype, sendTo[i], sendOwnerVec, sendOverlapSet, redistInf,
1648 neighbors);
1649 // +2, we need 2 integer more for the length of each part
1650 // (owner/overlap) of the array
1651 int buffersize=0;
1652 int tsize;
1653 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.communicator(), &buffersize);
1654 MPI_Pack_size(sendOwnerVec.size(), MPITraits<GI>::getType(), oocomm.communicator(), &tsize);
1655 buffersize +=tsize;
1656 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.communicator(), &tsize);
1657 buffersize +=tsize;
1658 MPI_Pack_size(sendOverlapSet.size(), MPITraits<GI>::getType(), oocomm.communicator(), &tsize);
1659 buffersize += tsize;
1660 MPI_Pack_size(1, MPITraits<std::size_t>::getType(), oocomm.communicator(), &tsize);
1661 buffersize += tsize;
1662 MPI_Pack_size(neighbors.size(), MPI_INT, oocomm.communicator(), &tsize);
1663 buffersize += tsize;
1664
1665 sendBuffers[i] = new char[buffersize];
1666
1667#ifdef DEBUG_REPART
1668 std::cout<<mype<<" sending "<<sendOwnerVec.size()<<" owner and "<<
1669 sendOverlapSet.size()<<" overlap to "<<sendTo[i]<<" buffersize="<<buffersize<<std::endl;
1670#endif
1671 createSendBuf(sendOwnerVec, sendOverlapSet, neighbors, sendBuffers[i], buffersize, oocomm.communicator());
1672 MPI_Issend(sendBuffers[i], buffersize, MPI_PACKED, sendTo[i], 99, oocomm.communicator(), requests+i);
1673 }
1674
1675 if(verbose) {
1676 oocomm.communicator().barrier();
1677 if(oocomm.communicator().rank()==0)
1678 std::cout<<" Creating sends took "<<
1679 time.elapsed()<<std::endl;
1680 }
1681 time.reset();
1682
1683 // Receive Messages
1684 int noRecv = recvFrom.size();
1685 int oldbuffersize=0;
1686 char* recvBuf = 0;
1687 while(noRecv>0) {
1688 // probe for an incoming message
1689 MPI_Status stat;
1690 MPI_Probe(MPI_ANY_SOURCE, 99, oocomm.communicator(), &stat);
1691 int buffersize;
1692 MPI_Get_count(&stat, MPI_PACKED, &buffersize);
1693
1694 if(oldbuffersize<buffersize) {
1695 // buffer too small, reallocate
1696 delete[] recvBuf;
1697 recvBuf = new char[buffersize];
1698 oldbuffersize = buffersize;
1699 }
1700 MPI_Recv(recvBuf, buffersize, MPI_PACKED, stat.MPI_SOURCE, 99, oocomm.communicator(), &stat);
1701 saveRecvBuf(recvBuf, buffersize, myOwnerVec, myOverlapSet, myNeighbors, redistInf,
1702 stat.MPI_SOURCE, oocomm.communicator());
1703 --noRecv;
1704 }
1705
1706 if(recvBuf)
1707 delete[] recvBuf;
1708
1709 time.reset();
1710 // Wait for sending messages to complete
1711 MPI_Status *statuses = new MPI_Status[noSendTo];
1712 int send = MPI_Waitall(noSendTo, requests, statuses);
1713
1714 // check for errors
1715 if(send==MPI_ERR_IN_STATUS) {
1716 std::cerr<<mype<<": Error in sending :"<<std::endl;
1717 // Search for the error
1718 for(int i=0; i< noSendTo; i++)
1719 if(statuses[i].MPI_ERROR!=MPI_SUCCESS) {
1720 char message[300];
1721 int messageLength;
1722 MPI_Error_string(statuses[i].MPI_ERROR, message, &messageLength);
1723 std::cerr<<" source="<<statuses[i].MPI_SOURCE<<" message: ";
1724 for(int j = 0; j < messageLength; j++)
1725 std::cout<<message[j];
1726 }
1727 std::cerr<<std::endl;
1728 }
1729
1730 if(verbose) {
1731 oocomm.communicator().barrier();
1732 if(oocomm.communicator().rank()==0)
1733 std::cout<<" Receiving and saving took "<<
1734 time.elapsed()<<std::endl;
1735 }
1736 time.reset();
1737
1738 for(int i=0; i < noSendTo; ++i)
1739 delete[] sendBuffers[i];
1740
1741 delete[] sendBuffers;
1742 delete[] statuses;
1743 delete[] requests;
1744
1745 redistInf.setCommunicator(oocomm.communicator());
1746
1747 //
1748 // 4.2) Create the IndexSet etc.
1749 //
1750
1751 // build the new outputIndexSet
1752
1753
1754 int color=0;
1755
1756 if (!existentOnNextLevel) {
1757 // this process is not used anymore
1758 color= MPI_UNDEFINED;
1759 }
1760 MPI_Comm outputComm;
1761
1762 MPI_Comm_split(oocomm.communicator(), color, oocomm.communicator().rank(), &outputComm);
1763 outcomm = std::make_shared<OOComm>(outputComm,SolverCategory::category(oocomm),true);
1764
1765 // translate neighbor ranks.
1766 int newrank=outcomm->communicator().rank();
1767 int *newranks=new int[oocomm.communicator().size()];
1768 std::vector<int> tneighbors;
1769 tneighbors.reserve(myNeighbors.size());
1770
1771 typename OOComm::ParallelIndexSet& outputIndexSet = outcomm->indexSet();
1772
1773 MPI_Allgather(&newrank, 1, MPI_INT, newranks, 1,
1774 MPI_INT, oocomm.communicator());
1775 typedef typename std::set<int>::const_iterator IIter;
1776
1777#ifdef DEBUG_REPART
1778 std::cout<<oocomm.communicator().rank()<<" ";
1779 for(IIter i=myNeighbors.begin(), end=myNeighbors.end();
1780 i!=end; ++i) {
1781 assert(newranks[*i]>=0);
1782 std::cout<<*i<<"->"<<newranks[*i]<<" ";
1783 tneighbors.push_back(newranks[*i]);
1784 }
1785 std::cout<<std::endl;
1786#else
1787 for(IIter i=myNeighbors.begin(), end=myNeighbors.end();
1788 i!=end; ++i) {
1789 tneighbors.push_back(newranks[*i]);
1790 }
1791#endif
1792 delete[] newranks;
1793 myNeighbors.clear();
1794
1795 if(verbose) {
1796 oocomm.communicator().barrier();
1797 if(oocomm.communicator().rank()==0)
1798 std::cout<<" Calculating new neighbours ("<<tneighbors.size()<<") took "<<
1799 time.elapsed()<<std::endl;
1800 }
1801 time.reset();
1802
1803
1804 outputIndexSet.beginResize();
1805 // 1) add the owner vertices
1806 // Sort the owners
1807 std::sort(myOwnerVec.begin(), myOwnerVec.end(), SortFirst());
1808 // The owners are sorted according to there global index
1809 // Therefore the entries of ownerVec are the same as the
1810 // ones in the resulting index set.
1811 typedef typename OOComm::ParallelIndexSet::LocalIndex LocalIndex;
1812 typedef typename std::vector<std::pair<GI,int> >::const_iterator VPIter;
1813 int i=0;
1814 for(VPIter g=myOwnerVec.begin(), end =myOwnerVec.end(); g!=end; ++g, ++i ) {
1815 outputIndexSet.add(g->first,LocalIndex(i, OwnerOverlapCopyAttributeSet::owner, true));
1816 redistInf.addReceiveIndex(g->second, i);
1817 }
1818
1819 if(verbose) {
1820 oocomm.communicator().barrier();
1821 if(oocomm.communicator().rank()==0)
1822 std::cout<<" Adding owner indices took "<<
1823 time.elapsed()<<std::endl;
1824 }
1825 time.reset();
1826
1827
1828 // After all the vertices are received, the vectors must
1829 // be "merged" together to create the final vectors.
1830 // Because some vertices that are sent as overlap could now
1831 // already included as owner vertiecs in the new partition
1832 mergeVec(myOwnerVec, myOverlapSet);
1833
1834 // Trick to free memory
1835 myOwnerVec.clear();
1836 myOwnerVec.swap(myOwnerVec);
1837
1838 if(verbose) {
1839 oocomm.communicator().barrier();
1840 if(oocomm.communicator().rank()==0)
1841 std::cout<<" Merging indices took "<<
1842 time.elapsed()<<std::endl;
1843 }
1844 time.reset();
1845
1846
1847 // 2) add the overlap vertices
1848 typedef typename std::set<GI>::const_iterator SIter;
1849 for(SIter g=myOverlapSet.begin(), end=myOverlapSet.end(); g!=end; ++g, i++) {
1850 outputIndexSet.add(*g,LocalIndex(i, OwnerOverlapCopyAttributeSet::copy, true));
1851 }
1852 myOverlapSet.clear();
1853 outputIndexSet.endResize();
1854
1855#ifdef DUNE_ISTL_WITH_CHECKING
1856 int numOfOwnVtx =0;
1857 typedef typename OOComm::ParallelIndexSet::const_iterator Iterator;
1858 Iterator end = outputIndexSet.end();
1859 for(Iterator index = outputIndexSet.begin(); index != end; ++index) {
1860 if (OwnerSet::contains(index->local().attribute())) {
1861 numOfOwnVtx++;
1862 }
1863 }
1864 numOfOwnVtx = oocomm.communicator().sum(numOfOwnVtx);
1865 // if(numOfOwnVtx!=indexMap.globalOwnerVertices)
1866 // {
1867 // std::cerr<<numOfOwnVtx<<"!="<<indexMap.globalOwnerVertices<<" owners missing or additional ones!"<<std::endl;
1868 // DUNE_THROW(ISTLError, numOfOwnVtx<<"!="<<indexMap.globalOwnerVertices<<" owners missing or additional ones"
1869 // <<" during repartitioning.");
1870 // }
1871 Iterator index=outputIndexSet.begin();
1872 if(index!=end) {
1873 ++index;
1874 for(Iterator old = outputIndexSet.begin(); index != end; old=index++) {
1875 if(old->global()>index->global())
1876 DUNE_THROW(ISTLError, "Index set's globalindex not sorted correctly");
1877 }
1878 }
1879#endif
1880 if(verbose) {
1881 oocomm.communicator().barrier();
1882 if(oocomm.communicator().rank()==0)
1883 std::cout<<" Adding overlap indices took "<<
1884 time.elapsed()<<std::endl;
1885 }
1886 time.reset();
1887
1888
1889 if(color != MPI_UNDEFINED) {
1890 outcomm->remoteIndices().setNeighbours(tneighbors);
1891 outcomm->remoteIndices().template rebuild<true>();
1892
1893 }
1894
1895 // release the memory
1896 delete[] sendTo;
1897
1898 if(verbose) {
1899 oocomm.communicator().barrier();
1900 if(oocomm.communicator().rank()==0)
1901 std::cout<<" Storing indexsets took "<<
1902 time.elapsed()<<std::endl;
1903 }
1904
1905#ifdef PERF_REPART
1906 // stop the time for step 4) and print the results
1907 t4=MPI_Wtime()-t4;
1908 tSum = t1 + t2 + t3 + t4;
1909 std::cout<<std::endl
1910 <<mype<<": WTime for step 1): "<<t1
1911 <<" 2): "<<t2
1912 <<" 3): "<<t3
1913 <<" 4): "<<t4
1914 <<" total: "<<tSum
1915 <<std::endl;
1916#endif
1917
1918 return color!=MPI_UNDEFINED;
1919
1920 }
1921#else
1922 template<class G, class P,class T1, class T2, class R>
1923 bool graphRepartition(const G& graph, P& oocomm, int nparts,
1924 std::shared_ptr<P>& outcomm,
1925 R& redistInf,
1926 bool v=false)
1927 {
1928 if(nparts!=oocomm.size())
1929 DUNE_THROW(NotImplemented, "only available for MPI programs");
1930 }
1931
1932
1933 template<class G, class P,class T1, class T2, class R>
1934 bool commGraphRepartition(const G& graph, P& oocomm, int nparts,
1935 std::shared_ptr<P>& outcomm,
1936 R& redistInf,
1937 bool v=false)
1938 {
1939 if(nparts!=oocomm.size())
1940 DUNE_THROW(NotImplemented, "only available for MPI programs");
1941 }
1942#endif // HAVE_MPI
1943} // end of namespace Dune
1944#endif
The (undirected) graph of a matrix.
Definition: graph.hh:49
Index Set Interface base class.
Definition: indexidset.hh:76
IndexType size(GeometryType type) const
Return total number of entities of given geometry type in entity set .
Definition: indexidset.hh:220
MPI_Comm communicator_
The MPI communicator we use.
Definition: interface.hh:314
size_type N() const
Return the number of rows.
Definition: matrix.hh:691
A class setting up standard communication for a two-valued attribute set with owner/overlap/copy sema...
Definition: owneroverlapcopy.hh:172
const ParallelIndexSet & indexSet() const
Get the underlying parallel index set.
Definition: owneroverlapcopy.hh:472
void copyCopyToAll(const T &source, T &dest) const
Communicate values from copy data points to all other data points.
Definition: owneroverlapcopy.hh:340
void copyOwnerToAll(const T &source, T &dest) const
Communicate values from owner data points to all other data points.
Definition: owneroverlapcopy.hh:323
const RemoteIndices & remoteIndices() const
Get the underlying remote indices.
Definition: owneroverlapcopy.hh:481
LocalIndex::Attribute Attribute
The type of the attribute.
Definition: remoteindices.hh:224
const_iterator end() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1535
int neighbours() const
Get the number of processors we share indices with.
Definition: remoteindices.hh:1452
const_iterator begin() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1528
A simple stop watch.
Definition: timer.hh:41
void reset() noexcept
Reset timer while keeping the running/stopped state.
Definition: timer.hh:55
double elapsed() const noexcept
Get elapsed user-time from last reset until now/last stop in seconds.
Definition: timer.hh:75
Provides utility classes for syncing distributed data via MPI communication.
Classes for building sets out of enumeration values.
Provides classes for building the matrix graph.
size_t size() const
Get the total number (public and nonpublic) indices.
void repairLocalIndexPointers(std::map< int, SLList< std::pair< typename T::GlobalIndex, typename T::LocalIndex::Attribute >, A > > &globalMap, RemoteIndices< T, A1 > &remoteIndices, const T &indexSet)
Repair the pointers to the local indices in the remote indices.
Definition: indicessyncer.hh:490
iterator begin()
Get an iterator over the indices positioned at the first index.
iterator end()
Get an iterator over the indices positioned after the last index.
const InformationMap & interfaces() const
Get information about the interfaces.
Definition: interface.hh:422
const IndexPair * pair(const std::size_t &local) const
Get the index pair corresponding to a local index.
void storeGlobalIndicesOfRemoteIndices(std::map< int, SLList< std::pair< typename T::GlobalIndex, typename T::LocalIndex::Attribute >, A > > &globalMap, const RemoteIndices< T, A1 > &remoteIndices)
Stores the corresponding global indices of the remote index information.
Definition: indicessyncer.hh:461
constexpr index_constant< 1 > _1
Compile time index with value 1.
Definition: indices.hh:54
typename FieldTraits< Type >::real_type real_t
Convenient access to FieldTraits<Type>::real_type.
Definition: typetraits.hh:382
#define DUNE_UNUSED
A macro for marking variables that the compiler mistakenly flags as unused, which sometimes happens d...
Definition: unused.hh:16
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
constexpr GeometryType vertex
GeometryType representing a vertex.
Definition: type.hh:797
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
DInfoType dinfo(std::cout)
Stream for informative output.
Definition: stdstreams.hh:138
DVerbType dverb(std::cout)
Singleton of verbose debug stream.
Definition: stdstreams.hh:114
Provides a map between global and local indices.
Class for adding missing indices of a distributed index set in a local communication.
Traits classes for mapping types onto MPI_Datatype.
Dune namespace.
Definition: alignedallocator.hh:14
void fillIndexSetHoles(const G &graph, Dune::OwnerOverlapCopyCommunication< T1, T2 > &oocomm)
Fills the holes in an index set.
Definition: repartition.hh:82
bool graphRepartition(const G &graph, Dune::OwnerOverlapCopyCommunication< T1, T2 > &oocomm, Metis::idx_t nparts, std::shared_ptr< Dune::OwnerOverlapCopyCommunication< T1, T2 > > &outcomm, RedistributeInterface &redistInf, bool verbose=false)
execute a graph repartition for a giving graph and indexset.
Definition: repartition.hh:1268
Classes providing communication interfaces for overlapping Schwarz methods.
Utilities for reduction like operations on ranges.
Classes describing a distributed indexset.
Standard Dune debug streams.
A traits class describing the mapping of types onto MPI_Datatypes.
Definition: mpitraits.hh:38
@ nonoverlapping
Category for non-overlapping solvers.
Definition: solvercategory.hh:25
static Category category(const OP &op, decltype(op.category()) *=nullptr)
Helperfunction to extract the solver category either from an enum, or from the newly introduced virtu...
Definition: solvercategory.hh:32
A simple timing class.
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)