Logo AND Algorithmique Numérique Distribuée

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