Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
alltoall implemented (almost opmpi algorithms)
[simgrid.git] / src / smpi / smpi_mpi.c
1
2
3 #include "private.h"
4 #include "smpi_coll_private.h"
5 #include "smpi_mpi_dt_private.h"
6
7 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_mpi, smpi,
8                                 "Logging specific to SMPI (mpi)");
9
10 int SMPI_MPI_Init(int *argc, char ***argv)
11 {
12   smpi_process_init(argc, argv);
13   smpi_bench_begin();
14   return MPI_SUCCESS;
15 }
16
17 int SMPI_MPI_Finalize()
18 {
19   smpi_bench_end();
20   smpi_process_finalize();
21   return MPI_SUCCESS;
22 }
23
24 // right now this just exits the current node, should send abort signal to all
25 // hosts in the communicator (TODO)
26 int SMPI_MPI_Abort(MPI_Comm comm, int errorcode)
27 {
28   smpi_exit(errorcode);
29   return 0;
30 }
31
32 int SMPI_MPI_Comm_size(MPI_Comm comm, int *size)
33 {
34   int retval = MPI_SUCCESS;
35
36   smpi_bench_end();
37
38   if (NULL == comm) {
39     retval = MPI_ERR_COMM;
40   } else if (NULL == size) {
41     retval = MPI_ERR_ARG;
42   } else {
43     *size = comm->size;
44   }
45
46   smpi_bench_begin();
47
48   return retval;
49 }
50
51 int SMPI_MPI_Comm_rank(MPI_Comm comm, int *rank)
52 {
53   int retval = MPI_SUCCESS;
54
55   smpi_bench_end();
56
57   if (NULL == comm) {
58     retval = MPI_ERR_COMM;
59   } else if (NULL == rank) {
60     retval = MPI_ERR_ARG;
61   } else {
62     *rank = smpi_mpi_comm_rank(comm);
63   }
64
65   smpi_bench_begin();
66
67   return retval;
68 }
69
70
71
72 /**
73  * Barrier
74  **/
75 int SMPI_MPI_Barrier(MPI_Comm comm)
76 {
77   int retval = MPI_SUCCESS;
78   int arity=4;
79
80   smpi_bench_end();
81
82   if (NULL == comm) {
83     retval = MPI_ERR_COMM;
84   } else {
85
86     /*
87      * original implemantation:
88      * retval = smpi_mpi_barrier(comm);
89      * this one is unrealistic: it just cond_waits, means no time.
90      */
91      retval = nary_tree_barrier( comm, arity );
92   }
93
94   smpi_bench_begin();
95
96   return retval;
97 }
98
99
100
101 int SMPI_MPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
102                    int tag, MPI_Comm comm, MPI_Request * request)
103 {
104   int retval = MPI_SUCCESS;
105
106   smpi_bench_end();
107
108   retval = smpi_create_request(buf, count, datatype, src, 0, tag, comm,
109                                request);
110   if (NULL != *request && MPI_SUCCESS == retval) {
111     retval = smpi_mpi_irecv(*request);
112   }
113
114   smpi_bench_begin();
115
116   return retval;
117 }
118
119 int SMPI_MPI_Recv(void *buf, int count, MPI_Datatype datatype, int src,
120                   int tag, MPI_Comm comm, MPI_Status * status)
121 {
122   int retval = MPI_SUCCESS;
123   smpi_mpi_request_t request;
124
125   smpi_bench_end();
126
127   retval = smpi_create_request(buf, count, datatype, src, 0, tag, comm,
128                                &request);
129   if (NULL != request && MPI_SUCCESS == retval) {
130     retval = smpi_mpi_irecv(request);
131     if (MPI_SUCCESS == retval) {
132       retval = smpi_mpi_wait(request, status);
133     }
134     xbt_mallocator_release(smpi_global->request_mallocator, request);
135   }
136
137   smpi_bench_begin();
138
139   return retval;
140 }
141
142 int SMPI_MPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
143                    int tag, MPI_Comm comm, MPI_Request * request)
144 {
145   int retval = MPI_SUCCESS;
146
147   smpi_bench_end();
148
149   retval = smpi_create_request(buf, count, datatype, 0, dst, tag, comm,
150                                request);
151   if (NULL != *request && MPI_SUCCESS == retval) {
152     retval = smpi_mpi_isend(*request);
153   }
154
155   smpi_bench_begin();
156
157   return retval;
158 }
159
160 /**
161  * MPI_Send user level
162  **/
163 int SMPI_MPI_Send(void *buf, int count, MPI_Datatype datatype, int dst,
164                   int tag, MPI_Comm comm)
165 {
166   int retval = MPI_SUCCESS;
167   smpi_mpi_request_t request;
168
169   smpi_bench_end();
170
171   retval = smpi_create_request(buf, count, datatype, 0, dst, tag, comm,
172                                &request);
173   if (NULL != request && MPI_SUCCESS == retval) {
174     retval = smpi_mpi_isend(request);
175     if (MPI_SUCCESS == retval) {
176       smpi_mpi_wait(request, MPI_STATUS_IGNORE);
177     }
178     xbt_mallocator_release(smpi_global->request_mallocator, request);
179   }
180
181   smpi_bench_begin();
182
183   return retval;
184 }
185
186
187 /**
188  * MPI_Sendrecv internal level 
189  **/
190 int smpi_mpi_sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dest, int sendtag, 
191                     void *recvbuf, int recvcount, MPI_Datatype recvtype, int source, int recvtag,
192                     MPI_Comm comm, MPI_Status *status)
193 {
194 int rank;
195 int retval = MPI_SUCCESS;
196 smpi_mpi_request_t srequest;
197 smpi_mpi_request_t rrequest;
198
199           rank = smpi_mpi_comm_rank(comm);
200
201           /* send */
202           retval = smpi_create_request(sendbuf, sendcount, sendtype, 
203                                 rank,dest,sendtag, 
204                                 comm, &srequest);
205           smpi_mpi_isend(srequest);
206         
207           /* recv */
208           retval = smpi_create_request(recvbuf, recvcount, recvtype, 
209                                 source, rank,recvtag, 
210                                 comm, &rrequest);
211           smpi_mpi_irecv(rrequest);
212
213           smpi_mpi_wait(srequest, MPI_STATUS_IGNORE);
214           //printf("[%d] isend request src=%d dst=%d tag=%d COMPLETED (retval=%d) \n",rank,rank,dest,sendtag,retval);
215           smpi_mpi_wait(rrequest, MPI_STATUS_IGNORE);
216           //printf("[%d] irecv request src=%d -> dst=%d tag=%d COMPLETED (retval=%d)\n",rank,source,rank,recvtag,retval);
217
218           return(retval);
219 }
220 /**
221  * MPI_Sendrecv user entry point
222  **/
223 int SMPI_MPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dest, int sendtag, 
224                     void *recvbuf, int recvcount, MPI_Datatype recvtype, int source, int recvtag,
225                     MPI_Comm comm, MPI_Status *status)
226 {
227 int retval = MPI_SUCCESS;
228
229   smpi_bench_end();
230   smpi_mpi_sendrecv( sendbuf, sendcount, sendtype, dest, sendtag, 
231                          recvbuf, recvcount, recvtype, source, recvtag,
232                          comm, status);
233   smpi_bench_begin();
234
235   return retval;
236
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   int rank = smpi_mpi_comm_rank(comm);
300
301   DEBUG1("<%d> entered smpi_mpi_bcast(). Calls nary_tree_bcast()",rank);
302   //retval = flat_tree_bcast(buf, count, datatype, root, comm);
303   retval = nary_tree_bcast(buf, count, datatype, root, comm, 2 );
304   return retval;
305 }
306
307 /**
308  * Bcast user entry point
309  **/
310 int SMPI_MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
311                    MPI_Comm comm)
312 {
313   int retval = MPI_SUCCESS;
314
315   smpi_bench_end();
316   smpi_mpi_bcast(buf,count,datatype,root,comm);
317   smpi_bench_begin();
318
319   return retval;
320 }
321
322
323
324 #ifdef DEBUG_REDUCE
325 /**
326  * debugging helper function
327  **/
328 static void print_buffer_int(void *buf, int len, char *msg, int rank)
329 {
330   int tmp, *v;
331   printf("**[%d] %s: ", rank, msg);
332   for (tmp = 0; tmp < len; tmp++) {
333     v = buf;
334     printf("[%d]", v[tmp]);
335   }
336   printf("\n");
337   free(msg);
338 }
339 static void print_buffer_double(void *buf, int len, char *msg, int rank)
340 {
341   int tmp;
342   double *v;
343   printf("**[%d] %s: ", rank, msg);
344   for (tmp = 0; tmp < len; tmp++) {
345     v = buf;
346     printf("[%lf]", v[tmp]);
347   }
348   printf("\n");
349   free(msg);
350 }
351
352
353 #endif
354 /**
355  * MPI_Reduce internal level 
356  **/
357 int smpi_mpi_reduce(void *sendbuf, void *recvbuf, int count,
358                 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
359 {
360         int retval = MPI_SUCCESS;
361         int rank;
362         int size;
363         int i;
364         int system_tag = 666;
365         smpi_mpi_request_t *requests;
366         smpi_mpi_request_t request;
367
368         smpi_bench_end();
369
370         rank = smpi_mpi_comm_rank(comm);
371         size = comm->size;
372         DEBUG1("<%d> entered smpi_mpi_reduce()",rank);
373
374         if (rank != root) {           // if i am not ROOT, simply send my buffer to root
375
376 #ifdef DEBUG_REDUCE
377                 print_buffer_int(sendbuf, count, xbt_strdup("sndbuf"), rank);
378 #endif
379                 retval = smpi_create_request(sendbuf, count, datatype, rank, root, system_tag, comm,
380                                         &request);
381                 smpi_mpi_isend(request);
382                 smpi_mpi_wait(request, MPI_STATUS_IGNORE);
383                 xbt_mallocator_release(smpi_global->request_mallocator, request);
384
385         } else {
386                 // i am the ROOT: wait for all buffers by creating one request by sender
387                 int src;
388                 requests = xbt_malloc((size-1) * sizeof(smpi_mpi_request_t));
389
390                 void **tmpbufs = xbt_malloc((size-1) * sizeof(void *));
391                 for (i = 0; i < size-1; i++) {
392                         // we need 1 buffer per request to store intermediate receptions
393                         tmpbufs[i] = xbt_malloc(count * datatype->size);
394                 }  
395                 // root: initiliaze recv buf with my own snd buf
396                 memcpy(recvbuf, sendbuf, count * datatype->size * sizeof(char));  
397
398                 // i can not use: 'request->forward = size-1;' (which would progagate size-1 receive reqs)
399                 // since we should op values as soon as one receiving request matches.
400                 for (i = 0; i < size-1; i++) {
401                         // reminder: for smpi_create_request() the src is always the process sending.
402                         src = i < root ? i : i + 1;
403                         retval = smpi_create_request(tmpbufs[i], count, datatype,
404                                         src, root, system_tag, comm, &(requests[i]));
405                         if (NULL != requests[i] && MPI_SUCCESS == retval) {
406                                 if (MPI_SUCCESS == retval) {
407                                         smpi_mpi_irecv(requests[i]);
408                                 }
409                         }
410                 }
411                 // now, wait for completion of all irecv's.
412                 for (i = 0; i < size-1; i++) {
413                         int index = MPI_UNDEFINED;
414                         smpi_mpi_waitany( size-1, requests, &index, MPI_STATUS_IGNORE);
415                         DEBUG3("<%d> waitany() unblocked by reception (completes request[%d]) (%d reqs remaining)",
416                                         rank,index,size-i-2);
417 #ifdef DEBUG_REDUCE
418                         print_buffer_int(tmpbufs[index], count, bprintf("tmpbufs[index=%d] (value received)", index),
419                                         rank);
420 #endif
421
422                         // arg 2 is modified
423                         op->func(tmpbufs[index], recvbuf, &count, &datatype);
424 #ifdef DEBUG_REDUCE
425                         print_buffer_int(recvbuf, count, xbt_strdup("rcvbuf"), rank);
426 #endif
427                         xbt_free(tmpbufs[index]);
428                         /* FIXME: with the following line, it  generates an
429                          * [xbt_ex/CRITICAL] Conditional list not empty 162518800.
430                          */
431                         // xbt_mallocator_release(smpi_global->request_mallocator, requests[index]);
432                 }
433                 xbt_free(requests);
434                 xbt_free(tmpbufs);
435         }
436         return retval;
437 }
438
439 /**
440  * MPI_Reduce user entry point
441  **/
442 int SMPI_MPI_Reduce(void *sendbuf, void *recvbuf, int count,
443                     MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
444 {
445 int retval = MPI_SUCCESS;
446
447           smpi_bench_end();
448
449           retval = smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
450
451           smpi_bench_begin();
452           return retval;
453 }
454
455
456
457 /**
458  * MPI_Allreduce
459  *
460  * Same as MPI_REDUCE except that the result appears in the receive buffer of all the group members.
461  **/
462 int SMPI_MPI_Allreduce( void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
463                     MPI_Op op, MPI_Comm comm )
464 {
465 int retval = MPI_SUCCESS;
466 int root=0;  // arbitrary choice
467
468           smpi_bench_end();
469
470           retval = smpi_mpi_reduce( sendbuf, recvbuf, count, datatype, op, root, comm);
471           if (MPI_SUCCESS != retval)
472                     return(retval);
473
474           retval = smpi_mpi_bcast( sendbuf, count, datatype, root, comm);
475
476           smpi_bench_end();
477           return( retval );
478 }
479
480
481 /**
482  * MPI_Scatter user entry point
483  **/
484 int SMPI_MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype datatype, 
485                          void *recvbuf, int recvcount, MPI_Datatype recvtype,
486                            int root, MPI_Comm comm)
487 {
488   int retval = MPI_SUCCESS;
489   int i;
490   int cnt=0;  
491   int rank;
492   int tag=0;
493   char *cptr;  // to manipulate the void * buffers
494   smpi_mpi_request_t *requests;
495   smpi_mpi_request_t request;
496   smpi_mpi_status_t status;
497
498
499   smpi_bench_end();
500
501   rank = smpi_mpi_comm_rank(comm);
502
503   requests = xbt_malloc((comm->size-1) * sizeof(smpi_mpi_request_t));
504   if (rank == root) {
505           // i am the root: distribute my sendbuf
506           //print_buffer_int(sendbuf, comm->size, xbt_strdup("rcvbuf"), rank);
507           cptr = sendbuf;
508           for (i=0; i < comm->size; i++) {
509                   if ( i!=root ) { // send to processes ...
510
511                           retval = smpi_create_request((void *)cptr, sendcount, 
512                                           datatype, root, i, tag, comm, &(requests[cnt]));
513                           if (NULL != requests[cnt] && MPI_SUCCESS == retval) {
514                                   if (MPI_SUCCESS == retval) {
515                                           smpi_mpi_isend(requests[cnt]);
516                                   }
517                                   }
518                                   cnt++;
519                         } 
520                         else { // ... except if it's me.
521                                   memcpy(recvbuf, (void *)cptr, recvcount*recvtype->size*sizeof(char));
522                         }
523                   cptr += sendcount*datatype->size;
524             }
525             for(i=0; i<cnt; i++) { // wait for send to complete
526                             /* FIXME: waitall() should be slightly better */
527                             smpi_mpi_wait(requests[i], &status);
528                             xbt_mallocator_release(smpi_global->request_mallocator, requests[i]);
529
530             }
531   } 
532   else {  // i am a non-root process: wait data from the root
533             retval = smpi_create_request(recvbuf,recvcount, 
534                                   recvtype, root, rank, tag, comm, &request);
535             if (NULL != request && MPI_SUCCESS == retval) {
536                         if (MPI_SUCCESS == retval) {
537                                   smpi_mpi_irecv(request);
538                         }
539             }
540             smpi_mpi_wait(request, &status);
541             xbt_mallocator_release(smpi_global->request_mallocator, request);
542   }
543   xbt_free(requests);
544
545   smpi_bench_begin();
546
547   return retval;
548 }
549
550
551 /**
552  * MPI_Alltoall user entry point
553  * 
554  * Uses the logic of OpenMPI (upto 1.2.7 or greater) for the optimizations
555  * ompi/mca/coll/tuned/coll_tuned_module.c
556  **/
557 int SMPI_MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype datatype, 
558                          void *recvbuf, int recvcount, MPI_Datatype recvtype,
559                            MPI_Comm comm)
560 {
561   int retval = MPI_SUCCESS;
562   int block_dsize;
563   int rank;
564
565   smpi_bench_end();
566
567   rank = smpi_mpi_comm_rank(comm);
568   block_dsize = datatype->size * sendcount;
569   DEBUG2("<%d> optimized alltoall() called. Block size sent to each rank: %d bytes.",rank,block_dsize);
570
571   if ((block_dsize < 200) && (comm->size > 12)) {
572             retval = smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, datatype,
573                                   recvbuf, recvcount, recvtype, comm);
574
575   } else if (block_dsize < 3000) {
576             retval = smpi_coll_tuned_alltoall_basic_linear(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