Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
581cd1dd3e984f339d0ee89d596c4b803d3504bc
[simgrid.git] / src / smpi / smpi_pmpi.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_pmpi, smpi,
12                                 "Logging specific to SMPI (pmpi)");
13
14 #ifdef HAVE_TRACING
15 //this function need to be here because of the calls to smpi_bench
16 void TRACE_smpi_set_category(const char *category)
17 {
18   //need to end bench otherwise categories for execution tasks are wrong
19   smpi_bench_end();
20   int ret;
21   if (!TRACE_is_enabled()){
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();
35 }
36 #endif
37
38 /* PMPI User level calls */
39
40 int PMPI_Init(int *argc, char ***argv)
41 {
42   smpi_process_init(argc, argv);
43 #ifdef HAVE_TRACING
44   TRACE_smpi_init(smpi_process_index());
45 #endif
46   smpi_bench_begin();
47   return MPI_SUCCESS;
48 }
49
50 int PMPI_Finalize(void)
51 {
52   smpi_bench_end();
53 #ifdef HAVE_TRACING
54   TRACE_smpi_finalize(smpi_process_index());
55 #endif
56   smpi_process_destroy();
57   return MPI_SUCCESS;
58 }
59
60 int PMPI_Init_thread(int *argc, char ***argv, int required, int *provided)
61 {
62   if (provided != NULL) {
63     *provided = MPI_THREAD_MULTIPLE;
64   }
65   return MPI_Init(argc, argv);
66 }
67
68 int PMPI_Query_thread(int *provided)
69 {
70   int retval;
71
72   smpi_bench_end();
73   if (provided == NULL) {
74     retval = MPI_ERR_ARG;
75   } else {
76     *provided = MPI_THREAD_MULTIPLE;
77     retval = MPI_SUCCESS;
78   }
79   smpi_bench_begin();
80   return retval;
81 }
82
83 int PMPI_Is_thread_main(int *flag)
84 {
85   int retval;
86
87   smpi_bench_end();
88   if (flag == NULL) {
89     retval = MPI_ERR_ARG;
90   } else {
91     *flag = smpi_process_index() == 0;
92     retval = MPI_SUCCESS;
93   }
94   smpi_bench_begin();
95   return retval;
96 }
97
98 int PMPI_Abort(MPI_Comm comm, int errorcode)
99 {
100   smpi_bench_end();
101   smpi_process_destroy();
102   // FIXME: should kill all processes in comm instead
103   SIMIX_req_process_kill(SIMIX_process_self());
104   return MPI_SUCCESS;
105 }
106
107 double PMPI_Wtime(void)
108 {
109   double time;
110
111   smpi_bench_end();
112   time = SIMIX_get_clock();
113   smpi_bench_begin();
114   return time;
115 }
116
117 int PMPI_Address(void *location, MPI_Aint * address)
118 {
119   int retval;
120
121   smpi_bench_end();
122   if (!address) {
123     retval = MPI_ERR_ARG;
124   } else {
125     *address = (MPI_Aint) location;
126   }
127   smpi_bench_begin();
128   return retval;
129 }
130
131 int PMPI_Type_free(MPI_Datatype * datatype)
132 {
133   int retval;
134
135   smpi_bench_end();
136   if (!datatype) {
137     retval = MPI_ERR_ARG;
138   } else {
139     // FIXME: always fail for now
140     retval = MPI_ERR_TYPE;
141   }
142   smpi_bench_begin();
143   return retval;
144 }
145
146 int PMPI_Type_size(MPI_Datatype datatype, int *size)
147 {
148   int retval;
149
150   smpi_bench_end();
151   if (datatype == MPI_DATATYPE_NULL) {
152     retval = MPI_ERR_TYPE;
153   } else if (size == NULL) {
154     retval = MPI_ERR_ARG;
155   } else {
156     *size = (int) smpi_datatype_size(datatype);
157     retval = MPI_SUCCESS;
158   }
159   smpi_bench_begin();
160   return retval;
161 }
162
163 int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
164 {
165   int retval;
166
167   smpi_bench_end();
168   if (datatype == MPI_DATATYPE_NULL) {
169     retval = MPI_ERR_TYPE;
170   } else if (lb == NULL || extent == NULL) {
171     retval = MPI_ERR_ARG;
172   } else {
173     retval = smpi_datatype_extent(datatype, lb, extent);
174   }
175   smpi_bench_begin();
176   return retval;
177 }
178
179 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
180 {
181   int retval;
182   MPI_Aint dummy;
183
184   smpi_bench_end();
185   if (datatype == MPI_DATATYPE_NULL) {
186     retval = MPI_ERR_TYPE;
187   } else if (extent == NULL) {
188     retval = MPI_ERR_ARG;
189   } else {
190     retval = smpi_datatype_extent(datatype, &dummy, extent);
191   }
192   smpi_bench_begin();
193   return retval;
194 }
195
196 int PMPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
197 {
198   int retval;
199
200   smpi_bench_end();
201   if (datatype == MPI_DATATYPE_NULL) {
202     retval = MPI_ERR_TYPE;
203   } else if (disp == NULL) {
204     retval = MPI_ERR_ARG;
205   } else {
206     *disp = smpi_datatype_lb(datatype);
207     retval = MPI_SUCCESS;
208   }
209   smpi_bench_begin();
210   return retval;
211 }
212
213 int PMPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
214 {
215   int retval;
216
217   smpi_bench_end();
218   if (datatype == MPI_DATATYPE_NULL) {
219     retval = MPI_ERR_TYPE;
220   } else if (disp == NULL) {
221     retval = MPI_ERR_ARG;
222   } else {
223     *disp = smpi_datatype_ub(datatype);
224     retval = MPI_SUCCESS;
225   }
226   smpi_bench_begin();
227   return retval;
228 }
229
230 int PMPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
231 {
232   int retval;
233
234   smpi_bench_end();
235   if (function == NULL || op == NULL) {
236     retval = MPI_ERR_ARG;
237   } else {
238     *op = smpi_op_new(function, commute);
239     retval = MPI_SUCCESS;
240   }
241   smpi_bench_begin();
242   return retval;
243 }
244
245 int PMPI_Op_free(MPI_Op * op)
246 {
247   int retval;
248
249   smpi_bench_end();
250   if (op == NULL) {
251     retval = MPI_ERR_ARG;
252   } else if (*op == MPI_OP_NULL) {
253     retval = MPI_ERR_OP;
254   } else {
255     smpi_op_destroy(*op);
256     *op = MPI_OP_NULL;
257     retval = MPI_SUCCESS;
258   }
259   smpi_bench_begin();
260   return retval;
261 }
262
263 int PMPI_Group_free(MPI_Group * group)
264 {
265   int retval;
266
267   smpi_bench_end();
268   if (group == NULL) {
269     retval = MPI_ERR_ARG;
270   } else {
271     smpi_group_destroy(*group);
272     *group = MPI_GROUP_NULL;
273     retval = MPI_SUCCESS;
274   }
275   smpi_bench_begin();
276   return retval;
277 }
278
279 int PMPI_Group_size(MPI_Group group, int *size)
280 {
281   int retval;
282
283   smpi_bench_end();
284   if (group == MPI_GROUP_NULL) {
285     retval = MPI_ERR_GROUP;
286   } else if (size == NULL) {
287     retval = MPI_ERR_ARG;
288   } else {
289     *size = smpi_group_size(group);
290     retval = MPI_SUCCESS;
291   }
292   smpi_bench_begin();
293   return retval;
294 }
295
296 int PMPI_Group_rank(MPI_Group group, int *rank)
297 {
298   int retval;
299
300   smpi_bench_end();
301   if (group == MPI_GROUP_NULL) {
302     retval = MPI_ERR_GROUP;
303   } else if (rank == NULL) {
304     retval = MPI_ERR_ARG;
305   } else {
306     *rank = smpi_group_rank(group, smpi_process_index());
307     retval = MPI_SUCCESS;
308   }
309   smpi_bench_begin();
310   return retval;
311 }
312
313 int PMPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1,
314                               MPI_Group group2, int *ranks2)
315 {
316   int retval, i, index;
317
318   smpi_bench_end();
319   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
320     retval = MPI_ERR_GROUP;
321   } else {
322     for (i = 0; i < n; i++) {
323       index = smpi_group_index(group1, ranks1[i]);
324       ranks2[i] = smpi_group_rank(group2, index);
325     }
326     retval = MPI_SUCCESS;
327   }
328   smpi_bench_begin();
329   return retval;
330 }
331
332 int PMPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
333 {
334   int retval;
335
336   smpi_bench_end();
337   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
338     retval = MPI_ERR_GROUP;
339   } else if (result == NULL) {
340     retval = MPI_ERR_ARG;
341   } else {
342     *result = smpi_group_compare(group1, group2);
343     retval = MPI_SUCCESS;
344   }
345   smpi_bench_begin();
346   return retval;
347 }
348
349 int PMPI_Group_union(MPI_Group group1, MPI_Group group2,
350                     MPI_Group * newgroup)
351 {
352   int retval, i, proc1, proc2, size, size2;
353
354   smpi_bench_end();
355   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
356     retval = MPI_ERR_GROUP;
357   } else if (newgroup == NULL) {
358     retval = MPI_ERR_ARG;
359   } else {
360     size = smpi_group_size(group1);
361     size2 = smpi_group_size(group2);
362     for (i = 0; i < size2; i++) {
363       proc2 = smpi_group_index(group2, i);
364       proc1 = smpi_group_rank(group1, proc2);
365       if (proc1 == MPI_UNDEFINED) {
366         size++;
367       }
368     }
369     if (size == 0) {
370       *newgroup = MPI_GROUP_EMPTY;
371     } else {
372       *newgroup = smpi_group_new(size);
373       size2 = smpi_group_size(group1);
374       for (i = 0; i < size2; i++) {
375         proc1 = smpi_group_index(group1, i);
376         smpi_group_set_mapping(*newgroup, proc1, i);
377       }
378       for (i = size2; i < size; i++) {
379         proc2 = smpi_group_index(group2, i - size2);
380         smpi_group_set_mapping(*newgroup, proc2, i);
381       }
382     }
383     smpi_group_use(*newgroup);
384     retval = MPI_SUCCESS;
385   }
386   smpi_bench_begin();
387   return retval;
388 }
389
390 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
391                            MPI_Group * newgroup)
392 {
393   int retval, i, proc1, proc2, size, size2;
394
395   smpi_bench_end();
396   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
397     retval = MPI_ERR_GROUP;
398   } else if (newgroup == NULL) {
399     retval = MPI_ERR_ARG;
400   } else {
401     size = smpi_group_size(group1);
402     size2 = smpi_group_size(group2);
403     for (i = 0; i < size2; i++) {
404       proc2 = smpi_group_index(group2, i);
405       proc1 = smpi_group_rank(group1, proc2);
406       if (proc1 == MPI_UNDEFINED) {
407         size--;
408       }
409     }
410     if (size == 0) {
411       *newgroup = MPI_GROUP_EMPTY;
412     } else {
413       *newgroup = smpi_group_new(size);
414       size2 = smpi_group_size(group1);
415       for (i = 0; i < size2; i++) {
416         proc1 = smpi_group_index(group1, i);
417         proc2 = smpi_group_rank(group2, proc1);
418         if (proc2 != MPI_UNDEFINED) {
419           smpi_group_set_mapping(*newgroup, proc1, i);
420         }
421       }
422     }
423     smpi_group_use(*newgroup);
424     retval = MPI_SUCCESS;
425   }
426   smpi_bench_begin();
427   return retval;
428 }
429
430 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
431 {
432   int retval, i, proc1, proc2, size, size2;
433
434   smpi_bench_end();
435   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
436     retval = MPI_ERR_GROUP;
437   } else if (newgroup == NULL) {
438     retval = MPI_ERR_ARG;
439   } else {
440     size = size2 = smpi_group_size(group1);
441     for (i = 0; i < size2; i++) {
442       proc1 = smpi_group_index(group1, i);
443       proc2 = smpi_group_rank(group2, proc1);
444       if (proc2 != MPI_UNDEFINED) {
445         size--;
446       }
447     }
448     if (size == 0) {
449       *newgroup = MPI_GROUP_EMPTY;
450     } else {
451       *newgroup = smpi_group_new(size);
452       for (i = 0; i < size2; i++) {
453         proc1 = smpi_group_index(group1, i);
454         proc2 = smpi_group_rank(group2, proc1);
455         if (proc2 == MPI_UNDEFINED) {
456           smpi_group_set_mapping(*newgroup, proc1, i);
457         }
458       }
459     }
460     smpi_group_use(*newgroup);
461     retval = MPI_SUCCESS;
462   }
463   smpi_bench_begin();
464   return retval;
465 }
466
467 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
468 {
469   int retval, i, index;
470
471   smpi_bench_end();
472   if (group == MPI_GROUP_NULL) {
473     retval = MPI_ERR_GROUP;
474   } else if (newgroup == NULL) {
475     retval = MPI_ERR_ARG;
476   } else {
477     if (n == 0) {
478       *newgroup = MPI_GROUP_EMPTY;
479     } else if (n == smpi_group_size(group)) {
480       *newgroup = group;
481     } else {
482       *newgroup = smpi_group_new(n);
483       for (i = 0; i < n; i++) {
484         index = smpi_group_index(group, ranks[i]);
485         smpi_group_set_mapping(*newgroup, index, i);
486       }
487     }
488     smpi_group_use(*newgroup);
489     retval = MPI_SUCCESS;
490   }
491   smpi_bench_begin();
492   return retval;
493 }
494
495 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
496 {
497   int retval, i, size, rank, index;
498
499   smpi_bench_end();
500   if (group == MPI_GROUP_NULL) {
501     retval = MPI_ERR_GROUP;
502   } else if (newgroup == NULL) {
503     retval = MPI_ERR_ARG;
504   } else {
505     if (n == 0) {
506       *newgroup = group;
507     } else if (n == smpi_group_size(group)) {
508       *newgroup = MPI_GROUP_EMPTY;
509     } else {
510       size = smpi_group_size(group) - n;
511       *newgroup = smpi_group_new(size);
512       rank = 0;
513       while (rank < size) {
514         for (i = 0; i < n; i++) {
515           if (ranks[i] == rank) {
516             break;
517           }
518         }
519         if (i >= n) {
520           index = smpi_group_index(group, rank);
521           smpi_group_set_mapping(*newgroup, index, rank);
522           rank++;
523         }
524       }
525     }
526     smpi_group_use(*newgroup);
527     retval = MPI_SUCCESS;
528   }
529   smpi_bench_begin();
530   return retval;
531 }
532
533 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
534                          MPI_Group * newgroup)
535 {
536   int retval, i, j, rank, size, index;
537
538   smpi_bench_end();
539   if (group == MPI_GROUP_NULL) {
540     retval = MPI_ERR_GROUP;
541   } else if (newgroup == NULL) {
542     retval = MPI_ERR_ARG;
543   } else {
544     if (n == 0) {
545       *newgroup = MPI_GROUP_EMPTY;
546     } else {
547       size = 0;
548       for (i = 0; i < n; i++) {
549         for (rank = ranges[i][0];       /* First */
550              rank >= 0 && rank <= ranges[i][1]; /* Last */
551              rank += ranges[i][2] /* Stride */ ) {
552           size++;
553         }
554       }
555       if (size == smpi_group_size(group)) {
556         *newgroup = group;
557       } else {
558         *newgroup = smpi_group_new(size);
559         j = 0;
560         for (i = 0; i < n; i++) {
561           for (rank = ranges[i][0];     /* First */
562                rank >= 0 && rank <= ranges[i][1];       /* Last */
563                rank += ranges[i][2] /* Stride */ ) {
564             index = smpi_group_index(group, rank);
565             smpi_group_set_mapping(*newgroup, index, j);
566             j++;
567           }
568         }
569       }
570     }
571     smpi_group_use(*newgroup);
572     retval = MPI_SUCCESS;
573   }
574   smpi_bench_begin();
575   return retval;
576 }
577
578 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
579                          MPI_Group * newgroup)
580 {
581   int retval, i, newrank, rank, size, index, add;
582
583   smpi_bench_end();
584   if (group == MPI_GROUP_NULL) {
585     retval = MPI_ERR_GROUP;
586   } else if (newgroup == NULL) {
587     retval = MPI_ERR_ARG;
588   } else {
589     if (n == 0) {
590       *newgroup = group;
591     } else {
592       size = smpi_group_size(group);
593       for (i = 0; i < n; i++) {
594         for (rank = ranges[i][0];       /* First */
595              rank >= 0 && rank <= ranges[i][1]; /* Last */
596              rank += ranges[i][2] /* Stride */ ) {
597           size--;
598         }
599       }
600       if (size == 0) {
601         *newgroup = MPI_GROUP_EMPTY;
602       } else {
603         *newgroup = smpi_group_new(size);
604         newrank = 0;
605         while (newrank < size) {
606           for (i = 0; i < n; i++) {
607             add = 1;
608             for (rank = ranges[i][0];   /* First */
609                  rank >= 0 && rank <= ranges[i][1];     /* Last */
610                  rank += ranges[i][2] /* Stride */ ) {
611               if (rank == newrank) {
612                 add = 0;
613                 break;
614               }
615             }
616             if (add == 1) {
617               index = smpi_group_index(group, newrank);
618               smpi_group_set_mapping(*newgroup, index, newrank);
619             }
620           }
621         }
622       }
623     }
624     smpi_group_use(*newgroup);
625     retval = MPI_SUCCESS;
626   }
627   smpi_bench_begin();
628   return retval;
629 }
630
631 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
632 {
633   int retval;
634
635   smpi_bench_end();
636   if (comm == MPI_COMM_NULL) {
637     retval = MPI_ERR_COMM;
638   } else {
639     *rank = smpi_comm_rank(comm);
640     retval = MPI_SUCCESS;
641   }
642   smpi_bench_begin();
643   return retval;
644 }
645
646 int PMPI_Comm_size(MPI_Comm comm, int *size)
647 {
648   int retval;
649
650   smpi_bench_end();
651   if (comm == MPI_COMM_NULL) {
652     retval = MPI_ERR_COMM;
653   } else if (size == NULL) {
654     retval = MPI_ERR_ARG;
655   } else {
656     *size = smpi_comm_size(comm);
657     retval = MPI_SUCCESS;
658   }
659   smpi_bench_begin();
660   return retval;
661 }
662
663 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
664 {
665   int retval;
666
667   smpi_bench_end();
668   if (comm == MPI_COMM_NULL) {
669     retval = MPI_ERR_COMM;
670   } else if (group == NULL) {
671     retval = MPI_ERR_ARG;
672   } else {
673     *group = smpi_comm_group(comm);
674     retval = MPI_SUCCESS;
675   }
676   smpi_bench_begin();
677   return retval;
678 }
679
680 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
681 {
682   int retval;
683
684   smpi_bench_end();
685   if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
686     retval = MPI_ERR_COMM;
687   } else if (result == NULL) {
688     retval = MPI_ERR_ARG;
689   } else {
690     if (comm1 == comm2) {       /* Same communicators means same groups */
691       *result = MPI_IDENT;
692     } else {
693       *result =
694           smpi_group_compare(smpi_comm_group(comm1),
695                              smpi_comm_group(comm2));
696       if (*result == MPI_IDENT) {
697         *result = MPI_CONGRUENT;
698       }
699     }
700     retval = MPI_SUCCESS;
701   }
702   smpi_bench_begin();
703   return retval;
704 }
705
706 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
707 {
708   int retval;
709
710   smpi_bench_end();
711   if (comm == MPI_COMM_NULL) {
712     retval = MPI_ERR_COMM;
713   } else if (newcomm == NULL) {
714     retval = MPI_ERR_ARG;
715   } else {
716     *newcomm = smpi_comm_new(smpi_comm_group(comm));
717     retval = MPI_SUCCESS;
718   }
719   smpi_bench_begin();
720   return retval;
721 }
722
723 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
724 {
725   int retval;
726
727   smpi_bench_end();
728   if (comm == MPI_COMM_NULL) {
729     retval = MPI_ERR_COMM;
730   } else if (group == MPI_GROUP_NULL) {
731     retval = MPI_ERR_GROUP;
732   } else if (newcomm == NULL) {
733     retval = MPI_ERR_ARG;
734   } else {
735     *newcomm = smpi_comm_new(group);
736     retval = MPI_SUCCESS;
737   }
738   smpi_bench_begin();
739   return retval;
740 }
741
742 int PMPI_Comm_free(MPI_Comm * comm)
743 {
744   int retval;
745
746   smpi_bench_end();
747   if (comm == NULL) {
748     retval = MPI_ERR_ARG;
749   } else if (*comm == MPI_COMM_NULL) {
750     retval = MPI_ERR_COMM;
751   } else {
752     smpi_comm_destroy(*comm);
753     *comm = MPI_COMM_NULL;
754     retval = MPI_SUCCESS;
755   }
756   smpi_bench_begin();
757   return retval;
758 }
759
760 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
761 {
762   int retval;
763
764   smpi_bench_end();
765   if (comm_out == NULL) {
766     retval = MPI_ERR_ARG;
767   } else if (comm == MPI_COMM_NULL) {
768     retval = MPI_ERR_COMM;
769   } else {
770     *comm_out = smpi_comm_split(comm, color, key);
771     retval = MPI_SUCCESS;
772   }
773   smpi_bench_begin();
774   return retval;
775 }
776
777 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
778                   int tag, MPI_Comm comm, MPI_Request * request)
779 {
780   int retval;
781
782   smpi_bench_end();
783   if (request == NULL) {
784     retval = MPI_ERR_ARG;
785   } else if (comm == MPI_COMM_NULL) {
786     retval = MPI_ERR_COMM;
787   } else {
788     *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
789     retval = MPI_SUCCESS;
790   }
791   smpi_bench_begin();
792   return retval;
793 }
794
795 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
796                   int tag, MPI_Comm comm, MPI_Request * request)
797 {
798   int retval;
799
800   smpi_bench_end();
801   if (request == NULL) {
802     retval = MPI_ERR_ARG;
803   } else if (comm == MPI_COMM_NULL) {
804     retval = MPI_ERR_COMM;
805   } else {
806     *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
807     retval = MPI_SUCCESS;
808   }
809   smpi_bench_begin();
810   return retval;
811 }
812
813 int PMPI_Start(MPI_Request * request)
814 {
815   int retval;
816
817   smpi_bench_end();
818   if (request == NULL) {
819     retval = MPI_ERR_ARG;
820   } else {
821     smpi_mpi_start(*request);
822     retval = MPI_SUCCESS;
823   }
824   smpi_bench_begin();
825   return retval;
826 }
827
828 int PMPI_Startall(int count, MPI_Request * requests)
829 {
830   int retval;
831
832   smpi_bench_end();
833   if (requests == NULL) {
834     retval = MPI_ERR_ARG;
835   } else {
836     smpi_mpi_startall(count, requests);
837     retval = MPI_SUCCESS;
838   }
839   smpi_bench_begin();
840   return retval;
841 }
842
843 int PMPI_Request_free(MPI_Request * request)
844 {
845   int retval;
846
847   smpi_bench_end();
848   if (request == NULL) {
849     retval = MPI_ERR_ARG;
850   } else {
851     smpi_mpi_request_free(request);
852     retval = MPI_SUCCESS;
853   }
854   smpi_bench_begin();
855   return retval;
856 }
857
858 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
859               int tag, MPI_Comm comm, MPI_Request * request)
860 {
861   int retval;
862
863   smpi_bench_end();
864 #ifdef HAVE_TRACING
865   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
866   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
867   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
868 #endif
869   if (request == NULL) {
870     retval = MPI_ERR_ARG;
871   } else if (comm == MPI_COMM_NULL) {
872     retval = MPI_ERR_COMM;
873   } else {
874     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
875     retval = MPI_SUCCESS;
876   }
877 #ifdef HAVE_TRACING
878   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
879   (*request)->recv = 1;
880 #endif
881   smpi_bench_begin();
882   return retval;
883 }
884
885 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
886               int tag, MPI_Comm comm, MPI_Request * request)
887 {
888   int retval;
889
890   smpi_bench_end();
891 #ifdef HAVE_TRACING
892   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
893   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
894   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
895   TRACE_smpi_send(rank, rank, dst_traced);
896 #endif
897   if (request == NULL) {
898     retval = MPI_ERR_ARG;
899   } else if (comm == MPI_COMM_NULL) {
900     retval = MPI_ERR_COMM;
901   } else {
902     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
903     retval = MPI_SUCCESS;
904   }
905 #ifdef HAVE_TRACING
906   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
907   (*request)->send = 1;
908 #endif
909   smpi_bench_begin();
910   return retval;
911 }
912
913 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
914              MPI_Comm comm, MPI_Status * status)
915 {
916   int retval;
917
918   smpi_bench_end();
919 #ifdef HAVE_TRACING
920   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
921   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
922   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
923 #endif
924   if (comm == MPI_COMM_NULL) {
925     retval = MPI_ERR_COMM;
926   } else {
927     smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
928     retval = MPI_SUCCESS;
929   }
930 #ifdef HAVE_TRACING
931   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
932   TRACE_smpi_recv(rank, src_traced, rank);
933 #endif
934   smpi_bench_begin();
935   return retval;
936 }
937
938 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
939              MPI_Comm comm)
940 {
941   int retval;
942
943   smpi_bench_end();
944 #ifdef HAVE_TRACING
945   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
946   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
947   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
948   TRACE_smpi_send(rank, rank, dst_traced);
949 #endif
950   if (comm == MPI_COMM_NULL) {
951     retval = MPI_ERR_COMM;
952   } else {
953     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
954     retval = MPI_SUCCESS;
955   }
956 #ifdef HAVE_TRACING
957   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
958 #endif
959   smpi_bench_begin();
960   return retval;
961 }
962
963 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
964                  int dst, int sendtag, void *recvbuf, int recvcount,
965                  MPI_Datatype recvtype, int src, int recvtag,
966                  MPI_Comm comm, MPI_Status * status)
967 {
968   int retval;
969
970   smpi_bench_end();
971 #ifdef HAVE_TRACING
972   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
973   int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
974   int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
975   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
976   TRACE_smpi_send(rank, rank, dst_traced);
977   TRACE_smpi_send(rank, src_traced, rank);
978 #endif
979   if (comm == MPI_COMM_NULL) {
980     retval = MPI_ERR_COMM;
981   } else if (sendtype == MPI_DATATYPE_NULL
982              || recvtype == MPI_DATATYPE_NULL) {
983     retval = MPI_ERR_TYPE;
984   } else {
985     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
986                       recvcount, recvtype, src, recvtag, comm, status);
987     retval = MPI_SUCCESS;
988   }
989 #ifdef HAVE_TRACING
990   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
991   TRACE_smpi_recv(rank, rank, dst_traced);
992   TRACE_smpi_recv(rank, src_traced, rank);
993 #endif
994   smpi_bench_begin();
995   return retval;
996 }
997
998 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
999                          int dst, int sendtag, int src, int recvtag,
1000                          MPI_Comm comm, MPI_Status * status)
1001 {
1002   //TODO: suboptimal implementation
1003   void *recvbuf;
1004   int retval, size;
1005
1006   size = smpi_datatype_size(datatype) * count;
1007   recvbuf = xbt_new(char, size);
1008   retval =
1009       MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1010                    datatype, src, recvtag, comm, status);
1011   memcpy(buf, recvbuf, size * sizeof(char));
1012   xbt_free(recvbuf);
1013   return retval;
1014 }
1015
1016 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1017 {
1018   int retval;
1019
1020   smpi_bench_end();
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();
1030   return retval;
1031 }
1032
1033 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1034                 MPI_Status * status)
1035 {
1036   int retval;
1037
1038   smpi_bench_end();
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();
1046   return retval;
1047 }
1048
1049 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1050 {
1051   int retval;
1052
1053   smpi_bench_end();
1054 #ifdef HAVE_TRACING
1055   int rank = request && (*request)->comm != MPI_COMM_NULL
1056       ? smpi_comm_rank((*request)->comm)
1057       : -1;
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();
1079   return retval;
1080 }
1081
1082 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1083 {
1084   int retval;
1085
1086   smpi_bench_end();
1087 #ifdef HAVE_TRACING
1088   //save requests information for tracing
1089   int i;
1090   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
1091   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
1092   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
1093   for (i = 0; i < count; i++) {
1094     MPI_Request req = requests[i];      //already received requests are no longer valid
1095     if (req) {
1096       int *asrc = xbt_new(int, 1);
1097       int *adst = xbt_new(int, 1);
1098       int *arecv = xbt_new(int, 1);
1099       *asrc = req->src;
1100       *adst = req->dst;
1101       *arecv = req->recv;
1102       xbt_dynar_insert_at(srcs, i, asrc);
1103       xbt_dynar_insert_at(dsts, i, adst);
1104       xbt_dynar_insert_at(recvs, i, arecv);
1105     } else {
1106       int *t = xbt_new(int, 1);
1107       xbt_dynar_insert_at(srcs, i, t);
1108       xbt_dynar_insert_at(dsts, i, t);
1109       xbt_dynar_insert_at(recvs, i, t);
1110     }
1111   }
1112   int rank_traced = smpi_comm_rank(MPI_COMM_WORLD);
1113   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1114 #endif
1115   if (index == NULL) {
1116     retval = MPI_ERR_ARG;
1117   } else {
1118     *index = smpi_mpi_waitany(count, requests, status);
1119     retval = MPI_SUCCESS;
1120   }
1121 #ifdef HAVE_TRACING
1122   int src_traced, dst_traced, is_wait_for_receive;
1123   xbt_dynar_get_cpy(srcs, *index, &src_traced);
1124   xbt_dynar_get_cpy(dsts, *index, &dst_traced);
1125   xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
1126   if (is_wait_for_receive) {
1127     TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1128   }
1129   TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1130   //clean-up of dynars
1131   xbt_free(srcs);
1132   xbt_free(dsts);
1133   xbt_free(recvs);
1134 #endif
1135   smpi_bench_begin();
1136   return retval;
1137 }
1138
1139 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1140 {
1141
1142   smpi_bench_end();
1143 #ifdef HAVE_TRACING
1144   //save information from requests
1145   int i;
1146   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), xbt_free);
1147   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), xbt_free);
1148   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), xbt_free);
1149   for (i = 0; i < count; i++) {
1150     MPI_Request req = requests[i];      //all req should be valid in Waitall
1151     int *asrc = xbt_new(int, 1);
1152     int *adst = xbt_new(int, 1);
1153     int *arecv = xbt_new(int, 1);
1154     *asrc = req->src;
1155     *adst = req->dst;
1156     *arecv = req->recv;
1157     xbt_dynar_insert_at(srcs, i, asrc);
1158     xbt_dynar_insert_at(dsts, i, adst);
1159     xbt_dynar_insert_at(recvs, i, arecv);
1160   }
1161   int rank_traced = smpi_comm_rank (MPI_COMM_WORLD);
1162   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1163 #endif
1164   smpi_mpi_waitall(count, requests, status);
1165 #ifdef HAVE_TRACING
1166   for (i = 0; i < count; i++) {
1167     int src_traced, dst_traced, is_wait_for_receive;
1168     xbt_dynar_get_cpy(srcs, i, &src_traced);
1169     xbt_dynar_get_cpy(dsts, i, &dst_traced);
1170     xbt_dynar_get_cpy(recvs, i, &is_wait_for_receive);
1171     if (is_wait_for_receive) {
1172       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1173     }
1174   }
1175   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1176   //clean-up of dynars
1177   xbt_free(srcs);
1178   xbt_free(dsts);
1179   xbt_free(recvs);
1180 #endif
1181   smpi_bench_begin();
1182   return MPI_SUCCESS;
1183 }
1184
1185 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1186                  int *indices, MPI_Status status[])
1187 {
1188   int retval;
1189
1190   smpi_bench_end();
1191   if (outcount == NULL || indices == NULL) {
1192     retval = MPI_ERR_ARG;
1193   } else {
1194     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1195     retval = MPI_SUCCESS;
1196   }
1197   smpi_bench_begin();
1198   return retval;
1199 }
1200
1201 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1202 {
1203   int retval;
1204
1205   smpi_bench_end();
1206 #ifdef HAVE_TRACING
1207   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1208   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1209   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1210 #endif
1211   if (comm == MPI_COMM_NULL) {
1212     retval = MPI_ERR_COMM;
1213   } else {
1214     smpi_mpi_bcast(buf, count, datatype, root, comm);
1215     retval = MPI_SUCCESS;
1216   }
1217 #ifdef HAVE_TRACING
1218   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1219 #endif
1220   smpi_bench_begin();
1221   return retval;
1222 }
1223
1224 int PMPI_Barrier(MPI_Comm comm)
1225 {
1226   int retval;
1227
1228   smpi_bench_end();
1229 #ifdef HAVE_TRACING
1230   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1231   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1232 #endif
1233   if (comm == MPI_COMM_NULL) {
1234     retval = MPI_ERR_COMM;
1235   } else {
1236     smpi_mpi_barrier(comm);
1237     retval = MPI_SUCCESS;
1238   }
1239 #ifdef HAVE_TRACING
1240   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1241 #endif
1242   smpi_bench_begin();
1243   return retval;
1244 }
1245
1246 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1247                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1248                int root, MPI_Comm comm)
1249 {
1250   int retval;
1251
1252   smpi_bench_end();
1253 #ifdef HAVE_TRACING
1254   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1255   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1256   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1257 #endif
1258   if (comm == MPI_COMM_NULL) {
1259     retval = MPI_ERR_COMM;
1260   } else if (sendtype == MPI_DATATYPE_NULL
1261              || recvtype == MPI_DATATYPE_NULL) {
1262     retval = MPI_ERR_TYPE;
1263   } else {
1264     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1265                     recvtype, root, comm);
1266     retval = MPI_SUCCESS;
1267   }
1268 #ifdef HAVE_TRACING
1269   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1270 #endif
1271   smpi_bench_begin();
1272   return retval;
1273 }
1274
1275 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1276                 void *recvbuf, int *recvcounts, int *displs,
1277                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1278 {
1279   int retval;
1280
1281   smpi_bench_end();
1282 #ifdef HAVE_TRACING
1283   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1284   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1285   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1286 #endif
1287   if (comm == MPI_COMM_NULL) {
1288     retval = MPI_ERR_COMM;
1289   } else if (sendtype == MPI_DATATYPE_NULL
1290              || recvtype == MPI_DATATYPE_NULL) {
1291     retval = MPI_ERR_TYPE;
1292   } else if (recvcounts == NULL || displs == NULL) {
1293     retval = MPI_ERR_ARG;
1294   } else {
1295     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1296                      displs, recvtype, root, comm);
1297     retval = MPI_SUCCESS;
1298   }
1299 #ifdef HAVE_TRACING
1300   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1301 #endif
1302   smpi_bench_begin();
1303   return retval;
1304 }
1305
1306 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1307                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1308                   MPI_Comm comm)
1309 {
1310   int retval;
1311
1312   smpi_bench_end();
1313 #ifdef HAVE_TRACING
1314   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1315   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1316 #endif
1317   if (comm == MPI_COMM_NULL) {
1318     retval = MPI_ERR_COMM;
1319   } else if (sendtype == MPI_DATATYPE_NULL
1320              || recvtype == MPI_DATATYPE_NULL) {
1321     retval = MPI_ERR_TYPE;
1322   } else {
1323     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1324                        recvtype, comm);
1325     retval = MPI_SUCCESS;
1326   }
1327 #ifdef HAVE_TRACING
1328   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1329 #endif
1330   smpi_bench_begin();
1331   return retval;
1332 }
1333
1334 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1335                    void *recvbuf, int *recvcounts, int *displs,
1336                    MPI_Datatype recvtype, MPI_Comm comm)
1337 {
1338   int retval;
1339
1340   smpi_bench_end();
1341 #ifdef HAVE_TRACING
1342   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1343   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1344 #endif
1345   if (comm == MPI_COMM_NULL) {
1346     retval = MPI_ERR_COMM;
1347   } else if (sendtype == MPI_DATATYPE_NULL
1348              || recvtype == MPI_DATATYPE_NULL) {
1349     retval = MPI_ERR_TYPE;
1350   } else if (recvcounts == NULL || displs == NULL) {
1351     retval = MPI_ERR_ARG;
1352   } else {
1353     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1354                         displs, recvtype, comm);
1355     retval = MPI_SUCCESS;
1356   }
1357 #ifdef HAVE_TRACING
1358   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1359 #endif
1360   smpi_bench_begin();
1361   return retval;
1362 }
1363
1364 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1365                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1366                 int root, MPI_Comm comm)
1367 {
1368   int retval;
1369
1370   smpi_bench_end();
1371 #ifdef HAVE_TRACING
1372   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1373   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1374   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1375 #endif
1376   if (comm == MPI_COMM_NULL) {
1377     retval = MPI_ERR_COMM;
1378   } else if (sendtype == MPI_DATATYPE_NULL
1379              || recvtype == MPI_DATATYPE_NULL) {
1380     retval = MPI_ERR_TYPE;
1381   } else {
1382     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1383                      recvtype, root, comm);
1384     retval = MPI_SUCCESS;
1385   }
1386 #ifdef HAVE_TRACING
1387   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1388 #endif
1389   smpi_bench_begin();
1390   return retval;
1391 }
1392
1393 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1394                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1395                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1396 {
1397   int retval;
1398
1399   smpi_bench_end();
1400 #ifdef HAVE_TRACING
1401   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1402   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1403   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1404 #endif
1405   if (comm == MPI_COMM_NULL) {
1406     retval = MPI_ERR_COMM;
1407   } else if (sendtype == MPI_DATATYPE_NULL
1408              || recvtype == MPI_DATATYPE_NULL) {
1409     retval = MPI_ERR_TYPE;
1410   } else if (sendcounts == NULL || displs == NULL) {
1411     retval = MPI_ERR_ARG;
1412   } else {
1413     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1414                       recvcount, recvtype, root, comm);
1415     retval = MPI_SUCCESS;
1416   }
1417 #ifdef HAVE_TRACING
1418   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1419 #endif
1420   smpi_bench_begin();
1421   return retval;
1422 }
1423
1424 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1425                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1426 {
1427   int retval;
1428
1429   smpi_bench_end();
1430 #ifdef HAVE_TRACING
1431   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1432   int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1433   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1434 #endif
1435   if (comm == MPI_COMM_NULL) {
1436     retval = MPI_ERR_COMM;
1437   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1438     retval = MPI_ERR_ARG;
1439   } else {
1440     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1441     retval = MPI_SUCCESS;
1442   }
1443 #ifdef HAVE_TRACING
1444   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1445 #endif
1446   smpi_bench_begin();
1447   return retval;
1448 }
1449
1450 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1451                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1452 {
1453   int retval;
1454
1455   smpi_bench_end();
1456 #ifdef HAVE_TRACING
1457   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1458   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1459 #endif
1460   if (comm == MPI_COMM_NULL) {
1461     retval = MPI_ERR_COMM;
1462   } else if (datatype == MPI_DATATYPE_NULL) {
1463     retval = MPI_ERR_TYPE;
1464   } else if (op == MPI_OP_NULL) {
1465     retval = MPI_ERR_OP;
1466   } else {
1467     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
1468     retval = MPI_SUCCESS;
1469   }
1470 #ifdef HAVE_TRACING
1471   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1472 #endif
1473   smpi_bench_begin();
1474   return retval;
1475 }
1476
1477 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
1478              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1479 {
1480   int retval;
1481
1482   smpi_bench_end();
1483 #ifdef HAVE_TRACING
1484   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1485   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1486 #endif
1487   if (comm == MPI_COMM_NULL) {
1488     retval = MPI_ERR_COMM;
1489   } else if (datatype == MPI_DATATYPE_NULL) {
1490     retval = MPI_ERR_TYPE;
1491   } else if (op == MPI_OP_NULL) {
1492     retval = MPI_ERR_OP;
1493   } else {
1494     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1495     retval = MPI_SUCCESS;
1496   }
1497 #ifdef HAVE_TRACING
1498   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1499 #endif
1500   smpi_bench_begin();
1501   return retval;
1502 }
1503
1504 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1505                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1506 {
1507   int retval, i, size, count;
1508   int *displs;
1509   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1510
1511   smpi_bench_end();
1512 #ifdef HAVE_TRACING
1513   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1514 #endif
1515   if (comm == MPI_COMM_NULL) {
1516     retval = MPI_ERR_COMM;
1517   } else if (datatype == MPI_DATATYPE_NULL) {
1518     retval = MPI_ERR_TYPE;
1519   } else if (op == MPI_OP_NULL) {
1520     retval = MPI_ERR_OP;
1521   } else if (recvcounts == NULL) {
1522     retval = MPI_ERR_ARG;
1523   } else {
1524     /* arbitrarily choose root as rank 0 */
1525     /* TODO: faster direct implementation ? */
1526     size = smpi_comm_size(comm);
1527     count = 0;
1528     displs = xbt_new(int, size);
1529     for (i = 0; i < size; i++) {
1530       count += recvcounts[i];
1531       displs[i] = 0;
1532     }
1533     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1534     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf,
1535                       recvcounts[rank], datatype, 0, comm);
1536     xbt_free(displs);
1537     retval = MPI_SUCCESS;
1538   }
1539 #ifdef HAVE_TRACING
1540   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1541 #endif
1542   smpi_bench_begin();
1543   return retval;
1544 }
1545
1546 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1547                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
1548                  MPI_Comm comm)
1549 {
1550   int retval, size, sendsize;
1551
1552   smpi_bench_end();
1553 #ifdef HAVE_TRACING
1554   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1555   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1556 #endif
1557   if (comm == MPI_COMM_NULL) {
1558     retval = MPI_ERR_COMM;
1559   } else if (sendtype == MPI_DATATYPE_NULL
1560              || recvtype == MPI_DATATYPE_NULL) {
1561     retval = MPI_ERR_TYPE;
1562   } else {
1563     size = smpi_comm_size(comm);
1564     sendsize = smpi_datatype_size(sendtype) * sendcount;
1565     if (sendsize < 200 && size > 12) {
1566       retval =
1567           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
1568                                          recvbuf, recvcount, recvtype,
1569                                          comm);
1570     } else if (sendsize < 3000) {
1571       retval =
1572           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
1573                                                 sendtype, recvbuf,
1574                                                 recvcount, recvtype, comm);
1575     } else {
1576       retval =
1577           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
1578                                             recvbuf, recvcount, recvtype,
1579                                             comm);
1580     }
1581   }
1582 #ifdef HAVE_TRACING
1583   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1584 #endif
1585   smpi_bench_begin();
1586   return retval;
1587 }
1588
1589 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
1590                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
1591                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
1592 {
1593   int retval;
1594
1595   smpi_bench_end();
1596 #ifdef HAVE_TRACING
1597   int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1598   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1599 #endif
1600   if (comm == MPI_COMM_NULL) {
1601     retval = MPI_ERR_COMM;
1602   } else if (sendtype == MPI_DATATYPE_NULL
1603              || recvtype == MPI_DATATYPE_NULL) {
1604     retval = MPI_ERR_TYPE;
1605   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
1606              || recvdisps == NULL) {
1607     retval = MPI_ERR_ARG;
1608   } else {
1609     retval =
1610         smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
1611                                   recvbuf, recvcounts, recvdisps, recvtype,
1612                                   comm);
1613   }
1614 #ifdef HAVE_TRACING
1615   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1616 #endif
1617   smpi_bench_begin();
1618   return retval;
1619 }
1620
1621
1622 int PMPI_Get_processor_name(char *name, int *resultlen)
1623 {
1624   int retval = MPI_SUCCESS;
1625
1626   smpi_bench_end();
1627   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
1628           MPI_MAX_PROCESSOR_NAME - 1);
1629   *resultlen =
1630       strlen(name) >
1631       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1632
1633   smpi_bench_begin();
1634   return retval;
1635 }
1636
1637 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
1638 {
1639   int retval = MPI_SUCCESS;
1640   size_t size;
1641
1642   smpi_bench_end();
1643   if (status == NULL || count == NULL) {
1644     retval = MPI_ERR_ARG;
1645   } else if (datatype == MPI_DATATYPE_NULL) {
1646     retval = MPI_ERR_TYPE;
1647   } else {
1648     size = smpi_datatype_size(datatype);
1649     if (size == 0) {
1650       *count = 0;
1651     } else if (status->count % size != 0) {
1652       retval = MPI_UNDEFINED;
1653     } else {
1654       *count = smpi_mpi_get_count(status, datatype);
1655     }
1656   }
1657   smpi_bench_begin();
1658   return retval;
1659 }
1660
1661 /* The following calls are not yet implemented and will fail at runtime. */
1662 /* Once implemented, please move them above this notice. */
1663
1664 static int not_yet_implemented(void) {
1665    xbt_die("Not yet implemented");
1666    return MPI_ERR_OTHER;
1667 }
1668
1669 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
1670    return not_yet_implemented();
1671 }
1672
1673 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
1674    return not_yet_implemented();
1675 }
1676
1677 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
1678    return not_yet_implemented();
1679 }
1680
1681 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
1682    return not_yet_implemented();
1683 }
1684
1685 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
1686    return not_yet_implemented();
1687 }
1688
1689 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
1690    return not_yet_implemented();
1691 }
1692
1693 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
1694    return not_yet_implemented();
1695 }
1696
1697 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
1698    return not_yet_implemented();
1699 }
1700
1701 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
1702    return not_yet_implemented();
1703 }
1704
1705 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
1706    return not_yet_implemented();
1707 }
1708
1709 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
1710    return not_yet_implemented();
1711 }
1712
1713 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
1714    return not_yet_implemented();
1715 }
1716
1717 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
1718    return not_yet_implemented();
1719 }
1720
1721 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
1722    return not_yet_implemented();
1723 }
1724
1725 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
1726    return not_yet_implemented();
1727 }
1728
1729 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
1730    return not_yet_implemented();
1731 }
1732
1733 int PMPI_Error_class(int errorcode, int* errorclass) {
1734    return not_yet_implemented();
1735 }
1736
1737 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
1738    return not_yet_implemented();
1739 }
1740
1741 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
1742    return not_yet_implemented();
1743 }
1744
1745 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
1746    return not_yet_implemented();
1747 }
1748
1749 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
1750    return not_yet_implemented();
1751 }
1752
1753 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
1754    return not_yet_implemented();
1755 }
1756
1757 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* newtype) {
1758    return not_yet_implemented();
1759 }
1760
1761 int PMPI_Cancel(MPI_Request* request) {
1762    return not_yet_implemented();
1763 }
1764
1765 int PMPI_Buffer_attach(void* buffer, int size) {
1766    return not_yet_implemented();
1767 }
1768
1769 int PMPI_Buffer_detach(void* buffer, int* size) {
1770    return not_yet_implemented();
1771 }
1772
1773 int PMPI_Testsome(int incount, MPI_Request* requests, int* outcount, int* indices, MPI_Status* statuses) {
1774    return not_yet_implemented();
1775 }
1776
1777 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
1778    return not_yet_implemented();
1779 }
1780
1781 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
1782    return not_yet_implemented();
1783 }
1784
1785 int PMPI_Type_commit(MPI_Datatype* datatype) {
1786    return not_yet_implemented();
1787 }
1788
1789 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
1790    return not_yet_implemented();
1791 }
1792
1793 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
1794    return not_yet_implemented();
1795 }
1796
1797 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* newtype) {
1798    return not_yet_implemented();
1799 }
1800
1801 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* newtype) {
1802    return not_yet_implemented();
1803 }
1804
1805 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* newtype) {
1806    return not_yet_implemented();
1807 }
1808
1809 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1810    return not_yet_implemented();
1811 }
1812
1813 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1814    return not_yet_implemented();
1815 }
1816
1817 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
1818    return not_yet_implemented();
1819 }
1820
1821 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
1822    return not_yet_implemented();
1823 }
1824
1825 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1826    return not_yet_implemented();
1827 }
1828
1829 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1830    return not_yet_implemented();
1831 }
1832
1833 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1834    return not_yet_implemented();
1835 }
1836
1837 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
1838    return not_yet_implemented();
1839 }
1840
1841 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
1842    return not_yet_implemented();
1843 }
1844
1845 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1846    return not_yet_implemented();
1847 }
1848
1849 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1850    return not_yet_implemented();
1851 }
1852
1853 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
1854    return not_yet_implemented();
1855 }
1856
1857 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
1858    return not_yet_implemented();
1859 }
1860
1861 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
1862    return not_yet_implemented();
1863 }
1864
1865 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
1866    return not_yet_implemented();
1867 }
1868
1869 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1870    return not_yet_implemented();
1871 }
1872
1873 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
1874    return not_yet_implemented();
1875 }
1876
1877 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
1878    return not_yet_implemented();
1879 }
1880
1881 int PMPI_Keyval_free(int* keyval) {
1882    return not_yet_implemented();
1883 }
1884
1885 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
1886    return not_yet_implemented();
1887 }
1888
1889 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
1890    return not_yet_implemented();
1891 }
1892
1893 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses) {
1894    return not_yet_implemented();
1895 }
1896
1897 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
1898    return not_yet_implemented();
1899 }
1900
1901 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
1902    return not_yet_implemented();
1903 }
1904
1905 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1906    return not_yet_implemented();
1907 }
1908
1909 int PMPI_Initialized(int* flag) {
1910    return not_yet_implemented();
1911 }