Logo AND Algorithmique Numérique Distribuée

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