Logo AND Algorithmique Numérique Distribuée

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