Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
39b292669e740e7817fe0bd59eb3b6ebf3c1c6f6
[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_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
766                   int tag, MPI_Comm comm, MPI_Request * request)
767 {
768   int retval;
769   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
770
771   smpi_bench_end(rank, "Send_init");
772   if (request == NULL) {
773     retval = MPI_ERR_ARG;
774   } else if (comm == MPI_COMM_NULL) {
775     retval = MPI_ERR_COMM;
776   } else {
777     *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
778     retval = MPI_SUCCESS;
779   }
780   smpi_bench_begin(rank, "Send_init");
781   return retval;
782 }
783
784 int MPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
785                   int tag, MPI_Comm comm, MPI_Request * request)
786 {
787   int retval;
788   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
789
790   smpi_bench_end(rank, "Recv_init");
791   if (request == NULL) {
792     retval = MPI_ERR_ARG;
793   } else if (comm == MPI_COMM_NULL) {
794     retval = MPI_ERR_COMM;
795   } else {
796     *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
797     retval = MPI_SUCCESS;
798   }
799   smpi_bench_begin(rank, "Recv_init");
800   return retval;
801 }
802
803 int MPI_Start(MPI_Request * request)
804 {
805   int retval;
806   MPI_Comm comm = (*request)->comm;
807   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
808
809   smpi_bench_end(rank, "Start");
810   if (request == NULL) {
811     retval = MPI_ERR_ARG;
812   } else {
813     smpi_mpi_start(*request);
814     retval = MPI_SUCCESS;
815   }
816   smpi_bench_begin(rank, "Start");
817   return retval;
818 }
819
820 int MPI_Startall(int count, MPI_Request * requests)
821 {
822   int retval;
823   MPI_Comm comm = count > 0
824       && requests ? requests[0]->comm : MPI_COMM_NULL;
825   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
826
827   smpi_bench_end(rank, "Startall");
828   if (requests == NULL) {
829     retval = MPI_ERR_ARG;
830   } else {
831     smpi_mpi_startall(count, requests);
832     retval = MPI_SUCCESS;
833   }
834   smpi_bench_begin(rank, "Startall");
835   return retval;
836 }
837
838 int MPI_Request_free(MPI_Request * request)
839 {
840   int retval;
841   MPI_Comm comm = (*request)->comm;
842   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
843
844   smpi_bench_end(rank, "Request_free");
845   if (request == NULL) {
846     retval = MPI_ERR_ARG;
847   } else {
848     smpi_mpi_request_free(request);
849     retval = MPI_SUCCESS;
850   }
851   smpi_bench_begin(rank, "Request_free");
852   return retval;
853 }
854
855 int MPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
856               int tag, MPI_Comm comm, MPI_Request * request)
857 {
858   int retval;
859   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
860
861   smpi_bench_end(rank, "Irecv");
862 #ifdef HAVE_TRACING
863   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
864   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
865 #endif
866   if (request == NULL) {
867     retval = MPI_ERR_ARG;
868   } else if (comm == MPI_COMM_NULL) {
869     retval = MPI_ERR_COMM;
870   } else {
871     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
872     retval = MPI_SUCCESS;
873   }
874 #ifdef HAVE_TRACING
875   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
876   (*request)->recv = 1;
877 #endif
878   smpi_bench_begin(rank, "Irecv");
879   return retval;
880 }
881
882 int MPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
883               int tag, MPI_Comm comm, MPI_Request * request)
884 {
885   int retval;
886   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
887
888   smpi_bench_end(rank, "Isend");
889 #ifdef HAVE_TRACING
890   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
891   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
892   TRACE_smpi_send(rank, rank, dst_traced);
893 #endif
894   if (request == NULL) {
895     retval = MPI_ERR_ARG;
896   } else if (comm == MPI_COMM_NULL) {
897     retval = MPI_ERR_COMM;
898   } else {
899     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
900     retval = MPI_SUCCESS;
901   }
902 #ifdef HAVE_TRACING
903   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
904   (*request)->send = 1;
905 #endif
906   smpi_bench_begin(rank, "Isend");
907   return retval;
908 }
909
910 int MPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
911              MPI_Comm comm, MPI_Status * status)
912 {
913   int retval;
914   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
915
916   smpi_bench_end(rank, "Recv");
917 #ifdef HAVE_TRACING
918   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
919   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
920 #endif
921   if (comm == MPI_COMM_NULL) {
922     retval = MPI_ERR_COMM;
923   } else {
924     smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
925     retval = MPI_SUCCESS;
926   }
927 #ifdef HAVE_TRACING
928   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
929   TRACE_smpi_recv(rank, src_traced, rank);
930 #endif
931   smpi_bench_begin(rank, "Recv");
932   return retval;
933 }
934
935 int MPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
936              MPI_Comm comm)
937 {
938   int retval;
939   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
940
941   smpi_bench_end(rank, "Send");
942 #ifdef HAVE_TRACING
943   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
944   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
945   TRACE_smpi_send(rank, rank, dst_traced);
946 #endif
947   if (comm == MPI_COMM_NULL) {
948     retval = MPI_ERR_COMM;
949   } else {
950     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
951     retval = MPI_SUCCESS;
952   }
953 #ifdef HAVE_TRACING
954   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
955 #endif
956   smpi_bench_begin(rank, "Send");
957   return retval;
958 }
959
960 int MPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
961                  int dst, int sendtag, void *recvbuf, int recvcount,
962                  MPI_Datatype recvtype, int src, int recvtag,
963                  MPI_Comm comm, MPI_Status * status)
964 {
965   int retval;
966   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
967
968   smpi_bench_end(rank, "Sendrecv");
969 #ifdef HAVE_TRACING
970   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
971   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
972   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
973   TRACE_smpi_send(rank, rank, dst_traced);
974   TRACE_smpi_send(rank, src_traced, rank);
975 #endif
976   if (comm == MPI_COMM_NULL) {
977     retval = MPI_ERR_COMM;
978   } else if (sendtype == MPI_DATATYPE_NULL
979              || recvtype == MPI_DATATYPE_NULL) {
980     retval = MPI_ERR_TYPE;
981   } else {
982     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
983                       recvcount, recvtype, src, recvtag, comm, status);
984     retval = MPI_SUCCESS;
985   }
986 #ifdef HAVE_TRACING
987   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
988   TRACE_smpi_recv(rank, rank, dst_traced);
989   TRACE_smpi_recv(rank, src_traced, rank);
990 #endif
991   smpi_bench_begin(rank, "Sendrecv");
992   return retval;
993 }
994
995 int MPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
996                          int dst, int sendtag, int src, int recvtag,
997                          MPI_Comm comm, MPI_Status * status)
998 {
999   //TODO: suboptimal implementation
1000   void *recvbuf;
1001   int retval, size;
1002
1003   size = smpi_datatype_size(datatype) * count;
1004   recvbuf = xbt_new(char, size);
1005   retval =
1006       MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1007                    datatype, src, recvtag, comm, status);
1008   memcpy(buf, recvbuf, size * sizeof(char));
1009   xbt_free(recvbuf);
1010   return retval;
1011 }
1012
1013 int MPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1014 {
1015   int retval;
1016   int rank = request && (*request)->comm != MPI_COMM_NULL
1017       ? smpi_comm_rank((*request)->comm)
1018       : -1;
1019
1020   smpi_bench_end(rank, "Test");
1021   if (request == NULL || flag == NULL) {
1022     retval = MPI_ERR_ARG;
1023   } else if (*request == MPI_REQUEST_NULL) {
1024     retval = MPI_ERR_REQUEST;
1025   } else {
1026     *flag = smpi_mpi_test(request, status);
1027     retval = MPI_SUCCESS;
1028   }
1029   smpi_bench_begin(rank, "Test");
1030   return retval;
1031 }
1032
1033 int MPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1034                 MPI_Status * status)
1035 {
1036   int retval;
1037
1038   smpi_bench_end(-1, NULL);     //FIXME
1039   if (index == NULL || flag == NULL) {
1040     retval = MPI_ERR_ARG;
1041   } else {
1042     *flag = smpi_mpi_testany(count, requests, index, status);
1043     retval = MPI_SUCCESS;
1044   }
1045   smpi_bench_begin(-1, NULL);
1046   return retval;
1047 }
1048
1049 int MPI_Wait(MPI_Request * request, MPI_Status * status)
1050 {
1051   int retval;
1052   int rank = request && (*request)->comm != MPI_COMM_NULL
1053       ? smpi_comm_rank((*request)->comm)
1054       : -1;
1055
1056   smpi_bench_end(rank, "Wait");
1057 #ifdef HAVE_TRACING
1058   MPI_Group group = smpi_comm_group((*request)->comm);
1059   int src_traced = smpi_group_rank(group, (*request)->src);
1060   int dst_traced = smpi_group_rank(group, (*request)->dst);
1061   int is_wait_for_receive = (*request)->recv;
1062   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1063 #endif
1064   if (request == NULL) {
1065     retval = MPI_ERR_ARG;
1066   } else if (*request == MPI_REQUEST_NULL) {
1067     retval = MPI_ERR_REQUEST;
1068   } else {
1069     smpi_mpi_wait(request, status);
1070     retval = MPI_SUCCESS;
1071   }
1072 #ifdef HAVE_TRACING
1073   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1074   if (is_wait_for_receive) {
1075     TRACE_smpi_recv(rank, src_traced, dst_traced);
1076   }
1077 #endif
1078   smpi_bench_begin(rank, "Wait");
1079   return retval;
1080 }
1081
1082 int MPI_Waitany(int count, MPI_Request requests[], int *index,
1083                 MPI_Status * status)
1084 {
1085   int retval;
1086
1087   smpi_bench_end(-1, NULL);     //FIXME
1088 #ifdef HAVE_TRACING
1089   //save requests information for tracing
1090   int i;
1091   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
1092   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
1093   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
1094   for (i = 0; i < count; i++) {
1095     MPI_Request req = requests[i];      //already received requests are no longer valid
1096     if (req) {
1097       int *asrc = xbt_new(int, 1);
1098       int *adst = xbt_new(int, 1);
1099       int *arecv = xbt_new(int, 1);
1100       *asrc = req->src;
1101       *adst = req->dst;
1102       *arecv = req->recv;
1103       xbt_dynar_insert_at(srcs, i, asrc);
1104       xbt_dynar_insert_at(dsts, i, adst);
1105       xbt_dynar_insert_at(recvs, i, arecv);
1106     } else {
1107       int *t = xbt_new(int, 1);
1108       xbt_dynar_insert_at(srcs, i, t);
1109       xbt_dynar_insert_at(dsts, i, t);
1110       xbt_dynar_insert_at(recvs, i, t);
1111     }
1112   }
1113   int rank_traced = smpi_comm_rank(MPI_COMM_WORLD);
1114   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1115 #endif
1116   if (index == NULL) {
1117     retval = MPI_ERR_ARG;
1118   } else {
1119     *index = smpi_mpi_waitany(count, requests, status);
1120     retval = MPI_SUCCESS;
1121   }
1122 #ifdef HAVE_TRACING
1123   int src_traced, dst_traced, is_wait_for_receive;
1124   xbt_dynar_get_cpy(srcs, *index, &src_traced);
1125   xbt_dynar_get_cpy(dsts, *index, &dst_traced);
1126   xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
1127   if (is_wait_for_receive) {
1128     TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1129   }
1130   TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1131   //clean-up of dynars
1132   xbt_free(srcs);
1133   xbt_free(dsts);
1134   xbt_free(recvs);
1135 #endif
1136   smpi_bench_begin(-1, NULL);
1137   return retval;
1138 }
1139
1140 int MPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1141 {
1142
1143   smpi_bench_end(-1, NULL);     //FIXME
1144 #ifdef HAVE_TRACING
1145   //save information from requests
1146   int i;
1147   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
1148   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
1149   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
1150   for (i = 0; i < count; i++) {
1151     MPI_Request req = requests[i];      //all req should be valid in Waitall
1152     int *asrc = xbt_new(int, 1);
1153     int *adst = xbt_new(int, 1);
1154     int *arecv = xbt_new(int, 1);
1155     *asrc = req->src;
1156     *adst = req->dst;
1157     *arecv = req->recv;
1158     xbt_dynar_insert_at(srcs, i, asrc);
1159     xbt_dynar_insert_at(dsts, i, adst);
1160     xbt_dynar_insert_at(recvs, i, arecv);
1161   }
1162   int rank_traced = smpi_comm_rank (MPI_COMM_WORLD);
1163   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1164 #endif
1165   smpi_mpi_waitall(count, requests, status);
1166 #ifdef HAVE_TRACING
1167   for (i = 0; i < count; i++) {
1168     int src_traced, dst_traced, is_wait_for_receive;
1169     xbt_dynar_get_cpy(srcs, i, &src_traced);
1170     xbt_dynar_get_cpy(dsts, i, &dst_traced);
1171     xbt_dynar_get_cpy(recvs, i, &is_wait_for_receive);
1172     if (is_wait_for_receive) {
1173       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1174     }
1175   }
1176   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1177   //clean-up of dynars
1178   xbt_free(srcs);
1179   xbt_free(dsts);
1180   xbt_free(recvs);
1181 #endif
1182   smpi_bench_begin(-1, NULL);
1183   return MPI_SUCCESS;
1184 }
1185
1186 int MPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1187                  int *indices, MPI_Status status[])
1188 {
1189   int retval;
1190
1191   smpi_bench_end(-1, NULL);     //FIXME
1192   if (outcount == NULL || indices == NULL) {
1193     retval = MPI_ERR_ARG;
1194   } else {
1195     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1196     retval = MPI_SUCCESS;
1197   }
1198   smpi_bench_begin(-1, NULL);
1199   return retval;
1200 }
1201
1202 int MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
1203               MPI_Comm comm)
1204 {
1205   int retval;
1206   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1207
1208   smpi_bench_end(rank, "Bcast");
1209 #ifdef HAVE_TRACING
1210   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1211   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1212 #endif
1213   if (comm == MPI_COMM_NULL) {
1214     retval = MPI_ERR_COMM;
1215   } else {
1216     smpi_mpi_bcast(buf, count, datatype, root, comm);
1217     retval = MPI_SUCCESS;
1218   }
1219 #ifdef HAVE_TRACING
1220   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1221 #endif
1222   smpi_bench_begin(rank, "Bcast");
1223   return retval;
1224 }
1225
1226 int MPI_Barrier(MPI_Comm comm)
1227 {
1228   int retval;
1229   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1230
1231   smpi_bench_end(rank, "Barrier");
1232 #ifdef HAVE_TRACING
1233   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1234 #endif
1235   if (comm == MPI_COMM_NULL) {
1236     retval = MPI_ERR_COMM;
1237   } else {
1238     smpi_mpi_barrier(comm);
1239     retval = MPI_SUCCESS;
1240   }
1241 #ifdef HAVE_TRACING
1242   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1243 #endif
1244   smpi_bench_begin(rank, "Barrier");
1245   return retval;
1246 }
1247
1248 int MPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1249                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1250                int root, MPI_Comm comm)
1251 {
1252   int retval;
1253   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1254
1255   smpi_bench_end(rank, "Gather");
1256 #ifdef HAVE_TRACING
1257   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1258   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1259 #endif
1260   if (comm == MPI_COMM_NULL) {
1261     retval = MPI_ERR_COMM;
1262   } else if (sendtype == MPI_DATATYPE_NULL
1263              || recvtype == MPI_DATATYPE_NULL) {
1264     retval = MPI_ERR_TYPE;
1265   } else {
1266     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1267                     recvtype, root, comm);
1268     retval = MPI_SUCCESS;
1269   }
1270 #ifdef HAVE_TRACING
1271   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1272 #endif
1273   smpi_bench_begin(rank, "Gather");
1274   return retval;
1275 }
1276
1277 int MPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1278                 void *recvbuf, int *recvcounts, int *displs,
1279                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1280 {
1281   int retval;
1282   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1283
1284   smpi_bench_end(rank, "Gatherv");
1285 #ifdef HAVE_TRACING
1286   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1287   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1288 #endif
1289   if (comm == MPI_COMM_NULL) {
1290     retval = MPI_ERR_COMM;
1291   } else if (sendtype == MPI_DATATYPE_NULL
1292              || recvtype == MPI_DATATYPE_NULL) {
1293     retval = MPI_ERR_TYPE;
1294   } else if (recvcounts == NULL || displs == NULL) {
1295     retval = MPI_ERR_ARG;
1296   } else {
1297     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1298                      displs, recvtype, root, comm);
1299     retval = MPI_SUCCESS;
1300   }
1301 #ifdef HAVE_TRACING
1302   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1303 #endif
1304   smpi_bench_begin(rank, "Gatherv");
1305   return retval;
1306 }
1307
1308 int MPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1309                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1310                   MPI_Comm comm)
1311 {
1312   int retval;
1313   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1314
1315   smpi_bench_end(rank, "Allgather");
1316 #ifdef HAVE_TRACING
1317   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1318 #endif
1319   if (comm == MPI_COMM_NULL) {
1320     retval = MPI_ERR_COMM;
1321   } else if (sendtype == MPI_DATATYPE_NULL
1322              || recvtype == MPI_DATATYPE_NULL) {
1323     retval = MPI_ERR_TYPE;
1324   } else {
1325     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1326                        recvtype, comm);
1327     retval = MPI_SUCCESS;
1328   }
1329 #ifdef HAVE_TRACING
1330   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1331 #endif
1332   smpi_bench_begin(rank, "Allgather");
1333   return retval;
1334 }
1335
1336 int MPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1337                    void *recvbuf, int *recvcounts, int *displs,
1338                    MPI_Datatype recvtype, MPI_Comm comm)
1339 {
1340   int retval;
1341   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1342
1343   smpi_bench_end(rank, "Allgatherv");
1344 #ifdef HAVE_TRACING
1345   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1346 #endif
1347   if (comm == MPI_COMM_NULL) {
1348     retval = MPI_ERR_COMM;
1349   } else if (sendtype == MPI_DATATYPE_NULL
1350              || recvtype == MPI_DATATYPE_NULL) {
1351     retval = MPI_ERR_TYPE;
1352   } else if (recvcounts == NULL || displs == NULL) {
1353     retval = MPI_ERR_ARG;
1354   } else {
1355     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1356                         displs, recvtype, comm);
1357     retval = MPI_SUCCESS;
1358   }
1359 #ifdef HAVE_TRACING
1360   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1361 #endif
1362   smpi_bench_begin(rank, "Allgatherv");
1363   return retval;
1364 }
1365
1366 int MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1367                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1368                 int root, MPI_Comm comm)
1369 {
1370   int retval;
1371   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1372
1373   smpi_bench_end(rank, "Scatter");
1374 #ifdef HAVE_TRACING
1375   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1376   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1377 #endif
1378   if (comm == MPI_COMM_NULL) {
1379     retval = MPI_ERR_COMM;
1380   } else if (sendtype == MPI_DATATYPE_NULL
1381              || recvtype == MPI_DATATYPE_NULL) {
1382     retval = MPI_ERR_TYPE;
1383   } else {
1384     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1385                      recvtype, root, comm);
1386     retval = MPI_SUCCESS;
1387   }
1388 #ifdef HAVE_TRACING
1389   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1390 #endif
1391   smpi_bench_begin(rank, "Scatter");
1392   return retval;
1393 }
1394
1395 int MPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1396                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1397                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1398 {
1399   int retval;
1400   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1401
1402   smpi_bench_end(rank, "Scatterv");
1403 #ifdef HAVE_TRACING
1404   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1405   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1406 #endif
1407   if (comm == MPI_COMM_NULL) {
1408     retval = MPI_ERR_COMM;
1409   } else if (sendtype == MPI_DATATYPE_NULL
1410              || recvtype == MPI_DATATYPE_NULL) {
1411     retval = MPI_ERR_TYPE;
1412   } else if (sendcounts == NULL || displs == NULL) {
1413     retval = MPI_ERR_ARG;
1414   } else {
1415     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1416                       recvcount, recvtype, root, comm);
1417     retval = MPI_SUCCESS;
1418   }
1419 #ifdef HAVE_TRACING
1420   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1421 #endif
1422   smpi_bench_begin(rank, "Scatterv");
1423   return retval;
1424 }
1425
1426 int MPI_Reduce(void *sendbuf, void *recvbuf, int count,
1427                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1428 {
1429   int retval;
1430   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1431
1432   smpi_bench_end(rank, "Reduce");
1433 #ifdef HAVE_TRACING
1434   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1435   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1436 #endif
1437   if (comm == MPI_COMM_NULL) {
1438     retval = MPI_ERR_COMM;
1439   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1440     retval = MPI_ERR_ARG;
1441   } else {
1442     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1443     retval = MPI_SUCCESS;
1444   }
1445 #ifdef HAVE_TRACING
1446   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1447 #endif
1448   smpi_bench_begin(rank, "Reduce");
1449   return retval;
1450 }
1451
1452 int MPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1453                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1454 {
1455   int retval;
1456   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1457
1458   smpi_bench_end(rank, "Allreduce");
1459 #ifdef HAVE_TRACING
1460   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1461 #endif
1462   if (comm == MPI_COMM_NULL) {
1463     retval = MPI_ERR_COMM;
1464   } else if (datatype == MPI_DATATYPE_NULL) {
1465     retval = MPI_ERR_TYPE;
1466   } else if (op == MPI_OP_NULL) {
1467     retval = MPI_ERR_OP;
1468   } else {
1469     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
1470     retval = MPI_SUCCESS;
1471   }
1472 #ifdef HAVE_TRACING
1473   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1474 #endif
1475   smpi_bench_begin(rank, "Allreduce");
1476   return retval;
1477 }
1478
1479 int MPI_Scan(void *sendbuf, void *recvbuf, int count,
1480              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1481 {
1482   int retval;
1483   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1484
1485   smpi_bench_end(rank, "Scan");
1486 #ifdef HAVE_TRACING
1487   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1488 #endif
1489   if (comm == MPI_COMM_NULL) {
1490     retval = MPI_ERR_COMM;
1491   } else if (datatype == MPI_DATATYPE_NULL) {
1492     retval = MPI_ERR_TYPE;
1493   } else if (op == MPI_OP_NULL) {
1494     retval = MPI_ERR_OP;
1495   } else {
1496     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1497     retval = MPI_SUCCESS;
1498   }
1499 #ifdef HAVE_TRACING
1500   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1501 #endif
1502   smpi_bench_begin(rank, "Scan");
1503   return retval;
1504 }
1505
1506 int MPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1507                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1508 {
1509   int retval, i, size, count;
1510   int *displs;
1511   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1512
1513   smpi_bench_end(rank, "Reduce_scatter");
1514 #ifdef HAVE_TRACING
1515   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1516 #endif
1517   if (comm == MPI_COMM_NULL) {
1518     retval = MPI_ERR_COMM;
1519   } else if (datatype == MPI_DATATYPE_NULL) {
1520     retval = MPI_ERR_TYPE;
1521   } else if (op == MPI_OP_NULL) {
1522     retval = MPI_ERR_OP;
1523   } else if (recvcounts == NULL) {
1524     retval = MPI_ERR_ARG;
1525   } else {
1526     /* arbitrarily choose root as rank 0 */
1527     /* TODO: faster direct implementation ? */
1528     size = smpi_comm_size(comm);
1529     count = 0;
1530     displs = xbt_new(int, size);
1531     for (i = 0; i < size; i++) {
1532       count += recvcounts[i];
1533       displs[i] = 0;
1534     }
1535     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1536     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf,
1537                       recvcounts[rank], datatype, 0, comm);
1538     xbt_free(displs);
1539     retval = MPI_SUCCESS;
1540   }
1541 #ifdef HAVE_TRACING
1542   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1543 #endif
1544   smpi_bench_begin(rank, "Reduce_scatter");
1545   return retval;
1546 }
1547
1548 int MPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1549                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
1550                  MPI_Comm comm)
1551 {
1552   int retval, size, sendsize;
1553   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1554
1555   smpi_bench_end(rank, "Alltoall");
1556 #ifdef HAVE_TRACING
1557   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1558 #endif
1559   if (comm == MPI_COMM_NULL) {
1560     retval = MPI_ERR_COMM;
1561   } else if (sendtype == MPI_DATATYPE_NULL
1562              || recvtype == MPI_DATATYPE_NULL) {
1563     retval = MPI_ERR_TYPE;
1564   } else {
1565     size = smpi_comm_size(comm);
1566     sendsize = smpi_datatype_size(sendtype) * sendcount;
1567     if (sendsize < 200 && size > 12) {
1568       retval =
1569           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
1570                                          recvbuf, recvcount, recvtype,
1571                                          comm);
1572     } else if (sendsize < 3000) {
1573       retval =
1574           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
1575                                                 sendtype, recvbuf,
1576                                                 recvcount, recvtype, comm);
1577     } else {
1578       retval =
1579           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
1580                                             recvbuf, recvcount, recvtype,
1581                                             comm);
1582     }
1583   }
1584 #ifdef HAVE_TRACING
1585   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1586 #endif
1587   smpi_bench_begin(rank, "Alltoall");
1588   return retval;
1589 }
1590
1591 int MPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
1592                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
1593                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
1594 {
1595   int retval;
1596   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1597
1598   smpi_bench_end(rank, "Alltoallv");
1599 #ifdef HAVE_TRACING
1600   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1601 #endif
1602   if (comm == MPI_COMM_NULL) {
1603     retval = MPI_ERR_COMM;
1604   } else if (sendtype == MPI_DATATYPE_NULL
1605              || recvtype == MPI_DATATYPE_NULL) {
1606     retval = MPI_ERR_TYPE;
1607   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
1608              || recvdisps == NULL) {
1609     retval = MPI_ERR_ARG;
1610   } else {
1611     retval =
1612         smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
1613                                   recvbuf, recvcounts, recvdisps, recvtype,
1614                                   comm);
1615   }
1616 #ifdef HAVE_TRACING
1617   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1618 #endif
1619   smpi_bench_begin(rank, "Alltoallv");
1620   return retval;
1621 }
1622
1623
1624 int MPI_Get_processor_name(char *name, int *resultlen)
1625 {
1626   int retval = MPI_SUCCESS;
1627
1628   smpi_bench_end(-1, NULL);
1629   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
1630           MPI_MAX_PROCESSOR_NAME - 1);
1631   *resultlen =
1632       strlen(name) >
1633       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1634
1635   smpi_bench_begin(-1, NULL);
1636   return retval;
1637 }
1638
1639 int MPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
1640 {
1641   int retval = MPI_SUCCESS;
1642   size_t size;
1643
1644   smpi_bench_end(-1, NULL);
1645   if (status == NULL || count == NULL) {
1646     retval = MPI_ERR_ARG;
1647   } else if (datatype == MPI_DATATYPE_NULL) {
1648     retval = MPI_ERR_TYPE;
1649   } else {
1650     size = smpi_datatype_size(datatype);
1651     if (size == 0) {
1652       *count = 0;
1653     } else if (status->count % size != 0) {
1654       retval = MPI_UNDEFINED;
1655     } else {
1656       *count = smpi_mpi_get_count(status, datatype);
1657     }
1658   }
1659   smpi_bench_begin(-1, NULL);
1660   return retval;
1661 }