Logo AND Algorithmique Numérique Distribuée

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