Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Remove the fork based checkpoint functionality.
[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 {
626     *rank = smpi_comm_rank(comm);
627     retval = MPI_SUCCESS;
628   }
629   smpi_bench_begin();
630   return retval;
631 }
632
633 int PMPI_Comm_size(MPI_Comm comm, int *size)
634 {
635   int retval;
636
637   smpi_bench_end();
638   if (comm == MPI_COMM_NULL) {
639     retval = MPI_ERR_COMM;
640   } else if (size == NULL) {
641     retval = MPI_ERR_ARG;
642   } else {
643     *size = smpi_comm_size(comm);
644     retval = MPI_SUCCESS;
645   }
646   smpi_bench_begin();
647   return retval;
648 }
649
650 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
651 {
652   int retval;
653
654   smpi_bench_end();
655   if (comm == MPI_COMM_NULL) {
656     retval = MPI_ERR_COMM;
657   } else if (group == NULL) {
658     retval = MPI_ERR_ARG;
659   } else {
660     *group = smpi_comm_group(comm);
661     retval = MPI_SUCCESS;
662   }
663   smpi_bench_begin();
664   return retval;
665 }
666
667 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
668 {
669   int retval;
670
671   smpi_bench_end();
672   if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
673     retval = MPI_ERR_COMM;
674   } else if (result == NULL) {
675     retval = MPI_ERR_ARG;
676   } else {
677     if (comm1 == comm2) {       /* Same communicators means same groups */
678       *result = MPI_IDENT;
679     } else {
680       *result =
681           smpi_group_compare(smpi_comm_group(comm1),
682                              smpi_comm_group(comm2));
683       if (*result == MPI_IDENT) {
684         *result = MPI_CONGRUENT;
685       }
686     }
687     retval = MPI_SUCCESS;
688   }
689   smpi_bench_begin();
690   return retval;
691 }
692
693 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
694 {
695   int retval;
696
697   smpi_bench_end();
698   if (comm == MPI_COMM_NULL) {
699     retval = MPI_ERR_COMM;
700   } else if (newcomm == NULL) {
701     retval = MPI_ERR_ARG;
702   } else {
703     *newcomm = smpi_comm_new(smpi_comm_group(comm));
704     retval = MPI_SUCCESS;
705   }
706   smpi_bench_begin();
707   return retval;
708 }
709
710 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
711 {
712   int retval;
713
714   smpi_bench_end();
715   if (comm == MPI_COMM_NULL) {
716     retval = MPI_ERR_COMM;
717   } else if (group == MPI_GROUP_NULL) {
718     retval = MPI_ERR_GROUP;
719   } else if (newcomm == NULL) {
720     retval = MPI_ERR_ARG;
721   } else {
722     *newcomm = smpi_comm_new(group);
723     retval = MPI_SUCCESS;
724   }
725   smpi_bench_begin();
726   return retval;
727 }
728
729 int PMPI_Comm_free(MPI_Comm * comm)
730 {
731   int retval;
732
733   smpi_bench_end();
734   if (comm == NULL) {
735     retval = MPI_ERR_ARG;
736   } else if (*comm == MPI_COMM_NULL) {
737     retval = MPI_ERR_COMM;
738   } else {
739     smpi_comm_destroy(*comm);
740     *comm = MPI_COMM_NULL;
741     retval = MPI_SUCCESS;
742   }
743   smpi_bench_begin();
744   return retval;
745 }
746
747 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
748 {
749   int retval;
750
751   smpi_bench_end();
752   if (comm_out == NULL) {
753     retval = MPI_ERR_ARG;
754   } else if (comm == MPI_COMM_NULL) {
755     retval = MPI_ERR_COMM;
756   } else {
757     *comm_out = smpi_comm_split(comm, color, key);
758     retval = MPI_SUCCESS;
759   }
760   smpi_bench_begin();
761   return retval;
762 }
763
764 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
765                   int tag, MPI_Comm comm, MPI_Request * request)
766 {
767   int retval;
768
769   smpi_bench_end();
770   if (request == NULL) {
771     retval = MPI_ERR_ARG;
772   } else if (comm == MPI_COMM_NULL) {
773     retval = MPI_ERR_COMM;
774   } else {
775     *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
776     retval = MPI_SUCCESS;
777   }
778   smpi_bench_begin();
779   return retval;
780 }
781
782 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
783                   int tag, MPI_Comm comm, MPI_Request * request)
784 {
785   int retval;
786
787   smpi_bench_end();
788   if (request == NULL) {
789     retval = MPI_ERR_ARG;
790   } else if (comm == MPI_COMM_NULL) {
791     retval = MPI_ERR_COMM;
792   } else {
793     *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
794     retval = MPI_SUCCESS;
795   }
796   smpi_bench_begin();
797   return retval;
798 }
799
800 int PMPI_Start(MPI_Request * request)
801 {
802   int retval;
803
804   smpi_bench_end();
805   if (request == NULL) {
806     retval = MPI_ERR_ARG;
807   } else {
808     smpi_mpi_start(*request);
809     retval = MPI_SUCCESS;
810   }
811   smpi_bench_begin();
812   return retval;
813 }
814
815 int PMPI_Startall(int count, MPI_Request * requests)
816 {
817   int retval;
818
819   smpi_bench_end();
820   if (requests == NULL) {
821     retval = MPI_ERR_ARG;
822   } else {
823     smpi_mpi_startall(count, requests);
824     retval = MPI_SUCCESS;
825   }
826   smpi_bench_begin();
827   return retval;
828 }
829
830 int PMPI_Request_free(MPI_Request * request)
831 {
832   int retval;
833
834   smpi_bench_end();
835   if (request == NULL) {
836     retval = MPI_ERR_ARG;
837   } else {
838     smpi_mpi_request_free(request);
839     retval = MPI_SUCCESS;
840   }
841   smpi_bench_begin();
842   return retval;
843 }
844
845 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
846               int tag, MPI_Comm comm, MPI_Request * request)
847 {
848   int retval;
849
850   smpi_bench_end();
851 #ifdef HAVE_TRACING
852   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
853   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
854   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
855 #endif
856   if (request == NULL) {
857     retval = MPI_ERR_ARG;
858   } else if (comm == MPI_COMM_NULL) {
859     retval = MPI_ERR_COMM;
860   } else {
861     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
862     retval = MPI_SUCCESS;
863   }
864 #ifdef HAVE_TRACING
865   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
866   (*request)->recv = 1;
867 #endif
868   smpi_bench_begin();
869   return retval;
870 }
871
872 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
873               int tag, MPI_Comm comm, MPI_Request * request)
874 {
875   int retval;
876
877   smpi_bench_end();
878 #ifdef HAVE_TRACING
879   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
880   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
881   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
882   TRACE_smpi_send(rank, rank, dst_traced);
883 #endif
884   if (request == NULL) {
885     retval = MPI_ERR_ARG;
886   } else if (comm == MPI_COMM_NULL) {
887     retval = MPI_ERR_COMM;
888   } else {
889     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
890     retval = MPI_SUCCESS;
891   }
892 #ifdef HAVE_TRACING
893   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
894   (*request)->send = 1;
895 #endif
896   smpi_bench_begin();
897   return retval;
898 }
899
900 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
901              MPI_Comm comm, MPI_Status * status)
902 {
903   int retval;
904
905   smpi_bench_end();
906 #ifdef HAVE_TRACING
907   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
908   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
909   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
910 #endif
911   if (comm == MPI_COMM_NULL) {
912     retval = MPI_ERR_COMM;
913   } else {
914     smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
915     retval = MPI_SUCCESS;
916   }
917 #ifdef HAVE_TRACING
918   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
919   TRACE_smpi_recv(rank, src_traced, rank);
920 #endif
921   smpi_bench_begin();
922   return retval;
923 }
924
925 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
926              MPI_Comm comm)
927 {
928   int retval;
929
930   smpi_bench_end();
931 #ifdef HAVE_TRACING
932   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
933   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
934   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
935   TRACE_smpi_send(rank, rank, dst_traced);
936 #endif
937   if (comm == MPI_COMM_NULL) {
938     retval = MPI_ERR_COMM;
939   } else {
940     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
941     retval = MPI_SUCCESS;
942   }
943 #ifdef HAVE_TRACING
944   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
945 #endif
946   smpi_bench_begin();
947   return retval;
948 }
949
950 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
951                  int dst, int sendtag, void *recvbuf, int recvcount,
952                  MPI_Datatype recvtype, int src, int recvtag,
953                  MPI_Comm comm, MPI_Status * status)
954 {
955   int retval;
956
957   smpi_bench_end();
958 #ifdef HAVE_TRACING
959   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
960   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
961   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
962   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
963   TRACE_smpi_send(rank, rank, dst_traced);
964   TRACE_smpi_send(rank, src_traced, rank);
965 #endif
966   if (comm == MPI_COMM_NULL) {
967     retval = MPI_ERR_COMM;
968   } else if (sendtype == MPI_DATATYPE_NULL
969              || recvtype == MPI_DATATYPE_NULL) {
970     retval = MPI_ERR_TYPE;
971   } else {
972     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
973                       recvcount, recvtype, src, recvtag, comm, status);
974     retval = MPI_SUCCESS;
975   }
976 #ifdef HAVE_TRACING
977   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
978   TRACE_smpi_recv(rank, rank, dst_traced);
979   TRACE_smpi_recv(rank, src_traced, rank);
980 #endif
981   smpi_bench_begin();
982   return retval;
983 }
984
985 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
986                          int dst, int sendtag, int src, int recvtag,
987                          MPI_Comm comm, MPI_Status * status)
988 {
989   //TODO: suboptimal implementation
990   void *recvbuf;
991   int retval, size;
992
993   size = smpi_datatype_size(datatype) * count;
994   recvbuf = xbt_new(char, size);
995   retval =
996       MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
997                    datatype, src, recvtag, comm, status);
998   memcpy(buf, recvbuf, size * sizeof(char));
999   xbt_free(recvbuf);
1000   return retval;
1001 }
1002
1003 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1004 {
1005   int retval;
1006
1007   smpi_bench_end();
1008   if (request == NULL || flag == NULL) {
1009     retval = MPI_ERR_ARG;
1010   } else if (*request == MPI_REQUEST_NULL) {
1011     retval = MPI_ERR_REQUEST;
1012   } else {
1013     *flag = smpi_mpi_test(request, status);
1014     retval = MPI_SUCCESS;
1015   }
1016   smpi_bench_begin();
1017   return retval;
1018 }
1019
1020 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1021                 MPI_Status * status)
1022 {
1023   int retval;
1024
1025   smpi_bench_end();
1026   if (index == NULL || flag == NULL) {
1027     retval = MPI_ERR_ARG;
1028   } else {
1029     *flag = smpi_mpi_testany(count, requests, index, status);
1030     retval = MPI_SUCCESS;
1031   }
1032   smpi_bench_begin();
1033   return retval;
1034 }
1035
1036 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1037 {
1038   int retval;
1039
1040   smpi_bench_end();
1041 #ifdef HAVE_TRACING
1042   int rank = request && (*request)->comm != MPI_COMM_NULL
1043       ? smpi_comm_rank((*request)->comm)
1044       : -1;
1045   MPI_Group group = smpi_comm_group((*request)->comm);
1046   int src_traced = smpi_group_rank(group, (*request)->src);
1047   int dst_traced = smpi_group_rank(group, (*request)->dst);
1048   int is_wait_for_receive = (*request)->recv;
1049   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1050 #endif
1051   if (request == NULL) {
1052     retval = MPI_ERR_ARG;
1053   } else if (*request == MPI_REQUEST_NULL) {
1054     retval = MPI_ERR_REQUEST;
1055   } else {
1056     smpi_mpi_wait(request, status);
1057     retval = MPI_SUCCESS;
1058   }
1059 #ifdef HAVE_TRACING
1060   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1061   if (is_wait_for_receive) {
1062     TRACE_smpi_recv(rank, src_traced, dst_traced);
1063   }
1064 #endif
1065   smpi_bench_begin();
1066   return retval;
1067 }
1068
1069 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1070 {
1071   int retval;
1072
1073   smpi_bench_end();
1074 #ifdef HAVE_TRACING
1075   //save requests information for tracing
1076   int i;
1077   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
1078   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
1079   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
1080   for (i = 0; i < count; i++) {
1081     MPI_Request req = requests[i];      //already received requests are no longer valid
1082     if (req) {
1083       int *asrc = xbt_new(int, 1);
1084       int *adst = xbt_new(int, 1);
1085       int *arecv = xbt_new(int, 1);
1086       *asrc = req->src;
1087       *adst = req->dst;
1088       *arecv = req->recv;
1089       xbt_dynar_insert_at(srcs, i, asrc);
1090       xbt_dynar_insert_at(dsts, i, adst);
1091       xbt_dynar_insert_at(recvs, i, arecv);
1092     } else {
1093       int *t = xbt_new(int, 1);
1094       xbt_dynar_insert_at(srcs, i, t);
1095       xbt_dynar_insert_at(dsts, i, t);
1096       xbt_dynar_insert_at(recvs, i, t);
1097     }
1098   }
1099   int rank_traced = smpi_comm_rank(MPI_COMM_WORLD);
1100   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1101 #endif
1102   if (index == NULL) {
1103     retval = MPI_ERR_ARG;
1104   } else {
1105     *index = smpi_mpi_waitany(count, requests, status);
1106     retval = MPI_SUCCESS;
1107   }
1108 #ifdef HAVE_TRACING
1109   int src_traced, dst_traced, is_wait_for_receive;
1110   xbt_dynar_get_cpy(srcs, *index, &src_traced);
1111   xbt_dynar_get_cpy(dsts, *index, &dst_traced);
1112   xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
1113   if (is_wait_for_receive) {
1114     TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1115   }
1116   TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1117   //clean-up of dynars
1118   xbt_free(srcs);
1119   xbt_free(dsts);
1120   xbt_free(recvs);
1121 #endif
1122   smpi_bench_begin();
1123   return retval;
1124 }
1125
1126 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1127 {
1128
1129   smpi_bench_end();
1130 #ifdef HAVE_TRACING
1131   //save information from requests
1132   int i;
1133   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
1134   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
1135   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
1136   for (i = 0; i < count; i++) {
1137     MPI_Request req = requests[i];      //all req should be valid in Waitall
1138     int *asrc = xbt_new(int, 1);
1139     int *adst = xbt_new(int, 1);
1140     int *arecv = xbt_new(int, 1);
1141     *asrc = req->src;
1142     *adst = req->dst;
1143     *arecv = req->recv;
1144     xbt_dynar_insert_at(srcs, i, asrc);
1145     xbt_dynar_insert_at(dsts, i, adst);
1146     xbt_dynar_insert_at(recvs, i, arecv);
1147   }
1148   int rank_traced = smpi_comm_rank (MPI_COMM_WORLD);
1149   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1150 #endif
1151   smpi_mpi_waitall(count, requests, status);
1152 #ifdef HAVE_TRACING
1153   for (i = 0; i < count; i++) {
1154     int src_traced, dst_traced, is_wait_for_receive;
1155     xbt_dynar_get_cpy(srcs, i, &src_traced);
1156     xbt_dynar_get_cpy(dsts, i, &dst_traced);
1157     xbt_dynar_get_cpy(recvs, i, &is_wait_for_receive);
1158     if (is_wait_for_receive) {
1159       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1160     }
1161   }
1162   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1163   //clean-up of dynars
1164   xbt_free(srcs);
1165   xbt_free(dsts);
1166   xbt_free(recvs);
1167 #endif
1168   smpi_bench_begin();
1169   return MPI_SUCCESS;
1170 }
1171
1172 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1173                  int *indices, MPI_Status status[])
1174 {
1175   int retval;
1176
1177   smpi_bench_end();
1178   if (outcount == NULL || indices == NULL) {
1179     retval = MPI_ERR_ARG;
1180   } else {
1181     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1182     retval = MPI_SUCCESS;
1183   }
1184   smpi_bench_begin();
1185   return retval;
1186 }
1187
1188 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1189 {
1190   int retval;
1191
1192   smpi_bench_end();
1193 #ifdef HAVE_TRACING
1194   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1195   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1196   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1197 #endif
1198   if (comm == MPI_COMM_NULL) {
1199     retval = MPI_ERR_COMM;
1200   } else {
1201     smpi_mpi_bcast(buf, count, datatype, root, comm);
1202     retval = MPI_SUCCESS;
1203   }
1204 #ifdef HAVE_TRACING
1205   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1206 #endif
1207   smpi_bench_begin();
1208   return retval;
1209 }
1210
1211 int PMPI_Barrier(MPI_Comm comm)
1212 {
1213   int retval;
1214
1215   smpi_bench_end();
1216 #ifdef HAVE_TRACING
1217   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1218   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1219 #endif
1220   if (comm == MPI_COMM_NULL) {
1221     retval = MPI_ERR_COMM;
1222   } else {
1223     smpi_mpi_barrier(comm);
1224     retval = MPI_SUCCESS;
1225   }
1226 #ifdef HAVE_TRACING
1227   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1228 #endif
1229   smpi_bench_begin();
1230   return retval;
1231 }
1232
1233 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1234                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1235                int root, MPI_Comm comm)
1236 {
1237   int retval;
1238
1239   smpi_bench_end();
1240 #ifdef HAVE_TRACING
1241   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1242   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1243   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1244 #endif
1245   if (comm == MPI_COMM_NULL) {
1246     retval = MPI_ERR_COMM;
1247   } else if (sendtype == MPI_DATATYPE_NULL
1248              || recvtype == MPI_DATATYPE_NULL) {
1249     retval = MPI_ERR_TYPE;
1250   } else {
1251     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1252                     recvtype, root, comm);
1253     retval = MPI_SUCCESS;
1254   }
1255 #ifdef HAVE_TRACING
1256   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1257 #endif
1258   smpi_bench_begin();
1259   return retval;
1260 }
1261
1262 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1263                 void *recvbuf, int *recvcounts, int *displs,
1264                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1265 {
1266   int retval;
1267
1268   smpi_bench_end();
1269 #ifdef HAVE_TRACING
1270   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1271   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1272   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1273 #endif
1274   if (comm == MPI_COMM_NULL) {
1275     retval = MPI_ERR_COMM;
1276   } else if (sendtype == MPI_DATATYPE_NULL
1277              || recvtype == MPI_DATATYPE_NULL) {
1278     retval = MPI_ERR_TYPE;
1279   } else if (recvcounts == NULL || displs == NULL) {
1280     retval = MPI_ERR_ARG;
1281   } else {
1282     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1283                      displs, recvtype, root, comm);
1284     retval = MPI_SUCCESS;
1285   }
1286 #ifdef HAVE_TRACING
1287   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1288 #endif
1289   smpi_bench_begin();
1290   return retval;
1291 }
1292
1293 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1294                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1295                   MPI_Comm comm)
1296 {
1297   int retval;
1298
1299   smpi_bench_end();
1300 #ifdef HAVE_TRACING
1301   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1302   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1303 #endif
1304   if (comm == MPI_COMM_NULL) {
1305     retval = MPI_ERR_COMM;
1306   } else if (sendtype == MPI_DATATYPE_NULL
1307              || recvtype == MPI_DATATYPE_NULL) {
1308     retval = MPI_ERR_TYPE;
1309   } else {
1310     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1311                        recvtype, comm);
1312     retval = MPI_SUCCESS;
1313   }
1314 #ifdef HAVE_TRACING
1315   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1316 #endif
1317   smpi_bench_begin();
1318   return retval;
1319 }
1320
1321 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1322                    void *recvbuf, int *recvcounts, int *displs,
1323                    MPI_Datatype recvtype, MPI_Comm comm)
1324 {
1325   int retval;
1326
1327   smpi_bench_end();
1328 #ifdef HAVE_TRACING
1329   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1330   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1331 #endif
1332   if (comm == MPI_COMM_NULL) {
1333     retval = MPI_ERR_COMM;
1334   } else if (sendtype == MPI_DATATYPE_NULL
1335              || recvtype == MPI_DATATYPE_NULL) {
1336     retval = MPI_ERR_TYPE;
1337   } else if (recvcounts == NULL || displs == NULL) {
1338     retval = MPI_ERR_ARG;
1339   } else {
1340     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1341                         displs, recvtype, comm);
1342     retval = MPI_SUCCESS;
1343   }
1344 #ifdef HAVE_TRACING
1345   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1346 #endif
1347   smpi_bench_begin();
1348   return retval;
1349 }
1350
1351 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1352                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1353                 int root, MPI_Comm comm)
1354 {
1355   int retval;
1356
1357   smpi_bench_end();
1358 #ifdef HAVE_TRACING
1359   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1360   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1361   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1362 #endif
1363   if (comm == MPI_COMM_NULL) {
1364     retval = MPI_ERR_COMM;
1365   } else if (sendtype == MPI_DATATYPE_NULL
1366              || recvtype == MPI_DATATYPE_NULL) {
1367     retval = MPI_ERR_TYPE;
1368   } else {
1369     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1370                      recvtype, root, comm);
1371     retval = MPI_SUCCESS;
1372   }
1373 #ifdef HAVE_TRACING
1374   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1375 #endif
1376   smpi_bench_begin();
1377   return retval;
1378 }
1379
1380 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1381                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1382                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1383 {
1384   int retval;
1385
1386   smpi_bench_end();
1387 #ifdef HAVE_TRACING
1388   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1389   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1390   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1391 #endif
1392   if (comm == MPI_COMM_NULL) {
1393     retval = MPI_ERR_COMM;
1394   } else if (sendtype == MPI_DATATYPE_NULL
1395              || recvtype == MPI_DATATYPE_NULL) {
1396     retval = MPI_ERR_TYPE;
1397   } else if (sendcounts == NULL || displs == NULL) {
1398     retval = MPI_ERR_ARG;
1399   } else {
1400     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1401                       recvcount, recvtype, root, comm);
1402     retval = MPI_SUCCESS;
1403   }
1404 #ifdef HAVE_TRACING
1405   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1406 #endif
1407   smpi_bench_begin();
1408   return retval;
1409 }
1410
1411 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1412                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1413 {
1414   int retval;
1415
1416   smpi_bench_end();
1417 #ifdef HAVE_TRACING
1418   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1419   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1420   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1421 #endif
1422   if (comm == MPI_COMM_NULL) {
1423     retval = MPI_ERR_COMM;
1424   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1425     retval = MPI_ERR_ARG;
1426   } else {
1427     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1428     retval = MPI_SUCCESS;
1429   }
1430 #ifdef HAVE_TRACING
1431   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1432 #endif
1433   smpi_bench_begin();
1434   return retval;
1435 }
1436
1437 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1438                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1439 {
1440   int retval;
1441
1442   smpi_bench_end();
1443 #ifdef HAVE_TRACING
1444   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1445   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1446 #endif
1447   if (comm == MPI_COMM_NULL) {
1448     retval = MPI_ERR_COMM;
1449   } else if (datatype == MPI_DATATYPE_NULL) {
1450     retval = MPI_ERR_TYPE;
1451   } else if (op == MPI_OP_NULL) {
1452     retval = MPI_ERR_OP;
1453   } else {
1454     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
1455     retval = MPI_SUCCESS;
1456   }
1457 #ifdef HAVE_TRACING
1458   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1459 #endif
1460   smpi_bench_begin();
1461   return retval;
1462 }
1463
1464 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
1465              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1466 {
1467   int retval;
1468
1469   smpi_bench_end();
1470 #ifdef HAVE_TRACING
1471   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1472   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1473 #endif
1474   if (comm == MPI_COMM_NULL) {
1475     retval = MPI_ERR_COMM;
1476   } else if (datatype == MPI_DATATYPE_NULL) {
1477     retval = MPI_ERR_TYPE;
1478   } else if (op == MPI_OP_NULL) {
1479     retval = MPI_ERR_OP;
1480   } else {
1481     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1482     retval = MPI_SUCCESS;
1483   }
1484 #ifdef HAVE_TRACING
1485   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1486 #endif
1487   smpi_bench_begin();
1488   return retval;
1489 }
1490
1491 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1492                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1493 {
1494   int retval, i, size, count;
1495   int *displs;
1496   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1497
1498   smpi_bench_end();
1499 #ifdef HAVE_TRACING
1500   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1501 #endif
1502   if (comm == MPI_COMM_NULL) {
1503     retval = MPI_ERR_COMM;
1504   } else if (datatype == MPI_DATATYPE_NULL) {
1505     retval = MPI_ERR_TYPE;
1506   } else if (op == MPI_OP_NULL) {
1507     retval = MPI_ERR_OP;
1508   } else if (recvcounts == NULL) {
1509     retval = MPI_ERR_ARG;
1510   } else {
1511     /* arbitrarily choose root as rank 0 */
1512     /* TODO: faster direct implementation ? */
1513     size = smpi_comm_size(comm);
1514     count = 0;
1515     displs = xbt_new(int, size);
1516     for (i = 0; i < size; i++) {
1517       count += recvcounts[i];
1518       displs[i] = 0;
1519     }
1520     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1521     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf,
1522                       recvcounts[rank], datatype, 0, comm);
1523     xbt_free(displs);
1524     retval = MPI_SUCCESS;
1525   }
1526 #ifdef HAVE_TRACING
1527   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1528 #endif
1529   smpi_bench_begin();
1530   return retval;
1531 }
1532
1533 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1534                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
1535                  MPI_Comm comm)
1536 {
1537   int retval, size, sendsize;
1538
1539   smpi_bench_end();
1540 #ifdef HAVE_TRACING
1541   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1542   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1543 #endif
1544   if (comm == MPI_COMM_NULL) {
1545     retval = MPI_ERR_COMM;
1546   } else if (sendtype == MPI_DATATYPE_NULL
1547              || recvtype == MPI_DATATYPE_NULL) {
1548     retval = MPI_ERR_TYPE;
1549   } else {
1550     size = smpi_comm_size(comm);
1551     sendsize = smpi_datatype_size(sendtype) * sendcount;
1552     if (sendsize < 200 && size > 12) {
1553       retval =
1554           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
1555                                          recvbuf, recvcount, recvtype,
1556                                          comm);
1557     } else if (sendsize < 3000) {
1558       retval =
1559           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
1560                                                 sendtype, recvbuf,
1561                                                 recvcount, recvtype, comm);
1562     } else {
1563       retval =
1564           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
1565                                             recvbuf, recvcount, recvtype,
1566                                             comm);
1567     }
1568   }
1569 #ifdef HAVE_TRACING
1570   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1571 #endif
1572   smpi_bench_begin();
1573   return retval;
1574 }
1575
1576 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
1577                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
1578                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
1579 {
1580   int retval;
1581
1582   smpi_bench_end();
1583 #ifdef HAVE_TRACING
1584   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
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 if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
1593              || recvdisps == NULL) {
1594     retval = MPI_ERR_ARG;
1595   } else {
1596     retval =
1597         smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
1598                                   recvbuf, recvcounts, recvdisps, recvtype,
1599                                   comm);
1600   }
1601 #ifdef HAVE_TRACING
1602   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1603 #endif
1604   smpi_bench_begin();
1605   return retval;
1606 }
1607
1608
1609 int PMPI_Get_processor_name(char *name, int *resultlen)
1610 {
1611   int retval = MPI_SUCCESS;
1612
1613   smpi_bench_end();
1614   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
1615           MPI_MAX_PROCESSOR_NAME - 1);
1616   *resultlen =
1617       strlen(name) >
1618       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1619
1620   smpi_bench_begin();
1621   return retval;
1622 }
1623
1624 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
1625 {
1626   int retval = MPI_SUCCESS;
1627   size_t size;
1628
1629   smpi_bench_end();
1630   if (status == NULL || count == NULL) {
1631     retval = MPI_ERR_ARG;
1632   } else if (datatype == MPI_DATATYPE_NULL) {
1633     retval = MPI_ERR_TYPE;
1634   } else {
1635     size = smpi_datatype_size(datatype);
1636     if (size == 0) {
1637       *count = 0;
1638     } else if (status->count % size != 0) {
1639       retval = MPI_UNDEFINED;
1640     } else {
1641       *count = smpi_mpi_get_count(status, datatype);
1642     }
1643   }
1644   smpi_bench_begin();
1645   return retval;
1646 }
1647
1648 /* The following calls are not yet implemented and will fail at runtime. */
1649 /* Once implemented, please move them above this notice. */
1650
1651 static int not_yet_implemented(void) {
1652    xbt_die("Not yet implemented");
1653    return MPI_ERR_OTHER;
1654 }
1655
1656 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
1657    return not_yet_implemented();
1658 }
1659
1660 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
1661    return not_yet_implemented();
1662 }
1663
1664 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
1665    return not_yet_implemented();
1666 }
1667
1668 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
1669    return not_yet_implemented();
1670 }
1671
1672 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
1673    return not_yet_implemented();
1674 }
1675
1676 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
1677    return not_yet_implemented();
1678 }
1679
1680 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
1681    return not_yet_implemented();
1682 }
1683
1684 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
1685    return not_yet_implemented();
1686 }
1687
1688 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
1689    return not_yet_implemented();
1690 }
1691
1692 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
1693    return not_yet_implemented();
1694 }
1695
1696 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
1697    return not_yet_implemented();
1698 }
1699
1700 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
1701    return not_yet_implemented();
1702 }
1703
1704 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
1705    return not_yet_implemented();
1706 }
1707
1708 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
1709    return not_yet_implemented();
1710 }
1711
1712 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
1713    return not_yet_implemented();
1714 }
1715
1716 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
1717    return not_yet_implemented();
1718 }
1719
1720 int PMPI_Error_class(int errorcode, int* errorclass) {
1721    return not_yet_implemented();
1722 }
1723
1724 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
1725    return not_yet_implemented();
1726 }
1727
1728 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
1729    return not_yet_implemented();
1730 }
1731
1732 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
1733    return not_yet_implemented();
1734 }
1735
1736 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
1737    return not_yet_implemented();
1738 }
1739
1740 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
1741    return not_yet_implemented();
1742 }
1743
1744 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* newtype) {
1745    return not_yet_implemented();
1746 }
1747
1748 int PMPI_Cancel(MPI_Request* request) {
1749    return not_yet_implemented();
1750 }
1751
1752 int PMPI_Buffer_attach(void* buffer, int size) {
1753    return not_yet_implemented();
1754 }
1755
1756 int PMPI_Buffer_detach(void* buffer, int* size) {
1757    return not_yet_implemented();
1758 }
1759
1760 int PMPI_Testsome(int incount, MPI_Request* requests, int* outcount, int* indices, MPI_Status* statuses) {
1761    return not_yet_implemented();
1762 }
1763
1764 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
1765    return not_yet_implemented();
1766 }
1767
1768 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
1769    return not_yet_implemented();
1770 }
1771
1772 int PMPI_Type_commit(MPI_Datatype* datatype) {
1773    return not_yet_implemented();
1774 }
1775
1776 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
1777    return not_yet_implemented();
1778 }
1779
1780 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
1781    return not_yet_implemented();
1782 }
1783
1784 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
1785    return not_yet_implemented();
1786 }
1787
1788 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* newtype) {
1789    return not_yet_implemented();
1790 }
1791
1792 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
1793    return not_yet_implemented();
1794 }
1795
1796 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1797    return not_yet_implemented();
1798 }
1799
1800 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1801    return not_yet_implemented();
1802 }
1803
1804 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
1805    return not_yet_implemented();
1806 }
1807
1808 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
1809    return not_yet_implemented();
1810 }
1811
1812 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1813    return not_yet_implemented();
1814 }
1815
1816 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1817    return not_yet_implemented();
1818 }
1819
1820 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1821    return not_yet_implemented();
1822 }
1823
1824 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
1825    return not_yet_implemented();
1826 }
1827
1828 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
1829    return not_yet_implemented();
1830 }
1831
1832 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1833    return not_yet_implemented();
1834 }
1835
1836 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1837    return not_yet_implemented();
1838 }
1839
1840 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
1841    return not_yet_implemented();
1842 }
1843
1844 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
1845    return not_yet_implemented();
1846 }
1847
1848 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
1849    return not_yet_implemented();
1850 }
1851
1852 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1853    return not_yet_implemented();
1854 }
1855
1856 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1857    return not_yet_implemented();
1858 }
1859
1860 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1861    return not_yet_implemented();
1862 }
1863
1864 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
1865    return not_yet_implemented();
1866 }
1867
1868 int PMPI_Keyval_free(int* keyval) {
1869    return not_yet_implemented();
1870 }
1871
1872 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
1873    return not_yet_implemented();
1874 }
1875
1876 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
1877    return not_yet_implemented();
1878 }
1879
1880 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses) {
1881    return not_yet_implemented();
1882 }
1883
1884 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
1885    return not_yet_implemented();
1886 }
1887
1888 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
1889    return not_yet_implemented();
1890 }
1891
1892 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1893    return not_yet_implemented();
1894 }
1895
1896 int PMPI_Initialized(int* flag) {
1897    return not_yet_implemented();
1898 }