Logo AND Algorithmique Numérique Distribuée

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