Logo AND Algorithmique Numérique Distribuée

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