Logo AND Algorithmique Numérique Distribuée

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