Logo AND Algorithmique Numérique Distribuée

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