Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
0f08ff8b356239a41d4a23eb81e7db6bf7ad8e86
[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
1089   //search for a suitable request to give the rank of current mpi proc
1090   MPI_Request req = NULL;
1091   for (i = 0; i < count && req == NULL; i++) {
1092     req = requests[i];
1093   }
1094   MPI_Comm comm = (req)->comm;
1095   int rank_traced = smpi_comm_rank(comm);
1096   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1097 #endif
1098   if (index == NULL) {
1099     retval = MPI_ERR_ARG;
1100   } else {
1101     *index = smpi_mpi_waitany(count, requests, status);
1102     retval = MPI_SUCCESS;
1103   }
1104 #ifdef HAVE_TRACING
1105   int src_traced, dst_traced, is_wait_for_receive;
1106   xbt_dynar_get_cpy(srcs, *index, &src_traced);
1107   xbt_dynar_get_cpy(dsts, *index, &dst_traced);
1108   xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
1109   if (is_wait_for_receive) {
1110     TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1111   }
1112   TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1113   //clean-up of dynars
1114   xbt_free(srcs);
1115   xbt_free(dsts);
1116   xbt_free(recvs);
1117 #endif
1118   smpi_bench_begin(-1, NULL);
1119   return retval;
1120 }
1121
1122 int MPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1123 {
1124
1125   smpi_bench_end(-1, NULL);     //FIXME
1126 #ifdef HAVE_TRACING
1127   //save information from requests
1128   int i;
1129   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
1130   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
1131   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
1132   for (i = 0; i < count; i++) {
1133     MPI_Request req = requests[i];      //all req should be valid in Waitall
1134     int *asrc = xbt_new(int, 1);
1135     int *adst = xbt_new(int, 1);
1136     int *arecv = xbt_new(int, 1);
1137     *asrc = req->src;
1138     *adst = req->dst;
1139     *arecv = req->recv;
1140     xbt_dynar_insert_at(srcs, i, asrc);
1141     xbt_dynar_insert_at(dsts, i, adst);
1142     xbt_dynar_insert_at(recvs, i, arecv);
1143   }
1144
1145 //  find my rank inside one of MPI_Comm's of the requests
1146   MPI_Request req = NULL;
1147   for (i = 0; i < count && req == NULL; i++) {
1148     req = requests[i];
1149   }
1150   MPI_Comm comm = (req)->comm;
1151   int rank_traced = smpi_comm_rank(comm);
1152   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1153 #endif
1154   smpi_mpi_waitall(count, requests, status);
1155 #ifdef HAVE_TRACING
1156   for (i = 0; i < count; i++) {
1157     int src_traced, dst_traced, is_wait_for_receive;
1158     xbt_dynar_get_cpy(srcs, i, &src_traced);
1159     xbt_dynar_get_cpy(dsts, i, &dst_traced);
1160     xbt_dynar_get_cpy(recvs, i, &is_wait_for_receive);
1161     if (is_wait_for_receive) {
1162       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1163     }
1164   }
1165   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1166   //clean-up of dynars
1167   xbt_free(srcs);
1168   xbt_free(dsts);
1169   xbt_free(recvs);
1170 #endif
1171   smpi_bench_begin(-1, NULL);
1172   return MPI_SUCCESS;
1173 }
1174
1175 int MPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1176                  int *indices, MPI_Status status[])
1177 {
1178   int retval;
1179
1180   smpi_bench_end(-1, NULL);     //FIXME
1181   if (outcount == NULL || indices == NULL) {
1182     retval = MPI_ERR_ARG;
1183   } else {
1184     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1185     retval = MPI_SUCCESS;
1186   }
1187   smpi_bench_begin(-1, NULL);
1188   return retval;
1189 }
1190
1191 int MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
1192               MPI_Comm comm)
1193 {
1194   int retval;
1195   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1196
1197   smpi_bench_end(rank, "Bcast");
1198 #ifdef HAVE_TRACING
1199   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1200   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1201 #endif
1202   if (comm == MPI_COMM_NULL) {
1203     retval = MPI_ERR_COMM;
1204   } else {
1205     smpi_mpi_bcast(buf, count, datatype, root, comm);
1206     retval = MPI_SUCCESS;
1207   }
1208 #ifdef HAVE_TRACING
1209   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1210 #endif
1211   smpi_bench_begin(rank, "Bcast");
1212   return retval;
1213 }
1214
1215 int MPI_Barrier(MPI_Comm comm)
1216 {
1217   int retval;
1218   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1219
1220   smpi_bench_end(rank, "Barrier");
1221 #ifdef HAVE_TRACING
1222   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1223 #endif
1224   if (comm == MPI_COMM_NULL) {
1225     retval = MPI_ERR_COMM;
1226   } else {
1227     smpi_mpi_barrier(comm);
1228     retval = MPI_SUCCESS;
1229   }
1230 #ifdef HAVE_TRACING
1231   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1232 #endif
1233   smpi_bench_begin(rank, "Barrier");
1234   return retval;
1235 }
1236
1237 int MPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1238                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1239                int root, MPI_Comm comm)
1240 {
1241   int retval;
1242   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1243
1244   smpi_bench_end(rank, "Gather");
1245 #ifdef HAVE_TRACING
1246   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1247   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1248 #endif
1249   if (comm == MPI_COMM_NULL) {
1250     retval = MPI_ERR_COMM;
1251   } else if (sendtype == MPI_DATATYPE_NULL
1252              || recvtype == MPI_DATATYPE_NULL) {
1253     retval = MPI_ERR_TYPE;
1254   } else {
1255     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1256                     recvtype, root, comm);
1257     retval = MPI_SUCCESS;
1258   }
1259 #ifdef HAVE_TRACING
1260   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1261 #endif
1262   smpi_bench_begin(rank, "Gather");
1263   return retval;
1264 }
1265
1266 int MPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1267                 void *recvbuf, int *recvcounts, int *displs,
1268                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1269 {
1270   int retval;
1271   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1272
1273   smpi_bench_end(rank, "Gatherv");
1274 #ifdef HAVE_TRACING
1275   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1276   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1277 #endif
1278   if (comm == MPI_COMM_NULL) {
1279     retval = MPI_ERR_COMM;
1280   } else if (sendtype == MPI_DATATYPE_NULL
1281              || recvtype == MPI_DATATYPE_NULL) {
1282     retval = MPI_ERR_TYPE;
1283   } else if (recvcounts == NULL || displs == NULL) {
1284     retval = MPI_ERR_ARG;
1285   } else {
1286     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1287                      displs, recvtype, root, comm);
1288     retval = MPI_SUCCESS;
1289   }
1290 #ifdef HAVE_TRACING
1291   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1292 #endif
1293   smpi_bench_begin(rank, "Gatherv");
1294   return retval;
1295 }
1296
1297 int MPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1298                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1299                   MPI_Comm comm)
1300 {
1301   int retval;
1302   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1303
1304   smpi_bench_end(rank, "Allgather");
1305 #ifdef HAVE_TRACING
1306   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1307 #endif
1308   if (comm == MPI_COMM_NULL) {
1309     retval = MPI_ERR_COMM;
1310   } else if (sendtype == MPI_DATATYPE_NULL
1311              || recvtype == MPI_DATATYPE_NULL) {
1312     retval = MPI_ERR_TYPE;
1313   } else {
1314     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1315                        recvtype, comm);
1316     retval = MPI_SUCCESS;
1317   }
1318 #ifdef HAVE_TRACING
1319   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1320 #endif
1321   smpi_bench_begin(rank, "Allgather");
1322   return retval;
1323 }
1324
1325 int MPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1326                    void *recvbuf, int *recvcounts, int *displs,
1327                    MPI_Datatype recvtype, MPI_Comm comm)
1328 {
1329   int retval;
1330   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1331
1332   smpi_bench_end(rank, "Allgatherv");
1333 #ifdef HAVE_TRACING
1334   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1335 #endif
1336   if (comm == MPI_COMM_NULL) {
1337     retval = MPI_ERR_COMM;
1338   } else if (sendtype == MPI_DATATYPE_NULL
1339              || recvtype == MPI_DATATYPE_NULL) {
1340     retval = MPI_ERR_TYPE;
1341   } else if (recvcounts == NULL || displs == NULL) {
1342     retval = MPI_ERR_ARG;
1343   } else {
1344     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1345                         displs, recvtype, comm);
1346     retval = MPI_SUCCESS;
1347   }
1348 #ifdef HAVE_TRACING
1349   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1350 #endif
1351   smpi_bench_begin(rank, "Allgatherv");
1352   return retval;
1353 }
1354
1355 int MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1356                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1357                 int root, MPI_Comm comm)
1358 {
1359   int retval;
1360   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1361
1362   smpi_bench_end(rank, "Scatter");
1363 #ifdef HAVE_TRACING
1364   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1365   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1366 #endif
1367   if (comm == MPI_COMM_NULL) {
1368     retval = MPI_ERR_COMM;
1369   } else if (sendtype == MPI_DATATYPE_NULL
1370              || recvtype == MPI_DATATYPE_NULL) {
1371     retval = MPI_ERR_TYPE;
1372   } else {
1373     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1374                      recvtype, root, comm);
1375     retval = MPI_SUCCESS;
1376   }
1377 #ifdef HAVE_TRACING
1378   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1379 #endif
1380   smpi_bench_begin(rank, "Scatter");
1381   return retval;
1382 }
1383
1384 int MPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1385                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1386                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1387 {
1388   int retval;
1389   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1390
1391   smpi_bench_end(rank, "Scatterv");
1392 #ifdef HAVE_TRACING
1393   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1394   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1395 #endif
1396   if (comm == MPI_COMM_NULL) {
1397     retval = MPI_ERR_COMM;
1398   } else if (sendtype == MPI_DATATYPE_NULL
1399              || recvtype == MPI_DATATYPE_NULL) {
1400     retval = MPI_ERR_TYPE;
1401   } else if (sendcounts == NULL || displs == NULL) {
1402     retval = MPI_ERR_ARG;
1403   } else {
1404     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1405                       recvcount, recvtype, root, comm);
1406     retval = MPI_SUCCESS;
1407   }
1408 #ifdef HAVE_TRACING
1409   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1410 #endif
1411   smpi_bench_begin(rank, "Scatterv");
1412   return retval;
1413 }
1414
1415 int MPI_Reduce(void *sendbuf, void *recvbuf, int count,
1416                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1417 {
1418   int retval;
1419   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1420
1421   smpi_bench_end(rank, "Reduce");
1422 #ifdef HAVE_TRACING
1423   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1424   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1425 #endif
1426   if (comm == MPI_COMM_NULL) {
1427     retval = MPI_ERR_COMM;
1428   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1429     retval = MPI_ERR_ARG;
1430   } else {
1431     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1432     retval = MPI_SUCCESS;
1433   }
1434 #ifdef HAVE_TRACING
1435   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1436 #endif
1437   smpi_bench_begin(rank, "Reduce");
1438   return retval;
1439 }
1440
1441 int MPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1442                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1443 {
1444   int retval;
1445   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1446
1447   smpi_bench_end(rank, "Allreduce");
1448 #ifdef HAVE_TRACING
1449   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1450 #endif
1451   if (comm == MPI_COMM_NULL) {
1452     retval = MPI_ERR_COMM;
1453   } else if (datatype == MPI_DATATYPE_NULL) {
1454     retval = MPI_ERR_TYPE;
1455   } else if (op == MPI_OP_NULL) {
1456     retval = MPI_ERR_OP;
1457   } else {
1458     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
1459     retval = MPI_SUCCESS;
1460   }
1461 #ifdef HAVE_TRACING
1462   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1463 #endif
1464   smpi_bench_begin(rank, "Allreduce");
1465   return retval;
1466 }
1467
1468 int MPI_Scan(void *sendbuf, void *recvbuf, int count,
1469              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1470 {
1471   int retval;
1472   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1473
1474   smpi_bench_end(rank, "Scan");
1475 #ifdef HAVE_TRACING
1476   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1477 #endif
1478   if (comm == MPI_COMM_NULL) {
1479     retval = MPI_ERR_COMM;
1480   } else if (datatype == MPI_DATATYPE_NULL) {
1481     retval = MPI_ERR_TYPE;
1482   } else if (op == MPI_OP_NULL) {
1483     retval = MPI_ERR_OP;
1484   } else {
1485     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1486     retval = MPI_SUCCESS;
1487   }
1488 #ifdef HAVE_TRACING
1489   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1490 #endif
1491   smpi_bench_begin(rank, "Scan");
1492   return retval;
1493 }
1494
1495 int MPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1496                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1497 {
1498   int retval, i, size, count;
1499   int *displs;
1500   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1501
1502   smpi_bench_end(rank, "Reduce_scatter");
1503 #ifdef HAVE_TRACING
1504   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1505 #endif
1506   if (comm == MPI_COMM_NULL) {
1507     retval = MPI_ERR_COMM;
1508   } else if (datatype == MPI_DATATYPE_NULL) {
1509     retval = MPI_ERR_TYPE;
1510   } else if (op == MPI_OP_NULL) {
1511     retval = MPI_ERR_OP;
1512   } else if (recvcounts == NULL) {
1513     retval = MPI_ERR_ARG;
1514   } else {
1515     /* arbitrarily choose root as rank 0 */
1516     /* TODO: faster direct implementation ? */
1517     size = smpi_comm_size(comm);
1518     count = 0;
1519     displs = xbt_new(int, size);
1520     for (i = 0; i < size; i++) {
1521       count += recvcounts[i];
1522       displs[i] = 0;
1523     }
1524     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1525     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf,
1526                       recvcounts[rank], datatype, 0, comm);
1527     xbt_free(displs);
1528     retval = MPI_SUCCESS;
1529   }
1530 #ifdef HAVE_TRACING
1531   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1532 #endif
1533   smpi_bench_begin(rank, "Reduce_scatter");
1534   return retval;
1535 }
1536
1537 int MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1538                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
1539                  MPI_Comm comm)
1540 {
1541   int retval, size, sendsize;
1542   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1543
1544   smpi_bench_end(rank, "Alltoall");
1545 #ifdef HAVE_TRACING
1546   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1547 #endif
1548   if (comm == MPI_COMM_NULL) {
1549     retval = MPI_ERR_COMM;
1550   } else if (sendtype == MPI_DATATYPE_NULL
1551              || recvtype == MPI_DATATYPE_NULL) {
1552     retval = MPI_ERR_TYPE;
1553   } else {
1554     size = smpi_comm_size(comm);
1555     sendsize = smpi_datatype_size(sendtype) * sendcount;
1556     if (sendsize < 200 && size > 12) {
1557       retval =
1558           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
1559                                          recvbuf, recvcount, recvtype,
1560                                          comm);
1561     } else if (sendsize < 3000) {
1562       retval =
1563           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
1564                                                 sendtype, recvbuf,
1565                                                 recvcount, recvtype, comm);
1566     } else {
1567       retval =
1568           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
1569                                             recvbuf, recvcount, recvtype,
1570                                             comm);
1571     }
1572   }
1573 #ifdef HAVE_TRACING
1574   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1575 #endif
1576   smpi_bench_begin(rank, "Alltoall");
1577   return retval;
1578 }
1579
1580 int MPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
1581                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
1582                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
1583 {
1584   int retval;
1585   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1586
1587   smpi_bench_end(rank, "Alltoallv");
1588 #ifdef HAVE_TRACING
1589   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1590 #endif
1591   if (comm == MPI_COMM_NULL) {
1592     retval = MPI_ERR_COMM;
1593   } else if (sendtype == MPI_DATATYPE_NULL
1594              || recvtype == MPI_DATATYPE_NULL) {
1595     retval = MPI_ERR_TYPE;
1596   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
1597              || recvdisps == NULL) {
1598     retval = MPI_ERR_ARG;
1599   } else {
1600     retval =
1601         smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
1602                                   recvbuf, recvcounts, recvdisps, recvtype,
1603                                   comm);
1604   }
1605 #ifdef HAVE_TRACING
1606   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1607 #endif
1608   smpi_bench_begin(rank, "Alltoallv");
1609   return retval;
1610 }
1611
1612
1613 int MPI_Get_processor_name(char *name, int *resultlen)
1614 {
1615   int retval = MPI_SUCCESS;
1616
1617   smpi_bench_end(-1, NULL);
1618   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
1619           MPI_MAX_PROCESSOR_NAME - 1);
1620   *resultlen =
1621       strlen(name) >
1622       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1623
1624   smpi_bench_begin(-1, NULL);
1625   return retval;
1626 }
1627
1628 int MPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
1629 {
1630   int retval = MPI_SUCCESS;
1631   size_t size;
1632
1633   smpi_bench_end(-1, NULL);
1634   if (status == NULL || count == NULL) {
1635     retval = MPI_ERR_ARG;
1636   } else if (datatype == MPI_DATATYPE_NULL) {
1637     retval = MPI_ERR_TYPE;
1638   } else {
1639     size = smpi_datatype_size(datatype);
1640     if (size == 0) {
1641       *count = 0;
1642     } else if (status->count % size != 0) {
1643       retval = MPI_UNDEFINED;
1644     } else {
1645       *count = smpi_mpi_get_count(status, datatype);
1646     }
1647   }
1648   smpi_bench_begin(-1, NULL);
1649   return retval;
1650 }