Logo AND Algorithmique Numérique Distribuée

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