Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
model-checker : need to switch between raw heap and std heap to remove ignore region
[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
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
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 =
124         parse_factor(sg_cfg_get_string("smpi/os"));
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_or(double size)
143 {
144   if (!smpi_or_values)
145     smpi_or_values =
146         parse_factor(sg_cfg_get_string("smpi/or"));
147
148   unsigned int iter = 0;
149   s_smpi_factor_t fact;
150   double current=0.0;
151   xbt_dynar_foreach(smpi_or_values, iter, fact) {
152     if (size <= fact.factor) {
153         XBT_DEBUG("or : %lf <= %ld return %f", size, fact.factor, current);
154       return current;
155     }else
156       current=fact.values[0]+fact.values[1]*size;
157   }
158   XBT_DEBUG("or : %lf > %ld return %f", size, fact.factor, current);
159
160   return current;
161 }
162
163 static MPI_Request build_request(void *buf, int count,
164                                  MPI_Datatype datatype, int src, int dst,
165                                  int tag, MPI_Comm comm, unsigned flags)
166 {
167   MPI_Request request;
168
169   void *old_buf = NULL;
170
171   request = xbt_new(s_smpi_mpi_request_t, 1);
172
173   s_smpi_subtype_t *subtype = datatype->substruct;
174
175   if(datatype->has_subtype == 1){
176     // This part handles the problem of non-contiguous memory
177     old_buf = buf;
178     buf = xbt_malloc(count*smpi_datatype_size(datatype));
179     if (flags & SEND) {
180       subtype->serialize(old_buf, buf, count, datatype->substruct);
181     }
182   }
183
184   request->buf = buf;
185   // This part handles the problem of non-contiguous memory (for the
186   // unserialisation at the reception)
187   request->old_buf = old_buf;
188   request->old_type = datatype;
189
190   request->size = smpi_datatype_size(datatype) * count;
191   request->src = src;
192   request->dst = dst;
193   request->tag = tag;
194   request->comm = comm;
195   request->action = NULL;
196   request->flags = flags;
197   request->detached = 0;
198   request->detached_sender = NULL;
199
200   request->truncated = 0;
201   request->real_size = 0;
202   request->real_tag = 0;
203
204   request->refcount=1;
205 #ifdef HAVE_TRACING
206   request->send = 0;
207   request->recv = 0;
208 #endif
209   if (flags & SEND) smpi_datatype_unuse(datatype);
210
211   return request;
212 }
213
214
215 void smpi_empty_status(MPI_Status * status) {
216   if(status != MPI_STATUS_IGNORE) {
217       status->MPI_SOURCE=MPI_ANY_SOURCE;
218       status->MPI_TAG=MPI_ANY_TAG;
219       status->count=0;
220   }
221 }
222
223 void smpi_action_trace_run(char *path)
224 {
225   char *name;
226   xbt_dynar_t todo;
227   xbt_dict_cursor_t cursor;
228
229   action_fp=NULL;
230   if (path) {
231     action_fp = fopen(path, "r");
232     xbt_assert(action_fp != NULL, "Cannot open %s: %s", path,
233                strerror(errno));
234   }
235
236   if (!xbt_dict_is_empty(action_queues)) {
237     XBT_WARN
238       ("Not all actions got consumed. If the simulation ended successfully (without deadlock), you may want to add new processes to your deployment file.");
239
240
241     xbt_dict_foreach(action_queues, cursor, name, todo) {
242       XBT_WARN("Still %lu actions for %s", xbt_dynar_length(todo), name);
243     }
244   }
245
246   if (path)
247     fclose(action_fp);
248   xbt_dict_free(&action_queues);
249   action_queues = xbt_dict_new_homogeneous(NULL);
250 }
251
252 static void smpi_mpi_request_free_voidp(void* request)
253 {
254   MPI_Request req = request;
255   smpi_mpi_request_free(&req);
256 }
257
258 /* MPI Low level calls */
259 MPI_Request smpi_mpi_send_init(void *buf, int count, MPI_Datatype datatype,
260                                int dst, int tag, MPI_Comm comm)
261 {
262   MPI_Request request =
263     build_request(buf, count, datatype, smpi_comm_rank(comm), dst, tag,
264                   comm, PERSISTENT | SEND);
265   request->refcount++;
266   return request;
267 }
268
269 MPI_Request smpi_mpi_recv_init(void *buf, int count, MPI_Datatype datatype,
270                                int src, int tag, MPI_Comm comm)
271 {
272   MPI_Request request =
273     build_request(buf, count, datatype, src, smpi_comm_rank(comm), tag,
274                   comm, PERSISTENT | RECV);
275   request->refcount++;
276   return request;
277 }
278
279 void smpi_mpi_start(MPI_Request request)
280 {
281   smx_rdv_t mailbox;
282
283   xbt_assert(!request->action,
284              "Cannot (re)start a non-finished communication");
285   if(request->flags & RECV) {
286     print_request("New recv", request);
287     if (request->size < sg_cfg_get_int("smpi/async_small_thres"))
288       mailbox = smpi_process_mailbox_small();
289     else
290       mailbox = smpi_process_mailbox();
291     // we make a copy here, as the size is modified by simix, and we may reuse the request in another receive later
292     request->real_size=request->size;
293     smpi_datatype_use(request->old_type);
294     request->action = simcall_comm_irecv(mailbox, request->buf, &request->real_size, &match_recv, request);
295
296     double sleeptime = smpi_or(request->size);
297     if(sleeptime!=0.0){
298         simcall_process_sleep(sleeptime);
299         XBT_DEBUG("receiving size of %ld : sleep %lf ", request->size, smpi_or(request->size));
300     }
301
302   } else {
303
304     int receiver = smpi_group_index(smpi_comm_group(request->comm), request->dst);
305 /*    if(receiver == MPI_UNDEFINED) {*/
306 /*      XBT_WARN("Trying to send a message to a wrong rank");*/
307 /*      return;*/
308 /*    }*/
309     print_request("New send", request);
310     if (request->size < sg_cfg_get_int("smpi/async_small_thres")) { // eager mode
311       mailbox = smpi_process_remote_mailbox_small(receiver);
312     }else{
313       XBT_DEBUG("Send request %p is not in the permanent receive mailbox (buf: %p)",request,request->buf);
314       mailbox = smpi_process_remote_mailbox(receiver);
315     }
316     if (request->size < 64*1024 ) { //(FIXME: this limit should be configurable)
317       void *oldbuf = NULL;
318       request->detached = 1;
319       request->refcount++;
320       if(request->old_type->has_subtype == 0){
321         oldbuf = request->buf;
322         if (oldbuf){
323           request->buf = xbt_malloc(request->size);
324           memcpy(request->buf,oldbuf,request->size);
325         }
326       }
327       XBT_DEBUG("Send request %p is detached; buf %p copied into %p",request,oldbuf,request->buf);
328     }
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     double sleeptime = smpi_os(request->size);
333     if(sleeptime!=0.0){
334         simcall_process_sleep(sleeptime);
335         XBT_DEBUG("sending size of %ld : sleep %lf ", request->size, smpi_os(request->size));
336     }
337     request->action =
338       simcall_comm_isend(mailbox, request->size, -1.0,
339                          request->buf, request->real_size,
340                          &match_send,
341                          &smpi_mpi_request_free_voidp, // how to free the userdata if a detached send fails
342                          request,
343                          // detach if msg size < eager/rdv switch limit
344                          request->detached);
345
346 #ifdef HAVE_TRACING
347     /* FIXME: detached sends are not traceable (request->action == NULL) */
348     if (request->action)
349       simcall_set_category(request->action, TRACE_internal_smpi_get_category());
350 #endif
351
352   }
353
354 }
355
356 void smpi_mpi_startall(int count, MPI_Request * requests)
357 {
358   int i;
359
360   for(i = 0; i < count; i++) {
361     smpi_mpi_start(requests[i]);
362   }
363 }
364
365 void smpi_mpi_request_free(MPI_Request * request)
366 {
367
368   if((*request) != MPI_REQUEST_NULL){
369     (*request)->refcount--;
370     if((*request)->refcount<0) xbt_die("wrong refcount");
371
372     if((*request)->refcount==0){
373         xbt_free(*request);
374         *request = MPI_REQUEST_NULL;
375     }
376   }else{
377       xbt_die("freeing an already free request");
378   }
379 }
380
381 MPI_Request smpi_isend_init(void *buf, int count, MPI_Datatype datatype,
382                             int dst, int tag, MPI_Comm comm)
383 {
384   MPI_Request request =
385     build_request(buf, count, datatype, smpi_comm_rank(comm), dst, tag,
386                   comm, NON_PERSISTENT | SEND);
387
388   return request;
389 }
390
391 MPI_Request smpi_mpi_isend(void *buf, int count, MPI_Datatype datatype,
392                            int dst, int tag, MPI_Comm comm)
393 {
394   MPI_Request request =
395       build_request(buf, count, datatype, smpi_comm_rank(comm), dst, tag,
396                     comm, NON_PERSISTENT | SEND);
397
398   smpi_mpi_start(request);
399   return request;
400 }
401
402 MPI_Request smpi_irecv_init(void *buf, int count, MPI_Datatype datatype,
403                             int src, int tag, MPI_Comm comm)
404 {
405   MPI_Request request =
406     build_request(buf, count, datatype, src, smpi_comm_rank(comm), tag,
407                   comm, NON_PERSISTENT | RECV);
408   return request;
409 }
410
411 MPI_Request smpi_mpi_irecv(void *buf, int count, MPI_Datatype datatype,
412                            int src, int tag, MPI_Comm comm)
413 {
414   MPI_Request request =
415       build_request(buf, count, datatype, src, smpi_comm_rank(comm), tag,
416                     comm, NON_PERSISTENT | RECV);
417
418   smpi_mpi_start(request);
419   return request;
420 }
421
422 void smpi_mpi_recv(void *buf, int count, MPI_Datatype datatype, int src,
423                    int tag, MPI_Comm comm, MPI_Status * status)
424 {
425   MPI_Request request;
426   request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
427   smpi_mpi_wait(&request, status);
428 }
429
430
431
432 void smpi_mpi_send(void *buf, int count, MPI_Datatype datatype, int dst,
433                    int tag, MPI_Comm comm)
434 {
435   MPI_Request request;
436   request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
437   smpi_mpi_wait(&request, MPI_STATUS_IGNORE);
438
439 }
440
441 void smpi_mpi_sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
442                        int dst, int sendtag, void *recvbuf, int recvcount,
443                        MPI_Datatype recvtype, int src, int recvtag,
444                        MPI_Comm comm, MPI_Status * status)
445 {
446   MPI_Request requests[2];
447   MPI_Status stats[2];
448
449   requests[0] =
450     smpi_isend_init(sendbuf, sendcount, sendtype, dst, sendtag, comm);
451   requests[1] =
452     smpi_irecv_init(recvbuf, recvcount, recvtype, src, recvtag, comm);
453   smpi_mpi_startall(2, requests);
454   smpi_mpi_waitall(2, requests, stats);
455   if(status != MPI_STATUS_IGNORE) {
456     // Copy receive status
457     memcpy(status, &stats[1], sizeof(MPI_Status));
458   }
459 }
460
461 int smpi_mpi_get_count(MPI_Status * status, MPI_Datatype datatype)
462 {
463   return status->count / smpi_datatype_size(datatype);
464 }
465
466 static void finish_wait(MPI_Request * request, MPI_Status * status)
467 {
468   MPI_Request req = *request;
469   if(!(req->detached && req->flags & SEND)){
470     if(status != MPI_STATUS_IGNORE) {
471       status->MPI_SOURCE = req->src == MPI_ANY_SOURCE ? req->real_src : req->src;
472       status->MPI_TAG = req->tag == MPI_ANY_TAG ? req->real_tag : req->tag;
473       if(req->truncated)
474       status->MPI_ERROR = MPI_ERR_TRUNCATE;
475       else status->MPI_ERROR = MPI_SUCCESS ;
476       // this handles the case were size in receive differs from size in send
477       // FIXME: really this should just contain the count of receive-type blocks,
478       // right?
479       status->count = req->real_size;
480     }
481
482     print_request("Finishing", req);
483     MPI_Datatype datatype = req->old_type;
484
485     if(datatype->has_subtype == 1){
486         // This part handles the problem of non-contignous memory
487         // the unserialization at the reception
488       s_smpi_subtype_t *subtype = datatype->substruct;
489       if(req->flags & RECV) {
490         subtype->unserialize(req->buf, req->old_buf, req->real_size/smpi_datatype_size(datatype) , datatype->substruct);
491       }
492       if(req->detached == 0) free(req->buf);
493     }
494     smpi_datatype_unuse(datatype);
495   }
496
497   if(req->detached_sender!=NULL){
498     smpi_mpi_request_free(&(req->detached_sender));
499   }
500
501   if(req->flags & NON_PERSISTENT) {
502     smpi_mpi_request_free(request);
503   } else {
504     req->action = NULL;
505   }
506 }
507
508 int smpi_mpi_test(MPI_Request * request, MPI_Status * status) {
509   int flag;
510
511   //assume that request is not MPI_REQUEST_NULL (filtered in PMPI_Test or smpi_mpi_testall before)
512   if ((*request)->action == NULL)
513     flag = 1;
514   else
515     flag = simcall_comm_test((*request)->action);
516   if(flag) {
517     (*request)->refcount++;
518     finish_wait(request, status);
519   }else{
520     smpi_empty_status(status);
521   }
522   return flag;
523 }
524
525 int smpi_mpi_testany(int count, MPI_Request requests[], int *index,
526                      MPI_Status * status)
527 {
528   xbt_dynar_t comms;
529   int i, flag, size;
530   int* map;
531
532   *index = MPI_UNDEFINED;
533   flag = 0;
534   if(count > 0) {
535     comms = xbt_dynar_new(sizeof(smx_action_t), NULL);
536     map = xbt_new(int, count);
537     size = 0;
538     for(i = 0; i < count; i++) {
539       if((requests[i]!=MPI_REQUEST_NULL) && requests[i]->action) {
540          xbt_dynar_push(comms, &requests[i]->action);
541          map[size] = i;
542          size++;
543       }
544     }
545     if(size > 0) {
546       i = simcall_comm_testany(comms);
547       // not MPI_UNDEFINED, as this is a simix return code
548       if(i != -1) {
549         *index = map[i];
550         finish_wait(&requests[*index], status);
551         flag = 1;
552       }
553     }else{
554         //all requests are null or inactive, return true
555         flag=1;
556         smpi_empty_status(status);
557     }
558     xbt_free(map);
559     xbt_dynar_free(&comms);
560   }
561
562   return flag;
563 }
564
565
566 int smpi_mpi_testall(int count, MPI_Request requests[],
567                      MPI_Status status[])
568 {
569   MPI_Status stat;
570   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
571   int flag=1;
572   int i;
573   for(i=0; i<count; i++){
574     if(requests[i]!= MPI_REQUEST_NULL){
575       if (smpi_mpi_test(&requests[i], pstat)!=1){
576         flag=0;
577       }
578     }else{
579       smpi_empty_status(pstat);
580     }
581     if(status != MPI_STATUSES_IGNORE) {
582       memcpy(&status[i], pstat, sizeof(*pstat));
583     }
584   }
585   return flag;
586 }
587
588 void smpi_mpi_probe(int source, int tag, MPI_Comm comm, MPI_Status* status){
589   int flag=0;
590   //FIXME find another wait to avoid busy waiting ?
591   // the issue here is that we have to wait on a nonexistent comm
592   while(flag==0){
593     smpi_mpi_iprobe(source, tag, comm, &flag, status);
594     XBT_DEBUG("Busy Waiting on probing : %d", flag);
595     if(!flag) {
596       simcall_process_sleep(0.0001);
597     }
598   }
599 }
600
601 void smpi_mpi_iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status){
602   MPI_Request request =build_request(NULL, 0, MPI_CHAR, source, smpi_comm_rank(comm), tag,
603             comm, NON_PERSISTENT | RECV);
604
605   // behave like a receive, but don't do it
606   smx_rdv_t mailbox;
607
608   print_request("New iprobe", request);
609   // We have to test both mailboxes as we don't know if we will receive one one or another
610     if (sg_cfg_get_int("smpi/async_small_thres")>0){
611         mailbox = smpi_process_mailbox_small();
612         XBT_DEBUG("trying to probe the perm recv mailbox");
613         request->action = simcall_comm_iprobe(mailbox, request->src, request->tag, &match_recv, (void*)request);
614     }
615     if (request->action==NULL){
616         mailbox = smpi_process_mailbox();
617         XBT_DEBUG("trying to probe the other mailbox");
618         request->action = simcall_comm_iprobe(mailbox, request->src, request->tag, &match_recv, (void*)request);
619     }
620
621   if(request->action){
622     MPI_Request req = (MPI_Request)SIMIX_comm_get_src_data(request->action);
623     *flag = 1;
624     if(status != MPI_STATUS_IGNORE) {
625       status->MPI_SOURCE = req->src;
626       status->MPI_TAG = req->tag;
627       status->MPI_ERROR = MPI_SUCCESS;
628       status->count = req->real_size;
629     }
630   }
631   else *flag = 0;
632   smpi_mpi_request_free(&request);
633
634   return;
635 }
636
637 void smpi_mpi_wait(MPI_Request * request, MPI_Status * status)
638 {
639   print_request("Waiting", *request);
640   if ((*request)->action != NULL) { // this is not a detached send
641     simcall_comm_wait((*request)->action, -1.0);
642   }
643   finish_wait(request, status);
644
645   // FIXME for a detached send, finish_wait is not called:
646 }
647
648 int smpi_mpi_waitany(int count, MPI_Request requests[],
649                      MPI_Status * status)
650 {
651   xbt_dynar_t comms;
652   int i, size, index;
653   int *map;
654
655   index = MPI_UNDEFINED;
656   if(count > 0) {
657     // Wait for a request to complete
658     comms = xbt_dynar_new(sizeof(smx_action_t), NULL);
659     map = xbt_new(int, count);
660     size = 0;
661     XBT_DEBUG("Wait for one of %d", count);
662     for(i = 0; i < count; i++) {
663       if(requests[i] != MPI_REQUEST_NULL) {
664         if (requests[i]->action != NULL) {
665           XBT_DEBUG("Waiting any %p ", requests[i]);
666           xbt_dynar_push(comms, &requests[i]->action);
667           map[size] = i;
668           size++;
669         }else{
670          //This is a finished detached request, let's return this one
671          size=0;//so we free the dynar but don't do the waitany call
672          index=i;
673          finish_wait(&requests[i], status);//cleanup if refcount = 0
674          requests[i]=MPI_REQUEST_NULL;//set to null
675          break;
676          }
677       }
678     }
679     if(size > 0) {
680       i = simcall_comm_waitany(comms);
681
682       // not MPI_UNDEFINED, as this is a simix return code
683       if (i != -1) {
684         index = map[i];
685         finish_wait(&requests[index], status);
686       }
687     }
688     xbt_free(map);
689     xbt_dynar_free(&comms);
690   }
691
692   if (index==MPI_UNDEFINED)
693     smpi_empty_status(status);
694
695   return index;
696 }
697
698 int smpi_mpi_waitall(int count, MPI_Request requests[],
699                       MPI_Status status[])
700 {
701   int  index, c;
702   MPI_Status stat;
703   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
704   int retvalue=MPI_SUCCESS;
705   //tag invalid requests in the set
706   for(c = 0; c < count; c++) {
707     if(requests[c]==MPI_REQUEST_NULL || requests[c]->dst == MPI_PROC_NULL ){
708       if(status != MPI_STATUSES_IGNORE)
709         smpi_empty_status(&status[c]);
710     }else if(requests[c]->src == MPI_PROC_NULL ){
711       if(status != MPI_STATUSES_IGNORE) {
712         smpi_empty_status(&status[c]);
713         status[c].MPI_SOURCE=MPI_PROC_NULL;
714       }
715     }
716   }
717   for(c = 0; c < count; c++) {
718       if(MC_is_active()) {
719         smpi_mpi_wait(&requests[c], pstat);
720         index = c;
721       } else {
722         index = smpi_mpi_waitany(count, requests, pstat);
723         if(index == MPI_UNDEFINED) {
724           break;
725        }
726       if(status != MPI_STATUSES_IGNORE) {
727         memcpy(&status[index], pstat, sizeof(*pstat));
728         if(status[index].MPI_ERROR==MPI_ERR_TRUNCATE)retvalue=MPI_ERR_IN_STATUS;
729
730       }
731     }
732   }
733
734   return retvalue;
735 }
736
737 int smpi_mpi_waitsome(int incount, MPI_Request requests[], int *indices,
738                       MPI_Status status[])
739 {
740   int i, count, index;
741   MPI_Status stat;
742   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
743
744   count = 0;
745   for(i = 0; i < incount; i++)
746   {
747     index=smpi_mpi_waitany(incount, requests, pstat);
748     if(index!=MPI_UNDEFINED){
749       indices[count] = index;
750       count++;
751       if(status != MPI_STATUSES_IGNORE) {
752         memcpy(&status[index], pstat, sizeof(*pstat));
753       }
754     }else{
755       return MPI_UNDEFINED;
756     }
757   }
758   return count;
759 }
760
761 int smpi_mpi_testsome(int incount, MPI_Request requests[], int *indices,
762                       MPI_Status status[])
763 {
764   int i, count, count_dead;
765   MPI_Status stat;
766   MPI_Status *pstat = status == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : &stat;
767
768   count = 0;
769   count_dead = 0;
770   for(i = 0; i < incount; i++) {
771     if((requests[i] != MPI_REQUEST_NULL)) {
772       if(smpi_mpi_test(&requests[i], pstat)) {
773          indices[count] = i;
774          count++;
775          if(status != MPI_STATUSES_IGNORE) {
776             memcpy(&status[i], pstat, sizeof(*pstat));
777          }
778       }
779     }else{
780       count_dead++;
781     }
782   }
783   if(count_dead==incount)return MPI_UNDEFINED;
784   else return count;
785 }
786
787 void smpi_mpi_bcast(void *buf, int count, MPI_Datatype datatype, int root,
788                     MPI_Comm comm)
789 {
790   // arity=2: a binary tree, arity=4 seem to be a good setting (see P2P-MPI))
791   nary_tree_bcast(buf, count, datatype, root, comm, 4);
792 }
793
794 void smpi_mpi_barrier(MPI_Comm comm)
795 {
796   // arity=2: a binary tree, arity=4 seem to be a good setting (see P2P-MPI))
797   nary_tree_barrier(comm, 4);
798 }
799
800 void smpi_mpi_gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
801                      void *recvbuf, int recvcount, MPI_Datatype recvtype,
802                      int root, MPI_Comm comm)
803 {
804   int system_tag = 666;
805   int rank, size, src, index;
806   MPI_Aint lb = 0, recvext = 0;
807   MPI_Request *requests;
808
809   rank = smpi_comm_rank(comm);
810   size = smpi_comm_size(comm);
811   if(rank != root) {
812     // Send buffer to root
813     smpi_mpi_send(sendbuf, sendcount, sendtype, root, system_tag, comm);
814   } else {
815     // FIXME: check for errors
816     smpi_datatype_extent(recvtype, &lb, &recvext);
817     // Local copy from root
818     smpi_datatype_copy(sendbuf, sendcount, sendtype,
819                        (char *)recvbuf + root * recvcount * recvext, recvcount, recvtype);
820     // Receive buffers from senders
821     requests = xbt_new(MPI_Request, size - 1);
822     index = 0;
823     for(src = 0; src < size; src++) {
824       if(src != root) {
825         requests[index] = smpi_irecv_init((char *)recvbuf + src * recvcount * recvext,
826                                           recvcount, recvtype,
827                                           src, system_tag, comm);
828         index++;
829       }
830     }
831     // Wait for completion of irecv's.
832     smpi_mpi_startall(size - 1, requests);
833     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
834     xbt_free(requests);
835   }
836 }
837
838 void smpi_mpi_gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
839                       void *recvbuf, int *recvcounts, int *displs,
840                       MPI_Datatype recvtype, int root, MPI_Comm comm)
841 {
842   int system_tag = 666;
843   int rank, size, src, index;
844   MPI_Aint lb = 0, recvext = 0;
845   MPI_Request *requests;
846
847   rank = smpi_comm_rank(comm);
848   size = smpi_comm_size(comm);
849   if(rank != root) {
850     // Send buffer to root
851     smpi_mpi_send(sendbuf, sendcount, sendtype, root, system_tag, comm);
852   } else {
853     // FIXME: check for errors
854     smpi_datatype_extent(recvtype, &lb, &recvext);
855     // Local copy from root
856     smpi_datatype_copy(sendbuf, sendcount, sendtype,
857                        (char *)recvbuf + displs[root] * recvext,
858                        recvcounts[root], recvtype);
859     // Receive buffers from senders
860     requests = xbt_new(MPI_Request, size - 1);
861     index = 0;
862     for(src = 0; src < size; src++) {
863       if(src != root) {
864         requests[index] =
865           smpi_irecv_init((char *)recvbuf + displs[src] * recvext,
866                           recvcounts[src], recvtype, src, system_tag, comm);
867         index++;
868       }
869     }
870     // Wait for completion of irecv's.
871     smpi_mpi_startall(size - 1, requests);
872     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
873     xbt_free(requests);
874   }
875 }
876
877 void smpi_mpi_allgather(void *sendbuf, int sendcount,
878                         MPI_Datatype sendtype, void *recvbuf,
879                         int recvcount, MPI_Datatype recvtype,
880                         MPI_Comm comm)
881 {
882   int system_tag = 666;
883   int rank, size, other, index;
884   MPI_Aint lb = 0, recvext = 0;
885   MPI_Request *requests;
886
887   rank = smpi_comm_rank(comm);
888   size = smpi_comm_size(comm);
889   // FIXME: check for errors
890   smpi_datatype_extent(recvtype, &lb, &recvext);
891   // Local copy from self
892   smpi_datatype_copy(sendbuf, sendcount, sendtype,
893                      (char *)recvbuf + rank * recvcount * recvext, recvcount,
894                      recvtype);
895   // Send/Recv buffers to/from others;
896   requests = xbt_new(MPI_Request, 2 * (size - 1));
897   index = 0;
898   for(other = 0; other < size; other++) {
899     if(other != rank) {
900       requests[index] =
901         smpi_isend_init(sendbuf, sendcount, sendtype, other, system_tag,
902                         comm);
903       index++;
904       requests[index] = smpi_irecv_init((char *)recvbuf + other * recvcount * recvext,
905                                         recvcount, recvtype, other,
906                                         system_tag, comm);
907       index++;
908     }
909   }
910   // Wait for completion of all comms.
911   smpi_mpi_startall(2 * (size - 1), requests);
912   smpi_mpi_waitall(2 * (size - 1), requests, MPI_STATUS_IGNORE);
913   xbt_free(requests);
914 }
915
916 void smpi_mpi_allgatherv(void *sendbuf, int sendcount,
917                          MPI_Datatype sendtype, void *recvbuf,
918                          int *recvcounts, int *displs,
919                          MPI_Datatype recvtype, MPI_Comm comm)
920 {
921   int system_tag = 666;
922   int rank, size, other, 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   // FIXME: check for errors
929   smpi_datatype_extent(recvtype, &lb, &recvext);
930   // Local copy from self
931   smpi_datatype_copy(sendbuf, sendcount, sendtype,
932                      (char *)recvbuf + displs[rank] * recvext,
933                      recvcounts[rank], recvtype);
934   // Send buffers to others;
935   requests = xbt_new(MPI_Request, 2 * (size - 1));
936   index = 0;
937   for(other = 0; other < size; other++) {
938     if(other != rank) {
939       requests[index] =
940         smpi_isend_init(sendbuf, sendcount, sendtype, other, system_tag,
941                         comm);
942       index++;
943       requests[index] =
944         smpi_irecv_init((char *)recvbuf + displs[other] * recvext, recvcounts[other],
945                         recvtype, other, system_tag, comm);
946       index++;
947     }
948   }
949   // Wait for completion of all comms.
950   smpi_mpi_startall(2 * (size - 1), requests);
951   smpi_mpi_waitall(2 * (size - 1), requests, MPI_STATUS_IGNORE);
952   xbt_free(requests);
953 }
954
955 void smpi_mpi_scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
956                       void *recvbuf, int recvcount, MPI_Datatype recvtype,
957                       int root, MPI_Comm comm)
958 {
959   int system_tag = 666;
960   int rank, size, dst, index;
961   MPI_Aint lb = 0, sendext = 0;
962   MPI_Request *requests;
963
964   rank = smpi_comm_rank(comm);
965   size = smpi_comm_size(comm);
966   if(rank != root) {
967     // Recv buffer from root
968     smpi_mpi_recv(recvbuf, recvcount, recvtype, root, system_tag, comm,
969                   MPI_STATUS_IGNORE);
970   } else {
971     // FIXME: check for errors
972     smpi_datatype_extent(sendtype, &lb, &sendext);
973     // Local copy from root
974     smpi_datatype_copy((char *)sendbuf + root * sendcount * sendext,
975                        sendcount, sendtype, recvbuf, recvcount, recvtype);
976     // Send buffers to receivers
977     requests = xbt_new(MPI_Request, size - 1);
978     index = 0;
979     for(dst = 0; dst < size; dst++) {
980       if(dst != root) {
981         requests[index] = smpi_isend_init((char *)sendbuf + dst * sendcount * sendext,
982                                           sendcount, sendtype, dst,
983                                           system_tag, comm);
984         index++;
985       }
986     }
987     // Wait for completion of isend's.
988     smpi_mpi_startall(size - 1, requests);
989     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
990     xbt_free(requests);
991   }
992 }
993
994 void smpi_mpi_scatterv(void *sendbuf, int *sendcounts, int *displs,
995                        MPI_Datatype sendtype, void *recvbuf, int recvcount,
996                        MPI_Datatype recvtype, int root, MPI_Comm comm)
997 {
998   int system_tag = 666;
999   int rank, size, dst, index;
1000   MPI_Aint lb = 0, sendext = 0;
1001   MPI_Request *requests;
1002
1003   rank = smpi_comm_rank(comm);
1004   size = smpi_comm_size(comm);
1005   if(rank != root) {
1006     // Recv buffer from root
1007     smpi_mpi_recv(recvbuf, recvcount, recvtype, root, system_tag, comm,
1008                   MPI_STATUS_IGNORE);
1009   } else {
1010     // FIXME: check for errors
1011     smpi_datatype_extent(sendtype, &lb, &sendext);
1012     // Local copy from root
1013     smpi_datatype_copy((char *)sendbuf + displs[root] * sendext, sendcounts[root],
1014                        sendtype, recvbuf, recvcount, recvtype);
1015     // Send buffers to receivers
1016     requests = xbt_new(MPI_Request, size - 1);
1017     index = 0;
1018     for(dst = 0; dst < size; dst++) {
1019       if(dst != root) {
1020         requests[index] =
1021           smpi_isend_init((char *)sendbuf + displs[dst] * sendext, sendcounts[dst],
1022                           sendtype, dst, system_tag, comm);
1023         index++;
1024       }
1025     }
1026     // Wait for completion of isend's.
1027     smpi_mpi_startall(size - 1, requests);
1028     smpi_mpi_waitall(size - 1, requests, MPI_STATUS_IGNORE);
1029     xbt_free(requests);
1030   }
1031 }
1032
1033 void smpi_mpi_reduce(void *sendbuf, void *recvbuf, int count,
1034                      MPI_Datatype datatype, MPI_Op op, int root,
1035                      MPI_Comm comm)
1036 {
1037   int system_tag = 666;
1038   int rank, size, src, index;
1039   MPI_Aint lb = 0, dataext = 0;
1040   MPI_Request *requests;
1041   void **tmpbufs;
1042
1043   rank = smpi_comm_rank(comm);
1044   size = smpi_comm_size(comm);
1045   if(rank != root) {
1046     // Send buffer to root
1047     smpi_mpi_send(sendbuf, count, datatype, root, system_tag, comm);
1048   } else {
1049     // FIXME: check for errors
1050     smpi_datatype_extent(datatype, &lb, &dataext);
1051     // Local copy from root
1052     if (sendbuf && recvbuf)
1053       smpi_datatype_copy(sendbuf, count, datatype, recvbuf, count, datatype);
1054     // Receive buffers from senders
1055     //TODO: make a MPI_barrier here ?
1056     requests = xbt_new(MPI_Request, size - 1);
1057     tmpbufs = xbt_new(void *, size - 1);
1058     index = 0;
1059     for(src = 0; src < size; src++) {
1060       if(src != root) {
1061         // FIXME: possibly overkill we we have contiguous/noncontiguous data
1062         //  mapping...
1063         tmpbufs[index] = xbt_malloc(count * dataext);
1064         requests[index] =
1065           smpi_irecv_init(tmpbufs[index], count, datatype, src,
1066                           system_tag, comm);
1067         index++;
1068       }
1069     }
1070     // Wait for completion of irecv's.
1071     smpi_mpi_startall(size - 1, requests);
1072     for(src = 0; src < size - 1; src++) {
1073       index = smpi_mpi_waitany(size - 1, requests, MPI_STATUS_IGNORE);
1074       XBT_DEBUG("finished waiting any request with index %d", index);
1075       if(index == MPI_UNDEFINED) {
1076         break;
1077       }
1078       if(op) /* op can be MPI_OP_NULL that does nothing */
1079         smpi_op_apply(op, tmpbufs[index], recvbuf, &count, &datatype);
1080     }
1081     for(index = 0; index < size - 1; index++) {
1082       xbt_free(tmpbufs[index]);
1083     }
1084     xbt_free(tmpbufs);
1085     xbt_free(requests);
1086   }
1087 }
1088
1089 void smpi_mpi_allreduce(void *sendbuf, void *recvbuf, int count,
1090                         MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1091 {
1092   smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1093   smpi_mpi_bcast(recvbuf, count, datatype, 0, comm);
1094 }
1095
1096 void smpi_mpi_scan(void *sendbuf, void *recvbuf, int count,
1097                    MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1098 {
1099   int system_tag = 666;
1100   int rank, size, other, index;
1101   MPI_Aint lb = 0, dataext = 0;
1102   MPI_Request *requests;
1103   void **tmpbufs;
1104
1105   rank = smpi_comm_rank(comm);
1106   size = smpi_comm_size(comm);
1107
1108   // FIXME: check for errors
1109   smpi_datatype_extent(datatype, &lb, &dataext);
1110
1111   // Local copy from self
1112   smpi_datatype_copy(sendbuf, count, datatype, recvbuf, count, datatype);
1113
1114   // Send/Recv buffers to/from others;
1115   requests = xbt_new(MPI_Request, size - 1);
1116   tmpbufs = xbt_new(void *, rank);
1117   index = 0;
1118   for(other = 0; other < rank; other++) {
1119     // FIXME: possibly overkill we we have contiguous/noncontiguous data
1120     // mapping...
1121     tmpbufs[index] = xbt_malloc(count * dataext);
1122     requests[index] =
1123       smpi_irecv_init(tmpbufs[index], count, datatype, other, system_tag,
1124                       comm);
1125     index++;
1126   }
1127   for(other = rank + 1; other < size; other++) {
1128     requests[index] =
1129       smpi_isend_init(sendbuf, count, datatype, other, system_tag, comm);
1130     index++;
1131   }
1132   // Wait for completion of all comms.
1133   smpi_mpi_startall(size - 1, requests);
1134   for(other = 0; other < size - 1; other++) {
1135     index = smpi_mpi_waitany(size - 1, requests, MPI_STATUS_IGNORE);
1136     if(index == MPI_UNDEFINED) {
1137       break;
1138     }
1139     if(index < rank) {
1140       // #Request is below rank: it's a irecv
1141       smpi_op_apply(op, tmpbufs[index], recvbuf, &count, &datatype);
1142     }
1143   }
1144   for(index = 0; index < rank; index++) {
1145     xbt_free(tmpbufs[index]);
1146   }
1147   xbt_free(tmpbufs);
1148   xbt_free(requests);
1149 }