Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
5079400097d06e15fd442c7decd50d8966770352
[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   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
976   TRACE_smpi_recv(rank, src_traced, rank);
977   TRACE_smpi_computing_in(rank);
978 #endif
979   smpi_bench_begin();
980   return retval;
981 }
982
983 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
984              MPI_Comm comm)
985 {
986   int retval;
987
988   smpi_bench_end();
989 #ifdef HAVE_TRACING
990   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
991   TRACE_smpi_computing_out(rank);
992   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
993   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
994   TRACE_smpi_send(rank, rank, dst_traced);
995 #endif
996   if (comm == MPI_COMM_NULL) {
997     retval = MPI_ERR_COMM;
998   } else {
999     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1000     retval = MPI_SUCCESS;
1001   }
1002 #ifdef HAVE_TRACING
1003   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1004   TRACE_smpi_computing_in(rank);
1005 #endif
1006   smpi_bench_begin();
1007   return retval;
1008 }
1009
1010 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1011                  int dst, int sendtag, void *recvbuf, int recvcount,
1012                  MPI_Datatype recvtype, int src, int recvtag,
1013                  MPI_Comm comm, MPI_Status * status)
1014 {
1015   int retval;
1016
1017   smpi_bench_end();
1018 #ifdef HAVE_TRACING
1019   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1020   TRACE_smpi_computing_out(rank);
1021   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
1022   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
1023   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1024   TRACE_smpi_send(rank, rank, dst_traced);
1025   TRACE_smpi_send(rank, src_traced, rank);
1026 #endif
1027   if (comm == MPI_COMM_NULL) {
1028     retval = MPI_ERR_COMM;
1029   } else if (sendtype == MPI_DATATYPE_NULL
1030              || recvtype == MPI_DATATYPE_NULL) {
1031     retval = MPI_ERR_TYPE;
1032   } else {
1033     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1034                       recvcount, recvtype, src, recvtag, comm, status);
1035     retval = MPI_SUCCESS;
1036   }
1037 #ifdef HAVE_TRACING
1038   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1039   TRACE_smpi_recv(rank, rank, dst_traced);
1040   TRACE_smpi_recv(rank, src_traced, rank);
1041   TRACE_smpi_computing_in(rank);
1042
1043 #endif
1044   smpi_bench_begin();
1045   return retval;
1046 }
1047
1048 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1049                          int dst, int sendtag, int src, int recvtag,
1050                          MPI_Comm comm, MPI_Status * status)
1051 {
1052   //TODO: suboptimal implementation
1053   void *recvbuf;
1054   int retval, size;
1055
1056   size = smpi_datatype_size(datatype) * count;
1057   recvbuf = xbt_new(char, size);
1058   retval =
1059       MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1060                    datatype, src, recvtag, comm, status);
1061   memcpy(buf, recvbuf, size * sizeof(char));
1062   xbt_free(recvbuf);
1063   return retval;
1064 }
1065
1066 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1067 {
1068   int retval;
1069
1070   smpi_bench_end();
1071   if (request == NULL || flag == NULL) {
1072     retval = MPI_ERR_ARG;
1073   } else if (*request == MPI_REQUEST_NULL) {
1074     retval = MPI_ERR_REQUEST;
1075   } else {
1076     *flag = smpi_mpi_test(request, status);
1077     retval = MPI_SUCCESS;
1078   }
1079   smpi_bench_begin();
1080   return retval;
1081 }
1082
1083 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1084                 MPI_Status * status)
1085 {
1086   int retval;
1087
1088   smpi_bench_end();
1089   if (index == NULL || flag == NULL) {
1090     retval = MPI_ERR_ARG;
1091   } else {
1092     *flag = smpi_mpi_testany(count, requests, index, status);
1093     retval = MPI_SUCCESS;
1094   }
1095   smpi_bench_begin();
1096   return retval;
1097 }
1098
1099
1100
1101 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1102   int retval;
1103   smpi_bench_end();
1104
1105   if (status == NULL) {
1106      retval = MPI_ERR_ARG;
1107   }else if (comm == MPI_COMM_NULL) {
1108         retval = MPI_ERR_COMM;
1109   } else {
1110         smpi_mpi_probe(source, tag, comm, status);
1111         retval = MPI_SUCCESS;
1112   }
1113   smpi_bench_begin();
1114   return retval;
1115 }
1116
1117
1118 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1119   int retval;
1120   smpi_bench_end();
1121
1122   if (flag == NULL) {
1123        retval = MPI_ERR_ARG;
1124     }else if (status == NULL) {
1125      retval = MPI_ERR_ARG;
1126   }else if (comm == MPI_COMM_NULL) {
1127         retval = MPI_ERR_COMM;
1128   } else {
1129         smpi_mpi_iprobe(source, tag, comm, flag, status);
1130         retval = MPI_SUCCESS;
1131   }
1132   smpi_bench_begin();
1133   return retval;
1134 }
1135
1136 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1137 {
1138   int retval;
1139
1140   smpi_bench_end();
1141 #ifdef HAVE_TRACING
1142   int rank = request && (*request)->comm != MPI_COMM_NULL
1143       ? smpi_comm_rank((*request)->comm)
1144       : -1;
1145   TRACE_smpi_computing_out(rank);
1146
1147   MPI_Group group = smpi_comm_group((*request)->comm);
1148   int src_traced = smpi_group_rank(group, (*request)->src);
1149   int dst_traced = smpi_group_rank(group, (*request)->dst);
1150   int is_wait_for_receive = (*request)->recv;
1151   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1152 #endif
1153   if (request == NULL) {
1154     retval = MPI_ERR_ARG;
1155   } else if (*request == MPI_REQUEST_NULL) {
1156     retval = MPI_ERR_REQUEST;
1157   } else {
1158     smpi_mpi_wait(request, status);
1159     retval = MPI_SUCCESS;
1160   }
1161 #ifdef HAVE_TRACING
1162   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1163   if (is_wait_for_receive) {
1164     TRACE_smpi_recv(rank, src_traced, dst_traced);
1165   }
1166   TRACE_smpi_computing_in(rank);
1167
1168 #endif
1169   smpi_bench_begin();
1170   return retval;
1171 }
1172
1173 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1174 {
1175   int retval;
1176
1177   smpi_bench_end();
1178 #ifdef HAVE_TRACING
1179   //save requests information for tracing
1180   int i;
1181   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), NULL);
1182   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), NULL);
1183   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), NULL);
1184   for (i = 0; i < count; i++) {
1185     MPI_Request req = requests[i];      //already received requests are no longer valid
1186     if (req) {
1187       int *asrc = xbt_new(int, 1);
1188       int *adst = xbt_new(int, 1);
1189       int *arecv = xbt_new(int, 1);
1190       *asrc = req->src;
1191       *adst = req->dst;
1192       *arecv = req->recv;
1193       xbt_dynar_insert_at(srcs, i, asrc);
1194       xbt_dynar_insert_at(dsts, i, adst);
1195       xbt_dynar_insert_at(recvs, i, arecv);
1196       xbt_free(asrc);
1197       xbt_free(adst);
1198       xbt_free(arecv);
1199     } else {
1200       int *t = xbt_new(int, 1);
1201       xbt_dynar_insert_at(srcs, i, t);
1202       xbt_dynar_insert_at(dsts, i, t);
1203       xbt_dynar_insert_at(recvs, i, t);
1204       xbt_free(t);
1205     }
1206   }
1207   int rank_traced = smpi_comm_rank(MPI_COMM_WORLD);
1208   TRACE_smpi_computing_out(rank_traced);
1209
1210   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1211
1212 #endif
1213   if (index == NULL) {
1214     retval = MPI_ERR_ARG;
1215   } else {
1216     *index = smpi_mpi_waitany(count, requests, status);
1217     retval = MPI_SUCCESS;
1218   }
1219 #ifdef HAVE_TRACING
1220   int src_traced, dst_traced, is_wait_for_receive;
1221   xbt_dynar_get_cpy(srcs, *index, &src_traced);
1222   xbt_dynar_get_cpy(dsts, *index, &dst_traced);
1223   xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
1224   if (is_wait_for_receive) {
1225     TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1226   }
1227   TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1228   //clean-up of dynars
1229   xbt_dynar_free(&srcs);
1230   xbt_dynar_free(&dsts);
1231   xbt_dynar_free(&recvs);
1232   TRACE_smpi_computing_in(rank_traced);
1233
1234 #endif
1235   smpi_bench_begin();
1236   return retval;
1237 }
1238
1239 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1240 {
1241
1242   smpi_bench_end();
1243 #ifdef HAVE_TRACING
1244   //save information from requests
1245   int i;
1246   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), NULL);
1247   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), NULL);
1248   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), NULL);
1249   for (i = 0; i < count; i++) {
1250     MPI_Request req = requests[i];      //all req should be valid in Waitall
1251     int *asrc = xbt_new(int, 1);
1252     int *adst = xbt_new(int, 1);
1253     int *arecv = xbt_new(int, 1);
1254     *asrc = req->src;
1255     *adst = req->dst;
1256     *arecv = req->recv;
1257     xbt_dynar_insert_at(srcs, i, asrc);
1258     xbt_dynar_insert_at(dsts, i, adst);
1259     xbt_dynar_insert_at(recvs, i, arecv);
1260     xbt_free(asrc);
1261     xbt_free(adst);
1262     xbt_free(arecv);
1263   }
1264   int rank_traced = smpi_comm_rank (MPI_COMM_WORLD);
1265   TRACE_smpi_computing_out(rank_traced);
1266
1267   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1268 #endif
1269   smpi_mpi_waitall(count, requests, status);
1270 #ifdef HAVE_TRACING
1271   for (i = 0; i < count; i++) {
1272     int src_traced, dst_traced, is_wait_for_receive;
1273     xbt_dynar_get_cpy(srcs, i, &src_traced);
1274     xbt_dynar_get_cpy(dsts, i, &dst_traced);
1275     xbt_dynar_get_cpy(recvs, i, &is_wait_for_receive);
1276     if (is_wait_for_receive) {
1277       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1278     }
1279   }
1280   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1281   //clean-up of dynars
1282   xbt_dynar_free(&srcs);
1283   xbt_dynar_free(&dsts);
1284   xbt_dynar_free(&recvs);
1285   TRACE_smpi_computing_in(rank_traced);
1286 #endif
1287   smpi_bench_begin();
1288   return MPI_SUCCESS;
1289 }
1290
1291 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1292                  int *indices, MPI_Status status[])
1293 {
1294   int retval;
1295
1296   smpi_bench_end();
1297   if (outcount == NULL || indices == NULL) {
1298     retval = MPI_ERR_ARG;
1299   } else {
1300     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1301     retval = MPI_SUCCESS;
1302   }
1303   smpi_bench_begin();
1304   return retval;
1305 }
1306
1307 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1308 {
1309   int retval;
1310
1311   smpi_bench_end();
1312 #ifdef HAVE_TRACING
1313   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1314   TRACE_smpi_computing_out(rank);
1315   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1316   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1317 #endif
1318   if (comm == MPI_COMM_NULL) {
1319     retval = MPI_ERR_COMM;
1320   } else {
1321     smpi_mpi_bcast(buf, count, datatype, root, comm);
1322     retval = MPI_SUCCESS;
1323   }
1324 #ifdef HAVE_TRACING
1325   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1326   TRACE_smpi_computing_in(rank);
1327 #endif
1328   smpi_bench_begin();
1329   return retval;
1330 }
1331
1332 int PMPI_Barrier(MPI_Comm comm)
1333 {
1334   int retval;
1335
1336   smpi_bench_end();
1337 #ifdef HAVE_TRACING
1338   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1339   TRACE_smpi_computing_out(rank);
1340   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1341 #endif
1342   if (comm == MPI_COMM_NULL) {
1343     retval = MPI_ERR_COMM;
1344   } else {
1345     smpi_mpi_barrier(comm);
1346     retval = MPI_SUCCESS;
1347   }
1348 #ifdef HAVE_TRACING
1349   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1350   TRACE_smpi_computing_in(rank);
1351 #endif
1352   smpi_bench_begin();
1353   return retval;
1354 }
1355
1356 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1357                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1358                int root, MPI_Comm comm)
1359 {
1360   int retval;
1361
1362   smpi_bench_end();
1363 #ifdef HAVE_TRACING
1364   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1365   TRACE_smpi_computing_out(rank);
1366   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1367   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1368 #endif
1369   if (comm == MPI_COMM_NULL) {
1370     retval = MPI_ERR_COMM;
1371   } else if (sendtype == MPI_DATATYPE_NULL
1372              || recvtype == MPI_DATATYPE_NULL) {
1373     retval = MPI_ERR_TYPE;
1374   } else {
1375     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1376                     recvtype, root, comm);
1377     retval = MPI_SUCCESS;
1378   }
1379 #ifdef HAVE_TRACING
1380   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1381   TRACE_smpi_computing_in(rank);
1382 #endif
1383   smpi_bench_begin();
1384   return retval;
1385 }
1386
1387 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1388                 void *recvbuf, int *recvcounts, int *displs,
1389                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1390 {
1391   int retval;
1392
1393   smpi_bench_end();
1394 #ifdef HAVE_TRACING
1395   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1396   TRACE_smpi_computing_out(rank);
1397   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1398   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1399 #endif
1400   if (comm == MPI_COMM_NULL) {
1401     retval = MPI_ERR_COMM;
1402   } else if (sendtype == MPI_DATATYPE_NULL
1403              || recvtype == MPI_DATATYPE_NULL) {
1404     retval = MPI_ERR_TYPE;
1405   } else if (recvcounts == NULL || displs == NULL) {
1406     retval = MPI_ERR_ARG;
1407   } else {
1408     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1409                      displs, recvtype, root, comm);
1410     retval = MPI_SUCCESS;
1411   }
1412 #ifdef HAVE_TRACING
1413   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1414   TRACE_smpi_computing_in(rank);
1415 #endif
1416   smpi_bench_begin();
1417   return retval;
1418 }
1419
1420 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1421                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1422                   MPI_Comm comm)
1423 {
1424   int retval;
1425
1426   smpi_bench_end();
1427 #ifdef HAVE_TRACING
1428   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1429   TRACE_smpi_computing_out(rank);
1430   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1431 #endif
1432   if (comm == MPI_COMM_NULL) {
1433     retval = MPI_ERR_COMM;
1434   } else if (sendtype == MPI_DATATYPE_NULL
1435              || recvtype == MPI_DATATYPE_NULL) {
1436     retval = MPI_ERR_TYPE;
1437   } else {
1438     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1439                        recvtype, comm);
1440     retval = MPI_SUCCESS;
1441   }
1442 #ifdef HAVE_TRACING
1443   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1444 #endif
1445   smpi_bench_begin();
1446   return retval;
1447 }
1448
1449 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1450                    void *recvbuf, int *recvcounts, int *displs,
1451                    MPI_Datatype recvtype, MPI_Comm comm)
1452 {
1453   int retval;
1454
1455   smpi_bench_end();
1456 #ifdef HAVE_TRACING
1457   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1458   TRACE_smpi_computing_out(rank);
1459   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1460 #endif
1461   if (comm == MPI_COMM_NULL) {
1462     retval = MPI_ERR_COMM;
1463   } else if (sendtype == MPI_DATATYPE_NULL
1464              || recvtype == MPI_DATATYPE_NULL) {
1465     retval = MPI_ERR_TYPE;
1466   } else if (recvcounts == NULL || displs == NULL) {
1467     retval = MPI_ERR_ARG;
1468   } else {
1469     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1470                         displs, recvtype, comm);
1471     retval = MPI_SUCCESS;
1472   }
1473 #ifdef HAVE_TRACING
1474   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1475   TRACE_smpi_computing_in(rank);
1476 #endif
1477   smpi_bench_begin();
1478   return retval;
1479 }
1480
1481 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1482                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1483                 int root, MPI_Comm comm)
1484 {
1485   int retval;
1486
1487   smpi_bench_end();
1488 #ifdef HAVE_TRACING
1489   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1490   TRACE_smpi_computing_out(rank);
1491   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1492   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1493 #endif
1494   if (comm == MPI_COMM_NULL) {
1495     retval = MPI_ERR_COMM;
1496   } else if (sendtype == MPI_DATATYPE_NULL
1497              || recvtype == MPI_DATATYPE_NULL) {
1498     retval = MPI_ERR_TYPE;
1499   } else {
1500     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1501                      recvtype, root, comm);
1502     retval = MPI_SUCCESS;
1503   }
1504 #ifdef HAVE_TRACING
1505   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1506   TRACE_smpi_computing_in(rank);
1507 #endif
1508   smpi_bench_begin();
1509   return retval;
1510 }
1511
1512 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1513                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1514                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1515 {
1516   int retval;
1517
1518   smpi_bench_end();
1519 #ifdef HAVE_TRACING
1520   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1521   TRACE_smpi_computing_out(rank);
1522   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1523   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1524 #endif
1525   if (comm == MPI_COMM_NULL) {
1526     retval = MPI_ERR_COMM;
1527   } else if (sendtype == MPI_DATATYPE_NULL
1528              || recvtype == MPI_DATATYPE_NULL) {
1529     retval = MPI_ERR_TYPE;
1530   } else if (sendcounts == NULL || displs == NULL) {
1531     retval = MPI_ERR_ARG;
1532   } else {
1533     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1534                       recvcount, recvtype, root, comm);
1535     retval = MPI_SUCCESS;
1536   }
1537 #ifdef HAVE_TRACING
1538   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1539   TRACE_smpi_computing_in(rank);
1540 #endif
1541   smpi_bench_begin();
1542   return retval;
1543 }
1544
1545 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1546                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1547 {
1548   int retval;
1549
1550   smpi_bench_end();
1551 #ifdef HAVE_TRACING
1552   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1553   TRACE_smpi_computing_out(rank);
1554   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1555   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1556 #endif
1557   if (comm == MPI_COMM_NULL) {
1558     retval = MPI_ERR_COMM;
1559   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1560     retval = MPI_ERR_ARG;
1561   } else {
1562     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1563     retval = MPI_SUCCESS;
1564   }
1565 #ifdef HAVE_TRACING
1566   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1567   TRACE_smpi_computing_in(rank);
1568 #endif
1569   smpi_bench_begin();
1570   return retval;
1571 }
1572
1573 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1574                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1575 {
1576   int retval;
1577
1578   smpi_bench_end();
1579 #ifdef HAVE_TRACING
1580   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1581   TRACE_smpi_computing_out(rank);
1582   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1583 #endif
1584   if (comm == MPI_COMM_NULL) {
1585     retval = MPI_ERR_COMM;
1586   } else if (datatype == MPI_DATATYPE_NULL) {
1587     retval = MPI_ERR_TYPE;
1588   } else if (op == MPI_OP_NULL) {
1589     retval = MPI_ERR_OP;
1590   } else {
1591     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
1592     retval = MPI_SUCCESS;
1593   }
1594 #ifdef HAVE_TRACING
1595   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1596   TRACE_smpi_computing_in(rank);
1597 #endif
1598   smpi_bench_begin();
1599   return retval;
1600 }
1601
1602 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
1603              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1604 {
1605   int retval;
1606
1607   smpi_bench_end();
1608 #ifdef HAVE_TRACING
1609   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1610   TRACE_smpi_computing_out(rank);
1611   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1612 #endif
1613   if (comm == MPI_COMM_NULL) {
1614     retval = MPI_ERR_COMM;
1615   } else if (datatype == MPI_DATATYPE_NULL) {
1616     retval = MPI_ERR_TYPE;
1617   } else if (op == MPI_OP_NULL) {
1618     retval = MPI_ERR_OP;
1619   } else {
1620     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1621     retval = MPI_SUCCESS;
1622   }
1623 #ifdef HAVE_TRACING
1624   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1625   TRACE_smpi_computing_in(rank);
1626 #endif
1627   smpi_bench_begin();
1628   return retval;
1629 }
1630
1631 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1632                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1633 {
1634   int retval, i, size, count;
1635   int *displs;
1636   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1637
1638   smpi_bench_end();
1639 #ifdef HAVE_TRACING
1640   TRACE_smpi_computing_out(rank);
1641   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1642 #endif
1643   if (comm == MPI_COMM_NULL) {
1644     retval = MPI_ERR_COMM;
1645   } else if (datatype == MPI_DATATYPE_NULL) {
1646     retval = MPI_ERR_TYPE;
1647   } else if (op == MPI_OP_NULL) {
1648     retval = MPI_ERR_OP;
1649   } else if (recvcounts == NULL) {
1650     retval = MPI_ERR_ARG;
1651   } else {
1652     /* arbitrarily choose root as rank 0 */
1653     /* TODO: faster direct implementation ? */
1654     size = smpi_comm_size(comm);
1655     count = 0;
1656     displs = xbt_new(int, size);
1657     for (i = 0; i < size; i++) {
1658       count += recvcounts[i];
1659       displs[i] = 0;
1660     }
1661     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1662     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf,
1663                       recvcounts[rank], datatype, 0, comm);
1664     xbt_free(displs);
1665     retval = MPI_SUCCESS;
1666   }
1667 #ifdef HAVE_TRACING
1668   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1669   TRACE_smpi_computing_in(rank);
1670 #endif
1671   smpi_bench_begin();
1672   return retval;
1673 }
1674
1675 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1676                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
1677                  MPI_Comm comm)
1678 {
1679   int retval, size, sendsize;
1680
1681   smpi_bench_end();
1682 #ifdef HAVE_TRACING
1683   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1684   TRACE_smpi_computing_out(rank);
1685   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1686 #endif
1687   if (comm == MPI_COMM_NULL) {
1688     retval = MPI_ERR_COMM;
1689   } else if (sendtype == MPI_DATATYPE_NULL
1690              || recvtype == MPI_DATATYPE_NULL) {
1691     retval = MPI_ERR_TYPE;
1692   } else {
1693     size = smpi_comm_size(comm);
1694     sendsize = smpi_datatype_size(sendtype) * sendcount;
1695     if (sendsize < 200 && size > 12) {
1696       retval =
1697           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
1698                                          recvbuf, recvcount, recvtype,
1699                                          comm);
1700     } else if (sendsize < 3000) {
1701       retval =
1702           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
1703                                                 sendtype, recvbuf,
1704                                                 recvcount, recvtype, comm);
1705     } else {
1706       retval =
1707           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
1708                                             recvbuf, recvcount, recvtype,
1709                                             comm);
1710     }
1711   }
1712 #ifdef HAVE_TRACING
1713   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1714   TRACE_smpi_computing_in(rank);
1715 #endif
1716   smpi_bench_begin();
1717   return retval;
1718 }
1719
1720 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
1721                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
1722                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
1723 {
1724   int retval;
1725
1726   smpi_bench_end();
1727 #ifdef HAVE_TRACING
1728   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1729   TRACE_smpi_computing_out(rank);
1730   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1731 #endif
1732   if (comm == MPI_COMM_NULL) {
1733     retval = MPI_ERR_COMM;
1734   } else if (sendtype == MPI_DATATYPE_NULL
1735              || recvtype == MPI_DATATYPE_NULL) {
1736     retval = MPI_ERR_TYPE;
1737   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
1738              || recvdisps == NULL) {
1739     retval = MPI_ERR_ARG;
1740   } else {
1741     retval =
1742         smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
1743                                   recvbuf, recvcounts, recvdisps, recvtype,
1744                                   comm);
1745   }
1746 #ifdef HAVE_TRACING
1747   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1748   TRACE_smpi_computing_in(rank);
1749 #endif
1750   smpi_bench_begin();
1751   return retval;
1752 }
1753
1754
1755 int PMPI_Get_processor_name(char *name, int *resultlen)
1756 {
1757   int retval = MPI_SUCCESS;
1758
1759   smpi_bench_end();
1760   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
1761           MPI_MAX_PROCESSOR_NAME - 1);
1762   *resultlen =
1763       strlen(name) >
1764       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1765
1766   smpi_bench_begin();
1767   return retval;
1768 }
1769
1770 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
1771 {
1772   int retval = MPI_SUCCESS;
1773   size_t size;
1774
1775   smpi_bench_end();
1776   if (status == NULL || count == NULL) {
1777     retval = MPI_ERR_ARG;
1778   } else if (datatype == MPI_DATATYPE_NULL) {
1779     retval = MPI_ERR_TYPE;
1780   } else {
1781     size = smpi_datatype_size(datatype);
1782     if (size == 0) {
1783       *count = 0;
1784     } else if (status->count % size != 0) {
1785       retval = MPI_UNDEFINED;
1786     } else {
1787       *count = smpi_mpi_get_count(status, datatype);
1788     }
1789   }
1790   smpi_bench_begin();
1791   return retval;
1792 }
1793
1794 /* The following calls are not yet implemented and will fail at runtime. */
1795 /* Once implemented, please move them above this notice. */
1796
1797 static int not_yet_implemented(void) {
1798    xbt_die("Not yet implemented");
1799    return MPI_ERR_OTHER;
1800 }
1801
1802 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
1803    return not_yet_implemented();
1804 }
1805
1806 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
1807    return not_yet_implemented();
1808 }
1809
1810 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
1811    return not_yet_implemented();
1812 }
1813
1814 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
1815    return not_yet_implemented();
1816 }
1817
1818 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
1819    return not_yet_implemented();
1820 }
1821
1822 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
1823    return not_yet_implemented();
1824 }
1825
1826 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
1827    return not_yet_implemented();
1828 }
1829
1830 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
1831    return not_yet_implemented();
1832 }
1833
1834 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
1835    return not_yet_implemented();
1836 }
1837
1838 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
1839    return not_yet_implemented();
1840 }
1841
1842 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
1843    return not_yet_implemented();
1844 }
1845
1846 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
1847    return not_yet_implemented();
1848 }
1849
1850 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
1851    return not_yet_implemented();
1852 }
1853
1854 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
1855    return not_yet_implemented();
1856 }
1857
1858 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
1859    return not_yet_implemented();
1860 }
1861
1862 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
1863    return not_yet_implemented();
1864 }
1865
1866 int PMPI_Error_class(int errorcode, int* errorclass) {
1867    return not_yet_implemented();
1868 }
1869
1870 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
1871    return not_yet_implemented();
1872 }
1873
1874 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
1875    return not_yet_implemented();
1876 }
1877
1878 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
1879    return not_yet_implemented();
1880 }
1881
1882 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
1883    return not_yet_implemented();
1884 }
1885
1886 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
1887    return not_yet_implemented();
1888 }
1889
1890 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* newtype) {
1891    return not_yet_implemented();
1892 }
1893
1894 int PMPI_Cancel(MPI_Request* request) {
1895    return not_yet_implemented();
1896 }
1897
1898 int PMPI_Buffer_attach(void* buffer, int size) {
1899    return not_yet_implemented();
1900 }
1901
1902 int PMPI_Buffer_detach(void* buffer, int* size) {
1903    return not_yet_implemented();
1904 }
1905
1906 int PMPI_Testsome(int incount, MPI_Request* requests, int* outcount, int* indices, MPI_Status* statuses) {
1907    return not_yet_implemented();
1908 }
1909
1910 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
1911    return not_yet_implemented();
1912 }
1913
1914 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
1915    return not_yet_implemented();
1916 }
1917
1918 int PMPI_Type_commit(MPI_Datatype* datatype) {
1919    return not_yet_implemented();
1920 }
1921
1922 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
1923    return not_yet_implemented();
1924 }
1925
1926 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
1927    return not_yet_implemented();
1928 }
1929
1930 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
1931    return not_yet_implemented();
1932 }
1933
1934 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* newtype) {
1935    return not_yet_implemented();
1936 }
1937
1938 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
1939    return not_yet_implemented();
1940 }
1941
1942 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1943    return not_yet_implemented();
1944 }
1945
1946 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1947    return not_yet_implemented();
1948 }
1949
1950 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
1951    return not_yet_implemented();
1952 }
1953
1954 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
1955    return not_yet_implemented();
1956 }
1957
1958 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1959    return not_yet_implemented();
1960 }
1961
1962 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1963    return not_yet_implemented();
1964 }
1965
1966 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1967    return not_yet_implemented();
1968 }
1969
1970 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
1971    return not_yet_implemented();
1972 }
1973
1974 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
1975    return not_yet_implemented();
1976 }
1977
1978 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1979    return not_yet_implemented();
1980 }
1981
1982
1983
1984 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
1985    return not_yet_implemented();
1986 }
1987
1988 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
1989    return not_yet_implemented();
1990 }
1991
1992 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
1993    return not_yet_implemented();
1994 }
1995
1996 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1997    return not_yet_implemented();
1998 }
1999
2000 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2001    return not_yet_implemented();
2002 }
2003
2004 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2005    return not_yet_implemented();
2006 }
2007
2008 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2009    return not_yet_implemented();
2010 }
2011
2012 int PMPI_Keyval_free(int* keyval) {
2013    return not_yet_implemented();
2014 }
2015
2016 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2017    return not_yet_implemented();
2018 }
2019
2020 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2021    return not_yet_implemented();
2022 }
2023
2024 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses) {
2025    return not_yet_implemented();
2026 }
2027
2028 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2029    return not_yet_implemented();
2030 }
2031
2032 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2033    return not_yet_implemented();
2034 }
2035
2036 int PMPI_Initialized(int* flag) {
2037    return not_yet_implemented();
2038 }