Logo AND Algorithmique Numérique Distribuée

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