Logo AND Algorithmique Numérique Distribuée

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