Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
model-checker : remove unused argument in functions for heap comparison algorithm
[simgrid.git] / src / smpi / smpi_pmpi.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_mpi_dt_private.h"
9
10 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_pmpi, smpi,
11                                 "Logging specific to SMPI (pmpi)");
12
13 #ifdef HAVE_TRACING
14 //this function need to be here because of the calls to smpi_bench
15 void TRACE_smpi_set_category(const char *category)
16 {
17   //need to end bench otherwise categories for execution tasks are wrong
18   smpi_bench_end();
19   TRACE_internal_smpi_set_category (category);
20   //begin bench after changing process's category
21   smpi_bench_begin();
22 }
23 #endif
24
25 /* PMPI User level calls */
26
27 int PMPI_Init(int *argc, char ***argv)
28 {
29   smpi_process_init(argc, argv);
30 #ifdef HAVE_TRACING
31   TRACE_smpi_init(smpi_process_index());
32 #endif
33   smpi_bench_begin();
34   return MPI_SUCCESS;
35 }
36
37 int PMPI_Finalize(void)
38 {
39   smpi_process_finalize();
40   smpi_bench_end();
41 #ifdef HAVE_TRACING
42   TRACE_smpi_finalize(smpi_process_index());
43 #endif
44   smpi_process_destroy();
45   return MPI_SUCCESS;
46 }
47
48 int PMPI_Init_thread(int *argc, char ***argv, int required, int *provided)
49 {
50   if (provided != NULL) {
51     *provided = MPI_THREAD_MULTIPLE;
52   }
53   return MPI_Init(argc, argv);
54 }
55
56 int PMPI_Query_thread(int *provided)
57 {
58   int retval;
59
60   smpi_bench_end();
61   if (provided == NULL) {
62     retval = MPI_ERR_ARG;
63   } else {
64     *provided = MPI_THREAD_MULTIPLE;
65     retval = MPI_SUCCESS;
66   }
67   smpi_bench_begin();
68   return retval;
69 }
70
71 int PMPI_Is_thread_main(int *flag)
72 {
73   int retval;
74
75   smpi_bench_end();
76   if (flag == NULL) {
77     retval = MPI_ERR_ARG;
78   } else {
79     *flag = smpi_process_index() == 0;
80     retval = MPI_SUCCESS;
81   }
82   smpi_bench_begin();
83   return retval;
84 }
85
86 int PMPI_Abort(MPI_Comm comm, int errorcode)
87 {
88   smpi_bench_end();
89   smpi_process_destroy();
90   // FIXME: should kill all processes in comm instead
91   simcall_process_kill(SIMIX_process_self());
92   return MPI_SUCCESS;
93 }
94
95 double PMPI_Wtime(void)
96 {
97   double time;
98
99   smpi_bench_end();
100   time = SIMIX_get_clock();
101   smpi_bench_begin();
102   return time;
103 }
104
105 int PMPI_Address(void *location, MPI_Aint * address)
106 {
107   int retval;
108
109   smpi_bench_end();
110   if (!address) {
111     retval = MPI_ERR_ARG;
112   } else {
113     *address = (MPI_Aint) location;
114   }
115   smpi_bench_begin();
116   return retval;
117 }
118
119 int PMPI_Type_free(MPI_Datatype * datatype)
120 {
121   int retval;
122
123   smpi_bench_end();
124   if (!datatype) {
125     retval = MPI_ERR_ARG;
126   } else {
127     // FIXME: always fail for now
128     retval = MPI_ERR_TYPE;
129   }
130   smpi_bench_begin();
131   return retval;
132 }
133
134 int PMPI_Type_size(MPI_Datatype datatype, int *size)
135 {
136   int retval;
137
138   smpi_bench_end();
139   if (datatype == MPI_DATATYPE_NULL) {
140     retval = MPI_ERR_TYPE;
141   } else if (size == NULL) {
142     retval = MPI_ERR_ARG;
143   } else {
144     *size = (int) smpi_datatype_size(datatype);
145     retval = MPI_SUCCESS;
146   }
147   smpi_bench_begin();
148   return retval;
149 }
150
151 int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
152 {
153   int retval;
154
155   smpi_bench_end();
156   if (datatype == MPI_DATATYPE_NULL) {
157     retval = MPI_ERR_TYPE;
158   } else if (lb == NULL || extent == NULL) {
159     retval = MPI_ERR_ARG;
160   } else {
161     retval = smpi_datatype_extent(datatype, lb, extent);
162   }
163   smpi_bench_begin();
164   return retval;
165 }
166
167 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
168 {
169   int retval;
170   MPI_Aint dummy;
171
172   smpi_bench_end();
173   if (datatype == MPI_DATATYPE_NULL) {
174     retval = MPI_ERR_TYPE;
175   } else if (extent == NULL) {
176     retval = MPI_ERR_ARG;
177   } else {
178     retval = smpi_datatype_extent(datatype, &dummy, extent);
179   }
180   smpi_bench_begin();
181   return retval;
182 }
183
184 int PMPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
185 {
186   int retval;
187
188   smpi_bench_end();
189   if (datatype == MPI_DATATYPE_NULL) {
190     retval = MPI_ERR_TYPE;
191   } else if (disp == NULL) {
192     retval = MPI_ERR_ARG;
193   } else {
194     *disp = smpi_datatype_lb(datatype);
195     retval = MPI_SUCCESS;
196   }
197   smpi_bench_begin();
198   return retval;
199 }
200
201 int PMPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
202 {
203   int retval;
204
205   smpi_bench_end();
206   if (datatype == MPI_DATATYPE_NULL) {
207     retval = MPI_ERR_TYPE;
208   } else if (disp == NULL) {
209     retval = MPI_ERR_ARG;
210   } else {
211     *disp = smpi_datatype_ub(datatype);
212     retval = MPI_SUCCESS;
213   }
214   smpi_bench_begin();
215   return retval;
216 }
217
218 int PMPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
219 {
220   int retval;
221
222   smpi_bench_end();
223   if (function == NULL || op == NULL) {
224     retval = MPI_ERR_ARG;
225   } else {
226     *op = smpi_op_new(function, commute);
227     retval = MPI_SUCCESS;
228   }
229   smpi_bench_begin();
230   return retval;
231 }
232
233 int PMPI_Op_free(MPI_Op * op)
234 {
235   int retval;
236
237   smpi_bench_end();
238   if (op == NULL) {
239     retval = MPI_ERR_ARG;
240   } else if (*op == MPI_OP_NULL) {
241     retval = MPI_ERR_OP;
242   } else {
243     smpi_op_destroy(*op);
244     *op = MPI_OP_NULL;
245     retval = MPI_SUCCESS;
246   }
247   smpi_bench_begin();
248   return retval;
249 }
250
251 int PMPI_Group_free(MPI_Group * group)
252 {
253   int retval;
254
255   smpi_bench_end();
256   if (group == NULL) {
257     retval = MPI_ERR_ARG;
258   } else {
259     smpi_group_destroy(*group);
260     *group = MPI_GROUP_NULL;
261     retval = MPI_SUCCESS;
262   }
263   smpi_bench_begin();
264   return retval;
265 }
266
267 int PMPI_Group_size(MPI_Group group, int *size)
268 {
269   int retval;
270
271   smpi_bench_end();
272   if (group == MPI_GROUP_NULL) {
273     retval = MPI_ERR_GROUP;
274   } else if (size == NULL) {
275     retval = MPI_ERR_ARG;
276   } else {
277     *size = smpi_group_size(group);
278     retval = MPI_SUCCESS;
279   }
280   smpi_bench_begin();
281   return retval;
282 }
283
284 int PMPI_Group_rank(MPI_Group group, int *rank)
285 {
286   int retval;
287
288   smpi_bench_end();
289   if (group == MPI_GROUP_NULL) {
290     retval = MPI_ERR_GROUP;
291   } else if (rank == NULL) {
292     retval = MPI_ERR_ARG;
293   } else {
294     *rank = smpi_group_rank(group, smpi_process_index());
295     retval = MPI_SUCCESS;
296   }
297   smpi_bench_begin();
298   return retval;
299 }
300
301 int PMPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1,
302                               MPI_Group group2, int *ranks2)
303 {
304   int retval, i, index;
305
306   smpi_bench_end();
307   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
308     retval = MPI_ERR_GROUP;
309   } else {
310     for (i = 0; i < n; i++) {
311       index = smpi_group_index(group1, ranks1[i]);
312       ranks2[i] = smpi_group_rank(group2, index);
313     }
314     retval = MPI_SUCCESS;
315   }
316   smpi_bench_begin();
317   return retval;
318 }
319
320 int PMPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
321 {
322   int retval;
323
324   smpi_bench_end();
325   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
326     retval = MPI_ERR_GROUP;
327   } else if (result == NULL) {
328     retval = MPI_ERR_ARG;
329   } else {
330     *result = smpi_group_compare(group1, group2);
331     retval = MPI_SUCCESS;
332   }
333   smpi_bench_begin();
334   return retval;
335 }
336
337 int PMPI_Group_union(MPI_Group group1, MPI_Group group2,
338                     MPI_Group * newgroup)
339 {
340   int retval, i, proc1, proc2, size, size2;
341
342   smpi_bench_end();
343   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
344     retval = MPI_ERR_GROUP;
345   } else if (newgroup == NULL) {
346     retval = MPI_ERR_ARG;
347   } else {
348     size = smpi_group_size(group1);
349     size2 = smpi_group_size(group2);
350     for (i = 0; i < size2; i++) {
351       proc2 = smpi_group_index(group2, i);
352       proc1 = smpi_group_rank(group1, proc2);
353       if (proc1 == MPI_UNDEFINED) {
354         size++;
355       }
356     }
357     if (size == 0) {
358       *newgroup = MPI_GROUP_EMPTY;
359     } else {
360       *newgroup = smpi_group_new(size);
361       size2 = smpi_group_size(group1);
362       for (i = 0; i < size2; i++) {
363         proc1 = smpi_group_index(group1, i);
364         smpi_group_set_mapping(*newgroup, proc1, i);
365       }
366       for (i = size2; i < size; i++) {
367         proc2 = smpi_group_index(group2, i - size2);
368         smpi_group_set_mapping(*newgroup, proc2, i);
369       }
370     }
371     smpi_group_use(*newgroup);
372     retval = MPI_SUCCESS;
373   }
374   smpi_bench_begin();
375   return retval;
376 }
377
378 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
379                            MPI_Group * newgroup)
380 {
381   int retval, i, proc1, proc2, size, size2;
382
383   smpi_bench_end();
384   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
385     retval = MPI_ERR_GROUP;
386   } else if (newgroup == NULL) {
387     retval = MPI_ERR_ARG;
388   } else {
389     size = smpi_group_size(group1);
390     size2 = smpi_group_size(group2);
391     for (i = 0; i < size2; i++) {
392       proc2 = smpi_group_index(group2, i);
393       proc1 = smpi_group_rank(group1, proc2);
394       if (proc1 == MPI_UNDEFINED) {
395         size--;
396       }
397     }
398     if (size == 0) {
399       *newgroup = MPI_GROUP_EMPTY;
400     } else {
401       *newgroup = smpi_group_new(size);
402       size2 = smpi_group_size(group1);
403       for (i = 0; i < size2; i++) {
404         proc1 = smpi_group_index(group1, i);
405         proc2 = smpi_group_rank(group2, proc1);
406         if (proc2 != MPI_UNDEFINED) {
407           smpi_group_set_mapping(*newgroup, proc1, i);
408         }
409       }
410     }
411     smpi_group_use(*newgroup);
412     retval = MPI_SUCCESS;
413   }
414   smpi_bench_begin();
415   return retval;
416 }
417
418 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
419 {
420   int retval, i, proc1, proc2, size, size2;
421
422   smpi_bench_end();
423   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
424     retval = MPI_ERR_GROUP;
425   } else if (newgroup == NULL) {
426     retval = MPI_ERR_ARG;
427   } else {
428     size = size2 = smpi_group_size(group1);
429     for (i = 0; i < size2; i++) {
430       proc1 = smpi_group_index(group1, i);
431       proc2 = smpi_group_rank(group2, proc1);
432       if (proc2 != MPI_UNDEFINED) {
433         size--;
434       }
435     }
436     if (size == 0) {
437       *newgroup = MPI_GROUP_EMPTY;
438     } else {
439       *newgroup = smpi_group_new(size);
440       for (i = 0; i < size2; i++) {
441         proc1 = smpi_group_index(group1, i);
442         proc2 = smpi_group_rank(group2, proc1);
443         if (proc2 == MPI_UNDEFINED) {
444           smpi_group_set_mapping(*newgroup, proc1, i);
445         }
446       }
447     }
448     smpi_group_use(*newgroup);
449     retval = MPI_SUCCESS;
450   }
451   smpi_bench_begin();
452   return retval;
453 }
454
455 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
456 {
457   int retval, i, index;
458
459   smpi_bench_end();
460   if (group == MPI_GROUP_NULL) {
461     retval = MPI_ERR_GROUP;
462   } else if (newgroup == NULL) {
463     retval = MPI_ERR_ARG;
464   } else {
465     if (n == 0) {
466       *newgroup = MPI_GROUP_EMPTY;
467     } else if (n == smpi_group_size(group)) {
468       *newgroup = group;
469     } else {
470       *newgroup = smpi_group_new(n);
471       for (i = 0; i < n; i++) {
472         index = smpi_group_index(group, ranks[i]);
473         smpi_group_set_mapping(*newgroup, index, i);
474       }
475     }
476     smpi_group_use(*newgroup);
477     retval = MPI_SUCCESS;
478   }
479   smpi_bench_begin();
480   return retval;
481 }
482
483 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
484 {
485   int retval, i, size, rank, index;
486
487   smpi_bench_end();
488   if (group == MPI_GROUP_NULL) {
489     retval = MPI_ERR_GROUP;
490   } else if (newgroup == NULL) {
491     retval = MPI_ERR_ARG;
492   } else {
493     if (n == 0) {
494       *newgroup = group;
495     } else if (n == smpi_group_size(group)) {
496       *newgroup = MPI_GROUP_EMPTY;
497     } else {
498       size = smpi_group_size(group) - n;
499       *newgroup = smpi_group_new(size);
500       rank = 0;
501       while (rank < size) {
502         for (i = 0; i < n; i++) {
503           if (ranks[i] == rank) {
504             break;
505           }
506         }
507         if (i >= n) {
508           index = smpi_group_index(group, rank);
509           smpi_group_set_mapping(*newgroup, index, rank);
510           rank++;
511         }
512       }
513     }
514     smpi_group_use(*newgroup);
515     retval = MPI_SUCCESS;
516   }
517   smpi_bench_begin();
518   return retval;
519 }
520
521 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
522                          MPI_Group * newgroup)
523 {
524   int retval, i, j, rank, size, index;
525
526   smpi_bench_end();
527   if (group == MPI_GROUP_NULL) {
528     retval = MPI_ERR_GROUP;
529   } else if (newgroup == NULL) {
530     retval = MPI_ERR_ARG;
531   } else {
532     if (n == 0) {
533       *newgroup = MPI_GROUP_EMPTY;
534     } else {
535       size = 0;
536       for (i = 0; i < n; i++) {
537         for (rank = ranges[i][0];       /* First */
538              rank >= 0 && rank <= ranges[i][1]; /* Last */
539              rank += ranges[i][2] /* Stride */ ) {
540           size++;
541         }
542       }
543       if (size == smpi_group_size(group)) {
544         *newgroup = group;
545       } else {
546         *newgroup = smpi_group_new(size);
547         j = 0;
548         for (i = 0; i < n; i++) {
549           for (rank = ranges[i][0];     /* First */
550                rank >= 0 && rank <= ranges[i][1];       /* Last */
551                rank += ranges[i][2] /* Stride */ ) {
552             index = smpi_group_index(group, rank);
553             smpi_group_set_mapping(*newgroup, index, j);
554             j++;
555           }
556         }
557       }
558     }
559     smpi_group_use(*newgroup);
560     retval = MPI_SUCCESS;
561   }
562   smpi_bench_begin();
563   return retval;
564 }
565
566 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
567                          MPI_Group * newgroup)
568 {
569   int retval, i, newrank, rank, size, index, add;
570
571   smpi_bench_end();
572   if (group == MPI_GROUP_NULL) {
573     retval = MPI_ERR_GROUP;
574   } else if (newgroup == NULL) {
575     retval = MPI_ERR_ARG;
576   } else {
577     if (n == 0) {
578       *newgroup = group;
579     } else {
580       size = smpi_group_size(group);
581       for (i = 0; i < n; i++) {
582         for (rank = ranges[i][0];       /* First */
583              rank >= 0 && rank <= ranges[i][1]; /* Last */
584              rank += ranges[i][2] /* Stride */ ) {
585           size--;
586         }
587       }
588       if (size == 0) {
589         *newgroup = MPI_GROUP_EMPTY;
590       } else {
591         *newgroup = smpi_group_new(size);
592         newrank = 0;
593         while (newrank < size) {
594           for (i = 0; i < n; i++) {
595             add = 1;
596             for (rank = ranges[i][0];   /* First */
597                  rank >= 0 && rank <= ranges[i][1];     /* Last */
598                  rank += ranges[i][2] /* Stride */ ) {
599               if (rank == newrank) {
600                 add = 0;
601                 break;
602               }
603             }
604             if (add == 1) {
605               index = smpi_group_index(group, newrank);
606               smpi_group_set_mapping(*newgroup, index, newrank);
607             }
608           }
609         }
610       }
611     }
612     smpi_group_use(*newgroup);
613     retval = MPI_SUCCESS;
614   }
615   smpi_bench_begin();
616   return retval;
617 }
618
619 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
620 {
621   int retval;
622
623   smpi_bench_end();
624   if (comm == MPI_COMM_NULL) {
625     retval = MPI_ERR_COMM;
626   } else if (rank == NULL) {
627     retval = MPI_ERR_ARG;
628   } else {
629     *rank = smpi_comm_rank(comm);
630     retval = MPI_SUCCESS;
631   }
632   smpi_bench_begin();
633   return retval;
634 }
635
636 int PMPI_Comm_size(MPI_Comm comm, int *size)
637 {
638   int retval;
639
640   smpi_bench_end();
641   if (comm == MPI_COMM_NULL) {
642     retval = MPI_ERR_COMM;
643   } else if (size == NULL) {
644     retval = MPI_ERR_ARG;
645   } else {
646     *size = smpi_comm_size(comm);
647     retval = MPI_SUCCESS;
648   }
649   smpi_bench_begin();
650   return retval;
651 }
652
653 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
654 {
655   int retval;
656
657   smpi_bench_end();
658   if (comm == MPI_COMM_NULL)  {
659     retval = MPI_ERR_COMM;
660   } else if (name == NULL || len == NULL)  {
661     retval = MPI_ERR_ARG;
662   } else {
663     smpi_comm_get_name(comm, name, len);
664     retval = MPI_SUCCESS;
665   }
666   smpi_bench_begin();
667   return retval;
668 }
669
670 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
671 {
672   int retval;
673
674   smpi_bench_end();
675   if (comm == MPI_COMM_NULL) {
676     retval = MPI_ERR_COMM;
677   } else if (group == NULL) {
678     retval = MPI_ERR_ARG;
679   } else {
680     *group = smpi_comm_group(comm);
681     retval = MPI_SUCCESS;
682   }
683   smpi_bench_begin();
684   return retval;
685 }
686
687 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
688 {
689   int retval;
690
691   smpi_bench_end();
692   if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
693     retval = MPI_ERR_COMM;
694   } else if (result == NULL) {
695     retval = MPI_ERR_ARG;
696   } else {
697     if (comm1 == comm2) {       /* Same communicators means same groups */
698       *result = MPI_IDENT;
699     } else {
700       *result =
701           smpi_group_compare(smpi_comm_group(comm1),
702                              smpi_comm_group(comm2));
703       if (*result == MPI_IDENT) {
704         *result = MPI_CONGRUENT;
705       }
706     }
707     retval = MPI_SUCCESS;
708   }
709   smpi_bench_begin();
710   return retval;
711 }
712
713 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
714 {
715   int retval;
716
717   smpi_bench_end();
718   if (comm == MPI_COMM_NULL) {
719     retval = MPI_ERR_COMM;
720   } else if (newcomm == NULL) {
721     retval = MPI_ERR_ARG;
722   } else {
723     *newcomm = smpi_comm_new(smpi_comm_group(comm));
724     retval = MPI_SUCCESS;
725   }
726   smpi_bench_begin();
727   return retval;
728 }
729
730 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
731 {
732   int retval;
733
734   smpi_bench_end();
735   if (comm == MPI_COMM_NULL) {
736     retval = MPI_ERR_COMM;
737   } else if (group == MPI_GROUP_NULL) {
738     retval = MPI_ERR_GROUP;
739   } else if (newcomm == NULL) {
740     retval = MPI_ERR_ARG;
741   } else {
742     *newcomm = smpi_comm_new(group);
743     retval = MPI_SUCCESS;
744   }
745   smpi_bench_begin();
746   return retval;
747 }
748
749 int PMPI_Comm_free(MPI_Comm * comm)
750 {
751   int retval;
752
753   smpi_bench_end();
754   if (comm == NULL) {
755     retval = MPI_ERR_ARG;
756   } else if (*comm == MPI_COMM_NULL) {
757     retval = MPI_ERR_COMM;
758   } else {
759     smpi_comm_destroy(*comm);
760     *comm = MPI_COMM_NULL;
761     retval = MPI_SUCCESS;
762   }
763   smpi_bench_begin();
764   return retval;
765 }
766
767 int PMPI_Comm_disconnect(MPI_Comm * comm)
768 {
769   /* TODO: wait until all communication in comm are done */
770   int retval;
771
772   smpi_bench_end();
773   if (comm == NULL) {
774     retval = MPI_ERR_ARG;
775   } else if (*comm == MPI_COMM_NULL) {
776     retval = MPI_ERR_COMM;
777   } else {
778     smpi_comm_destroy(*comm);
779     *comm = MPI_COMM_NULL;
780     retval = MPI_SUCCESS;
781   }
782   smpi_bench_begin();
783   return retval;
784 }
785
786 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
787 {
788   int retval;
789
790   smpi_bench_end();
791   if (comm_out == NULL) {
792     retval = MPI_ERR_ARG;
793   } else if (comm == MPI_COMM_NULL) {
794     retval = MPI_ERR_COMM;
795   } else {
796     *comm_out = smpi_comm_split(comm, color, key);
797     retval = MPI_SUCCESS;
798   }
799   smpi_bench_begin();
800   return retval;
801 }
802
803 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
804                   int tag, MPI_Comm comm, MPI_Request * request)
805 {
806   int retval;
807
808   smpi_bench_end();
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_send_init(buf, count, datatype, dst, tag, comm);
815     retval = MPI_SUCCESS;
816   }
817   smpi_bench_begin();
818   return retval;
819 }
820
821 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
822                   int tag, MPI_Comm comm, MPI_Request * request)
823 {
824   int retval;
825
826   smpi_bench_end();
827   if (request == NULL) {
828     retval = MPI_ERR_ARG;
829   } else if (comm == MPI_COMM_NULL) {
830     retval = MPI_ERR_COMM;
831   } else {
832     *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
833     retval = MPI_SUCCESS;
834   }
835   smpi_bench_begin();
836   return retval;
837 }
838
839 int PMPI_Start(MPI_Request * request)
840 {
841   int retval;
842
843   smpi_bench_end();
844   if (request == NULL) {
845     retval = MPI_ERR_ARG;
846   } else {
847     smpi_mpi_start(*request);
848     retval = MPI_SUCCESS;
849   }
850   smpi_bench_begin();
851   return retval;
852 }
853
854 int PMPI_Startall(int count, MPI_Request * requests)
855 {
856   int retval;
857
858   smpi_bench_end();
859   if (requests == NULL) {
860     retval = MPI_ERR_ARG;
861   } else {
862     smpi_mpi_startall(count, requests);
863     retval = MPI_SUCCESS;
864   }
865   smpi_bench_begin();
866   return retval;
867 }
868
869 int PMPI_Request_free(MPI_Request * request)
870 {
871   int retval;
872
873   smpi_bench_end();
874   if (request == NULL) {
875     retval = MPI_ERR_ARG;
876   } else {
877     smpi_mpi_request_free(request);
878     retval = MPI_SUCCESS;
879   }
880   smpi_bench_begin();
881   return retval;
882 }
883
884 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
885               int tag, MPI_Comm comm, MPI_Request * request)
886 {
887   int retval;
888
889   smpi_bench_end();
890 #ifdef HAVE_TRACING
891   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
892   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
893   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
894 #endif
895   if (request == NULL) {
896     retval = MPI_ERR_ARG;
897   } else if (comm == MPI_COMM_NULL) {
898     retval = MPI_ERR_COMM;
899   } else {
900     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
901     retval = MPI_SUCCESS;
902   }
903 #ifdef HAVE_TRACING
904   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
905   (*request)->recv = 1;
906 #endif
907   smpi_bench_begin();
908   return retval;
909 }
910
911 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
912               int tag, MPI_Comm comm, MPI_Request * request)
913 {
914   int retval;
915
916   smpi_bench_end();
917 #ifdef HAVE_TRACING
918   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
919   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
920   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
921   TRACE_smpi_send(rank, rank, dst_traced);
922 #endif
923   if (request == NULL) {
924     retval = MPI_ERR_ARG;
925   } else if (comm == MPI_COMM_NULL) {
926     retval = MPI_ERR_COMM;
927   } else {
928     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
929     retval = MPI_SUCCESS;
930   }
931 #ifdef HAVE_TRACING
932   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
933   (*request)->send = 1;
934 #endif
935   smpi_bench_begin();
936   return retval;
937 }
938
939 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
940              MPI_Comm comm, MPI_Status * status)
941 {
942   int retval;
943
944   smpi_bench_end();
945 #ifdef HAVE_TRACING
946   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
947   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
948   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
949 #endif
950   if (comm == MPI_COMM_NULL) {
951     retval = MPI_ERR_COMM;
952   } else {
953     smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
954     retval = MPI_SUCCESS;
955   }
956 #ifdef HAVE_TRACING
957   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
958   TRACE_smpi_recv(rank, src_traced, rank);
959 #endif
960   smpi_bench_begin();
961   return retval;
962 }
963
964 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
965              MPI_Comm comm)
966 {
967   int retval;
968
969   smpi_bench_end();
970 #ifdef HAVE_TRACING
971   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
972   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
973   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
974   TRACE_smpi_send(rank, rank, dst_traced);
975 #endif
976   if (comm == MPI_COMM_NULL) {
977     retval = MPI_ERR_COMM;
978   } else {
979     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
980     retval = MPI_SUCCESS;
981   }
982 #ifdef HAVE_TRACING
983   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
984 #endif
985   smpi_bench_begin();
986   return retval;
987 }
988
989 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
990                  int dst, int sendtag, void *recvbuf, int recvcount,
991                  MPI_Datatype recvtype, int src, int recvtag,
992                  MPI_Comm comm, MPI_Status * status)
993 {
994   int retval;
995
996   smpi_bench_end();
997 #ifdef HAVE_TRACING
998   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
999   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
1000   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
1001   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1002   TRACE_smpi_send(rank, rank, dst_traced);
1003   TRACE_smpi_send(rank, src_traced, rank);
1004 #endif
1005   if (comm == MPI_COMM_NULL) {
1006     retval = MPI_ERR_COMM;
1007   } else if (sendtype == MPI_DATATYPE_NULL
1008              || recvtype == MPI_DATATYPE_NULL) {
1009     retval = MPI_ERR_TYPE;
1010   } else {
1011     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1012                       recvcount, recvtype, src, recvtag, comm, status);
1013     retval = MPI_SUCCESS;
1014   }
1015 #ifdef HAVE_TRACING
1016   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1017   TRACE_smpi_recv(rank, rank, dst_traced);
1018   TRACE_smpi_recv(rank, src_traced, rank);
1019 #endif
1020   smpi_bench_begin();
1021   return retval;
1022 }
1023
1024 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1025                          int dst, int sendtag, int src, int recvtag,
1026                          MPI_Comm comm, MPI_Status * status)
1027 {
1028   //TODO: suboptimal implementation
1029   void *recvbuf;
1030   int retval, size;
1031
1032   size = smpi_datatype_size(datatype) * count;
1033   recvbuf = xbt_new(char, size);
1034   retval =
1035       MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1036                    datatype, src, recvtag, comm, status);
1037   memcpy(buf, recvbuf, size * sizeof(char));
1038   xbt_free(recvbuf);
1039   return retval;
1040 }
1041
1042 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1043 {
1044   int retval;
1045
1046   smpi_bench_end();
1047   if (request == NULL || flag == NULL) {
1048     retval = MPI_ERR_ARG;
1049   } else if (*request == MPI_REQUEST_NULL) {
1050     retval = MPI_ERR_REQUEST;
1051   } else {
1052     *flag = smpi_mpi_test(request, status);
1053     retval = MPI_SUCCESS;
1054   }
1055   smpi_bench_begin();
1056   return retval;
1057 }
1058
1059 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1060                 MPI_Status * status)
1061 {
1062   int retval;
1063
1064   smpi_bench_end();
1065   if (index == NULL || flag == NULL) {
1066     retval = MPI_ERR_ARG;
1067   } else {
1068     *flag = smpi_mpi_testany(count, requests, index, status);
1069     retval = MPI_SUCCESS;
1070   }
1071   smpi_bench_begin();
1072   return retval;
1073 }
1074
1075 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1076 {
1077   int retval;
1078
1079   smpi_bench_end();
1080 #ifdef HAVE_TRACING
1081   int rank = request && (*request)->comm != MPI_COMM_NULL
1082       ? smpi_comm_rank((*request)->comm)
1083       : -1;
1084   MPI_Group group = smpi_comm_group((*request)->comm);
1085   int src_traced = smpi_group_rank(group, (*request)->src);
1086   int dst_traced = smpi_group_rank(group, (*request)->dst);
1087   int is_wait_for_receive = (*request)->recv;
1088   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1089 #endif
1090   if (request == NULL) {
1091     retval = MPI_ERR_ARG;
1092   } else if (*request == MPI_REQUEST_NULL) {
1093     retval = MPI_ERR_REQUEST;
1094   } else {
1095     smpi_mpi_wait(request, status);
1096     retval = MPI_SUCCESS;
1097   }
1098 #ifdef HAVE_TRACING
1099   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1100   if (is_wait_for_receive) {
1101     TRACE_smpi_recv(rank, src_traced, dst_traced);
1102   }
1103 #endif
1104   smpi_bench_begin();
1105   return retval;
1106 }
1107
1108 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1109 {
1110   int retval;
1111
1112   smpi_bench_end();
1113 #ifdef HAVE_TRACING
1114   //save requests information for tracing
1115   int i;
1116   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), NULL);
1117   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), NULL);
1118   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), NULL);
1119   for (i = 0; i < count; i++) {
1120     MPI_Request req = requests[i];      //already received requests are no longer valid
1121     if (req) {
1122       int *asrc = xbt_new(int, 1);
1123       int *adst = xbt_new(int, 1);
1124       int *arecv = xbt_new(int, 1);
1125       *asrc = req->src;
1126       *adst = req->dst;
1127       *arecv = req->recv;
1128       xbt_dynar_insert_at(srcs, i, asrc);
1129       xbt_dynar_insert_at(dsts, i, adst);
1130       xbt_dynar_insert_at(recvs, i, arecv);
1131       xbt_free(asrc);
1132       xbt_free(adst);
1133       xbt_free(arecv);
1134     } else {
1135       int *t = xbt_new(int, 1);
1136       xbt_dynar_insert_at(srcs, i, t);
1137       xbt_dynar_insert_at(dsts, i, t);
1138       xbt_dynar_insert_at(recvs, i, t);
1139       xbt_free(t);
1140     }
1141   }
1142   int rank_traced = smpi_comm_rank(MPI_COMM_WORLD);
1143   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1144 #endif
1145   if (index == NULL) {
1146     retval = MPI_ERR_ARG;
1147   } else {
1148     *index = smpi_mpi_waitany(count, requests, status);
1149     retval = MPI_SUCCESS;
1150   }
1151 #ifdef HAVE_TRACING
1152   int src_traced, dst_traced, is_wait_for_receive;
1153   xbt_dynar_get_cpy(srcs, *index, &src_traced);
1154   xbt_dynar_get_cpy(dsts, *index, &dst_traced);
1155   xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
1156   if (is_wait_for_receive) {
1157     TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1158   }
1159   TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1160   //clean-up of dynars
1161   xbt_dynar_free(&srcs);
1162   xbt_dynar_free(&dsts);
1163   xbt_dynar_free(&recvs);
1164 #endif
1165   smpi_bench_begin();
1166   return retval;
1167 }
1168
1169 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1170 {
1171
1172   smpi_bench_end();
1173 #ifdef HAVE_TRACING
1174   //save information from requests
1175   int i;
1176   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), NULL);
1177   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), NULL);
1178   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), NULL);
1179   for (i = 0; i < count; i++) {
1180     MPI_Request req = requests[i];      //all req should be valid in Waitall
1181     int *asrc = xbt_new(int, 1);
1182     int *adst = xbt_new(int, 1);
1183     int *arecv = xbt_new(int, 1);
1184     *asrc = req->src;
1185     *adst = req->dst;
1186     *arecv = req->recv;
1187     xbt_dynar_insert_at(srcs, i, asrc);
1188     xbt_dynar_insert_at(dsts, i, adst);
1189     xbt_dynar_insert_at(recvs, i, arecv);
1190     xbt_free(asrc);
1191     xbt_free(adst);
1192     xbt_free(arecv);
1193   }
1194   int rank_traced = smpi_comm_rank (MPI_COMM_WORLD);
1195   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1196 #endif
1197   smpi_mpi_waitall(count, requests, status);
1198 #ifdef HAVE_TRACING
1199   for (i = 0; i < count; i++) {
1200     int src_traced, dst_traced, is_wait_for_receive;
1201     xbt_dynar_get_cpy(srcs, i, &src_traced);
1202     xbt_dynar_get_cpy(dsts, i, &dst_traced);
1203     xbt_dynar_get_cpy(recvs, i, &is_wait_for_receive);
1204     if (is_wait_for_receive) {
1205       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1206     }
1207   }
1208   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1209   //clean-up of dynars
1210   xbt_dynar_free(&srcs);
1211   xbt_dynar_free(&dsts);
1212   xbt_dynar_free(&recvs);
1213 #endif
1214   smpi_bench_begin();
1215   return MPI_SUCCESS;
1216 }
1217
1218 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1219                  int *indices, MPI_Status status[])
1220 {
1221   int retval;
1222
1223   smpi_bench_end();
1224   if (outcount == NULL || indices == NULL) {
1225     retval = MPI_ERR_ARG;
1226   } else {
1227     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1228     retval = MPI_SUCCESS;
1229   }
1230   smpi_bench_begin();
1231   return retval;
1232 }
1233
1234 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1235 {
1236   int retval;
1237
1238   smpi_bench_end();
1239 #ifdef HAVE_TRACING
1240   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1241   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1242   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1243 #endif
1244   if (comm == MPI_COMM_NULL) {
1245     retval = MPI_ERR_COMM;
1246   } else {
1247     smpi_mpi_bcast(buf, count, datatype, root, comm);
1248     retval = MPI_SUCCESS;
1249   }
1250 #ifdef HAVE_TRACING
1251   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1252 #endif
1253   smpi_bench_begin();
1254   return retval;
1255 }
1256
1257 int PMPI_Barrier(MPI_Comm comm)
1258 {
1259   int retval;
1260
1261   smpi_bench_end();
1262 #ifdef HAVE_TRACING
1263   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1264   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1265 #endif
1266   if (comm == MPI_COMM_NULL) {
1267     retval = MPI_ERR_COMM;
1268   } else {
1269     smpi_mpi_barrier(comm);
1270     retval = MPI_SUCCESS;
1271   }
1272 #ifdef HAVE_TRACING
1273   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1274 #endif
1275   smpi_bench_begin();
1276   return retval;
1277 }
1278
1279 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1280                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1281                int root, MPI_Comm comm)
1282 {
1283   int retval;
1284
1285   smpi_bench_end();
1286 #ifdef HAVE_TRACING
1287   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1288   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1289   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1290 #endif
1291   if (comm == MPI_COMM_NULL) {
1292     retval = MPI_ERR_COMM;
1293   } else if (sendtype == MPI_DATATYPE_NULL
1294              || recvtype == MPI_DATATYPE_NULL) {
1295     retval = MPI_ERR_TYPE;
1296   } else {
1297     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1298                     recvtype, root, comm);
1299     retval = MPI_SUCCESS;
1300   }
1301 #ifdef HAVE_TRACING
1302   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1303 #endif
1304   smpi_bench_begin();
1305   return retval;
1306 }
1307
1308 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1309                 void *recvbuf, int *recvcounts, int *displs,
1310                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1311 {
1312   int retval;
1313
1314   smpi_bench_end();
1315 #ifdef HAVE_TRACING
1316   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1317   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1318   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1319 #endif
1320   if (comm == MPI_COMM_NULL) {
1321     retval = MPI_ERR_COMM;
1322   } else if (sendtype == MPI_DATATYPE_NULL
1323              || recvtype == MPI_DATATYPE_NULL) {
1324     retval = MPI_ERR_TYPE;
1325   } else if (recvcounts == NULL || displs == NULL) {
1326     retval = MPI_ERR_ARG;
1327   } else {
1328     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1329                      displs, recvtype, root, comm);
1330     retval = MPI_SUCCESS;
1331   }
1332 #ifdef HAVE_TRACING
1333   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1334 #endif
1335   smpi_bench_begin();
1336   return retval;
1337 }
1338
1339 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1340                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1341                   MPI_Comm comm)
1342 {
1343   int retval;
1344
1345   smpi_bench_end();
1346 #ifdef HAVE_TRACING
1347   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1348   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1349 #endif
1350   if (comm == MPI_COMM_NULL) {
1351     retval = MPI_ERR_COMM;
1352   } else if (sendtype == MPI_DATATYPE_NULL
1353              || recvtype == MPI_DATATYPE_NULL) {
1354     retval = MPI_ERR_TYPE;
1355   } else {
1356     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1357                        recvtype, comm);
1358     retval = MPI_SUCCESS;
1359   }
1360 #ifdef HAVE_TRACING
1361   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1362 #endif
1363   smpi_bench_begin();
1364   return retval;
1365 }
1366
1367 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1368                    void *recvbuf, int *recvcounts, int *displs,
1369                    MPI_Datatype recvtype, MPI_Comm comm)
1370 {
1371   int retval;
1372
1373   smpi_bench_end();
1374 #ifdef HAVE_TRACING
1375   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1376   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1377 #endif
1378   if (comm == MPI_COMM_NULL) {
1379     retval = MPI_ERR_COMM;
1380   } else if (sendtype == MPI_DATATYPE_NULL
1381              || recvtype == MPI_DATATYPE_NULL) {
1382     retval = MPI_ERR_TYPE;
1383   } else if (recvcounts == NULL || displs == NULL) {
1384     retval = MPI_ERR_ARG;
1385   } else {
1386     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1387                         displs, recvtype, comm);
1388     retval = MPI_SUCCESS;
1389   }
1390 #ifdef HAVE_TRACING
1391   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1392 #endif
1393   smpi_bench_begin();
1394   return retval;
1395 }
1396
1397 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1398                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1399                 int root, MPI_Comm comm)
1400 {
1401   int retval;
1402
1403   smpi_bench_end();
1404 #ifdef HAVE_TRACING
1405   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1406   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1407   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1408 #endif
1409   if (comm == MPI_COMM_NULL) {
1410     retval = MPI_ERR_COMM;
1411   } else if (sendtype == MPI_DATATYPE_NULL
1412              || recvtype == MPI_DATATYPE_NULL) {
1413     retval = MPI_ERR_TYPE;
1414   } else {
1415     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1416                      recvtype, root, comm);
1417     retval = MPI_SUCCESS;
1418   }
1419 #ifdef HAVE_TRACING
1420   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1421 #endif
1422   smpi_bench_begin();
1423   return retval;
1424 }
1425
1426 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1427                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1428                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1429 {
1430   int retval;
1431
1432   smpi_bench_end();
1433 #ifdef HAVE_TRACING
1434   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1435   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1436   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1437 #endif
1438   if (comm == MPI_COMM_NULL) {
1439     retval = MPI_ERR_COMM;
1440   } else if (sendtype == MPI_DATATYPE_NULL
1441              || recvtype == MPI_DATATYPE_NULL) {
1442     retval = MPI_ERR_TYPE;
1443   } else if (sendcounts == NULL || displs == NULL) {
1444     retval = MPI_ERR_ARG;
1445   } else {
1446     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1447                       recvcount, recvtype, root, comm);
1448     retval = MPI_SUCCESS;
1449   }
1450 #ifdef HAVE_TRACING
1451   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1452 #endif
1453   smpi_bench_begin();
1454   return retval;
1455 }
1456
1457 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1458                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1459 {
1460   int retval;
1461
1462   smpi_bench_end();
1463 #ifdef HAVE_TRACING
1464   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1465   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1466   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1467 #endif
1468   if (comm == MPI_COMM_NULL) {
1469     retval = MPI_ERR_COMM;
1470   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1471     retval = MPI_ERR_ARG;
1472   } else {
1473     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1474     retval = MPI_SUCCESS;
1475   }
1476 #ifdef HAVE_TRACING
1477   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1478 #endif
1479   smpi_bench_begin();
1480   return retval;
1481 }
1482
1483 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1484                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1485 {
1486   int retval;
1487
1488   smpi_bench_end();
1489 #ifdef HAVE_TRACING
1490   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1491   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1492 #endif
1493   if (comm == MPI_COMM_NULL) {
1494     retval = MPI_ERR_COMM;
1495   } else if (datatype == MPI_DATATYPE_NULL) {
1496     retval = MPI_ERR_TYPE;
1497   } else if (op == MPI_OP_NULL) {
1498     retval = MPI_ERR_OP;
1499   } else {
1500     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
1501     retval = MPI_SUCCESS;
1502   }
1503 #ifdef HAVE_TRACING
1504   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1505 #endif
1506   smpi_bench_begin();
1507   return retval;
1508 }
1509
1510 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
1511              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1512 {
1513   int retval;
1514
1515   smpi_bench_end();
1516 #ifdef HAVE_TRACING
1517   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1518   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1519 #endif
1520   if (comm == MPI_COMM_NULL) {
1521     retval = MPI_ERR_COMM;
1522   } else if (datatype == MPI_DATATYPE_NULL) {
1523     retval = MPI_ERR_TYPE;
1524   } else if (op == MPI_OP_NULL) {
1525     retval = MPI_ERR_OP;
1526   } else {
1527     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1528     retval = MPI_SUCCESS;
1529   }
1530 #ifdef HAVE_TRACING
1531   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1532 #endif
1533   smpi_bench_begin();
1534   return retval;
1535 }
1536
1537 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1538                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1539 {
1540   int retval, i, size, count;
1541   int *displs;
1542   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1543
1544   smpi_bench_end();
1545 #ifdef HAVE_TRACING
1546   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1547 #endif
1548   if (comm == MPI_COMM_NULL) {
1549     retval = MPI_ERR_COMM;
1550   } else if (datatype == MPI_DATATYPE_NULL) {
1551     retval = MPI_ERR_TYPE;
1552   } else if (op == MPI_OP_NULL) {
1553     retval = MPI_ERR_OP;
1554   } else if (recvcounts == NULL) {
1555     retval = MPI_ERR_ARG;
1556   } else {
1557     /* arbitrarily choose root as rank 0 */
1558     /* TODO: faster direct implementation ? */
1559     size = smpi_comm_size(comm);
1560     count = 0;
1561     displs = xbt_new(int, size);
1562     for (i = 0; i < size; i++) {
1563       count += recvcounts[i];
1564       displs[i] = 0;
1565     }
1566     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1567     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf,
1568                       recvcounts[rank], datatype, 0, comm);
1569     xbt_free(displs);
1570     retval = MPI_SUCCESS;
1571   }
1572 #ifdef HAVE_TRACING
1573   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1574 #endif
1575   smpi_bench_begin();
1576   return retval;
1577 }
1578
1579 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1580                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
1581                  MPI_Comm comm)
1582 {
1583   int retval, size, sendsize;
1584
1585   smpi_bench_end();
1586 #ifdef HAVE_TRACING
1587   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1588   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1589 #endif
1590   if (comm == MPI_COMM_NULL) {
1591     retval = MPI_ERR_COMM;
1592   } else if (sendtype == MPI_DATATYPE_NULL
1593              || recvtype == MPI_DATATYPE_NULL) {
1594     retval = MPI_ERR_TYPE;
1595   } else {
1596     size = smpi_comm_size(comm);
1597     sendsize = smpi_datatype_size(sendtype) * sendcount;
1598     if (sendsize < 200 && size > 12) {
1599       retval =
1600           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
1601                                          recvbuf, recvcount, recvtype,
1602                                          comm);
1603     } else if (sendsize < 3000) {
1604       retval =
1605           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
1606                                                 sendtype, recvbuf,
1607                                                 recvcount, recvtype, comm);
1608     } else {
1609       retval =
1610           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
1611                                             recvbuf, recvcount, recvtype,
1612                                             comm);
1613     }
1614   }
1615 #ifdef HAVE_TRACING
1616   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1617 #endif
1618   smpi_bench_begin();
1619   return retval;
1620 }
1621
1622 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
1623                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
1624                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
1625 {
1626   int retval;
1627
1628   smpi_bench_end();
1629 #ifdef HAVE_TRACING
1630   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1631   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1632 #endif
1633   if (comm == MPI_COMM_NULL) {
1634     retval = MPI_ERR_COMM;
1635   } else if (sendtype == MPI_DATATYPE_NULL
1636              || recvtype == MPI_DATATYPE_NULL) {
1637     retval = MPI_ERR_TYPE;
1638   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
1639              || recvdisps == NULL) {
1640     retval = MPI_ERR_ARG;
1641   } else {
1642     retval =
1643         smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
1644                                   recvbuf, recvcounts, recvdisps, recvtype,
1645                                   comm);
1646   }
1647 #ifdef HAVE_TRACING
1648   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1649 #endif
1650   smpi_bench_begin();
1651   return retval;
1652 }
1653
1654
1655 int PMPI_Get_processor_name(char *name, int *resultlen)
1656 {
1657   int retval = MPI_SUCCESS;
1658
1659   smpi_bench_end();
1660   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
1661           MPI_MAX_PROCESSOR_NAME - 1);
1662   *resultlen =
1663       strlen(name) >
1664       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1665
1666   smpi_bench_begin();
1667   return retval;
1668 }
1669
1670 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
1671 {
1672   int retval = MPI_SUCCESS;
1673   size_t size;
1674
1675   smpi_bench_end();
1676   if (status == NULL || count == NULL) {
1677     retval = MPI_ERR_ARG;
1678   } else if (datatype == MPI_DATATYPE_NULL) {
1679     retval = MPI_ERR_TYPE;
1680   } else {
1681     size = smpi_datatype_size(datatype);
1682     if (size == 0) {
1683       *count = 0;
1684     } else if (status->count % size != 0) {
1685       retval = MPI_UNDEFINED;
1686     } else {
1687       *count = smpi_mpi_get_count(status, datatype);
1688     }
1689   }
1690   smpi_bench_begin();
1691   return retval;
1692 }
1693
1694 /* The following calls are not yet implemented and will fail at runtime. */
1695 /* Once implemented, please move them above this notice. */
1696
1697 static int not_yet_implemented(void) {
1698    xbt_die("Not yet implemented");
1699    return MPI_ERR_OTHER;
1700 }
1701
1702 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
1703    return not_yet_implemented();
1704 }
1705
1706 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
1707    return not_yet_implemented();
1708 }
1709
1710 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
1711    return not_yet_implemented();
1712 }
1713
1714 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
1715    return not_yet_implemented();
1716 }
1717
1718 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
1719    return not_yet_implemented();
1720 }
1721
1722 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
1723    return not_yet_implemented();
1724 }
1725
1726 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
1727    return not_yet_implemented();
1728 }
1729
1730 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
1731    return not_yet_implemented();
1732 }
1733
1734 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
1735    return not_yet_implemented();
1736 }
1737
1738 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
1739    return not_yet_implemented();
1740 }
1741
1742 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
1743    return not_yet_implemented();
1744 }
1745
1746 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
1747    return not_yet_implemented();
1748 }
1749
1750 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
1751    return not_yet_implemented();
1752 }
1753
1754 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
1755    return not_yet_implemented();
1756 }
1757
1758 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
1759    return not_yet_implemented();
1760 }
1761
1762 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
1763    return not_yet_implemented();
1764 }
1765
1766 int PMPI_Error_class(int errorcode, int* errorclass) {
1767    return not_yet_implemented();
1768 }
1769
1770 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
1771    return not_yet_implemented();
1772 }
1773
1774 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
1775    return not_yet_implemented();
1776 }
1777
1778 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
1779    return not_yet_implemented();
1780 }
1781
1782 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
1783    return not_yet_implemented();
1784 }
1785
1786 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
1787    return not_yet_implemented();
1788 }
1789
1790 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* newtype) {
1791    return not_yet_implemented();
1792 }
1793
1794 int PMPI_Cancel(MPI_Request* request) {
1795    return not_yet_implemented();
1796 }
1797
1798 int PMPI_Buffer_attach(void* buffer, int size) {
1799    return not_yet_implemented();
1800 }
1801
1802 int PMPI_Buffer_detach(void* buffer, int* size) {
1803    return not_yet_implemented();
1804 }
1805
1806 int PMPI_Testsome(int incount, MPI_Request* requests, int* outcount, int* indices, MPI_Status* statuses) {
1807    return not_yet_implemented();
1808 }
1809
1810 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
1811    return not_yet_implemented();
1812 }
1813
1814 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
1815    return not_yet_implemented();
1816 }
1817
1818 int PMPI_Type_commit(MPI_Datatype* datatype) {
1819    return not_yet_implemented();
1820 }
1821
1822 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
1823    return not_yet_implemented();
1824 }
1825
1826 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
1827    return not_yet_implemented();
1828 }
1829
1830 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
1831    return not_yet_implemented();
1832 }
1833
1834 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* newtype) {
1835    return not_yet_implemented();
1836 }
1837
1838 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
1839    return not_yet_implemented();
1840 }
1841
1842 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1843    return not_yet_implemented();
1844 }
1845
1846 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1847    return not_yet_implemented();
1848 }
1849
1850 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
1851    return not_yet_implemented();
1852 }
1853
1854 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
1855    return not_yet_implemented();
1856 }
1857
1858 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1859    return not_yet_implemented();
1860 }
1861
1862 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1863    return not_yet_implemented();
1864 }
1865
1866 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1867    return not_yet_implemented();
1868 }
1869
1870 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
1871    return not_yet_implemented();
1872 }
1873
1874 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
1875    return not_yet_implemented();
1876 }
1877
1878 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1879    return not_yet_implemented();
1880 }
1881
1882 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1883    return not_yet_implemented();
1884 }
1885
1886 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
1887    return not_yet_implemented();
1888 }
1889
1890 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
1891    return not_yet_implemented();
1892 }
1893
1894 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
1895    return not_yet_implemented();
1896 }
1897
1898 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1899    return not_yet_implemented();
1900 }
1901
1902 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1903    return not_yet_implemented();
1904 }
1905
1906 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1907    return not_yet_implemented();
1908 }
1909
1910 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
1911    return not_yet_implemented();
1912 }
1913
1914 int PMPI_Keyval_free(int* keyval) {
1915    return not_yet_implemented();
1916 }
1917
1918 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
1919    return not_yet_implemented();
1920 }
1921
1922 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
1923    return not_yet_implemented();
1924 }
1925
1926 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses) {
1927    return not_yet_implemented();
1928 }
1929
1930 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
1931    return not_yet_implemented();
1932 }
1933
1934 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
1935    return not_yet_implemented();
1936 }
1937
1938 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1939    return not_yet_implemented();
1940 }
1941
1942 int PMPI_Initialized(int* flag) {
1943    return not_yet_implemented();
1944 }