Logo AND Algorithmique Numérique Distribuée

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