Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
separate user-level and internal-level for SMPI_MPI_{Reduce,Bcast}
[simgrid.git] / src / smpi / smpi_mpi.c
1
2
3 #include "private.h"
4 #include "smpi_coll_private.h"
5
6 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_mpi, smpi,
7                                 "Logging specific to SMPI (mpi)");
8
9 int SMPI_MPI_Init(int *argc, char ***argv)
10 {
11   smpi_process_init(argc, argv);
12   smpi_bench_begin();
13   return MPI_SUCCESS;
14 }
15
16 int SMPI_MPI_Finalize()
17 {
18   smpi_bench_end();
19   smpi_process_finalize();
20   return MPI_SUCCESS;
21 }
22
23 // right now this just exits the current node, should send abort signal to all
24 // hosts in the communicator (TODO)
25 int SMPI_MPI_Abort(MPI_Comm comm, int errorcode)
26 {
27   smpi_exit(errorcode);
28   return 0;
29 }
30
31 int SMPI_MPI_Comm_size(MPI_Comm comm, int *size)
32 {
33   int retval = MPI_SUCCESS;
34
35   smpi_bench_end();
36
37   if (NULL == comm) {
38     retval = MPI_ERR_COMM;
39   } else if (NULL == size) {
40     retval = MPI_ERR_ARG;
41   } else {
42     *size = comm->size;
43   }
44
45   smpi_bench_begin();
46
47   return retval;
48 }
49
50 int SMPI_MPI_Comm_rank(MPI_Comm comm, int *rank)
51 {
52   int retval = MPI_SUCCESS;
53
54   smpi_bench_end();
55
56   if (NULL == comm) {
57     retval = MPI_ERR_COMM;
58   } else if (NULL == rank) {
59     retval = MPI_ERR_ARG;
60   } else {
61     *rank = smpi_mpi_comm_rank(comm);
62   }
63
64   smpi_bench_begin();
65
66   return retval;
67 }
68
69 int SMPI_MPI_Type_size(MPI_Datatype datatype, size_t * size)
70 {
71   int retval = MPI_SUCCESS;
72
73   smpi_bench_end();
74
75   if (NULL == datatype) {
76     retval = MPI_ERR_TYPE;
77   } else if (NULL == size) {
78     retval = MPI_ERR_ARG;
79   } else {
80     *size = datatype->size;
81   }
82
83   smpi_bench_begin();
84
85   return retval;
86 }
87
88 int SMPI_MPI_Barrier(MPI_Comm comm)
89 {
90   int retval = MPI_SUCCESS;
91   int arity=4;
92
93   smpi_bench_end();
94
95   if (NULL == comm) {
96     retval = MPI_ERR_COMM;
97   } else {
98
99     /*
100      * original implemantation:
101      * retval = smpi_mpi_barrier(comm);
102      * this one is unrealistic: it just cond_waits, means no time.
103      */
104      retval = nary_tree_barrier( comm, arity );
105   }
106
107   smpi_bench_begin();
108
109   return retval;
110 }
111
112 int SMPI_MPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
113                    int tag, MPI_Comm comm, MPI_Request * request)
114 {
115   int retval = MPI_SUCCESS;
116
117   smpi_bench_end();
118
119   retval = smpi_create_request(buf, count, datatype, src, 0, tag, comm,
120                                request);
121   if (NULL != *request && MPI_SUCCESS == retval) {
122     retval = smpi_mpi_irecv(*request);
123   }
124
125   smpi_bench_begin();
126
127   return retval;
128 }
129
130 int SMPI_MPI_Recv(void *buf, int count, MPI_Datatype datatype, int src,
131                   int tag, MPI_Comm comm, MPI_Status * status)
132 {
133   int retval = MPI_SUCCESS;
134   smpi_mpi_request_t request;
135
136   smpi_bench_end();
137
138   retval = smpi_create_request(buf, count, datatype, src, 0, tag, comm,
139                                &request);
140   if (NULL != request && MPI_SUCCESS == retval) {
141     retval = smpi_mpi_irecv(request);
142     if (MPI_SUCCESS == retval) {
143       retval = smpi_mpi_wait(request, status);
144     }
145     xbt_mallocator_release(smpi_global->request_mallocator, request);
146   }
147
148   smpi_bench_begin();
149
150   return retval;
151 }
152
153 int SMPI_MPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
154                    int tag, MPI_Comm comm, MPI_Request * request)
155 {
156   int retval = MPI_SUCCESS;
157
158   smpi_bench_end();
159
160   retval = smpi_create_request(buf, count, datatype, 0, dst, tag, comm,
161                                request);
162   if (NULL != *request && MPI_SUCCESS == retval) {
163     retval = smpi_mpi_isend(*request);
164   }
165
166   smpi_bench_begin();
167
168   return retval;
169 }
170
171 int SMPI_MPI_Send(void *buf, int count, MPI_Datatype datatype, int dst,
172                   int tag, MPI_Comm comm)
173 {
174   int retval = MPI_SUCCESS;
175   smpi_mpi_request_t request;
176
177   smpi_bench_end();
178
179   retval = smpi_create_request(buf, count, datatype, 0, dst, tag, comm,
180                                &request);
181   if (NULL != request && MPI_SUCCESS == retval) {
182     retval = smpi_mpi_isend(request);
183     if (MPI_SUCCESS == retval) {
184       smpi_mpi_wait(request, MPI_STATUS_IGNORE);
185     }
186     xbt_mallocator_release(smpi_global->request_mallocator, request);
187   }
188
189   smpi_bench_begin();
190
191   return retval;
192 }
193
194 /**
195  * MPI_Sendrecv
196  **/
197 int SMPI_MPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dest, int sendtag, 
198                     void *recvbuf, int recvcount, MPI_Datatype recvtype, int source, int recvtag,
199                     MPI_Comm comm, MPI_Status *status)
200 {
201 int rank;
202 int retval = MPI_SUCCESS;
203 smpi_mpi_request_t srequest;
204 smpi_mpi_request_t rrequest;
205
206           rank = smpi_mpi_comm_rank(comm);
207
208           /* send */
209           /* -------------*/
210             retval = smpi_create_request(sendbuf, sendcount, sendtype, 
211                                 rank,dest,sendtag, 
212                                 comm, &srequest);
213           printf("[%d] isend request src=%d -> dst=%d (retval=%d)\n",rank,rank,dest,retval);
214           smpi_mpi_isend(srequest);
215         
216           
217           //retval = MPI_Isend( sendbuf, sendcount, sendtype, dest, sendtag, MPI_COMM_WORLD, &srequest);
218
219
220           /* recv */
221           retval = smpi_create_request(recvbuf, recvcount, recvtype, 
222                                 source, rank,recvtag, 
223                                 comm, &rrequest);
224           printf("[%d] irecv request src=%d -> dst=%d (retval=%d)\n",rank,source,rank,retval);
225           smpi_mpi_irecv(rrequest);
226
227           //retval = MPI_Irecv( recvbuf, recvcount, recvtype, source, recvtag, MPI_COMM_WORLD, &rrequest);
228
229
230           smpi_mpi_wait(srequest, MPI_STATUS_IGNORE);
231           printf("[%d] isend request src=%d dst=%d tag=%d COMPLETED (retval=%d) \n",rank,rank,dest,sendtag,retval);
232
233           smpi_mpi_wait(rrequest, MPI_STATUS_IGNORE);
234           printf("[%d] irecv request src=%d -> dst=%d tag=%d COMPLETED (retval=%d)\n",rank,source,rank,recvtag,retval);
235
236           return(retval);
237 }
238
239
240 /**
241  * MPI_Wait and friends
242  **/
243 int SMPI_MPI_Wait(MPI_Request * request, MPI_Status * status)
244 {
245   return smpi_mpi_wait(*request, status);
246 }
247
248 int SMPI_MPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
249 {
250   return smpi_mpi_waitall(count, requests, status);
251 }
252
253 int SMPI_MPI_Waitany(int count, MPI_Request requests[], int *index,
254                      MPI_Status status[])
255 {
256   return smpi_mpi_waitany(count, requests, index, status);
257 }
258
259 /**
260  * MPI_Bcast
261  **/
262
263 /**
264  * flat bcast 
265  **/
266 int flat_tree_bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
267 int flat_tree_bcast(void *buf, int count, MPI_Datatype datatype, int root,
268                 MPI_Comm comm)
269 {
270         int rank;
271         int retval = MPI_SUCCESS;
272         smpi_mpi_request_t request;
273
274         rank = smpi_mpi_comm_rank(comm);
275         if (rank == root) {
276                 retval = smpi_create_request(buf, count, datatype, root,
277                                 (root + 1) % comm->size, 0, comm, &request);
278                 request->forward = comm->size - 1;
279                 smpi_mpi_isend(request);
280         } else {
281                 retval = smpi_create_request(buf, count, datatype, MPI_ANY_SOURCE, rank,
282                                 0, comm, &request);
283                 smpi_mpi_irecv(request);
284         }
285
286         smpi_mpi_wait(request, MPI_STATUS_IGNORE);
287         xbt_mallocator_release(smpi_global->request_mallocator, request);
288
289         return(retval);
290
291 }
292 /**
293  * Bcast internal level 
294  **/
295 int smpi_mpi_bcast(void *buf, int count, MPI_Datatype datatype, int root,
296                    MPI_Comm comm)
297 {
298   int retval = MPI_SUCCESS;
299   //retval = flat_tree_bcast(buf, count, datatype, root, comm);
300   retval = nary_tree_bcast(buf, count, datatype, root, comm, 2 );
301   return retval;
302 }
303
304 /**
305  * Bcast user entry point
306  **/
307 int SMPI_MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
308                    MPI_Comm comm)
309 {
310   int retval = MPI_SUCCESS;
311
312   smpi_bench_end();
313   smpi_mpi_bcast(buf,count,datatype,root,comm);
314   smpi_bench_begin();
315
316   return retval;
317 }
318
319
320
321 #ifdef DEBUG_REDUCE
322 /**
323  * debugging helper function
324  **/
325 static void print_buffer_int(void *buf, int len, char *msg, int rank)
326 {
327   int tmp, *v;
328   printf("**[%d] %s: ", rank, msg);
329   for (tmp = 0; tmp < len; tmp++) {
330     v = buf;
331     printf("[%d]", v[tmp]);
332   }
333   printf("\n");
334   free(msg);
335 }
336 static void print_buffer_double(void *buf, int len, char *msg, int rank)
337 {
338   int tmp;
339   double *v;
340   printf("**[%d] %s: ", rank, msg);
341   for (tmp = 0; tmp < len; tmp++) {
342     v = buf;
343     printf("[%lf]", v[tmp]);
344   }
345   printf("\n");
346   free(msg);
347 }
348
349
350 #endif
351 /**
352  * MPI_Reduce
353  **/
354 int smpi_mpi_reduce(void *sendbuf, void *recvbuf, int count,
355                 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
356 {
357         int retval = MPI_SUCCESS;
358         int rank;
359         int size;
360         int i;
361         int tag = 0;
362         smpi_mpi_request_t *requests;
363         smpi_mpi_request_t request;
364
365         smpi_bench_end();
366
367         rank = smpi_mpi_comm_rank(comm);
368         size = comm->size;
369
370         if (rank != root) {           // if i am not ROOT, simply send my buffer to root
371
372 #ifdef DEBUG_REDUCE
373                 print_buffer_int(sendbuf, count, xbt_strdup("sndbuf"), rank);
374 #endif
375                 retval = smpi_create_request(sendbuf, count, datatype, rank, root, tag, comm,
376                                         &request);
377                 smpi_mpi_isend(request);
378                 smpi_mpi_wait(request, MPI_STATUS_IGNORE);
379                 xbt_mallocator_release(smpi_global->request_mallocator, request);
380
381         } else {
382                 // i am the ROOT: wait for all buffers by creating one request by sender
383                 int src;
384                 requests = xbt_malloc((size-1) * sizeof(smpi_mpi_request_t));
385
386                 void **tmpbufs = xbt_malloc((size-1) * sizeof(void *));
387                 for (i = 0; i < size-1; i++) {
388                         // we need 1 buffer per request to store intermediate receptions
389                         tmpbufs[i] = xbt_malloc(count * datatype->size);
390                 }  
391                 // root: initiliaze recv buf with my own snd buf
392                 memcpy(recvbuf, sendbuf, count * datatype->size * sizeof(char));  
393
394                 // i can not use: 'request->forward = size-1;' (which would progagate size-1 receive reqs)
395                 // since we should op values as soon as one receiving request matches.
396                 for (i = 0; i < size-1; i++) {
397                         // reminder: for smpi_create_request() the src is always the process sending.
398                         src = i < root ? i : i + 1;
399                         retval = smpi_create_request(tmpbufs[i], count, datatype,
400                                         src, root, tag, comm, &(requests[i]));
401                         if (NULL != requests[i] && MPI_SUCCESS == retval) {
402                                 if (MPI_SUCCESS == retval) {
403                                         smpi_mpi_irecv(requests[i]);
404                                 }
405                         }
406                 }
407                 // now, wait for completion of all irecv's.
408                 for (i = 0; i < size-1; i++) {
409                         int index = MPI_UNDEFINED;
410                         smpi_mpi_waitany( size-1, requests, &index, MPI_STATUS_IGNORE);
411 #ifdef DEBUG_REDUCE
412                         printf ("MPI_Waitany() unblocked: root received (completes req[index=%d])\n",index);
413                         print_buffer_int(tmpbufs[index], count, bprintf("tmpbufs[index=%d] (value received)", index),
414                                         rank);
415 #endif
416
417                         // arg 2 is modified
418                         op->func(tmpbufs[index], recvbuf, &count, &datatype);
419 #ifdef DEBUG_REDUCE
420                         print_buffer_int(recvbuf, count, xbt_strdup("rcvbuf"), rank);
421 #endif
422                         xbt_free(tmpbufs[index]);
423                         /* FIXME: with the following line, it  generates an
424                          * [xbt_ex/CRITICAL] Conditional list not empty 162518800.
425                          */
426                         // xbt_mallocator_release(smpi_global->request_mallocator, requests[index]);
427                 }
428                 xbt_free(requests);
429                 xbt_free(tmpbufs);
430         }
431         return retval;
432 }
433
434 /**
435  * MPI_Reduce user entry point
436  **/
437 int SMPI_MPI_Reduce(void *sendbuf, void *recvbuf, int count,
438                     MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
439 {
440 int retval = MPI_SUCCESS;
441
442           smpi_bench_end();
443
444           retval = smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
445
446           smpi_bench_begin();
447           return retval;
448 }
449
450
451
452 /**
453  * MPI_Allreduce
454  *
455  * Same as MPI_REDUCE except that the result appears in the receive buffer of all the group members.
456  **/
457 int SMPI_MPI_Allreduce( void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
458                     MPI_Op op, MPI_Comm comm )
459 {
460 int retval = MPI_SUCCESS;
461 int root=0;  // arbitrary choice
462
463           smpi_bench_end();
464
465           retval = smpi_mpi_reduce( sendbuf, recvbuf, count, datatype, op, root, comm);
466           if (MPI_SUCCESS != retval)
467                     return(retval);
468
469           retval = smpi_mpi_bcast( sendbuf, count, datatype, root, comm);
470
471           smpi_bench_end();
472           return( retval );
473 }
474
475
476 /**
477  * MPI_Scatter user entry point
478  **/
479 //int SMPI_MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype datatype, 
480 //                       void *recvbuf, int recvcount, MPI_Datatype recvtype,int root,
481 //                      MPI_Comm comm);
482 int SMPI_MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype datatype, 
483                          void *recvbuf, int recvcount, MPI_Datatype recvtype,
484                            int root, MPI_Comm comm)
485 {
486   int retval = MPI_SUCCESS;
487   int i;
488   int cnt=0;  
489   int rank;
490   int tag=0;
491   char *cptr;  // to manipulate the void * buffers
492   smpi_mpi_request_t *requests;
493   smpi_mpi_request_t request;
494   smpi_mpi_status_t status;
495
496
497   smpi_bench_end();
498
499   rank = smpi_mpi_comm_rank(comm);
500
501   requests = xbt_malloc((comm->size-1) * sizeof(smpi_mpi_request_t));
502   if (rank == root) {
503           // i am the root: distribute my sendbuf
504           //print_buffer_int(sendbuf, comm->size, xbt_strdup("rcvbuf"), rank);
505           cptr = sendbuf;
506           for (i=0; i < comm->size; i++) {
507                   if ( i!=root ) { // send to processes ...
508
509                           retval = smpi_create_request((void *)cptr, sendcount, 
510                                           datatype, root, i, tag, comm, &(requests[cnt]));
511                           if (NULL != requests[cnt] && MPI_SUCCESS == retval) {
512                                   if (MPI_SUCCESS == retval) {
513                                           smpi_mpi_isend(requests[cnt]);
514                                   }
515                                   }
516                                   cnt++;
517                         } 
518                         else { // ... except if it's me.
519                                   memcpy(recvbuf, (void *)cptr, recvcount*recvtype->size*sizeof(char));
520                         }
521                   cptr += sendcount*datatype->size;
522             }
523             for(i=0; i<cnt; i++) { // wait for send to complete
524                             /* FIXME: waitall() should be slightly better */
525                             smpi_mpi_wait(requests[i], &status);
526                             xbt_mallocator_release(smpi_global->request_mallocator, requests[i]);
527
528             }
529   } 
530   else {  // i am a non-root process: wait data from the root
531             retval = smpi_create_request(recvbuf,recvcount, 
532                                   recvtype, root, rank, tag, comm, &request);
533             if (NULL != request && MPI_SUCCESS == retval) {
534                         if (MPI_SUCCESS == retval) {
535                                   smpi_mpi_irecv(request);
536                         }
537             }
538             smpi_mpi_wait(request, &status);
539             xbt_mallocator_release(smpi_global->request_mallocator, request);
540   }
541   xbt_free(requests);
542
543   smpi_bench_begin();
544
545   return retval;
546 }
547
548
549 /**
550  * MPI_Alltoall user entry point
551  * 
552  * Uses the logic of OpenMPI (upto 1.2.7 or greater) for the optimizations
553  * ompi/mca/coll/tuned/coll_tuned_module.c
554  **/
555 int SMPI_MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype datatype, 
556                          void *recvbuf, int recvcount, MPI_Datatype recvtype,
557                            MPI_Comm comm)
558 {
559   int retval = MPI_SUCCESS;
560   int block_dsize;
561   int rank;
562
563   smpi_bench_end();
564
565   rank = smpi_mpi_comm_rank(comm);
566   block_dsize = datatype->size * sendcount;
567
568   if ((block_dsize < 200) && (comm->size > 12)) {
569             retval = smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, datatype,
570                                   recvbuf, recvcount, recvtype, comm);
571
572   } else if (block_dsize < 3000) {
573 /* use this one !!          retval = smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount, datatype,
574                                   recvbuf, recvcount, recvtype, comm);
575                                   */
576   retval = smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, datatype,
577                                   recvbuf, recvcount, recvtype, comm);
578   } else {
579
580   retval = smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, datatype,
581                                   recvbuf, recvcount, recvtype, comm);
582   }
583
584   smpi_bench_begin();
585
586   return retval;
587 }
588
589
590
591
592 // used by comm_split to sort ranks based on key values
593 int smpi_compare_rankkeys(const void *a, const void *b);
594 int smpi_compare_rankkeys(const void *a, const void *b)
595 {
596   int *x = (int *) a;
597   int *y = (int *) b;
598
599   if (x[1] < y[1])
600     return -1;
601
602   if (x[1] == y[1]) {
603     if (x[0] < y[0])
604       return -1;
605     if (x[0] == y[0])
606       return 0;
607     return 1;
608   }
609
610   return 1;
611 }
612
613 int SMPI_MPI_Comm_split(MPI_Comm comm, int color, int key,
614                         MPI_Comm * comm_out)
615 {
616   int retval = MPI_SUCCESS;
617
618   int index, rank;
619   smpi_mpi_request_t request;
620   int colorkey[2];
621   smpi_mpi_status_t status;
622
623   smpi_bench_end();
624
625   // FIXME: need to test parameters
626
627   index = smpi_process_index();
628   rank = comm->index_to_rank_map[index];
629
630   // default output
631   comm_out = NULL;
632
633   // root node does most of the real work
634   if (0 == rank) {
635     int colormap[comm->size];
636     int keymap[comm->size];
637     int rankkeymap[comm->size * 2];
638     int i, j;
639     smpi_mpi_communicator_t tempcomm = NULL;
640     int count;
641     int indextmp;
642
643     colormap[0] = color;
644     keymap[0] = key;
645
646     // FIXME: use scatter/gather or similar instead of individual comms
647     for (i = 1; i < comm->size; i++) {
648       retval = smpi_create_request(colorkey, 2, MPI_INT, MPI_ANY_SOURCE,
649                                    rank, MPI_ANY_TAG, comm, &request);
650       smpi_mpi_irecv(request);
651       smpi_mpi_wait(request, &status);
652       colormap[status.MPI_SOURCE] = colorkey[0];
653       keymap[status.MPI_SOURCE] = colorkey[1];
654       xbt_mallocator_release(smpi_global->request_mallocator, request);
655     }
656
657     for (i = 0; i < comm->size; i++) {
658       if (MPI_UNDEFINED == colormap[i]) {
659         continue;
660       }
661       // make a list of nodes with current color and sort by keys
662       count = 0;
663       for (j = i; j < comm->size; j++) {
664         if (colormap[i] == colormap[j]) {
665           colormap[j] = MPI_UNDEFINED;
666           rankkeymap[count * 2] = j;
667           rankkeymap[count * 2 + 1] = keymap[j];
668           count++;
669         }
670       }
671       qsort(rankkeymap, count, sizeof(int) * 2, &smpi_compare_rankkeys);
672
673       // new communicator
674       tempcomm = xbt_new(s_smpi_mpi_communicator_t, 1);
675       tempcomm->barrier_count = 0;
676       tempcomm->size = count;
677       tempcomm->barrier_mutex = SIMIX_mutex_init();
678       tempcomm->barrier_cond = SIMIX_cond_init();
679       tempcomm->rank_to_index_map = xbt_new(int, count);
680       tempcomm->index_to_rank_map = xbt_new(int, smpi_global->process_count);
681       for (j = 0; j < smpi_global->process_count; j++) {
682         tempcomm->index_to_rank_map[j] = -1;
683       }
684       for (j = 0; j < count; j++) {
685         indextmp = comm->rank_to_index_map[rankkeymap[j * 2]];
686         tempcomm->rank_to_index_map[j] = indextmp;
687         tempcomm->index_to_rank_map[indextmp] = j;
688       }
689       for (j = 0; j < count; j++) {
690         if (rankkeymap[j * 2]) {
691           retval = smpi_create_request(&j, 1, MPI_INT, 0,
692                                        rankkeymap[j * 2], 0, comm, &request);
693           request->data = tempcomm;
694           smpi_mpi_isend(request);
695           smpi_mpi_wait(request, &status);
696           xbt_mallocator_release(smpi_global->request_mallocator, request);
697         } else {
698           *comm_out = tempcomm;
699         }
700       }
701     }
702   } else {
703     colorkey[0] = color;
704     colorkey[1] = key;
705     retval = smpi_create_request(colorkey, 2, MPI_INT, rank, 0, 0, comm,
706                                  &request);
707     smpi_mpi_isend(request);
708     smpi_mpi_wait(request, &status);
709     xbt_mallocator_release(smpi_global->request_mallocator, request);
710     if (MPI_UNDEFINED != color) {
711       retval = smpi_create_request(colorkey, 1, MPI_INT, 0, rank, 0, comm,
712                                    &request);
713       smpi_mpi_irecv(request);
714       smpi_mpi_wait(request, &status);
715       *comm_out = request->data;
716     }
717   }
718
719   smpi_bench_begin();
720
721   return retval;
722 }
723
724 double SMPI_MPI_Wtime(void)
725 {
726   return (SIMIX_get_clock());
727 }
728
729