Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
tiny doc formating improvement
[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   long 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 int factor_cmp(const s_smpi_factor_multival_t& pa, const s_smpi_factor_multival_t& pb)
94 {
95   return (pa.factor > pb.factor) ? 1 :
96          (pa.factor < pb.factor) ? -1 : 0;
97 }
98
99 static std::vector<s_smpi_factor_multival_t> parse_factor(const char *smpi_coef_string)
100 {
101   std::vector<s_smpi_factor_multival_t> smpi_factor;
102   s_smpi_factor_multival_t fact;
103
104   /** Setup the tokenizer that parses the string **/
105   typedef boost::tokenizer<boost::char_separator<char>> Tokenizer;
106   boost::char_separator<char> sep(";");
107   boost::char_separator<char> factor_separator(":");
108   std::string tmp_string(smpi_coef_string);
109   Tokenizer tokens(tmp_string, sep);
110
111   /** 
112    * Iterate over patterns like A:B:C:D;E:F;G:H
113    * These will be broken down into:
114    * A --> B, C, D
115    * E --> F
116    * G --> H
117    */
118   for (Tokenizer::iterator token_iter = tokens.begin();
119          token_iter != tokens.end(); token_iter++) {
120     Tokenizer factor_values(*token_iter, factor_separator);
121
122     if (factor_values.begin() == factor_values.end()) {
123       xbt_die("Malformed radical for smpi factor: '%s'", smpi_coef_string);
124     }
125     unsigned int iteration = 0;
126     for (Tokenizer::iterator factor_iter = factor_values.begin();
127          factor_iter != factor_values.end(); factor_iter++, iteration++) {
128       char *errmsg;
129
130       if (factor_iter == factor_values.begin()) { /* first element */
131         errmsg = bprintf("Invalid factor in chunk #%zu: %%s", smpi_factor.size()+1);
132         fact.factor = xbt_str_parse_int(factor_iter->c_str(), errmsg);
133       }
134       else {
135         errmsg = bprintf("Invalid factor value %d in chunk #%zu: %%s", iteration, smpi_factor.size()+1);
136         fact.values.push_back(xbt_str_parse_double((*factor_iter).c_str(), errmsg));
137       }
138       xbt_free(errmsg);
139     }
140
141     smpi_factor.push_back(fact);
142     XBT_DEBUG("smpi_factor:\t%ld : %zu values, first: %f", fact.factor, smpi_factor.size(), fact.values[0]);
143   }
144   std::sort(smpi_factor.begin(), smpi_factor.end(), &factor_cmp);
145   for (auto& fact : smpi_factor) {
146     XBT_DEBUG("smpi_factor:\t%ld : %zu values, first: %f", fact.factor, smpi_factor.size() ,fact.values[0]);
147   }
148
149   return smpi_factor;
150 }
151
152 static double smpi_os(double size)
153 {
154   if (smpi_os_values.empty()) {
155     smpi_os_values = parse_factor(xbt_cfg_get_string("smpi/os"));
156   }
157   double current=0.0;
158   // Iterate over all the sections that were specified and find the right
159   // value. (fact.factor represents the interval sizes; we want to find the
160   // section that has fact.factor <= size and no other such fact.factor <= size)
161   // Note: parse_factor() (used before) already sorts the dynar we iterate over!
162   for (auto& fact : smpi_os_values) {
163     if (size <= fact.factor) { // Values already too large, use the previously
164                                // computed value of current!
165         XBT_DEBUG("os : %f <= %ld return %f", size, fact.factor, current);
166       return current;
167     }else{
168       // If the next section is too large, the current section must be used.
169       // Hence, save the cost, as we might have to use it.
170       current = fact.values[0]+fact.values[1]*size;
171     }
172   }
173   XBT_DEBUG("Searching for smpi/os: %f is larger than the largest boundary, return %f", size, current);
174
175   return current;
176 }
177
178 static double smpi_ois(double size)
179 {
180   if (smpi_ois_values.empty()) {
181     smpi_ois_values = parse_factor(xbt_cfg_get_string("smpi/ois"));
182   }
183   double current=0.0;
184   // Iterate over all the sections that were specified and find the right value. (fact.factor represents the interval
185   // sizes; we want to find the section that has fact.factor <= size and no other such fact.factor <= size)
186   // Note: parse_factor() (used before) already sorts the dynar we iterate over!
187   for (auto& fact : smpi_ois_values) {
188     if (size <= fact.factor) { // Values already too large, use the previously  computed value of current!
189         XBT_DEBUG("ois : %f <= %ld return %f", size, fact.factor, current);
190       return current;
191     }else{
192       // If the next section is too large, the current section must be used.
193       // Hence, save the cost, as we might have to use it.
194       current = fact.values[0]+fact.values[1]*size;
195     }
196   }
197   XBT_DEBUG("Searching for smpi/ois: %f is larger than the largest boundary, return %f", size, current);
198
199   return current;
200 }
201
202 static double smpi_or(double size)
203 {
204   if (smpi_or_values.empty()) {
205     smpi_or_values = parse_factor(xbt_cfg_get_string("smpi/or"));
206   }
207   double current=0.0;
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 : %f <= %ld return %f", 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: %f is larger than largest boundary, return %f", 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   xbt_dynar_t comms;
777   int i;
778   int flag = 0;
779
780   *index = MPI_UNDEFINED;
781   comms = xbt_dynar_new(sizeof(smx_synchro_t), nullptr);
782   std::vector<int> map; /** Maps all matching comms back to their location in requests **/
783   for(i = 0; i < count; i++) {
784     if ((requests[i] != MPI_REQUEST_NULL) && requests[i]->action && !(requests[i]->flags & PREPARED)) {
785        xbt_dynar_push(comms, &requests[i]->action);
786        map.push_back(i);
787     }
788   }
789   if(!map.empty()) {
790     //multiplier to the sleeptime, to increase speed of execution, each failed testany will increase it
791     static int nsleeps = 1;
792     if(smpi_test_sleep > 0) 
793       simcall_process_sleep(nsleeps*smpi_test_sleep);
794
795     i = simcall_comm_testany(comms); // The i-th element in comms matches!
796     if (i != -1) { // -1 is not MPI_UNDEFINED but a SIMIX return code. (nothing matches)
797       *index = map[i]; 
798       finish_wait(&requests[*index], status);
799       flag             = 1;
800       nsleeps          = 1;
801       if (requests[*index] != MPI_REQUEST_NULL && (requests[*index]->flags & NON_PERSISTENT)) {
802         requests[*index] = MPI_REQUEST_NULL;
803       }
804     } else {
805       nsleeps++;
806     }
807   } else {
808       flag = 0; // all requests are null or inactive, return false
809       smpi_empty_status(status);
810   }
811   xbt_dynar_free(&comms);
812
813   return flag;
814 }
815
816 int smpi_mpi_testall(int count, MPI_Request requests[], MPI_Status status[])
817 {
818   MPI_Status stat;
819   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
820   int flag=1;
821   int i;
822   for(i=0; i<count; i++){
823     if (requests[i] != MPI_REQUEST_NULL && !(requests[i]->flags & PREPARED)) {
824       if (smpi_mpi_test(&requests[i], pstat)!=1){
825         flag=0;
826       }else{
827           requests[i]=MPI_REQUEST_NULL;
828       }
829     }else{
830       smpi_empty_status(pstat);
831     }
832     if(status != MPI_STATUSES_IGNORE) {
833       status[i] = *pstat;
834     }
835   }
836   return flag;
837 }
838
839 void smpi_mpi_probe(int source, int tag, MPI_Comm comm, MPI_Status* status){
840   int flag=0;
841   //FIXME find another way to avoid busy waiting ?
842   // the issue here is that we have to wait on a nonexistent comm
843   while(flag==0){
844     smpi_mpi_iprobe(source, tag, comm, &flag, status);
845     XBT_DEBUG("Busy Waiting on probing : %d", flag);
846   }
847 }
848
849 void smpi_mpi_iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status){
850
851   MPI_Request request = build_request(nullptr, 0, MPI_CHAR, source == MPI_ANY_SOURCE ? MPI_ANY_SOURCE :
852                  smpi_group_index(smpi_comm_group(comm), source), smpi_comm_rank(comm), tag, comm, PERSISTENT | RECV);
853
854   // to avoid deadlock, we have to sleep some time here, or the timer won't advance and we will only do iprobe simcalls
855   // (especially when used as a break condition, such as while(MPI_Iprobe(...)) ... )
856   // multiplier to the sleeptime, to increase speed of execution, each failed iprobe will increase it
857   static int nsleeps = 1;
858   if(smpi_iprobe_sleep > 0)  
859     simcall_process_sleep(nsleeps*smpi_iprobe_sleep);
860   // behave like a receive, but don't do it
861   smx_mailbox_t mailbox;
862
863   print_request("New iprobe", request);
864   // We have to test both mailboxes as we don't know if we will receive one one or another
865   if (xbt_cfg_get_int("smpi/async-small-thresh") > 0){
866       mailbox = smpi_process_mailbox_small();
867       XBT_DEBUG("Trying to probe the perm recv mailbox");
868       request->action = simcall_comm_iprobe(mailbox, 0, request->src, request->tag, &match_recv, static_cast<void*>(request));
869   }
870
871   if (request->action == nullptr){
872     mailbox = smpi_process_mailbox();
873     XBT_DEBUG("trying to probe the other mailbox");
874     request->action = simcall_comm_iprobe(mailbox, 0, request->src,request->tag, &match_recv, static_cast<void*>(request));
875   }
876
877   if (request->action != nullptr){
878     simgrid::simix::Comm *sync_comm = static_cast<simgrid::simix::Comm*>(request->action);
879     MPI_Request req                 = static_cast<MPI_Request>(sync_comm->src_data);
880     *flag = 1;
881     if(status != MPI_STATUS_IGNORE && (req->flags & PREPARED) == 0) {
882       status->MPI_SOURCE = smpi_group_rank(smpi_comm_group(comm), req->src);
883       status->MPI_TAG    = req->tag;
884       status->MPI_ERROR  = MPI_SUCCESS;
885       status->count      = req->real_size;
886     }
887     nsleeps = 1;//reset the number of sleeps we will do next time
888   }
889   else {
890     *flag = 0;
891     nsleeps++;
892   }
893   smpi_mpi_request_free(&request);
894
895   return;
896 }
897
898 void smpi_mpi_wait(MPI_Request * request, MPI_Status * status)
899 {
900   print_request("Waiting", *request);
901   if ((*request)->flags & PREPARED) {
902     smpi_empty_status(status);
903     return;
904   }
905
906   if ((*request)->action != nullptr)
907     // this is not a detached send
908     simcall_comm_wait((*request)->action, -1.0);
909
910   finish_wait(request, status);
911   if (*request != MPI_REQUEST_NULL && (((*request)->flags & NON_PERSISTENT)!=0))
912       *request = MPI_REQUEST_NULL;
913 }
914
915 int smpi_mpi_waitany(int count, MPI_Request requests[], MPI_Status * status)
916 {
917   xbt_dynar_t comms;
918   int i;
919   int size = 0;
920   int index = MPI_UNDEFINED;
921   int *map;
922
923   if(count > 0) {
924     // Wait for a request to complete
925     comms = xbt_dynar_new(sizeof(smx_synchro_t), nullptr);
926     map = xbt_new(int, count);
927     XBT_DEBUG("Wait for one of %d", count);
928     for(i = 0; i < count; i++) {
929       if (requests[i] != MPI_REQUEST_NULL && !(requests[i]->flags & PREPARED) && !(requests[i]->flags & FINISHED)) {
930         if (requests[i]->action != nullptr) {
931           XBT_DEBUG("Waiting any %p ", requests[i]);
932           xbt_dynar_push(comms, &requests[i]->action);
933           map[size] = i;
934           size++;
935         }else{
936          //This is a finished detached request, let's return this one
937          size=0;//so we free the dynar but don't do the waitany call
938          index=i;
939          finish_wait(&requests[i], status);//cleanup if refcount = 0
940          if (requests[i] != MPI_REQUEST_NULL && (requests[i]->flags & NON_PERSISTENT))
941          requests[i]=MPI_REQUEST_NULL;//set to null
942          break;
943          }
944       }
945     }
946     if(size > 0) {
947       i = simcall_comm_waitany(comms);
948
949       // not MPI_UNDEFINED, as this is a simix return code
950       if (i != -1) {
951         index = map[i];
952         finish_wait(&requests[index], status);
953         if (requests[i] != MPI_REQUEST_NULL && (requests[i]->flags & NON_PERSISTENT))
954         requests[index] = MPI_REQUEST_NULL;
955       }
956     }
957     xbt_free(map);
958     xbt_dynar_free(&comms);
959   }
960
961   if (index==MPI_UNDEFINED)
962     smpi_empty_status(status);
963
964   return index;
965 }
966
967 int smpi_mpi_waitall(int count, MPI_Request requests[], MPI_Status status[])
968 {
969   int  index, c;
970   MPI_Status stat;
971   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
972   int retvalue = MPI_SUCCESS;
973   //tag invalid requests in the set
974   if (status != MPI_STATUSES_IGNORE) {
975     for (c = 0; c < count; c++) {
976       if (requests[c] == MPI_REQUEST_NULL || requests[c]->dst == MPI_PROC_NULL || (requests[c]->flags & PREPARED)) {
977         smpi_empty_status(&status[c]);
978       } else if (requests[c]->src == MPI_PROC_NULL) {
979         smpi_empty_status(&status[c]);
980         status[c].MPI_SOURCE = MPI_PROC_NULL;
981       }
982     }
983   }
984   for(c = 0; c < count; c++) {
985
986     if (MC_is_active() || MC_record_replay_is_active()) {
987       smpi_mpi_wait(&requests[c], pstat);
988       index = c;
989     } else {
990       index = smpi_mpi_waitany(count, requests, pstat);
991       if (index == MPI_UNDEFINED)
992         break;
993       if (requests[index] != MPI_REQUEST_NULL && (requests[index]->flags & NON_PERSISTENT))
994       requests[index]=MPI_REQUEST_NULL;
995     }
996     if (status != MPI_STATUSES_IGNORE) {
997       status[index] = *pstat;
998       if (status[index].MPI_ERROR == MPI_ERR_TRUNCATE)
999         retvalue = MPI_ERR_IN_STATUS;
1000     }
1001   }
1002
1003   return retvalue;
1004 }
1005
1006 int smpi_mpi_waitsome(int incount, MPI_Request requests[], int *indices, MPI_Status status[])
1007 {
1008   int i;
1009   int count = 0;
1010   int index;
1011   MPI_Status stat;
1012   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
1013
1014   for(i = 0; i < incount; i++)
1015   {
1016     index=smpi_mpi_waitany(incount, requests, pstat);
1017     if(index!=MPI_UNDEFINED){
1018       indices[count] = index;
1019       count++;
1020       if(status != MPI_STATUSES_IGNORE) {
1021         status[index] = *pstat;
1022       }
1023      if (requests[index] != MPI_REQUEST_NULL && (requests[index]->flags & NON_PERSISTENT))
1024      requests[index]=MPI_REQUEST_NULL;
1025     }else{
1026       return MPI_UNDEFINED;
1027     }
1028   }
1029   return count;
1030 }
1031
1032 int smpi_mpi_testsome(int incount, MPI_Request requests[], int *indices, MPI_Status status[])
1033 {
1034   int i;
1035   int count = 0;
1036   int count_dead = 0;
1037   MPI_Status stat;
1038   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
1039
1040   for(i = 0; i < incount; i++) {
1041     if((requests[i] != MPI_REQUEST_NULL)) {
1042       if(smpi_mpi_test(&requests[i], pstat)) {
1043          indices[i] = 1;
1044          count++;
1045          if(status != MPI_STATUSES_IGNORE) {
1046            status[i] = *pstat;
1047          }
1048          if ((requests[i] != MPI_REQUEST_NULL) && requests[i]->flags & NON_PERSISTENT)
1049          requests[i]=MPI_REQUEST_NULL;
1050       }
1051     }else{
1052       count_dead++;
1053     }
1054   }
1055   if(count_dead==incount)
1056     return MPI_UNDEFINED;
1057   else return count;
1058 }
1059
1060 void smpi_mpi_bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1061 {
1062     smpi_coll_tuned_bcast_binomial_tree(buf, count, datatype, root, comm);
1063 }
1064
1065 void smpi_mpi_barrier(MPI_Comm comm)
1066 {
1067     smpi_coll_tuned_barrier_ompi_basic_linear(comm);
1068 }
1069
1070 void smpi_mpi_gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1071                      void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
1072 {
1073   int system_tag = COLL_TAG_GATHER;
1074   int rank, size, src, index;
1075   MPI_Aint lb = 0, recvext = 0;
1076   MPI_Request *requests;
1077
1078   rank = smpi_comm_rank(comm);
1079   size = smpi_comm_size(comm);
1080   if(rank != root) {
1081     // Send buffer to root
1082     smpi_mpi_send(sendbuf, sendcount, sendtype, root, system_tag, comm);
1083   } else {
1084     smpi_datatype_extent(recvtype, &lb, &recvext);
1085     // Local copy from root
1086     smpi_datatype_copy(sendbuf, sendcount, sendtype, static_cast<char*>(recvbuf) + root * recvcount * recvext, recvcount, recvtype);
1087     // Receive buffers from senders
1088     requests = xbt_new(MPI_Request, size - 1);
1089     index = 0;
1090     for(src = 0; src < size; src++) {
1091       if(src != root) {
1092         requests[index] = smpi_irecv_init(static_cast<char*>(recvbuf) + src * recvcount * recvext, recvcount, recvtype,
1093                                           src, system_tag, comm);
1094         index++;
1095       }
1096     }
1097     // Wait for completion of irecv's.
1098     smpi_mpi_startall(size - 1, requests);
1099     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
1100     for(src = 0; src < size-1; src++) {
1101       smpi_mpi_request_free(&requests[src]);
1102     }
1103     xbt_free(requests);
1104   }
1105 }
1106
1107 void smpi_mpi_reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts, MPI_Datatype datatype, MPI_Op op,
1108                              MPI_Comm comm)
1109 {
1110     int i, size, count;
1111     int *displs;
1112     int rank = smpi_comm_rank(comm);
1113     void *tmpbuf;
1114
1115     /* arbitrarily choose root as rank 0 */
1116     size = smpi_comm_size(comm);
1117     count = 0;
1118     displs = xbt_new(int, size);
1119     for (i = 0; i < size; i++) {
1120       displs[i] = count;
1121       count += recvcounts[i];
1122     }
1123     tmpbuf=static_cast<void*>(smpi_get_tmp_sendbuffer(count*smpi_datatype_get_extent(datatype)));
1124
1125     mpi_coll_reduce_fun(sendbuf, tmpbuf, count, datatype, op, 0, comm);
1126     smpi_mpi_scatterv(tmpbuf, recvcounts, displs, datatype, recvbuf, recvcounts[rank], datatype, 0, comm);
1127     xbt_free(displs);
1128     smpi_free_tmp_buffer(tmpbuf);
1129 }
1130
1131 void smpi_mpi_gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int *recvcounts, int *displs,
1132                       MPI_Datatype recvtype, int root, MPI_Comm comm)
1133 {
1134   int system_tag = COLL_TAG_GATHERV;
1135   int rank, size, src, index;
1136   MPI_Aint lb = 0, recvext = 0;
1137   MPI_Request *requests;
1138
1139   rank = smpi_comm_rank(comm);
1140   size = smpi_comm_size(comm);
1141   if(rank != root) {
1142     // Send buffer to root
1143     smpi_mpi_send(sendbuf, sendcount, sendtype, root, system_tag, comm);
1144   } else {
1145     smpi_datatype_extent(recvtype, &lb, &recvext);
1146     // Local copy from root
1147     smpi_datatype_copy(sendbuf, sendcount, sendtype, static_cast<char*>(recvbuf) + displs[root] * recvext,
1148                        recvcounts[root], recvtype);
1149     // Receive buffers from senders
1150     requests = xbt_new(MPI_Request, size - 1);
1151     index = 0;
1152     for(src = 0; src < size; src++) {
1153       if(src != root) {
1154         requests[index] = smpi_irecv_init(static_cast<char*>(recvbuf) + displs[src] * recvext,
1155                           recvcounts[src], recvtype, src, system_tag, comm);
1156         index++;
1157       }
1158     }
1159     // Wait for completion of irecv's.
1160     smpi_mpi_startall(size - 1, requests);
1161     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
1162     for(src = 0; src < size-1; src++) {
1163       smpi_mpi_request_free(&requests[src]);
1164     }
1165     xbt_free(requests);
1166   }
1167 }
1168
1169 void smpi_mpi_allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1170                         void *recvbuf,int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
1171 {
1172   int system_tag = COLL_TAG_ALLGATHER;
1173   int rank, size, other, index;
1174   MPI_Aint lb = 0, recvext = 0;
1175   MPI_Request *requests;
1176
1177   rank = smpi_comm_rank(comm);
1178   size = smpi_comm_size(comm);
1179   // FIXME: check for errors
1180   smpi_datatype_extent(recvtype, &lb, &recvext);
1181   // Local copy from self
1182   smpi_datatype_copy(sendbuf, sendcount, sendtype, static_cast<char *>(recvbuf) + rank * recvcount * recvext, recvcount, recvtype);
1183   // Send/Recv buffers to/from others;
1184   requests = xbt_new(MPI_Request, 2 * (size - 1));
1185   index = 0;
1186   for(other = 0; other < size; other++) {
1187     if(other != rank) {
1188       requests[index] = smpi_isend_init(sendbuf, sendcount, sendtype, other, system_tag,comm);
1189       index++;
1190       requests[index] = smpi_irecv_init(static_cast<char *>(recvbuf) + other * recvcount * recvext, recvcount, recvtype, other,
1191                                         system_tag, comm);
1192       index++;
1193     }
1194   }
1195   // Wait for completion of all comms.
1196   smpi_mpi_startall(2 * (size - 1), requests);
1197   smpi_mpi_waitall(2 * (size - 1), requests, MPI_STATUS_IGNORE);
1198   for(other = 0; other < 2*(size-1); other++) {
1199     smpi_mpi_request_free(&requests[other]);
1200   }
1201   xbt_free(requests);
1202 }
1203
1204 void smpi_mpi_allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf,
1205                          int *recvcounts, int *displs, MPI_Datatype recvtype, MPI_Comm comm)
1206 {
1207   int system_tag = COLL_TAG_ALLGATHERV;
1208   int rank, size, other, index;
1209   MPI_Aint lb = 0, recvext = 0;
1210   MPI_Request *requests;
1211
1212   rank = smpi_comm_rank(comm);
1213   size = smpi_comm_size(comm);
1214   smpi_datatype_extent(recvtype, &lb, &recvext);
1215   // Local copy from self
1216   smpi_datatype_copy(sendbuf, sendcount, sendtype, static_cast<char *>(recvbuf) + displs[rank] * recvext,recvcounts[rank], recvtype);
1217   // Send buffers to others;
1218   requests = xbt_new(MPI_Request, 2 * (size - 1));
1219   index = 0;
1220   for(other = 0; other < size; other++) {
1221     if(other != rank) {
1222       requests[index] =
1223         smpi_isend_init(sendbuf, sendcount, sendtype, other, system_tag, comm);
1224       index++;
1225       requests[index] = smpi_irecv_init(static_cast<char *>(recvbuf) + displs[other] * recvext, recvcounts[other],
1226                           recvtype, other, system_tag, comm);
1227       index++;
1228     }
1229   }
1230   // Wait for completion of all comms.
1231   smpi_mpi_startall(2 * (size - 1), requests);
1232   smpi_mpi_waitall(2 * (size - 1), requests, MPI_STATUS_IGNORE);
1233   for(other = 0; other < 2*(size-1); other++) {
1234     smpi_mpi_request_free(&requests[other]);
1235   }
1236   xbt_free(requests);
1237 }
1238
1239 void smpi_mpi_scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1240                       void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
1241 {
1242   int system_tag = COLL_TAG_SCATTER;
1243   int rank, size, dst, index;
1244   MPI_Aint lb = 0, sendext = 0;
1245   MPI_Request *requests;
1246
1247   rank = smpi_comm_rank(comm);
1248   size = smpi_comm_size(comm);
1249   if(rank != root) {
1250     // Recv buffer from root
1251     smpi_mpi_recv(recvbuf, recvcount, recvtype, root, system_tag, comm, MPI_STATUS_IGNORE);
1252   } else {
1253     smpi_datatype_extent(sendtype, &lb, &sendext);
1254     // Local copy from root
1255     if(recvbuf!=MPI_IN_PLACE){
1256         smpi_datatype_copy(static_cast<char *>(sendbuf) + root * sendcount * sendext,
1257                            sendcount, sendtype, recvbuf, recvcount, recvtype);
1258     }
1259     // Send buffers to receivers
1260     requests = xbt_new(MPI_Request, size - 1);
1261     index = 0;
1262     for(dst = 0; dst < size; dst++) {
1263       if(dst != root) {
1264         requests[index] = smpi_isend_init(static_cast<char *>(sendbuf) + dst * sendcount * sendext, sendcount, sendtype, dst,
1265                                           system_tag, comm);
1266         index++;
1267       }
1268     }
1269     // Wait for completion of isend's.
1270     smpi_mpi_startall(size - 1, requests);
1271     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
1272     for(dst = 0; dst < size-1; dst++) {
1273       smpi_mpi_request_free(&requests[dst]);
1274     }
1275     xbt_free(requests);
1276   }
1277 }
1278
1279 void smpi_mpi_scatterv(void *sendbuf, int *sendcounts, int *displs, MPI_Datatype sendtype, void *recvbuf, int recvcount,
1280                        MPI_Datatype recvtype, int root, MPI_Comm comm)
1281 {
1282   int system_tag = COLL_TAG_SCATTERV;
1283   int rank, size, dst, index;
1284   MPI_Aint lb = 0, sendext = 0;
1285   MPI_Request *requests;
1286
1287   rank = smpi_comm_rank(comm);
1288   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     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 rank, size, src, index;
1324   MPI_Aint lb = 0, dataext = 0;
1325   MPI_Request *requests;
1326   void **tmpbufs;
1327
1328   char* sendtmpbuf = static_cast<char *>(sendbuf);
1329
1330
1331   rank = smpi_comm_rank(comm);
1332   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 rank, size, other, index;
1402   MPI_Aint lb = 0, dataext = 0;
1403   MPI_Request *requests;
1404   void **tmpbufs;
1405
1406   rank = smpi_comm_rank(comm);
1407   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 rank, size, other, index;
1464   MPI_Aint lb = 0, dataext = 0;
1465   MPI_Request *requests;
1466   void **tmpbufs;
1467   int recvbuf_is_empty=1;
1468   rank = smpi_comm_rank(comm);
1469   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] =
1480       smpi_irecv_init(tmpbufs[index], count, datatype, other, system_tag, comm);
1481     index++;
1482   }
1483   for(other = rank + 1; other < size; other++) {
1484     requests[index] =
1485       smpi_isend_init(sendbuf, count, datatype, other, system_tag, comm);
1486     index++;
1487   }
1488   // Wait for completion of all comms.
1489   smpi_mpi_startall(size - 1, requests);
1490   if(smpi_op_is_commute(op)){
1491     for(other = 0; other < size - 1; other++) {
1492       index = smpi_mpi_waitany(size - 1, requests, MPI_STATUS_IGNORE);
1493       if(index == MPI_UNDEFINED) {
1494         break;
1495       }
1496       if(index < rank) {
1497         if(recvbuf_is_empty){
1498           smpi_datatype_copy(tmpbufs[index], count, datatype, recvbuf, count, datatype);
1499           recvbuf_is_empty=0;
1500         }else
1501         // #Request is below rank: it's a irecv
1502         smpi_op_apply(op, tmpbufs[index], recvbuf, &count, &datatype);
1503       }
1504     }
1505   }else{
1506     //non commutative case, wait in order
1507     for(other = 0; other < size - 1; other++) {
1508       smpi_mpi_wait(&(requests[other]), MPI_STATUS_IGNORE);
1509       if(index < rank) {
1510           if(recvbuf_is_empty){
1511             smpi_datatype_copy(tmpbufs[other], count, datatype, recvbuf, count, datatype);
1512             recvbuf_is_empty=0;
1513           }else smpi_op_apply(op, tmpbufs[other], recvbuf, &count, &datatype);
1514       }
1515     }
1516   }
1517   for(index = 0; index < rank; index++) {
1518     smpi_free_tmp_buffer(tmpbufs[index]);
1519   }
1520   for(index = 0; index < size-1; index++) {
1521     smpi_mpi_request_free(&requests[index]);
1522   }
1523   xbt_free(tmpbufs);
1524   xbt_free(requests);
1525 }