Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
7bb89d2505693a9d5bf4a8370130bfc37acab54b
[simgrid.git] / src / smpi / smpi_mpi.c
1 /* $Id$tag */
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 /* MPI User level calls */
11
12 int MPI_Init(int* argc, char*** argv) {
13   smpi_process_init(argc, argv);
14   smpi_bench_begin(-1, NULL);
15   return MPI_SUCCESS;
16 }
17
18 int MPI_Finalize(void) {
19   smpi_bench_end(-1, NULL);
20   smpi_process_destroy();
21   return MPI_SUCCESS;
22 }
23
24 int MPI_Init_thread(int* argc, char*** argv, int required, int* provided) {
25   if(provided != NULL) {
26     *provided = MPI_THREAD_MULTIPLE;
27   }
28   return MPI_Init(argc, argv);
29 }
30
31 int MPI_Query_thread(int* provided) {
32   int retval;
33
34   smpi_bench_end(-1, NULL);
35   if(provided == NULL) {
36     retval = MPI_ERR_ARG;
37   } else {
38     *provided = MPI_THREAD_MULTIPLE;
39     retval = MPI_SUCCESS;
40   }
41   smpi_bench_begin(-1, NULL);
42   return retval;
43 }
44
45 int MPI_Is_thread_main(int* flag) {
46   int retval;
47
48   smpi_bench_end(-1, NULL);
49   if(flag == NULL) {
50     retval = MPI_ERR_ARG;
51   } else {
52     *flag = smpi_process_index() == 0;
53     retval = MPI_SUCCESS;
54   }
55   smpi_bench_begin(-1, NULL);
56   return retval;
57 }
58
59 int MPI_Abort(MPI_Comm comm, int errorcode) {
60   smpi_bench_end(-1, NULL);
61   smpi_process_destroy();
62   // FIXME: should kill all processes in comm instead
63   SIMIX_process_kill(SIMIX_process_self());
64   return MPI_SUCCESS;
65 }
66
67 double MPI_Wtime(void) {
68   double time;
69
70   smpi_bench_end(-1, NULL);
71   time = SIMIX_get_clock();
72   smpi_bench_begin(-1, NULL);
73   return time;
74 }
75
76 int MPI_Type_size(MPI_Datatype datatype, size_t* size) {
77   int retval;
78
79   smpi_bench_end(-1, NULL);
80   if(datatype == MPI_DATATYPE_NULL) {
81     retval = MPI_ERR_TYPE;
82   } else if(size == NULL) {
83     retval = MPI_ERR_ARG;
84   } else {
85     *size = smpi_datatype_size(datatype);
86     retval = MPI_SUCCESS;
87   }
88   smpi_bench_begin(-1, NULL);
89   return retval;
90 }
91
92 int MPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint* lb, MPI_Aint* extent) {
93   int retval;
94
95   smpi_bench_end(-1, NULL);
96   if(datatype == MPI_DATATYPE_NULL) {
97     retval = MPI_ERR_TYPE;
98   } else if(lb == NULL || extent == NULL) {
99     retval = MPI_ERR_ARG;
100   } else {
101     retval = smpi_datatype_extent(datatype, lb, extent);
102   }
103   smpi_bench_begin(-1, NULL);
104   return retval;
105 }
106
107 int MPI_Type_lb(MPI_Datatype datatype, MPI_Aint* disp) {
108   int retval;
109
110   smpi_bench_end(-1, NULL);
111   if(datatype == MPI_DATATYPE_NULL) {
112     retval = MPI_ERR_TYPE;
113   } else if(disp == NULL) {
114     retval = MPI_ERR_ARG;
115   } else {
116     *disp = smpi_datatype_lb(datatype);
117     retval = MPI_SUCCESS;
118   }
119   smpi_bench_begin(-1, NULL);
120   return retval;
121 }
122
123 int MPI_Type_ub(MPI_Datatype datatype, MPI_Aint* disp) {
124   int retval;
125
126   smpi_bench_end(-1, NULL);
127   if(datatype == MPI_DATATYPE_NULL) {
128     retval = MPI_ERR_TYPE;
129   } else if(disp == NULL) {
130     retval = MPI_ERR_ARG;
131   } else {
132     *disp = smpi_datatype_ub(datatype);
133     retval = MPI_SUCCESS;
134   }
135   smpi_bench_begin(-1, NULL);
136   return retval;
137 }
138
139 int MPI_Op_create(MPI_User_function* function, int commute, MPI_Op* op) {
140   int retval;
141
142   smpi_bench_end(-1, NULL);
143   if(function == NULL || op == NULL) {
144     retval = MPI_ERR_ARG;
145   } else {
146     *op = smpi_op_new(function, commute);
147     retval = MPI_SUCCESS;
148   }
149   smpi_bench_begin(-1, NULL);
150   return retval;
151 }
152
153 int MPI_Op_free(MPI_Op* op) {
154   int retval;
155
156   smpi_bench_end(-1, NULL);
157   if(op == NULL) {
158     retval = MPI_ERR_ARG;
159   } else if(*op == MPI_OP_NULL) {
160     retval = MPI_ERR_OP;
161   } else {
162     smpi_op_destroy(*op);
163     *op = MPI_OP_NULL;
164     retval = MPI_SUCCESS;
165   }
166   smpi_bench_begin(-1, NULL);
167   return retval;
168 }
169
170 int MPI_Group_free(MPI_Group *group) {
171   int retval;
172
173   smpi_bench_end(-1, NULL);
174   if(group == NULL) {
175     retval = MPI_ERR_ARG;
176   } else {
177     smpi_group_destroy(*group);
178     *group = MPI_GROUP_NULL;
179     retval = MPI_SUCCESS;
180   }
181   smpi_bench_begin(-1, NULL);
182   return retval;
183 }
184
185 int MPI_Group_size(MPI_Group group, int* size) {
186   int retval;
187
188   smpi_bench_end(-1, NULL);
189   if(group == MPI_GROUP_NULL) {
190     retval = MPI_ERR_GROUP;
191   } else if(size == NULL) {
192     retval = MPI_ERR_ARG;
193   } else {
194     *size = smpi_group_size(group);
195     retval = MPI_SUCCESS;
196   }
197   smpi_bench_begin(-1, NULL);
198   return retval;
199 }
200
201 int MPI_Group_rank(MPI_Group group, int* rank) {
202   int retval;
203
204   smpi_bench_end(-1, NULL);
205   if(group == MPI_GROUP_NULL) {
206     retval = MPI_ERR_GROUP;
207   } else if(rank == NULL) {
208     retval = MPI_ERR_ARG;
209   } else {
210     *rank = smpi_group_rank(group, smpi_process_index());
211     retval = MPI_SUCCESS;
212   }
213   smpi_bench_begin(-1, NULL);
214   return retval;
215 }
216
217 int MPI_Group_translate_ranks (MPI_Group group1, int n, int* ranks1, MPI_Group group2, int* ranks2) {
218   int retval, i, index;
219
220   smpi_bench_end(-1, NULL);
221   if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
222     retval = MPI_ERR_GROUP;
223   } else {
224     for(i = 0; i < n; i++) {
225       index = smpi_group_index(group1, ranks1[i]);
226       ranks2[i] = smpi_group_rank(group2, index);
227     }
228     retval = MPI_SUCCESS;
229   }
230   smpi_bench_begin(-1, NULL);
231   return retval;
232 }
233
234 int MPI_Group_compare(MPI_Group group1, MPI_Group group2, int* result) {
235   int retval;
236
237   smpi_bench_end(-1, NULL);
238   if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
239     retval = MPI_ERR_GROUP;
240   } else if(result == NULL) {
241     retval = MPI_ERR_ARG;
242   } else {
243     *result = smpi_group_compare(group1, group2);
244     retval = MPI_SUCCESS;
245   }
246   smpi_bench_begin(-1, NULL);
247   return retval;
248 }
249
250 int MPI_Group_union(MPI_Group group1, MPI_Group group2, MPI_Group* newgroup) {
251   int retval, i, proc1, proc2, size, size2;
252
253   smpi_bench_end(-1, NULL);
254   if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
255     retval = MPI_ERR_GROUP;
256   } else if(newgroup == NULL) {
257     retval = MPI_ERR_ARG;
258   } else {
259     size = smpi_group_size(group1);
260     size2 = smpi_group_size(group2);
261     for(i = 0; i < size2; i++) {
262       proc2 = smpi_group_index(group2, i);
263       proc1 = smpi_group_rank(group1, proc2);
264       if(proc1 == MPI_UNDEFINED) {
265         size++;
266       }
267     }
268     if(size == 0) {
269       *newgroup = MPI_GROUP_EMPTY;
270     } else {
271       *newgroup = smpi_group_new(size);
272       size2 = smpi_group_size(group1);
273       for(i = 0; i < size2; i++) {
274         proc1 = smpi_group_index(group1, i);
275         smpi_group_set_mapping(*newgroup, proc1, i);
276       }
277       for(i = size2; i < size; i++) {
278         proc2 = smpi_group_index(group2, i - size2);
279         smpi_group_set_mapping(*newgroup, proc2, i);
280       }
281     }
282     smpi_group_use(*newgroup);
283     retval = MPI_SUCCESS;
284   }
285   smpi_bench_begin(-1, NULL);
286   return retval;
287 }
288
289 int MPI_Group_intersection(MPI_Group group1, MPI_Group group2, MPI_Group* newgroup) {
290    int retval, i, proc1, proc2, size, size2;
291
292   smpi_bench_end(-1, NULL);
293   if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
294     retval = MPI_ERR_GROUP;
295   } else if(newgroup == NULL) {
296     retval = MPI_ERR_ARG;
297   } else {
298     size = smpi_group_size(group1);
299     size2 = smpi_group_size(group2);
300     for(i = 0; i < size2; i++) {
301       proc2 = smpi_group_index(group2, i);
302       proc1 = smpi_group_rank(group1, proc2);
303       if(proc1 == MPI_UNDEFINED) {
304         size--;
305       }
306     }
307     if(size == 0) {
308       *newgroup = MPI_GROUP_EMPTY;
309     } else {
310       *newgroup = smpi_group_new(size);
311       size2 = smpi_group_size(group1);
312       for(i = 0; i < size2; i++) {
313         proc1 = smpi_group_index(group1, i);
314         proc2 = smpi_group_rank(group2, proc1);
315         if(proc2 != MPI_UNDEFINED) {
316           smpi_group_set_mapping(*newgroup, proc1, i);
317         }
318       }
319     }
320     smpi_group_use(*newgroup);
321     retval = MPI_SUCCESS;
322   }
323   smpi_bench_begin(-1, NULL);
324   return retval;
325 }
326
327 int MPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group* newgroup) {
328   int retval, i, proc1, proc2, size, size2;
329
330   smpi_bench_end(-1, NULL);
331   if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
332     retval = MPI_ERR_GROUP;
333   } else if(newgroup == NULL) {
334     retval = MPI_ERR_ARG;
335   } else {
336     size = size2 = smpi_group_size(group1);
337     for(i = 0; i < size2; i++) {
338       proc1 = smpi_group_index(group1, i);
339       proc2 = smpi_group_rank(group2, proc1);
340       if(proc2 != MPI_UNDEFINED) {
341         size--;
342       }
343     }
344     if(size == 0) {
345       *newgroup = MPI_GROUP_EMPTY;
346     } else {
347       *newgroup = smpi_group_new(size);
348       for(i = 0; i < size2; i++) {
349         proc1 = smpi_group_index(group1, i);
350         proc2 = smpi_group_rank(group2, proc1);
351         if(proc2 == MPI_UNDEFINED) {
352           smpi_group_set_mapping(*newgroup, proc1, i);
353         }
354       }
355     }
356     smpi_group_use(*newgroup);
357     retval = MPI_SUCCESS;
358   }
359   smpi_bench_begin(-1, NULL);
360   return retval;
361 }
362
363 int MPI_Group_incl(MPI_Group group, int n, int* ranks, MPI_Group* newgroup) {
364   int retval, i, index;
365
366   smpi_bench_end(-1, NULL);
367   if(group == MPI_GROUP_NULL) {
368     retval = MPI_ERR_GROUP;
369   } else if(newgroup == NULL) {
370     retval = MPI_ERR_ARG;
371   } else {
372     if(n == 0) {
373       *newgroup = MPI_GROUP_EMPTY;
374     } else if(n == smpi_group_size(group)) {
375       *newgroup = group;
376     } else {
377       *newgroup = smpi_group_new(n);
378       for(i = 0; i < n; i++) {
379         index = smpi_group_index(group, ranks[i]);
380         smpi_group_set_mapping(*newgroup, index, i);
381       }
382     }
383     smpi_group_use(*newgroup);
384     retval = MPI_SUCCESS;
385   }
386   smpi_bench_begin(-1, NULL);
387   return retval;
388 }
389
390 int MPI_Group_excl(MPI_Group group, int n, int* ranks, MPI_Group* newgroup) {
391   int retval, i, size, rank, index;
392
393   smpi_bench_end(-1, NULL);
394   if(group == MPI_GROUP_NULL) {
395     retval = MPI_ERR_GROUP;
396   } else if(newgroup == NULL) {
397     retval = MPI_ERR_ARG;
398   } else {
399     if(n == 0) {
400       *newgroup = group;
401     } else if(n == smpi_group_size(group)) {
402       *newgroup = MPI_GROUP_EMPTY;
403     } else {
404       size = smpi_group_size(group) - n;
405       *newgroup = smpi_group_new(size);
406       rank = 0;
407       while(rank < size) {
408         for(i = 0; i < n; i++) {
409           if(ranks[i] == rank) {
410             break;
411           }
412         }
413         if(i >= n) {
414           index = smpi_group_index(group, rank);
415           smpi_group_set_mapping(*newgroup, index, rank);
416           rank++;
417         }
418       }
419     }
420     smpi_group_use(*newgroup);
421     retval = MPI_SUCCESS;
422   }
423   smpi_bench_begin(-1, NULL);
424   return retval;
425 }
426
427 int MPI_Group_range_incl(MPI_Group group, int n, int ranges[][3], MPI_Group* newgroup) {
428   int retval, i, j, rank, size, index;
429
430   smpi_bench_end(-1, NULL);
431   if(group == MPI_GROUP_NULL) {
432     retval = MPI_ERR_GROUP;
433   } else if(newgroup == NULL) {
434     retval = MPI_ERR_ARG;
435   } else {
436     if(n == 0) {
437       *newgroup = MPI_GROUP_EMPTY;
438     } else {
439       size = 0;
440       for(i = 0; i < n; i++) {
441         for(rank = ranges[i][0]; /* First */
442             rank >= 0 && rank <= ranges[i][1]; /* Last */
443             rank += ranges[i][2] /* Stride */) {
444           size++;
445         }
446       }
447       if(size == smpi_group_size(group)) {
448         *newgroup = group;
449       } else {
450         *newgroup = smpi_group_new(size);
451         j = 0;
452         for(i = 0; i < n; i++) {
453           for(rank = ranges[i][0]; /* First */
454               rank >= 0 && rank <= ranges[i][1]; /* Last */
455               rank += ranges[i][2] /* Stride */) {
456             index = smpi_group_index(group, rank);
457             smpi_group_set_mapping(*newgroup, index, j);
458             j++;
459           }
460         }
461       }
462     }
463     smpi_group_use(*newgroup);
464     retval = MPI_SUCCESS;
465   }
466   smpi_bench_begin(-1, NULL);
467   return retval;
468 }
469
470 int MPI_Group_range_excl(MPI_Group group, int n, int ranges[][3], MPI_Group* newgroup) {
471   int retval, i, newrank, rank, size, index, add;
472
473   smpi_bench_end(-1, NULL);
474   if(group == MPI_GROUP_NULL) {
475     retval = MPI_ERR_GROUP;
476   } else if(newgroup == NULL) {
477     retval = MPI_ERR_ARG;
478   } else {
479     if(n == 0) {
480       *newgroup = group;
481     } else {
482       size = smpi_group_size(group);
483       for(i = 0; i < n; i++) {
484         for(rank = ranges[i][0]; /* First */
485             rank >= 0 && rank <= ranges[i][1]; /* Last */
486             rank += ranges[i][2] /* Stride */) {
487           size--;
488         }
489       }
490       if(size == 0) {
491         *newgroup = MPI_GROUP_EMPTY;
492       } else {
493         *newgroup = smpi_group_new(size);
494         newrank = 0;
495         while(newrank < size) {
496           for(i = 0; i < n; i++) {
497             add = 1;
498             for(rank = ranges[i][0]; /* First */
499                 rank >= 0 && rank <= ranges[i][1]; /* Last */
500                 rank += ranges[i][2] /* Stride */) {
501               if(rank == newrank) {
502                 add = 0;
503                 break;
504               }
505             }
506             if(add == 1) {
507               index = smpi_group_index(group, newrank);
508               smpi_group_set_mapping(*newgroup, index, newrank);
509             }
510           }
511         }
512       }
513     }
514     smpi_group_use(*newgroup);
515     retval = MPI_SUCCESS;
516   }
517   smpi_bench_begin(-1, NULL);
518   return retval;
519 }
520
521 int MPI_Comm_rank(MPI_Comm comm, int* rank) {
522   int retval;
523
524   smpi_bench_end(-1, NULL);
525   if(comm == MPI_COMM_NULL) {
526     retval = MPI_ERR_COMM;
527   } else {
528     *rank = smpi_comm_rank(comm);
529     retval = MPI_SUCCESS;
530   }
531   smpi_bench_begin(-1, NULL);
532   return retval;
533 }
534
535 int MPI_Comm_size(MPI_Comm comm, int* size) {
536   int retval;
537
538   smpi_bench_end(-1, NULL);
539   if(comm == MPI_COMM_NULL) {
540     retval = MPI_ERR_COMM;
541   } else if(size == NULL) {
542     retval = MPI_ERR_ARG;
543   } else {
544     *size = smpi_comm_size(comm);
545     retval = MPI_SUCCESS;
546   }
547   smpi_bench_begin(-1, NULL);
548   return retval;
549 }
550
551 int MPI_Comm_group(MPI_Comm comm, MPI_Group* group) {
552   int retval;
553
554   smpi_bench_end(-1, NULL);
555   if(comm == MPI_COMM_NULL) {
556     retval = MPI_ERR_COMM;
557   } else if(group == NULL) {
558     retval = MPI_ERR_ARG;
559   } else {
560     *group = smpi_comm_group(comm);
561     retval = MPI_SUCCESS;
562   }
563   smpi_bench_begin(-1, NULL);
564   return retval;
565 }
566
567 int MPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int* result) {
568   int retval;
569
570   smpi_bench_end(-1, NULL);
571   if(comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
572     retval = MPI_ERR_COMM;
573   } else if(result == NULL) {
574     retval = MPI_ERR_ARG;
575   } else {
576     if(comm1 == comm2) { /* Same communicators means same groups */
577       *result = MPI_IDENT;
578     } else {
579       *result = smpi_group_compare(smpi_comm_group(comm1), smpi_comm_group(comm2));
580       if(*result == MPI_IDENT) {
581         *result = MPI_CONGRUENT;
582       }
583     }
584     retval = MPI_SUCCESS;
585   }
586   smpi_bench_begin(-1, NULL);
587   return retval;
588 }
589
590 int MPI_Comm_dup(MPI_Comm comm, MPI_Comm* newcomm) {
591   int retval;
592
593   smpi_bench_end(-1, NULL);
594   if(comm == MPI_COMM_NULL) {
595     retval = MPI_ERR_COMM;
596   } else if(newcomm == NULL) {
597     retval = MPI_ERR_ARG;
598   } else {
599     *newcomm = smpi_comm_new(smpi_comm_group(comm));
600     retval = MPI_SUCCESS;
601   }
602   smpi_bench_begin(-1, NULL);
603   return retval;
604 }
605
606 int MPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm* newcomm) {
607   int retval;
608
609   smpi_bench_end(-1, NULL);
610   if(comm == MPI_COMM_NULL) {
611     retval = MPI_ERR_COMM;
612   } else if(group == MPI_GROUP_NULL) {
613     retval = MPI_ERR_GROUP;
614   } else if(newcomm == NULL) {
615     retval = MPI_ERR_ARG;
616   } else {
617     *newcomm = smpi_comm_new(group);
618     retval = MPI_SUCCESS;
619   }
620   smpi_bench_begin(-1, NULL);
621   return retval;
622 }
623
624 int MPI_Comm_free(MPI_Comm* comm) {
625   int retval;
626
627   smpi_bench_end(-1, NULL);
628   if(comm == NULL) {
629     retval = MPI_ERR_ARG;
630   } else if(*comm == MPI_COMM_NULL) {
631     retval = MPI_ERR_COMM;
632   } else {
633     smpi_comm_destroy(*comm);
634     *comm = MPI_COMM_NULL;
635     retval = MPI_SUCCESS;
636   }
637   smpi_bench_begin(-1, NULL);
638   return retval;
639 }
640
641 int MPI_Irecv(void* buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request* request) {
642   int retval;
643   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
644
645   smpi_bench_end(rank, "Irecv");
646   if(request == NULL) {
647     retval = MPI_ERR_ARG;
648   } else if (comm == MPI_COMM_NULL) {
649     retval = MPI_ERR_COMM;
650   } else {
651     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
652     retval = MPI_SUCCESS;
653   }
654   smpi_bench_begin(rank, "Irecv");
655   return retval;
656 }
657
658 int MPI_Isend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
659   int retval;
660   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
661
662   smpi_bench_end(rank, "Isend");
663   if(request == NULL) {
664     retval = MPI_ERR_ARG;
665   } else if (comm == MPI_COMM_NULL) {
666     retval = MPI_ERR_COMM;
667   } else {
668     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
669     retval = MPI_SUCCESS;
670   }
671   smpi_bench_begin(rank, "Isend");
672   return retval;
673 }
674
675 int MPI_Recv(void* buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status* status) {
676   int retval;
677   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
678
679   smpi_bench_end(rank, "Recv");
680   if (comm == MPI_COMM_NULL) {
681     retval = MPI_ERR_COMM;
682   } else {
683     smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
684     retval = MPI_SUCCESS;
685   }
686   smpi_bench_begin(rank, "Recv");
687   return retval;
688 }
689
690 int MPI_Send(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
691   int retval;
692   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
693
694   smpi_bench_end(rank, "Send");
695   if (comm == MPI_COMM_NULL) {
696     retval = MPI_ERR_COMM;
697   } else {
698     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
699     retval = MPI_SUCCESS;
700   }
701   smpi_bench_begin(rank, "Send");
702   return retval;
703 }
704
705 int MPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf, int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status) {
706   int retval;
707   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
708
709   smpi_bench_end(rank, "Sendrecv");
710   if (comm == MPI_COMM_NULL) {
711     retval = MPI_ERR_COMM;
712   } else if (sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
713     retval = MPI_ERR_TYPE;
714   } else {
715     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src, recvtag, comm, status);
716     retval = MPI_SUCCESS;
717   }
718   smpi_bench_begin(rank, "Sendrecv");
719   return retval;
720 }
721
722 int MPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag, MPI_Comm comm, MPI_Status* status) {
723   //TODO: suboptimal implementation
724   void* recvbuf;
725   int retval, size;
726
727   size = smpi_datatype_size(datatype) * count;
728   recvbuf = xbt_new(char, size);
729   retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
730   memcpy(buf, recvbuf, size * sizeof(char));
731   xbt_free(recvbuf);
732   return retval;
733 }
734
735 int MPI_Test(MPI_Request* request, int* flag, MPI_Status* status) {
736   int retval;
737   int rank = request && (*request)->comm != MPI_COMM_NULL
738              ? smpi_comm_rank((*request)->comm)
739              : -1;
740
741   smpi_bench_end(rank, "Test");
742   if(request == NULL || flag == NULL) {
743     retval = MPI_ERR_ARG;
744   } else if(*request == MPI_REQUEST_NULL) {
745     retval = MPI_ERR_REQUEST;
746   } else {
747     *flag = smpi_mpi_test(request, status);
748     retval = MPI_SUCCESS;
749   }
750   smpi_bench_begin(rank, "Test");
751   return retval;
752 }
753
754 int MPI_Testany(int count, MPI_Request requests[], int* index, int* flag, MPI_Status* status) {
755   int retval;
756
757   smpi_bench_end(-1, NULL); //FIXME
758   if(index == NULL || flag == NULL) {
759     retval = MPI_ERR_ARG;
760   } else {
761     *flag = smpi_mpi_testany(count, requests, index, status);
762     retval = MPI_SUCCESS;
763   }
764   smpi_bench_begin(-1, NULL);
765   return retval;
766 }
767
768 int MPI_Wait(MPI_Request* request, MPI_Status* status) {
769   int retval;
770   int rank = request && (*request)->comm != MPI_COMM_NULL
771              ? smpi_comm_rank((*request)->comm)
772              : -1;
773
774   smpi_bench_end(rank, "Wait");
775   if(request == NULL) {
776     retval = MPI_ERR_ARG;
777   } else if(*request == MPI_REQUEST_NULL) {
778     retval = MPI_ERR_REQUEST;
779   } else {
780     smpi_mpi_wait(request, status);
781     retval = MPI_SUCCESS;
782   }
783   smpi_bench_begin(rank, "Wait");
784   return retval;
785 }
786
787 int MPI_Waitany(int count, MPI_Request requests[], int* index, MPI_Status* status) {
788   int retval;
789
790   smpi_bench_end(-1, NULL); //FIXME
791   if(index == NULL) {
792     retval = MPI_ERR_ARG;
793   } else {
794     *index = smpi_mpi_waitany(count, requests, status);
795     retval = MPI_SUCCESS;
796   }
797   smpi_bench_begin(-1, NULL);
798   return retval;
799 }
800
801 int MPI_Waitall(int count, MPI_Request requests[],  MPI_Status status[]) {
802   smpi_bench_end(-1, NULL); //FIXME
803   smpi_mpi_waitall(count, requests, status);
804   smpi_bench_begin(-1, NULL);
805   return MPI_SUCCESS;
806 }
807
808 int MPI_Waitsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[]) {
809   int retval;
810
811   smpi_bench_end(-1, NULL); //FIXME
812   if(outcount == NULL || indices == NULL) {
813     retval = MPI_ERR_ARG;
814   } else {
815     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
816     retval = MPI_SUCCESS;
817   }
818   smpi_bench_begin(-1, NULL);
819   return retval;
820 }
821
822 int MPI_Bcast(void* buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm) {
823   int retval;
824   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
825
826   smpi_bench_end(rank, "Bcast");
827   if(comm == MPI_COMM_NULL) {
828     retval = MPI_ERR_COMM;
829   } else {
830     smpi_mpi_bcast(buf, count, datatype, root, comm);
831     retval = MPI_SUCCESS;
832   }
833   smpi_bench_begin(rank, "Bcast");
834   return retval;
835 }
836
837 int MPI_Barrier(MPI_Comm comm) {
838   int retval;
839   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
840
841   smpi_bench_end(rank, "Barrier");
842   if(comm == MPI_COMM_NULL) {
843     retval = MPI_ERR_COMM;
844   } else {
845     smpi_mpi_barrier(comm);
846     retval = MPI_SUCCESS;
847   }
848   smpi_bench_begin(rank, "Barrier");
849   return retval;
850 }
851
852 int MPI_Gather(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm) {
853   int retval;
854   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
855
856   smpi_bench_end(rank, "Gather");
857   if(comm == MPI_COMM_NULL) {
858     retval = MPI_ERR_COMM;
859   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
860     retval = MPI_ERR_TYPE;
861   } else {
862     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
863     retval = MPI_SUCCESS;
864   }
865   smpi_bench_begin(rank, "Gather");
866   return retval;
867 }
868
869 int MPI_Gatherv(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int* recvcounts, int* displs, MPI_Datatype recvtype, int root, MPI_Comm comm) {
870   int retval;
871   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
872
873   smpi_bench_end(rank, "Gatherv");
874   if(comm == MPI_COMM_NULL) {
875     retval = MPI_ERR_COMM;
876   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
877     retval = MPI_ERR_TYPE;
878   } else if(recvcounts == NULL || displs == NULL) {
879     retval = MPI_ERR_ARG;
880   } else {
881     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm);
882     retval = MPI_SUCCESS;
883   }
884   smpi_bench_begin(rank, "Gatherv");
885   return retval;
886 }
887
888 int MPI_Allgather(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm) {
889   int retval;
890   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
891
892   smpi_bench_end(rank, "Allgather");
893   if(comm == MPI_COMM_NULL) {
894     retval = MPI_ERR_COMM;
895   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
896     retval = MPI_ERR_TYPE;
897   } else {
898     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
899     retval = MPI_SUCCESS;
900   }
901   smpi_bench_begin(rank, "Allgather");
902   return retval;
903 }
904
905 int MPI_Allgatherv(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int* recvcounts, int* displs, MPI_Datatype recvtype, MPI_Comm comm) {
906   int retval;
907   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
908
909   smpi_bench_end(rank, "Allgatherv");
910   if(comm == MPI_COMM_NULL) {
911     retval = MPI_ERR_COMM;
912   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
913     retval = MPI_ERR_TYPE;
914   } else if(recvcounts == NULL || displs == NULL) {
915     retval = MPI_ERR_ARG;
916   } else {
917     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
918     retval = MPI_SUCCESS;
919   }
920   smpi_bench_begin(rank, "Allgatherv");
921   return retval;
922 }
923
924 int MPI_Scatter(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm) {
925   int retval;
926   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
927
928   smpi_bench_end(rank, "Scatter");
929   if(comm == MPI_COMM_NULL) {
930     retval = MPI_ERR_COMM;
931   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
932     retval = MPI_ERR_TYPE;
933   } else {
934     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
935     retval = MPI_SUCCESS;
936   }
937   smpi_bench_begin(rank, "Scatter");
938   return retval;
939 }
940
941 int MPI_Scatterv(void* sendbuf, int* sendcounts, int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm) {
942   int retval;
943   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
944
945   smpi_bench_end(rank, "Scatterv");
946   if(comm == MPI_COMM_NULL) {
947     retval = MPI_ERR_COMM;
948   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
949     retval = MPI_ERR_TYPE;
950   } else if(sendcounts == NULL || displs == NULL) {
951     retval = MPI_ERR_ARG;
952   } else {
953     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
954     retval = MPI_SUCCESS;
955   }
956   smpi_bench_begin(rank, "Scatterv");
957   return retval;
958 }
959
960 int MPI_Reduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm) {
961   int retval;
962   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
963
964   smpi_bench_end(rank, "Reduce");
965   if(comm == MPI_COMM_NULL) {
966     retval = MPI_ERR_COMM;
967   } else if(datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
968     retval = MPI_ERR_ARG;
969   } else {
970     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
971     retval = MPI_SUCCESS;
972   }
973   smpi_bench_begin(rank, "Reduce");
974   return retval;
975 }
976
977 int MPI_Allreduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) {
978   int retval;
979   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
980
981   smpi_bench_end(rank, "Allreduce");
982   if(comm == MPI_COMM_NULL) {
983     retval = MPI_ERR_COMM;
984   } else if(datatype == MPI_DATATYPE_NULL) {
985     retval = MPI_ERR_TYPE;
986   } else if(op == MPI_OP_NULL) {
987     retval = MPI_ERR_OP;
988   } else {
989     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
990     retval = MPI_SUCCESS;
991   }
992   smpi_bench_begin(rank, "Allreduce");
993   return retval;
994 }
995
996 int MPI_Reduce_scatter(void* sendbuf, void* recvbuf, int* recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) {
997   int retval, i, size, count;
998   int* displs;
999   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1000
1001   smpi_bench_end(rank, "Reduce_scatter");
1002   if(comm == MPI_COMM_NULL) {
1003     retval = MPI_ERR_COMM;
1004   } else if(datatype == MPI_DATATYPE_NULL) {
1005     retval = MPI_ERR_TYPE;
1006   } else if(op == MPI_OP_NULL) {
1007     retval = MPI_ERR_OP;
1008   } else if(recvcounts == NULL) {
1009     retval = MPI_ERR_ARG;
1010   } else {
1011     /* arbitrarily choose root as rank 0 */
1012     /* TODO: faster direct implementation ? */
1013     size = smpi_comm_size(comm);
1014     count = 0;
1015     displs = xbt_new(int, size);
1016     for(i = 0; i < size; i++) {
1017       count += recvcounts[i];
1018       displs[i] = 0;
1019     }
1020     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1021     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf, recvcounts[rank], datatype, 0, comm);
1022     xbt_free(displs);
1023     retval = MPI_SUCCESS;
1024   }
1025   smpi_bench_begin(rank, "Reduce_scatter");
1026   return retval;
1027 }
1028
1029 int MPI_Alltoall(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm) {
1030   int retval, size, sendsize;
1031   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1032
1033   smpi_bench_end(rank, "Alltoall");
1034   if(comm == MPI_COMM_NULL) {
1035     retval = MPI_ERR_COMM;
1036   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
1037     retval = MPI_ERR_TYPE;
1038   } else {
1039     size = smpi_comm_size(comm);
1040     sendsize = smpi_datatype_size(sendtype) * sendcount;
1041     if(sendsize < 200 && size > 12) {
1042       retval = smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1043     } else if(sendsize < 3000) {
1044       retval = smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1045     } else {
1046       retval = smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1047     }
1048   }
1049   smpi_bench_begin(rank, "Alltoall");
1050   return retval;
1051 }
1052
1053 int MPI_Alltoallv(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype sendtype, void* recvbuf, int *recvcounts, int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm) {
1054   int retval;
1055   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1056
1057   smpi_bench_end(rank, "Alltoallv");
1058   if(comm == MPI_COMM_NULL) {
1059     retval = MPI_ERR_COMM;
1060   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
1061     retval = MPI_ERR_TYPE;
1062   } else if(sendcounts == NULL || senddisps == NULL || recvcounts == NULL || recvdisps == NULL) {
1063     retval = MPI_ERR_ARG;
1064   } else {
1065     retval = smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype, recvbuf, recvcounts, recvdisps, recvtype, comm); 
1066   }
1067   smpi_bench_begin(rank, "Alltoallv");
1068   return retval;
1069 }
1070
1071
1072 int MPI_Get_processor_name( char *name, int *resultlen ) {
1073   int retval = MPI_SUCCESS;
1074   smpi_bench_end(-1, NULL);
1075   strcpy( name , SIMIX_host_get_name(SIMIX_host_self()));
1076   name[MPI_MAX_PROCESSOR_NAME-1]='\0';
1077   *resultlen= strlen(name) > MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1078
1079   smpi_bench_begin(-1, NULL);
1080   return retval;
1081 }
1082