Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Make bprintf abort on error, and define bvprintf accordingly.
[simgrid.git] / src / smpi / smpi_mpi.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_coll_private.h"
9 #include "smpi_mpi_dt_private.h"
10
11 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_mpi, smpi,
12                                 "Logging specific to SMPI (mpi)");
13
14 /* MPI User level calls */
15
16 int MPI_Init(int *argc, char ***argv)
17 {
18   smpi_process_init(argc, argv);
19 #ifdef HAVE_TRACING
20   TRACE_smpi_init(smpi_process_index());
21 #endif
22   smpi_bench_begin(-1, NULL);
23   return MPI_SUCCESS;
24 }
25
26 int MPI_Finalize(void)
27 {
28   smpi_bench_end(-1, NULL);
29 #ifdef HAVE_TRACING
30   TRACE_smpi_finalize(smpi_process_index());
31 #endif
32   smpi_process_destroy();
33   return MPI_SUCCESS;
34 }
35
36 int MPI_Init_thread(int *argc, char ***argv, int required, int *provided)
37 {
38   if (provided != NULL) {
39     *provided = MPI_THREAD_MULTIPLE;
40   }
41   return MPI_Init(argc, argv);
42 }
43
44 int MPI_Query_thread(int *provided)
45 {
46   int retval;
47
48   smpi_bench_end(-1, NULL);
49   if (provided == NULL) {
50     retval = MPI_ERR_ARG;
51   } else {
52     *provided = MPI_THREAD_MULTIPLE;
53     retval = MPI_SUCCESS;
54   }
55   smpi_bench_begin(-1, NULL);
56   return retval;
57 }
58
59 int MPI_Is_thread_main(int *flag)
60 {
61   int retval;
62
63   smpi_bench_end(-1, NULL);
64   if (flag == NULL) {
65     retval = MPI_ERR_ARG;
66   } else {
67     *flag = smpi_process_index() == 0;
68     retval = MPI_SUCCESS;
69   }
70   smpi_bench_begin(-1, NULL);
71   return retval;
72 }
73
74 int MPI_Abort(MPI_Comm comm, int errorcode)
75 {
76   smpi_bench_end(-1, NULL);
77   smpi_process_destroy();
78   // FIXME: should kill all processes in comm instead
79   SIMIX_process_kill(SIMIX_process_self());
80   return MPI_SUCCESS;
81 }
82
83 double MPI_Wtime(void)
84 {
85   double time;
86
87   smpi_bench_end(-1, NULL);
88   time = SIMIX_get_clock();
89   smpi_bench_begin(-1, NULL);
90   return time;
91 }
92
93 int MPI_Address(void *location, MPI_Aint * address)
94 {
95   int retval;
96
97   smpi_bench_end(-1, NULL);
98   if (!address) {
99     retval = MPI_ERR_ARG;
100   } else {
101     *address = (MPI_Aint) location;
102   }
103   smpi_bench_begin(-1, NULL);
104   return retval;
105 }
106
107 int MPI_Type_free(MPI_Datatype * datatype)
108 {
109   int retval;
110
111   smpi_bench_end(-1, NULL);
112   if (!datatype) {
113     retval = MPI_ERR_ARG;
114   } else {
115     // FIXME: always fail for now
116     retval = MPI_ERR_TYPE;
117   }
118   smpi_bench_begin(-1, NULL);
119   return retval;
120 }
121
122 int MPI_Type_size(MPI_Datatype datatype, int *size)
123 {
124   int retval;
125
126   smpi_bench_end(-1, NULL);
127   if (datatype == MPI_DATATYPE_NULL) {
128     retval = MPI_ERR_TYPE;
129   } else if (size == NULL) {
130     retval = MPI_ERR_ARG;
131   } else {
132     *size = (int) smpi_datatype_size(datatype);
133     retval = MPI_SUCCESS;
134   }
135   smpi_bench_begin(-1, NULL);
136   return retval;
137 }
138
139 int MPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb,
140                         MPI_Aint * extent)
141 {
142   int retval;
143
144   smpi_bench_end(-1, NULL);
145   if (datatype == MPI_DATATYPE_NULL) {
146     retval = MPI_ERR_TYPE;
147   } else if (lb == NULL || extent == NULL) {
148     retval = MPI_ERR_ARG;
149   } else {
150     retval = smpi_datatype_extent(datatype, lb, extent);
151   }
152   smpi_bench_begin(-1, NULL);
153   return retval;
154 }
155
156 int MPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
157 {
158   int retval;
159   MPI_Aint dummy;
160
161   smpi_bench_end(-1, NULL);
162   if (datatype == MPI_DATATYPE_NULL) {
163     retval = MPI_ERR_TYPE;
164   } else if (extent == NULL) {
165     retval = MPI_ERR_ARG;
166   } else {
167     retval = smpi_datatype_extent(datatype, &dummy, extent);
168   }
169   smpi_bench_begin(-1, NULL);
170   return retval;
171 }
172
173 int MPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
174 {
175   int retval;
176
177   smpi_bench_end(-1, NULL);
178   if (datatype == MPI_DATATYPE_NULL) {
179     retval = MPI_ERR_TYPE;
180   } else if (disp == NULL) {
181     retval = MPI_ERR_ARG;
182   } else {
183     *disp = smpi_datatype_lb(datatype);
184     retval = MPI_SUCCESS;
185   }
186   smpi_bench_begin(-1, NULL);
187   return retval;
188 }
189
190 int MPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
191 {
192   int retval;
193
194   smpi_bench_end(-1, NULL);
195   if (datatype == MPI_DATATYPE_NULL) {
196     retval = MPI_ERR_TYPE;
197   } else if (disp == NULL) {
198     retval = MPI_ERR_ARG;
199   } else {
200     *disp = smpi_datatype_ub(datatype);
201     retval = MPI_SUCCESS;
202   }
203   smpi_bench_begin(-1, NULL);
204   return retval;
205 }
206
207 int MPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
208 {
209   int retval;
210
211   smpi_bench_end(-1, NULL);
212   if (function == NULL || op == NULL) {
213     retval = MPI_ERR_ARG;
214   } else {
215     *op = smpi_op_new(function, commute);
216     retval = MPI_SUCCESS;
217   }
218   smpi_bench_begin(-1, NULL);
219   return retval;
220 }
221
222 int MPI_Op_free(MPI_Op * op)
223 {
224   int retval;
225
226   smpi_bench_end(-1, NULL);
227   if (op == NULL) {
228     retval = MPI_ERR_ARG;
229   } else if (*op == MPI_OP_NULL) {
230     retval = MPI_ERR_OP;
231   } else {
232     smpi_op_destroy(*op);
233     *op = MPI_OP_NULL;
234     retval = MPI_SUCCESS;
235   }
236   smpi_bench_begin(-1, NULL);
237   return retval;
238 }
239
240 int MPI_Group_free(MPI_Group * group)
241 {
242   int retval;
243
244   smpi_bench_end(-1, NULL);
245   if (group == NULL) {
246     retval = MPI_ERR_ARG;
247   } else {
248     smpi_group_destroy(*group);
249     *group = MPI_GROUP_NULL;
250     retval = MPI_SUCCESS;
251   }
252   smpi_bench_begin(-1, NULL);
253   return retval;
254 }
255
256 int MPI_Group_size(MPI_Group group, int *size)
257 {
258   int retval;
259
260   smpi_bench_end(-1, NULL);
261   if (group == MPI_GROUP_NULL) {
262     retval = MPI_ERR_GROUP;
263   } else if (size == NULL) {
264     retval = MPI_ERR_ARG;
265   } else {
266     *size = smpi_group_size(group);
267     retval = MPI_SUCCESS;
268   }
269   smpi_bench_begin(-1, NULL);
270   return retval;
271 }
272
273 int MPI_Group_rank(MPI_Group group, int *rank)
274 {
275   int retval;
276
277   smpi_bench_end(-1, NULL);
278   if (group == MPI_GROUP_NULL) {
279     retval = MPI_ERR_GROUP;
280   } else if (rank == NULL) {
281     retval = MPI_ERR_ARG;
282   } else {
283     *rank = smpi_group_rank(group, smpi_process_index());
284     retval = MPI_SUCCESS;
285   }
286   smpi_bench_begin(-1, NULL);
287   return retval;
288 }
289
290 int MPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1,
291                               MPI_Group group2, int *ranks2)
292 {
293   int retval, i, index;
294
295   smpi_bench_end(-1, NULL);
296   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
297     retval = MPI_ERR_GROUP;
298   } else {
299     for (i = 0; i < n; i++) {
300       index = smpi_group_index(group1, ranks1[i]);
301       ranks2[i] = smpi_group_rank(group2, index);
302     }
303     retval = MPI_SUCCESS;
304   }
305   smpi_bench_begin(-1, NULL);
306   return retval;
307 }
308
309 int MPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
310 {
311   int retval;
312
313   smpi_bench_end(-1, NULL);
314   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
315     retval = MPI_ERR_GROUP;
316   } else if (result == NULL) {
317     retval = MPI_ERR_ARG;
318   } else {
319     *result = smpi_group_compare(group1, group2);
320     retval = MPI_SUCCESS;
321   }
322   smpi_bench_begin(-1, NULL);
323   return retval;
324 }
325
326 int MPI_Group_union(MPI_Group group1, MPI_Group group2,
327                     MPI_Group * newgroup)
328 {
329   int retval, i, proc1, proc2, size, size2;
330
331   smpi_bench_end(-1, NULL);
332   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
333     retval = MPI_ERR_GROUP;
334   } else if (newgroup == NULL) {
335     retval = MPI_ERR_ARG;
336   } else {
337     size = smpi_group_size(group1);
338     size2 = smpi_group_size(group2);
339     for (i = 0; i < size2; i++) {
340       proc2 = smpi_group_index(group2, i);
341       proc1 = smpi_group_rank(group1, proc2);
342       if (proc1 == MPI_UNDEFINED) {
343         size++;
344       }
345     }
346     if (size == 0) {
347       *newgroup = MPI_GROUP_EMPTY;
348     } else {
349       *newgroup = smpi_group_new(size);
350       size2 = smpi_group_size(group1);
351       for (i = 0; i < size2; i++) {
352         proc1 = smpi_group_index(group1, i);
353         smpi_group_set_mapping(*newgroup, proc1, i);
354       }
355       for (i = size2; i < size; i++) {
356         proc2 = smpi_group_index(group2, i - size2);
357         smpi_group_set_mapping(*newgroup, proc2, i);
358       }
359     }
360     smpi_group_use(*newgroup);
361     retval = MPI_SUCCESS;
362   }
363   smpi_bench_begin(-1, NULL);
364   return retval;
365 }
366
367 int MPI_Group_intersection(MPI_Group group1, MPI_Group group2,
368                            MPI_Group * newgroup)
369 {
370   int retval, i, proc1, proc2, size, size2;
371
372   smpi_bench_end(-1, NULL);
373   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
374     retval = MPI_ERR_GROUP;
375   } else if (newgroup == NULL) {
376     retval = MPI_ERR_ARG;
377   } else {
378     size = smpi_group_size(group1);
379     size2 = smpi_group_size(group2);
380     for (i = 0; i < size2; i++) {
381       proc2 = smpi_group_index(group2, i);
382       proc1 = smpi_group_rank(group1, proc2);
383       if (proc1 == MPI_UNDEFINED) {
384         size--;
385       }
386     }
387     if (size == 0) {
388       *newgroup = MPI_GROUP_EMPTY;
389     } else {
390       *newgroup = smpi_group_new(size);
391       size2 = smpi_group_size(group1);
392       for (i = 0; i < size2; i++) {
393         proc1 = smpi_group_index(group1, i);
394         proc2 = smpi_group_rank(group2, proc1);
395         if (proc2 != MPI_UNDEFINED) {
396           smpi_group_set_mapping(*newgroup, proc1, i);
397         }
398       }
399     }
400     smpi_group_use(*newgroup);
401     retval = MPI_SUCCESS;
402   }
403   smpi_bench_begin(-1, NULL);
404   return retval;
405 }
406
407 int MPI_Group_difference(MPI_Group group1, MPI_Group group2,
408                          MPI_Group * newgroup)
409 {
410   int retval, i, proc1, proc2, size, size2;
411
412   smpi_bench_end(-1, NULL);
413   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
414     retval = MPI_ERR_GROUP;
415   } else if (newgroup == NULL) {
416     retval = MPI_ERR_ARG;
417   } else {
418     size = size2 = smpi_group_size(group1);
419     for (i = 0; i < size2; i++) {
420       proc1 = smpi_group_index(group1, i);
421       proc2 = smpi_group_rank(group2, proc1);
422       if (proc2 != MPI_UNDEFINED) {
423         size--;
424       }
425     }
426     if (size == 0) {
427       *newgroup = MPI_GROUP_EMPTY;
428     } else {
429       *newgroup = smpi_group_new(size);
430       for (i = 0; i < size2; i++) {
431         proc1 = smpi_group_index(group1, i);
432         proc2 = smpi_group_rank(group2, proc1);
433         if (proc2 == MPI_UNDEFINED) {
434           smpi_group_set_mapping(*newgroup, proc1, i);
435         }
436       }
437     }
438     smpi_group_use(*newgroup);
439     retval = MPI_SUCCESS;
440   }
441   smpi_bench_begin(-1, NULL);
442   return retval;
443 }
444
445 int MPI_Group_incl(MPI_Group group, int n, int *ranks,
446                    MPI_Group * newgroup)
447 {
448   int retval, i, index;
449
450   smpi_bench_end(-1, NULL);
451   if (group == MPI_GROUP_NULL) {
452     retval = MPI_ERR_GROUP;
453   } else if (newgroup == NULL) {
454     retval = MPI_ERR_ARG;
455   } else {
456     if (n == 0) {
457       *newgroup = MPI_GROUP_EMPTY;
458     } else if (n == smpi_group_size(group)) {
459       *newgroup = group;
460     } else {
461       *newgroup = smpi_group_new(n);
462       for (i = 0; i < n; i++) {
463         index = smpi_group_index(group, ranks[i]);
464         smpi_group_set_mapping(*newgroup, index, i);
465       }
466     }
467     smpi_group_use(*newgroup);
468     retval = MPI_SUCCESS;
469   }
470   smpi_bench_begin(-1, NULL);
471   return retval;
472 }
473
474 int MPI_Group_excl(MPI_Group group, int n, int *ranks,
475                    MPI_Group * newgroup)
476 {
477   int retval, i, size, rank, index;
478
479   smpi_bench_end(-1, NULL);
480   if (group == MPI_GROUP_NULL) {
481     retval = MPI_ERR_GROUP;
482   } else if (newgroup == NULL) {
483     retval = MPI_ERR_ARG;
484   } else {
485     if (n == 0) {
486       *newgroup = group;
487     } else if (n == smpi_group_size(group)) {
488       *newgroup = MPI_GROUP_EMPTY;
489     } else {
490       size = smpi_group_size(group) - n;
491       *newgroup = smpi_group_new(size);
492       rank = 0;
493       while (rank < size) {
494         for (i = 0; i < n; i++) {
495           if (ranks[i] == rank) {
496             break;
497           }
498         }
499         if (i >= n) {
500           index = smpi_group_index(group, rank);
501           smpi_group_set_mapping(*newgroup, index, rank);
502           rank++;
503         }
504       }
505     }
506     smpi_group_use(*newgroup);
507     retval = MPI_SUCCESS;
508   }
509   smpi_bench_begin(-1, NULL);
510   return retval;
511 }
512
513 int MPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
514                          MPI_Group * newgroup)
515 {
516   int retval, i, j, rank, size, index;
517
518   smpi_bench_end(-1, NULL);
519   if (group == MPI_GROUP_NULL) {
520     retval = MPI_ERR_GROUP;
521   } else if (newgroup == NULL) {
522     retval = MPI_ERR_ARG;
523   } else {
524     if (n == 0) {
525       *newgroup = MPI_GROUP_EMPTY;
526     } else {
527       size = 0;
528       for (i = 0; i < n; i++) {
529         for (rank = ranges[i][0];       /* First */
530              rank >= 0 && rank <= ranges[i][1]; /* Last */
531              rank += ranges[i][2] /* Stride */ ) {
532           size++;
533         }
534       }
535       if (size == smpi_group_size(group)) {
536         *newgroup = group;
537       } else {
538         *newgroup = smpi_group_new(size);
539         j = 0;
540         for (i = 0; i < n; i++) {
541           for (rank = ranges[i][0];     /* First */
542                rank >= 0 && rank <= ranges[i][1];       /* Last */
543                rank += ranges[i][2] /* Stride */ ) {
544             index = smpi_group_index(group, rank);
545             smpi_group_set_mapping(*newgroup, index, j);
546             j++;
547           }
548         }
549       }
550     }
551     smpi_group_use(*newgroup);
552     retval = MPI_SUCCESS;
553   }
554   smpi_bench_begin(-1, NULL);
555   return retval;
556 }
557
558 int MPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
559                          MPI_Group * newgroup)
560 {
561   int retval, i, newrank, rank, size, index, add;
562
563   smpi_bench_end(-1, NULL);
564   if (group == MPI_GROUP_NULL) {
565     retval = MPI_ERR_GROUP;
566   } else if (newgroup == NULL) {
567     retval = MPI_ERR_ARG;
568   } else {
569     if (n == 0) {
570       *newgroup = group;
571     } else {
572       size = smpi_group_size(group);
573       for (i = 0; i < n; i++) {
574         for (rank = ranges[i][0];       /* First */
575              rank >= 0 && rank <= ranges[i][1]; /* Last */
576              rank += ranges[i][2] /* Stride */ ) {
577           size--;
578         }
579       }
580       if (size == 0) {
581         *newgroup = MPI_GROUP_EMPTY;
582       } else {
583         *newgroup = smpi_group_new(size);
584         newrank = 0;
585         while (newrank < size) {
586           for (i = 0; i < n; i++) {
587             add = 1;
588             for (rank = ranges[i][0];   /* First */
589                  rank >= 0 && rank <= ranges[i][1];     /* Last */
590                  rank += ranges[i][2] /* Stride */ ) {
591               if (rank == newrank) {
592                 add = 0;
593                 break;
594               }
595             }
596             if (add == 1) {
597               index = smpi_group_index(group, newrank);
598               smpi_group_set_mapping(*newgroup, index, newrank);
599             }
600           }
601         }
602       }
603     }
604     smpi_group_use(*newgroup);
605     retval = MPI_SUCCESS;
606   }
607   smpi_bench_begin(-1, NULL);
608   return retval;
609 }
610
611 int MPI_Comm_rank(MPI_Comm comm, int *rank)
612 {
613   int retval;
614
615   smpi_bench_end(-1, NULL);
616   if (comm == MPI_COMM_NULL) {
617     retval = MPI_ERR_COMM;
618   } else {
619     *rank = smpi_comm_rank(comm);
620     retval = MPI_SUCCESS;
621   }
622   smpi_bench_begin(-1, NULL);
623   return retval;
624 }
625
626 int MPI_Comm_size(MPI_Comm comm, int *size)
627 {
628   int retval;
629
630   smpi_bench_end(-1, NULL);
631   if (comm == MPI_COMM_NULL) {
632     retval = MPI_ERR_COMM;
633   } else if (size == NULL) {
634     retval = MPI_ERR_ARG;
635   } else {
636     *size = smpi_comm_size(comm);
637     retval = MPI_SUCCESS;
638   }
639   smpi_bench_begin(-1, NULL);
640   return retval;
641 }
642
643 int MPI_Comm_group(MPI_Comm comm, MPI_Group * group)
644 {
645   int retval;
646
647   smpi_bench_end(-1, NULL);
648   if (comm == MPI_COMM_NULL) {
649     retval = MPI_ERR_COMM;
650   } else if (group == NULL) {
651     retval = MPI_ERR_ARG;
652   } else {
653     *group = smpi_comm_group(comm);
654     retval = MPI_SUCCESS;
655   }
656   smpi_bench_begin(-1, NULL);
657   return retval;
658 }
659
660 int MPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
661 {
662   int retval;
663
664   smpi_bench_end(-1, NULL);
665   if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
666     retval = MPI_ERR_COMM;
667   } else if (result == NULL) {
668     retval = MPI_ERR_ARG;
669   } else {
670     if (comm1 == comm2) {       /* Same communicators means same groups */
671       *result = MPI_IDENT;
672     } else {
673       *result =
674           smpi_group_compare(smpi_comm_group(comm1),
675                              smpi_comm_group(comm2));
676       if (*result == MPI_IDENT) {
677         *result = MPI_CONGRUENT;
678       }
679     }
680     retval = MPI_SUCCESS;
681   }
682   smpi_bench_begin(-1, NULL);
683   return retval;
684 }
685
686 int MPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
687 {
688   int retval;
689
690   smpi_bench_end(-1, NULL);
691   if (comm == MPI_COMM_NULL) {
692     retval = MPI_ERR_COMM;
693   } else if (newcomm == NULL) {
694     retval = MPI_ERR_ARG;
695   } else {
696     *newcomm = smpi_comm_new(smpi_comm_group(comm));
697     retval = MPI_SUCCESS;
698   }
699   smpi_bench_begin(-1, NULL);
700   return retval;
701 }
702
703 int MPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
704 {
705   int retval;
706
707   smpi_bench_end(-1, NULL);
708   if (comm == MPI_COMM_NULL) {
709     retval = MPI_ERR_COMM;
710   } else if (group == MPI_GROUP_NULL) {
711     retval = MPI_ERR_GROUP;
712   } else if (newcomm == NULL) {
713     retval = MPI_ERR_ARG;
714   } else {
715     *newcomm = smpi_comm_new(group);
716     retval = MPI_SUCCESS;
717   }
718   smpi_bench_begin(-1, NULL);
719   return retval;
720 }
721
722 int MPI_Comm_free(MPI_Comm * comm)
723 {
724   int retval;
725
726   smpi_bench_end(-1, NULL);
727   if (comm == NULL) {
728     retval = MPI_ERR_ARG;
729   } else if (*comm == MPI_COMM_NULL) {
730     retval = MPI_ERR_COMM;
731   } else {
732     smpi_comm_destroy(*comm);
733     *comm = MPI_COMM_NULL;
734     retval = MPI_SUCCESS;
735   }
736   smpi_bench_begin(-1, NULL);
737   return retval;
738 }
739
740 int MPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
741                   int tag, MPI_Comm comm, MPI_Request * request)
742 {
743   int retval;
744   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
745
746   smpi_bench_end(rank, "Send_init");
747   if (request == NULL) {
748     retval = MPI_ERR_ARG;
749   } else if (comm == MPI_COMM_NULL) {
750     retval = MPI_ERR_COMM;
751   } else {
752     *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
753     retval = MPI_SUCCESS;
754   }
755   smpi_bench_begin(rank, "Send_init");
756   return retval;
757 }
758
759 int MPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
760                   int tag, MPI_Comm comm, MPI_Request * request)
761 {
762   int retval;
763   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
764
765   smpi_bench_end(rank, "Recv_init");
766   if (request == NULL) {
767     retval = MPI_ERR_ARG;
768   } else if (comm == MPI_COMM_NULL) {
769     retval = MPI_ERR_COMM;
770   } else {
771     *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
772     retval = MPI_SUCCESS;
773   }
774   smpi_bench_begin(rank, "Recv_init");
775   return retval;
776 }
777
778 int MPI_Start(MPI_Request * request)
779 {
780   int retval;
781   MPI_Comm comm = (*request)->comm;
782   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
783
784   smpi_bench_end(rank, "Start");
785   if (request == NULL) {
786     retval = MPI_ERR_ARG;
787   } else {
788     smpi_mpi_start(*request);
789     retval = MPI_SUCCESS;
790   }
791   smpi_bench_begin(rank, "Start");
792   return retval;
793 }
794
795 int MPI_Startall(int count, MPI_Request * requests)
796 {
797   int retval;
798   MPI_Comm comm = count > 0
799       && requests ? requests[0]->comm : MPI_COMM_NULL;
800   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
801
802   smpi_bench_end(rank, "Startall");
803   if (requests == NULL) {
804     retval = MPI_ERR_ARG;
805   } else {
806     smpi_mpi_startall(count, requests);
807     retval = MPI_SUCCESS;
808   }
809   smpi_bench_begin(rank, "Startall");
810   return retval;
811 }
812
813 int MPI_Request_free(MPI_Request * request)
814 {
815   int retval;
816   MPI_Comm comm = (*request)->comm;
817   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
818
819   smpi_bench_end(rank, "Request_free");
820   if (request == NULL) {
821     retval = MPI_ERR_ARG;
822   } else {
823     smpi_mpi_request_free(request);
824     retval = MPI_SUCCESS;
825   }
826   smpi_bench_begin(rank, "Request_free");
827   return retval;
828 }
829
830 int MPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
831               int tag, MPI_Comm comm, MPI_Request * request)
832 {
833   int retval;
834   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
835
836   smpi_bench_end(rank, "Irecv");
837 #ifdef HAVE_TRACING
838   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
839   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
840 #endif
841   if (request == NULL) {
842     retval = MPI_ERR_ARG;
843   } else if (comm == MPI_COMM_NULL) {
844     retval = MPI_ERR_COMM;
845   } else {
846     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
847     retval = MPI_SUCCESS;
848   }
849 #ifdef HAVE_TRACING
850   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
851   (*request)->recv = 1;
852 #endif
853   smpi_bench_begin(rank, "Irecv");
854   return retval;
855 }
856
857 int MPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
858               int tag, MPI_Comm comm, MPI_Request * request)
859 {
860   int retval;
861   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
862
863   smpi_bench_end(rank, "Isend");
864 #ifdef HAVE_TRACING
865   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
866   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
867   TRACE_smpi_send(rank, rank, dst_traced);
868 #endif
869   if (request == NULL) {
870     retval = MPI_ERR_ARG;
871   } else if (comm == MPI_COMM_NULL) {
872     retval = MPI_ERR_COMM;
873   } else {
874     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
875     retval = MPI_SUCCESS;
876   }
877 #ifdef HAVE_TRACING
878   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
879   (*request)->send = 1;
880 #endif
881   smpi_bench_begin(rank, "Isend");
882   return retval;
883 }
884
885 int MPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
886              MPI_Comm comm, MPI_Status * status)
887 {
888   int retval;
889   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
890
891   smpi_bench_end(rank, "Recv");
892 #ifdef HAVE_TRACING
893   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
894   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
895 #endif
896   if (comm == MPI_COMM_NULL) {
897     retval = MPI_ERR_COMM;
898   } else {
899     smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
900     retval = MPI_SUCCESS;
901   }
902 #ifdef HAVE_TRACING
903   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
904   TRACE_smpi_recv(rank, src_traced, rank);
905 #endif
906   smpi_bench_begin(rank, "Recv");
907   return retval;
908 }
909
910 int MPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
911              MPI_Comm comm)
912 {
913   int retval;
914   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
915
916   smpi_bench_end(rank, "Send");
917 #ifdef HAVE_TRACING
918   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
919   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
920   TRACE_smpi_send(rank, rank, dst_traced);
921 #endif
922   if (comm == MPI_COMM_NULL) {
923     retval = MPI_ERR_COMM;
924   } else {
925     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
926     retval = MPI_SUCCESS;
927   }
928 #ifdef HAVE_TRACING
929   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
930 #endif
931   smpi_bench_begin(rank, "Send");
932   return retval;
933 }
934
935 int MPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
936                  int dst, int sendtag, void *recvbuf, int recvcount,
937                  MPI_Datatype recvtype, int src, int recvtag,
938                  MPI_Comm comm, MPI_Status * status)
939 {
940   int retval;
941   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
942
943   smpi_bench_end(rank, "Sendrecv");
944 #ifdef HAVE_TRACING
945   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
946   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
947   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
948   TRACE_smpi_send(rank, rank, dst_traced);
949   TRACE_smpi_send(rank, src_traced, rank);
950 #endif
951   if (comm == MPI_COMM_NULL) {
952     retval = MPI_ERR_COMM;
953   } else if (sendtype == MPI_DATATYPE_NULL
954              || recvtype == MPI_DATATYPE_NULL) {
955     retval = MPI_ERR_TYPE;
956   } else {
957     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
958                       recvcount, recvtype, src, recvtag, comm, status);
959     retval = MPI_SUCCESS;
960   }
961 #ifdef HAVE_TRACING
962   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
963   TRACE_smpi_recv(rank, rank, dst_traced);
964   TRACE_smpi_recv(rank, src_traced, rank);
965 #endif
966   smpi_bench_begin(rank, "Sendrecv");
967   return retval;
968 }
969
970 int MPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
971                          int dst, int sendtag, int src, int recvtag,
972                          MPI_Comm comm, MPI_Status * status)
973 {
974   //TODO: suboptimal implementation
975   void *recvbuf;
976   int retval, size;
977
978   size = smpi_datatype_size(datatype) * count;
979   recvbuf = xbt_new(char, size);
980   retval =
981       MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
982                    datatype, src, recvtag, comm, status);
983   memcpy(buf, recvbuf, size * sizeof(char));
984   xbt_free(recvbuf);
985   return retval;
986 }
987
988 int MPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
989 {
990   int retval;
991   int rank = request && (*request)->comm != MPI_COMM_NULL
992       ? smpi_comm_rank((*request)->comm)
993       : -1;
994
995   smpi_bench_end(rank, "Test");
996   if (request == NULL || flag == NULL) {
997     retval = MPI_ERR_ARG;
998   } else if (*request == MPI_REQUEST_NULL) {
999     retval = MPI_ERR_REQUEST;
1000   } else {
1001     *flag = smpi_mpi_test(request, status);
1002     retval = MPI_SUCCESS;
1003   }
1004   smpi_bench_begin(rank, "Test");
1005   return retval;
1006 }
1007
1008 int MPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1009                 MPI_Status * status)
1010 {
1011   int retval;
1012
1013   smpi_bench_end(-1, NULL);     //FIXME
1014   if (index == NULL || flag == NULL) {
1015     retval = MPI_ERR_ARG;
1016   } else {
1017     *flag = smpi_mpi_testany(count, requests, index, status);
1018     retval = MPI_SUCCESS;
1019   }
1020   smpi_bench_begin(-1, NULL);
1021   return retval;
1022 }
1023
1024 int MPI_Wait(MPI_Request * request, MPI_Status * status)
1025 {
1026   int retval;
1027   int rank = request && (*request)->comm != MPI_COMM_NULL
1028       ? smpi_comm_rank((*request)->comm)
1029       : -1;
1030
1031   smpi_bench_end(rank, "Wait");
1032 #ifdef HAVE_TRACING
1033   MPI_Group group = smpi_comm_group((*request)->comm);
1034   int src_traced = smpi_group_rank(group, (*request)->src);
1035   int dst_traced = smpi_group_rank(group, (*request)->dst);
1036   int is_wait_for_receive = (*request)->recv;
1037   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1038 #endif
1039   if (request == NULL) {
1040     retval = MPI_ERR_ARG;
1041   } else if (*request == MPI_REQUEST_NULL) {
1042     retval = MPI_ERR_REQUEST;
1043   } else {
1044     smpi_mpi_wait(request, status);
1045     retval = MPI_SUCCESS;
1046   }
1047 #ifdef HAVE_TRACING
1048   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1049   if (is_wait_for_receive) {
1050     TRACE_smpi_recv(rank, src_traced, dst_traced);
1051   }
1052 #endif
1053   smpi_bench_begin(rank, "Wait");
1054   return retval;
1055 }
1056
1057 int MPI_Waitany(int count, MPI_Request requests[], int *index,
1058                 MPI_Status * status)
1059 {
1060   int retval;
1061
1062   smpi_bench_end(-1, NULL);     //FIXME
1063 #ifdef HAVE_TRACING
1064   //save requests information for tracing
1065   int i;
1066   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
1067   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
1068   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
1069   for (i = 0; i < count; i++) {
1070     MPI_Request req = requests[i];      //already received requests are no longer valid
1071     if (req) {
1072       int *asrc = xbt_new(int, 1);
1073       int *adst = xbt_new(int, 1);
1074       int *arecv = xbt_new(int, 1);
1075       *asrc = req->src;
1076       *adst = req->dst;
1077       *arecv = req->recv;
1078       xbt_dynar_insert_at(srcs, i, asrc);
1079       xbt_dynar_insert_at(dsts, i, adst);
1080       xbt_dynar_insert_at(recvs, i, arecv);
1081     } else {
1082       int *t = xbt_new(int, 1);
1083       xbt_dynar_insert_at(srcs, i, t);
1084       xbt_dynar_insert_at(dsts, i, t);
1085       xbt_dynar_insert_at(recvs, i, t);
1086     }
1087   }
1088   int rank_traced = smpi_comm_rank(MPI_COMM_WORLD);
1089   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1090 #endif
1091   if (index == NULL) {
1092     retval = MPI_ERR_ARG;
1093   } else {
1094     *index = smpi_mpi_waitany(count, requests, status);
1095     retval = MPI_SUCCESS;
1096   }
1097 #ifdef HAVE_TRACING
1098   int src_traced, dst_traced, is_wait_for_receive;
1099   xbt_dynar_get_cpy(srcs, *index, &src_traced);
1100   xbt_dynar_get_cpy(dsts, *index, &dst_traced);
1101   xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
1102   if (is_wait_for_receive) {
1103     TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1104   }
1105   TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1106   //clean-up of dynars
1107   xbt_free(srcs);
1108   xbt_free(dsts);
1109   xbt_free(recvs);
1110 #endif
1111   smpi_bench_begin(-1, NULL);
1112   return retval;
1113 }
1114
1115 int MPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1116 {
1117
1118   smpi_bench_end(-1, NULL);     //FIXME
1119 #ifdef HAVE_TRACING
1120   //save information from requests
1121   int i;
1122   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
1123   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
1124   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
1125   for (i = 0; i < count; i++) {
1126     MPI_Request req = requests[i];      //all req should be valid in Waitall
1127     int *asrc = xbt_new(int, 1);
1128     int *adst = xbt_new(int, 1);
1129     int *arecv = xbt_new(int, 1);
1130     *asrc = req->src;
1131     *adst = req->dst;
1132     *arecv = req->recv;
1133     xbt_dynar_insert_at(srcs, i, asrc);
1134     xbt_dynar_insert_at(dsts, i, adst);
1135     xbt_dynar_insert_at(recvs, i, arecv);
1136   }
1137   int rank_traced = smpi_comm_rank (MPI_COMM_WORLD);
1138   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1139 #endif
1140   smpi_mpi_waitall(count, requests, status);
1141 #ifdef HAVE_TRACING
1142   for (i = 0; i < count; i++) {
1143     int src_traced, dst_traced, is_wait_for_receive;
1144     xbt_dynar_get_cpy(srcs, i, &src_traced);
1145     xbt_dynar_get_cpy(dsts, i, &dst_traced);
1146     xbt_dynar_get_cpy(recvs, i, &is_wait_for_receive);
1147     if (is_wait_for_receive) {
1148       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1149     }
1150   }
1151   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1152   //clean-up of dynars
1153   xbt_free(srcs);
1154   xbt_free(dsts);
1155   xbt_free(recvs);
1156 #endif
1157   smpi_bench_begin(-1, NULL);
1158   return MPI_SUCCESS;
1159 }
1160
1161 int MPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1162                  int *indices, MPI_Status status[])
1163 {
1164   int retval;
1165
1166   smpi_bench_end(-1, NULL);     //FIXME
1167   if (outcount == NULL || indices == NULL) {
1168     retval = MPI_ERR_ARG;
1169   } else {
1170     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1171     retval = MPI_SUCCESS;
1172   }
1173   smpi_bench_begin(-1, NULL);
1174   return retval;
1175 }
1176
1177 int MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
1178               MPI_Comm comm)
1179 {
1180   int retval;
1181   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1182
1183   smpi_bench_end(rank, "Bcast");
1184 #ifdef HAVE_TRACING
1185   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1186   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1187 #endif
1188   if (comm == MPI_COMM_NULL) {
1189     retval = MPI_ERR_COMM;
1190   } else {
1191     smpi_mpi_bcast(buf, count, datatype, root, comm);
1192     retval = MPI_SUCCESS;
1193   }
1194 #ifdef HAVE_TRACING
1195   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1196 #endif
1197   smpi_bench_begin(rank, "Bcast");
1198   return retval;
1199 }
1200
1201 int MPI_Barrier(MPI_Comm comm)
1202 {
1203   int retval;
1204   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1205
1206   smpi_bench_end(rank, "Barrier");
1207 #ifdef HAVE_TRACING
1208   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1209 #endif
1210   if (comm == MPI_COMM_NULL) {
1211     retval = MPI_ERR_COMM;
1212   } else {
1213     smpi_mpi_barrier(comm);
1214     retval = MPI_SUCCESS;
1215   }
1216 #ifdef HAVE_TRACING
1217   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1218 #endif
1219   smpi_bench_begin(rank, "Barrier");
1220   return retval;
1221 }
1222
1223 int MPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1224                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1225                int root, MPI_Comm comm)
1226 {
1227   int retval;
1228   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1229
1230   smpi_bench_end(rank, "Gather");
1231 #ifdef HAVE_TRACING
1232   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1233   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1234 #endif
1235   if (comm == MPI_COMM_NULL) {
1236     retval = MPI_ERR_COMM;
1237   } else if (sendtype == MPI_DATATYPE_NULL
1238              || recvtype == MPI_DATATYPE_NULL) {
1239     retval = MPI_ERR_TYPE;
1240   } else {
1241     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1242                     recvtype, root, comm);
1243     retval = MPI_SUCCESS;
1244   }
1245 #ifdef HAVE_TRACING
1246   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1247 #endif
1248   smpi_bench_begin(rank, "Gather");
1249   return retval;
1250 }
1251
1252 int MPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1253                 void *recvbuf, int *recvcounts, int *displs,
1254                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1255 {
1256   int retval;
1257   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1258
1259   smpi_bench_end(rank, "Gatherv");
1260 #ifdef HAVE_TRACING
1261   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1262   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1263 #endif
1264   if (comm == MPI_COMM_NULL) {
1265     retval = MPI_ERR_COMM;
1266   } else if (sendtype == MPI_DATATYPE_NULL
1267              || recvtype == MPI_DATATYPE_NULL) {
1268     retval = MPI_ERR_TYPE;
1269   } else if (recvcounts == NULL || displs == NULL) {
1270     retval = MPI_ERR_ARG;
1271   } else {
1272     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1273                      displs, recvtype, root, comm);
1274     retval = MPI_SUCCESS;
1275   }
1276 #ifdef HAVE_TRACING
1277   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1278 #endif
1279   smpi_bench_begin(rank, "Gatherv");
1280   return retval;
1281 }
1282
1283 int MPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1284                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1285                   MPI_Comm comm)
1286 {
1287   int retval;
1288   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1289
1290   smpi_bench_end(rank, "Allgather");
1291 #ifdef HAVE_TRACING
1292   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1293 #endif
1294   if (comm == MPI_COMM_NULL) {
1295     retval = MPI_ERR_COMM;
1296   } else if (sendtype == MPI_DATATYPE_NULL
1297              || recvtype == MPI_DATATYPE_NULL) {
1298     retval = MPI_ERR_TYPE;
1299   } else {
1300     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1301                        recvtype, comm);
1302     retval = MPI_SUCCESS;
1303   }
1304 #ifdef HAVE_TRACING
1305   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1306 #endif
1307   smpi_bench_begin(rank, "Allgather");
1308   return retval;
1309 }
1310
1311 int MPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1312                    void *recvbuf, int *recvcounts, int *displs,
1313                    MPI_Datatype recvtype, MPI_Comm comm)
1314 {
1315   int retval;
1316   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1317
1318   smpi_bench_end(rank, "Allgatherv");
1319 #ifdef HAVE_TRACING
1320   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1321 #endif
1322   if (comm == MPI_COMM_NULL) {
1323     retval = MPI_ERR_COMM;
1324   } else if (sendtype == MPI_DATATYPE_NULL
1325              || recvtype == MPI_DATATYPE_NULL) {
1326     retval = MPI_ERR_TYPE;
1327   } else if (recvcounts == NULL || displs == NULL) {
1328     retval = MPI_ERR_ARG;
1329   } else {
1330     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1331                         displs, recvtype, comm);
1332     retval = MPI_SUCCESS;
1333   }
1334 #ifdef HAVE_TRACING
1335   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1336 #endif
1337   smpi_bench_begin(rank, "Allgatherv");
1338   return retval;
1339 }
1340
1341 int MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1342                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1343                 int root, MPI_Comm comm)
1344 {
1345   int retval;
1346   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1347
1348   smpi_bench_end(rank, "Scatter");
1349 #ifdef HAVE_TRACING
1350   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1351   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1352 #endif
1353   if (comm == MPI_COMM_NULL) {
1354     retval = MPI_ERR_COMM;
1355   } else if (sendtype == MPI_DATATYPE_NULL
1356              || recvtype == MPI_DATATYPE_NULL) {
1357     retval = MPI_ERR_TYPE;
1358   } else {
1359     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1360                      recvtype, root, comm);
1361     retval = MPI_SUCCESS;
1362   }
1363 #ifdef HAVE_TRACING
1364   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1365 #endif
1366   smpi_bench_begin(rank, "Scatter");
1367   return retval;
1368 }
1369
1370 int MPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1371                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1372                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1373 {
1374   int retval;
1375   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1376
1377   smpi_bench_end(rank, "Scatterv");
1378 #ifdef HAVE_TRACING
1379   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1380   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1381 #endif
1382   if (comm == MPI_COMM_NULL) {
1383     retval = MPI_ERR_COMM;
1384   } else if (sendtype == MPI_DATATYPE_NULL
1385              || recvtype == MPI_DATATYPE_NULL) {
1386     retval = MPI_ERR_TYPE;
1387   } else if (sendcounts == NULL || displs == NULL) {
1388     retval = MPI_ERR_ARG;
1389   } else {
1390     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1391                       recvcount, recvtype, root, comm);
1392     retval = MPI_SUCCESS;
1393   }
1394 #ifdef HAVE_TRACING
1395   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1396 #endif
1397   smpi_bench_begin(rank, "Scatterv");
1398   return retval;
1399 }
1400
1401 int MPI_Reduce(void *sendbuf, void *recvbuf, int count,
1402                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1403 {
1404   int retval;
1405   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1406
1407   smpi_bench_end(rank, "Reduce");
1408 #ifdef HAVE_TRACING
1409   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1410   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1411 #endif
1412   if (comm == MPI_COMM_NULL) {
1413     retval = MPI_ERR_COMM;
1414   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1415     retval = MPI_ERR_ARG;
1416   } else {
1417     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1418     retval = MPI_SUCCESS;
1419   }
1420 #ifdef HAVE_TRACING
1421   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1422 #endif
1423   smpi_bench_begin(rank, "Reduce");
1424   return retval;
1425 }
1426
1427 int MPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1428                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1429 {
1430   int retval;
1431   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1432
1433   smpi_bench_end(rank, "Allreduce");
1434 #ifdef HAVE_TRACING
1435   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1436 #endif
1437   if (comm == MPI_COMM_NULL) {
1438     retval = MPI_ERR_COMM;
1439   } else if (datatype == MPI_DATATYPE_NULL) {
1440     retval = MPI_ERR_TYPE;
1441   } else if (op == MPI_OP_NULL) {
1442     retval = MPI_ERR_OP;
1443   } else {
1444     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
1445     retval = MPI_SUCCESS;
1446   }
1447 #ifdef HAVE_TRACING
1448   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1449 #endif
1450   smpi_bench_begin(rank, "Allreduce");
1451   return retval;
1452 }
1453
1454 int MPI_Scan(void *sendbuf, void *recvbuf, int count,
1455              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1456 {
1457   int retval;
1458   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1459
1460   smpi_bench_end(rank, "Scan");
1461 #ifdef HAVE_TRACING
1462   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1463 #endif
1464   if (comm == MPI_COMM_NULL) {
1465     retval = MPI_ERR_COMM;
1466   } else if (datatype == MPI_DATATYPE_NULL) {
1467     retval = MPI_ERR_TYPE;
1468   } else if (op == MPI_OP_NULL) {
1469     retval = MPI_ERR_OP;
1470   } else {
1471     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1472     retval = MPI_SUCCESS;
1473   }
1474 #ifdef HAVE_TRACING
1475   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1476 #endif
1477   smpi_bench_begin(rank, "Scan");
1478   return retval;
1479 }
1480
1481 int MPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1482                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1483 {
1484   int retval, i, size, count;
1485   int *displs;
1486   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1487
1488   smpi_bench_end(rank, "Reduce_scatter");
1489 #ifdef HAVE_TRACING
1490   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1491 #endif
1492   if (comm == MPI_COMM_NULL) {
1493     retval = MPI_ERR_COMM;
1494   } else if (datatype == MPI_DATATYPE_NULL) {
1495     retval = MPI_ERR_TYPE;
1496   } else if (op == MPI_OP_NULL) {
1497     retval = MPI_ERR_OP;
1498   } else if (recvcounts == NULL) {
1499     retval = MPI_ERR_ARG;
1500   } else {
1501     /* arbitrarily choose root as rank 0 */
1502     /* TODO: faster direct implementation ? */
1503     size = smpi_comm_size(comm);
1504     count = 0;
1505     displs = xbt_new(int, size);
1506     for (i = 0; i < size; i++) {
1507       count += recvcounts[i];
1508       displs[i] = 0;
1509     }
1510     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1511     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf,
1512                       recvcounts[rank], datatype, 0, comm);
1513     xbt_free(displs);
1514     retval = MPI_SUCCESS;
1515   }
1516 #ifdef HAVE_TRACING
1517   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1518 #endif
1519   smpi_bench_begin(rank, "Reduce_scatter");
1520   return retval;
1521 }
1522
1523 int MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1524                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
1525                  MPI_Comm comm)
1526 {
1527   int retval, size, sendsize;
1528   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1529
1530   smpi_bench_end(rank, "Alltoall");
1531 #ifdef HAVE_TRACING
1532   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1533 #endif
1534   if (comm == MPI_COMM_NULL) {
1535     retval = MPI_ERR_COMM;
1536   } else if (sendtype == MPI_DATATYPE_NULL
1537              || recvtype == MPI_DATATYPE_NULL) {
1538     retval = MPI_ERR_TYPE;
1539   } else {
1540     size = smpi_comm_size(comm);
1541     sendsize = smpi_datatype_size(sendtype) * sendcount;
1542     if (sendsize < 200 && size > 12) {
1543       retval =
1544           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
1545                                          recvbuf, recvcount, recvtype,
1546                                          comm);
1547     } else if (sendsize < 3000) {
1548       retval =
1549           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
1550                                                 sendtype, recvbuf,
1551                                                 recvcount, recvtype, comm);
1552     } else {
1553       retval =
1554           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
1555                                             recvbuf, recvcount, recvtype,
1556                                             comm);
1557     }
1558   }
1559 #ifdef HAVE_TRACING
1560   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1561 #endif
1562   smpi_bench_begin(rank, "Alltoall");
1563   return retval;
1564 }
1565
1566 int MPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
1567                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
1568                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
1569 {
1570   int retval;
1571   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1572
1573   smpi_bench_end(rank, "Alltoallv");
1574 #ifdef HAVE_TRACING
1575   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1576 #endif
1577   if (comm == MPI_COMM_NULL) {
1578     retval = MPI_ERR_COMM;
1579   } else if (sendtype == MPI_DATATYPE_NULL
1580              || recvtype == MPI_DATATYPE_NULL) {
1581     retval = MPI_ERR_TYPE;
1582   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
1583              || recvdisps == NULL) {
1584     retval = MPI_ERR_ARG;
1585   } else {
1586     retval =
1587         smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
1588                                   recvbuf, recvcounts, recvdisps, recvtype,
1589                                   comm);
1590   }
1591 #ifdef HAVE_TRACING
1592   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1593 #endif
1594   smpi_bench_begin(rank, "Alltoallv");
1595   return retval;
1596 }
1597
1598
1599 int MPI_Get_processor_name(char *name, int *resultlen)
1600 {
1601   int retval = MPI_SUCCESS;
1602
1603   smpi_bench_end(-1, NULL);
1604   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
1605           MPI_MAX_PROCESSOR_NAME - 1);
1606   *resultlen =
1607       strlen(name) >
1608       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1609
1610   smpi_bench_begin(-1, NULL);
1611   return retval;
1612 }
1613
1614 int MPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
1615 {
1616   int retval = MPI_SUCCESS;
1617   size_t size;
1618
1619   smpi_bench_end(-1, NULL);
1620   if (status == NULL || count == NULL) {
1621     retval = MPI_ERR_ARG;
1622   } else if (datatype == MPI_DATATYPE_NULL) {
1623     retval = MPI_ERR_TYPE;
1624   } else {
1625     size = smpi_datatype_size(datatype);
1626     if (size == 0) {
1627       *count = 0;
1628     } else if (status->count % size != 0) {
1629       retval = MPI_UNDEFINED;
1630     } else {
1631       *count = smpi_mpi_get_count(status, datatype);
1632     }
1633   }
1634   smpi_bench_begin(-1, NULL);
1635   return retval;
1636 }