Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
d7b67c078ea37d8a50949aaa042362c367acc147
[simgrid.git] / src / smpi / smpi_base.c
1 /* Copyright (c) 2007, 2008, 2009, 2010. 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 "private.h"
8 #include "xbt/virtu.h"
9 #include "mc/mc.h"
10 #include "xbt/replay.h"
11 #include <errno.h>
12 #include "simix/smx_private.h"
13 #include "surf/surf.h"
14 #include "simgrid/sg_config.h"
15 #include "colls/colls.h"
16
17 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_base, smpi, "Logging specific to SMPI (base)");
18
19
20 static int match_recv(void* a, void* b, smx_action_t ignored) {
21    MPI_Request ref = (MPI_Request)a;
22    MPI_Request req = (MPI_Request)b;
23    XBT_DEBUG("Trying to match a recv of src %d against %d, tag %d against %d",ref->src,req->src, ref->tag, req->tag);
24
25   xbt_assert(ref, "Cannot match recv against null reference");
26   xbt_assert(req, "Cannot match recv against null request");
27   if((ref->src == MPI_ANY_SOURCE || req->src == ref->src)
28     && (ref->tag == MPI_ANY_TAG || req->tag == ref->tag)){
29     //we match, we can transfer some values
30     // FIXME : move this to the copy function ?
31     if(ref->src == MPI_ANY_SOURCE)ref->real_src = req->src;
32     if(ref->tag == MPI_ANY_TAG)ref->real_tag = req->tag;
33     if(ref->real_size < req->real_size) ref->truncated = 1;
34     if(req->detached==1){
35         ref->detached_sender=req; //tie the sender to the receiver, as it is detached and has to be freed in the receiver
36     }
37     return 1;
38   }else return 0;
39 }
40
41 static int match_send(void* a, void* b,smx_action_t ignored) {
42    MPI_Request ref = (MPI_Request)a;
43    MPI_Request req = (MPI_Request)b;
44    XBT_DEBUG("Trying to match a send of src %d against %d, tag %d against %d",ref->src,req->src, ref->tag, req->tag);
45    xbt_assert(ref, "Cannot match send against null reference");
46    xbt_assert(req, "Cannot match send against null request");
47
48    if((req->src == MPI_ANY_SOURCE || req->src == ref->src)
49              && (req->tag == MPI_ANY_TAG || req->tag == ref->tag))
50    {
51      if(req->src == MPI_ANY_SOURCE)req->real_src = ref->src;
52      if(req->tag == MPI_ANY_TAG)req->real_tag = ref->tag;
53      if(req->real_size < ref->real_size) req->truncated = 1;
54      if(ref->detached==1){
55          req->detached_sender=ref; //tie the sender to the receiver, as it is detached and has to be freed in the receiver
56      }
57
58      return 1;
59    } else return 0;
60 }
61
62
63 typedef struct s_smpi_factor *smpi_factor_t;
64 typedef struct s_smpi_factor {
65   long factor;
66   int nb_values;
67   double values[4];//arbitrary set to 4
68 } s_smpi_factor_t;
69 xbt_dynar_t smpi_os_values = NULL;
70 xbt_dynar_t smpi_or_values = NULL;
71 xbt_dynar_t smpi_ois_values = NULL;
72
73 // Methods used to parse and store the values for timing injections in smpi
74 // These are taken from surf/network.c and generalized to have more factors
75 // These methods should be merged with those in surf/network.c (moved somewhere in xbt ?)
76
77 static int factor_cmp(const void *pa, const void *pb)
78 {
79   return (((s_smpi_factor_t*)pa)->factor > ((s_smpi_factor_t*)pb)->factor);
80 }
81
82
83 static xbt_dynar_t parse_factor(const char *smpi_coef_string)
84 {
85   char *value = NULL;
86   unsigned int iter = 0;
87   s_smpi_factor_t fact;
88   int i=0;
89   xbt_dynar_t smpi_factor, radical_elements, radical_elements2 = NULL;
90
91   smpi_factor = xbt_dynar_new(sizeof(s_smpi_factor_t), NULL);
92   radical_elements = xbt_str_split(smpi_coef_string, ";");
93   xbt_dynar_foreach(radical_elements, iter, value) {
94     fact.nb_values=0;
95     radical_elements2 = xbt_str_split(value, ":");
96     if (xbt_dynar_length(radical_elements2) <2 || xbt_dynar_length(radical_elements2) > 5)
97       xbt_die("Malformed radical for smpi factor!");
98     for(i =0; i<xbt_dynar_length(radical_elements2);i++ ){
99         if (i==0){
100            fact.factor = atol(xbt_dynar_get_as(radical_elements2, i, char *));
101         }else{
102            fact.values[fact.nb_values] = atof(xbt_dynar_get_as(radical_elements2, i, char *));
103            fact.nb_values++;
104         }
105     }
106
107     xbt_dynar_push_as(smpi_factor, s_smpi_factor_t, fact);
108     XBT_DEBUG("smpi_factor:\t%ld : %d values, first: %f", fact.factor, fact.nb_values ,fact.values[0]);
109     xbt_dynar_free(&radical_elements2);
110   }
111   xbt_dynar_free(&radical_elements);
112   iter=0;
113   xbt_dynar_sort(smpi_factor, &factor_cmp);
114   xbt_dynar_foreach(smpi_factor, iter, fact) {
115     XBT_DEBUG("smpi_factor:\t%ld : %d values, first: %f", fact.factor, fact.nb_values ,fact.values[0]);
116   }
117   return smpi_factor;
118 }
119
120 static double smpi_os(double size)
121 {
122   if (!smpi_os_values) {
123     smpi_os_values = parse_factor(sg_cfg_get_string("smpi/os"));
124     smpi_register_static(smpi_os_values, xbt_dynar_free_voidp);
125   }
126   unsigned int iter = 0;
127   s_smpi_factor_t fact;
128   double current=0.0;
129   xbt_dynar_foreach(smpi_os_values, iter, fact) {
130     if (size <= fact.factor) {
131         XBT_DEBUG("os : %lf <= %ld return %f", size, fact.factor, current);
132       return current;
133     }else{
134       current=fact.values[0]+fact.values[1]*size;
135     }
136   }
137   XBT_DEBUG("os : %lf > %ld return %f", size, fact.factor, current);
138
139   return current;
140 }
141
142 static double smpi_ois(double size)
143 {
144   if (!smpi_ois_values) {
145     smpi_ois_values = parse_factor(sg_cfg_get_string("smpi/ois"));
146     smpi_register_static(smpi_ois_values, xbt_dynar_free_voidp);
147   }
148   unsigned int iter = 0;
149   s_smpi_factor_t fact;
150   double current=0.0;
151   xbt_dynar_foreach(smpi_ois_values, iter, fact) {
152     if (size <= fact.factor) {
153         XBT_DEBUG("ois : %lf <= %ld return %f", size, fact.factor, current);
154       return current;
155     }else{
156       current=fact.values[0]+fact.values[1]*size;
157     }
158   }
159   XBT_DEBUG("ois : %lf > %ld return %f", size, fact.factor, current);
160
161   return current;
162 }
163
164 static double smpi_or(double size)
165 {
166   if (!smpi_or_values) {
167     smpi_or_values = parse_factor(sg_cfg_get_string("smpi/or"));
168     smpi_register_static(smpi_or_values, xbt_dynar_free_voidp);
169   }
170   unsigned int iter = 0;
171   s_smpi_factor_t fact;
172   double current=0.0;
173   xbt_dynar_foreach(smpi_or_values, iter, fact) {
174     if (size <= fact.factor) {
175         XBT_DEBUG("or : %lf <= %ld return %f", size, fact.factor, current);
176       return current;
177     }else
178       current=fact.values[0]+fact.values[1]*size;
179   }
180   XBT_DEBUG("or : %lf > %ld return %f", size, fact.factor, current);
181
182   return current;
183 }
184
185 static MPI_Request build_request(void *buf, int count,
186                                  MPI_Datatype datatype, int src, int dst,
187                                  int tag, MPI_Comm comm, unsigned flags)
188 {
189   MPI_Request request;
190
191   void *old_buf = NULL;
192
193   request = xbt_new(s_smpi_mpi_request_t, 1);
194
195   s_smpi_subtype_t *subtype = datatype->substruct;
196
197   if(datatype->has_subtype == 1){
198     // This part handles the problem of non-contiguous memory
199     old_buf = buf;
200     buf = count==0 ? NULL : xbt_malloc(count*smpi_datatype_size(datatype));
201     if (flags & SEND) {
202       subtype->serialize(old_buf, buf, count, datatype->substruct);
203     }
204   }
205
206   request->buf = buf;
207   // This part handles the problem of non-contiguous memory (for the
208   // unserialisation at the reception)
209   request->old_buf = old_buf;
210   request->old_type = datatype;
211
212   request->size = smpi_datatype_size(datatype) * count;
213   request->src = src;
214   request->dst = dst;
215   request->tag = tag;
216   request->comm = comm;
217   request->action = NULL;
218   request->flags = flags;
219   request->detached = 0;
220   request->detached_sender = NULL;
221
222   request->truncated = 0;
223   request->real_size = 0;
224   request->real_tag = 0;
225
226   request->refcount=1;
227 #ifdef HAVE_TRACING
228   request->send = 0;
229   request->recv = 0;
230 #endif
231   if (flags & SEND) smpi_datatype_unuse(datatype);
232
233   return request;
234 }
235
236
237 void smpi_empty_status(MPI_Status * status)
238 {
239   if(status != MPI_STATUS_IGNORE) {
240     status->MPI_SOURCE = MPI_ANY_SOURCE;
241     status->MPI_TAG = MPI_ANY_TAG;
242     status->MPI_ERROR = MPI_SUCCESS;
243     status->count=0;
244   }
245 }
246
247 void smpi_action_trace_run(char *path)
248 {
249   char *name;
250   xbt_dynar_t todo;
251   xbt_dict_cursor_t cursor;
252
253   action_fp=NULL;
254   if (path) {
255     action_fp = fopen(path, "r");
256     xbt_assert(action_fp != NULL, "Cannot open %s: %s", path,
257                strerror(errno));
258   }
259
260   if (!xbt_dict_is_empty(action_queues)) {
261     XBT_WARN
262       ("Not all actions got consumed. If the simulation ended successfully (without deadlock), you may want to add new processes to your deployment file.");
263
264
265     xbt_dict_foreach(action_queues, cursor, name, todo) {
266       XBT_WARN("Still %lu actions for %s", xbt_dynar_length(todo), name);
267     }
268   }
269
270   if (path)
271     fclose(action_fp);
272   xbt_dict_free(&action_queues);
273   action_queues = xbt_dict_new_homogeneous(NULL);
274 }
275
276 static void smpi_mpi_request_free_voidp(void* request)
277 {
278   MPI_Request req = request;
279   smpi_mpi_request_free(&req);
280 }
281
282 /* MPI Low level calls */
283 MPI_Request smpi_mpi_send_init(void *buf, int count, MPI_Datatype datatype,
284                                int dst, int tag, MPI_Comm comm)
285 {
286   MPI_Request request =
287     build_request(buf, count, datatype, smpi_process_index(), smpi_group_index(smpi_comm_group(comm), dst), tag,
288                   comm, PERSISTENT | SEND);
289   request->refcount++;
290   return request;
291 }
292
293 MPI_Request smpi_mpi_ssend_init(void *buf, int count, MPI_Datatype datatype,
294                                int dst, int tag, MPI_Comm comm)
295 {
296   MPI_Request request =
297     build_request(buf, count, datatype, smpi_process_index(), smpi_group_index(smpi_comm_group(comm), dst), tag,
298                   comm, PERSISTENT | SSEND | SEND);
299   request->refcount++;
300   return request;
301 }
302
303 MPI_Request smpi_mpi_recv_init(void *buf, int count, MPI_Datatype datatype,
304                                int src, int tag, MPI_Comm comm)
305 {
306   MPI_Request request =
307     build_request(buf, count, datatype, src == MPI_ANY_SOURCE ? MPI_ANY_SOURCE : smpi_group_index(smpi_comm_group(comm), src), smpi_process_index(), tag,
308                   comm, PERSISTENT | RECV);
309   request->refcount++;
310   return request;
311 }
312
313 void smpi_mpi_start(MPI_Request request)
314 {
315   smx_rdv_t mailbox;
316
317   xbt_assert(!request->action,
318              "Cannot (re)start a non-finished communication");
319   if(request->flags & RECV) {
320     print_request("New recv", request);
321     if (request->size < sg_cfg_get_int("smpi/async_small_thres"))
322       mailbox = smpi_process_mailbox_small();
323     else
324       mailbox = smpi_process_mailbox();
325     // we make a copy here, as the size is modified by simix, and we may reuse the request in another receive later
326     request->real_size=request->size;
327     smpi_datatype_use(request->old_type);
328     request->action = simcall_comm_irecv(mailbox, request->buf, &request->real_size, &match_recv, request);
329
330     //integrate pseudo-timing for buffering of small messages, do not bother to execute the simcall if 0
331     double sleeptime = request->detached ? smpi_or(request->size) : 0.0;
332     if(sleeptime!=0.0){
333         simcall_process_sleep(sleeptime);
334         XBT_DEBUG("receiving size of %zu : sleep %lf ", request->size, smpi_or(request->size));
335     }
336
337   } else {
338
339
340     int receiver = request->dst;//smpi_group_index(smpi_comm_group(request->comm), request->dst);
341
342     #ifdef HAVE_TRACING
343       int rank = smpi_process_index();
344       if (TRACE_smpi_view_internals()) {
345         TRACE_smpi_send(rank, rank, receiver);
346       }
347     #endif
348 /*    if(receiver == MPI_UNDEFINED) {*/
349 /*      XBT_WARN("Trying to send a message to a wrong rank");*/
350 /*      return;*/
351 /*    }*/
352     print_request("New send", request);
353     if (request->size < sg_cfg_get_int("smpi/async_small_thres")) { // eager mode
354       mailbox = smpi_process_remote_mailbox_small(receiver);
355     }else{
356       XBT_DEBUG("Send request %p is not in the permanent receive mailbox (buf: %p)",request,request->buf);
357       mailbox = smpi_process_remote_mailbox(receiver);
358     }
359     if ( (! (request->flags & SSEND)) && (request->size < sg_cfg_get_int("smpi/send_is_detached_thres"))) {
360       void *oldbuf = NULL;
361       request->detached = 1;
362       request->refcount++;
363       if(request->old_type->has_subtype == 0){
364         oldbuf = request->buf;
365         if (oldbuf && request->size!=0){
366           request->buf = xbt_malloc(request->size);
367           memcpy(request->buf,oldbuf,request->size);
368         }
369       }
370       XBT_DEBUG("Send request %p is detached; buf %p copied into %p",request,oldbuf,request->buf);
371     }
372
373     // we make a copy here, as the size is modified by simix, and we may reuse the request in another receive later
374     request->real_size=request->size;
375     smpi_datatype_use(request->old_type);
376
377     //if we are giving back the control to the user without waiting for completion, we have to inject timings
378     double sleeptime =0.0;
379     if(request->detached || (request->flags & (ISEND|SSEND))){// issend should be treated as isend
380       //isend and send timings may be different
381       sleeptime = (request->flags & ISEND)? smpi_ois(request->size) : smpi_os(request->size);
382     }
383
384     if(sleeptime!=0.0){
385         simcall_process_sleep(sleeptime);
386         XBT_DEBUG("sending size of %zu : sleep %lf ", request->size, smpi_os(request->size));
387     }
388
389     request->action =
390       simcall_comm_isend(mailbox, request->size, -1.0,
391                          request->buf, request->real_size,
392                          &match_send,
393                          &smpi_mpi_request_free_voidp, // how to free the userdata if a detached send fails
394                          request,
395                          // detach if msg size < eager/rdv switch limit
396                          request->detached);
397
398 #ifdef HAVE_TRACING
399     /* FIXME: detached sends are not traceable (request->action == NULL) */
400     if (request->action)
401       simcall_set_category(request->action, TRACE_internal_smpi_get_category());
402
403 #endif
404
405   }
406
407 }
408
409 void smpi_mpi_startall(int count, MPI_Request * requests)
410 {
411   int i;
412
413   for(i = 0; i < count; i++) {
414     smpi_mpi_start(requests[i]);
415   }
416 }
417
418 void smpi_mpi_request_free(MPI_Request * request)
419 {
420
421   if((*request) != MPI_REQUEST_NULL){
422     (*request)->refcount--;
423     if((*request)->refcount<0) xbt_die("wrong refcount");
424
425     if((*request)->refcount==0){
426         print_request("Destroying", (*request));
427         xbt_free(*request);
428         *request = MPI_REQUEST_NULL;
429     }else{
430         print_request("Decrementing", (*request));
431
432     }
433   }else{
434       xbt_die("freeing an already free request");
435   }
436 }
437
438 MPI_Request smpi_isend_init(void *buf, int count, MPI_Datatype datatype,
439                             int dst, int tag, MPI_Comm comm)
440 {
441   MPI_Request request =
442     build_request(buf, count, datatype, smpi_process_index(), smpi_group_index(smpi_comm_group(comm), dst), tag,
443                   comm, NON_PERSISTENT | SEND);
444
445   return request;
446 }
447
448 MPI_Request smpi_mpi_isend(void *buf, int count, MPI_Datatype datatype,
449                            int dst, int tag, MPI_Comm comm)
450 {
451   MPI_Request request =
452       build_request(buf, count, datatype, smpi_process_index(), smpi_group_index(smpi_comm_group(comm), dst), tag,
453                     comm, NON_PERSISTENT | ISEND | SEND);
454
455   smpi_mpi_start(request);
456   return request;
457 }
458
459 MPI_Request smpi_mpi_issend(void *buf, int count, MPI_Datatype datatype,
460                            int dst, int tag, MPI_Comm comm)
461 {
462   MPI_Request request =
463       build_request(buf, count, datatype, smpi_process_index(), smpi_group_index(smpi_comm_group(comm), dst), tag,
464                     comm, NON_PERSISTENT | ISEND | SSEND | SEND);
465   smpi_mpi_start(request);
466   return request;
467 }
468
469
470
471 MPI_Request smpi_irecv_init(void *buf, int count, MPI_Datatype datatype,
472                             int src, int tag, MPI_Comm comm)
473 {
474   MPI_Request request =
475     build_request(buf, count, datatype, src == MPI_ANY_SOURCE ? MPI_ANY_SOURCE : smpi_group_index(smpi_comm_group(comm), src), smpi_process_index(), tag,
476                   comm, NON_PERSISTENT | RECV);
477   return request;
478 }
479
480 MPI_Request smpi_mpi_irecv(void *buf, int count, MPI_Datatype datatype,
481                            int src, int tag, MPI_Comm comm)
482 {
483   MPI_Request request =
484       build_request(buf, count, datatype, src == MPI_ANY_SOURCE ? MPI_ANY_SOURCE : smpi_group_index(smpi_comm_group(comm), src), smpi_process_index(), tag,
485                     comm, NON_PERSISTENT | RECV);
486
487   smpi_mpi_start(request);
488   return request;
489 }
490
491 void smpi_mpi_recv(void *buf, int count, MPI_Datatype datatype, int src,
492                    int tag, MPI_Comm comm, MPI_Status * status)
493 {
494   MPI_Request request;
495   request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
496   smpi_mpi_wait(&request, status);
497 }
498
499
500
501 void smpi_mpi_send(void *buf, int count, MPI_Datatype datatype, int dst,
502                    int tag, MPI_Comm comm)
503 {
504   MPI_Request request =
505       build_request(buf, count, datatype, smpi_process_index(), smpi_group_index(smpi_comm_group(comm), dst), tag,
506                     comm, NON_PERSISTENT | SEND);
507   smpi_mpi_start(request);
508   smpi_mpi_wait(&request, MPI_STATUS_IGNORE);
509
510 }
511
512 void smpi_mpi_ssend(void *buf, int count, MPI_Datatype datatype,
513                            int dst, int tag, MPI_Comm comm)
514 {
515   MPI_Request request =
516       build_request(buf, count, datatype, smpi_process_index(), smpi_group_index(smpi_comm_group(comm), dst), tag,
517                     comm, NON_PERSISTENT | SSEND | SEND);
518
519   smpi_mpi_start(request);
520   smpi_mpi_wait(&request, MPI_STATUS_IGNORE);
521 }
522
523 void smpi_mpi_sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
524                        int dst, int sendtag, void *recvbuf, int recvcount,
525                        MPI_Datatype recvtype, int src, int recvtag,
526                        MPI_Comm comm, MPI_Status * status)
527 {
528   MPI_Request requests[2];
529   MPI_Status stats[2];
530   int myid=smpi_process_index();
531   if ((smpi_group_index(smpi_comm_group(comm), dst) == myid) && (smpi_group_index(smpi_comm_group(comm), src) == myid)) {
532       smpi_datatype_copy(sendbuf, sendcount, sendtype,
533                                      recvbuf, recvcount, recvtype);
534       return;
535   }
536   requests[0] =
537     smpi_isend_init(sendbuf, sendcount, sendtype, dst, sendtag, comm);
538   requests[1] =
539     smpi_irecv_init(recvbuf, recvcount, recvtype, src, recvtag, comm);
540   smpi_mpi_startall(2, requests);
541   smpi_mpi_waitall(2, requests, stats);
542   if(status != MPI_STATUS_IGNORE) {
543     // Copy receive status
544     *status = stats[1];
545   }
546 }
547
548 int smpi_mpi_get_count(MPI_Status * status, MPI_Datatype datatype)
549 {
550   return status->count / smpi_datatype_size(datatype);
551 }
552
553 static void finish_wait(MPI_Request * request, MPI_Status * status)
554 {
555   MPI_Request req = *request;
556   if(status != MPI_STATUS_IGNORE)
557     smpi_empty_status(status);
558
559   if(!(req->detached && req->flags & SEND)){
560     if(status != MPI_STATUS_IGNORE) {
561       int src = req->src == MPI_ANY_SOURCE ? req->real_src : req->src;
562       status->MPI_SOURCE = smpi_group_rank(smpi_comm_group(req->comm), src);
563       status->MPI_TAG = req->tag == MPI_ANY_TAG ? req->real_tag : req->tag;
564       status->MPI_ERROR = req->truncated ? MPI_ERR_TRUNCATE : MPI_SUCCESS;
565       // this handles the case were size in receive differs from size in send
566       // FIXME: really this should just contain the count of receive-type blocks,
567       // right?
568       status->count = req->real_size;
569     }
570
571     print_request("Finishing", req);
572     MPI_Datatype datatype = req->old_type;
573
574     if(datatype->has_subtype == 1){
575         // This part handles the problem of non-contignous memory
576         // the unserialization at the reception
577       s_smpi_subtype_t *subtype = datatype->substruct;
578       if(req->flags & RECV) {
579         subtype->unserialize(req->buf, req->old_buf, req->real_size/smpi_datatype_size(datatype) , datatype->substruct);
580       }
581       if(req->detached == 0) free(req->buf);
582     }
583     smpi_datatype_unuse(datatype);
584
585   }
586
587 #ifdef HAVE_TRACING
588     if (TRACE_smpi_view_internals()) {
589       if(req->flags & RECV){
590         int rank = smpi_process_index();
591         int  src_traced = smpi_group_index(smpi_comm_group(req->comm), req->src == MPI_ANY_SOURCE ? req->real_src : req->src);
592         TRACE_smpi_recv(rank, src_traced, rank);
593       }
594     }
595 #endif
596
597   if(req->detached_sender!=NULL){
598     smpi_mpi_request_free(&(req->detached_sender));
599   }
600
601   if(req->flags & NON_PERSISTENT) {
602     smpi_mpi_request_free(request);
603   } else {
604     req->action = NULL;
605   }
606 }
607
608 int smpi_mpi_test(MPI_Request * request, MPI_Status * status) {
609   int flag;
610
611   //assume that request is not MPI_REQUEST_NULL (filtered in PMPI_Test or smpi_mpi_testall before)
612   if ((*request)->action == NULL)
613     flag = 1;
614   else
615     flag = simcall_comm_test((*request)->action);
616   if(flag) {
617     finish_wait(request, status);
618     request=MPI_REQUEST_NULL;
619   }else{
620     smpi_empty_status(status);
621   }
622   return flag;
623 }
624
625 int smpi_mpi_testany(int count, MPI_Request requests[], int *index,
626                      MPI_Status * status)
627 {
628   xbt_dynar_t comms;
629   int i, flag, size;
630   int* map;
631
632   *index = MPI_UNDEFINED;
633   flag = 0;
634   if(count > 0) {
635     comms = xbt_dynar_new(sizeof(smx_action_t), NULL);
636     map = xbt_new(int, count);
637     size = 0;
638     for(i = 0; i < count; i++) {
639       if((requests[i]!=MPI_REQUEST_NULL) && requests[i]->action) {
640          xbt_dynar_push(comms, &requests[i]->action);
641          map[size] = i;
642          size++;
643       }
644     }
645     if(size > 0) {
646       i = simcall_comm_testany(comms);
647       // not MPI_UNDEFINED, as this is a simix return code
648       if(i != -1) {
649         *index = map[i];
650         finish_wait(&requests[*index], status);
651         flag = 1;
652       }
653     }else{
654         //all requests are null or inactive, return true
655         flag=1;
656         smpi_empty_status(status);
657     }
658     xbt_free(map);
659     xbt_dynar_free(&comms);
660   }
661
662   return flag;
663 }
664
665
666 int smpi_mpi_testall(int count, MPI_Request requests[],
667                      MPI_Status status[])
668 {
669   MPI_Status stat;
670   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
671   int flag=1;
672   int i;
673   for(i=0; i<count; i++){
674     if(requests[i]!= MPI_REQUEST_NULL){
675       if (smpi_mpi_test(&requests[i], pstat)!=1){
676         flag=0;
677       }
678     }else{
679       smpi_empty_status(pstat);
680     }
681     if(status != MPI_STATUSES_IGNORE) {
682       status[i] = *pstat;
683     }
684   }
685   return flag;
686 }
687
688 void smpi_mpi_probe(int source, int tag, MPI_Comm comm, MPI_Status* status){
689   int flag=0;
690   //FIXME find another wait to avoid busy waiting ?
691   // the issue here is that we have to wait on a nonexistent comm
692   while(flag==0){
693     smpi_mpi_iprobe(source, tag, comm, &flag, status);
694     XBT_DEBUG("Busy Waiting on probing : %d", flag);
695   }
696 }
697
698 void smpi_mpi_iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status){
699
700   MPI_Request request =build_request(NULL, 0, MPI_CHAR, source == MPI_ANY_SOURCE ? MPI_ANY_SOURCE : smpi_group_index(smpi_comm_group(comm), source), smpi_comm_rank(comm), tag,
701             comm, NON_PERSISTENT | RECV);
702
703   //to avoid deadlock, we have to sleep some time here, or the timer won't advance and we will only do iprobe simcalls
704   double sleeptime= sg_cfg_get_double("smpi/iprobe");
705   //multiplier to the sleeptime, to increase speed of execution, each failed iprobe will increase it
706   static int nsleeps = 1;
707
708   simcall_process_sleep(sleeptime);
709
710   // behave like a receive, but don't do it
711   smx_rdv_t mailbox;
712
713   print_request("New iprobe", request);
714   // We have to test both mailboxes as we don't know if we will receive one one or another
715     if (sg_cfg_get_int("smpi/async_small_thres")>0){
716         mailbox = smpi_process_mailbox_small();
717         XBT_DEBUG("trying to probe the perm recv mailbox");
718         request->action = simcall_comm_iprobe(mailbox, request->src, request->tag, &match_recv, (void*)request);
719     }
720     if (request->action==NULL){
721         mailbox = smpi_process_mailbox();
722         XBT_DEBUG("trying to probe the other mailbox");
723         request->action = simcall_comm_iprobe(mailbox, request->src, request->tag, &match_recv, (void*)request);
724     }
725
726   if(request->action){
727     MPI_Request req = (MPI_Request)SIMIX_comm_get_src_data(request->action);
728     *flag = 1;
729     if(status != MPI_STATUS_IGNORE) {
730       status->MPI_SOURCE = smpi_group_rank(smpi_comm_group(comm), req->src);
731       status->MPI_TAG = req->tag;
732       status->MPI_ERROR = MPI_SUCCESS;
733       status->count = req->real_size;
734     }
735     nsleeps=1;//reset the number of sleeps we will do next time
736   }
737   else {
738       *flag = 0;
739       nsleeps++;
740   }
741   smpi_mpi_request_free(&request);
742
743   return;
744 }
745
746 void smpi_mpi_wait(MPI_Request * request, MPI_Status * status)
747 {
748   print_request("Waiting", *request);
749   if ((*request)->action != NULL) { // this is not a detached send
750     simcall_comm_wait((*request)->action, -1.0);
751   }
752   finish_wait(request, status);
753
754   // FIXME for a detached send, finish_wait is not called:
755 }
756
757 int smpi_mpi_waitany(int count, MPI_Request requests[],
758                      MPI_Status * status)
759 {
760   xbt_dynar_t comms;
761   int i, size, index;
762   int *map;
763
764   index = MPI_UNDEFINED;
765   if(count > 0) {
766     // Wait for a request to complete
767     comms = xbt_dynar_new(sizeof(smx_action_t), NULL);
768     map = xbt_new(int, count);
769     size = 0;
770     XBT_DEBUG("Wait for one of %d", count);
771     for(i = 0; i < count; i++) {
772       if(requests[i] != MPI_REQUEST_NULL) {
773         if (requests[i]->action != NULL) {
774           XBT_DEBUG("Waiting any %p ", requests[i]);
775           xbt_dynar_push(comms, &requests[i]->action);
776           map[size] = i;
777           size++;
778         }else{
779          //This is a finished detached request, let's return this one
780          size=0;//so we free the dynar but don't do the waitany call
781          index=i;
782          finish_wait(&requests[i], status);//cleanup if refcount = 0
783          requests[i]=MPI_REQUEST_NULL;//set to null
784          break;
785          }
786       }
787     }
788     if(size > 0) {
789       i = simcall_comm_waitany(comms);
790
791       // not MPI_UNDEFINED, as this is a simix return code
792       if (i != -1) {
793         index = map[i];
794         finish_wait(&requests[index], status);
795       }
796     }
797     xbt_free(map);
798     xbt_dynar_free(&comms);
799   }
800
801   if (index==MPI_UNDEFINED)
802     smpi_empty_status(status);
803
804   return index;
805 }
806
807 int smpi_mpi_waitall(int count, MPI_Request requests[],
808                       MPI_Status status[])
809 {
810   int  index, c;
811   MPI_Status stat;
812   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
813   int retvalue = MPI_SUCCESS;
814   //tag invalid requests in the set
815   if (status != MPI_STATUSES_IGNORE) {
816     for (c = 0; c < count; c++) {
817       if (requests[c] == MPI_REQUEST_NULL || requests[c]->dst == MPI_PROC_NULL) {
818         smpi_empty_status(&status[c]);
819       } else if (requests[c]->src == MPI_PROC_NULL) {
820         smpi_empty_status(&status[c]);
821         status[c].MPI_SOURCE = MPI_PROC_NULL;
822       }
823     }
824   }
825   for(c = 0; c < count; c++) {
826     if (MC_is_active()) {
827       smpi_mpi_wait(&requests[c], pstat);
828       index = c;
829     } else {
830       index = smpi_mpi_waitany(count, requests, pstat);
831       if (index == MPI_UNDEFINED)
832         break;
833       requests[index]=MPI_REQUEST_NULL;
834     }
835     if (status != MPI_STATUSES_IGNORE) {
836       status[index] = *pstat;
837       if (status[index].MPI_ERROR == MPI_ERR_TRUNCATE)
838         retvalue = MPI_ERR_IN_STATUS;
839     }
840   }
841
842   return retvalue;
843 }
844
845 int smpi_mpi_waitsome(int incount, MPI_Request requests[], int *indices,
846                       MPI_Status status[])
847 {
848   int i, count, index;
849   MPI_Status stat;
850   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
851
852   count = 0;
853   for(i = 0; i < incount; i++)
854   {
855     index=smpi_mpi_waitany(incount, requests, pstat);
856     if(index!=MPI_UNDEFINED){
857       indices[count] = index;
858       count++;
859       if(status != MPI_STATUSES_IGNORE) {
860         status[index] = *pstat;
861       }
862      requests[index]=MPI_REQUEST_NULL;
863     }else{
864       return MPI_UNDEFINED;
865     }
866   }
867   return count;
868 }
869
870 int smpi_mpi_testsome(int incount, MPI_Request requests[], int *indices,
871                       MPI_Status status[])
872 {
873   int i, count, count_dead;
874   MPI_Status stat;
875   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
876
877   count = 0;
878   count_dead = 0;
879   for(i = 0; i < incount; i++) {
880     if((requests[i] != MPI_REQUEST_NULL)) {
881       if(smpi_mpi_test(&requests[i], pstat)) {
882          indices[count] = i;
883          count++;
884          if(status != MPI_STATUSES_IGNORE) {
885            status[i] = *pstat;
886          }
887          requests[i]=MPI_REQUEST_NULL;
888
889       }
890     }else{
891       count_dead++;
892     }
893   }
894   if(count_dead==incount)return MPI_UNDEFINED;
895   else return count;
896 }
897
898 void smpi_mpi_bcast(void *buf, int count, MPI_Datatype datatype, int root,
899                     MPI_Comm comm)
900 {
901   // arity=2: a binary tree, arity=4 seem to be a good setting (see P2P-MPI))
902   nary_tree_bcast(buf, count, datatype, root, comm, 4);
903 }
904
905 void smpi_mpi_barrier(MPI_Comm comm)
906 {
907   // arity=2: a binary tree, arity=4 seem to be a good setting (see P2P-MPI))
908   nary_tree_barrier(comm, 4);
909 }
910
911 void smpi_mpi_gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
912                      void *recvbuf, int recvcount, MPI_Datatype recvtype,
913                      int root, MPI_Comm comm)
914 {
915   int system_tag = COLL_TAG_GATHER;
916   int rank, size, src, index;
917   MPI_Aint lb = 0, recvext = 0;
918   MPI_Request *requests;
919
920   rank = smpi_comm_rank(comm);
921   size = smpi_comm_size(comm);
922   if(rank != root) {
923     // Send buffer to root
924     smpi_mpi_send(sendbuf, sendcount, sendtype, root, system_tag, comm);
925   } else {
926     // FIXME: check for errors
927     smpi_datatype_extent(recvtype, &lb, &recvext);
928     // Local copy from root
929     smpi_datatype_copy(sendbuf, sendcount, sendtype,
930                        (char *)recvbuf + root * recvcount * recvext, recvcount, recvtype);
931     // Receive buffers from senders
932     requests = xbt_new(MPI_Request, size - 1);
933     index = 0;
934     for(src = 0; src < size; src++) {
935       if(src != root) {
936         requests[index] = smpi_irecv_init((char *)recvbuf + src * recvcount * recvext,
937                                           recvcount, recvtype,
938                                           src, system_tag, comm);
939         index++;
940       }
941     }
942     // Wait for completion of irecv's.
943     smpi_mpi_startall(size - 1, requests);
944     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
945     xbt_free(requests);
946   }
947 }
948
949
950 void smpi_mpi_reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
951                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
952 {
953     int i, size, count;
954     int *displs;
955     int rank = smpi_process_index();
956     void *tmpbuf;
957
958     /* arbitrarily choose root as rank 0 */
959     size = smpi_comm_size(comm);
960     count = 0;
961     displs = xbt_new(int, size);
962     for (i = 0; i < size; i++) {
963       displs[i] = count;
964       count += recvcounts[i];
965     }
966     tmpbuf=(void*)xbt_malloc(count*smpi_datatype_get_extent(datatype));
967     mpi_coll_reduce_fun(sendbuf, tmpbuf, count, datatype, op, 0, comm);
968     smpi_mpi_scatterv(tmpbuf, recvcounts, displs, datatype, recvbuf,
969                       recvcounts[rank], datatype, 0, comm);
970     xbt_free(displs);
971     xbt_free(tmpbuf);
972 }
973
974 void smpi_mpi_gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
975                       void *recvbuf, int *recvcounts, int *displs,
976                       MPI_Datatype recvtype, int root, MPI_Comm comm)
977 {
978   int system_tag = COLL_TAG_GATHERV;
979   int rank, size, src, index;
980   MPI_Aint lb = 0, recvext = 0;
981   MPI_Request *requests;
982
983   rank = smpi_comm_rank(comm);
984   size = smpi_comm_size(comm);
985   if(rank != root) {
986     // Send buffer to root
987     smpi_mpi_send(sendbuf, sendcount, sendtype, root, system_tag, comm);
988   } else {
989     // FIXME: check for errors
990     smpi_datatype_extent(recvtype, &lb, &recvext);
991     // Local copy from root
992     smpi_datatype_copy(sendbuf, sendcount, sendtype,
993                        (char *)recvbuf + displs[root] * recvext,
994                        recvcounts[root], recvtype);
995     // Receive buffers from senders
996     requests = xbt_new(MPI_Request, size - 1);
997     index = 0;
998     for(src = 0; src < size; src++) {
999       if(src != root) {
1000         requests[index] =
1001           smpi_irecv_init((char *)recvbuf + displs[src] * recvext,
1002                           recvcounts[src], recvtype, src, system_tag, comm);
1003         index++;
1004       }
1005     }
1006     // Wait for completion of irecv's.
1007     smpi_mpi_startall(size - 1, requests);
1008     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
1009     xbt_free(requests);
1010   }
1011 }
1012
1013 void smpi_mpi_allgather(void *sendbuf, int sendcount,
1014                         MPI_Datatype sendtype, void *recvbuf,
1015                         int recvcount, MPI_Datatype recvtype,
1016                         MPI_Comm comm)
1017 {
1018   int system_tag = COLL_TAG_ALLGATHER;
1019   int rank, size, other, index;
1020   MPI_Aint lb = 0, recvext = 0;
1021   MPI_Request *requests;
1022
1023   rank = smpi_comm_rank(comm);
1024   size = smpi_comm_size(comm);
1025   // FIXME: check for errors
1026   smpi_datatype_extent(recvtype, &lb, &recvext);
1027   // Local copy from self
1028   smpi_datatype_copy(sendbuf, sendcount, sendtype,
1029                      (char *)recvbuf + rank * recvcount * recvext, recvcount,
1030                      recvtype);
1031   // Send/Recv buffers to/from others;
1032   requests = xbt_new(MPI_Request, 2 * (size - 1));
1033   index = 0;
1034   for(other = 0; other < size; other++) {
1035     if(other != rank) {
1036       requests[index] =
1037         smpi_isend_init(sendbuf, sendcount, sendtype, other, system_tag,
1038                         comm);
1039       index++;
1040       requests[index] = smpi_irecv_init((char *)recvbuf + other * recvcount * recvext,
1041                                         recvcount, recvtype, other,
1042                                         system_tag, comm);
1043       index++;
1044     }
1045   }
1046   // Wait for completion of all comms.
1047   smpi_mpi_startall(2 * (size - 1), requests);
1048   smpi_mpi_waitall(2 * (size - 1), requests, MPI_STATUS_IGNORE);
1049   xbt_free(requests);
1050 }
1051
1052 void smpi_mpi_allgatherv(void *sendbuf, int sendcount,
1053                          MPI_Datatype sendtype, void *recvbuf,
1054                          int *recvcounts, int *displs,
1055                          MPI_Datatype recvtype, MPI_Comm comm)
1056 {
1057   int system_tag = COLL_TAG_ALLGATHERV;
1058   int rank, size, other, index;
1059   MPI_Aint lb = 0, recvext = 0;
1060   MPI_Request *requests;
1061
1062   rank = smpi_comm_rank(comm);
1063   size = smpi_comm_size(comm);
1064   // FIXME: check for errors
1065   smpi_datatype_extent(recvtype, &lb, &recvext);
1066   // Local copy from self
1067   smpi_datatype_copy(sendbuf, sendcount, sendtype,
1068                      (char *)recvbuf + displs[rank] * recvext,
1069                      recvcounts[rank], recvtype);
1070   // Send buffers to others;
1071   requests = xbt_new(MPI_Request, 2 * (size - 1));
1072   index = 0;
1073   for(other = 0; other < size; other++) {
1074     if(other != rank) {
1075       requests[index] =
1076         smpi_isend_init(sendbuf, sendcount, sendtype, other, system_tag,
1077                         comm);
1078       index++;
1079       requests[index] =
1080         smpi_irecv_init((char *)recvbuf + displs[other] * recvext, recvcounts[other],
1081                         recvtype, other, system_tag, comm);
1082       index++;
1083     }
1084   }
1085   // Wait for completion of all comms.
1086   smpi_mpi_startall(2 * (size - 1), requests);
1087   smpi_mpi_waitall(2 * (size - 1), requests, MPI_STATUS_IGNORE);
1088   xbt_free(requests);
1089 }
1090
1091 void smpi_mpi_scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1092                       void *recvbuf, int recvcount, MPI_Datatype recvtype,
1093                       int root, MPI_Comm comm)
1094 {
1095   int system_tag = COLL_TAG_SCATTER;
1096   int rank, size, dst, index;
1097   MPI_Aint lb = 0, sendext = 0;
1098   MPI_Request *requests;
1099
1100   rank = smpi_comm_rank(comm);
1101   size = smpi_comm_size(comm);
1102   if(rank != root) {
1103     // Recv buffer from root
1104     smpi_mpi_recv(recvbuf, recvcount, recvtype, root, system_tag, comm,
1105                   MPI_STATUS_IGNORE);
1106   } else {
1107     // FIXME: check for errors
1108     smpi_datatype_extent(sendtype, &lb, &sendext);
1109     // Local copy from root
1110     smpi_datatype_copy((char *)sendbuf + root * sendcount * sendext,
1111                        sendcount, sendtype, recvbuf, recvcount, recvtype);
1112     // Send buffers to receivers
1113     requests = xbt_new(MPI_Request, size - 1);
1114     index = 0;
1115     for(dst = 0; dst < size; dst++) {
1116       if(dst != root) {
1117         requests[index] = smpi_isend_init((char *)sendbuf + dst * sendcount * sendext,
1118                                           sendcount, sendtype, dst,
1119                                           system_tag, comm);
1120         index++;
1121       }
1122     }
1123     // Wait for completion of isend's.
1124     smpi_mpi_startall(size - 1, requests);
1125     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
1126     xbt_free(requests);
1127   }
1128 }
1129
1130 void smpi_mpi_scatterv(void *sendbuf, int *sendcounts, int *displs,
1131                        MPI_Datatype sendtype, void *recvbuf, int recvcount,
1132                        MPI_Datatype recvtype, int root, MPI_Comm comm)
1133 {
1134   int system_tag = COLL_TAG_SCATTERV;
1135   int rank, size, dst, index;
1136   MPI_Aint lb = 0, sendext = 0;
1137   MPI_Request *requests;
1138
1139   rank = smpi_comm_rank(comm);
1140   size = smpi_comm_size(comm);
1141   if(rank != root) {
1142     // Recv buffer from root
1143     smpi_mpi_recv(recvbuf, recvcount, recvtype, root, system_tag, comm,
1144                   MPI_STATUS_IGNORE);
1145   } else {
1146     // FIXME: check for errors
1147     smpi_datatype_extent(sendtype, &lb, &sendext);
1148     // Local copy from root
1149     smpi_datatype_copy((char *)sendbuf + displs[root] * sendext, sendcounts[root],
1150                        sendtype, recvbuf, recvcount, recvtype);
1151     // Send buffers to receivers
1152     requests = xbt_new(MPI_Request, size - 1);
1153     index = 0;
1154     for(dst = 0; dst < size; dst++) {
1155       if(dst != root) {
1156         requests[index] =
1157           smpi_isend_init((char *)sendbuf + displs[dst] * sendext, sendcounts[dst],
1158                           sendtype, dst, system_tag, comm);
1159         index++;
1160       }
1161     }
1162     // Wait for completion of isend's.
1163     smpi_mpi_startall(size - 1, requests);
1164     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
1165     xbt_free(requests);
1166   }
1167 }
1168
1169 void smpi_mpi_reduce(void *sendbuf, void *recvbuf, int count,
1170                      MPI_Datatype datatype, MPI_Op op, int root,
1171                      MPI_Comm comm)
1172 {
1173   int system_tag = COLL_TAG_REDUCE;
1174   int rank, size, src, index;
1175   MPI_Aint lb = 0, dataext = 0;
1176   MPI_Request *requests;
1177   void **tmpbufs;
1178
1179   char* sendtmpbuf = (char*) sendbuf;
1180   if( sendbuf == MPI_IN_PLACE ) {
1181       sendtmpbuf = (char *)recvbuf;
1182   }
1183
1184   rank = smpi_comm_rank(comm);
1185   size = smpi_comm_size(comm);
1186   //non commutative case, use a working algo from openmpi
1187   if(!smpi_op_is_commute(op)){
1188     smpi_coll_tuned_reduce_ompi_basic_linear(sendtmpbuf, recvbuf, count,
1189                      datatype, op, root, comm);
1190     return;
1191   }
1192   
1193   if(rank != root) {
1194     // Send buffer to root
1195     smpi_mpi_send(sendtmpbuf, count, datatype, root, system_tag, comm);
1196   } else {
1197     // FIXME: check for errors
1198     smpi_datatype_extent(datatype, &lb, &dataext);
1199     // Local copy from root
1200     if (sendtmpbuf && recvbuf)
1201       smpi_datatype_copy(sendtmpbuf, count, datatype, recvbuf, count, datatype);
1202     // Receive buffers from senders
1203     //TODO: make a MPI_barrier here ?
1204     requests = xbt_new(MPI_Request, size - 1);
1205     tmpbufs = xbt_new(void *, size - 1);
1206     index = 0;
1207     for(src = 0; src < size; src++) {
1208       if(src != root) {
1209         // FIXME: possibly overkill we we have contiguous/noncontiguous data
1210         //  mapping...
1211         tmpbufs[index] = xbt_malloc(count * dataext);
1212         requests[index] =
1213           smpi_irecv_init(tmpbufs[index], count, datatype, src,
1214                           system_tag, comm);
1215         index++;
1216       }
1217     }
1218     // Wait for completion of irecv's.
1219     smpi_mpi_startall(size - 1, requests);
1220     for(src = 0; src < size - 1; src++) {
1221       index = smpi_mpi_waitany(size - 1, requests, MPI_STATUS_IGNORE);
1222       XBT_DEBUG("finished waiting any request with index %d", index);
1223       if(index == MPI_UNDEFINED) {
1224         break;
1225       }
1226       if(op) /* op can be MPI_OP_NULL that does nothing */
1227         smpi_op_apply(op, tmpbufs[index], recvbuf, &count, &datatype);
1228     }
1229     for(index = 0; index < size - 1; index++) {
1230       xbt_free(tmpbufs[index]);
1231     }
1232     xbt_free(tmpbufs);
1233     xbt_free(requests);
1234   }
1235 }
1236
1237 void smpi_mpi_allreduce(void *sendbuf, void *recvbuf, int count,
1238                         MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1239 {
1240   smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1241   smpi_mpi_bcast(recvbuf, count, datatype, 0, comm);
1242 }
1243
1244 void smpi_mpi_scan(void *sendbuf, void *recvbuf, int count,
1245                    MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1246 {
1247   int system_tag = 888;
1248   int rank, size, other, index;
1249   MPI_Aint lb = 0, dataext = 0;
1250   MPI_Request *requests;
1251   void **tmpbufs;
1252
1253   rank = smpi_comm_rank(comm);
1254   size = smpi_comm_size(comm);
1255
1256   // FIXME: check for errors
1257   smpi_datatype_extent(datatype, &lb, &dataext);
1258
1259   // Local copy from self
1260   smpi_datatype_copy(sendbuf, count, datatype, recvbuf, count, datatype);
1261
1262   // Send/Recv buffers to/from others;
1263   requests = xbt_new(MPI_Request, size - 1);
1264   tmpbufs = xbt_new(void *, rank);
1265   index = 0;
1266   for(other = 0; other < rank; other++) {
1267     // FIXME: possibly overkill we we have contiguous/noncontiguous data
1268     // mapping...
1269     tmpbufs[index] = xbt_malloc(count * dataext);
1270     requests[index] =
1271       smpi_irecv_init(tmpbufs[index], count, datatype, other, system_tag,
1272                       comm);
1273     index++;
1274   }
1275   for(other = rank + 1; other < size; other++) {
1276     requests[index] =
1277       smpi_isend_init(sendbuf, count, datatype, other, system_tag, comm);
1278     index++;
1279   }
1280   // Wait for completion of all comms.
1281   smpi_mpi_startall(size - 1, requests);
1282   for(other = 0; other < size - 1; other++) {
1283     index = smpi_mpi_waitany(size - 1, requests, MPI_STATUS_IGNORE);
1284     if(index == MPI_UNDEFINED) {
1285       break;
1286     }
1287     if(index < rank) {
1288       // #Request is below rank: it's a irecv
1289       smpi_op_apply(op, tmpbufs[index], recvbuf, &count, &datatype);
1290     }
1291   }
1292   for(index = 0; index < rank; index++) {
1293     xbt_free(tmpbufs[index]);
1294   }
1295   xbt_free(tmpbufs);
1296   xbt_free(requests);
1297 }