Logo AND Algorithmique Numérique Distribuée

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