Logo AND Algorithmique Numérique Distribuée

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