Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
67499281b204cea9e0670a2aa17dbe503da49828
[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_coll_private.h"
9 #include "smpi_mpi_dt_private.h"
10
11 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_pmpi, smpi,
12                                 "Logging specific to SMPI (pmpi)");
13
14 #ifdef HAVE_TRACING
15 //this function need to be here because of the calls to smpi_bench
16 void TRACE_smpi_set_category(const char *category)
17 {
18   //need to end bench otherwise categories for execution tasks are wrong
19   smpi_bench_end();
20   TRACE_internal_smpi_set_category (category);
21   //begin bench after changing process's category
22   smpi_bench_begin();
23 }
24 #endif
25
26 /* PMPI User level calls */
27
28 int PMPI_Init(int *argc, char ***argv)
29 {
30   smpi_process_init(argc, argv);
31 #ifdef HAVE_TRACING
32   TRACE_smpi_init(smpi_process_index());
33 #endif
34   smpi_bench_begin();
35   return MPI_SUCCESS;
36 }
37
38 int PMPI_Finalize(void)
39 {
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   SIMIX_req_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 {
627     *rank = smpi_comm_rank(comm);
628     retval = MPI_SUCCESS;
629   }
630   smpi_bench_begin();
631   return retval;
632 }
633
634 int PMPI_Comm_size(MPI_Comm comm, int *size)
635 {
636   int retval;
637
638   smpi_bench_end();
639   if (comm == MPI_COMM_NULL) {
640     retval = MPI_ERR_COMM;
641   } else if (size == NULL) {
642     retval = MPI_ERR_ARG;
643   } else {
644     *size = smpi_comm_size(comm);
645     retval = MPI_SUCCESS;
646   }
647   smpi_bench_begin();
648   return retval;
649 }
650
651 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
652 {
653   int retval;
654
655   smpi_bench_end();
656   if (comm == MPI_COMM_NULL) {
657     retval = MPI_ERR_COMM;
658   } else if (group == NULL) {
659     retval = MPI_ERR_ARG;
660   } else {
661     *group = smpi_comm_group(comm);
662     retval = MPI_SUCCESS;
663   }
664   smpi_bench_begin();
665   return retval;
666 }
667
668 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
669 {
670   int retval;
671
672   smpi_bench_end();
673   if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
674     retval = MPI_ERR_COMM;
675   } else if (result == NULL) {
676     retval = MPI_ERR_ARG;
677   } else {
678     if (comm1 == comm2) {       /* Same communicators means same groups */
679       *result = MPI_IDENT;
680     } else {
681       *result =
682           smpi_group_compare(smpi_comm_group(comm1),
683                              smpi_comm_group(comm2));
684       if (*result == MPI_IDENT) {
685         *result = MPI_CONGRUENT;
686       }
687     }
688     retval = MPI_SUCCESS;
689   }
690   smpi_bench_begin();
691   return retval;
692 }
693
694 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
695 {
696   int retval;
697
698   smpi_bench_end();
699   if (comm == MPI_COMM_NULL) {
700     retval = MPI_ERR_COMM;
701   } else if (newcomm == NULL) {
702     retval = MPI_ERR_ARG;
703   } else {
704     *newcomm = smpi_comm_new(smpi_comm_group(comm));
705     retval = MPI_SUCCESS;
706   }
707   smpi_bench_begin();
708   return retval;
709 }
710
711 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
712 {
713   int retval;
714
715   smpi_bench_end();
716   if (comm == MPI_COMM_NULL) {
717     retval = MPI_ERR_COMM;
718   } else if (group == MPI_GROUP_NULL) {
719     retval = MPI_ERR_GROUP;
720   } else if (newcomm == NULL) {
721     retval = MPI_ERR_ARG;
722   } else {
723     *newcomm = smpi_comm_new(group);
724     retval = MPI_SUCCESS;
725   }
726   smpi_bench_begin();
727   return retval;
728 }
729
730 int PMPI_Comm_free(MPI_Comm * comm)
731 {
732   int retval;
733
734   smpi_bench_end();
735   if (comm == NULL) {
736     retval = MPI_ERR_ARG;
737   } else if (*comm == MPI_COMM_NULL) {
738     retval = MPI_ERR_COMM;
739   } else {
740     smpi_comm_destroy(*comm);
741     *comm = MPI_COMM_NULL;
742     retval = MPI_SUCCESS;
743   }
744   smpi_bench_begin();
745   return retval;
746 }
747
748 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
749 {
750   int retval;
751
752   smpi_bench_end();
753   if (comm_out == NULL) {
754     retval = MPI_ERR_ARG;
755   } else if (comm == MPI_COMM_NULL) {
756     retval = MPI_ERR_COMM;
757   } else {
758     *comm_out = smpi_comm_split(comm, color, key);
759     retval = MPI_SUCCESS;
760   }
761   smpi_bench_begin();
762   return retval;
763 }
764
765 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
766                   int tag, MPI_Comm comm, MPI_Request * request)
767 {
768   int retval;
769
770   smpi_bench_end();
771   if (request == NULL) {
772     retval = MPI_ERR_ARG;
773   } else if (comm == MPI_COMM_NULL) {
774     retval = MPI_ERR_COMM;
775   } else {
776     *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
777     retval = MPI_SUCCESS;
778   }
779   smpi_bench_begin();
780   return retval;
781 }
782
783 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
784                   int tag, MPI_Comm comm, MPI_Request * request)
785 {
786   int retval;
787
788   smpi_bench_end();
789   if (request == NULL) {
790     retval = MPI_ERR_ARG;
791   } else if (comm == MPI_COMM_NULL) {
792     retval = MPI_ERR_COMM;
793   } else {
794     *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
795     retval = MPI_SUCCESS;
796   }
797   smpi_bench_begin();
798   return retval;
799 }
800
801 int PMPI_Start(MPI_Request * request)
802 {
803   int retval;
804
805   smpi_bench_end();
806   if (request == NULL) {
807     retval = MPI_ERR_ARG;
808   } else {
809     smpi_mpi_start(*request);
810     retval = MPI_SUCCESS;
811   }
812   smpi_bench_begin();
813   return retval;
814 }
815
816 int PMPI_Startall(int count, MPI_Request * requests)
817 {
818   int retval;
819
820   smpi_bench_end();
821   if (requests == NULL) {
822     retval = MPI_ERR_ARG;
823   } else {
824     smpi_mpi_startall(count, requests);
825     retval = MPI_SUCCESS;
826   }
827   smpi_bench_begin();
828   return retval;
829 }
830
831 int PMPI_Request_free(MPI_Request * request)
832 {
833   int retval;
834
835   smpi_bench_end();
836   if (request == NULL) {
837     retval = MPI_ERR_ARG;
838   } else {
839     smpi_mpi_request_free(request);
840     retval = MPI_SUCCESS;
841   }
842   smpi_bench_begin();
843   return retval;
844 }
845
846 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
847               int tag, MPI_Comm comm, MPI_Request * request)
848 {
849   int retval;
850
851   smpi_bench_end();
852 #ifdef HAVE_TRACING
853   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
854   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
855   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
856 #endif
857   if (request == NULL) {
858     retval = MPI_ERR_ARG;
859   } else if (comm == MPI_COMM_NULL) {
860     retval = MPI_ERR_COMM;
861   } else {
862     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
863     retval = MPI_SUCCESS;
864   }
865 #ifdef HAVE_TRACING
866   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
867   (*request)->recv = 1;
868 #endif
869   smpi_bench_begin();
870   return retval;
871 }
872
873 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
874               int tag, MPI_Comm comm, MPI_Request * request)
875 {
876   int retval;
877
878   smpi_bench_end();
879 #ifdef HAVE_TRACING
880   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
881   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
882   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
883   TRACE_smpi_send(rank, rank, dst_traced);
884 #endif
885   if (request == NULL) {
886     retval = MPI_ERR_ARG;
887   } else if (comm == MPI_COMM_NULL) {
888     retval = MPI_ERR_COMM;
889   } else {
890     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
891     retval = MPI_SUCCESS;
892   }
893 #ifdef HAVE_TRACING
894   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
895   (*request)->send = 1;
896 #endif
897   smpi_bench_begin();
898   return retval;
899 }
900
901 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
902              MPI_Comm comm, MPI_Status * status)
903 {
904   int retval;
905
906   smpi_bench_end();
907 #ifdef HAVE_TRACING
908   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
909   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
910   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
911 #endif
912   if (comm == MPI_COMM_NULL) {
913     retval = MPI_ERR_COMM;
914   } else {
915     smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
916     retval = MPI_SUCCESS;
917   }
918 #ifdef HAVE_TRACING
919   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
920   TRACE_smpi_recv(rank, src_traced, rank);
921 #endif
922   smpi_bench_begin();
923   return retval;
924 }
925
926 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
927              MPI_Comm comm)
928 {
929   int retval;
930
931   smpi_bench_end();
932 #ifdef HAVE_TRACING
933   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
934   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
935   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
936   TRACE_smpi_send(rank, rank, dst_traced);
937 #endif
938   if (comm == MPI_COMM_NULL) {
939     retval = MPI_ERR_COMM;
940   } else {
941     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
942     retval = MPI_SUCCESS;
943   }
944 #ifdef HAVE_TRACING
945   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
946 #endif
947   smpi_bench_begin();
948   return retval;
949 }
950
951 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
952                  int dst, int sendtag, void *recvbuf, int recvcount,
953                  MPI_Datatype recvtype, int src, int recvtag,
954                  MPI_Comm comm, MPI_Status * status)
955 {
956   int retval;
957
958   smpi_bench_end();
959 #ifdef HAVE_TRACING
960   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
961   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
962   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
963   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
964   TRACE_smpi_send(rank, rank, dst_traced);
965   TRACE_smpi_send(rank, src_traced, rank);
966 #endif
967   if (comm == MPI_COMM_NULL) {
968     retval = MPI_ERR_COMM;
969   } else if (sendtype == MPI_DATATYPE_NULL
970              || recvtype == MPI_DATATYPE_NULL) {
971     retval = MPI_ERR_TYPE;
972   } else {
973     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
974                       recvcount, recvtype, src, recvtag, comm, status);
975     retval = MPI_SUCCESS;
976   }
977 #ifdef HAVE_TRACING
978   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
979   TRACE_smpi_recv(rank, rank, dst_traced);
980   TRACE_smpi_recv(rank, src_traced, rank);
981 #endif
982   smpi_bench_begin();
983   return retval;
984 }
985
986 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
987                          int dst, int sendtag, int src, int recvtag,
988                          MPI_Comm comm, MPI_Status * status)
989 {
990   //TODO: suboptimal implementation
991   void *recvbuf;
992   int retval, size;
993
994   size = smpi_datatype_size(datatype) * count;
995   recvbuf = xbt_new(char, size);
996   retval =
997       MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
998                    datatype, src, recvtag, comm, status);
999   memcpy(buf, recvbuf, size * sizeof(char));
1000   xbt_free(recvbuf);
1001   return retval;
1002 }
1003
1004 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1005 {
1006   int retval;
1007
1008   smpi_bench_end();
1009   if (request == NULL || flag == NULL) {
1010     retval = MPI_ERR_ARG;
1011   } else if (*request == MPI_REQUEST_NULL) {
1012     retval = MPI_ERR_REQUEST;
1013   } else {
1014     *flag = smpi_mpi_test(request, status);
1015     retval = MPI_SUCCESS;
1016   }
1017   smpi_bench_begin();
1018   return retval;
1019 }
1020
1021 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1022                 MPI_Status * status)
1023 {
1024   int retval;
1025
1026   smpi_bench_end();
1027   if (index == NULL || flag == NULL) {
1028     retval = MPI_ERR_ARG;
1029   } else {
1030     *flag = smpi_mpi_testany(count, requests, index, status);
1031     retval = MPI_SUCCESS;
1032   }
1033   smpi_bench_begin();
1034   return retval;
1035 }
1036
1037 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1038 {
1039   int retval;
1040
1041   smpi_bench_end();
1042 #ifdef HAVE_TRACING
1043   int rank = request && (*request)->comm != MPI_COMM_NULL
1044       ? smpi_comm_rank((*request)->comm)
1045       : -1;
1046   MPI_Group group = smpi_comm_group((*request)->comm);
1047   int src_traced = smpi_group_rank(group, (*request)->src);
1048   int dst_traced = smpi_group_rank(group, (*request)->dst);
1049   int is_wait_for_receive = (*request)->recv;
1050   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1051 #endif
1052   if (request == NULL) {
1053     retval = MPI_ERR_ARG;
1054   } else if (*request == MPI_REQUEST_NULL) {
1055     retval = MPI_ERR_REQUEST;
1056   } else {
1057     smpi_mpi_wait(request, status);
1058     retval = MPI_SUCCESS;
1059   }
1060 #ifdef HAVE_TRACING
1061   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1062   if (is_wait_for_receive) {
1063     TRACE_smpi_recv(rank, src_traced, dst_traced);
1064   }
1065 #endif
1066   smpi_bench_begin();
1067   return retval;
1068 }
1069
1070 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1071 {
1072   int retval;
1073
1074   smpi_bench_end();
1075 #ifdef HAVE_TRACING
1076   //save requests information for tracing
1077   int i;
1078   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
1079   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
1080   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
1081   for (i = 0; i < count; i++) {
1082     MPI_Request req = requests[i];      //already received requests are no longer valid
1083     if (req) {
1084       int *asrc = xbt_new(int, 1);
1085       int *adst = xbt_new(int, 1);
1086       int *arecv = xbt_new(int, 1);
1087       *asrc = req->src;
1088       *adst = req->dst;
1089       *arecv = req->recv;
1090       xbt_dynar_insert_at(srcs, i, asrc);
1091       xbt_dynar_insert_at(dsts, i, adst);
1092       xbt_dynar_insert_at(recvs, i, arecv);
1093     } else {
1094       int *t = xbt_new(int, 1);
1095       xbt_dynar_insert_at(srcs, i, t);
1096       xbt_dynar_insert_at(dsts, i, t);
1097       xbt_dynar_insert_at(recvs, i, t);
1098     }
1099   }
1100   int rank_traced = smpi_comm_rank(MPI_COMM_WORLD);
1101   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1102 #endif
1103   if (index == NULL) {
1104     retval = MPI_ERR_ARG;
1105   } else {
1106     *index = smpi_mpi_waitany(count, requests, status);
1107     retval = MPI_SUCCESS;
1108   }
1109 #ifdef HAVE_TRACING
1110   int src_traced, dst_traced, is_wait_for_receive;
1111   xbt_dynar_get_cpy(srcs, *index, &src_traced);
1112   xbt_dynar_get_cpy(dsts, *index, &dst_traced);
1113   xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
1114   if (is_wait_for_receive) {
1115     TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1116   }
1117   TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1118   //clean-up of dynars
1119   xbt_free(srcs);
1120   xbt_free(dsts);
1121   xbt_free(recvs);
1122 #endif
1123   smpi_bench_begin();
1124   return retval;
1125 }
1126
1127 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1128 {
1129
1130   smpi_bench_end();
1131 #ifdef HAVE_TRACING
1132   //save information from requests
1133   int i;
1134   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
1135   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
1136   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
1137   for (i = 0; i < count; i++) {
1138     MPI_Request req = requests[i];      //all req should be valid in Waitall
1139     int *asrc = xbt_new(int, 1);
1140     int *adst = xbt_new(int, 1);
1141     int *arecv = xbt_new(int, 1);
1142     *asrc = req->src;
1143     *adst = req->dst;
1144     *arecv = req->recv;
1145     xbt_dynar_insert_at(srcs, i, asrc);
1146     xbt_dynar_insert_at(dsts, i, adst);
1147     xbt_dynar_insert_at(recvs, i, arecv);
1148   }
1149   int rank_traced = smpi_comm_rank (MPI_COMM_WORLD);
1150   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1151 #endif
1152   smpi_mpi_waitall(count, requests, status);
1153 #ifdef HAVE_TRACING
1154   for (i = 0; i < count; i++) {
1155     int src_traced, dst_traced, is_wait_for_receive;
1156     xbt_dynar_get_cpy(srcs, i, &src_traced);
1157     xbt_dynar_get_cpy(dsts, i, &dst_traced);
1158     xbt_dynar_get_cpy(recvs, i, &is_wait_for_receive);
1159     if (is_wait_for_receive) {
1160       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1161     }
1162   }
1163   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1164   //clean-up of dynars
1165   xbt_free(srcs);
1166   xbt_free(dsts);
1167   xbt_free(recvs);
1168 #endif
1169   smpi_bench_begin();
1170   return MPI_SUCCESS;
1171 }
1172
1173 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1174                  int *indices, MPI_Status status[])
1175 {
1176   int retval;
1177
1178   smpi_bench_end();
1179   if (outcount == NULL || indices == NULL) {
1180     retval = MPI_ERR_ARG;
1181   } else {
1182     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1183     retval = MPI_SUCCESS;
1184   }
1185   smpi_bench_begin();
1186   return retval;
1187 }
1188
1189 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1190 {
1191   int retval;
1192
1193   smpi_bench_end();
1194 #ifdef HAVE_TRACING
1195   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1196   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1197   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1198 #endif
1199   if (comm == MPI_COMM_NULL) {
1200     retval = MPI_ERR_COMM;
1201   } else {
1202     smpi_mpi_bcast(buf, count, datatype, root, comm);
1203     retval = MPI_SUCCESS;
1204   }
1205 #ifdef HAVE_TRACING
1206   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1207 #endif
1208   smpi_bench_begin();
1209   return retval;
1210 }
1211
1212 int PMPI_Barrier(MPI_Comm comm)
1213 {
1214   int retval;
1215
1216   smpi_bench_end();
1217 #ifdef HAVE_TRACING
1218   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1219   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1220 #endif
1221   if (comm == MPI_COMM_NULL) {
1222     retval = MPI_ERR_COMM;
1223   } else {
1224     smpi_mpi_barrier(comm);
1225     retval = MPI_SUCCESS;
1226   }
1227 #ifdef HAVE_TRACING
1228   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1229 #endif
1230   smpi_bench_begin();
1231   return retval;
1232 }
1233
1234 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1235                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1236                int root, MPI_Comm comm)
1237 {
1238   int retval;
1239
1240   smpi_bench_end();
1241 #ifdef HAVE_TRACING
1242   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1243   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1244   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1245 #endif
1246   if (comm == MPI_COMM_NULL) {
1247     retval = MPI_ERR_COMM;
1248   } else if (sendtype == MPI_DATATYPE_NULL
1249              || recvtype == MPI_DATATYPE_NULL) {
1250     retval = MPI_ERR_TYPE;
1251   } else {
1252     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1253                     recvtype, root, comm);
1254     retval = MPI_SUCCESS;
1255   }
1256 #ifdef HAVE_TRACING
1257   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1258 #endif
1259   smpi_bench_begin();
1260   return retval;
1261 }
1262
1263 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1264                 void *recvbuf, int *recvcounts, int *displs,
1265                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1266 {
1267   int retval;
1268
1269   smpi_bench_end();
1270 #ifdef HAVE_TRACING
1271   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1272   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1273   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1274 #endif
1275   if (comm == MPI_COMM_NULL) {
1276     retval = MPI_ERR_COMM;
1277   } else if (sendtype == MPI_DATATYPE_NULL
1278              || recvtype == MPI_DATATYPE_NULL) {
1279     retval = MPI_ERR_TYPE;
1280   } else if (recvcounts == NULL || displs == NULL) {
1281     retval = MPI_ERR_ARG;
1282   } else {
1283     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1284                      displs, recvtype, root, comm);
1285     retval = MPI_SUCCESS;
1286   }
1287 #ifdef HAVE_TRACING
1288   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1289 #endif
1290   smpi_bench_begin();
1291   return retval;
1292 }
1293
1294 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1295                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1296                   MPI_Comm comm)
1297 {
1298   int retval;
1299
1300   smpi_bench_end();
1301 #ifdef HAVE_TRACING
1302   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1303   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1304 #endif
1305   if (comm == MPI_COMM_NULL) {
1306     retval = MPI_ERR_COMM;
1307   } else if (sendtype == MPI_DATATYPE_NULL
1308              || recvtype == MPI_DATATYPE_NULL) {
1309     retval = MPI_ERR_TYPE;
1310   } else {
1311     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1312                        recvtype, comm);
1313     retval = MPI_SUCCESS;
1314   }
1315 #ifdef HAVE_TRACING
1316   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1317 #endif
1318   smpi_bench_begin();
1319   return retval;
1320 }
1321
1322 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1323                    void *recvbuf, int *recvcounts, int *displs,
1324                    MPI_Datatype recvtype, MPI_Comm comm)
1325 {
1326   int retval;
1327
1328   smpi_bench_end();
1329 #ifdef HAVE_TRACING
1330   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1331   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1332 #endif
1333   if (comm == MPI_COMM_NULL) {
1334     retval = MPI_ERR_COMM;
1335   } else if (sendtype == MPI_DATATYPE_NULL
1336              || recvtype == MPI_DATATYPE_NULL) {
1337     retval = MPI_ERR_TYPE;
1338   } else if (recvcounts == NULL || displs == NULL) {
1339     retval = MPI_ERR_ARG;
1340   } else {
1341     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1342                         displs, recvtype, comm);
1343     retval = MPI_SUCCESS;
1344   }
1345 #ifdef HAVE_TRACING
1346   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1347 #endif
1348   smpi_bench_begin();
1349   return retval;
1350 }
1351
1352 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1353                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1354                 int root, MPI_Comm comm)
1355 {
1356   int retval;
1357
1358   smpi_bench_end();
1359 #ifdef HAVE_TRACING
1360   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1361   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1362   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1363 #endif
1364   if (comm == MPI_COMM_NULL) {
1365     retval = MPI_ERR_COMM;
1366   } else if (sendtype == MPI_DATATYPE_NULL
1367              || recvtype == MPI_DATATYPE_NULL) {
1368     retval = MPI_ERR_TYPE;
1369   } else {
1370     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1371                      recvtype, root, comm);
1372     retval = MPI_SUCCESS;
1373   }
1374 #ifdef HAVE_TRACING
1375   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1376 #endif
1377   smpi_bench_begin();
1378   return retval;
1379 }
1380
1381 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1382                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1383                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1384 {
1385   int retval;
1386
1387   smpi_bench_end();
1388 #ifdef HAVE_TRACING
1389   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1390   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1391   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1392 #endif
1393   if (comm == MPI_COMM_NULL) {
1394     retval = MPI_ERR_COMM;
1395   } else if (sendtype == MPI_DATATYPE_NULL
1396              || recvtype == MPI_DATATYPE_NULL) {
1397     retval = MPI_ERR_TYPE;
1398   } else if (sendcounts == NULL || displs == NULL) {
1399     retval = MPI_ERR_ARG;
1400   } else {
1401     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1402                       recvcount, recvtype, root, comm);
1403     retval = MPI_SUCCESS;
1404   }
1405 #ifdef HAVE_TRACING
1406   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1407 #endif
1408   smpi_bench_begin();
1409   return retval;
1410 }
1411
1412 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1413                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1414 {
1415   int retval;
1416
1417   smpi_bench_end();
1418 #ifdef HAVE_TRACING
1419   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1420   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1421   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1422 #endif
1423   if (comm == MPI_COMM_NULL) {
1424     retval = MPI_ERR_COMM;
1425   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1426     retval = MPI_ERR_ARG;
1427   } else {
1428     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1429     retval = MPI_SUCCESS;
1430   }
1431 #ifdef HAVE_TRACING
1432   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1433 #endif
1434   smpi_bench_begin();
1435   return retval;
1436 }
1437
1438 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1439                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1440 {
1441   int retval;
1442
1443   smpi_bench_end();
1444 #ifdef HAVE_TRACING
1445   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1446   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1447 #endif
1448   if (comm == MPI_COMM_NULL) {
1449     retval = MPI_ERR_COMM;
1450   } else if (datatype == MPI_DATATYPE_NULL) {
1451     retval = MPI_ERR_TYPE;
1452   } else if (op == MPI_OP_NULL) {
1453     retval = MPI_ERR_OP;
1454   } else {
1455     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
1456     retval = MPI_SUCCESS;
1457   }
1458 #ifdef HAVE_TRACING
1459   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1460 #endif
1461   smpi_bench_begin();
1462   return retval;
1463 }
1464
1465 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
1466              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1467 {
1468   int retval;
1469
1470   smpi_bench_end();
1471 #ifdef HAVE_TRACING
1472   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1473   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1474 #endif
1475   if (comm == MPI_COMM_NULL) {
1476     retval = MPI_ERR_COMM;
1477   } else if (datatype == MPI_DATATYPE_NULL) {
1478     retval = MPI_ERR_TYPE;
1479   } else if (op == MPI_OP_NULL) {
1480     retval = MPI_ERR_OP;
1481   } else {
1482     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1483     retval = MPI_SUCCESS;
1484   }
1485 #ifdef HAVE_TRACING
1486   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1487 #endif
1488   smpi_bench_begin();
1489   return retval;
1490 }
1491
1492 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1493                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1494 {
1495   int retval, i, size, count;
1496   int *displs;
1497   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1498
1499   smpi_bench_end();
1500 #ifdef HAVE_TRACING
1501   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1502 #endif
1503   if (comm == MPI_COMM_NULL) {
1504     retval = MPI_ERR_COMM;
1505   } else if (datatype == MPI_DATATYPE_NULL) {
1506     retval = MPI_ERR_TYPE;
1507   } else if (op == MPI_OP_NULL) {
1508     retval = MPI_ERR_OP;
1509   } else if (recvcounts == NULL) {
1510     retval = MPI_ERR_ARG;
1511   } else {
1512     /* arbitrarily choose root as rank 0 */
1513     /* TODO: faster direct implementation ? */
1514     size = smpi_comm_size(comm);
1515     count = 0;
1516     displs = xbt_new(int, size);
1517     for (i = 0; i < size; i++) {
1518       count += recvcounts[i];
1519       displs[i] = 0;
1520     }
1521     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1522     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf,
1523                       recvcounts[rank], datatype, 0, comm);
1524     xbt_free(displs);
1525     retval = MPI_SUCCESS;
1526   }
1527 #ifdef HAVE_TRACING
1528   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1529 #endif
1530   smpi_bench_begin();
1531   return retval;
1532 }
1533
1534 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1535                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
1536                  MPI_Comm comm)
1537 {
1538   int retval, size, sendsize;
1539
1540   smpi_bench_end();
1541 #ifdef HAVE_TRACING
1542   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1543   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1544 #endif
1545   if (comm == MPI_COMM_NULL) {
1546     retval = MPI_ERR_COMM;
1547   } else if (sendtype == MPI_DATATYPE_NULL
1548              || recvtype == MPI_DATATYPE_NULL) {
1549     retval = MPI_ERR_TYPE;
1550   } else {
1551     size = smpi_comm_size(comm);
1552     sendsize = smpi_datatype_size(sendtype) * sendcount;
1553     if (sendsize < 200 && size > 12) {
1554       retval =
1555           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
1556                                          recvbuf, recvcount, recvtype,
1557                                          comm);
1558     } else if (sendsize < 3000) {
1559       retval =
1560           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
1561                                                 sendtype, recvbuf,
1562                                                 recvcount, recvtype, comm);
1563     } else {
1564       retval =
1565           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
1566                                             recvbuf, recvcount, recvtype,
1567                                             comm);
1568     }
1569   }
1570 #ifdef HAVE_TRACING
1571   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1572 #endif
1573   smpi_bench_begin();
1574   return retval;
1575 }
1576
1577 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
1578                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
1579                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
1580 {
1581   int retval;
1582
1583   smpi_bench_end();
1584 #ifdef HAVE_TRACING
1585   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1586   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1587 #endif
1588   if (comm == MPI_COMM_NULL) {
1589     retval = MPI_ERR_COMM;
1590   } else if (sendtype == MPI_DATATYPE_NULL
1591              || recvtype == MPI_DATATYPE_NULL) {
1592     retval = MPI_ERR_TYPE;
1593   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
1594              || recvdisps == NULL) {
1595     retval = MPI_ERR_ARG;
1596   } else {
1597     retval =
1598         smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
1599                                   recvbuf, recvcounts, recvdisps, recvtype,
1600                                   comm);
1601   }
1602 #ifdef HAVE_TRACING
1603   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1604 #endif
1605   smpi_bench_begin();
1606   return retval;
1607 }
1608
1609
1610 int PMPI_Get_processor_name(char *name, int *resultlen)
1611 {
1612   int retval = MPI_SUCCESS;
1613
1614   smpi_bench_end();
1615   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
1616           MPI_MAX_PROCESSOR_NAME - 1);
1617   *resultlen =
1618       strlen(name) >
1619       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1620
1621   smpi_bench_begin();
1622   return retval;
1623 }
1624
1625 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
1626 {
1627   int retval = MPI_SUCCESS;
1628   size_t size;
1629
1630   smpi_bench_end();
1631   if (status == NULL || count == NULL) {
1632     retval = MPI_ERR_ARG;
1633   } else if (datatype == MPI_DATATYPE_NULL) {
1634     retval = MPI_ERR_TYPE;
1635   } else {
1636     size = smpi_datatype_size(datatype);
1637     if (size == 0) {
1638       *count = 0;
1639     } else if (status->count % size != 0) {
1640       retval = MPI_UNDEFINED;
1641     } else {
1642       *count = smpi_mpi_get_count(status, datatype);
1643     }
1644   }
1645   smpi_bench_begin();
1646   return retval;
1647 }
1648
1649 /* The following calls are not yet implemented and will fail at runtime. */
1650 /* Once implemented, please move them above this notice. */
1651
1652 static int not_yet_implemented(void) {
1653    xbt_die("Not yet implemented");
1654    return MPI_ERR_OTHER;
1655 }
1656
1657 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
1658    return not_yet_implemented();
1659 }
1660
1661 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
1662    return not_yet_implemented();
1663 }
1664
1665 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
1666    return not_yet_implemented();
1667 }
1668
1669 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
1670    return not_yet_implemented();
1671 }
1672
1673 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
1674    return not_yet_implemented();
1675 }
1676
1677 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
1678    return not_yet_implemented();
1679 }
1680
1681 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
1682    return not_yet_implemented();
1683 }
1684
1685 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
1686    return not_yet_implemented();
1687 }
1688
1689 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
1690    return not_yet_implemented();
1691 }
1692
1693 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
1694    return not_yet_implemented();
1695 }
1696
1697 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
1698    return not_yet_implemented();
1699 }
1700
1701 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
1702    return not_yet_implemented();
1703 }
1704
1705 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
1706    return not_yet_implemented();
1707 }
1708
1709 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
1710    return not_yet_implemented();
1711 }
1712
1713 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
1714    return not_yet_implemented();
1715 }
1716
1717 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
1718    return not_yet_implemented();
1719 }
1720
1721 int PMPI_Error_class(int errorcode, int* errorclass) {
1722    return not_yet_implemented();
1723 }
1724
1725 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
1726    return not_yet_implemented();
1727 }
1728
1729 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
1730    return not_yet_implemented();
1731 }
1732
1733 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
1734    return not_yet_implemented();
1735 }
1736
1737 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
1738    return not_yet_implemented();
1739 }
1740
1741 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
1742    return not_yet_implemented();
1743 }
1744
1745 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* newtype) {
1746    return not_yet_implemented();
1747 }
1748
1749 int PMPI_Cancel(MPI_Request* request) {
1750    return not_yet_implemented();
1751 }
1752
1753 int PMPI_Buffer_attach(void* buffer, int size) {
1754    return not_yet_implemented();
1755 }
1756
1757 int PMPI_Buffer_detach(void* buffer, int* size) {
1758    return not_yet_implemented();
1759 }
1760
1761 int PMPI_Testsome(int incount, MPI_Request* requests, int* outcount, int* indices, MPI_Status* statuses) {
1762    return not_yet_implemented();
1763 }
1764
1765 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
1766    return not_yet_implemented();
1767 }
1768
1769 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
1770    return not_yet_implemented();
1771 }
1772
1773 int PMPI_Type_commit(MPI_Datatype* datatype) {
1774    return not_yet_implemented();
1775 }
1776
1777 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
1778    return not_yet_implemented();
1779 }
1780
1781 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
1782    return not_yet_implemented();
1783 }
1784
1785 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
1786    return not_yet_implemented();
1787 }
1788
1789 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* newtype) {
1790    return not_yet_implemented();
1791 }
1792
1793 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
1794    return not_yet_implemented();
1795 }
1796
1797 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1798    return not_yet_implemented();
1799 }
1800
1801 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1802    return not_yet_implemented();
1803 }
1804
1805 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
1806    return not_yet_implemented();
1807 }
1808
1809 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
1810    return not_yet_implemented();
1811 }
1812
1813 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1814    return not_yet_implemented();
1815 }
1816
1817 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1818    return not_yet_implemented();
1819 }
1820
1821 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1822    return not_yet_implemented();
1823 }
1824
1825 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
1826    return not_yet_implemented();
1827 }
1828
1829 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
1830    return not_yet_implemented();
1831 }
1832
1833 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1834    return not_yet_implemented();
1835 }
1836
1837 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1838    return not_yet_implemented();
1839 }
1840
1841 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
1842    return not_yet_implemented();
1843 }
1844
1845 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
1846    return not_yet_implemented();
1847 }
1848
1849 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
1850    return not_yet_implemented();
1851 }
1852
1853 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1854    return not_yet_implemented();
1855 }
1856
1857 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1858    return not_yet_implemented();
1859 }
1860
1861 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1862    return not_yet_implemented();
1863 }
1864
1865 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
1866    return not_yet_implemented();
1867 }
1868
1869 int PMPI_Keyval_free(int* keyval) {
1870    return not_yet_implemented();
1871 }
1872
1873 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
1874    return not_yet_implemented();
1875 }
1876
1877 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
1878    return not_yet_implemented();
1879 }
1880
1881 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses) {
1882    return not_yet_implemented();
1883 }
1884
1885 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
1886    return not_yet_implemented();
1887 }
1888
1889 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
1890    return not_yet_implemented();
1891 }
1892
1893 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1894    return not_yet_implemented();
1895 }
1896
1897 int PMPI_Initialized(int* flag) {
1898    return not_yet_implemented();
1899 }