Dune Core Modules (2.7.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#ifndef DUNE_INTERFACE_HH
4#define DUNE_INTERFACE_HH
5
6#if HAVE_MPI
7
8#include "remoteindices.hh"
10
11namespace Dune
12{
33 {
34 public:
35 class RemoteIndicesStateError : public InvalidStateException
36 {};
37
38 virtual ~InterfaceBuilder()
39 {}
40
41 protected:
46 {}
47
85 template<class R, class T1, class T2, class Op, bool send>
86 void buildInterface (const R& remoteIndices,
87 const T1& sourceFlags, const T2& destFlags,
88 Op& functor) const;
89 };
90
99 {
100
101 public:
102
106 size_t size() const
107 {
108 return size_;
109 }
114 std::size_t& operator[](size_t i)
115 {
116 assert(i<size_);
117 return indices_[i];
118 }
123 std::size_t operator[](size_t i) const
124 {
125 assert(i<size_);
126 return indices_[i];
127 }
132 void reserve(size_t size)
133 {
134 indices_ = new std::size_t[size];
135 maxSize_ = size;
136
137 }
141 void free()
142 {
143 if(indices_)
144 delete[] indices_;
145 maxSize_ = 0;
146 size_=0;
147 indices_=0;
148 }
152 void add(std::size_t index)
153 {
154 assert(size_<maxSize_);
155 indices_[size_++]=index;
156 }
157
159 : size_(0), maxSize_(0), indices_(0)
160 {}
161
162 virtual ~InterfaceInformation()
163 {}
164
165 bool operator!=(const InterfaceInformation& o) const
166 {
167 return !operator==(o);
168 }
169
170 bool operator==(const InterfaceInformation& o) const
171 {
172 if(size_!=o.size_)
173 return false;
174 for(std::size_t i=0; i< size_; ++i)
175 if(indices_[i]!=o.indices_[i])
176 return false;
177 return true;
178 }
179
180 private:
184 size_t size_;
188 size_t maxSize_;
192 std::size_t* indices_;
193 };
194
207 {
208
209 public:
214 typedef std::map<int,std::pair<InterfaceInformation,InterfaceInformation> > InformationMap;
215
232 template<typename R, typename T1, typename T2>
233 void build(const R& remoteIndices, const T1& sourceFlags,
234 const T2& destFlags);
235
239 void free();
240
244 MPI_Comm communicator() const;
245
254 const InformationMap& interfaces() const;
255
256 Interface(MPI_Comm comm)
257 : communicator_(comm), interfaces_()
258 {}
259
260 Interface()
261 : communicator_(MPI_COMM_NULL), interfaces_()
262 {}
263
267 void print() const;
268
269 bool operator!=(const Interface& o) const
270 {
271 return ! operator==(o);
272 }
273
274 bool operator==(const Interface& o) const
275 {
276 if(communicator_!=o.communicator_)
277 return false;
278 if(interfaces_.size()!=o.interfaces_.size())
279 return false;
280 typedef InformationMap::const_iterator MIter;
281
282 for(MIter m=interfaces_.begin(), om=o.interfaces_.begin();
283 m!=interfaces_.end(); ++m, ++om)
284 {
285 if(om->first!=m->first)
286 return false;
287 if(om->second.first!=om->second.first)
288 return false;
289 if(om->second.second!=om->second.second)
290 return false;
291 }
292 return true;
293 }
294
298 virtual ~Interface();
299
300 void strip();
301 protected:
302
312
315
316 private:
324 InformationMap interfaces_;
325
326 template<bool send>
327 class InformationBuilder
328 {
329 public:
330 InformationBuilder(InformationMap& interfaces)
331 : interfaces_(interfaces)
332 {}
333
334 void reserve(int proc, int size)
335 {
336 if(send)
337 interfaces_[proc].first.reserve(size);
338 else
339 interfaces_[proc].second.reserve(size);
340 }
341 void add(int proc, std::size_t local)
342 {
343 if(send) {
344 interfaces_[proc].first.add(local);
345 }else{
346 interfaces_[proc].second.add(local);
347 }
348 }
349
350 private:
351 InformationMap& interfaces_;
352 };
353 };
354
355 template<class R, class T1, class T2, class Op, bool send>
356 void InterfaceBuilder::buildInterface(const R& remoteIndices, const T1& sourceFlags, const T2& destFlags, Op& interfaceInformation) const
357 {
358
359 if(!remoteIndices.isSynced())
360 DUNE_THROW(RemoteIndicesStateError,"RemoteIndices is not in sync with the index set. Call RemoteIndices::rebuild first!");
361 // Allocate the memory for the data type construction.
362 typedef R RemoteIndices;
363 typedef typename RemoteIndices::RemoteIndexMap::const_iterator const_iterator;
364
365 const const_iterator end=remoteIndices.end();
366
367 int rank;
368
369 MPI_Comm_rank(remoteIndices.communicator(), &rank);
370
371 // Allocate memory for the type construction.
372 for(const_iterator process=remoteIndices.begin(); process != end; ++process) {
373 // Messure the number of indices send to the remote process first
374 int size=0;
375 typedef typename RemoteIndices::RemoteIndexList::const_iterator RemoteIterator;
376 const RemoteIterator remoteEnd = send ? process->second.first->end() :
377 process->second.second->end();
378 RemoteIterator remote = send ? process->second.first->begin() : process->second.second->begin();
379
380 while(remote!=remoteEnd) {
381 if( send ? destFlags.contains(remote->attribute()) :
382 sourceFlags.contains(remote->attribute())) {
383
384 // do we send the index?
385 if( send ? sourceFlags.contains(remote->localIndexPair().local().attribute()) :
386 destFlags.contains(remote->localIndexPair().local().attribute()))
387 ++size;
388 }
389 ++remote;
390 }
391 interfaceInformation.reserve(process->first, size);
392 }
393
394 // compare the local and remote indices and set up the types
395
396 for(const_iterator process=remoteIndices.begin(); process != end; ++process) {
397 typedef typename RemoteIndices::RemoteIndexList::const_iterator RemoteIterator;
398 const RemoteIterator remoteEnd = send ? process->second.first->end() :
399 process->second.second->end();
400 RemoteIterator remote = send ? process->second.first->begin() : process->second.second->begin();
401
402 while(remote!=remoteEnd) {
403 if( send ? destFlags.contains(remote->attribute()) :
404 sourceFlags.contains(remote->attribute())) {
405 // do we send the index?
406 if( send ? sourceFlags.contains(remote->localIndexPair().local().attribute()) :
407 destFlags.contains(remote->localIndexPair().local().attribute()))
408 interfaceInformation.add(process->first,remote->localIndexPair().local().local());
409 }
410 ++remote;
411 }
412 }
413 }
414
415 inline MPI_Comm Interface::communicator() const
416 {
417 return communicator_;
418
419 }
420
421
422 inline const std::map<int,std::pair<InterfaceInformation,InterfaceInformation> >& Interface::interfaces() const
423 {
424 return interfaces_;
425 }
426
427 inline std::map<int,std::pair<InterfaceInformation,InterfaceInformation> >& Interface::interfaces()
428 {
429 return interfaces_;
430 }
431
432 inline void Interface::print() const
433 {
434 typedef InformationMap::const_iterator const_iterator;
435 const const_iterator end=interfaces_.end();
436 int rank;
437 MPI_Comm_rank(communicator(), &rank);
438
439 for(const_iterator infoPair=interfaces_.begin(); infoPair!=end; ++infoPair) {
440 {
441 std::cout<<rank<<": send for process "<<infoPair->first<<": ";
442 const InterfaceInformation& info(infoPair->second.first);
443 for(size_t i=0; i < info.size(); i++)
444 std::cout<<info[i]<<" ";
445 std::cout<<std::endl;
446 } {
447
448 std::cout<<rank<<": receive for process "<<infoPair->first<<": ";
449 const InterfaceInformation& info(infoPair->second.second);
450 for(size_t i=0; i < info.size(); i++)
451 std::cout<<info[i]<<" ";
452 std::cout<<std::endl;
453 }
454
455 }
456 }
457
458 template<typename R, typename T1, typename T2>
459 inline void Interface::build(const R& remoteIndices, const T1& sourceFlags,
460 const T2& destFlags)
461 {
462 communicator_=remoteIndices.communicator();
463
464 assert(interfaces_.empty());
465
466 // Build the send interface
467 InformationBuilder<true> sendInformation(interfaces_);
468 this->template buildInterface<R,T1,T2,InformationBuilder<true>,true>(remoteIndices, sourceFlags,
469 destFlags, sendInformation);
470
471 // Build the receive interface
472 InformationBuilder<false> recvInformation(interfaces_);
473 this->template buildInterface<R,T1,T2,InformationBuilder<false>,false>(remoteIndices,sourceFlags,
474 destFlags, recvInformation);
475 strip();
476 }
477 inline void Interface::strip()
478 {
479 typedef InformationMap::iterator const_iterator;
480 for(const_iterator interfacePair = interfaces_.begin(); interfacePair != interfaces_.end();)
481 if(interfacePair->second.first.size()==0 && interfacePair->second.second.size()==0) {
482 interfacePair->second.first.free();
483 interfacePair->second.second.free();
484 const_iterator toerase=interfacePair++;
485 interfaces_.erase(toerase);
486 }else
487 ++interfacePair;
488 }
489
490 inline void Interface::free()
491 {
492 typedef InformationMap::iterator iterator;
493 typedef InformationMap::const_iterator const_iterator;
494 const const_iterator end = interfaces_.end();
495 for(iterator interfacePair = interfaces_.begin(); interfacePair != end; ++interfacePair) {
496 interfacePair->second.first.free();
497 interfacePair->second.second.free();
498 }
499 interfaces_.clear();
500 }
501
503 {
504 free();
505 }
508 inline std::ostream& operator<<(std::ostream& os, const Interface& interface)
509 {
510 typedef Interface::InformationMap InfoMap;
511 typedef InfoMap::const_iterator Iter;
512 for(Iter i=interface.interfaces().begin(), end = interface.interfaces().end();
513 i!=end; ++i)
514 {
515 os<<i->first<<": [ source=[";
516 for(std::size_t j=0; j < i->second.first.size(); ++j)
517 os<<i->second.first[j]<<" ";
518 os<<"] size="<<i->second.first.size()<<", target=[";
519 for(std::size_t j=0; j < i->second.second.size(); ++j)
520 os<<i->second.second[j]<<" ";
521 os<<"] size="<<i->second.second.size()<<"\n";
522 }
523 return os;
524 }
525}
526#endif // HAVE_MPI
527
528#endif
Base class of all classes representing a communication interface.
Definition: interface.hh:33
InterfaceBuilder()
Not for public use.
Definition: interface.hh:45
Information describing an interface.
Definition: interface.hh:99
void free()
Definition: interface.hh:141
void add(std::size_t index)
Add a new index to the interface.
Definition: interface.hh:152
std::size_t & operator[](size_t i)
Get the local index for an entry.
Definition: interface.hh:114
size_t size() const
Get the number of entries in the interface.
Definition: interface.hh:106
void reserve(size_t size)
Reserve space for a number of entries.
Definition: interface.hh:132
std::size_t operator[](size_t i) const
Get the local index for an entry.
Definition: interface.hh:123
Communication interface between remote and local indices.
Definition: interface.hh:207
MPI_Comm communicator_
The MPI communicator we use.
Definition: interface.hh:314
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:214
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:279
The indices present on remote processes.
Definition: remoteindices.hh:187
const_iterator end() const
Get an iterator over all remote index lists.
Definition: remoteindices.hh:1535
A constant iterator for the SLList.
Definition: sllist.hh:369
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:356
void build(const R &remoteIndices, const T1 &sourceFlags, const T2 &destFlags)
Builds the interface.
Definition: interface.hh:459
void print() const
Print the interface to std::out for debugging.
Definition: interface.hh:432
const InformationMap & interfaces() const
Get information about the interfaces.
Definition: interface.hh:422
MPI_Comm communicator() const
Get the MPI Communicator.
Definition: interface.hh:415
virtual ~Interface()
Destructor.
Definition: interface.hh:502
void free()
Frees memory allocated during the build.
Definition: interface.hh:490
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
Dune namespace.
Definition: alignedallocator.hh:14
Classes describing a distributed indexset.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Oct 13, 22:30, 2024)