Dune Core Modules (2.3.1)

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