Dune Core Modules (2.10.0)

interface.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// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_COMMON_PARALLEL_INTERFACE_HH
6#define DUNE_COMMON_PARALLEL_INTERFACE_HH
7
8#if HAVE_MPI
9
10#include <cassert>
11#include <cstddef>
12#include <iostream>
13#include <map>
14#include <tuple>
15
16#include <mpi.h>
17
21
22namespace Dune
23{
44 {
45 public:
46 class RemoteIndicesStateError : public InvalidStateException
47 {};
48
49 virtual ~InterfaceBuilder()
50 {}
51
52 protected:
57 {}
58
96 template<class R, class T1, class T2, class Op, bool send>
97 void buildInterface (const R& remoteIndices,
98 const T1& sourceFlags, const T2& destFlags,
99 Op& functor) const;
100 };
101
110 {
111
112 public:
113
117 size_t size() const
118 {
119 return size_;
120 }
125 std::size_t& operator[](size_t i)
126 {
127 assert(i<size_);
128 return indices_[i];
129 }
134 std::size_t operator[](size_t i) const
135 {
136 assert(i<size_);
137 return indices_[i];
138 }
143 void reserve(size_t size)
144 {
145 indices_ = new std::size_t[size];
146 maxSize_ = size;
147
148 }
152 void free()
153 {
154 if(indices_)
155 delete[] indices_;
156 maxSize_ = 0;
157 size_=0;
158 indices_=0;
159 }
163 void add(std::size_t index)
164 {
165 assert(size_<maxSize_);
166 indices_[size_++]=index;
167 }
168
170 : size_(0), maxSize_(0), indices_(0)
171 {}
172
173 virtual ~InterfaceInformation()
174 {}
175
176 bool operator!=(const InterfaceInformation& o) const
177 {
178 return !operator==(o);
179 }
180
181 bool operator==(const InterfaceInformation& o) const
182 {
183 if(size_!=o.size_)
184 return false;
185 for(std::size_t i=0; i< size_; ++i)
186 if(indices_[i]!=o.indices_[i])
187 return false;
188 return true;
189 }
190
191 private:
195 size_t size_;
199 size_t maxSize_;
203 std::size_t* indices_;
204 };
205
218 {
219
220 public:
225 typedef std::map<int,std::pair<InterfaceInformation,InterfaceInformation> > InformationMap;
226
243 template<typename R, typename T1, typename T2>
244 void build(const R& remoteIndices, const T1& sourceFlags,
245 const T2& destFlags);
246
250 void free();
251
255 MPI_Comm communicator() const;
256
265 const InformationMap& interfaces() const;
266
267 Interface(MPI_Comm comm)
268 : communicator_(comm), interfaces_()
269 {}
270
271 Interface()
272 : communicator_(MPI_COMM_NULL), interfaces_()
273 {}
274
278 void print() const;
279
280 bool operator!=(const Interface& o) const
281 {
282 return ! operator==(o);
283 }
284
285 bool operator==(const Interface& o) const
286 {
287 if(communicator_!=o.communicator_)
288 return false;
289 if(interfaces_.size()!=o.interfaces_.size())
290 return false;
291 typedef InformationMap::const_iterator MIter;
292
293 for(MIter m=interfaces_.begin(), om=o.interfaces_.begin();
294 m!=interfaces_.end(); ++m, ++om)
295 {
296 if(om->first!=m->first)
297 return false;
298 if(om->second.first!=om->second.first)
299 return false;
300 if(om->second.second!=om->second.second)
301 return false;
302 }
303 return true;
304 }
305
309 virtual ~Interface();
310
311 void strip();
312 protected:
313
323
326
327 private:
335 InformationMap interfaces_;
336
337 template<bool send>
338 class InformationBuilder
339 {
340 public:
341 InformationBuilder(InformationMap& interfaces)
342 : interfaces_(interfaces)
343 {}
344
345 void reserve(int proc, int size)
346 {
347 if(send)
348 interfaces_[proc].first.reserve(size);
349 else
350 interfaces_[proc].second.reserve(size);
351 }
352 void add(int proc, std::size_t local)
353 {
354 if(send) {
355 interfaces_[proc].first.add(local);
356 }else{
357 interfaces_[proc].second.add(local);
358 }
359 }
360
361 private:
362 InformationMap& interfaces_;
363 };
364 };
365
366 template<class R, class T1, class T2, class Op, bool send>
367 void InterfaceBuilder::buildInterface(const R& remoteIndices, const T1& sourceFlags, const T2& destFlags, Op& interfaceInformation) const
368 {
369
370 if(!remoteIndices.isSynced())
371 DUNE_THROW(RemoteIndicesStateError,"RemoteIndices is not in sync with the index set. Call RemoteIndices::rebuild first!");
372 // Allocate the memory for the data type construction.
373 typedef R RemoteIndices;
374 typedef typename RemoteIndices::RemoteIndexMap::const_iterator const_iterator;
375
376 const const_iterator end=remoteIndices.end();
377
378 int rank;
379
380 MPI_Comm_rank(remoteIndices.communicator(), &rank);
381
382 // Allocate memory for the type construction.
383 for(const_iterator process=remoteIndices.begin(); process != end; ++process) {
384 // Measure the number of indices sent to the remote process first
385 int size=0;
386 typedef typename RemoteIndices::RemoteIndexList::const_iterator RemoteIterator;
387 const RemoteIterator remoteEnd = send ? process->second.first->end() :
388 process->second.second->end();
389 RemoteIterator remote = send ? process->second.first->begin() : process->second.second->begin();
390
391 while(remote!=remoteEnd) {
392 if( send ? destFlags.contains(remote->attribute()) :
393 sourceFlags.contains(remote->attribute())) {
394
395 // do we send the index?
396 if( send ? sourceFlags.contains(remote->localIndexPair().local().attribute()) :
397 destFlags.contains(remote->localIndexPair().local().attribute()))
398 ++size;
399 }
400 ++remote;
401 }
402 interfaceInformation.reserve(process->first, size);
403 }
404
405 // compare the local and remote indices and set up the types
406
407 for(const_iterator process=remoteIndices.begin(); process != end; ++process) {
408 typedef typename RemoteIndices::RemoteIndexList::const_iterator RemoteIterator;
409 const RemoteIterator remoteEnd = send ? process->second.first->end() :
410 process->second.second->end();
411 RemoteIterator remote = send ? process->second.first->begin() : process->second.second->begin();
412
413 while(remote!=remoteEnd) {
414 if( send ? destFlags.contains(remote->attribute()) :
415 sourceFlags.contains(remote->attribute())) {
416 // do we send the index?
417 if( send ? sourceFlags.contains(remote->localIndexPair().local().attribute()) :
418 destFlags.contains(remote->localIndexPair().local().attribute()))
419 interfaceInformation.add(process->first,remote->localIndexPair().local().local());
420 }
421 ++remote;
422 }
423 }
424 }
425
426 inline MPI_Comm Interface::communicator() const
427 {
428 return communicator_;
429
430 }
431
432
433 inline const std::map<int,std::pair<InterfaceInformation,InterfaceInformation> >& Interface::interfaces() const
434 {
435 return interfaces_;
436 }
437
438 inline std::map<int,std::pair<InterfaceInformation,InterfaceInformation> >& Interface::interfaces()
439 {
440 return interfaces_;
441 }
442
443 inline void Interface::print() const
444 {
445 typedef InformationMap::const_iterator const_iterator;
446 const const_iterator end=interfaces_.end();
447 int rank;
448 MPI_Comm_rank(communicator(), &rank);
449
450 for(const_iterator infoPair=interfaces_.begin(); infoPair!=end; ++infoPair) {
451 {
452 std::cout<<rank<<": send for process "<<infoPair->first<<": ";
453 const InterfaceInformation& info(infoPair->second.first);
454 for(size_t i=0; i < info.size(); i++)
455 std::cout<<info[i]<<" ";
456 std::cout<<std::endl;
457 } {
458
459 std::cout<<rank<<": receive for process "<<infoPair->first<<": ";
460 const InterfaceInformation& info(infoPair->second.second);
461 for(size_t i=0; i < info.size(); i++)
462 std::cout<<info[i]<<" ";
463 std::cout<<std::endl;
464 }
465
466 }
467 }
468
469 template<typename R, typename T1, typename T2>
470 inline void Interface::build(const R& remoteIndices, const T1& sourceFlags,
471 const T2& destFlags)
472 {
473 communicator_=remoteIndices.communicator();
474
475 assert(interfaces_.empty());
476
477 // Build the send interface
478 InformationBuilder<true> sendInformation(interfaces_);
479 this->template buildInterface<R,T1,T2,InformationBuilder<true>,true>(remoteIndices, sourceFlags,
480 destFlags, sendInformation);
481
482 // Build the receive interface
483 InformationBuilder<false> recvInformation(interfaces_);
484 this->template buildInterface<R,T1,T2,InformationBuilder<false>,false>(remoteIndices,sourceFlags,
485 destFlags, recvInformation);
486 strip();
487 }
488 inline void Interface::strip()
489 {
490 typedef InformationMap::iterator const_iterator;
491 for(const_iterator interfacePair = interfaces_.begin(); interfacePair != interfaces_.end();)
492 if(interfacePair->second.first.size()==0 && interfacePair->second.second.size()==0) {
493 interfacePair->second.first.free();
494 interfacePair->second.second.free();
495 const_iterator toerase=interfacePair++;
496 interfaces_.erase(toerase);
497 }else
498 ++interfacePair;
499 }
500
501 inline void Interface::free()
502 {
503 typedef InformationMap::iterator iterator;
504 typedef InformationMap::const_iterator const_iterator;
505 const const_iterator end = interfaces_.end();
506 for(iterator interfacePair = interfaces_.begin(); interfacePair != end; ++interfacePair) {
507 interfacePair->second.first.free();
508 interfacePair->second.second.free();
509 }
510 interfaces_.clear();
511 }
512
514 {
515 free();
516 }
519 inline std::ostream& operator<<(std::ostream& os, const Interface& interface)
520 {
521 typedef Interface::InformationMap InfoMap;
522 typedef InfoMap::const_iterator Iter;
523 for(Iter i=interface.interfaces().begin(), end = interface.interfaces().end();
524 i!=end; ++i)
525 {
526 os<<i->first<<": [ source=[";
527 for(std::size_t j=0; j < i->second.first.size(); ++j)
528 os<<i->second.first[j]<<" ";
529 os<<"] size="<<i->second.first.size()<<", target=[";
530 for(std::size_t j=0; j < i->second.second.size(); ++j)
531 os<<i->second.second[j]<<" ";
532 os<<"] size="<<i->second.second.size()<<"\n";
533 }
534 return os;
535 }
536}
537#endif // HAVE_MPI
538
539#endif // DUNE_COMMON_PARALLEL_INTERFACE_HH
Base class of all classes representing a communication interface.
Definition: interface.hh:44
InterfaceBuilder()
Not for public use.
Definition: interface.hh:56
Information describing an interface.
Definition: interface.hh:110
void free()
Definition: interface.hh:152
void add(std::size_t index)
Add a new index to the interface.
Definition: interface.hh:163
std::size_t & operator[](size_t i)
Get the local index for an entry.
Definition: interface.hh:125
size_t size() const
Get the number of entries in the interface.
Definition: interface.hh:117
void reserve(size_t size)
Reserve space for a number of entries.
Definition: interface.hh:143
std::size_t operator[](size_t i) const
Get the local index for an entry.
Definition: interface.hh:134
Communication interface between remote and local indices.
Definition: interface.hh:218
MPI_Comm communicator_
The MPI communicator we use.
Definition: interface.hh:325
std::map< int, std::pair< InterfaceInformation, InterfaceInformation > > InformationMap
The type of the map form process number to InterfaceInformation for sending and receiving to and from...
Definition: interface.hh:225
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:281
The indices present on remote processes.
Definition: remoteindices.hh:190
A constant iterator for the SLList.
Definition: sllist.hh:371
A few common exception classes.
Classes for building sets out of enumeration values.
void buildInterface(const R &remoteIndices, const T1 &sourceFlags, const T2 &destFlags, Op &functor) const
Builds the interface between remote processes.
Definition: interface.hh:367
void build(const R &remoteIndices, const T1 &sourceFlags, const T2 &destFlags)
Builds the interface.
Definition: interface.hh:470
void print() const
Print the interface to std::out for debugging.
Definition: interface.hh:443
const_iterator end() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1528
const InformationMap & interfaces() const
Get information about the interfaces.
Definition: interface.hh:433
MPI_Comm communicator() const
Get the MPI Communicator.
Definition: interface.hh:426
virtual ~Interface()
Destructor.
Definition: interface.hh:513
void free()
Frees memory allocated during the build.
Definition: interface.hh:501
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
Classes describing a distributed indexset.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 26, 23:30, 2024)