Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge pull request #112 from glesserd/teshpy
[simgrid.git] / src / smpi / smpi_base.cpp
1 /* Copyright (c) 2007-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 <xbt/config.hpp>
8 #include <boost/tokenizer.hpp>
9 #include <algorithm>
10
11 #include "private.h"
12 #include "xbt/virtu.h"
13 #include "mc/mc.h"
14 #include "src/mc/mc_replay.h"
15 #include "xbt/replay.h"
16 #include <errno.h>
17 #include "src/simix/smx_private.h"
18 #include "surf/surf.h"
19 #include "simgrid/sg_config.h"
20 #include "colls/colls.h"
21
22 #include "src/kernel/activity/SynchroComm.hpp"
23
24 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_base, smpi, "Logging specific to SMPI (base)");
25
26 static int match_recv(void* a, void* b, smx_activity_t ignored) {
27    MPI_Request ref = static_cast<MPI_Request>(a);
28    MPI_Request req = static_cast<MPI_Request>(b);
29    XBT_DEBUG("Trying to match a recv of src %d against %d, tag %d against %d",ref->src,req->src, ref->tag, req->tag);
30
31   xbt_assert(ref, "Cannot match recv against null reference");
32   xbt_assert(req, "Cannot match recv against null request");
33   if((ref->src == MPI_ANY_SOURCE || req->src == ref->src)
34     && ((ref->tag == MPI_ANY_TAG && req->tag >=0) || req->tag == ref->tag)){
35     //we match, we can transfer some values
36     if(ref->src == MPI_ANY_SOURCE)
37         ref->real_src = req->src;
38     if(ref->tag == MPI_ANY_TAG)
39         ref->real_tag = req->tag;
40     if(ref->real_size < req->real_size) 
41         ref->truncated = 1;
42     if(req->detached==1)
43         ref->detached_sender=req; //tie the sender to the receiver, as it is detached and has to be freed in the receiver
44     XBT_DEBUG("match succeeded");
45     return 1;
46   }else return 0;
47 }
48
49 static int match_send(void* a, void* b,smx_activity_t ignored) {
50    MPI_Request ref = static_cast<MPI_Request>(a);
51    MPI_Request req = static_cast<MPI_Request>(b);
52    XBT_DEBUG("Trying to match a send of src %d against %d, tag %d against %d",ref->src,req->src, ref->tag, req->tag);
53    xbt_assert(ref, "Cannot match send against null reference");
54    xbt_assert(req, "Cannot match send against null request");
55
56    if((req->src == MPI_ANY_SOURCE || req->src == ref->src)
57              && ((req->tag == MPI_ANY_TAG && ref->tag >=0)|| req->tag == ref->tag))
58    {
59      if(req->src == MPI_ANY_SOURCE)
60         req->real_src = ref->src;
61      if(req->tag == MPI_ANY_TAG)
62         req->real_tag = ref->tag;
63      if(req->real_size < ref->real_size) 
64         req->truncated = 1;
65      if(ref->detached==1)
66          req->detached_sender=ref; //tie the sender to the receiver, as it is detached and has to be freed in the receiver
67      XBT_DEBUG("match succeeded");
68      return 1;
69    } else return 0;
70 }
71
72 // Methods used to parse and store the values for timing injections in smpi
73 // These are taken from surf/network.c and generalized to have more values for each factor
74 typedef struct s_smpi_factor_multival *smpi_os_factor_multival_t;
75 typedef struct s_smpi_factor_multival { // FIXME: this should be merged (deduplicated) with s_smpi_factor defined in network_smpi.c
76   size_t factor=0;
77   std::vector<double> values;
78 } s_smpi_factor_multival_t;
79
80 std::vector<s_smpi_factor_multival_t> smpi_os_values;
81 std::vector<s_smpi_factor_multival_t> smpi_or_values;
82 std::vector<s_smpi_factor_multival_t> smpi_ois_values;
83
84 static simgrid::config::Flag<double> smpi_wtime_sleep(
85   "smpi/wtime", "Minimum time to inject inside a call to MPI_Wtime", 0.0);
86 static simgrid::config::Flag<double> smpi_init_sleep(
87   "smpi/init", "Time to inject inside a call to MPI_Init", 0.0);
88 static simgrid::config::Flag<double> smpi_iprobe_sleep(
89   "smpi/iprobe", "Minimum time to inject inside a call to MPI_Iprobe", 1e-4);
90 static simgrid::config::Flag<double> smpi_test_sleep(
91   "smpi/test", "Minimum time to inject inside a call to MPI_Test", 1e-4);
92
93 static std::vector<s_smpi_factor_multival_t> parse_factor(const char *smpi_coef_string)
94 {
95   std::vector<s_smpi_factor_multival_t> smpi_factor;
96
97   /** Setup the tokenizer that parses the string **/
98   typedef boost::tokenizer<boost::char_separator<char>> Tokenizer;
99   boost::char_separator<char> sep(";");
100   boost::char_separator<char> factor_separator(":");
101   std::string tmp_string(smpi_coef_string);
102   Tokenizer tokens(tmp_string, sep);
103
104   /** 
105    * Iterate over patterns like A:B:C:D;E:F;G:H
106    * These will be broken down into:
107    * A --> B, C, D
108    * E --> F
109    * G --> H
110    */
111   for (Tokenizer::iterator token_iter = tokens.begin();
112          token_iter != tokens.end(); token_iter++) {
113 XBT_DEBUG("token : %s", token_iter->c_str());
114     Tokenizer factor_values(*token_iter, factor_separator);
115     s_smpi_factor_multival_t fact;
116     if (factor_values.begin() == factor_values.end()) {
117       xbt_die("Malformed radical for smpi factor: '%s'", smpi_coef_string);
118     }
119     unsigned int iteration = 0;
120     for (Tokenizer::iterator factor_iter = factor_values.begin();
121          factor_iter != factor_values.end(); factor_iter++, iteration++) {
122       char *errmsg;
123
124       if (factor_iter == factor_values.begin()) { /* first element */
125         errmsg = bprintf("Invalid factor in chunk #%zu: %%s", smpi_factor.size()+1);
126         fact.factor = xbt_str_parse_int(factor_iter->c_str(), errmsg);
127       }
128       else {
129         errmsg = bprintf("Invalid factor value %d in chunk #%zu: %%s", iteration, smpi_factor.size()+1);
130         fact.values.push_back(xbt_str_parse_double(factor_iter->c_str(), errmsg));
131       }
132       xbt_free(errmsg);
133     }
134
135     smpi_factor.push_back(fact);
136     XBT_DEBUG("smpi_factor:\t%zu : %zu values, first: %f", fact.factor, smpi_factor.size(), fact.values[0]);
137   }
138   std::sort(smpi_factor.begin(), smpi_factor.end(),
139             [](const s_smpi_factor_multival_t &pa,
140                const s_smpi_factor_multival_t &pb) {
141               return (pa.factor < pb.factor);
142             });
143   for (auto& fact : smpi_factor) {
144     XBT_DEBUG("smpi_factor:\t%zu : %zu values, first: %f", fact.factor, smpi_factor.size() ,fact.values[0]);
145   }
146
147   return smpi_factor;
148 }
149
150 static double smpi_os(size_t size)
151 {
152   if (smpi_os_values.empty()) {
153     smpi_os_values = parse_factor(xbt_cfg_get_string("smpi/os"));
154   }
155   double current=smpi_os_values.empty()?0.0:smpi_os_values[0].values[0]+smpi_os_values[0].values[1]*size;
156   // Iterate over all the sections that were specified and find the right
157   // value. (fact.factor represents the interval sizes; we want to find the
158   // section that has fact.factor <= size and no other such fact.factor <= size)
159   // Note: parse_factor() (used before) already sorts the vector we iterate over!
160   for (auto& fact : smpi_os_values) {
161     if (size <= fact.factor) { // Values already too large, use the previously
162                                // computed value of current!
163         XBT_DEBUG("os : %zu <= %zu return %.10f", size, fact.factor, current);
164       return current;
165     }else{
166       // If the next section is too large, the current section must be used.
167       // Hence, save the cost, as we might have to use it.
168       current = fact.values[0]+fact.values[1]*size;
169     }
170   }
171   XBT_DEBUG("Searching for smpi/os: %zu is larger than the largest boundary, return %.10f", size, current);
172
173   return current;
174 }
175
176 static double smpi_ois(size_t size)
177 {
178   if (smpi_ois_values.empty()) {
179     smpi_ois_values = parse_factor(xbt_cfg_get_string("smpi/ois"));
180   }
181   double current=smpi_ois_values.empty()?0.0:smpi_ois_values[0].values[0]+smpi_ois_values[0].values[1]*size;
182   // Iterate over all the sections that were specified and find the right value. (fact.factor represents the interval
183   // sizes; we want to find the section that has fact.factor <= size and no other such fact.factor <= size)
184   // Note: parse_factor() (used before) already sorts the vector we iterate over!
185   for (auto& fact : smpi_ois_values) {
186     if (size <= fact.factor) { // Values already too large, use the previously  computed value of current!
187         XBT_DEBUG("ois : %zu <= %zu return %.10f", size, fact.factor, current);
188       return current;
189     }else{
190       // If the next section is too large, the current section must be used.
191       // Hence, save the cost, as we might have to use it.
192       current = fact.values[0]+fact.values[1]*size;
193     }
194   }
195   XBT_DEBUG("Searching for smpi/ois: %zu is larger than the largest boundary, return %.10f", size, current);
196
197   return current;
198 }
199
200 static double smpi_or(size_t size)
201 {
202   if (smpi_or_values.empty()) {
203     smpi_or_values = parse_factor(xbt_cfg_get_string("smpi/or"));
204   }
205   
206   double current=smpi_or_values.empty()?0.0:smpi_or_values.front().values[0]+smpi_or_values.front().values[1]*size;
207
208   // Iterate over all the sections that were specified and find the right value. (fact.factor represents the interval
209   // sizes; we want to find the section that has fact.factor <= size and no other such fact.factor <= size)
210   // Note: parse_factor() (used before) already sorts the vector we iterate over!
211   for (auto fact : smpi_or_values) {
212     if (size <= fact.factor) { // Values already too large, use the previously computed value of current!
213         XBT_DEBUG("or : %zu <= %zu return %.10f", size, fact.factor, current);
214       return current;
215     } else {
216       // If the next section is too large, the current section must be used.
217       // Hence, save the cost, as we might have to use it.
218       current=fact.values[0]+fact.values[1]*size;
219     }
220   }
221   XBT_DEBUG("smpi_or: %zu is larger than largest boundary, return %.10f", size, current);
222
223   return current;
224 }
225
226 void smpi_mpi_init() {
227   if(smpi_init_sleep > 0) 
228     simcall_process_sleep(smpi_init_sleep);
229 }
230
231 double smpi_mpi_wtime(){
232   double time;
233   if (smpi_process_initialized() != 0 && smpi_process_finalized() == 0 && smpi_process_get_sampling() == 0) {
234     smpi_bench_end();
235     time = SIMIX_get_clock();
236     // to avoid deadlocks if used as a break condition, such as
237     //     while (MPI_Wtime(...) < time_limit) {
238     //       ....
239     //     }
240     // because the time will not normally advance when only calls to MPI_Wtime
241     // are made -> deadlock (MPI_Wtime never reaches the time limit)
242     if(smpi_wtime_sleep > 0) 
243       simcall_process_sleep(smpi_wtime_sleep);
244     smpi_bench_begin();
245   } else {
246     time = SIMIX_get_clock();
247   }
248   return time;
249 }
250
251 static MPI_Request build_request(void *buf, int count, MPI_Datatype datatype, int src, int dst, int tag, MPI_Comm comm,
252                                  unsigned flags)
253 {
254   MPI_Request request = nullptr;
255
256   void *old_buf = nullptr;
257
258   request = xbt_new(s_smpi_mpi_request_t, 1);
259
260   s_smpi_subtype_t *subtype = static_cast<s_smpi_subtype_t*>(datatype->substruct);
261
262   if((((flags & RECV) != 0) && ((flags & ACCUMULATE) !=0)) || (datatype->sizeof_substruct != 0)){
263     // This part handles the problem of non-contiguous memory
264     old_buf = buf;
265     buf = count==0 ? nullptr : xbt_malloc(count*smpi_datatype_size(datatype));
266     if ((datatype->sizeof_substruct != 0) && ((flags & SEND) != 0)) {
267       subtype->serialize(old_buf, buf, count, datatype->substruct);
268     }
269   }
270
271   request->buf      = buf;
272   // This part handles the problem of non-contiguous memory (for the unserialisation at the reception)
273   request->old_buf  = old_buf;
274   request->old_type = datatype;
275
276   request->size = smpi_datatype_size(datatype) * count;
277   smpi_datatype_use(datatype);
278   request->src  = src;
279   request->dst  = dst;
280   request->tag  = tag;
281   request->comm = comm;
282   smpi_comm_use(request->comm);
283   request->action          = nullptr;
284   request->flags           = flags;
285   request->detached        = 0;
286   request->detached_sender = nullptr;
287   request->real_src        = 0;
288   request->truncated       = 0;
289   request->real_size       = 0;
290   request->real_tag        = 0;
291   if (flags & PERSISTENT)
292     request->refcount = 1;
293   else
294     request->refcount = 0;
295   request->op   = MPI_REPLACE;
296   request->send = 0;
297   request->recv = 0;
298
299   return request;
300 }
301
302 void smpi_empty_status(MPI_Status * status)
303 {
304   if(status != MPI_STATUS_IGNORE) {
305     status->MPI_SOURCE = MPI_ANY_SOURCE;
306     status->MPI_TAG = MPI_ANY_TAG;
307     status->MPI_ERROR = MPI_SUCCESS;
308     status->count=0;
309   }
310 }
311
312 static void smpi_mpi_request_free_voidp(void* request)
313 {
314   MPI_Request req = static_cast<MPI_Request>(request);
315   smpi_mpi_request_free(&req);
316 }
317
318 /* MPI Low level calls */
319 MPI_Request smpi_mpi_send_init(void *buf, int count, MPI_Datatype datatype,
320                                int dst, int tag, MPI_Comm comm)
321 {
322   MPI_Request request = nullptr; /* MC needs the comm to be set to nullptr during the call */
323   request = build_request(buf==MPI_BOTTOM ? nullptr : buf, count, datatype, smpi_process_index(),
324                           smpi_group_index(smpi_comm_group(comm), dst), tag, comm, PERSISTENT | SEND | PREPARED);
325   return request;
326 }
327
328 MPI_Request smpi_mpi_ssend_init(void *buf, int count, MPI_Datatype datatype,
329                                int dst, int tag, MPI_Comm comm)
330 {
331   MPI_Request request = nullptr; /* MC needs the comm to be set to nullptr during the call */
332   request = build_request(buf==MPI_BOTTOM ? nullptr : buf, count, datatype, smpi_process_index(),
333                         smpi_group_index(smpi_comm_group(comm), dst), tag, comm, PERSISTENT | SSEND | SEND | PREPARED);
334   return request;
335 }
336
337 MPI_Request smpi_mpi_recv_init(void *buf, int count, MPI_Datatype datatype,
338                                int src, int tag, MPI_Comm comm)
339 {
340   MPI_Request request = nullptr; /* MC needs the comm to be set to nullptr during the call */
341   request = build_request(buf==MPI_BOTTOM ? nullptr : buf, count, datatype,
342                           src == MPI_ANY_SOURCE ? MPI_ANY_SOURCE : smpi_group_index(smpi_comm_group(comm), src),
343                           smpi_process_index(), tag, comm, PERSISTENT | RECV | PREPARED);
344   return request;
345 }
346
347 void smpi_mpi_start(MPI_Request request)
348 {
349   smx_mailbox_t mailbox;
350
351   xbt_assert(request->action == nullptr, "Cannot (re-)start unfinished communication");
352   request->flags &= ~PREPARED;
353   request->flags &= ~FINISHED;
354   request->refcount++;
355
356   if ((request->flags & RECV) != 0) {
357     print_request("New recv", request);
358
359     int async_small_thresh = xbt_cfg_get_int("smpi/async-small-thresh");
360
361     xbt_mutex_t mut = smpi_process_mailboxes_mutex();
362     if (async_small_thresh != 0 || (request->flags & RMA) != 0)
363       xbt_mutex_acquire(mut);
364
365     if (async_small_thresh == 0 && (request->flags & RMA) == 0 ) {
366       mailbox = smpi_process_mailbox();
367     } 
368     else if (((request->flags & RMA) != 0) || static_cast<int>(request->size) < async_small_thresh) {
369       //We have to check both mailboxes (because SSEND messages are sent to the large mbox).
370       //begin with the more appropriate one : the small one.
371       mailbox = smpi_process_mailbox_small();
372       XBT_DEBUG("Is there a corresponding send already posted in the small mailbox %p (in case of SSEND)?", mailbox);
373       smx_activity_t action = simcall_comm_iprobe(mailbox, 0, request->src,request->tag, &match_recv, static_cast<void*>(request));
374     
375       if (action == nullptr) {
376         mailbox = smpi_process_mailbox();
377         XBT_DEBUG("No, nothing in the small mailbox test the other one : %p", mailbox);
378         action = simcall_comm_iprobe(mailbox, 0, request->src,request->tag, &match_recv, static_cast<void*>(request));
379         if (action == nullptr) {
380           XBT_DEBUG("Still nothing, switch back to the small mailbox : %p", mailbox);
381           mailbox = smpi_process_mailbox_small();
382         }
383       }
384       else {
385         XBT_DEBUG("yes there was something for us in the large mailbox");
386       }
387     }
388     else {
389       mailbox = smpi_process_mailbox_small();
390       XBT_DEBUG("Is there a corresponding send already posted the small mailbox?");
391       smx_activity_t action = simcall_comm_iprobe(mailbox, 0, request->src,request->tag, &match_recv, (void*)request);
392     
393       if (action == nullptr) {
394         XBT_DEBUG("No, nothing in the permanent receive mailbox");
395         mailbox = smpi_process_mailbox();
396       }
397       else {
398         XBT_DEBUG("yes there was something for us in the small mailbox");
399       }
400     }
401
402     // we make a copy here, as the size is modified by simix, and we may reuse the request in another receive later
403     request->real_size=request->size;
404     request->action = simcall_comm_irecv(SIMIX_process_self(), mailbox, request->buf, &request->real_size, &match_recv,
405                                          ! smpi_process_get_replaying()? &smpi_comm_copy_buffer_callback
406                                          : &smpi_comm_null_copy_buffer_callback, request, -1.0);
407         XBT_DEBUG("recv simcall posted");
408
409     if (async_small_thresh != 0 || (request->flags & RMA) != 0 )
410       xbt_mutex_release(mut);
411   }
412   else { /* the RECV flag was not set, so this is a send */
413     int receiver = request->dst;
414
415     int rank = request->src;
416     if (TRACE_smpi_view_internals()) {
417       TRACE_smpi_send(rank, rank, receiver,request->size);
418     }
419     print_request("New send", request);
420
421     void* buf = request->buf;
422     if ( (request->flags & SSEND) == 0 
423         && ( (request->flags & RMA) != 0 || static_cast<int>(request->size) < xbt_cfg_get_int("smpi/send-is-detached-thresh") ) ) {
424       void *oldbuf = nullptr;
425       request->detached = 1;
426       XBT_DEBUG("Send request %p is detached", request);
427       request->refcount++;
428       if(request->old_type->sizeof_substruct == 0){
429         oldbuf = request->buf;
430         if (!smpi_process_get_replaying() && oldbuf != nullptr && request->size!=0){
431           if((smpi_privatize_global_variables != 0)
432             && (static_cast<char*>(request->buf) >= smpi_start_data_exe)
433             && (static_cast<char*>(request->buf) < smpi_start_data_exe + smpi_size_data_exe )){
434             XBT_DEBUG("Privatization : We are sending from a zone inside global memory. Switch data segment ");
435             smpi_switch_data_segment(request->src);
436           }
437           buf = xbt_malloc(request->size);
438           memcpy(buf,oldbuf,request->size);
439           XBT_DEBUG("buf %p copied into %p",oldbuf,buf);
440         }
441       }
442     }
443
444     //if we are giving back the control to the user without waiting for completion, we have to inject timings
445     double sleeptime = 0.0;
446     if(request->detached != 0 || ((request->flags & (ISEND|SSEND)) != 0)){// issend should be treated as isend
447       //isend and send timings may be different
448       sleeptime = ((request->flags & ISEND) != 0) ? smpi_ois(request->size) : smpi_os(request->size);
449     }
450
451     if(sleeptime > 0.0){
452         simcall_process_sleep(sleeptime);
453         XBT_DEBUG("sending size of %zu : sleep %f ", request->size, sleeptime);
454     } 
455
456     int async_small_thresh = xbt_cfg_get_int("smpi/async-small-thresh");
457
458     xbt_mutex_t mut=smpi_process_remote_mailboxes_mutex(receiver);
459
460     if (async_small_thresh != 0 || (request->flags & RMA) != 0)
461       xbt_mutex_acquire(mut);
462
463     if (!(async_small_thresh != 0 || (request->flags & RMA) !=0)) {
464       mailbox = smpi_process_remote_mailbox(receiver);
465     }
466     else if (((request->flags & RMA) != 0) || static_cast<int>(request->size) < async_small_thresh) { // eager mode
467       mailbox = smpi_process_remote_mailbox(receiver);
468       XBT_DEBUG("Is there a corresponding recv already posted in the large mailbox %p?", mailbox);
469       smx_activity_t action = simcall_comm_iprobe(mailbox, 1,request->dst, request->tag, &match_send, static_cast<void*>(request));
470       if (action == nullptr) {
471         if ((request->flags & SSEND) == 0){
472           mailbox = smpi_process_remote_mailbox_small(receiver);
473           XBT_DEBUG("No, nothing in the large mailbox, message is to be sent on the small one %p", mailbox);
474         } 
475         else {
476           mailbox = smpi_process_remote_mailbox_small(receiver);
477           XBT_DEBUG("SSEND : Is there a corresponding recv already posted in the small mailbox %p?", mailbox);
478           action = simcall_comm_iprobe(mailbox, 1,request->dst, request->tag, &match_send, static_cast<void*>(request));
479           if (action == nullptr) {
480             XBT_DEBUG("No, we are first, send to large mailbox");
481             mailbox = smpi_process_remote_mailbox(receiver);
482           }
483         }
484       }
485       else {
486         XBT_DEBUG("Yes there was something for us in the large mailbox");
487       }
488     }
489     else {
490       mailbox = smpi_process_remote_mailbox(receiver);
491       XBT_DEBUG("Send request %p is in the large mailbox %p (buf: %p)",mailbox, request,request->buf);
492     }
493
494     // we make a copy here, as the size is modified by simix, and we may reuse the request in another receive later
495     request->real_size=request->size;
496     request->action = simcall_comm_isend(SIMIX_process_from_PID(request->src+1), mailbox, request->size, -1.0,
497                                          buf, request->real_size, &match_send,
498                          &xbt_free_f, // how to free the userdata if a detached send fails
499                          !smpi_process_get_replaying() ? &smpi_comm_copy_buffer_callback
500                          : &smpi_comm_null_copy_buffer_callback, request,
501                          // detach if msg size < eager/rdv switch limit
502                          request->detached);
503     XBT_DEBUG("send simcall posted");
504
505     /* FIXME: detached sends are not traceable (request->action == nullptr) */
506     if (request->action != nullptr)
507       simcall_set_category(request->action, TRACE_internal_smpi_get_category());
508
509     if (async_small_thresh != 0 || ((request->flags & RMA)!=0))
510       xbt_mutex_release(mut);
511   }
512 }
513
514 void smpi_mpi_startall(int count, MPI_Request * requests)
515 {
516   if(requests== nullptr) 
517     return;
518
519   for(int i = 0; i < count; i++) {
520     smpi_mpi_start(requests[i]);
521   }
522 }
523
524 void smpi_mpi_request_free(MPI_Request * request)
525 {
526   if((*request) != MPI_REQUEST_NULL){
527     (*request)->refcount--;
528     if((*request)->refcount<0) xbt_die("wrong refcount");
529
530     if((*request)->refcount==0){
531         smpi_datatype_unuse((*request)->old_type);
532         smpi_comm_unuse((*request)->comm);
533         print_request("Destroying", (*request));
534         xbt_free(*request);
535         *request = MPI_REQUEST_NULL;
536     }else{
537         print_request("Decrementing", (*request));
538     }
539   }else{
540       xbt_die("freeing an already free request");
541   }
542 }
543
544 MPI_Request smpi_rma_send_init(void *buf, int count, MPI_Datatype datatype, int src, int dst, int tag, MPI_Comm comm,
545                                MPI_Op op)
546 {
547   MPI_Request request = nullptr; /* MC needs the comm to be set to nullptr during the call */
548   if(op==MPI_OP_NULL){
549     request = build_request(buf==MPI_BOTTOM ? nullptr : buf , count, datatype, src, dst, tag,
550                             comm, RMA | NON_PERSISTENT | ISEND | SEND | PREPARED);
551   }else{
552     request = build_request(buf==MPI_BOTTOM ? nullptr : buf, count, datatype,  src, dst, tag,
553                             comm, RMA | NON_PERSISTENT | ISEND | SEND | PREPARED | ACCUMULATE);
554     request->op = op;
555   }
556   return request;
557 }
558
559 MPI_Request smpi_rma_recv_init(void *buf, int count, MPI_Datatype datatype, int src, int dst, int tag, MPI_Comm comm,
560                                MPI_Op op)
561 {
562   MPI_Request request = nullptr; /* MC needs the comm to be set to nullptr during the call */
563   if(op==MPI_OP_NULL){
564     request = build_request(buf==MPI_BOTTOM ? nullptr : buf, count, datatype,  src, dst, tag,
565                             comm, RMA | NON_PERSISTENT | RECV | PREPARED);
566   }else{
567     request = build_request(buf==MPI_BOTTOM ? nullptr : buf, count, datatype,  src, dst, tag,
568                             comm, RMA | NON_PERSISTENT | RECV | PREPARED | ACCUMULATE);
569     request->op = op;
570   }
571   return request;
572 }
573
574 MPI_Request smpi_isend_init(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
575 {
576   MPI_Request request = nullptr; /* MC needs the comm to be set to nullptr during the call */
577   request = build_request(buf==MPI_BOTTOM ? nullptr : buf , count, datatype, smpi_process_index(),
578                           smpi_group_index(smpi_comm_group(comm), dst), tag,comm, PERSISTENT | ISEND | SEND | PREPARED);
579   return request;
580 }
581
582 MPI_Request smpi_mpi_isend(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
583 {
584   MPI_Request request = nullptr; /* MC needs the comm to be set to nullptr during the call */
585   request =  build_request(buf==MPI_BOTTOM ? nullptr : buf, count, datatype, smpi_process_index(),
586                            smpi_group_index(smpi_comm_group(comm), dst), tag, comm, NON_PERSISTENT | ISEND | SEND);
587   smpi_mpi_start(request);
588   return request;
589 }
590
591 MPI_Request smpi_mpi_issend(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
592 {
593   MPI_Request request = nullptr; /* MC needs the comm to be set to nullptr during the call */
594   request = build_request(buf==MPI_BOTTOM ? nullptr : buf, count, datatype, smpi_process_index(),
595                         smpi_group_index(smpi_comm_group(comm), dst), tag,comm, NON_PERSISTENT | ISEND | SSEND | SEND);
596   smpi_mpi_start(request);
597   return request;
598 }
599
600 MPI_Request smpi_irecv_init(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm)
601 {
602   MPI_Request request = nullptr; /* MC needs the comm to be set to nullptr during the call */
603   request = build_request(buf==MPI_BOTTOM ? nullptr : buf, count, datatype, src == MPI_ANY_SOURCE ? MPI_ANY_SOURCE :
604                           smpi_group_index(smpi_comm_group(comm), src), smpi_process_index(), tag,
605                           comm, PERSISTENT | RECV | PREPARED);
606   return request;
607 }
608
609 MPI_Request smpi_mpi_irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm)
610 {
611   MPI_Request request = nullptr; /* MC needs the comm to be set to nullptr during the call */
612   request = build_request(buf==MPI_BOTTOM ? nullptr : buf, count, datatype, src == MPI_ANY_SOURCE ? MPI_ANY_SOURCE :
613                           smpi_group_index(smpi_comm_group(comm), src), smpi_process_index(), tag, comm,
614                           NON_PERSISTENT | RECV);
615   smpi_mpi_start(request);
616   return request;
617 }
618
619 void smpi_mpi_recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
620 {
621   MPI_Request request = nullptr; /* MC needs the comm to be set to nullptr during the call */
622   request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
623   smpi_mpi_wait(&request, status);
624   request = nullptr;
625 }
626
627 void smpi_mpi_send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
628 {
629   MPI_Request request = nullptr; /* MC needs the comm to be set to nullptr during the call */
630   request = build_request(buf==MPI_BOTTOM ? nullptr : buf, count, datatype, smpi_process_index(),
631                           smpi_group_index(smpi_comm_group(comm), dst), tag, comm, NON_PERSISTENT | SEND);
632
633   smpi_mpi_start(request);
634   smpi_mpi_wait(&request, MPI_STATUS_IGNORE);
635   request = nullptr;
636 }
637
638 void smpi_mpi_ssend(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
639 {
640   MPI_Request request = nullptr; /* MC needs the comm to be set to nullptr during the call */
641   request = build_request(buf==MPI_BOTTOM ? nullptr : buf, count, datatype, smpi_process_index(),
642                           smpi_group_index(smpi_comm_group(comm), dst), tag, comm, NON_PERSISTENT | SSEND | SEND);
643
644   smpi_mpi_start(request);
645   smpi_mpi_wait(&request, MPI_STATUS_IGNORE);
646   request = nullptr;
647 }
648
649 void smpi_mpi_sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,int dst, int sendtag,
650                        void *recvbuf, int recvcount, MPI_Datatype recvtype, int src, int recvtag,
651                        MPI_Comm comm, MPI_Status * status)
652 {
653   MPI_Request requests[2];
654   MPI_Status stats[2];
655   int myid=smpi_process_index();
656   if ((smpi_group_index(smpi_comm_group(comm), dst) == myid) && (smpi_group_index(smpi_comm_group(comm), src) == myid)){
657       smpi_datatype_copy(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype);
658       return;
659   }
660   requests[0] = smpi_isend_init(sendbuf, sendcount, sendtype, dst, sendtag, comm);
661   requests[1] = smpi_irecv_init(recvbuf, recvcount, recvtype, src, recvtag, comm);
662   smpi_mpi_startall(2, requests);
663   smpi_mpi_waitall(2, requests, stats);
664   smpi_mpi_request_free(&requests[0]);
665   smpi_mpi_request_free(&requests[1]);
666   if(status != MPI_STATUS_IGNORE) {
667     // Copy receive status
668     *status = stats[1];
669   }
670 }
671
672 int smpi_mpi_get_count(MPI_Status * status, MPI_Datatype datatype)
673 {
674   return status->count / smpi_datatype_size(datatype);
675 }
676
677 static void finish_wait(MPI_Request * request, MPI_Status * status)
678 {
679   MPI_Request req = *request;
680   smpi_empty_status(status);
681
682   if(!((req->detached != 0) && ((req->flags & SEND) != 0)) && ((req->flags & PREPARED) == 0)){
683     if(status != MPI_STATUS_IGNORE) {
684       int src = req->src == MPI_ANY_SOURCE ? req->real_src : req->src;
685       status->MPI_SOURCE = smpi_group_rank(smpi_comm_group(req->comm), src);
686       status->MPI_TAG = req->tag == MPI_ANY_TAG ? req->real_tag : req->tag;
687       status->MPI_ERROR = req->truncated != 0 ? MPI_ERR_TRUNCATE : MPI_SUCCESS;
688       // this handles the case were size in receive differs from size in send
689       status->count = req->real_size;
690     }
691
692     print_request("Finishing", req);
693     MPI_Datatype datatype = req->old_type;
694
695     if(((req->flags & ACCUMULATE) != 0) || (datatype->sizeof_substruct != 0)){
696       if (!smpi_process_get_replaying()){
697         if( smpi_privatize_global_variables != 0 && (static_cast<char*>(req->old_buf) >= smpi_start_data_exe)
698             && ((char*)req->old_buf < smpi_start_data_exe + smpi_size_data_exe )){
699             XBT_VERB("Privatization : We are unserializing to a zone in global memory - Switch data segment ");
700             smpi_switch_data_segment(smpi_process_index());
701         }
702       }
703
704       if(datatype->sizeof_substruct != 0){
705         // This part handles the problem of non-contignous memory the unserialization at the reception
706         s_smpi_subtype_t *subtype = static_cast<s_smpi_subtype_t*>(datatype->substruct);
707         if(req->flags & RECV)
708           subtype->unserialize(req->buf, req->old_buf, req->real_size/smpi_datatype_size(datatype) ,
709                                datatype->substruct, req->op);
710         xbt_free(req->buf);
711       }else if(req->flags & RECV){//apply op on contiguous buffer for accumulate
712           int n =req->real_size/smpi_datatype_size(datatype);
713           smpi_op_apply(req->op, req->buf, req->old_buf, &n, &datatype);
714           xbt_free(req->buf);
715       }
716     }
717   }
718
719   if (TRACE_smpi_view_internals() && ((req->flags & RECV) != 0)){
720     int rank = smpi_process_index();
721     int src_traced = (req->src == MPI_ANY_SOURCE ? req->real_src : req->src);
722     TRACE_smpi_recv(rank, src_traced, rank);
723   }
724
725   if(req->detached_sender != nullptr){
726
727     //integrate pseudo-timing for buffering of small messages, do not bother to execute the simcall if 0
728     double sleeptime = smpi_or(req->real_size);
729     if(sleeptime > 0.0){
730         simcall_process_sleep(sleeptime);
731         XBT_DEBUG("receiving size of %zu : sleep %f ", req->real_size, sleeptime);
732     }
733     smpi_mpi_request_free(&(req->detached_sender));
734   }
735   if(req->flags & PERSISTENT)
736     req->action = nullptr;
737   req->flags |= FINISHED;
738
739   smpi_mpi_request_free(request);
740 }
741
742 int smpi_mpi_test(MPI_Request * request, MPI_Status * status) {
743   //assume that request is not MPI_REQUEST_NULL (filtered in PMPI_Test or smpi_mpi_testall before)
744
745   // to avoid deadlocks if used as a break condition, such as
746   //     while (MPI_Test(request, flag, status) && flag) {
747   //     }
748   // because the time will not normally advance when only calls to MPI_Test are made -> deadlock
749   // multiplier to the sleeptime, to increase speed of execution, each failed test will increase it
750   static int nsleeps = 1;
751   if(smpi_test_sleep > 0)  
752     simcall_process_sleep(nsleeps*smpi_test_sleep);
753
754   smpi_empty_status(status);
755   int flag = 1;
756   if (((*request)->flags & PREPARED) == 0) {
757     if ((*request)->action != nullptr)
758       flag = simcall_comm_test((*request)->action);
759     if (flag) {
760       finish_wait(request, status);
761       nsleeps=1;//reset the number of sleeps we will do next time
762       if (*request != MPI_REQUEST_NULL && ((*request)->flags & PERSISTENT)==0)
763       *request = MPI_REQUEST_NULL;
764     }else{
765       nsleeps++;
766     }
767   }
768   return flag;
769 }
770
771 int smpi_mpi_testany(int count, MPI_Request requests[], int *index, MPI_Status * status)
772 {
773   std::vector<simgrid::kernel::activity::ActivityImpl*> comms;
774   comms.reserve(count);
775
776   int i;
777   int flag = 0;
778
779   *index = MPI_UNDEFINED;
780
781   std::vector<int> map; /** Maps all matching comms back to their location in requests **/
782   for(i = 0; i < count; i++) {
783     if ((requests[i] != MPI_REQUEST_NULL) && requests[i]->action && !(requests[i]->flags & PREPARED)) {
784        comms.push_back(requests[i]->action);
785        map.push_back(i);
786     }
787   }
788   if(!map.empty()) {
789     //multiplier to the sleeptime, to increase speed of execution, each failed testany will increase it
790     static int nsleeps = 1;
791     if(smpi_test_sleep > 0) 
792       simcall_process_sleep(nsleeps*smpi_test_sleep);
793
794     i = simcall_comm_testany(comms.data(), comms.size()); // The i-th element in comms matches!
795     if (i != -1) { // -1 is not MPI_UNDEFINED but a SIMIX return code. (nothing matches)
796       *index = map[i]; 
797       finish_wait(&requests[*index], status);
798       flag             = 1;
799       nsleeps          = 1;
800       if (requests[*index] != MPI_REQUEST_NULL && (requests[*index]->flags & NON_PERSISTENT)) {
801         requests[*index] = MPI_REQUEST_NULL;
802       }
803     } else {
804       nsleeps++;
805     }
806   } else {
807       //all requests are null or inactive, return true
808       flag = 1;
809       smpi_empty_status(status);
810   }
811
812   return flag;
813 }
814
815 int smpi_mpi_testall(int count, MPI_Request requests[], MPI_Status status[])
816 {
817   MPI_Status stat;
818   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
819   int flag=1;
820   for(int i=0; i<count; i++){
821     if (requests[i] != MPI_REQUEST_NULL && !(requests[i]->flags & PREPARED)) {
822       if (smpi_mpi_test(&requests[i], pstat)!=1){
823         flag=0;
824       }else{
825           requests[i]=MPI_REQUEST_NULL;
826       }
827     }else{
828       smpi_empty_status(pstat);
829     }
830     if(status != MPI_STATUSES_IGNORE) {
831       status[i] = *pstat;
832     }
833   }
834   return flag;
835 }
836
837 void smpi_mpi_probe(int source, int tag, MPI_Comm comm, MPI_Status* status){
838   int flag=0;
839   //FIXME find another way to avoid busy waiting ?
840   // the issue here is that we have to wait on a nonexistent comm
841   while(flag==0){
842     smpi_mpi_iprobe(source, tag, comm, &flag, status);
843     XBT_DEBUG("Busy Waiting on probing : %d", flag);
844   }
845 }
846
847 void smpi_mpi_iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status){
848
849   MPI_Request request = build_request(nullptr, 0, MPI_CHAR, source == MPI_ANY_SOURCE ? MPI_ANY_SOURCE :
850                  smpi_group_index(smpi_comm_group(comm), source), smpi_comm_rank(comm), tag, comm, PERSISTENT | RECV);
851
852   // to avoid deadlock, we have to sleep some time here, or the timer won't advance and we will only do iprobe simcalls
853   // (especially when used as a break condition, such as while(MPI_Iprobe(...)) ... )
854   // multiplier to the sleeptime, to increase speed of execution, each failed iprobe will increase it
855   static int nsleeps = 1;
856   if(smpi_iprobe_sleep > 0)  
857     simcall_process_sleep(nsleeps*smpi_iprobe_sleep);
858   // behave like a receive, but don't do it
859   smx_mailbox_t mailbox;
860
861   print_request("New iprobe", request);
862   // We have to test both mailboxes as we don't know if we will receive one one or another
863   if (xbt_cfg_get_int("smpi/async-small-thresh") > 0){
864       mailbox = smpi_process_mailbox_small();
865       XBT_DEBUG("Trying to probe the perm recv mailbox");
866       request->action = simcall_comm_iprobe(mailbox, 0, request->src, request->tag, &match_recv, static_cast<void*>(request));
867   }
868
869   if (request->action == nullptr){
870     mailbox = smpi_process_mailbox();
871     XBT_DEBUG("trying to probe the other mailbox");
872     request->action = simcall_comm_iprobe(mailbox, 0, request->src,request->tag, &match_recv, static_cast<void*>(request));
873   }
874
875   if (request->action != nullptr){
876     simgrid::kernel::activity::Comm *sync_comm = static_cast<simgrid::kernel::activity::Comm*>(request->action);
877     MPI_Request req                            = static_cast<MPI_Request>(sync_comm->src_data);
878     *flag = 1;
879     if(status != MPI_STATUS_IGNORE && (req->flags & PREPARED) == 0) {
880       status->MPI_SOURCE = smpi_group_rank(smpi_comm_group(comm), req->src);
881       status->MPI_TAG    = req->tag;
882       status->MPI_ERROR  = MPI_SUCCESS;
883       status->count      = req->real_size;
884     }
885     nsleeps = 1;//reset the number of sleeps we will do next time
886   }
887   else {
888     *flag = 0;
889     nsleeps++;
890   }
891   smpi_mpi_request_free(&request);
892
893   return;
894 }
895
896 void smpi_mpi_wait(MPI_Request * request, MPI_Status * status)
897 {
898   print_request("Waiting", *request);
899   if ((*request)->flags & PREPARED) {
900     smpi_empty_status(status);
901     return;
902   }
903
904   if ((*request)->action != nullptr)
905     // this is not a detached send
906     simcall_comm_wait((*request)->action, -1.0);
907
908   finish_wait(request, status);
909   if (*request != MPI_REQUEST_NULL && (((*request)->flags & NON_PERSISTENT)!=0))
910       *request = MPI_REQUEST_NULL;
911 }
912
913 int smpi_mpi_waitany(int count, MPI_Request requests[], MPI_Status * status)
914 {
915   xbt_dynar_t comms;
916   int i;
917   int size = 0;
918   int index = MPI_UNDEFINED;
919   int *map;
920
921   if(count > 0) {
922     // Wait for a request to complete
923     comms = xbt_dynar_new(sizeof(smx_activity_t), nullptr);
924     map = xbt_new(int, count);
925     XBT_DEBUG("Wait for one of %d", count);
926     for(i = 0; i < count; i++) {
927       if (requests[i] != MPI_REQUEST_NULL && !(requests[i]->flags & PREPARED) && !(requests[i]->flags & FINISHED)) {
928         if (requests[i]->action != nullptr) {
929           XBT_DEBUG("Waiting any %p ", requests[i]);
930           xbt_dynar_push(comms, &requests[i]->action);
931           map[size] = i;
932           size++;
933         }else{
934          //This is a finished detached request, let's return this one
935          size=0;//so we free the dynar but don't do the waitany call
936          index=i;
937          finish_wait(&requests[i], status);//cleanup if refcount = 0
938          if (requests[i] != MPI_REQUEST_NULL && (requests[i]->flags & NON_PERSISTENT))
939          requests[i]=MPI_REQUEST_NULL;//set to null
940          break;
941          }
942       }
943     }
944     if(size > 0) {
945       i = simcall_comm_waitany(comms, -1);
946
947       // not MPI_UNDEFINED, as this is a simix return code
948       if (i != -1) {
949         index = map[i];
950         finish_wait(&requests[index], status);
951         if (requests[i] != MPI_REQUEST_NULL && (requests[i]->flags & NON_PERSISTENT))
952         requests[index] = MPI_REQUEST_NULL;
953       }
954     }
955     xbt_free(map);
956     xbt_dynar_free(&comms);
957   }
958
959   if (index==MPI_UNDEFINED)
960     smpi_empty_status(status);
961
962   return index;
963 }
964
965 int smpi_mpi_waitall(int count, MPI_Request requests[], MPI_Status status[])
966 {
967   int  index, c;
968   MPI_Status stat;
969   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
970   int retvalue = MPI_SUCCESS;
971   //tag invalid requests in the set
972   if (status != MPI_STATUSES_IGNORE) {
973     for (c = 0; c < count; c++) {
974       if (requests[c] == MPI_REQUEST_NULL || requests[c]->dst == MPI_PROC_NULL || (requests[c]->flags & PREPARED)) {
975         smpi_empty_status(&status[c]);
976       } else if (requests[c]->src == MPI_PROC_NULL) {
977         smpi_empty_status(&status[c]);
978         status[c].MPI_SOURCE = MPI_PROC_NULL;
979       }
980     }
981   }
982   for(c = 0; c < count; c++) {
983
984     if (MC_is_active() || MC_record_replay_is_active()) {
985       smpi_mpi_wait(&requests[c], pstat);
986       index = c;
987     } else {
988       index = smpi_mpi_waitany(count, requests, pstat);
989       if (index == MPI_UNDEFINED)
990         break;
991       if (requests[index] != MPI_REQUEST_NULL && (requests[index]->flags & NON_PERSISTENT))
992       requests[index]=MPI_REQUEST_NULL;
993     }
994     if (status != MPI_STATUSES_IGNORE) {
995       status[index] = *pstat;
996       if (status[index].MPI_ERROR == MPI_ERR_TRUNCATE)
997         retvalue = MPI_ERR_IN_STATUS;
998     }
999   }
1000
1001   return retvalue;
1002 }
1003
1004 int smpi_mpi_waitsome(int incount, MPI_Request requests[], int *indices, MPI_Status status[])
1005 {
1006   int i;
1007   int count = 0;
1008   int index;
1009   MPI_Status stat;
1010   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
1011
1012   for(i = 0; i < incount; i++)
1013   {
1014     index=smpi_mpi_waitany(incount, requests, pstat);
1015     if(index!=MPI_UNDEFINED){
1016       indices[count] = index;
1017       count++;
1018       if(status != MPI_STATUSES_IGNORE) {
1019         status[index] = *pstat;
1020       }
1021      if (requests[index] != MPI_REQUEST_NULL && (requests[index]->flags & NON_PERSISTENT))
1022      requests[index]=MPI_REQUEST_NULL;
1023     }else{
1024       return MPI_UNDEFINED;
1025     }
1026   }
1027   return count;
1028 }
1029
1030 int smpi_mpi_testsome(int incount, MPI_Request requests[], int *indices, MPI_Status status[])
1031 {
1032   int i;
1033   int count = 0;
1034   int count_dead = 0;
1035   MPI_Status stat;
1036   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
1037
1038   for(i = 0; i < incount; i++) {
1039     if((requests[i] != MPI_REQUEST_NULL)) {
1040       if(smpi_mpi_test(&requests[i], pstat)) {
1041          indices[i] = 1;
1042          count++;
1043          if(status != MPI_STATUSES_IGNORE) {
1044            status[i] = *pstat;
1045          }
1046          if ((requests[i] != MPI_REQUEST_NULL) && requests[i]->flags & NON_PERSISTENT)
1047          requests[i]=MPI_REQUEST_NULL;
1048       }
1049     }else{
1050       count_dead++;
1051     }
1052   }
1053   if(count_dead==incount)
1054     return MPI_UNDEFINED;
1055   else return count;
1056 }
1057
1058 void smpi_mpi_bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1059 {
1060     smpi_coll_tuned_bcast_binomial_tree(buf, count, datatype, root, comm);
1061 }
1062
1063 void smpi_mpi_barrier(MPI_Comm comm)
1064 {
1065     smpi_coll_tuned_barrier_ompi_basic_linear(comm);
1066 }
1067
1068 void smpi_mpi_gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1069                      void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
1070 {
1071   int system_tag = COLL_TAG_GATHER;
1072   int rank, size, src, index;
1073   MPI_Aint lb = 0, recvext = 0;
1074   MPI_Request *requests;
1075
1076   rank = smpi_comm_rank(comm);
1077   size = smpi_comm_size(comm);
1078   if(rank != root) {
1079     // Send buffer to root
1080     smpi_mpi_send(sendbuf, sendcount, sendtype, root, system_tag, comm);
1081   } else {
1082     smpi_datatype_extent(recvtype, &lb, &recvext);
1083     // Local copy from root
1084     smpi_datatype_copy(sendbuf, sendcount, sendtype, static_cast<char*>(recvbuf) + root * recvcount * recvext, recvcount, recvtype);
1085     // Receive buffers from senders
1086     requests = xbt_new(MPI_Request, size - 1);
1087     index = 0;
1088     for(src = 0; src < size; src++) {
1089       if(src != root) {
1090         requests[index] = smpi_irecv_init(static_cast<char*>(recvbuf) + src * recvcount * recvext, recvcount, recvtype,
1091                                           src, system_tag, comm);
1092         index++;
1093       }
1094     }
1095     // Wait for completion of irecv's.
1096     smpi_mpi_startall(size - 1, requests);
1097     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
1098     for(src = 0; src < size-1; src++) {
1099       smpi_mpi_request_free(&requests[src]);
1100     }
1101     xbt_free(requests);
1102   }
1103 }
1104
1105 void smpi_mpi_reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts, MPI_Datatype datatype, MPI_Op op,
1106                              MPI_Comm comm)
1107 {
1108     int i, size, count;
1109     int *displs;
1110     int rank = smpi_comm_rank(comm);
1111     void *tmpbuf;
1112
1113     /* arbitrarily choose root as rank 0 */
1114     size = smpi_comm_size(comm);
1115     count = 0;
1116     displs = xbt_new(int, size);
1117     for (i = 0; i < size; i++) {
1118       displs[i] = count;
1119       count += recvcounts[i];
1120     }
1121     tmpbuf=static_cast<void*>(smpi_get_tmp_sendbuffer(count*smpi_datatype_get_extent(datatype)));
1122
1123     mpi_coll_reduce_fun(sendbuf, tmpbuf, count, datatype, op, 0, comm);
1124     smpi_mpi_scatterv(tmpbuf, recvcounts, displs, datatype, recvbuf, recvcounts[rank], datatype, 0, comm);
1125     xbt_free(displs);
1126     smpi_free_tmp_buffer(tmpbuf);
1127 }
1128
1129 void smpi_mpi_gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int *recvcounts, int *displs,
1130                       MPI_Datatype recvtype, int root, MPI_Comm comm)
1131 {
1132   int system_tag = COLL_TAG_GATHERV;
1133   int rank, size, src, index;
1134   MPI_Aint lb = 0, recvext = 0;
1135   MPI_Request *requests;
1136
1137   rank = smpi_comm_rank(comm);
1138   size = smpi_comm_size(comm);
1139   if(rank != root) {
1140     // Send buffer to root
1141     smpi_mpi_send(sendbuf, sendcount, sendtype, root, system_tag, comm);
1142   } else {
1143     smpi_datatype_extent(recvtype, &lb, &recvext);
1144     // Local copy from root
1145     smpi_datatype_copy(sendbuf, sendcount, sendtype, static_cast<char*>(recvbuf) + displs[root] * recvext,
1146                        recvcounts[root], recvtype);
1147     // Receive buffers from senders
1148     requests = xbt_new(MPI_Request, size - 1);
1149     index = 0;
1150     for(src = 0; src < size; src++) {
1151       if(src != root) {
1152         requests[index] = smpi_irecv_init(static_cast<char*>(recvbuf) + displs[src] * recvext,
1153                           recvcounts[src], recvtype, src, system_tag, comm);
1154         index++;
1155       }
1156     }
1157     // Wait for completion of irecv's.
1158     smpi_mpi_startall(size - 1, requests);
1159     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
1160     for(src = 0; src < size-1; src++) {
1161       smpi_mpi_request_free(&requests[src]);
1162     }
1163     xbt_free(requests);
1164   }
1165 }
1166
1167 void smpi_mpi_allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1168                         void *recvbuf,int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
1169 {
1170   int system_tag = COLL_TAG_ALLGATHER;
1171   int rank, size, other, index;
1172   MPI_Aint lb = 0, recvext = 0;
1173   MPI_Request *requests;
1174
1175   rank = smpi_comm_rank(comm);
1176   size = smpi_comm_size(comm);
1177   // FIXME: check for errors
1178   smpi_datatype_extent(recvtype, &lb, &recvext);
1179   // Local copy from self
1180   smpi_datatype_copy(sendbuf, sendcount, sendtype, static_cast<char *>(recvbuf) + rank * recvcount * recvext, recvcount, recvtype);
1181   // Send/Recv buffers to/from others;
1182   requests = xbt_new(MPI_Request, 2 * (size - 1));
1183   index = 0;
1184   for(other = 0; other < size; other++) {
1185     if(other != rank) {
1186       requests[index] = smpi_isend_init(sendbuf, sendcount, sendtype, other, system_tag,comm);
1187       index++;
1188       requests[index] = smpi_irecv_init(static_cast<char *>(recvbuf) + other * recvcount * recvext, recvcount, recvtype, other,
1189                                         system_tag, comm);
1190       index++;
1191     }
1192   }
1193   // Wait for completion of all comms.
1194   smpi_mpi_startall(2 * (size - 1), requests);
1195   smpi_mpi_waitall(2 * (size - 1), requests, MPI_STATUS_IGNORE);
1196   for(other = 0; other < 2*(size-1); other++) {
1197     smpi_mpi_request_free(&requests[other]);
1198   }
1199   xbt_free(requests);
1200 }
1201
1202 void smpi_mpi_allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf,
1203                          int *recvcounts, int *displs, MPI_Datatype recvtype, MPI_Comm comm)
1204 {
1205   int system_tag = COLL_TAG_ALLGATHERV;
1206   int rank, size, other, index;
1207   MPI_Aint lb = 0, recvext = 0;
1208   MPI_Request *requests;
1209
1210   rank = smpi_comm_rank(comm);
1211   size = smpi_comm_size(comm);
1212   smpi_datatype_extent(recvtype, &lb, &recvext);
1213   // Local copy from self
1214   smpi_datatype_copy(sendbuf, sendcount, sendtype, static_cast<char *>(recvbuf) + displs[rank] * recvext,recvcounts[rank], recvtype);
1215   // Send buffers to others;
1216   requests = xbt_new(MPI_Request, 2 * (size - 1));
1217   index = 0;
1218   for(other = 0; other < size; other++) {
1219     if(other != rank) {
1220       requests[index] =
1221         smpi_isend_init(sendbuf, sendcount, sendtype, other, system_tag, comm);
1222       index++;
1223       requests[index] = smpi_irecv_init(static_cast<char *>(recvbuf) + displs[other] * recvext, recvcounts[other],
1224                           recvtype, other, system_tag, comm);
1225       index++;
1226     }
1227   }
1228   // Wait for completion of all comms.
1229   smpi_mpi_startall(2 * (size - 1), requests);
1230   smpi_mpi_waitall(2 * (size - 1), requests, MPI_STATUS_IGNORE);
1231   for(other = 0; other < 2*(size-1); other++) {
1232     smpi_mpi_request_free(&requests[other]);
1233   }
1234   xbt_free(requests);
1235 }
1236
1237 void smpi_mpi_scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1238                       void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
1239 {
1240   int system_tag = COLL_TAG_SCATTER;
1241   int dst;
1242   MPI_Aint lb = 0;
1243   MPI_Aint sendext = 0;
1244   MPI_Request *requests;
1245
1246   int rank = smpi_comm_rank(comm);
1247   int size = smpi_comm_size(comm);
1248   if(rank != root) {
1249     // Recv buffer from root
1250     smpi_mpi_recv(recvbuf, recvcount, recvtype, root, system_tag, comm, MPI_STATUS_IGNORE);
1251   } else {
1252     smpi_datatype_extent(sendtype, &lb, &sendext);
1253     // Local copy from root
1254     if(recvbuf!=MPI_IN_PLACE){
1255         smpi_datatype_copy(static_cast<char *>(sendbuf) + root * sendcount * sendext,
1256                            sendcount, sendtype, recvbuf, recvcount, recvtype);
1257     }
1258     // Send buffers to receivers
1259     requests = xbt_new(MPI_Request, size - 1);
1260     int index = 0;
1261     for(dst = 0; dst < size; dst++) {
1262       if(dst != root) {
1263         requests[index] = smpi_isend_init(static_cast<char *>(sendbuf) + dst * sendcount * sendext, sendcount, sendtype, dst,
1264                                           system_tag, comm);
1265         index++;
1266       }
1267     }
1268     // Wait for completion of isend's.
1269     smpi_mpi_startall(size - 1, requests);
1270     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
1271     for(dst = 0; dst < size-1; dst++) {
1272       smpi_mpi_request_free(&requests[dst]);
1273     }
1274     xbt_free(requests);
1275   }
1276 }
1277
1278 void smpi_mpi_scatterv(void *sendbuf, int *sendcounts, int *displs, MPI_Datatype sendtype, void *recvbuf, int recvcount,
1279                        MPI_Datatype recvtype, int root, MPI_Comm comm)
1280 {
1281   int system_tag = COLL_TAG_SCATTERV;
1282   int dst;
1283   MPI_Aint lb = 0;
1284   MPI_Aint sendext = 0;
1285   MPI_Request *requests;
1286
1287   int rank = smpi_comm_rank(comm);
1288   int size = smpi_comm_size(comm);
1289   if(rank != root) {
1290     // Recv buffer from root
1291     smpi_mpi_recv(recvbuf, recvcount, recvtype, root, system_tag, comm, MPI_STATUS_IGNORE);
1292   } else {
1293     smpi_datatype_extent(sendtype, &lb, &sendext);
1294     // Local copy from root
1295     if(recvbuf!=MPI_IN_PLACE){
1296       smpi_datatype_copy(static_cast<char *>(sendbuf) + displs[root] * sendext, sendcounts[root],
1297                        sendtype, recvbuf, recvcount, recvtype);
1298     }
1299     // Send buffers to receivers
1300     requests = xbt_new(MPI_Request, size - 1);
1301     int index = 0;
1302     for(dst = 0; dst < size; dst++) {
1303       if(dst != root) {
1304         requests[index] = smpi_isend_init(static_cast<char *>(sendbuf) + displs[dst] * sendext, sendcounts[dst],
1305                             sendtype, dst, system_tag, comm);
1306         index++;
1307       }
1308     }
1309     // Wait for completion of isend's.
1310     smpi_mpi_startall(size - 1, requests);
1311     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
1312     for(dst = 0; dst < size-1; dst++) {
1313       smpi_mpi_request_free(&requests[dst]);
1314     }
1315     xbt_free(requests);
1316   }
1317 }
1318
1319 void smpi_mpi_reduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root,
1320                      MPI_Comm comm)
1321 {
1322   int system_tag = COLL_TAG_REDUCE;
1323   int src, index;
1324   MPI_Aint lb = 0;
1325   MPI_Aint dataext = 0;
1326   MPI_Request *requests;
1327   void **tmpbufs;
1328
1329   char* sendtmpbuf = static_cast<char *>(sendbuf);
1330
1331   int rank = smpi_comm_rank(comm);
1332   int size = smpi_comm_size(comm);
1333   //non commutative case, use a working algo from openmpi
1334   if(!smpi_op_is_commute(op)){
1335     smpi_coll_tuned_reduce_ompi_basic_linear(sendtmpbuf, recvbuf, count, datatype, op, root, comm);
1336     return;
1337   }
1338
1339   if( sendbuf == MPI_IN_PLACE ) {
1340     sendtmpbuf = static_cast<char *>(smpi_get_tmp_sendbuffer(count*smpi_datatype_get_extent(datatype)));
1341     smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
1342   }
1343   
1344   if(rank != root) {
1345     // Send buffer to root
1346     smpi_mpi_send(sendtmpbuf, count, datatype, root, system_tag, comm);
1347   } else {
1348     smpi_datatype_extent(datatype, &lb, &dataext);
1349     // Local copy from root
1350     if (sendtmpbuf != nullptr && recvbuf != nullptr)
1351       smpi_datatype_copy(sendtmpbuf, count, datatype, recvbuf, count, datatype);
1352     // Receive buffers from senders
1353     requests = xbt_new(MPI_Request, size - 1);
1354     tmpbufs = xbt_new(void *, size - 1);
1355     index = 0;
1356     for(src = 0; src < size; src++) {
1357       if(src != root) {
1358          if (!smpi_process_get_replaying())
1359           tmpbufs[index] = xbt_malloc(count * dataext);
1360          else
1361            tmpbufs[index] = smpi_get_tmp_sendbuffer(count * dataext);
1362         requests[index] =
1363           smpi_irecv_init(tmpbufs[index], count, datatype, src, system_tag, comm);
1364         index++;
1365       }
1366     }
1367     // Wait for completion of irecv's.
1368     smpi_mpi_startall(size - 1, requests);
1369     for(src = 0; src < size - 1; src++) {
1370       index = smpi_mpi_waitany(size - 1, requests, MPI_STATUS_IGNORE);
1371       XBT_DEBUG("finished waiting any request with index %d", index);
1372       if(index == MPI_UNDEFINED) {
1373         break;
1374       }else{
1375         smpi_mpi_request_free(&requests[index]);
1376       }
1377       if(op) /* op can be MPI_OP_NULL that does nothing */
1378         smpi_op_apply(op, tmpbufs[index], recvbuf, &count, &datatype);
1379     }
1380       for(index = 0; index < size - 1; index++) {
1381         smpi_free_tmp_buffer(tmpbufs[index]);
1382       }
1383     xbt_free(tmpbufs);
1384     xbt_free(requests);
1385
1386   }
1387   if( sendbuf == MPI_IN_PLACE ) {
1388     smpi_free_tmp_buffer(sendtmpbuf);
1389   }
1390 }
1391
1392 void smpi_mpi_allreduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1393 {
1394   smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1395   smpi_mpi_bcast(recvbuf, count, datatype, 0, comm);
1396 }
1397
1398 void smpi_mpi_scan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1399 {
1400   int system_tag = -888;
1401   int other, index;
1402   MPI_Aint lb = 0, dataext = 0;
1403   MPI_Request *requests;
1404   void **tmpbufs;
1405
1406   int rank = smpi_comm_rank(comm);
1407   int size = smpi_comm_size(comm);
1408
1409   smpi_datatype_extent(datatype, &lb, &dataext);
1410
1411   // Local copy from self
1412   smpi_datatype_copy(sendbuf, count, datatype, recvbuf, count, datatype);
1413
1414   // Send/Recv buffers to/from others;
1415   requests = xbt_new(MPI_Request, size - 1);
1416   tmpbufs = xbt_new(void *, rank);
1417   index = 0;
1418   for(other = 0; other < rank; other++) {
1419     tmpbufs[index] = smpi_get_tmp_sendbuffer(count * dataext);
1420     requests[index] = smpi_irecv_init(tmpbufs[index], count, datatype, other, system_tag, comm);
1421     index++;
1422   }
1423   for(other = rank + 1; other < size; other++) {
1424     requests[index] = smpi_isend_init(sendbuf, count, datatype, other, system_tag, comm);
1425     index++;
1426   }
1427   // Wait for completion of all comms.
1428   smpi_mpi_startall(size - 1, requests);
1429
1430   if(smpi_op_is_commute(op)){
1431     for(other = 0; other < size - 1; other++) {
1432       index = smpi_mpi_waitany(size - 1, requests, MPI_STATUS_IGNORE);
1433       if(index == MPI_UNDEFINED) {
1434         break;
1435       }
1436       if(index < rank) {
1437         // #Request is below rank: it's a irecv
1438         smpi_op_apply(op, tmpbufs[index], recvbuf, &count, &datatype);
1439       }
1440     }
1441   }else{
1442     //non commutative case, wait in order
1443     for(other = 0; other < size - 1; other++) {
1444       smpi_mpi_wait(&(requests[other]), MPI_STATUS_IGNORE);
1445       if(index < rank) {
1446         smpi_op_apply(op, tmpbufs[other], recvbuf, &count, &datatype);
1447       }
1448     }
1449   }
1450   for(index = 0; index < rank; index++) {
1451     smpi_free_tmp_buffer(tmpbufs[index]);
1452   }
1453   for(index = 0; index < size-1; index++) {
1454     smpi_mpi_request_free(&requests[index]);
1455   }
1456   xbt_free(tmpbufs);
1457   xbt_free(requests);
1458 }
1459
1460 void smpi_mpi_exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1461 {
1462   int system_tag = -888;
1463   int other, index;
1464   MPI_Aint lb = 0, dataext = 0;
1465   MPI_Request *requests;
1466   void **tmpbufs;
1467   int recvbuf_is_empty=1;
1468   int rank = smpi_comm_rank(comm);
1469   int size = smpi_comm_size(comm);
1470
1471   smpi_datatype_extent(datatype, &lb, &dataext);
1472
1473   // Send/Recv buffers to/from others;
1474   requests = xbt_new(MPI_Request, size - 1);
1475   tmpbufs = xbt_new(void *, rank);
1476   index = 0;
1477   for(other = 0; other < rank; other++) {
1478     tmpbufs[index] = smpi_get_tmp_sendbuffer(count * dataext);
1479     requests[index] = smpi_irecv_init(tmpbufs[index], count, datatype, other, system_tag, comm);
1480     index++;
1481   }
1482   for(other = rank + 1; other < size; other++) {
1483     requests[index] = smpi_isend_init(sendbuf, count, datatype, other, system_tag, comm);
1484     index++;
1485   }
1486   // Wait for completion of all comms.
1487   smpi_mpi_startall(size - 1, requests);
1488   if(smpi_op_is_commute(op)){
1489     for(other = 0; other < size - 1; other++) {
1490       index = smpi_mpi_waitany(size - 1, requests, MPI_STATUS_IGNORE);
1491       if(index == MPI_UNDEFINED) {
1492         break;
1493       }
1494       if(index < rank) {
1495         if(recvbuf_is_empty){
1496           smpi_datatype_copy(tmpbufs[index], count, datatype, recvbuf, count, datatype);
1497           recvbuf_is_empty=0;
1498         } else
1499           // #Request is below rank: it's a irecv
1500           smpi_op_apply(op, tmpbufs[index], recvbuf, &count, &datatype);
1501       }
1502     }
1503   }else{
1504     //non commutative case, wait in order
1505     for(other = 0; other < size - 1; other++) {
1506       smpi_mpi_wait(&(requests[other]), MPI_STATUS_IGNORE);
1507       if(index < rank) {
1508           if(recvbuf_is_empty){
1509             smpi_datatype_copy(tmpbufs[other], count, datatype, recvbuf, count, datatype);
1510             recvbuf_is_empty=0;
1511           } else
1512             smpi_op_apply(op, tmpbufs[other], recvbuf, &count, &datatype);
1513       }
1514     }
1515   }
1516   for(index = 0; index < rank; index++) {
1517     smpi_free_tmp_buffer(tmpbufs[index]);
1518   }
1519   for(index = 0; index < size-1; index++) {
1520     smpi_mpi_request_free(&requests[index]);
1521   }
1522   xbt_free(tmpbufs);
1523   xbt_free(requests);
1524 }