Logo AND Algorithmique Numérique Distribuée

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