Logo AND Algorithmique Numérique Distribuée

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