Logo AND Algorithmique Numérique Distribuée

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