Logo AND Algorithmique Numérique Distribuée

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