Logo AND Algorithmique Numérique Distribuée

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