Logo AND Algorithmique Numérique Distribuée

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