Geant4  10.02.p01
G4MPIutils.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 #include "G4MPIutils.hh"
27 #include <iostream>
28 #include <cstdlib>
29 #include <algorithm>
30 #include <functional>
31 #include <assert.h>
32 #include <utility>
33 #include "globals.hh"
34 
35 
36 
38  std::vector<G4mpi::rank_t>& input ) {
39  using namespace G4mpi;
40  //Check validity of input
41  std::sort(input.begin(),input.end());
42  if ( input.size() < 1 || input[0] != 0 ) {
43  G4Exception("G4mpi::buildCommunicationMap(...)","G4mpi001",FatalException,
44  "Empty input or cannot find rank 0 in input.");
45  }
46  //Requested that no duplicates!
47  std::vector<rank_t> copy(input.size());
48  std::copy(input.begin(),input.end(),copy.begin());
49  copy.erase( std::unique(copy.begin(),copy.end()),copy.end());
50  if ( copy != input )
51  {
52  G4Exception("G4mpi::buildCommunicationMap(...)","G4mpi001",FatalException,
53  "There are duplicates in list of input ranks.");
54  }
55 
56  //The final communication map
57  commMap_t mymap;
58  //The communication map key
59  int cycle = 0;
60  //The communication map value
61  std::vector<couple_t> couples;
62  //Start a loop (on cycles) that will break
63  do {
64  //An helper container
65  std::vector<rank_t> receiving;
66  couples.clear();
67  //Loop on all input until there is
68  //at least a couple
69  while ( input.size() > 1 ) {
70  //Sort input in ascending order
71  std::sort( input.begin(),input.end() );
72  //Pop from back of the input a couple
73  const auto& send_ = input.back();
74  input.pop_back();
75  const auto& rec_ = input.back();
76  input.pop_back();
77  //Actually the receiving is not empty,
78  //remember it because we have to add it back
79  receiving.push_back( rec_ );
80  //This is a couple for this cycle
81  couples.push_back( std::make_pair(send_,rec_) );
82  }
83  //Populate final map for this cycle
84  mymap[cycle++]=couples;
85  //Let's put back in the input container the receivers
86  input.insert( input.end() , receiving.begin() , receiving.end() );
87  //Let's continue until there is ony one rank in input (number 0)
88  } while ( input.size()!=1 );
89  return mymap;
90 }
91 
92 //Test function
93 int _testMe(int argc,char** argv) {
94  using namespace G4mpi;
95  unsigned int worldSize = 10;
96  if ( argc > 1 ) worldSize = atoi(argv[1]);
97  unsigned int myRank = worldSize-1;
98  if ( argc > 2 ) myRank = atoi(argv[2]);
99  std::cout<<"World size: "<<worldSize<<std::endl;
100 
101  assert( myRank < worldSize);
102  //MPI function stubs
103  auto MPI_Receive = [](const rank_t& s, const rank_t& r) {
104  std::cout<<"MPI_Receive from: "<<s<<" to "<<r<<std::endl;
105  return 0;
106  };
107  auto MPI_Send = [](const rank_t& s, const rank_t& r) {
108  std::cout<<"MPI_Send from: "<<s<<" to "<<r<<std::endl;;
109  return 0;
110  };
111  auto MPI_Barrier = [] {
112  std::cout<<"MPI_Barrier"<<std::endl;
113  return 0;
114  };
115 
116  //Build the initial network of ranks
117  std::vector<rank_t> ranks(worldSize);
118  for ( unsigned int i = 0 ; i<worldSize ; ++i ) {
119  if ( i != 2 )
120  { ranks.push_back(i);
121  } else {
122  ranks.push_back(i);
123  ranks.push_back(i);
124  }
125  }
126  //Remove duplicates
127  ranks.erase( std::unique(ranks.begin(),ranks.end()),ranks.end());
128 
129  //Optimize network trafic
130  if ( ranks.size() == 1 ) {
131  std::cout<<"only one rank, nothing to do"<<std::endl;
132  return 0;
133  }
134  auto comms = G4mpi::buildCommunicationMap( ranks );
135  assert( ranks.size() == 1 && ranks[0] == 0 );
136 
137  std::cout<<"Communiction Map (size: "<<comms.size()<<"):"<<std::endl;
138  for ( const auto& x : comms ) {
139  std::cout<<"Cycle "<<x.first<<": ";
140  for ( const auto& y : x.second ) {
141  std::cout<<y.first<<"->"<<y.second<<", ";
142  }
143  std::cout<<std::endl;
144  }
145 
146  std::cout<<"Simulate communication pattern for rank: "<<myRank<<std::endl;
147  for (const auto& x: comms ) {
148  std::cout<<"Cycle "<<x.first<<std::endl;
149  for ( const auto& y : x.second ) {
150  if ( myRank == y.first ) { MPI_Send(y.first,y.second); }
151  else if ( myRank == y.second ) { MPI_Receive(y.first,y.second); }
152  }
153  //Important: Wait for this cycle to end before going to the next, even if
154  //this rank did not do anything
155  //This is needed to be sure that the redcutions are done correctly
156  MPI_Barrier();
157  }
158  return 0;
159 }
160 
161 void G4mpi::Merge( std::function<void(unsigned int)> senderF ,
162  std::function<void(unsigned int)> receiverF,
163  std::function<void(void)> barrierF ,
164  unsigned int commSize ,
165  unsigned int myrank) {
166  //Optimize communications between ranks
167  std::vector<G4mpi::rank_t> ranks(commSize);
168  std::iota(ranks.begin(),ranks.end(),0); //{0,1,2,3,...}
169  auto comms = G4mpi::buildCommunicationMap(ranks);
170  //Loop on cycles of communications
171  for ( const auto& cycle : comms ) {
172  //Each cycle is a set of communications between ranks, it is guarantted that
173  //each rank participate in one and only one communication for each cycle
174  for (const auto& pattern : cycle.second ) {
175  //pattern is a couple: sender,receiver
176  if ( myrank == pattern.first ) {
177  //Send to destination
178  senderF(pattern.second);
179  }
180  else if ( myrank == pattern.second ) {
181  //Receive from source
182  receiverF(pattern.first);
183  }
184  }
185  //Important: Wait for this cycle to end before going to the next, even if this rank
186  //did not do anything
187  //This is needed to be sure that the redcutions are done correctly
188  barrierF();
189  }
190 }
commMap_t buildCommunicationMap(std::vector< rank_t > &input)
Definition: G4MPIutils.cc:37
void Merge(std::function< void(unsigned int)> senderF, std::function< void(unsigned int)> receiverF, std::function< void(void)> barrierF, unsigned int commSize, unsigned int myrank)
Definition: G4MPIutils.cc:161
G4double(* function)(G4double)
static const double s
Definition: G4SIunits.hh:168
unsigned int rank_t
Definition: G4MPIutils.hh:45
std::map< int, std::vector< couple_t > > commMap_t
Definition: G4MPIutils.hh:50
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4double x[NPOINTSGL]
int _testMe(int argc, char **argv)
Definition: G4MPIutils.cc:93