Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
c1f52536d6244626fda4ade49d1675c7e4ea801a
[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 /**
196  * MPI_Sendrecv internal level 
197  **/
198 int smpi_mpi_sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dest, int sendtag, 
199                     void *recvbuf, int recvcount, MPI_Datatype recvtype, int source, int recvtag,
200                     MPI_Comm comm, MPI_Status *status)
201 {
202 int rank;
203 int retval = MPI_SUCCESS;
204 smpi_mpi_request_t srequest;
205 smpi_mpi_request_t rrequest;
206
207           rank = smpi_mpi_comm_rank(comm);
208
209           /* send */
210           retval = smpi_create_request(sendbuf, sendcount, sendtype, 
211                                 rank,dest,sendtag, 
212                                 comm, &srequest);
213           smpi_mpi_isend(srequest);
214         
215           /* recv */
216           retval = smpi_create_request(recvbuf, recvcount, recvtype, 
217                                 source, rank,recvtag, 
218                                 comm, &rrequest);
219           smpi_mpi_irecv(rrequest);
220
221           smpi_mpi_wait(srequest, MPI_STATUS_IGNORE);
222           //printf("[%d] isend request src=%d dst=%d tag=%d COMPLETED (retval=%d) \n",rank,rank,dest,sendtag,retval);
223           smpi_mpi_wait(rrequest, MPI_STATUS_IGNORE);
224           //printf("[%d] irecv request src=%d -> dst=%d tag=%d COMPLETED (retval=%d)\n",rank,source,rank,recvtag,retval);
225
226           return(retval);
227 }
228 /**
229  * MPI_Sendrecv user entry point
230  **/
231 int SMPI_MPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dest, int sendtag, 
232                     void *recvbuf, int recvcount, MPI_Datatype recvtype, int source, int recvtag,
233                     MPI_Comm comm, MPI_Status *status)
234 {
235 int retval = MPI_SUCCESS;
236
237   smpi_bench_end();
238   smpi_mpi_sendrecv( sendbuf, sendcount, sendtype, dest, sendtag, 
239                          recvbuf, recvcount, recvtype, source, recvtag,
240                          comm, status);
241   smpi_bench_begin();
242
243   return retval;
244
245         
246 }
247
248 /**
249  * MPI_Wait and friends
250  **/
251 int SMPI_MPI_Wait(MPI_Request * request, MPI_Status * status)
252 {
253   return smpi_mpi_wait(*request, status);
254 }
255
256 int SMPI_MPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
257 {
258   return smpi_mpi_waitall(count, requests, status);
259 }
260
261 int SMPI_MPI_Waitany(int count, MPI_Request requests[], int *index,
262                      MPI_Status status[])
263 {
264   return smpi_mpi_waitany(count, requests, index, status);
265 }
266
267 /**
268  * MPI_Bcast
269  **/
270
271 /**
272  * flat bcast 
273  **/
274 int flat_tree_bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm);
275 int flat_tree_bcast(void *buf, int count, MPI_Datatype datatype, int root,
276                 MPI_Comm comm)
277 {
278         int rank;
279         int retval = MPI_SUCCESS;
280         smpi_mpi_request_t request;
281
282         rank = smpi_mpi_comm_rank(comm);
283         if (rank == root) {
284                 retval = smpi_create_request(buf, count, datatype, root,
285                                 (root + 1) % comm->size, 0, comm, &request);
286                 request->forward = comm->size - 1;
287                 smpi_mpi_isend(request);
288         } else {
289                 retval = smpi_create_request(buf, count, datatype, MPI_ANY_SOURCE, rank,
290                                 0, comm, &request);
291                 smpi_mpi_irecv(request);
292         }
293
294         smpi_mpi_wait(request, MPI_STATUS_IGNORE);
295         xbt_mallocator_release(smpi_global->request_mallocator, request);
296
297         return(retval);
298
299 }
300 /**
301  * Bcast internal level 
302  **/
303 int smpi_mpi_bcast(void *buf, int count, MPI_Datatype datatype, int root,
304                    MPI_Comm comm)
305 {
306   int retval = MPI_SUCCESS;
307   //retval = flat_tree_bcast(buf, count, datatype, root, comm);
308   retval = nary_tree_bcast(buf, count, datatype, root, comm, 2 );
309   return retval;
310 }
311
312 /**
313  * Bcast user entry point
314  **/
315 int SMPI_MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
316                    MPI_Comm comm)
317 {
318   int retval = MPI_SUCCESS;
319
320   smpi_bench_end();
321   smpi_mpi_bcast(buf,count,datatype,root,comm);
322   smpi_bench_begin();
323
324   return retval;
325 }
326
327
328
329 #ifdef DEBUG_REDUCE
330 /**
331  * debugging helper function
332  **/
333 static void print_buffer_int(void *buf, int len, char *msg, int rank)
334 {
335   int tmp, *v;
336   printf("**[%d] %s: ", rank, msg);
337   for (tmp = 0; tmp < len; tmp++) {
338     v = buf;
339     printf("[%d]", v[tmp]);
340   }
341   printf("\n");
342   free(msg);
343 }
344 static void print_buffer_double(void *buf, int len, char *msg, int rank)
345 {
346   int tmp;
347   double *v;
348   printf("**[%d] %s: ", rank, msg);
349   for (tmp = 0; tmp < len; tmp++) {
350     v = buf;
351     printf("[%lf]", v[tmp]);
352   }
353   printf("\n");
354   free(msg);
355 }
356
357
358 #endif
359 /**
360  * MPI_Reduce internal level 
361  **/
362 int smpi_mpi_reduce(void *sendbuf, void *recvbuf, int count,
363                 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
364 {
365         int retval = MPI_SUCCESS;
366         int rank;
367         int size;
368         int i;
369         int tag = 0;
370         smpi_mpi_request_t *requests;
371         smpi_mpi_request_t request;
372
373         smpi_bench_end();
374
375         rank = smpi_mpi_comm_rank(comm);
376         size = comm->size;
377
378         if (rank != root) {           // if i am not ROOT, simply send my buffer to root
379
380 #ifdef DEBUG_REDUCE
381                 print_buffer_int(sendbuf, count, xbt_strdup("sndbuf"), rank);
382 #endif
383                 retval = smpi_create_request(sendbuf, count, datatype, rank, root, tag, comm,
384                                         &request);
385                 smpi_mpi_isend(request);
386                 smpi_mpi_wait(request, MPI_STATUS_IGNORE);
387                 xbt_mallocator_release(smpi_global->request_mallocator, request);
388
389         } else {
390                 // i am the ROOT: wait for all buffers by creating one request by sender
391                 int src;
392                 requests = xbt_malloc((size-1) * sizeof(smpi_mpi_request_t));
393
394                 void **tmpbufs = xbt_malloc((size-1) * sizeof(void *));
395                 for (i = 0; i < size-1; i++) {
396                         // we need 1 buffer per request to store intermediate receptions
397                         tmpbufs[i] = xbt_malloc(count * datatype->size);
398                 }  
399                 // root: initiliaze recv buf with my own snd buf
400                 memcpy(recvbuf, sendbuf, count * datatype->size * sizeof(char));  
401
402                 // i can not use: 'request->forward = size-1;' (which would progagate size-1 receive reqs)
403                 // since we should op values as soon as one receiving request matches.
404                 for (i = 0; i < size-1; i++) {
405                         // reminder: for smpi_create_request() the src is always the process sending.
406                         src = i < root ? i : i + 1;
407                         retval = smpi_create_request(tmpbufs[i], count, datatype,
408                                         src, root, tag, comm, &(requests[i]));
409                         if (NULL != requests[i] && MPI_SUCCESS == retval) {
410                                 if (MPI_SUCCESS == retval) {
411                                         smpi_mpi_irecv(requests[i]);
412                                 }
413                         }
414                 }
415                 // now, wait for completion of all irecv's.
416                 for (i = 0; i < size-1; i++) {
417                         int index = MPI_UNDEFINED;
418                         smpi_mpi_waitany( size-1, requests, &index, MPI_STATUS_IGNORE);
419 #ifdef DEBUG_REDUCE
420                         printf ("MPI_Waitany() unblocked: root received (completes req[index=%d])\n",index);
421                         print_buffer_int(tmpbufs[index], count, bprintf("tmpbufs[index=%d] (value received)", index),
422                                         rank);
423 #endif
424
425                         // arg 2 is modified
426                         op->func(tmpbufs[index], recvbuf, &count, &datatype);
427 #ifdef DEBUG_REDUCE
428                         print_buffer_int(recvbuf, count, xbt_strdup("rcvbuf"), rank);
429 #endif
430                         xbt_free(tmpbufs[index]);
431                         /* FIXME: with the following line, it  generates an
432                          * [xbt_ex/CRITICAL] Conditional list not empty 162518800.
433                          */
434                         // xbt_mallocator_release(smpi_global->request_mallocator, requests[index]);
435                 }
436                 xbt_free(requests);
437                 xbt_free(tmpbufs);
438         }
439         return retval;
440 }
441
442 /**
443  * MPI_Reduce user entry point
444  **/
445 int SMPI_MPI_Reduce(void *sendbuf, void *recvbuf, int count,
446                     MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
447 {
448 int retval = MPI_SUCCESS;
449
450           smpi_bench_end();
451
452           retval = smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
453
454           smpi_bench_begin();
455           return retval;
456 }
457
458
459
460 /**
461  * MPI_Allreduce
462  *
463  * Same as MPI_REDUCE except that the result appears in the receive buffer of all the group members.
464  **/
465 int SMPI_MPI_Allreduce( void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
466                     MPI_Op op, MPI_Comm comm )
467 {
468 int retval = MPI_SUCCESS;
469 int root=0;  // arbitrary choice
470
471           smpi_bench_end();
472
473           retval = smpi_mpi_reduce( sendbuf, recvbuf, count, datatype, op, root, comm);
474           if (MPI_SUCCESS != retval)
475                     return(retval);
476
477           retval = smpi_mpi_bcast( sendbuf, count, datatype, root, comm);
478
479           smpi_bench_end();
480           return( retval );
481 }
482
483
484 /**
485  * MPI_Scatter user entry point
486  **/
487 int SMPI_MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype datatype, 
488                          void *recvbuf, int recvcount, MPI_Datatype recvtype,
489                            int root, MPI_Comm comm)
490 {
491   int retval = MPI_SUCCESS;
492   int i;
493   int cnt=0;  
494   int rank;
495   int tag=0;
496   char *cptr;  // to manipulate the void * buffers
497   smpi_mpi_request_t *requests;
498   smpi_mpi_request_t request;
499   smpi_mpi_status_t status;
500
501
502   smpi_bench_end();
503
504   rank = smpi_mpi_comm_rank(comm);
505
506   requests = xbt_malloc((comm->size-1) * sizeof(smpi_mpi_request_t));
507   if (rank == root) {
508           // i am the root: distribute my sendbuf
509           //print_buffer_int(sendbuf, comm->size, xbt_strdup("rcvbuf"), rank);
510           cptr = sendbuf;
511           for (i=0; i < comm->size; i++) {
512                   if ( i!=root ) { // send to processes ...
513
514                           retval = smpi_create_request((void *)cptr, sendcount, 
515                                           datatype, root, i, tag, comm, &(requests[cnt]));
516                           if (NULL != requests[cnt] && MPI_SUCCESS == retval) {
517                                   if (MPI_SUCCESS == retval) {
518                                           smpi_mpi_isend(requests[cnt]);
519                                   }
520                                   }
521                                   cnt++;
522                         } 
523                         else { // ... except if it's me.
524                                   memcpy(recvbuf, (void *)cptr, recvcount*recvtype->size*sizeof(char));
525                         }
526                   cptr += sendcount*datatype->size;
527             }
528             for(i=0; i<cnt; i++) { // wait for send to complete
529                             /* FIXME: waitall() should be slightly better */
530                             smpi_mpi_wait(requests[i], &status);
531                             xbt_mallocator_release(smpi_global->request_mallocator, requests[i]);
532
533             }
534   } 
535   else {  // i am a non-root process: wait data from the root
536             retval = smpi_create_request(recvbuf,recvcount, 
537                                   recvtype, root, rank, tag, comm, &request);
538             if (NULL != request && MPI_SUCCESS == retval) {
539                         if (MPI_SUCCESS == retval) {
540                                   smpi_mpi_irecv(request);
541                         }
542             }
543             smpi_mpi_wait(request, &status);
544             xbt_mallocator_release(smpi_global->request_mallocator, request);
545   }
546   xbt_free(requests);
547
548   smpi_bench_begin();
549
550   return retval;
551 }
552
553
554 /**
555  * MPI_Alltoall user entry point
556  * 
557  * Uses the logic of OpenMPI (upto 1.2.7 or greater) for the optimizations
558  * ompi/mca/coll/tuned/coll_tuned_module.c
559  **/
560 int SMPI_MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype datatype, 
561                          void *recvbuf, int recvcount, MPI_Datatype recvtype,
562                            MPI_Comm comm)
563 {
564   int retval = MPI_SUCCESS;
565   int block_dsize;
566   int rank;
567
568   smpi_bench_end();
569
570   rank = smpi_mpi_comm_rank(comm);
571   block_dsize = datatype->size * sendcount;
572
573   if ((block_dsize < 200) && (comm->size > 12)) {
574             retval = smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, datatype,
575                                   recvbuf, recvcount, recvtype, comm);
576
577   } else if (block_dsize < 3000) {
578 /* use this one !!          retval = smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount, datatype,
579                                   recvbuf, recvcount, recvtype, comm);
580                                   */
581   retval = smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, datatype,
582                                   recvbuf, recvcount, recvtype, comm);
583   } else {
584
585   retval = smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, datatype,
586                                   recvbuf, recvcount, recvtype, comm);
587   }
588
589   smpi_bench_begin();
590
591   return retval;
592 }
593
594
595
596
597 // used by comm_split to sort ranks based on key values
598 int smpi_compare_rankkeys(const void *a, const void *b);
599 int smpi_compare_rankkeys(const void *a, const void *b)
600 {
601   int *x = (int *) a;
602   int *y = (int *) b;
603
604   if (x[1] < y[1])
605     return -1;
606
607   if (x[1] == y[1]) {
608     if (x[0] < y[0])
609       return -1;
610     if (x[0] == y[0])
611       return 0;
612     return 1;
613   }
614
615   return 1;
616 }
617
618 int SMPI_MPI_Comm_split(MPI_Comm comm, int color, int key,
619                         MPI_Comm * comm_out)
620 {
621   int retval = MPI_SUCCESS;
622
623   int index, rank;
624   smpi_mpi_request_t request;
625   int colorkey[2];
626   smpi_mpi_status_t status;
627
628   smpi_bench_end();
629
630   // FIXME: need to test parameters
631
632   index = smpi_process_index();
633   rank = comm->index_to_rank_map[index];
634
635   // default output
636   comm_out = NULL;
637
638   // root node does most of the real work
639   if (0 == rank) {
640     int colormap[comm->size];
641     int keymap[comm->size];
642     int rankkeymap[comm->size * 2];
643     int i, j;
644     smpi_mpi_communicator_t tempcomm = NULL;
645     int count;
646     int indextmp;
647
648     colormap[0] = color;
649     keymap[0] = key;
650
651     // FIXME: use scatter/gather or similar instead of individual comms
652     for (i = 1; i < comm->size; i++) {
653       retval = smpi_create_request(colorkey, 2, MPI_INT, MPI_ANY_SOURCE,
654                                    rank, MPI_ANY_TAG, comm, &request);
655       smpi_mpi_irecv(request);
656       smpi_mpi_wait(request, &status);
657       colormap[status.MPI_SOURCE] = colorkey[0];
658       keymap[status.MPI_SOURCE] = colorkey[1];
659       xbt_mallocator_release(smpi_global->request_mallocator, request);
660     }
661
662     for (i = 0; i < comm->size; i++) {
663       if (MPI_UNDEFINED == colormap[i]) {
664         continue;
665       }
666       // make a list of nodes with current color and sort by keys
667       count = 0;
668       for (j = i; j < comm->size; j++) {
669         if (colormap[i] == colormap[j]) {
670           colormap[j] = MPI_UNDEFINED;
671           rankkeymap[count * 2] = j;
672           rankkeymap[count * 2 + 1] = keymap[j];
673           count++;
674         }
675       }
676       qsort(rankkeymap, count, sizeof(int) * 2, &smpi_compare_rankkeys);
677
678       // new communicator
679       tempcomm = xbt_new(s_smpi_mpi_communicator_t, 1);
680       tempcomm->barrier_count = 0;
681       tempcomm->size = count;
682       tempcomm->barrier_mutex = SIMIX_mutex_init();
683       tempcomm->barrier_cond = SIMIX_cond_init();
684       tempcomm->rank_to_index_map = xbt_new(int, count);
685       tempcomm->index_to_rank_map = xbt_new(int, smpi_global->process_count);
686       for (j = 0; j < smpi_global->process_count; j++) {
687         tempcomm->index_to_rank_map[j] = -1;
688       }
689       for (j = 0; j < count; j++) {
690         indextmp = comm->rank_to_index_map[rankkeymap[j * 2]];
691         tempcomm->rank_to_index_map[j] = indextmp;
692         tempcomm->index_to_rank_map[indextmp] = j;
693       }
694       for (j = 0; j < count; j++) {
695         if (rankkeymap[j * 2]) {
696           retval = smpi_create_request(&j, 1, MPI_INT, 0,
697                                        rankkeymap[j * 2], 0, comm, &request);
698           request->data = tempcomm;
699           smpi_mpi_isend(request);
700           smpi_mpi_wait(request, &status);
701           xbt_mallocator_release(smpi_global->request_mallocator, request);
702         } else {
703           *comm_out = tempcomm;
704         }
705       }
706     }
707   } else {
708     colorkey[0] = color;
709     colorkey[1] = key;
710     retval = smpi_create_request(colorkey, 2, MPI_INT, rank, 0, 0, comm,
711                                  &request);
712     smpi_mpi_isend(request);
713     smpi_mpi_wait(request, &status);
714     xbt_mallocator_release(smpi_global->request_mallocator, request);
715     if (MPI_UNDEFINED != color) {
716       retval = smpi_create_request(colorkey, 1, MPI_INT, 0, rank, 0, comm,
717                                    &request);
718       smpi_mpi_irecv(request);
719       smpi_mpi_wait(request, &status);
720       *comm_out = request->data;
721     }
722   }
723
724   smpi_bench_begin();
725
726   return retval;
727 }
728
729 double SMPI_MPI_Wtime(void)
730 {
731   return (SIMIX_get_clock());
732 }
733
734