Logo AND Algorithmique Numérique Distribuée

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