Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
5436c0a4ccbae37f60893d0890395bba3773eb65
[simgrid.git] / src / smpi / smpi_mpi.c
1 /* $Id$tag */
2
3 #include "private.h"
4 #include "smpi_coll_private.h"
5 #include "smpi_mpi_dt_private.h"
6
7 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_mpi, smpi,
8                                 "Logging specific to SMPI (mpi)");
9
10 /* MPI User level calls */
11
12 int MPI_Init(int* argc, char*** argv) {
13   smpi_process_init(argc, argv);
14   smpi_bench_begin();
15   return MPI_SUCCESS;
16 }
17
18 int MPI_Finalize(void) {
19   smpi_bench_end();
20   smpi_process_destroy();
21   return MPI_SUCCESS;
22 }
23
24 int MPI_Init_thread(int* argc, char*** argv, int required, int* provided) {
25   if(provided != NULL) {
26     *provided = MPI_THREAD_MULTIPLE;
27   }
28   return MPI_Init(argc, argv);
29 }
30
31 int MPI_Query_thread(int* provided) {
32   int retval;
33
34   smpi_bench_end();
35   if(provided == NULL) {
36     retval = MPI_ERR_ARG;
37   } else {
38     *provided = MPI_THREAD_MULTIPLE;
39     retval = MPI_SUCCESS;
40   }
41   smpi_bench_begin();
42   return retval;
43 }
44
45 int MPI_Is_thread_main(int* flag) {
46   int retval;
47
48   smpi_bench_end();
49   if(flag == NULL) {
50     retval = MPI_ERR_ARG;
51   } else {
52     *flag = smpi_process_index() == 0;
53     retval = MPI_SUCCESS;
54   }
55   smpi_bench_begin();
56   return retval;
57 }
58
59 int MPI_Abort(MPI_Comm comm, int errorcode) {
60   smpi_bench_end();
61   smpi_process_destroy();
62   // FIXME: should kill all processes in comm instead
63   SIMIX_process_kill(SIMIX_process_self());
64   return MPI_SUCCESS;
65 }
66
67 double MPI_Wtime(void) {
68   double time;
69
70   smpi_bench_end();
71   time = SIMIX_get_clock();
72   smpi_bench_begin();
73   return time;
74 }
75
76 int MPI_Type_size(MPI_Datatype datatype, size_t* size) {
77   int retval;
78
79   smpi_bench_end();
80   if(datatype == MPI_DATATYPE_NULL) {
81     retval = MPI_ERR_TYPE;
82   } else if(size == NULL) {
83     retval = MPI_ERR_ARG;
84   } else {
85     *size = smpi_datatype_size(datatype);
86     retval = MPI_SUCCESS;
87   }
88   smpi_bench_begin();
89   return retval;
90 }
91
92 int MPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint* lb, MPI_Aint* extent) {
93   int retval;
94
95   smpi_bench_end();
96   if(datatype == MPI_DATATYPE_NULL) {
97     retval = MPI_ERR_TYPE;
98   } else if(lb == NULL || extent == NULL) {
99     retval = MPI_ERR_ARG;
100   } else {
101     retval = smpi_datatype_extent(datatype, lb, extent);
102   }
103   smpi_bench_begin();
104   return retval;
105 }
106
107 int MPI_Type_lb(MPI_Datatype datatype, MPI_Aint* disp) {
108   int retval;
109
110   smpi_bench_end();
111   if(datatype == MPI_DATATYPE_NULL) {
112     retval = MPI_ERR_TYPE;
113   } else if(disp == NULL) {
114     retval = MPI_ERR_ARG;
115   } else {
116     *disp = smpi_datatype_lb(datatype);
117     retval = MPI_SUCCESS;
118   }
119   smpi_bench_begin();
120   return retval;
121 }
122
123 int MPI_Type_ub(MPI_Datatype datatype, MPI_Aint* disp) {
124   int retval;
125
126   smpi_bench_end();
127   if(datatype == MPI_DATATYPE_NULL) {
128     retval = MPI_ERR_TYPE;
129   } else if(disp == NULL) {
130     retval = MPI_ERR_ARG;
131   } else {
132     *disp = smpi_datatype_ub(datatype);
133     retval = MPI_SUCCESS;
134   }
135   smpi_bench_begin();
136   return retval;
137 }
138
139 int MPI_Op_create(MPI_User_function* function, int commute, MPI_Op* op) {
140   int retval;
141
142   smpi_bench_end();
143   if(function == NULL || op == NULL) {
144     retval = MPI_ERR_ARG;
145   } else {
146     *op = smpi_op_new(function, commute);
147     retval = MPI_SUCCESS;
148   }
149   smpi_bench_begin();
150   return retval;
151 }
152
153 int MPI_Op_free(MPI_Op* op) {
154   int retval;
155
156   smpi_bench_end();
157   if(op == NULL) {
158     retval = MPI_ERR_ARG;
159   } else if(*op == MPI_OP_NULL) {
160     retval = MPI_ERR_OP;
161   } else {
162     smpi_op_destroy(*op);
163     *op = MPI_OP_NULL;
164     retval = MPI_SUCCESS;
165   }
166   smpi_bench_begin();
167   return retval;
168 }
169
170 int MPI_Group_free(MPI_Group *group) {
171   int retval;
172
173   smpi_bench_end();
174   if(group == NULL) {
175     retval = MPI_ERR_ARG;
176   } else {
177     smpi_group_destroy(*group);
178     *group = MPI_GROUP_NULL;
179     retval = MPI_SUCCESS;
180   }
181   smpi_bench_begin();
182   return retval;
183 }
184
185 int MPI_Group_size(MPI_Group group, int* size) {
186   int retval;
187
188   smpi_bench_end();
189   if(group == MPI_GROUP_NULL) {
190     retval = MPI_ERR_GROUP;
191   } else if(size == NULL) {
192     retval = MPI_ERR_ARG;
193   } else {
194     *size = smpi_group_size(group);
195     retval = MPI_SUCCESS;
196   }
197   smpi_bench_begin();
198   return retval;
199 }
200
201 int MPI_Group_rank(MPI_Group group, int* rank) {
202   int retval;
203
204   smpi_bench_end();
205   if(group == MPI_GROUP_NULL) {
206     retval = MPI_ERR_GROUP;
207   } else if(rank == NULL) {
208     retval = MPI_ERR_ARG;
209   } else {
210     *rank = smpi_group_rank(group, smpi_process_index());
211     retval = MPI_SUCCESS;
212   }
213   smpi_bench_begin();
214   return retval;
215 }
216
217 int MPI_Group_translate_ranks (MPI_Group group1, int n, int* ranks1, MPI_Group group2, int* ranks2) {
218   int retval, i, index;
219
220   smpi_bench_end();
221   if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
222     retval = MPI_ERR_GROUP;
223   } else {
224     for(i = 0; i < n; i++) {
225       index = smpi_group_index(group1, ranks1[i]);
226       ranks2[i] = smpi_group_rank(group2, index);
227     }
228     retval = MPI_SUCCESS;
229   }
230   smpi_bench_begin();
231   return retval;
232 }
233
234 int MPI_Group_compare(MPI_Group group1, MPI_Group group2, int* result) {
235   int retval;
236
237   smpi_bench_end();
238   if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
239     retval = MPI_ERR_GROUP;
240   } else if(result == NULL) {
241     retval = MPI_ERR_ARG;
242   } else {
243     *result = smpi_group_compare(group1, group2);
244     retval = MPI_SUCCESS;
245   }
246   smpi_bench_begin();
247   return retval;
248 }
249
250 int MPI_Group_union(MPI_Group group1, MPI_Group group2, MPI_Group* newgroup) {
251   int retval, i, proc1, proc2, size, size2;
252
253   smpi_bench_end();
254   if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
255     retval = MPI_ERR_GROUP;
256   } else if(newgroup == NULL) {
257     retval = MPI_ERR_ARG;
258   } else {
259     size = smpi_group_size(group1);
260     size2 = smpi_group_size(group2);
261     for(i = 0; i < size2; i++) {
262       proc2 = smpi_group_index(group2, i);
263       proc1 = smpi_group_rank(group1, proc2);
264       if(proc1 == MPI_UNDEFINED) {
265         size++;
266       }
267     }
268     if(size == 0) {
269       *newgroup = MPI_GROUP_EMPTY;
270     } else {
271       *newgroup = smpi_group_new(size);
272       size2 = smpi_group_size(group1);
273       for(i = 0; i < size2; i++) {
274         proc1 = smpi_group_index(group1, i);
275         smpi_group_set_mapping(*newgroup, proc1, i);
276       }
277       for(i = size2; i < size; i++) {
278         proc2 = smpi_group_index(group2, i - size2);
279         smpi_group_set_mapping(*newgroup, proc2, i);
280       }
281     }
282     smpi_group_use(*newgroup);
283     retval = MPI_SUCCESS;
284   }
285   smpi_bench_begin();
286   return retval;
287 }
288
289 int MPI_Group_intersection(MPI_Group group1, MPI_Group group2, MPI_Group* newgroup) {
290    int retval, i, proc1, proc2, size, size2;
291
292   smpi_bench_end();
293   if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
294     retval = MPI_ERR_GROUP;
295   } else if(newgroup == NULL) {
296     retval = MPI_ERR_ARG;
297   } else {
298     size = smpi_group_size(group1);
299     size2 = smpi_group_size(group2);
300     for(i = 0; i < size2; i++) {
301       proc2 = smpi_group_index(group2, i);
302       proc1 = smpi_group_rank(group1, proc2);
303       if(proc1 == MPI_UNDEFINED) {
304         size--;
305       }
306     }
307     if(size == 0) {
308       *newgroup = MPI_GROUP_EMPTY;
309     } else {
310       *newgroup = smpi_group_new(size);
311       size2 = smpi_group_size(group1);
312       for(i = 0; i < size2; i++) {
313         proc1 = smpi_group_index(group1, i);
314         proc2 = smpi_group_rank(group2, proc1);
315         if(proc2 != MPI_UNDEFINED) {
316           smpi_group_set_mapping(*newgroup, proc1, i);
317         }
318       }
319     }
320     smpi_group_use(*newgroup);
321     retval = MPI_SUCCESS;
322   }
323   smpi_bench_begin();
324   return retval;
325 }
326
327 int MPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group* newgroup) {
328   int retval, i, proc1, proc2, size, size2;
329
330   smpi_bench_end();
331   if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
332     retval = MPI_ERR_GROUP;
333   } else if(newgroup == NULL) {
334     retval = MPI_ERR_ARG;
335   } else {
336     size = size2 = smpi_group_size(group1);
337     for(i = 0; i < size2; i++) {
338       proc1 = smpi_group_index(group1, i);
339       proc2 = smpi_group_rank(group2, proc1);
340       if(proc2 != MPI_UNDEFINED) {
341         size--;
342       }
343     }
344     if(size == 0) {
345       *newgroup = MPI_GROUP_EMPTY;
346     } else {
347       *newgroup = smpi_group_new(size);
348       for(i = 0; i < size2; i++) {
349         proc1 = smpi_group_index(group1, i);
350         proc2 = smpi_group_rank(group2, proc1);
351         if(proc2 == MPI_UNDEFINED) {
352           smpi_group_set_mapping(*newgroup, proc1, i);
353         }
354       }
355     }
356     smpi_group_use(*newgroup);
357     retval = MPI_SUCCESS;
358   }
359   smpi_bench_begin();
360   return retval;
361 }
362
363 int MPI_Group_incl(MPI_Group group, int n, int* ranks, MPI_Group* newgroup) {
364   int retval, i, index;
365
366   smpi_bench_end();
367   if(group == MPI_GROUP_NULL) {
368     retval = MPI_ERR_GROUP;
369   } else if(newgroup == NULL) {
370     retval = MPI_ERR_ARG;
371   } else {
372     if(n == 0) {
373       *newgroup = MPI_GROUP_EMPTY;
374     } else if(n == smpi_group_size(group)) {
375       *newgroup = group;
376     } else {
377       *newgroup = smpi_group_new(n);
378       for(i = 0; i < n; i++) {
379         index = smpi_group_index(group, ranks[i]);
380         smpi_group_set_mapping(*newgroup, index, i);
381       }
382     }
383     smpi_group_use(*newgroup);
384     retval = MPI_SUCCESS;
385   }
386   smpi_bench_begin();
387   return retval;
388 }
389
390 int MPI_Group_excl(MPI_Group group, int n, int* ranks, MPI_Group* newgroup) {
391   int retval, i, size, rank, index;
392
393   smpi_bench_end();
394   if(group == MPI_GROUP_NULL) {
395     retval = MPI_ERR_GROUP;
396   } else if(newgroup == NULL) {
397     retval = MPI_ERR_ARG;
398   } else {
399     if(n == 0) {
400       *newgroup = group;
401     } else if(n == smpi_group_size(group)) {
402       *newgroup = MPI_GROUP_EMPTY;
403     } else {
404       size = smpi_group_size(group) - n;
405       *newgroup = smpi_group_new(size);
406       rank = 0;
407       while(rank < size) {
408         for(i = 0; i < n; i++) {
409           if(ranks[i] == rank) {
410             break;
411           }
412         }
413         if(i >= n) {
414           index = smpi_group_index(group, rank);
415           smpi_group_set_mapping(*newgroup, index, rank);
416           rank++;
417         }
418       }
419     }
420     smpi_group_use(*newgroup);
421     retval = MPI_SUCCESS;
422   }
423   smpi_bench_begin();
424   return retval;
425 }
426
427 int MPI_Group_range_incl(MPI_Group group, int n, int ranges[][3], MPI_Group* newgroup) {
428   int retval, i, j, rank, size, index;
429
430   smpi_bench_end();
431   if(group == MPI_GROUP_NULL) {
432     retval = MPI_ERR_GROUP;
433   } else if(newgroup == NULL) {
434     retval = MPI_ERR_ARG;
435   } else {
436     if(n == 0) {
437       *newgroup = MPI_GROUP_EMPTY;
438     } else {
439       size = 0;
440       for(i = 0; i < n; i++) {
441         for(rank = ranges[i][0]; /* First */
442             rank >= 0 && rank <= ranges[i][1]; /* Last */
443             rank += ranges[i][2] /* Stride */) {
444           size++;
445         }
446       }
447       if(size == smpi_group_size(group)) {
448         *newgroup = group;
449       } else {
450         *newgroup = smpi_group_new(size);
451         j = 0;
452         for(i = 0; i < n; i++) {
453           for(rank = ranges[i][0]; /* First */
454               rank >= 0 && rank <= ranges[i][1]; /* Last */
455               rank += ranges[i][2] /* Stride */) {
456             index = smpi_group_index(group, rank);
457             smpi_group_set_mapping(*newgroup, index, j);
458             j++;
459           }
460         }
461       }
462     }
463     smpi_group_use(*newgroup);
464     retval = MPI_SUCCESS;
465   }
466   smpi_bench_begin();
467   return retval;
468 }
469
470 int MPI_Group_range_excl(MPI_Group group, int n, int ranges[][3], MPI_Group* newgroup) {
471   int retval, i, newrank, rank, size, index, add;
472
473   smpi_bench_end();
474   if(group == MPI_GROUP_NULL) {
475     retval = MPI_ERR_GROUP;
476   } else if(newgroup == NULL) {
477     retval = MPI_ERR_ARG;
478   } else {
479     if(n == 0) {
480       *newgroup = group;
481     } else {
482       size = smpi_group_size(group);
483       for(i = 0; i < n; i++) {
484         for(rank = ranges[i][0]; /* First */
485             rank >= 0 && rank <= ranges[i][1]; /* Last */
486             rank += ranges[i][2] /* Stride */) {
487           size--;
488         }
489       }
490       if(size == 0) {
491         *newgroup = MPI_GROUP_EMPTY;
492       } else {
493         *newgroup = smpi_group_new(size);
494         newrank = 0;
495         while(newrank < size) {
496           for(i = 0; i < n; i++) {
497             add = 1;
498             for(rank = ranges[i][0]; /* First */
499                 rank >= 0 && rank <= ranges[i][1]; /* Last */
500                 rank += ranges[i][2] /* Stride */) {
501               if(rank == newrank) {
502                 add = 0;
503                 break;
504               }
505             }
506             if(add == 1) {
507               index = smpi_group_index(group, newrank);
508               smpi_group_set_mapping(*newgroup, index, newrank);
509             }
510           }
511         }
512       }
513     }
514     smpi_group_use(*newgroup);
515     retval = MPI_SUCCESS;
516   }
517   smpi_bench_begin();
518   return retval;
519 }
520
521 int MPI_Comm_rank(MPI_Comm comm, int* rank) {
522   int retval;
523
524   smpi_bench_end();
525   if(comm == MPI_COMM_NULL) {
526     retval = MPI_ERR_COMM;
527   } else {
528     *rank = smpi_comm_rank(comm);
529     retval = MPI_SUCCESS;
530   }
531   smpi_bench_begin();
532   return retval;
533 }
534
535 int MPI_Comm_size(MPI_Comm comm, int* size) {
536   int retval;
537
538   smpi_bench_end();
539   if(comm == MPI_COMM_NULL) {
540     retval = MPI_ERR_COMM;
541   } else if(size == NULL) {
542     retval = MPI_ERR_ARG;
543   } else {
544     *size = smpi_comm_size(comm);
545     retval = MPI_SUCCESS;
546   }
547   smpi_bench_begin();
548   return retval;
549 }
550
551 int MPI_Comm_group(MPI_Comm comm, MPI_Group* group) {
552   int retval;
553
554   smpi_bench_end();
555   if(comm == MPI_COMM_NULL) {
556     retval = MPI_ERR_COMM;
557   } else if(group == NULL) {
558     retval = MPI_ERR_ARG;
559   } else {
560     *group = smpi_comm_group(comm);
561     retval = MPI_SUCCESS;
562   }
563   smpi_bench_begin();
564   return retval;
565 }
566
567 int MPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int* result) {
568   int retval;
569
570   smpi_bench_end();
571   if(comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
572     retval = MPI_ERR_COMM;
573   } else if(result == NULL) {
574     retval = MPI_ERR_ARG;
575   } else {
576     if(comm1 == comm2) { /* Same communicators means same groups */
577       *result = MPI_IDENT;
578     } else {
579       *result = smpi_group_compare(smpi_comm_group(comm1), smpi_comm_group(comm2));
580       if(*result == MPI_IDENT) {
581         *result = MPI_CONGRUENT;
582       }
583     }
584     retval = MPI_SUCCESS;
585   }
586   smpi_bench_begin();
587   return retval;
588 }
589
590 int MPI_Comm_dup(MPI_Comm comm, MPI_Comm* newcomm) {
591   int retval;
592
593   smpi_bench_end();
594   if(comm == MPI_COMM_NULL) {
595     retval = MPI_ERR_COMM;
596   } else if(newcomm == NULL) {
597     retval = MPI_ERR_ARG;
598   } else {
599     *newcomm = smpi_comm_new(smpi_comm_group(comm));
600     retval = MPI_SUCCESS;
601   }
602   smpi_bench_begin();
603   return retval;
604 }
605
606 int MPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm* newcomm) {
607   int retval;
608
609   smpi_bench_end();
610   if(comm == MPI_COMM_NULL) {
611     retval = MPI_ERR_COMM;
612   } else if(group == MPI_GROUP_NULL) {
613     retval = MPI_ERR_GROUP;
614   } else if(newcomm == NULL) {
615     retval = MPI_ERR_ARG;
616   } else {
617     *newcomm = smpi_comm_new(group);
618     retval = MPI_SUCCESS;
619   }
620   smpi_bench_begin();
621   return retval;
622 }
623
624 int MPI_Comm_free(MPI_Comm* comm) {
625   int retval;
626
627   smpi_bench_end();
628   if(comm == NULL) {
629     retval = MPI_ERR_ARG;
630   } else if(*comm == MPI_COMM_NULL) {
631     retval = MPI_ERR_COMM;
632   } else {
633     smpi_comm_destroy(*comm);
634     *comm = MPI_COMM_NULL;
635     retval = MPI_SUCCESS;
636   }
637   smpi_bench_begin();
638   return retval;
639 }
640
641 int MPI_Irecv(void* buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request* request) {
642   int retval;
643
644   smpi_bench_end();
645   if(request == NULL) {
646     retval = MPI_ERR_ARG;
647   } else {
648     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
649     retval = MPI_SUCCESS;
650   }
651   smpi_bench_begin();
652   return retval;
653 }
654
655 int MPI_Isend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
656   int retval;
657
658   smpi_bench_end();
659   if(request == NULL) {
660     retval = MPI_ERR_ARG;
661   } else {
662     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
663     retval = MPI_SUCCESS;
664   }
665   smpi_bench_begin();
666   return retval;
667 }
668
669 int MPI_Recv(void* buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status* status) {
670   smpi_bench_end();
671   smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
672   smpi_bench_begin();
673   return MPI_SUCCESS;
674 }
675
676 int MPI_Send(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
677   smpi_bench_end();
678   smpi_mpi_send(buf, count, datatype, dst, tag, comm);
679   smpi_bench_begin();
680   return MPI_SUCCESS;
681 }
682
683 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) {
684   smpi_bench_end();
685   smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src, recvtag, comm, status);
686   smpi_bench_begin();
687   return MPI_SUCCESS;
688 }
689
690 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) {
691   //TODO: suboptimal implementation
692   void* recvbuf;
693   int retval, size;
694
695   size = smpi_datatype_size(datatype) * count;
696   recvbuf = xbt_new(char, size);
697   retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
698   memcpy(buf, recvbuf, size * sizeof(char));
699   xbt_free(recvbuf);
700   return retval;
701 }
702
703 int MPI_Test(MPI_Request* request, int* flag, MPI_Status* status) {
704   int retval;
705
706   smpi_bench_end();
707   if(request == NULL || flag == NULL) {
708     retval = MPI_ERR_ARG;
709   } else if(*request == MPI_REQUEST_NULL) {
710     retval = MPI_ERR_REQUEST;
711   } else {
712     *flag = smpi_mpi_test(request, status);
713     retval = MPI_SUCCESS;
714   }
715   smpi_bench_begin();
716   return retval;
717 }
718
719 int MPI_Testany(int count, MPI_Request requests[], int* index, int* flag, MPI_Status* status) {
720   int retval;
721
722   smpi_bench_end();
723   if(index == NULL || flag == NULL) {
724     retval = MPI_ERR_ARG;
725   } else {
726     *flag = smpi_mpi_testany(count, requests, index, status);
727     retval = MPI_SUCCESS;
728   }
729   smpi_bench_begin();
730   return retval;
731 }
732
733 int MPI_Wait(MPI_Request* request, MPI_Status* status) {
734   int retval;
735
736   smpi_bench_end();
737   if(request == NULL) {
738     retval = MPI_ERR_ARG;
739   } else if(*request == MPI_REQUEST_NULL) {
740     retval = MPI_ERR_REQUEST;
741   } else {
742     smpi_mpi_wait(request, status);
743     retval = MPI_SUCCESS;
744   }
745   smpi_bench_begin();
746   return retval;
747 }
748
749 int MPI_Waitany(int count, MPI_Request requests[], int* index, MPI_Status* status) {
750   int retval;
751
752   smpi_bench_end();
753   if(index == NULL) {
754     retval = MPI_ERR_ARG;
755   } else {
756     *index = smpi_mpi_waitany(count, requests, status);
757     retval = MPI_SUCCESS;
758   }
759   smpi_bench_begin();
760   return retval;
761 }
762
763 int MPI_Waitall(int count, MPI_Request requests[],  MPI_Status status[]) {
764   smpi_bench_end();
765   smpi_mpi_waitall(count, requests, status);
766   smpi_bench_begin();
767   return MPI_SUCCESS;
768 }
769
770 int MPI_Waitsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[]) {
771   int retval;
772
773   smpi_bench_end();
774   if(outcount == NULL || indices == NULL) {
775     retval = MPI_ERR_ARG;
776   } else {
777     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
778     retval = MPI_SUCCESS;
779   }
780   smpi_bench_begin();
781   return retval;
782 }
783
784 int MPI_Bcast(void* buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm) {
785   int retval;
786
787   smpi_bench_end();
788   if(comm == MPI_COMM_NULL) {
789     retval = MPI_ERR_COMM;
790   } else {
791     smpi_mpi_bcast(buf, count, datatype, root, comm);
792     retval = MPI_SUCCESS;
793   }
794   smpi_bench_begin();
795   return retval;
796 }
797
798 int MPI_Barrier(MPI_Comm comm) {
799   int retval;
800
801   smpi_bench_end();
802   if(comm == MPI_COMM_NULL) {
803     retval = MPI_ERR_COMM;
804   } else {
805     smpi_mpi_barrier(comm);
806     retval = MPI_SUCCESS;
807   }
808   smpi_bench_begin();
809   return retval;
810 }
811
812 int MPI_Gather(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm) {
813   int retval;
814
815   smpi_bench_end();
816   if(comm == MPI_COMM_NULL) {
817     retval = MPI_ERR_COMM;
818   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
819     retval = MPI_ERR_TYPE;
820   } else {
821     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
822     retval = MPI_SUCCESS;
823   }
824   smpi_bench_begin();
825   return retval;
826 }
827
828 int MPI_Gatherv(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int* recvcounts, int* displs, MPI_Datatype recvtype, int root, MPI_Comm comm) {
829   int retval;
830
831   smpi_bench_end();
832   if(comm == MPI_COMM_NULL) {
833     retval = MPI_ERR_COMM;
834   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
835     retval = MPI_ERR_TYPE;
836   } else if(recvcounts == NULL || displs == NULL) {
837     retval = MPI_ERR_ARG;
838   } else {
839     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm);
840     retval = MPI_SUCCESS;
841   }
842   smpi_bench_begin();
843   return retval;
844 }
845
846 int MPI_Allgather(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm) {
847   int retval;
848
849   smpi_bench_end();
850   if(comm == MPI_COMM_NULL) {
851     retval = MPI_ERR_COMM;
852   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
853     retval = MPI_ERR_TYPE;
854   } else {
855     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
856     retval = MPI_SUCCESS;
857   }
858   smpi_bench_begin();
859   return retval;
860 }
861
862 int MPI_Allgatherv(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int* recvcounts, int* displs, MPI_Datatype recvtype, MPI_Comm comm) {
863   int retval;
864
865   smpi_bench_end();
866   if(comm == MPI_COMM_NULL) {
867     retval = MPI_ERR_COMM;
868   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
869     retval = MPI_ERR_TYPE;
870   } else if(recvcounts == NULL || displs == NULL) {
871     retval = MPI_ERR_ARG;
872   } else {
873     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
874     retval = MPI_SUCCESS;
875   }
876   smpi_bench_begin();
877   return retval;
878 }
879
880 int MPI_Scatter(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm) {
881   int retval;
882
883   smpi_bench_end();
884   if(comm == MPI_COMM_NULL) {
885     retval = MPI_ERR_COMM;
886   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
887     retval = MPI_ERR_TYPE;
888   } else {
889     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
890     retval = MPI_SUCCESS;
891   }
892   smpi_bench_begin();
893   return retval;
894 }
895
896 int MPI_Scatterv(void* sendbuf, int* sendcounts, int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm) {
897   int retval;
898
899   smpi_bench_end();
900   if(comm == MPI_COMM_NULL) {
901     retval = MPI_ERR_COMM;
902   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
903     retval = MPI_ERR_TYPE;
904   } else if(sendcounts == NULL || displs == NULL) {
905     retval = MPI_ERR_ARG;
906   } else {
907     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
908     retval = MPI_SUCCESS;
909   }
910   smpi_bench_begin();
911   return retval;
912 }
913
914 int MPI_Reduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm) {
915   int retval;
916
917   smpi_bench_end();
918   if(comm == MPI_COMM_NULL) {
919     retval = MPI_ERR_COMM;
920   } else if(datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
921     retval = MPI_ERR_ARG;
922   } else {
923     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
924     retval = MPI_SUCCESS;
925   }
926   smpi_bench_begin();
927   return retval;
928 }
929
930 int MPI_Allreduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) {
931   int retval;
932
933   smpi_bench_end();
934   if(comm == MPI_COMM_NULL) {
935     retval = MPI_ERR_COMM;
936   } else if(datatype == MPI_DATATYPE_NULL) {
937     retval = MPI_ERR_TYPE;
938   } else if(op == MPI_OP_NULL) {
939     retval = MPI_ERR_OP;
940   } else {
941     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
942     retval = MPI_SUCCESS;
943   }
944   smpi_bench_begin();
945   return retval;
946 }
947
948 int MPI_Reduce_scatter(void* sendbuf, void* recvbuf, int* recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) {
949   int retval, i, rank, size, count;
950   int* displs;
951
952   smpi_bench_end();
953   if(comm == MPI_COMM_NULL) {
954     retval = MPI_ERR_COMM;
955   } else if(datatype == MPI_DATATYPE_NULL) {
956     retval = MPI_ERR_TYPE;
957   } else if(op == MPI_OP_NULL) {
958     retval = MPI_ERR_OP;
959   } else if(recvcounts == NULL) {
960     retval = MPI_ERR_ARG;
961   } else {
962     /* arbitrarily choose root as rank 0 */
963     /* TODO: faster direct implementation ? */
964     rank = smpi_comm_rank(comm);
965     size = smpi_comm_size(comm);
966     count = 0;
967     displs = xbt_new(int, size);
968     for(i = 0; i < size; i++) {
969       count += recvcounts[i];
970       displs[i] = 0;
971     }
972     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
973     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf, recvcounts[rank], datatype, 0, comm);
974     xbt_free(displs);
975     retval = MPI_SUCCESS;
976   }
977   smpi_bench_begin();
978   return retval;
979 }
980
981 /**
982  * MPI_Alltoall user entry point
983  * 
984  * Uses the logic of OpenMPI (upto 1.2.7 or greater) for the optimizations
985  * ompi/mca/coll/tuned/coll_tuned_module.c
986  **/
987
988 int MPI_Alltoall(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm) {
989   int retval, size, sendsize;
990
991   smpi_bench_end();
992   if(comm == MPI_COMM_NULL) {
993     retval = MPI_ERR_COMM;
994   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
995     retval = MPI_ERR_TYPE;
996   } else {
997     size = smpi_comm_size(comm);
998     sendsize = smpi_datatype_size(sendtype) * sendcount;
999     if(sendsize < 200 && size > 12) {
1000       retval = smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1001     } else if(sendsize < 3000) {
1002       retval = smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1003     } else {
1004       retval = smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1005     }
1006   }
1007   smpi_bench_begin();
1008   return retval;
1009 }
1010
1011 int MPI_Alltoallv(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype sendtype, void* recvbuf, int *recvcounts, int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm) {
1012   int retval;
1013
1014   smpi_bench_end();
1015   if(comm == MPI_COMM_NULL) {
1016     retval = MPI_ERR_COMM;
1017   } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
1018     retval = MPI_ERR_TYPE;
1019   } else if(sendcounts == NULL || senddisps == NULL || recvcounts == NULL || recvdisps == NULL) {
1020     retval = MPI_ERR_ARG;
1021   } else {
1022     retval = smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype, recvbuf, recvcounts, recvdisps, recvtype, comm); 
1023   }
1024   smpi_bench_begin();
1025   return retval;
1026 }