Logo AND Algorithmique Numérique Distribuée

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