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