Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
a3b9f5f0fe6d5e414e8b678fdeb927a0b78c8a1e
[simgrid.git] / src / smpi / smpi_comm.cpp
1 /* Copyright (c) 2010-2017. The SimGrid Team. All rights reserved.          */
2
3 /* This program is free software; you can redistribute it and/or modify it
4  * under the terms of the license (GNU LGPL) which comes with this package. */
5
6 #include "simgrid/s4u/Host.hpp"
7
8 #include "private.h"
9 #include "src/simix/smx_private.h"
10
11 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_comm, smpi, "Logging specific to SMPI (comm)");
12
13  simgrid::smpi::Comm mpi_MPI_COMM_UNINITIALIZED;
14 MPI_Comm MPI_COMM_UNINITIALIZED=&mpi_MPI_COMM_UNINITIALIZED;
15
16 /* Support for cartesian topology was added, but there are 2 other types of topology, graph et dist graph. In order to
17  * support them, we have to add a field MPIR_Topo_type, and replace the MPI_Topology field by an union. */
18
19 static int smpi_compare_rankmap(const void *a, const void *b)
20 {
21   const int* x = static_cast<const int*>(a);
22   const int* y = static_cast<const int*>(b);
23
24   if (x[1] < y[1]) {
25     return -1;
26   }
27   if (x[1] == y[1]) {
28     if (x[0] < y[0]) {
29       return -1;
30     }
31     if (x[0] == y[0]) {
32       return 0;
33     }
34     return 1;
35   }
36   return 1;
37 }
38
39 namespace simgrid{
40 namespace smpi{
41
42 std::unordered_map<int, smpi_key_elem> Comm::keyvals_;
43 int Comm::keyval_id_=0;
44
45 Comm::Comm(MPI_Group group, MPI_Topology topo) : group_(group), topo_(topo)
46 {
47   refcount_=1;
48   topoType_ = MPI_INVALID_TOPO;
49   intra_comm_ = MPI_COMM_NULL;
50   leaders_comm_ = MPI_COMM_NULL;
51   is_uniform_=1;
52   non_uniform_map_ = nullptr;
53   leaders_map_ = nullptr;
54   is_blocked_=0;
55 }
56
57 void Comm::destroy(Comm* comm)
58 {
59   if (comm == MPI_COMM_UNINITIALIZED){
60     Comm::destroy(smpi_process()->comm_world());
61     return;
62   }
63   delete comm->topo_; // there's no use count on topos
64   Comm::unref(comm);
65 }
66
67 int Comm::dup(MPI_Comm* newcomm){
68   if(smpi_privatize_global_variables == SMPI_PRIVATIZE_MMAP){ //we need to switch as the called function may silently touch global variables
69      smpi_switch_data_segment(smpi_process()->index());
70    }
71   MPI_Group cp = new  Group(this->group());
72   (*newcomm) = new  Comm(cp, this->topo());
73   int ret = MPI_SUCCESS;
74
75   if(!attributes()->empty()){
76     int flag;
77     void* value_out;
78     for(auto it : *attributes()){
79       smpi_key_elem elem = keyvals_.at(it.first);
80       if (elem != nullptr && elem->copy_fn.comm_copy_fn != MPI_NULL_COPY_FN) {
81         ret = elem->copy_fn.comm_copy_fn(this, it.first, nullptr, it.second, &value_out, &flag);
82         if (ret != MPI_SUCCESS) {
83           Comm::destroy(*newcomm);
84           *newcomm = MPI_COMM_NULL;
85           return ret;
86         }
87         if (flag){
88           elem->refcount++;
89           (*newcomm)->attributes()->insert({it.first, value_out});
90         }
91       }
92       }
93     }
94   return ret;
95 }
96
97 MPI_Group Comm::group()
98 {
99   if (this == MPI_COMM_UNINITIALIZED)
100     return smpi_process()->comm_world()->group();
101   return group_;
102 }
103
104 MPI_Topology Comm::topo() {
105   return topo_;
106 }
107
108 int Comm::size()
109 {
110   if (this == MPI_COMM_UNINITIALIZED)
111     return smpi_process()->comm_world()->size();
112   return group_->size();
113 }
114
115 int Comm::rank()
116 {
117   if (this == MPI_COMM_UNINITIALIZED)
118     return smpi_process()->comm_world()->rank();
119   return group_->rank(smpi_process()->index());
120 }
121
122 void Comm::get_name (char* name, int* len)
123 {
124   if (this == MPI_COMM_UNINITIALIZED){
125     smpi_process()->comm_world()->get_name(name, len);
126     return;
127   }
128   if(this == MPI_COMM_WORLD) {
129     strncpy(name, "WORLD",5);
130     *len = 5;
131   } else {
132     *len = snprintf(name, MPI_MAX_NAME_STRING, "%p", this);
133   }
134 }
135
136 void Comm::set_leaders_comm(MPI_Comm leaders){
137   if (this == MPI_COMM_UNINITIALIZED){
138     smpi_process()->comm_world()->set_leaders_comm(leaders);
139     return;
140   }
141   leaders_comm_=leaders;
142 }
143
144 void Comm::set_intra_comm(MPI_Comm leaders){
145   intra_comm_=leaders;
146 }
147
148 int* Comm::get_non_uniform_map(){
149   if (this == MPI_COMM_UNINITIALIZED)
150     return smpi_process()->comm_world()->get_non_uniform_map();
151   return non_uniform_map_;
152 }
153
154 int* Comm::get_leaders_map(){
155   if (this == MPI_COMM_UNINITIALIZED)
156     return smpi_process()->comm_world()->get_leaders_map();
157   return leaders_map_;
158 }
159
160 MPI_Comm Comm::get_leaders_comm(){
161   if (this == MPI_COMM_UNINITIALIZED)
162     return smpi_process()->comm_world()->get_leaders_comm();
163   return leaders_comm_;
164 }
165
166 MPI_Comm Comm::get_intra_comm(){
167   if (this == MPI_COMM_UNINITIALIZED || this==MPI_COMM_WORLD) 
168     return smpi_process()->comm_intra();
169   else return intra_comm_;
170 }
171
172 int Comm::is_uniform(){
173   if (this == MPI_COMM_UNINITIALIZED)
174     return smpi_process()->comm_world()->is_uniform();
175   return is_uniform_;
176 }
177
178 int Comm::is_blocked(){
179   if (this == MPI_COMM_UNINITIALIZED)
180     return smpi_process()->comm_world()->is_blocked();
181   return is_blocked_;
182 }
183
184 MPI_Comm Comm::split(int color, int key)
185 {
186   if (this == MPI_COMM_UNINITIALIZED)
187     return smpi_process()->comm_world()->split(color, key);
188   int system_tag = 123;
189   int* recvbuf;
190
191   MPI_Group group_root = nullptr;
192   MPI_Group group_out  = nullptr;
193   MPI_Group group      = this->group();
194   int rank             = this->rank();
195   int size             = this->size();
196   /* Gather all colors and keys on rank 0 */
197   int* sendbuf = xbt_new(int, 2);
198   sendbuf[0] = color;
199   sendbuf[1] = key;
200   if(rank == 0) {
201     recvbuf = xbt_new(int, 2 * size);
202   } else {
203     recvbuf = nullptr;
204   }
205   Coll_gather_default::gather(sendbuf, 2, MPI_INT, recvbuf, 2, MPI_INT, 0, this);
206   xbt_free(sendbuf);
207   /* Do the actual job */
208   if(rank == 0) {
209     MPI_Group* group_snd = xbt_new(MPI_Group, size);
210     int* rankmap         = xbt_new(int, 2 * size);
211     for (int i = 0; i < size; i++) {
212       if (recvbuf[2 * i] != MPI_UNDEFINED) {
213         int count = 0;
214         for (int j = i + 1; j < size; j++) {
215           if(recvbuf[2 * i] == recvbuf[2 * j]) {
216             recvbuf[2 * j] = MPI_UNDEFINED;
217             rankmap[2 * count] = j;
218             rankmap[2 * count + 1] = recvbuf[2 * j + 1];
219             count++;
220           }
221         }
222         /* Add self in the group */
223         recvbuf[2 * i] = MPI_UNDEFINED;
224         rankmap[2 * count] = i;
225         rankmap[2 * count + 1] = recvbuf[2 * i + 1];
226         count++;
227         qsort(rankmap, count, 2 * sizeof(int), &smpi_compare_rankmap);
228         group_out = new  Group(count);
229         if (i == 0) {
230           group_root = group_out; /* Save root's group */
231         }
232         for (int j = 0; j < count; j++) {
233           int index = group->index(rankmap[2 * j]);
234           group_out->set_mapping(index, j);
235         }
236         MPI_Request* requests = xbt_new(MPI_Request, count);
237         int reqs              = 0;
238         for (int j = 0; j < count; j++) {
239           if(rankmap[2 * j] != 0) {
240             group_snd[reqs]=new  Group(group_out);
241             requests[reqs] = Request::isend(&(group_snd[reqs]), 1, MPI_PTR, rankmap[2 * j], system_tag, this);
242             reqs++;
243           }
244         }
245         if(i != 0 && group_out != MPI_COMM_WORLD->group() && group_out != MPI_GROUP_EMPTY)
246           Group::unref(group_out);
247         
248         Request::waitall(reqs, requests, MPI_STATUS_IGNORE);
249         xbt_free(requests);
250       }
251     }
252     xbt_free(recvbuf);
253     xbt_free(rankmap);
254     xbt_free(group_snd);
255     group_out = group_root; /* exit with root's group */
256   } else {
257     if(color != MPI_UNDEFINED) {
258       Request::recv(&group_out, 1, MPI_PTR, 0, system_tag, this, MPI_STATUS_IGNORE);
259     } /* otherwise, exit with group_out == nullptr */
260   }
261   return group_out!=nullptr ? new  Comm(group_out, nullptr) : MPI_COMM_NULL;
262 }
263
264 void Comm::ref(){
265   if (this == MPI_COMM_UNINITIALIZED){
266     smpi_process()->comm_world()->ref();
267     return;
268   }
269   group_->ref();
270   refcount_++;
271 }
272
273 void Comm::cleanup_smp(){
274   if (intra_comm_ != MPI_COMM_NULL)
275     Comm::unref(intra_comm_);
276   if (leaders_comm_ != MPI_COMM_NULL)
277     Comm::unref(leaders_comm_);
278   if (non_uniform_map_ != nullptr)
279     xbt_free(non_uniform_map_);
280   if (leaders_map_ != nullptr)
281     xbt_free(leaders_map_);
282 }
283
284 void Comm::unref(Comm* comm){
285   if (comm == MPI_COMM_UNINITIALIZED){
286     Comm::unref(smpi_process()->comm_world());
287     return;
288   }
289   comm->refcount_--;
290   Group::unref(comm->group_);
291
292   if(comm->refcount_==0){
293     comm->cleanup_smp();
294     comm->cleanup_attr<Comm>();
295     delete comm;
296   }
297 }
298
299 static int compare_ints (const void *a, const void *b)
300 {
301   const int *da = static_cast<const int *>(a);
302   const int *db = static_cast<const int *>(b);
303
304   return static_cast<int>(*da > *db) - static_cast<int>(*da < *db);
305 }
306
307 void Comm::init_smp(){
308   int leader = -1;
309
310   if (this == MPI_COMM_UNINITIALIZED)
311     smpi_process()->comm_world()->init_smp();
312
313   int comm_size = this->size();
314   
315   // If we are in replay - perform an ugly hack  
316   // tell SimGrid we are not in replay for a while, because we need the buffers to be copied for the following calls
317   bool replaying = false; //cache data to set it back again after
318   if(smpi_process()->replaying()){
319    replaying=true;
320    smpi_process()->set_replaying(false);
321   }
322
323   if(smpi_privatize_global_variables == SMPI_PRIVATIZE_MMAP){ //we need to switch as the called function may silently touch global variables
324      smpi_switch_data_segment(smpi_process()->index());
325    }
326   //identify neighbours in comm
327   //get the indexes of all processes sharing the same simix host
328   xbt_swag_t process_list = SIMIX_host_self()->extension<simgrid::simix::Host>()->process_list;
329   int intra_comm_size     = 0;
330   int min_index           = INT_MAX;//the minimum index will be the leader
331   smx_actor_t actor       = nullptr;
332   xbt_swag_foreach(actor, process_list) {
333     int index = actor->pid -1;
334
335     if(this->group()->rank(index)!=MPI_UNDEFINED){
336       intra_comm_size++;
337       //the process is in the comm
338       if(index < min_index)
339         min_index=index;
340     }
341   }
342   XBT_DEBUG("number of processes deployed on my node : %d", intra_comm_size);
343   MPI_Group group_intra = new  Group(intra_comm_size);
344   int i = 0;
345   actor = nullptr;
346   xbt_swag_foreach(actor, process_list) {
347     int index = actor->pid -1;
348     if(this->group()->rank(index)!=MPI_UNDEFINED){
349       group_intra->set_mapping(index, i);
350       i++;
351     }
352   }
353
354   MPI_Comm comm_intra = new  Comm(group_intra, nullptr);
355   leader=min_index;
356
357   int * leaders_map= static_cast<int*>(xbt_malloc0(sizeof(int)*comm_size));
358   int * leader_list= static_cast<int*>(xbt_malloc0(sizeof(int)*comm_size));
359   for(i=0; i<comm_size; i++){
360       leader_list[i]=-1;
361   }
362
363   Coll_allgather_mpich::allgather(&leader, 1, MPI_INT , leaders_map, 1, MPI_INT, this);
364
365   if(smpi_privatize_global_variables == SMPI_PRIVATIZE_MMAP){ //we need to switch as the called function may silently touch global variables
366      smpi_switch_data_segment(smpi_process()->index());
367    }
368
369   if(leaders_map_==nullptr){
370     leaders_map_= leaders_map;
371   }else{
372     xbt_free(leaders_map);
373   }
374   int j=0;
375   int leader_group_size = 0;
376   for(i=0; i<comm_size; i++){
377       int already_done=0;
378       for(j=0;j<leader_group_size; j++){
379         if(leaders_map_[i]==leader_list[j]){
380             already_done=1;
381         }
382       }
383       if(already_done==0){
384         leader_list[leader_group_size]=leaders_map_[i];
385         leader_group_size++;
386       }
387   }
388   qsort(leader_list, leader_group_size, sizeof(int),compare_ints);
389
390   MPI_Group leaders_group = new  Group(leader_group_size);
391
392   MPI_Comm leader_comm = MPI_COMM_NULL;
393   if(MPI_COMM_WORLD!=MPI_COMM_UNINITIALIZED && this!=MPI_COMM_WORLD){
394     //create leader_communicator
395     for (i=0; i< leader_group_size;i++)
396       leaders_group->set_mapping(leader_list[i], i);
397     leader_comm = new  Comm(leaders_group, nullptr);
398     this->set_leaders_comm(leader_comm);
399     this->set_intra_comm(comm_intra);
400
401    //create intracommunicator
402   }else{
403     for (i=0; i< leader_group_size;i++)
404       leaders_group->set_mapping(leader_list[i], i);
405
406     if(this->get_leaders_comm()==MPI_COMM_NULL){
407       leader_comm = new  Comm(leaders_group, nullptr);
408       this->set_leaders_comm(leader_comm);
409     }else{
410       leader_comm=this->get_leaders_comm();
411       Group::unref(leaders_group);
412     }
413     smpi_process()->set_comm_intra(comm_intra);
414   }
415
416   int is_uniform = 1;
417
418   // Are the nodes uniform ? = same number of process/node
419   int my_local_size=comm_intra->size();
420   if(comm_intra->rank()==0) {
421     int* non_uniform_map = xbt_new0(int,leader_group_size);
422     Coll_allgather_mpich::allgather(&my_local_size, 1, MPI_INT,
423         non_uniform_map, 1, MPI_INT, leader_comm);
424     for(i=0; i < leader_group_size; i++) {
425       if(non_uniform_map[0] != non_uniform_map[i]) {
426         is_uniform = 0;
427         break;
428       }
429     }
430     if(is_uniform==0 && this->is_uniform()!=0){
431         non_uniform_map_= non_uniform_map;
432     }else{
433         xbt_free(non_uniform_map);
434     }
435     is_uniform_=is_uniform;
436   }
437   Coll_bcast_mpich::bcast(&(is_uniform_),1, MPI_INT, 0, comm_intra );
438
439   if(smpi_privatize_global_variables == SMPI_PRIVATIZE_MMAP){ //we need to switch as the called function may silently touch global variables
440      smpi_switch_data_segment(smpi_process()->index());
441    }
442   // Are the ranks blocked ? = allocated contiguously on the SMP nodes
443   int is_blocked=1;
444   int prev=this->group()->rank(comm_intra->group()->index(0));
445     for (i=1; i<my_local_size; i++){
446       int that=this->group()->rank(comm_intra->group()->index(i));
447       if(that!=prev+1){
448         is_blocked=0;
449         break;
450       }
451       prev = that;
452   }
453
454   int global_blocked;
455   Coll_allreduce_default::allreduce(&is_blocked, &(global_blocked), 1, MPI_INT, MPI_LAND, this);
456
457   if(MPI_COMM_WORLD==MPI_COMM_UNINITIALIZED || this==MPI_COMM_WORLD){
458     if(this->rank()==0){
459         is_blocked_=global_blocked;
460     }
461   }else{
462     is_blocked_=global_blocked;
463   }
464   xbt_free(leader_list);
465   
466   if(replaying)
467     smpi_process()->set_replaying(true); 
468 }
469
470 MPI_Comm Comm::f2c(int id) {
471   if(id == -2) {
472     return MPI_COMM_SELF;
473   } else if(id==0){
474     return MPI_COMM_WORLD;
475   } else if(F2C::f2c_lookup() != nullptr && id >= 0) {
476       char key[KEY_SIZE];
477       MPI_Comm tmp =  static_cast<MPI_Comm>(xbt_dict_get_or_null(F2C::f2c_lookup(),get_key_id(key, id)));
478       return tmp != nullptr ? tmp : MPI_COMM_NULL ;
479   } else {
480     return MPI_COMM_NULL;
481   }
482 }
483
484 void Comm::free_f(int id) {
485   char key[KEY_SIZE];
486   xbt_dict_remove(F2C::f2c_lookup(), id==0? get_key(key, id) : get_key_id(key, id));
487 }
488
489 int Comm::add_f() {
490   if(F2C::f2c_lookup()==nullptr){
491     F2C::set_f2c_lookup(xbt_dict_new_homogeneous(nullptr));
492   }
493   char key[KEY_SIZE];
494   xbt_dict_set(F2C::f2c_lookup(), this==MPI_COMM_WORLD? get_key(key, F2C::f2c_id()) : get_key_id(key,F2C::f2c_id()), this, nullptr);
495   f2c_id_increment();
496   return F2C::f2c_id()-1;
497 }
498
499
500 void Comm::add_rma_win(MPI_Win win){
501   rma_wins_.push_back(win);
502 }
503
504 void Comm::remove_rma_win(MPI_Win win){
505   rma_wins_.remove(win);
506 }
507
508 void Comm::finish_rma_calls(){
509   for(auto it : rma_wins_){
510     if(it->rank()==this->rank()){//is it ours (for MPI_COMM_WORLD)?
511       int finished = it->finish_comms();
512       XBT_DEBUG("Barrier for rank %d - Finished %d RMA calls",this->rank(), finished);
513     }
514   }
515 }
516
517
518 }
519 }
520
521