Logo AND Algorithmique Numérique Distribuée

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