Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
64f0161b41c108d176ec692c73fcd813720c7b7a
[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     *flag=TRUE;
1081     retval = MPI_ERR_REQUEST;
1082   } else {
1083     *flag = smpi_mpi_test(request, status);
1084     retval = MPI_SUCCESS;
1085   }
1086   smpi_bench_begin();
1087   return retval;
1088 }
1089
1090 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1091                 MPI_Status * status)
1092 {
1093   int retval;
1094
1095   smpi_bench_end();
1096   if (index == NULL || flag == NULL) {
1097     retval = MPI_ERR_ARG;
1098   } else {
1099     *flag = smpi_mpi_testany(count, requests, index, status);
1100     retval = MPI_SUCCESS;
1101   }
1102   smpi_bench_begin();
1103   return retval;
1104 }
1105
1106 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1107 {
1108   int retval;
1109
1110   smpi_bench_end();
1111   if (flag == NULL) {
1112     retval = MPI_ERR_ARG;
1113   } else {
1114     *flag = smpi_mpi_testall(count, requests, statuses);
1115     retval = MPI_SUCCESS;
1116   }
1117   smpi_bench_begin();
1118   return retval;
1119 }
1120
1121 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1122   int retval;
1123   smpi_bench_end();
1124
1125   if (status == NULL) {
1126      retval = MPI_ERR_ARG;
1127   }else if (comm == MPI_COMM_NULL) {
1128         retval = MPI_ERR_COMM;
1129   } else {
1130         smpi_mpi_probe(source, tag, comm, status);
1131         retval = MPI_SUCCESS;
1132   }
1133   smpi_bench_begin();
1134   return retval;
1135 }
1136
1137
1138 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1139   int retval;
1140   smpi_bench_end();
1141
1142   if (flag == NULL) {
1143        retval = MPI_ERR_ARG;
1144     }else if (status == NULL) {
1145      retval = MPI_ERR_ARG;
1146   }else if (comm == MPI_COMM_NULL) {
1147         retval = MPI_ERR_COMM;
1148   } else {
1149         smpi_mpi_iprobe(source, tag, comm, flag, status);
1150         retval = MPI_SUCCESS;
1151   }
1152   smpi_bench_begin();
1153   return retval;
1154 }
1155
1156 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1157 {
1158   int retval;
1159
1160   smpi_bench_end();
1161 #ifdef HAVE_TRACING
1162   int rank = request && (*request)->comm != MPI_COMM_NULL
1163       ? smpi_process_index()
1164       : -1;
1165   TRACE_smpi_computing_out(rank);
1166
1167   MPI_Group group = smpi_comm_group((*request)->comm);
1168   int src_traced = smpi_group_index(group, (*request)->src);
1169   int dst_traced = smpi_group_index(group, (*request)->dst);
1170   int is_wait_for_receive = (*request)->recv;
1171   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1172 #endif
1173   if (request == NULL) {
1174     retval = MPI_ERR_ARG;
1175   } else if (*request == MPI_REQUEST_NULL) {
1176     retval = MPI_ERR_REQUEST;
1177   } else {
1178     smpi_mpi_wait(request, status);
1179     retval = MPI_SUCCESS;
1180   }
1181 #ifdef HAVE_TRACING
1182   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1183   if (is_wait_for_receive) {
1184     TRACE_smpi_recv(rank, src_traced, dst_traced);
1185   }
1186   TRACE_smpi_computing_in(rank);
1187
1188 #endif
1189   smpi_bench_begin();
1190   return retval;
1191 }
1192
1193 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1194 {
1195   int retval;
1196
1197   smpi_bench_end();
1198 #ifdef HAVE_TRACING
1199   //save requests information for tracing
1200   int i;
1201   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), NULL);
1202   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), NULL);
1203   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), NULL);
1204   for (i = 0; i < count; i++) {
1205     MPI_Request req = requests[i];      //already received requests are no longer valid
1206     if (req) {
1207       int *asrc = xbt_new(int, 1);
1208       int *adst = xbt_new(int, 1);
1209       int *arecv = xbt_new(int, 1);
1210       *asrc = req->src;
1211       *adst = req->dst;
1212       *arecv = req->recv;
1213       xbt_dynar_insert_at(srcs, i, asrc);
1214       xbt_dynar_insert_at(dsts, i, adst);
1215       xbt_dynar_insert_at(recvs, i, arecv);
1216       xbt_free(asrc);
1217       xbt_free(adst);
1218       xbt_free(arecv);
1219     } else {
1220       int *t = xbt_new(int, 1);
1221       xbt_dynar_insert_at(srcs, i, t);
1222       xbt_dynar_insert_at(dsts, i, t);
1223       xbt_dynar_insert_at(recvs, i, t);
1224       xbt_free(t);
1225     }
1226   }
1227   int rank_traced = smpi_process_index();
1228   TRACE_smpi_computing_out(rank_traced);
1229
1230   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1231
1232 #endif
1233   if (index == NULL) {
1234     retval = MPI_ERR_ARG;
1235   } else {
1236     *index = smpi_mpi_waitany(count, requests, status);
1237     retval = MPI_SUCCESS;
1238   }
1239 #ifdef HAVE_TRACING
1240   int src_traced, dst_traced, is_wait_for_receive;
1241   xbt_dynar_get_cpy(srcs, *index, &src_traced);
1242   xbt_dynar_get_cpy(dsts, *index, &dst_traced);
1243   xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
1244   if (is_wait_for_receive) {
1245     TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1246   }
1247   TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1248   //clean-up of dynars
1249   xbt_dynar_free(&srcs);
1250   xbt_dynar_free(&dsts);
1251   xbt_dynar_free(&recvs);
1252   TRACE_smpi_computing_in(rank_traced);
1253
1254 #endif
1255   smpi_bench_begin();
1256   return retval;
1257 }
1258
1259 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1260 {
1261
1262   smpi_bench_end();
1263 #ifdef HAVE_TRACING
1264   //save information from requests
1265   int i;
1266   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), NULL);
1267   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), NULL);
1268   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), NULL);
1269   for (i = 0; i < count; i++) {
1270     MPI_Request req = requests[i];      //all req should be valid in Waitall
1271     int *asrc = xbt_new(int, 1);
1272     int *adst = xbt_new(int, 1);
1273     int *arecv = xbt_new(int, 1);
1274     *asrc = req->src;
1275     *adst = req->dst;
1276     *arecv = req->recv;
1277     xbt_dynar_insert_at(srcs, i, asrc);
1278     xbt_dynar_insert_at(dsts, i, adst);
1279     xbt_dynar_insert_at(recvs, i, arecv);
1280     xbt_free(asrc);
1281     xbt_free(adst);
1282     xbt_free(arecv);
1283   }
1284   int rank_traced = smpi_process_index();
1285   TRACE_smpi_computing_out(rank_traced);
1286
1287   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1288 #endif
1289   smpi_mpi_waitall(count, requests, status);
1290 #ifdef HAVE_TRACING
1291   for (i = 0; i < count; i++) {
1292     int src_traced, dst_traced, is_wait_for_receive;
1293     xbt_dynar_get_cpy(srcs, i, &src_traced);
1294     xbt_dynar_get_cpy(dsts, i, &dst_traced);
1295     xbt_dynar_get_cpy(recvs, i, &is_wait_for_receive);
1296     if (is_wait_for_receive) {
1297       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1298     }
1299   }
1300   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1301   //clean-up of dynars
1302   xbt_dynar_free(&srcs);
1303   xbt_dynar_free(&dsts);
1304   xbt_dynar_free(&recvs);
1305   TRACE_smpi_computing_in(rank_traced);
1306 #endif
1307   smpi_bench_begin();
1308   return MPI_SUCCESS;
1309 }
1310
1311 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1312                  int *indices, MPI_Status status[])
1313 {
1314   int retval;
1315
1316   smpi_bench_end();
1317   if (outcount == NULL || indices == NULL) {
1318     retval = MPI_ERR_ARG;
1319   } else {
1320     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1321     retval = MPI_SUCCESS;
1322   }
1323   smpi_bench_begin();
1324   return retval;
1325 }
1326
1327 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1328                  int* indices, MPI_Status status[])
1329 {
1330   int retval;
1331
1332    smpi_bench_end();
1333    if (outcount == NULL || indices == NULL) {
1334      retval = MPI_ERR_ARG;
1335    } else {
1336      *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1337      retval = MPI_SUCCESS;
1338    }
1339    smpi_bench_begin();
1340    return retval;
1341 }
1342
1343
1344 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1345 {
1346   int retval;
1347
1348   smpi_bench_end();
1349 #ifdef HAVE_TRACING
1350   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1351   TRACE_smpi_computing_out(rank);
1352   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1353   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1354 #endif
1355   if (comm == MPI_COMM_NULL) {
1356     retval = MPI_ERR_COMM;
1357   } else {
1358     smpi_mpi_bcast(buf, count, datatype, root, comm);
1359     retval = MPI_SUCCESS;
1360   }
1361 #ifdef HAVE_TRACING
1362   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1363   TRACE_smpi_computing_in(rank);
1364 #endif
1365   smpi_bench_begin();
1366   return retval;
1367 }
1368
1369 int PMPI_Barrier(MPI_Comm comm)
1370 {
1371   int retval;
1372
1373   smpi_bench_end();
1374 #ifdef HAVE_TRACING
1375   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1376   TRACE_smpi_computing_out(rank);
1377   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1378 #endif
1379   if (comm == MPI_COMM_NULL) {
1380     retval = MPI_ERR_COMM;
1381   } else {
1382     smpi_mpi_barrier(comm);
1383     retval = MPI_SUCCESS;
1384   }
1385 #ifdef HAVE_TRACING
1386   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1387   TRACE_smpi_computing_in(rank);
1388 #endif
1389   smpi_bench_begin();
1390   return retval;
1391 }
1392
1393 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1394                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1395                int root, MPI_Comm comm)
1396 {
1397   int retval;
1398
1399   smpi_bench_end();
1400 #ifdef HAVE_TRACING
1401   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1402   TRACE_smpi_computing_out(rank);
1403   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1404   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1405 #endif
1406   if (comm == MPI_COMM_NULL) {
1407     retval = MPI_ERR_COMM;
1408   } else if (sendtype == MPI_DATATYPE_NULL
1409              || recvtype == MPI_DATATYPE_NULL) {
1410     retval = MPI_ERR_TYPE;
1411   } else {
1412     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1413                     recvtype, root, comm);
1414     retval = MPI_SUCCESS;
1415   }
1416 #ifdef HAVE_TRACING
1417   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1418   TRACE_smpi_computing_in(rank);
1419 #endif
1420   smpi_bench_begin();
1421   return retval;
1422 }
1423
1424 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1425                 void *recvbuf, int *recvcounts, int *displs,
1426                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1427 {
1428   int retval;
1429
1430   smpi_bench_end();
1431 #ifdef HAVE_TRACING
1432   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1433   TRACE_smpi_computing_out(rank);
1434   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1435   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1436 #endif
1437   if (comm == MPI_COMM_NULL) {
1438     retval = MPI_ERR_COMM;
1439   } else if (sendtype == MPI_DATATYPE_NULL
1440              || recvtype == MPI_DATATYPE_NULL) {
1441     retval = MPI_ERR_TYPE;
1442   } else if (recvcounts == NULL || displs == NULL) {
1443     retval = MPI_ERR_ARG;
1444   } else {
1445     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1446                      displs, recvtype, root, comm);
1447     retval = MPI_SUCCESS;
1448   }
1449 #ifdef HAVE_TRACING
1450   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1451   TRACE_smpi_computing_in(rank);
1452 #endif
1453   smpi_bench_begin();
1454   return retval;
1455 }
1456
1457 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1458                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1459                   MPI_Comm comm)
1460 {
1461   int retval;
1462
1463   smpi_bench_end();
1464 #ifdef HAVE_TRACING
1465   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1466   TRACE_smpi_computing_out(rank);
1467   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1468 #endif
1469   if (comm == MPI_COMM_NULL) {
1470     retval = MPI_ERR_COMM;
1471   } else if (sendtype == MPI_DATATYPE_NULL
1472              || recvtype == MPI_DATATYPE_NULL) {
1473     retval = MPI_ERR_TYPE;
1474   } else {
1475     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1476                        recvtype, comm);
1477     retval = MPI_SUCCESS;
1478   }
1479 #ifdef HAVE_TRACING
1480   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1481 #endif
1482   smpi_bench_begin();
1483   return retval;
1484 }
1485
1486 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1487                    void *recvbuf, int *recvcounts, int *displs,
1488                    MPI_Datatype recvtype, MPI_Comm comm)
1489 {
1490   int retval;
1491
1492   smpi_bench_end();
1493 #ifdef HAVE_TRACING
1494   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1495   TRACE_smpi_computing_out(rank);
1496   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1497 #endif
1498   if (comm == MPI_COMM_NULL) {
1499     retval = MPI_ERR_COMM;
1500   } else if (sendtype == MPI_DATATYPE_NULL
1501              || recvtype == MPI_DATATYPE_NULL) {
1502     retval = MPI_ERR_TYPE;
1503   } else if (recvcounts == NULL || displs == NULL) {
1504     retval = MPI_ERR_ARG;
1505   } else {
1506     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1507                         displs, recvtype, comm);
1508     retval = MPI_SUCCESS;
1509   }
1510 #ifdef HAVE_TRACING
1511   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1512   TRACE_smpi_computing_in(rank);
1513 #endif
1514   smpi_bench_begin();
1515   return retval;
1516 }
1517
1518 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1519                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1520                 int root, MPI_Comm comm)
1521 {
1522   int retval;
1523
1524   smpi_bench_end();
1525 #ifdef HAVE_TRACING
1526   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1527   TRACE_smpi_computing_out(rank);
1528   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1529
1530   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1531 #endif
1532   if (comm == MPI_COMM_NULL) {
1533     retval = MPI_ERR_COMM;
1534   } else if (sendtype == MPI_DATATYPE_NULL
1535              || recvtype == MPI_DATATYPE_NULL) {
1536     retval = MPI_ERR_TYPE;
1537   } else {
1538     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1539                      recvtype, root, comm);
1540     retval = MPI_SUCCESS;
1541   }
1542 #ifdef HAVE_TRACING
1543   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1544   TRACE_smpi_computing_in(rank);
1545 #endif
1546   smpi_bench_begin();
1547   return retval;
1548 }
1549
1550 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1551                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1552                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1553 {
1554   int retval;
1555
1556   smpi_bench_end();
1557 #ifdef HAVE_TRACING
1558   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1559   TRACE_smpi_computing_out(rank);
1560   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1561   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1562 #endif
1563   if (comm == MPI_COMM_NULL) {
1564     retval = MPI_ERR_COMM;
1565   } else if (sendtype == MPI_DATATYPE_NULL
1566              || recvtype == MPI_DATATYPE_NULL) {
1567     retval = MPI_ERR_TYPE;
1568   } else if (sendcounts == NULL || displs == NULL) {
1569     retval = MPI_ERR_ARG;
1570   } else {
1571     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1572                       recvcount, recvtype, root, comm);
1573     retval = MPI_SUCCESS;
1574   }
1575 #ifdef HAVE_TRACING
1576   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1577   TRACE_smpi_computing_in(rank);
1578 #endif
1579   smpi_bench_begin();
1580   return retval;
1581 }
1582
1583 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1584                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1585 {
1586   int retval;
1587
1588   smpi_bench_end();
1589 #ifdef HAVE_TRACING
1590   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1591   TRACE_smpi_computing_out(rank);
1592   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1593   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1594 #endif
1595   if (comm == MPI_COMM_NULL) {
1596     retval = MPI_ERR_COMM;
1597   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1598     retval = MPI_ERR_ARG;
1599   } else {
1600     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1601     retval = MPI_SUCCESS;
1602   }
1603 #ifdef HAVE_TRACING
1604   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1605   TRACE_smpi_computing_in(rank);
1606 #endif
1607   smpi_bench_begin();
1608   return retval;
1609 }
1610
1611 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1612                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1613 {
1614   int retval;
1615
1616   smpi_bench_end();
1617 #ifdef HAVE_TRACING
1618   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1619   TRACE_smpi_computing_out(rank);
1620   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1621 #endif
1622   if (comm == MPI_COMM_NULL) {
1623     retval = MPI_ERR_COMM;
1624   } else if (datatype == MPI_DATATYPE_NULL) {
1625     retval = MPI_ERR_TYPE;
1626   } else if (op == MPI_OP_NULL) {
1627     retval = MPI_ERR_OP;
1628   } else {
1629     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
1630     retval = MPI_SUCCESS;
1631   }
1632 #ifdef HAVE_TRACING
1633   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1634   TRACE_smpi_computing_in(rank);
1635 #endif
1636   smpi_bench_begin();
1637   return retval;
1638 }
1639
1640 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
1641              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1642 {
1643   int retval;
1644
1645   smpi_bench_end();
1646 #ifdef HAVE_TRACING
1647   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1648   TRACE_smpi_computing_out(rank);
1649   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1650 #endif
1651   if (comm == MPI_COMM_NULL) {
1652     retval = MPI_ERR_COMM;
1653   } else if (datatype == MPI_DATATYPE_NULL) {
1654     retval = MPI_ERR_TYPE;
1655   } else if (op == MPI_OP_NULL) {
1656     retval = MPI_ERR_OP;
1657   } else {
1658     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1659     retval = MPI_SUCCESS;
1660   }
1661 #ifdef HAVE_TRACING
1662   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1663   TRACE_smpi_computing_in(rank);
1664 #endif
1665   smpi_bench_begin();
1666   return retval;
1667 }
1668
1669 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1670                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1671 {
1672   int retval, i, size, count;
1673   int *displs;
1674   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1675
1676   smpi_bench_end();
1677 #ifdef HAVE_TRACING
1678   TRACE_smpi_computing_out(rank);
1679   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1680 #endif
1681   if (comm == MPI_COMM_NULL) {
1682     retval = MPI_ERR_COMM;
1683   } else if (datatype == MPI_DATATYPE_NULL) {
1684     retval = MPI_ERR_TYPE;
1685   } else if (op == MPI_OP_NULL) {
1686     retval = MPI_ERR_OP;
1687   } else if (recvcounts == NULL) {
1688     retval = MPI_ERR_ARG;
1689   } else {
1690     /* arbitrarily choose root as rank 0 */
1691     /* TODO: faster direct implementation ? */
1692     size = smpi_comm_size(comm);
1693     count = 0;
1694     displs = xbt_new(int, size);
1695     for (i = 0; i < size; i++) {
1696       count += recvcounts[i];
1697       displs[i] = 0;
1698     }
1699     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1700     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf,
1701                       recvcounts[rank], datatype, 0, comm);
1702     xbt_free(displs);
1703     retval = MPI_SUCCESS;
1704   }
1705 #ifdef HAVE_TRACING
1706   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1707   TRACE_smpi_computing_in(rank);
1708 #endif
1709   smpi_bench_begin();
1710   return retval;
1711 }
1712
1713 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1714                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
1715                  MPI_Comm comm)
1716 {
1717   int retval, size, sendsize;
1718
1719   smpi_bench_end();
1720 #ifdef HAVE_TRACING
1721   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1722   TRACE_smpi_computing_out(rank);
1723   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1724 #endif
1725   if (comm == MPI_COMM_NULL) {
1726     retval = MPI_ERR_COMM;
1727   } else if (sendtype == MPI_DATATYPE_NULL
1728              || recvtype == MPI_DATATYPE_NULL) {
1729     retval = MPI_ERR_TYPE;
1730   } else {
1731     size = smpi_comm_size(comm);
1732     sendsize = smpi_datatype_size(sendtype) * sendcount;
1733     if (sendsize < 200 && size > 12) {
1734       retval =
1735           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
1736                                          recvbuf, recvcount, recvtype,
1737                                          comm);
1738     } else if (sendsize < 3000) {
1739       retval =
1740           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
1741                                                 sendtype, recvbuf,
1742                                                 recvcount, recvtype, comm);
1743     } else {
1744       retval =
1745           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
1746                                             recvbuf, recvcount, recvtype,
1747                                             comm);
1748     }
1749   }
1750 #ifdef HAVE_TRACING
1751   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1752   TRACE_smpi_computing_in(rank);
1753 #endif
1754   smpi_bench_begin();
1755   return retval;
1756 }
1757
1758 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
1759                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
1760                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
1761 {
1762   int retval;
1763
1764   smpi_bench_end();
1765 #ifdef HAVE_TRACING
1766   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1767   TRACE_smpi_computing_out(rank);
1768   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1769 #endif
1770   if (comm == MPI_COMM_NULL) {
1771     retval = MPI_ERR_COMM;
1772   } else if (sendtype == MPI_DATATYPE_NULL
1773              || recvtype == MPI_DATATYPE_NULL) {
1774     retval = MPI_ERR_TYPE;
1775   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
1776              || recvdisps == NULL) {
1777     retval = MPI_ERR_ARG;
1778   } else {
1779     retval =
1780         smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
1781                                   recvbuf, recvcounts, recvdisps, recvtype,
1782                                   comm);
1783   }
1784 #ifdef HAVE_TRACING
1785   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1786   TRACE_smpi_computing_in(rank);
1787 #endif
1788   smpi_bench_begin();
1789   return retval;
1790 }
1791
1792
1793 int PMPI_Get_processor_name(char *name, int *resultlen)
1794 {
1795   int retval = MPI_SUCCESS;
1796
1797   smpi_bench_end();
1798   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
1799           MPI_MAX_PROCESSOR_NAME - 1);
1800   *resultlen =
1801       strlen(name) >
1802       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1803
1804   smpi_bench_begin();
1805   return retval;
1806 }
1807
1808 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
1809 {
1810   int retval = MPI_SUCCESS;
1811   size_t size;
1812
1813   smpi_bench_end();
1814   if (status == NULL || count == NULL) {
1815     retval = MPI_ERR_ARG;
1816   } else if (datatype == MPI_DATATYPE_NULL) {
1817     retval = MPI_ERR_TYPE;
1818   } else {
1819     size = smpi_datatype_size(datatype);
1820     if (size == 0) {
1821       *count = 0;
1822     } else if (status->count % size != 0) {
1823       retval = MPI_UNDEFINED;
1824     } else {
1825       *count = smpi_mpi_get_count(status, datatype);
1826     }
1827   }
1828   smpi_bench_begin();
1829   return retval;
1830 }
1831
1832 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
1833   int retval;
1834
1835   smpi_bench_end();
1836   if (old_type == MPI_DATATYPE_NULL) {
1837     retval = MPI_ERR_TYPE;
1838   } else if (count<=0){
1839     retval = MPI_ERR_COUNT;
1840   } else {
1841     retval = smpi_datatype_contiguous(count, old_type, new_type);
1842   }
1843   smpi_bench_begin();
1844   return retval;
1845 }
1846
1847 int PMPI_Type_commit(MPI_Datatype* datatype) {
1848   int retval;
1849
1850   smpi_bench_end();
1851   if (datatype == MPI_DATATYPE_NULL) {
1852     retval = MPI_ERR_TYPE;
1853   } else {
1854     smpi_datatype_commit(datatype);
1855     retval = MPI_SUCCESS;
1856   }
1857   smpi_bench_begin();
1858   return retval;
1859 }
1860
1861
1862 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
1863   int retval;
1864
1865   smpi_bench_end();
1866   if (old_type == MPI_DATATYPE_NULL) {
1867     retval = MPI_ERR_TYPE;
1868   } else if (count<=0 || blocklen<=0){
1869     retval = MPI_ERR_COUNT;
1870   } else {
1871     retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
1872   }
1873   smpi_bench_begin();
1874   return retval;
1875 }
1876
1877 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
1878   int retval;
1879
1880   smpi_bench_end();
1881   if (old_type == MPI_DATATYPE_NULL) {
1882     retval = MPI_ERR_TYPE;
1883   } else if (count<=0 || blocklen<=0){
1884     retval = MPI_ERR_COUNT;
1885   } else {
1886     retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
1887   }
1888   smpi_bench_begin();
1889   return retval;
1890 }
1891
1892
1893 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
1894   int retval;
1895
1896   smpi_bench_end();
1897   if (old_type == MPI_DATATYPE_NULL) {
1898     retval = MPI_ERR_TYPE;
1899   } else if (count<=0){
1900     retval = MPI_ERR_COUNT;
1901   } else {
1902     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
1903   }
1904   smpi_bench_begin();
1905   return retval;
1906 }
1907
1908 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
1909   int retval;
1910
1911   smpi_bench_end();
1912   if (old_type == MPI_DATATYPE_NULL) {
1913     retval = MPI_ERR_TYPE;
1914   } else if (count<=0){
1915     retval = MPI_ERR_COUNT;
1916   } else {
1917     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
1918   }
1919   smpi_bench_begin();
1920   return retval;
1921 }
1922
1923
1924 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
1925   int retval;
1926
1927   smpi_bench_end();
1928   if (count<=0){
1929     retval = MPI_ERR_COUNT;
1930   } else {
1931     retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
1932   }
1933   smpi_bench_begin();
1934   return retval;}
1935
1936
1937 /* The following calls are not yet implemented and will fail at runtime. */
1938 /* Once implemented, please move them above this notice. */
1939
1940 static int not_yet_implemented(void) {
1941           XBT_WARN("Not yet implemented");
1942    return MPI_SUCCESS;
1943 }
1944
1945 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
1946    return not_yet_implemented();
1947 }
1948
1949 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
1950    return not_yet_implemented();
1951 }
1952
1953 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
1954    return not_yet_implemented();
1955 }
1956
1957 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
1958    return not_yet_implemented();
1959 }
1960
1961 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
1962    return not_yet_implemented();
1963 }
1964
1965 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
1966    return not_yet_implemented();
1967 }
1968
1969 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
1970    return not_yet_implemented();
1971 }
1972
1973 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
1974    return not_yet_implemented();
1975 }
1976
1977 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
1978    return not_yet_implemented();
1979 }
1980
1981 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
1982    return not_yet_implemented();
1983 }
1984
1985 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
1986    return not_yet_implemented();
1987 }
1988
1989 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
1990    return not_yet_implemented();
1991 }
1992
1993 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
1994    return not_yet_implemented();
1995 }
1996
1997 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
1998    return not_yet_implemented();
1999 }
2000
2001 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2002    return not_yet_implemented();
2003 }
2004
2005 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2006    return not_yet_implemented();
2007 }
2008
2009 int PMPI_Error_class(int errorcode, int* errorclass) {
2010    return not_yet_implemented();
2011 }
2012
2013 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2014    return not_yet_implemented();
2015 }
2016
2017 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2018    return not_yet_implemented();
2019 }
2020
2021 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2022    return not_yet_implemented();
2023 }
2024
2025 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2026    return not_yet_implemented();
2027 }
2028
2029 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2030    return not_yet_implemented();
2031 }
2032
2033
2034 int PMPI_Cancel(MPI_Request* request) {
2035    return not_yet_implemented();
2036 }
2037
2038 int PMPI_Buffer_attach(void* buffer, int size) {
2039    return not_yet_implemented();
2040 }
2041
2042 int PMPI_Buffer_detach(void* buffer, int* size) {
2043    return not_yet_implemented();
2044 }
2045
2046 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2047    return not_yet_implemented();
2048 }
2049
2050 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2051 {
2052    return not_yet_implemented();
2053 }
2054
2055 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2056    return not_yet_implemented();
2057 }
2058
2059 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2060    return not_yet_implemented();
2061 }
2062
2063 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2064    return not_yet_implemented();
2065 }
2066
2067 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2068    return not_yet_implemented();
2069 }
2070
2071 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2072    return not_yet_implemented();
2073 }
2074
2075 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2076    return not_yet_implemented();
2077 }
2078
2079 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2080    return not_yet_implemented();
2081 }
2082
2083 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2084    return not_yet_implemented();
2085 }
2086
2087 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2088    return not_yet_implemented();
2089 }
2090
2091 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2092    return not_yet_implemented();
2093 }
2094
2095 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2096    return not_yet_implemented();
2097 }
2098
2099
2100 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2101    return not_yet_implemented();
2102 }
2103
2104 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2105    return not_yet_implemented();
2106 }
2107
2108 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2109    return not_yet_implemented();
2110 }
2111
2112 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2113    return not_yet_implemented();
2114 }
2115
2116 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2117    return not_yet_implemented();
2118 }
2119
2120 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2121    return not_yet_implemented();
2122 }
2123
2124 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2125    return not_yet_implemented();
2126 }
2127
2128 int PMPI_Keyval_free(int* keyval) {
2129    return not_yet_implemented();
2130 }
2131
2132 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2133    return not_yet_implemented();
2134 }
2135
2136 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2137    return not_yet_implemented();
2138 }
2139
2140 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2141    return not_yet_implemented();
2142 }
2143
2144 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2145    return not_yet_implemented();
2146 }
2147
2148 int PMPI_Initialized(int* flag) {
2149    return not_yet_implemented();
2150 }