Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
MPI_Comm_split is back to life.
[simgrid.git] / src / smpi / smpi_mpi.c
1 /* Copyright (c) 2007, 2008, 2009, 2010. The SimGrid Team.
2  * All rights reserved.                                                     */
3
4 /* This program is free software; you can redistribute it and/or modify it
5   * under the terms of the license (GNU LGPL) which comes with this package. */
6
7 #include "private.h"
8 #include "smpi_coll_private.h"
9 #include "smpi_mpi_dt_private.h"
10
11 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_mpi, smpi,
12                                 "Logging specific to SMPI (mpi)");
13
14 #ifdef HAVE_TRACING
15 //this function need to be here because of the calls to smpi_bench
16 int TRACE_smpi_set_category(const char *category)
17 {
18   //need to end bench otherwise categories for execution tasks are wrong
19   smpi_bench_end (-1, NULL);
20   int ret;
21   if (!IS_TRACING){
22     ret = 1;
23   }else{
24     if (category != NULL) {
25       ret = TRACE_category(category);
26       TRACE_category_set(SIMIX_process_self(), category);
27     }else{
28       //if category is NULL, trace of platform is disabled for this process
29       TRACE_category_unset(SIMIX_process_self());
30       ret = 0;
31     }
32   }
33   //begin bench after changing process's category
34   smpi_bench_begin (-1, NULL);
35   return ret;
36 }
37 #endif
38
39 /* MPI User level calls */
40
41 int MPI_Init(int *argc, char ***argv)
42 {
43   smpi_process_init(argc, argv);
44 #ifdef HAVE_TRACING
45   TRACE_smpi_init(smpi_process_index());
46 #endif
47   smpi_bench_begin(-1, NULL);
48   return MPI_SUCCESS;
49 }
50
51 int MPI_Finalize(void)
52 {
53   smpi_bench_end(-1, NULL);
54 #ifdef HAVE_TRACING
55   TRACE_smpi_finalize(smpi_process_index());
56 #endif
57   smpi_process_destroy();
58   return MPI_SUCCESS;
59 }
60
61 int MPI_Init_thread(int *argc, char ***argv, int required, int *provided)
62 {
63   if (provided != NULL) {
64     *provided = MPI_THREAD_MULTIPLE;
65   }
66   return MPI_Init(argc, argv);
67 }
68
69 int MPI_Query_thread(int *provided)
70 {
71   int retval;
72
73   smpi_bench_end(-1, NULL);
74   if (provided == NULL) {
75     retval = MPI_ERR_ARG;
76   } else {
77     *provided = MPI_THREAD_MULTIPLE;
78     retval = MPI_SUCCESS;
79   }
80   smpi_bench_begin(-1, NULL);
81   return retval;
82 }
83
84 int MPI_Is_thread_main(int *flag)
85 {
86   int retval;
87
88   smpi_bench_end(-1, NULL);
89   if (flag == NULL) {
90     retval = MPI_ERR_ARG;
91   } else {
92     *flag = smpi_process_index() == 0;
93     retval = MPI_SUCCESS;
94   }
95   smpi_bench_begin(-1, NULL);
96   return retval;
97 }
98
99 int MPI_Abort(MPI_Comm comm, int errorcode)
100 {
101   smpi_bench_end(-1, NULL);
102   smpi_process_destroy();
103   // FIXME: should kill all processes in comm instead
104   SIMIX_process_kill(SIMIX_process_self());
105   return MPI_SUCCESS;
106 }
107
108 double MPI_Wtime(void)
109 {
110   double time;
111
112   smpi_bench_end(-1, NULL);
113   time = SIMIX_get_clock();
114   smpi_bench_begin(-1, NULL);
115   return time;
116 }
117
118 int MPI_Address(void *location, MPI_Aint * address)
119 {
120   int retval;
121
122   smpi_bench_end(-1, NULL);
123   if (!address) {
124     retval = MPI_ERR_ARG;
125   } else {
126     *address = (MPI_Aint) location;
127   }
128   smpi_bench_begin(-1, NULL);
129   return retval;
130 }
131
132 int MPI_Type_free(MPI_Datatype * datatype)
133 {
134   int retval;
135
136   smpi_bench_end(-1, NULL);
137   if (!datatype) {
138     retval = MPI_ERR_ARG;
139   } else {
140     // FIXME: always fail for now
141     retval = MPI_ERR_TYPE;
142   }
143   smpi_bench_begin(-1, NULL);
144   return retval;
145 }
146
147 int MPI_Type_size(MPI_Datatype datatype, int *size)
148 {
149   int retval;
150
151   smpi_bench_end(-1, NULL);
152   if (datatype == MPI_DATATYPE_NULL) {
153     retval = MPI_ERR_TYPE;
154   } else if (size == NULL) {
155     retval = MPI_ERR_ARG;
156   } else {
157     *size = (int) smpi_datatype_size(datatype);
158     retval = MPI_SUCCESS;
159   }
160   smpi_bench_begin(-1, NULL);
161   return retval;
162 }
163
164 int MPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb,
165                         MPI_Aint * extent)
166 {
167   int retval;
168
169   smpi_bench_end(-1, NULL);
170   if (datatype == MPI_DATATYPE_NULL) {
171     retval = MPI_ERR_TYPE;
172   } else if (lb == NULL || extent == NULL) {
173     retval = MPI_ERR_ARG;
174   } else {
175     retval = smpi_datatype_extent(datatype, lb, extent);
176   }
177   smpi_bench_begin(-1, NULL);
178   return retval;
179 }
180
181 int MPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
182 {
183   int retval;
184   MPI_Aint dummy;
185
186   smpi_bench_end(-1, NULL);
187   if (datatype == MPI_DATATYPE_NULL) {
188     retval = MPI_ERR_TYPE;
189   } else if (extent == NULL) {
190     retval = MPI_ERR_ARG;
191   } else {
192     retval = smpi_datatype_extent(datatype, &dummy, extent);
193   }
194   smpi_bench_begin(-1, NULL);
195   return retval;
196 }
197
198 int MPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
199 {
200   int retval;
201
202   smpi_bench_end(-1, NULL);
203   if (datatype == MPI_DATATYPE_NULL) {
204     retval = MPI_ERR_TYPE;
205   } else if (disp == NULL) {
206     retval = MPI_ERR_ARG;
207   } else {
208     *disp = smpi_datatype_lb(datatype);
209     retval = MPI_SUCCESS;
210   }
211   smpi_bench_begin(-1, NULL);
212   return retval;
213 }
214
215 int MPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
216 {
217   int retval;
218
219   smpi_bench_end(-1, NULL);
220   if (datatype == MPI_DATATYPE_NULL) {
221     retval = MPI_ERR_TYPE;
222   } else if (disp == NULL) {
223     retval = MPI_ERR_ARG;
224   } else {
225     *disp = smpi_datatype_ub(datatype);
226     retval = MPI_SUCCESS;
227   }
228   smpi_bench_begin(-1, NULL);
229   return retval;
230 }
231
232 int MPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
233 {
234   int retval;
235
236   smpi_bench_end(-1, NULL);
237   if (function == NULL || op == NULL) {
238     retval = MPI_ERR_ARG;
239   } else {
240     *op = smpi_op_new(function, commute);
241     retval = MPI_SUCCESS;
242   }
243   smpi_bench_begin(-1, NULL);
244   return retval;
245 }
246
247 int MPI_Op_free(MPI_Op * op)
248 {
249   int retval;
250
251   smpi_bench_end(-1, NULL);
252   if (op == NULL) {
253     retval = MPI_ERR_ARG;
254   } else if (*op == MPI_OP_NULL) {
255     retval = MPI_ERR_OP;
256   } else {
257     smpi_op_destroy(*op);
258     *op = MPI_OP_NULL;
259     retval = MPI_SUCCESS;
260   }
261   smpi_bench_begin(-1, NULL);
262   return retval;
263 }
264
265 int MPI_Group_free(MPI_Group * group)
266 {
267   int retval;
268
269   smpi_bench_end(-1, NULL);
270   if (group == NULL) {
271     retval = MPI_ERR_ARG;
272   } else {
273     smpi_group_destroy(*group);
274     *group = MPI_GROUP_NULL;
275     retval = MPI_SUCCESS;
276   }
277   smpi_bench_begin(-1, NULL);
278   return retval;
279 }
280
281 int MPI_Group_size(MPI_Group group, int *size)
282 {
283   int retval;
284
285   smpi_bench_end(-1, NULL);
286   if (group == MPI_GROUP_NULL) {
287     retval = MPI_ERR_GROUP;
288   } else if (size == NULL) {
289     retval = MPI_ERR_ARG;
290   } else {
291     *size = smpi_group_size(group);
292     retval = MPI_SUCCESS;
293   }
294   smpi_bench_begin(-1, NULL);
295   return retval;
296 }
297
298 int MPI_Group_rank(MPI_Group group, int *rank)
299 {
300   int retval;
301
302   smpi_bench_end(-1, NULL);
303   if (group == MPI_GROUP_NULL) {
304     retval = MPI_ERR_GROUP;
305   } else if (rank == NULL) {
306     retval = MPI_ERR_ARG;
307   } else {
308     *rank = smpi_group_rank(group, smpi_process_index());
309     retval = MPI_SUCCESS;
310   }
311   smpi_bench_begin(-1, NULL);
312   return retval;
313 }
314
315 int MPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1,
316                               MPI_Group group2, int *ranks2)
317 {
318   int retval, i, index;
319
320   smpi_bench_end(-1, NULL);
321   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
322     retval = MPI_ERR_GROUP;
323   } else {
324     for (i = 0; i < n; i++) {
325       index = smpi_group_index(group1, ranks1[i]);
326       ranks2[i] = smpi_group_rank(group2, index);
327     }
328     retval = MPI_SUCCESS;
329   }
330   smpi_bench_begin(-1, NULL);
331   return retval;
332 }
333
334 int MPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
335 {
336   int retval;
337
338   smpi_bench_end(-1, NULL);
339   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
340     retval = MPI_ERR_GROUP;
341   } else if (result == NULL) {
342     retval = MPI_ERR_ARG;
343   } else {
344     *result = smpi_group_compare(group1, group2);
345     retval = MPI_SUCCESS;
346   }
347   smpi_bench_begin(-1, NULL);
348   return retval;
349 }
350
351 int MPI_Group_union(MPI_Group group1, MPI_Group group2,
352                     MPI_Group * newgroup)
353 {
354   int retval, i, proc1, proc2, size, size2;
355
356   smpi_bench_end(-1, NULL);
357   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
358     retval = MPI_ERR_GROUP;
359   } else if (newgroup == NULL) {
360     retval = MPI_ERR_ARG;
361   } else {
362     size = smpi_group_size(group1);
363     size2 = smpi_group_size(group2);
364     for (i = 0; i < size2; i++) {
365       proc2 = smpi_group_index(group2, i);
366       proc1 = smpi_group_rank(group1, proc2);
367       if (proc1 == MPI_UNDEFINED) {
368         size++;
369       }
370     }
371     if (size == 0) {
372       *newgroup = MPI_GROUP_EMPTY;
373     } else {
374       *newgroup = smpi_group_new(size);
375       size2 = smpi_group_size(group1);
376       for (i = 0; i < size2; i++) {
377         proc1 = smpi_group_index(group1, i);
378         smpi_group_set_mapping(*newgroup, proc1, i);
379       }
380       for (i = size2; i < size; i++) {
381         proc2 = smpi_group_index(group2, i - size2);
382         smpi_group_set_mapping(*newgroup, proc2, i);
383       }
384     }
385     smpi_group_use(*newgroup);
386     retval = MPI_SUCCESS;
387   }
388   smpi_bench_begin(-1, NULL);
389   return retval;
390 }
391
392 int MPI_Group_intersection(MPI_Group group1, MPI_Group group2,
393                            MPI_Group * newgroup)
394 {
395   int retval, i, proc1, proc2, size, size2;
396
397   smpi_bench_end(-1, NULL);
398   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
399     retval = MPI_ERR_GROUP;
400   } else if (newgroup == NULL) {
401     retval = MPI_ERR_ARG;
402   } else {
403     size = smpi_group_size(group1);
404     size2 = smpi_group_size(group2);
405     for (i = 0; i < size2; i++) {
406       proc2 = smpi_group_index(group2, i);
407       proc1 = smpi_group_rank(group1, proc2);
408       if (proc1 == MPI_UNDEFINED) {
409         size--;
410       }
411     }
412     if (size == 0) {
413       *newgroup = MPI_GROUP_EMPTY;
414     } else {
415       *newgroup = smpi_group_new(size);
416       size2 = smpi_group_size(group1);
417       for (i = 0; i < size2; i++) {
418         proc1 = smpi_group_index(group1, i);
419         proc2 = smpi_group_rank(group2, proc1);
420         if (proc2 != MPI_UNDEFINED) {
421           smpi_group_set_mapping(*newgroup, proc1, i);
422         }
423       }
424     }
425     smpi_group_use(*newgroup);
426     retval = MPI_SUCCESS;
427   }
428   smpi_bench_begin(-1, NULL);
429   return retval;
430 }
431
432 int MPI_Group_difference(MPI_Group group1, MPI_Group group2,
433                          MPI_Group * newgroup)
434 {
435   int retval, i, proc1, proc2, size, size2;
436
437   smpi_bench_end(-1, NULL);
438   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
439     retval = MPI_ERR_GROUP;
440   } else if (newgroup == NULL) {
441     retval = MPI_ERR_ARG;
442   } else {
443     size = size2 = smpi_group_size(group1);
444     for (i = 0; i < size2; i++) {
445       proc1 = smpi_group_index(group1, i);
446       proc2 = smpi_group_rank(group2, proc1);
447       if (proc2 != MPI_UNDEFINED) {
448         size--;
449       }
450     }
451     if (size == 0) {
452       *newgroup = MPI_GROUP_EMPTY;
453     } else {
454       *newgroup = smpi_group_new(size);
455       for (i = 0; i < size2; i++) {
456         proc1 = smpi_group_index(group1, i);
457         proc2 = smpi_group_rank(group2, proc1);
458         if (proc2 == MPI_UNDEFINED) {
459           smpi_group_set_mapping(*newgroup, proc1, i);
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_incl(MPI_Group group, int n, int *ranks,
471                    MPI_Group * newgroup)
472 {
473   int retval, i, index;
474
475   smpi_bench_end(-1, NULL);
476   if (group == MPI_GROUP_NULL) {
477     retval = MPI_ERR_GROUP;
478   } else if (newgroup == NULL) {
479     retval = MPI_ERR_ARG;
480   } else {
481     if (n == 0) {
482       *newgroup = MPI_GROUP_EMPTY;
483     } else if (n == smpi_group_size(group)) {
484       *newgroup = group;
485     } else {
486       *newgroup = smpi_group_new(n);
487       for (i = 0; i < n; i++) {
488         index = smpi_group_index(group, ranks[i]);
489         smpi_group_set_mapping(*newgroup, index, i);
490       }
491     }
492     smpi_group_use(*newgroup);
493     retval = MPI_SUCCESS;
494   }
495   smpi_bench_begin(-1, NULL);
496   return retval;
497 }
498
499 int MPI_Group_excl(MPI_Group group, int n, int *ranks,
500                    MPI_Group * newgroup)
501 {
502   int retval, i, size, rank, index;
503
504   smpi_bench_end(-1, NULL);
505   if (group == MPI_GROUP_NULL) {
506     retval = MPI_ERR_GROUP;
507   } else if (newgroup == NULL) {
508     retval = MPI_ERR_ARG;
509   } else {
510     if (n == 0) {
511       *newgroup = group;
512     } else if (n == smpi_group_size(group)) {
513       *newgroup = MPI_GROUP_EMPTY;
514     } else {
515       size = smpi_group_size(group) - n;
516       *newgroup = smpi_group_new(size);
517       rank = 0;
518       while (rank < size) {
519         for (i = 0; i < n; i++) {
520           if (ranks[i] == rank) {
521             break;
522           }
523         }
524         if (i >= n) {
525           index = smpi_group_index(group, rank);
526           smpi_group_set_mapping(*newgroup, index, rank);
527           rank++;
528         }
529       }
530     }
531     smpi_group_use(*newgroup);
532     retval = MPI_SUCCESS;
533   }
534   smpi_bench_begin(-1, NULL);
535   return retval;
536 }
537
538 int MPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
539                          MPI_Group * newgroup)
540 {
541   int retval, i, j, rank, size, index;
542
543   smpi_bench_end(-1, NULL);
544   if (group == MPI_GROUP_NULL) {
545     retval = MPI_ERR_GROUP;
546   } else if (newgroup == NULL) {
547     retval = MPI_ERR_ARG;
548   } else {
549     if (n == 0) {
550       *newgroup = MPI_GROUP_EMPTY;
551     } else {
552       size = 0;
553       for (i = 0; i < n; i++) {
554         for (rank = ranges[i][0];       /* First */
555              rank >= 0 && rank <= ranges[i][1]; /* Last */
556              rank += ranges[i][2] /* Stride */ ) {
557           size++;
558         }
559       }
560       if (size == smpi_group_size(group)) {
561         *newgroup = group;
562       } else {
563         *newgroup = smpi_group_new(size);
564         j = 0;
565         for (i = 0; i < n; i++) {
566           for (rank = ranges[i][0];     /* First */
567                rank >= 0 && rank <= ranges[i][1];       /* Last */
568                rank += ranges[i][2] /* Stride */ ) {
569             index = smpi_group_index(group, rank);
570             smpi_group_set_mapping(*newgroup, index, j);
571             j++;
572           }
573         }
574       }
575     }
576     smpi_group_use(*newgroup);
577     retval = MPI_SUCCESS;
578   }
579   smpi_bench_begin(-1, NULL);
580   return retval;
581 }
582
583 int MPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
584                          MPI_Group * newgroup)
585 {
586   int retval, i, newrank, rank, size, index, add;
587
588   smpi_bench_end(-1, NULL);
589   if (group == MPI_GROUP_NULL) {
590     retval = MPI_ERR_GROUP;
591   } else if (newgroup == NULL) {
592     retval = MPI_ERR_ARG;
593   } else {
594     if (n == 0) {
595       *newgroup = group;
596     } else {
597       size = smpi_group_size(group);
598       for (i = 0; i < n; i++) {
599         for (rank = ranges[i][0];       /* First */
600              rank >= 0 && rank <= ranges[i][1]; /* Last */
601              rank += ranges[i][2] /* Stride */ ) {
602           size--;
603         }
604       }
605       if (size == 0) {
606         *newgroup = MPI_GROUP_EMPTY;
607       } else {
608         *newgroup = smpi_group_new(size);
609         newrank = 0;
610         while (newrank < size) {
611           for (i = 0; i < n; i++) {
612             add = 1;
613             for (rank = ranges[i][0];   /* First */
614                  rank >= 0 && rank <= ranges[i][1];     /* Last */
615                  rank += ranges[i][2] /* Stride */ ) {
616               if (rank == newrank) {
617                 add = 0;
618                 break;
619               }
620             }
621             if (add == 1) {
622               index = smpi_group_index(group, newrank);
623               smpi_group_set_mapping(*newgroup, index, newrank);
624             }
625           }
626         }
627       }
628     }
629     smpi_group_use(*newgroup);
630     retval = MPI_SUCCESS;
631   }
632   smpi_bench_begin(-1, NULL);
633   return retval;
634 }
635
636 int MPI_Comm_rank(MPI_Comm comm, int *rank)
637 {
638   int retval;
639
640   smpi_bench_end(-1, NULL);
641   if (comm == MPI_COMM_NULL) {
642     retval = MPI_ERR_COMM;
643   } else {
644     *rank = smpi_comm_rank(comm);
645     retval = MPI_SUCCESS;
646   }
647   smpi_bench_begin(-1, NULL);
648   return retval;
649 }
650
651 int MPI_Comm_size(MPI_Comm comm, int *size)
652 {
653   int retval;
654
655   smpi_bench_end(-1, NULL);
656   if (comm == MPI_COMM_NULL) {
657     retval = MPI_ERR_COMM;
658   } else if (size == NULL) {
659     retval = MPI_ERR_ARG;
660   } else {
661     *size = smpi_comm_size(comm);
662     retval = MPI_SUCCESS;
663   }
664   smpi_bench_begin(-1, NULL);
665   return retval;
666 }
667
668 int MPI_Comm_group(MPI_Comm comm, MPI_Group * group)
669 {
670   int retval;
671
672   smpi_bench_end(-1, NULL);
673   if (comm == MPI_COMM_NULL) {
674     retval = MPI_ERR_COMM;
675   } else if (group == NULL) {
676     retval = MPI_ERR_ARG;
677   } else {
678     *group = smpi_comm_group(comm);
679     retval = MPI_SUCCESS;
680   }
681   smpi_bench_begin(-1, NULL);
682   return retval;
683 }
684
685 int MPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
686 {
687   int retval;
688
689   smpi_bench_end(-1, NULL);
690   if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
691     retval = MPI_ERR_COMM;
692   } else if (result == NULL) {
693     retval = MPI_ERR_ARG;
694   } else {
695     if (comm1 == comm2) {       /* Same communicators means same groups */
696       *result = MPI_IDENT;
697     } else {
698       *result =
699           smpi_group_compare(smpi_comm_group(comm1),
700                              smpi_comm_group(comm2));
701       if (*result == MPI_IDENT) {
702         *result = MPI_CONGRUENT;
703       }
704     }
705     retval = MPI_SUCCESS;
706   }
707   smpi_bench_begin(-1, NULL);
708   return retval;
709 }
710
711 int MPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
712 {
713   int retval;
714
715   smpi_bench_end(-1, NULL);
716   if (comm == MPI_COMM_NULL) {
717     retval = MPI_ERR_COMM;
718   } else if (newcomm == NULL) {
719     retval = MPI_ERR_ARG;
720   } else {
721     *newcomm = smpi_comm_new(smpi_comm_group(comm));
722     retval = MPI_SUCCESS;
723   }
724   smpi_bench_begin(-1, NULL);
725   return retval;
726 }
727
728 int MPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
729 {
730   int retval;
731
732   smpi_bench_end(-1, NULL);
733   if (comm == MPI_COMM_NULL) {
734     retval = MPI_ERR_COMM;
735   } else if (group == MPI_GROUP_NULL) {
736     retval = MPI_ERR_GROUP;
737   } else if (newcomm == NULL) {
738     retval = MPI_ERR_ARG;
739   } else {
740     *newcomm = smpi_comm_new(group);
741     retval = MPI_SUCCESS;
742   }
743   smpi_bench_begin(-1, NULL);
744   return retval;
745 }
746
747 int MPI_Comm_free(MPI_Comm * comm)
748 {
749   int retval;
750
751   smpi_bench_end(-1, NULL);
752   if (comm == NULL) {
753     retval = MPI_ERR_ARG;
754   } else if (*comm == MPI_COMM_NULL) {
755     retval = MPI_ERR_COMM;
756   } else {
757     smpi_comm_destroy(*comm);
758     *comm = MPI_COMM_NULL;
759     retval = MPI_SUCCESS;
760   }
761   smpi_bench_begin(-1, NULL);
762   return retval;
763 }
764
765 int MPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
766 {
767   int retval;
768
769   smpi_bench_end(-1, NULL);
770   if (comm_out == NULL) {
771     retval = MPI_ERR_ARG;
772   } else if (comm == MPI_COMM_NULL) {
773     retval = MPI_ERR_COMM;
774   } else {
775     printf("[%d] gives %d, %d\n", smpi_comm_rank(comm), color, key);
776     *comm_out = smpi_comm_split(comm, color, key);
777     retval = MPI_SUCCESS;
778   }
779   smpi_bench_begin(-1, NULL);
780   return retval;
781 }
782
783 int MPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
784                   int tag, MPI_Comm comm, MPI_Request * request)
785 {
786   int retval;
787   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
788
789   smpi_bench_end(rank, "Send_init");
790   if (request == NULL) {
791     retval = MPI_ERR_ARG;
792   } else if (comm == MPI_COMM_NULL) {
793     retval = MPI_ERR_COMM;
794   } else {
795     *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
796     retval = MPI_SUCCESS;
797   }
798   smpi_bench_begin(rank, "Send_init");
799   return retval;
800 }
801
802 int MPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
803                   int tag, MPI_Comm comm, MPI_Request * request)
804 {
805   int retval;
806   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
807
808   smpi_bench_end(rank, "Recv_init");
809   if (request == NULL) {
810     retval = MPI_ERR_ARG;
811   } else if (comm == MPI_COMM_NULL) {
812     retval = MPI_ERR_COMM;
813   } else {
814     *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
815     retval = MPI_SUCCESS;
816   }
817   smpi_bench_begin(rank, "Recv_init");
818   return retval;
819 }
820
821 int MPI_Start(MPI_Request * request)
822 {
823   int retval;
824   MPI_Comm comm = (*request)->comm;
825   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
826
827   smpi_bench_end(rank, "Start");
828   if (request == NULL) {
829     retval = MPI_ERR_ARG;
830   } else {
831     smpi_mpi_start(*request);
832     retval = MPI_SUCCESS;
833   }
834   smpi_bench_begin(rank, "Start");
835   return retval;
836 }
837
838 int MPI_Startall(int count, MPI_Request * requests)
839 {
840   int retval;
841   MPI_Comm comm = count > 0
842       && requests ? requests[0]->comm : MPI_COMM_NULL;
843   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
844
845   smpi_bench_end(rank, "Startall");
846   if (requests == NULL) {
847     retval = MPI_ERR_ARG;
848   } else {
849     smpi_mpi_startall(count, requests);
850     retval = MPI_SUCCESS;
851   }
852   smpi_bench_begin(rank, "Startall");
853   return retval;
854 }
855
856 int MPI_Request_free(MPI_Request * request)
857 {
858   int retval;
859   MPI_Comm comm = (*request)->comm;
860   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
861
862   smpi_bench_end(rank, "Request_free");
863   if (request == NULL) {
864     retval = MPI_ERR_ARG;
865   } else {
866     smpi_mpi_request_free(request);
867     retval = MPI_SUCCESS;
868   }
869   smpi_bench_begin(rank, "Request_free");
870   return retval;
871 }
872
873 int MPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
874               int tag, MPI_Comm comm, MPI_Request * request)
875 {
876   int retval;
877   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
878
879   smpi_bench_end(rank, "Irecv");
880 #ifdef HAVE_TRACING
881   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
882   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
883 #endif
884   if (request == NULL) {
885     retval = MPI_ERR_ARG;
886   } else if (comm == MPI_COMM_NULL) {
887     retval = MPI_ERR_COMM;
888   } else {
889     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
890     retval = MPI_SUCCESS;
891   }
892 #ifdef HAVE_TRACING
893   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
894   (*request)->recv = 1;
895 #endif
896   smpi_bench_begin(rank, "Irecv");
897   return retval;
898 }
899
900 int MPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
901               int tag, MPI_Comm comm, MPI_Request * request)
902 {
903   int retval;
904   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
905
906   smpi_bench_end(rank, "Isend");
907 #ifdef HAVE_TRACING
908   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
909   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
910   TRACE_smpi_send(rank, rank, dst_traced);
911 #endif
912   if (request == NULL) {
913     retval = MPI_ERR_ARG;
914   } else if (comm == MPI_COMM_NULL) {
915     retval = MPI_ERR_COMM;
916   } else {
917     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
918     retval = MPI_SUCCESS;
919   }
920 #ifdef HAVE_TRACING
921   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
922   (*request)->send = 1;
923 #endif
924   smpi_bench_begin(rank, "Isend");
925   return retval;
926 }
927
928 int MPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
929              MPI_Comm comm, MPI_Status * status)
930 {
931   int retval;
932   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
933
934   smpi_bench_end(rank, "Recv");
935 #ifdef HAVE_TRACING
936   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
937   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
938 #endif
939   if (comm == MPI_COMM_NULL) {
940     retval = MPI_ERR_COMM;
941   } else {
942     smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
943     retval = MPI_SUCCESS;
944   }
945 #ifdef HAVE_TRACING
946   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
947   TRACE_smpi_recv(rank, src_traced, rank);
948 #endif
949   smpi_bench_begin(rank, "Recv");
950   return retval;
951 }
952
953 int MPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
954              MPI_Comm comm)
955 {
956   int retval;
957   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
958
959   smpi_bench_end(rank, "Send");
960 #ifdef HAVE_TRACING
961   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
962   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
963   TRACE_smpi_send(rank, rank, dst_traced);
964 #endif
965   if (comm == MPI_COMM_NULL) {
966     retval = MPI_ERR_COMM;
967   } else {
968     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
969     retval = MPI_SUCCESS;
970   }
971 #ifdef HAVE_TRACING
972   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
973 #endif
974   smpi_bench_begin(rank, "Send");
975   return retval;
976 }
977
978 int MPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
979                  int dst, int sendtag, void *recvbuf, int recvcount,
980                  MPI_Datatype recvtype, int src, int recvtag,
981                  MPI_Comm comm, MPI_Status * status)
982 {
983   int retval;
984   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
985
986   smpi_bench_end(rank, "Sendrecv");
987 #ifdef HAVE_TRACING
988   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
989   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
990   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
991   TRACE_smpi_send(rank, rank, dst_traced);
992   TRACE_smpi_send(rank, src_traced, rank);
993 #endif
994   if (comm == MPI_COMM_NULL) {
995     retval = MPI_ERR_COMM;
996   } else if (sendtype == MPI_DATATYPE_NULL
997              || recvtype == MPI_DATATYPE_NULL) {
998     retval = MPI_ERR_TYPE;
999   } else {
1000     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1001                       recvcount, recvtype, src, recvtag, comm, status);
1002     retval = MPI_SUCCESS;
1003   }
1004 #ifdef HAVE_TRACING
1005   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1006   TRACE_smpi_recv(rank, rank, dst_traced);
1007   TRACE_smpi_recv(rank, src_traced, rank);
1008 #endif
1009   smpi_bench_begin(rank, "Sendrecv");
1010   return retval;
1011 }
1012
1013 int MPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1014                          int dst, int sendtag, int src, int recvtag,
1015                          MPI_Comm comm, MPI_Status * status)
1016 {
1017   //TODO: suboptimal implementation
1018   void *recvbuf;
1019   int retval, size;
1020
1021   size = smpi_datatype_size(datatype) * count;
1022   recvbuf = xbt_new(char, size);
1023   retval =
1024       MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1025                    datatype, src, recvtag, comm, status);
1026   memcpy(buf, recvbuf, size * sizeof(char));
1027   xbt_free(recvbuf);
1028   return retval;
1029 }
1030
1031 int MPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1032 {
1033   int retval;
1034   int rank = request && (*request)->comm != MPI_COMM_NULL
1035       ? smpi_comm_rank((*request)->comm)
1036       : -1;
1037
1038   smpi_bench_end(rank, "Test");
1039   if (request == NULL || flag == NULL) {
1040     retval = MPI_ERR_ARG;
1041   } else if (*request == MPI_REQUEST_NULL) {
1042     retval = MPI_ERR_REQUEST;
1043   } else {
1044     *flag = smpi_mpi_test(request, status);
1045     retval = MPI_SUCCESS;
1046   }
1047   smpi_bench_begin(rank, "Test");
1048   return retval;
1049 }
1050
1051 int MPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1052                 MPI_Status * status)
1053 {
1054   int retval;
1055
1056   smpi_bench_end(-1, NULL);     //FIXME
1057   if (index == NULL || flag == NULL) {
1058     retval = MPI_ERR_ARG;
1059   } else {
1060     *flag = smpi_mpi_testany(count, requests, index, status);
1061     retval = MPI_SUCCESS;
1062   }
1063   smpi_bench_begin(-1, NULL);
1064   return retval;
1065 }
1066
1067 int MPI_Wait(MPI_Request * request, MPI_Status * status)
1068 {
1069   int retval;
1070   int rank = request && (*request)->comm != MPI_COMM_NULL
1071       ? smpi_comm_rank((*request)->comm)
1072       : -1;
1073
1074   smpi_bench_end(rank, "Wait");
1075 #ifdef HAVE_TRACING
1076   MPI_Group group = smpi_comm_group((*request)->comm);
1077   int src_traced = smpi_group_rank(group, (*request)->src);
1078   int dst_traced = smpi_group_rank(group, (*request)->dst);
1079   int is_wait_for_receive = (*request)->recv;
1080   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1081 #endif
1082   if (request == NULL) {
1083     retval = MPI_ERR_ARG;
1084   } else if (*request == MPI_REQUEST_NULL) {
1085     retval = MPI_ERR_REQUEST;
1086   } else {
1087     smpi_mpi_wait(request, status);
1088     retval = MPI_SUCCESS;
1089   }
1090 #ifdef HAVE_TRACING
1091   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1092   if (is_wait_for_receive) {
1093     TRACE_smpi_recv(rank, src_traced, dst_traced);
1094   }
1095 #endif
1096   smpi_bench_begin(rank, "Wait");
1097   return retval;
1098 }
1099
1100 int MPI_Waitany(int count, MPI_Request requests[], int *index,
1101                 MPI_Status * status)
1102 {
1103   int retval;
1104
1105   smpi_bench_end(-1, NULL);     //FIXME
1106 #ifdef HAVE_TRACING
1107   //save requests information for tracing
1108   int i;
1109   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
1110   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
1111   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
1112   for (i = 0; i < count; i++) {
1113     MPI_Request req = requests[i];      //already received requests are no longer valid
1114     if (req) {
1115       int *asrc = xbt_new(int, 1);
1116       int *adst = xbt_new(int, 1);
1117       int *arecv = xbt_new(int, 1);
1118       *asrc = req->src;
1119       *adst = req->dst;
1120       *arecv = req->recv;
1121       xbt_dynar_insert_at(srcs, i, asrc);
1122       xbt_dynar_insert_at(dsts, i, adst);
1123       xbt_dynar_insert_at(recvs, i, arecv);
1124     } else {
1125       int *t = xbt_new(int, 1);
1126       xbt_dynar_insert_at(srcs, i, t);
1127       xbt_dynar_insert_at(dsts, i, t);
1128       xbt_dynar_insert_at(recvs, i, t);
1129     }
1130   }
1131   int rank_traced = smpi_comm_rank(MPI_COMM_WORLD);
1132   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1133 #endif
1134   if (index == NULL) {
1135     retval = MPI_ERR_ARG;
1136   } else {
1137     *index = smpi_mpi_waitany(count, requests, status);
1138     retval = MPI_SUCCESS;
1139   }
1140 #ifdef HAVE_TRACING
1141   int src_traced, dst_traced, is_wait_for_receive;
1142   xbt_dynar_get_cpy(srcs, *index, &src_traced);
1143   xbt_dynar_get_cpy(dsts, *index, &dst_traced);
1144   xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
1145   if (is_wait_for_receive) {
1146     TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1147   }
1148   TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1149   //clean-up of dynars
1150   xbt_free(srcs);
1151   xbt_free(dsts);
1152   xbt_free(recvs);
1153 #endif
1154   smpi_bench_begin(-1, NULL);
1155   return retval;
1156 }
1157
1158 int MPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1159 {
1160
1161   smpi_bench_end(-1, NULL);     //FIXME
1162 #ifdef HAVE_TRACING
1163   //save information from requests
1164   int i;
1165   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
1166   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
1167   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
1168   for (i = 0; i < count; i++) {
1169     MPI_Request req = requests[i];      //all req should be valid in Waitall
1170     int *asrc = xbt_new(int, 1);
1171     int *adst = xbt_new(int, 1);
1172     int *arecv = xbt_new(int, 1);
1173     *asrc = req->src;
1174     *adst = req->dst;
1175     *arecv = req->recv;
1176     xbt_dynar_insert_at(srcs, i, asrc);
1177     xbt_dynar_insert_at(dsts, i, adst);
1178     xbt_dynar_insert_at(recvs, i, arecv);
1179   }
1180   int rank_traced = smpi_comm_rank (MPI_COMM_WORLD);
1181   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1182 #endif
1183   smpi_mpi_waitall(count, requests, status);
1184 #ifdef HAVE_TRACING
1185   for (i = 0; i < count; i++) {
1186     int src_traced, dst_traced, is_wait_for_receive;
1187     xbt_dynar_get_cpy(srcs, i, &src_traced);
1188     xbt_dynar_get_cpy(dsts, i, &dst_traced);
1189     xbt_dynar_get_cpy(recvs, i, &is_wait_for_receive);
1190     if (is_wait_for_receive) {
1191       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1192     }
1193   }
1194   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1195   //clean-up of dynars
1196   xbt_free(srcs);
1197   xbt_free(dsts);
1198   xbt_free(recvs);
1199 #endif
1200   smpi_bench_begin(-1, NULL);
1201   return MPI_SUCCESS;
1202 }
1203
1204 int MPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1205                  int *indices, MPI_Status status[])
1206 {
1207   int retval;
1208
1209   smpi_bench_end(-1, NULL);     //FIXME
1210   if (outcount == NULL || indices == NULL) {
1211     retval = MPI_ERR_ARG;
1212   } else {
1213     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1214     retval = MPI_SUCCESS;
1215   }
1216   smpi_bench_begin(-1, NULL);
1217   return retval;
1218 }
1219
1220 int MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
1221               MPI_Comm comm)
1222 {
1223   int retval;
1224   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1225
1226   smpi_bench_end(rank, "Bcast");
1227 #ifdef HAVE_TRACING
1228   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1229   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1230 #endif
1231   if (comm == MPI_COMM_NULL) {
1232     retval = MPI_ERR_COMM;
1233   } else {
1234     smpi_mpi_bcast(buf, count, datatype, root, comm);
1235     retval = MPI_SUCCESS;
1236   }
1237 #ifdef HAVE_TRACING
1238   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1239 #endif
1240   smpi_bench_begin(rank, "Bcast");
1241   return retval;
1242 }
1243
1244 int MPI_Barrier(MPI_Comm comm)
1245 {
1246   int retval;
1247   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1248
1249   smpi_bench_end(rank, "Barrier");
1250 #ifdef HAVE_TRACING
1251   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1252 #endif
1253   if (comm == MPI_COMM_NULL) {
1254     retval = MPI_ERR_COMM;
1255   } else {
1256     smpi_mpi_barrier(comm);
1257     retval = MPI_SUCCESS;
1258   }
1259 #ifdef HAVE_TRACING
1260   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1261 #endif
1262   smpi_bench_begin(rank, "Barrier");
1263   return retval;
1264 }
1265
1266 int MPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1267                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1268                int root, MPI_Comm comm)
1269 {
1270   int retval;
1271   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1272
1273   smpi_bench_end(rank, "Gather");
1274 #ifdef HAVE_TRACING
1275   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1276   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1277 #endif
1278   if (comm == MPI_COMM_NULL) {
1279     retval = MPI_ERR_COMM;
1280   } else if (sendtype == MPI_DATATYPE_NULL
1281              || recvtype == MPI_DATATYPE_NULL) {
1282     retval = MPI_ERR_TYPE;
1283   } else {
1284     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1285                     recvtype, root, comm);
1286     retval = MPI_SUCCESS;
1287   }
1288 #ifdef HAVE_TRACING
1289   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1290 #endif
1291   smpi_bench_begin(rank, "Gather");
1292   return retval;
1293 }
1294
1295 int MPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1296                 void *recvbuf, int *recvcounts, int *displs,
1297                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1298 {
1299   int retval;
1300   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1301
1302   smpi_bench_end(rank, "Gatherv");
1303 #ifdef HAVE_TRACING
1304   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1305   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1306 #endif
1307   if (comm == MPI_COMM_NULL) {
1308     retval = MPI_ERR_COMM;
1309   } else if (sendtype == MPI_DATATYPE_NULL
1310              || recvtype == MPI_DATATYPE_NULL) {
1311     retval = MPI_ERR_TYPE;
1312   } else if (recvcounts == NULL || displs == NULL) {
1313     retval = MPI_ERR_ARG;
1314   } else {
1315     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1316                      displs, recvtype, root, comm);
1317     retval = MPI_SUCCESS;
1318   }
1319 #ifdef HAVE_TRACING
1320   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1321 #endif
1322   smpi_bench_begin(rank, "Gatherv");
1323   return retval;
1324 }
1325
1326 int MPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1327                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1328                   MPI_Comm comm)
1329 {
1330   int retval;
1331   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1332
1333   smpi_bench_end(rank, "Allgather");
1334 #ifdef HAVE_TRACING
1335   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1336 #endif
1337   if (comm == MPI_COMM_NULL) {
1338     retval = MPI_ERR_COMM;
1339   } else if (sendtype == MPI_DATATYPE_NULL
1340              || recvtype == MPI_DATATYPE_NULL) {
1341     retval = MPI_ERR_TYPE;
1342   } else {
1343     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1344                        recvtype, comm);
1345     retval = MPI_SUCCESS;
1346   }
1347 #ifdef HAVE_TRACING
1348   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1349 #endif
1350   smpi_bench_begin(rank, "Allgather");
1351   return retval;
1352 }
1353
1354 int MPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1355                    void *recvbuf, int *recvcounts, int *displs,
1356                    MPI_Datatype recvtype, MPI_Comm comm)
1357 {
1358   int retval;
1359   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1360
1361   smpi_bench_end(rank, "Allgatherv");
1362 #ifdef HAVE_TRACING
1363   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1364 #endif
1365   if (comm == MPI_COMM_NULL) {
1366     retval = MPI_ERR_COMM;
1367   } else if (sendtype == MPI_DATATYPE_NULL
1368              || recvtype == MPI_DATATYPE_NULL) {
1369     retval = MPI_ERR_TYPE;
1370   } else if (recvcounts == NULL || displs == NULL) {
1371     retval = MPI_ERR_ARG;
1372   } else {
1373     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1374                         displs, recvtype, comm);
1375     retval = MPI_SUCCESS;
1376   }
1377 #ifdef HAVE_TRACING
1378   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1379 #endif
1380   smpi_bench_begin(rank, "Allgatherv");
1381   return retval;
1382 }
1383
1384 int MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1385                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1386                 int root, MPI_Comm comm)
1387 {
1388   int retval;
1389   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1390
1391   smpi_bench_end(rank, "Scatter");
1392 #ifdef HAVE_TRACING
1393   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1394   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1395 #endif
1396   if (comm == MPI_COMM_NULL) {
1397     retval = MPI_ERR_COMM;
1398   } else if (sendtype == MPI_DATATYPE_NULL
1399              || recvtype == MPI_DATATYPE_NULL) {
1400     retval = MPI_ERR_TYPE;
1401   } else {
1402     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1403                      recvtype, root, comm);
1404     retval = MPI_SUCCESS;
1405   }
1406 #ifdef HAVE_TRACING
1407   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1408 #endif
1409   smpi_bench_begin(rank, "Scatter");
1410   return retval;
1411 }
1412
1413 int MPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1414                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1415                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1416 {
1417   int retval;
1418   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1419
1420   smpi_bench_end(rank, "Scatterv");
1421 #ifdef HAVE_TRACING
1422   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1423   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1424 #endif
1425   if (comm == MPI_COMM_NULL) {
1426     retval = MPI_ERR_COMM;
1427   } else if (sendtype == MPI_DATATYPE_NULL
1428              || recvtype == MPI_DATATYPE_NULL) {
1429     retval = MPI_ERR_TYPE;
1430   } else if (sendcounts == NULL || displs == NULL) {
1431     retval = MPI_ERR_ARG;
1432   } else {
1433     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1434                       recvcount, recvtype, root, comm);
1435     retval = MPI_SUCCESS;
1436   }
1437 #ifdef HAVE_TRACING
1438   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1439 #endif
1440   smpi_bench_begin(rank, "Scatterv");
1441   return retval;
1442 }
1443
1444 int MPI_Reduce(void *sendbuf, void *recvbuf, int count,
1445                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1446 {
1447   int retval;
1448   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1449
1450   smpi_bench_end(rank, "Reduce");
1451 #ifdef HAVE_TRACING
1452   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1453   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1454 #endif
1455   if (comm == MPI_COMM_NULL) {
1456     retval = MPI_ERR_COMM;
1457   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1458     retval = MPI_ERR_ARG;
1459   } else {
1460     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1461     retval = MPI_SUCCESS;
1462   }
1463 #ifdef HAVE_TRACING
1464   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1465 #endif
1466   smpi_bench_begin(rank, "Reduce");
1467   return retval;
1468 }
1469
1470 int MPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1471                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1472 {
1473   int retval;
1474   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1475
1476   smpi_bench_end(rank, "Allreduce");
1477 #ifdef HAVE_TRACING
1478   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1479 #endif
1480   if (comm == MPI_COMM_NULL) {
1481     retval = MPI_ERR_COMM;
1482   } else if (datatype == MPI_DATATYPE_NULL) {
1483     retval = MPI_ERR_TYPE;
1484   } else if (op == MPI_OP_NULL) {
1485     retval = MPI_ERR_OP;
1486   } else {
1487     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
1488     retval = MPI_SUCCESS;
1489   }
1490 #ifdef HAVE_TRACING
1491   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1492 #endif
1493   smpi_bench_begin(rank, "Allreduce");
1494   return retval;
1495 }
1496
1497 int MPI_Scan(void *sendbuf, void *recvbuf, int count,
1498              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1499 {
1500   int retval;
1501   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1502
1503   smpi_bench_end(rank, "Scan");
1504 #ifdef HAVE_TRACING
1505   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1506 #endif
1507   if (comm == MPI_COMM_NULL) {
1508     retval = MPI_ERR_COMM;
1509   } else if (datatype == MPI_DATATYPE_NULL) {
1510     retval = MPI_ERR_TYPE;
1511   } else if (op == MPI_OP_NULL) {
1512     retval = MPI_ERR_OP;
1513   } else {
1514     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1515     retval = MPI_SUCCESS;
1516   }
1517 #ifdef HAVE_TRACING
1518   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1519 #endif
1520   smpi_bench_begin(rank, "Scan");
1521   return retval;
1522 }
1523
1524 int MPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1525                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1526 {
1527   int retval, i, size, count;
1528   int *displs;
1529   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1530
1531   smpi_bench_end(rank, "Reduce_scatter");
1532 #ifdef HAVE_TRACING
1533   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1534 #endif
1535   if (comm == MPI_COMM_NULL) {
1536     retval = MPI_ERR_COMM;
1537   } else if (datatype == MPI_DATATYPE_NULL) {
1538     retval = MPI_ERR_TYPE;
1539   } else if (op == MPI_OP_NULL) {
1540     retval = MPI_ERR_OP;
1541   } else if (recvcounts == NULL) {
1542     retval = MPI_ERR_ARG;
1543   } else {
1544     /* arbitrarily choose root as rank 0 */
1545     /* TODO: faster direct implementation ? */
1546     size = smpi_comm_size(comm);
1547     count = 0;
1548     displs = xbt_new(int, size);
1549     for (i = 0; i < size; i++) {
1550       count += recvcounts[i];
1551       displs[i] = 0;
1552     }
1553     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1554     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf,
1555                       recvcounts[rank], datatype, 0, comm);
1556     xbt_free(displs);
1557     retval = MPI_SUCCESS;
1558   }
1559 #ifdef HAVE_TRACING
1560   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1561 #endif
1562   smpi_bench_begin(rank, "Reduce_scatter");
1563   return retval;
1564 }
1565
1566 int MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1567                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
1568                  MPI_Comm comm)
1569 {
1570   int retval, size, sendsize;
1571   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1572
1573   smpi_bench_end(rank, "Alltoall");
1574 #ifdef HAVE_TRACING
1575   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1576 #endif
1577   if (comm == MPI_COMM_NULL) {
1578     retval = MPI_ERR_COMM;
1579   } else if (sendtype == MPI_DATATYPE_NULL
1580              || recvtype == MPI_DATATYPE_NULL) {
1581     retval = MPI_ERR_TYPE;
1582   } else {
1583     size = smpi_comm_size(comm);
1584     sendsize = smpi_datatype_size(sendtype) * sendcount;
1585     if (sendsize < 200 && size > 12) {
1586       retval =
1587           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
1588                                          recvbuf, recvcount, recvtype,
1589                                          comm);
1590     } else if (sendsize < 3000) {
1591       retval =
1592           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
1593                                                 sendtype, recvbuf,
1594                                                 recvcount, recvtype, comm);
1595     } else {
1596       retval =
1597           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
1598                                             recvbuf, recvcount, recvtype,
1599                                             comm);
1600     }
1601   }
1602 #ifdef HAVE_TRACING
1603   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1604 #endif
1605   smpi_bench_begin(rank, "Alltoall");
1606   return retval;
1607 }
1608
1609 int MPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
1610                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
1611                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
1612 {
1613   int retval;
1614   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1615
1616   smpi_bench_end(rank, "Alltoallv");
1617 #ifdef HAVE_TRACING
1618   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1619 #endif
1620   if (comm == MPI_COMM_NULL) {
1621     retval = MPI_ERR_COMM;
1622   } else if (sendtype == MPI_DATATYPE_NULL
1623              || recvtype == MPI_DATATYPE_NULL) {
1624     retval = MPI_ERR_TYPE;
1625   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
1626              || recvdisps == NULL) {
1627     retval = MPI_ERR_ARG;
1628   } else {
1629     retval =
1630         smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
1631                                   recvbuf, recvcounts, recvdisps, recvtype,
1632                                   comm);
1633   }
1634 #ifdef HAVE_TRACING
1635   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1636 #endif
1637   smpi_bench_begin(rank, "Alltoallv");
1638   return retval;
1639 }
1640
1641
1642 int MPI_Get_processor_name(char *name, int *resultlen)
1643 {
1644   int retval = MPI_SUCCESS;
1645
1646   smpi_bench_end(-1, NULL);
1647   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
1648           MPI_MAX_PROCESSOR_NAME - 1);
1649   *resultlen =
1650       strlen(name) >
1651       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1652
1653   smpi_bench_begin(-1, NULL);
1654   return retval;
1655 }
1656
1657 int MPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
1658 {
1659   int retval = MPI_SUCCESS;
1660   size_t size;
1661
1662   smpi_bench_end(-1, NULL);
1663   if (status == NULL || count == NULL) {
1664     retval = MPI_ERR_ARG;
1665   } else if (datatype == MPI_DATATYPE_NULL) {
1666     retval = MPI_ERR_TYPE;
1667   } else {
1668     size = smpi_datatype_size(datatype);
1669     if (size == 0) {
1670       *count = 0;
1671     } else if (status->count % size != 0) {
1672       retval = MPI_UNDEFINED;
1673     } else {
1674       *count = smpi_mpi_get_count(status, datatype);
1675     }
1676   }
1677   smpi_bench_begin(-1, NULL);
1678   return retval;
1679 }