Dune Core Modules (2.9.0)

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