Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
8765cf70f20b632d141b425642659f700f91b0ee
[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_Finalized(int* flag)
54 {
55   *flag=smpi_process_finalized();
56   return MPI_SUCCESS;
57 }
58
59 int PMPI_Get_version (int *version,int *subversion){
60   *version = MPI_VERSION;
61   *subversion= MPI_SUBVERSION;
62   return MPI_SUCCESS;
63 }
64
65 int PMPI_Get_library_version (char *version,int *len){
66   int retval = MPI_SUCCESS;
67   smpi_bench_end();
68   snprintf(version,MPI_MAX_LIBRARY_VERSION_STRING,"SMPI Version %d.%d. Copyright The Simgrid Team 2007-2013",SIMGRID_VERSION_MAJOR,
69           SIMGRID_VERSION_MINOR);
70   *len = strlen(version) > MPI_MAX_LIBRARY_VERSION_STRING ? MPI_MAX_LIBRARY_VERSION_STRING : strlen(version);
71   smpi_bench_begin();
72   return retval;
73 }
74
75 int PMPI_Init_thread(int *argc, char ***argv, int required, int *provided)
76 {
77   if (provided != NULL) {
78     *provided = MPI_THREAD_MULTIPLE;
79   }
80   return MPI_Init(argc, argv);
81 }
82
83 int PMPI_Query_thread(int *provided)
84 {
85   int retval;
86
87   smpi_bench_end();
88   if (provided == NULL) {
89     retval = MPI_ERR_ARG;
90   } else {
91     *provided = MPI_THREAD_MULTIPLE;
92     retval = MPI_SUCCESS;
93   }
94   smpi_bench_begin();
95   return retval;
96 }
97
98 int PMPI_Is_thread_main(int *flag)
99 {
100   int retval;
101
102   smpi_bench_end();
103   if (flag == NULL) {
104     retval = MPI_ERR_ARG;
105   } else {
106     *flag = smpi_process_index() == 0;
107     retval = MPI_SUCCESS;
108   }
109   smpi_bench_begin();
110   return retval;
111 }
112
113 int PMPI_Abort(MPI_Comm comm, int errorcode)
114 {
115   smpi_bench_end();
116   smpi_process_destroy();
117 #ifdef HAVE_TRACING
118   int rank = smpi_process_index();
119   TRACE_smpi_computing_out(rank);
120 #endif
121   // FIXME: should kill all processes in comm instead
122   simcall_process_kill(SIMIX_process_self());
123   return MPI_SUCCESS;
124 }
125
126 double PMPI_Wtime(void)
127 {
128   double time;
129
130   smpi_bench_end();
131   time = SIMIX_get_clock();
132   smpi_bench_begin();
133   return time;
134 }
135 extern double sg_maxmin_precision;
136 double PMPI_Wtick(void)
137 {
138   return sg_maxmin_precision;
139 }
140
141 int PMPI_Address(void *location, MPI_Aint * address)
142 {
143   int retval;
144
145   smpi_bench_end();
146   if (!address) {
147     retval = MPI_ERR_ARG;
148   } else {
149     *address = (MPI_Aint) location;
150     retval = MPI_SUCCESS;
151   }
152   smpi_bench_begin();
153   return retval;
154 }
155
156 int PMPI_Get_address(void *location, MPI_Aint * address)
157 {
158   return PMPI_Address(location, address);
159 }
160
161 int PMPI_Type_free(MPI_Datatype * datatype)
162 {
163   int retval;
164
165   smpi_bench_end();
166   if (!datatype) {
167     retval = MPI_ERR_ARG;
168   } else {
169     smpi_datatype_free(datatype);
170     retval = MPI_SUCCESS;
171   }
172   smpi_bench_begin();
173   return retval;
174 }
175
176 int PMPI_Type_size(MPI_Datatype datatype, int *size)
177 {
178   int retval;
179
180   smpi_bench_end();
181   if (datatype == MPI_DATATYPE_NULL) {
182     retval = MPI_ERR_TYPE;
183   } else if (size == NULL) {
184     retval = MPI_ERR_ARG;
185   } else {
186     *size = (int) smpi_datatype_size(datatype);
187     retval = MPI_SUCCESS;
188   }
189   smpi_bench_begin();
190   return retval;
191 }
192
193 int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
194 {
195   int retval;
196
197   smpi_bench_end();
198   if (datatype == MPI_DATATYPE_NULL) {
199     retval = MPI_ERR_TYPE;
200   } else if (lb == NULL || extent == NULL) {
201     retval = MPI_ERR_ARG;
202   } else {
203     retval = smpi_datatype_extent(datatype, lb, extent);
204   }
205   smpi_bench_begin();
206   return retval;
207 }
208
209 int PMPI_Type_get_true_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
210 {
211   return PMPI_Type_get_extent(datatype, lb, extent);
212 }
213
214 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
215 {
216   int retval;
217
218   smpi_bench_end();
219   if (datatype == MPI_DATATYPE_NULL) {
220     retval = MPI_ERR_TYPE;
221   } else if (extent == NULL) {
222     retval = MPI_ERR_ARG;
223   } else {
224     *extent = smpi_datatype_get_extent(datatype);
225     retval = MPI_SUCCESS;
226   }
227   smpi_bench_begin();
228   return retval;
229 }
230
231 int PMPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
232 {
233   int retval;
234
235   smpi_bench_end();
236   if (datatype == MPI_DATATYPE_NULL) {
237     retval = MPI_ERR_TYPE;
238   } else if (disp == NULL) {
239     retval = MPI_ERR_ARG;
240   } else {
241     *disp = smpi_datatype_lb(datatype);
242     retval = MPI_SUCCESS;
243   }
244   smpi_bench_begin();
245   return retval;
246 }
247
248 int PMPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
249 {
250   int retval;
251
252   smpi_bench_end();
253   if (datatype == MPI_DATATYPE_NULL) {
254     retval = MPI_ERR_TYPE;
255   } else if (disp == NULL) {
256     retval = MPI_ERR_ARG;
257   } else {
258     *disp = smpi_datatype_ub(datatype);
259     retval = MPI_SUCCESS;
260   }
261   smpi_bench_begin();
262   return retval;
263 }
264
265 int PMPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
266 {
267   int retval;
268
269   smpi_bench_end();
270   if (function == NULL || op == NULL) {
271     retval = MPI_ERR_ARG;
272   } else {
273     *op = smpi_op_new(function, commute);
274     retval = MPI_SUCCESS;
275   }
276   smpi_bench_begin();
277   return retval;
278 }
279
280 int PMPI_Op_free(MPI_Op * op)
281 {
282   int retval;
283
284   smpi_bench_end();
285   if (op == NULL) {
286     retval = MPI_ERR_ARG;
287   } else if (*op == MPI_OP_NULL) {
288     retval = MPI_ERR_OP;
289   } else {
290     smpi_op_destroy(*op);
291     *op = MPI_OP_NULL;
292     retval = MPI_SUCCESS;
293   }
294   smpi_bench_begin();
295   return retval;
296 }
297
298 int PMPI_Group_free(MPI_Group * group)
299 {
300   int retval;
301
302   smpi_bench_end();
303   if (group == NULL) {
304     retval = MPI_ERR_ARG;
305   } else {
306     if(*group!= smpi_comm_group(MPI_COMM_WORLD))// do not free the group of the comm_world
307     smpi_group_destroy(*group);
308     *group = MPI_GROUP_NULL;
309     retval = MPI_SUCCESS;
310   }
311   smpi_bench_begin();
312   return retval;
313 }
314
315 int PMPI_Group_size(MPI_Group group, int *size)
316 {
317   int retval;
318
319   smpi_bench_end();
320   if (group == MPI_GROUP_NULL) {
321     retval = MPI_ERR_GROUP;
322   } else if (size == NULL) {
323     retval = MPI_ERR_ARG;
324   } else {
325     *size = smpi_group_size(group);
326     retval = MPI_SUCCESS;
327   }
328   smpi_bench_begin();
329   return retval;
330 }
331
332 int PMPI_Group_rank(MPI_Group group, int *rank)
333 {
334   int retval;
335
336   smpi_bench_end();
337   if (group == MPI_GROUP_NULL) {
338     retval = MPI_ERR_GROUP;
339   } else if (rank == NULL) {
340     retval = MPI_ERR_ARG;
341   } else {
342     *rank = smpi_group_rank(group, smpi_process_index());
343     retval = MPI_SUCCESS;
344   }
345   smpi_bench_begin();
346   return retval;
347 }
348
349 int PMPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1,
350                               MPI_Group group2, int *ranks2)
351 {
352   int retval, i, index;
353
354   smpi_bench_end();
355   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
356     retval = MPI_ERR_GROUP;
357   } else {
358     for (i = 0; i < n; i++) {
359       index = smpi_group_index(group1, ranks1[i]);
360       ranks2[i] = smpi_group_rank(group2, index);
361     }
362     retval = MPI_SUCCESS;
363   }
364   smpi_bench_begin();
365   return retval;
366 }
367
368 int PMPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
369 {
370   int retval;
371
372   smpi_bench_end();
373   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
374     retval = MPI_ERR_GROUP;
375   } else if (result == NULL) {
376     retval = MPI_ERR_ARG;
377   } else {
378     *result = smpi_group_compare(group1, group2);
379     retval = MPI_SUCCESS;
380   }
381   smpi_bench_begin();
382   return retval;
383 }
384
385 int PMPI_Group_union(MPI_Group group1, MPI_Group group2,
386                     MPI_Group * newgroup)
387 {
388   int retval, i, proc1, proc2, size, size2;
389
390   smpi_bench_end();
391   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
392     retval = MPI_ERR_GROUP;
393   } else if (newgroup == NULL) {
394     retval = MPI_ERR_ARG;
395   } else {
396     size = smpi_group_size(group1);
397     size2 = smpi_group_size(group2);
398     for (i = 0; i < size2; i++) {
399       proc2 = smpi_group_index(group2, i);
400       proc1 = smpi_group_rank(group1, proc2);
401       if (proc1 == MPI_UNDEFINED) {
402         size++;
403       }
404     }
405     if (size == 0) {
406       *newgroup = MPI_GROUP_EMPTY;
407     } else {
408       *newgroup = smpi_group_new(size);
409       size2 = smpi_group_size(group1);
410       for (i = 0; i < size2; i++) {
411         proc1 = smpi_group_index(group1, i);
412         smpi_group_set_mapping(*newgroup, proc1, i);
413       }
414       for (i = size2; i < size; i++) {
415         proc2 = smpi_group_index(group2, i - size2);
416         smpi_group_set_mapping(*newgroup, proc2, i);
417       }
418     }
419     smpi_group_use(*newgroup);
420     retval = MPI_SUCCESS;
421   }
422   smpi_bench_begin();
423   return retval;
424 }
425
426 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
427                            MPI_Group * newgroup)
428 {
429   int retval, i, proc1, proc2, size;
430
431   smpi_bench_end();
432   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
433     retval = MPI_ERR_GROUP;
434   } else if (newgroup == NULL) {
435     retval = MPI_ERR_ARG;
436   } else {
437     size = smpi_group_size(group2);
438     for (i = 0; i < size; i++) {
439       proc2 = smpi_group_index(group2, i);
440       proc1 = smpi_group_rank(group1, proc2);
441       if (proc1 == MPI_UNDEFINED) {
442         size--;
443       }
444     }
445     if (size == 0) {
446       *newgroup = MPI_GROUP_EMPTY;
447     } else {
448       *newgroup = smpi_group_new(size);
449       int j=0;
450       for (i = 0; i < smpi_group_size(group2); i++) {
451         proc2 = smpi_group_index(group2, i);
452         proc1 = smpi_group_rank(group1, proc2);
453         if (proc1 != MPI_UNDEFINED) {
454           smpi_group_set_mapping(*newgroup, proc2, j);
455           j++;
456         }
457       }
458     }
459     smpi_group_use(*newgroup);
460     retval = MPI_SUCCESS;
461   }
462   smpi_bench_begin();
463   return retval;
464 }
465
466 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
467 {
468   int retval, i, proc1, proc2, size, size2;
469
470   smpi_bench_end();
471   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
472     retval = MPI_ERR_GROUP;
473   } else if (newgroup == NULL) {
474     retval = MPI_ERR_ARG;
475   } else {
476     size = size2 = smpi_group_size(group1);
477     for (i = 0; i < size2; i++) {
478       proc1 = smpi_group_index(group1, i);
479       proc2 = smpi_group_rank(group2, proc1);
480       if (proc2 != MPI_UNDEFINED) {
481         size--;
482       }
483     }
484     if (size == 0) {
485       *newgroup = MPI_GROUP_EMPTY;
486     } else {
487       *newgroup = smpi_group_new(size);
488       for (i = 0; i < size2; i++) {
489         proc1 = smpi_group_index(group1, i);
490         proc2 = smpi_group_rank(group2, proc1);
491         if (proc2 == MPI_UNDEFINED) {
492           smpi_group_set_mapping(*newgroup, proc1, i);
493         }
494       }
495     }
496     smpi_group_use(*newgroup);
497     retval = MPI_SUCCESS;
498   }
499   smpi_bench_begin();
500   return retval;
501 }
502
503 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
504 {
505   int retval, i, index;
506
507   smpi_bench_end();
508   if (group == MPI_GROUP_NULL) {
509     retval = MPI_ERR_GROUP;
510   } else if (newgroup == NULL) {
511     retval = MPI_ERR_ARG;
512   } else {
513     if (n == 0) {
514       *newgroup = MPI_GROUP_EMPTY;
515     } else if (n == smpi_group_size(group)) {
516       *newgroup = group;
517     } else {
518       *newgroup = smpi_group_new(n);
519       for (i = 0; i < n; i++) {
520         index = smpi_group_index(group, ranks[i]);
521         smpi_group_set_mapping(*newgroup, index, i);
522       }
523     }
524     smpi_group_use(*newgroup);
525     retval = MPI_SUCCESS;
526   }
527   smpi_bench_begin();
528   return retval;
529 }
530
531 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
532 {
533   int retval, i, j, newsize, oldsize, index;
534
535   smpi_bench_end();
536   if (group == MPI_GROUP_NULL) {
537     retval = MPI_ERR_GROUP;
538   } else if (newgroup == NULL) {
539     retval = MPI_ERR_ARG;
540   } else {
541     if (n == 0) {
542       *newgroup = group;
543     } else if (n == smpi_group_size(group)) {
544       *newgroup = MPI_GROUP_EMPTY;
545     } else {
546       oldsize=smpi_group_size(group);
547       newsize = oldsize - n;
548       *newgroup = smpi_group_new(newsize);
549
550       int* to_exclude=xbt_new(int, smpi_group_size(group));
551       for(i=0; i<oldsize; i++)
552         to_exclude[i]=0;
553       for(i=0; i<n; i++)
554         to_exclude[ranks[i]]=1;
555
556       j=0;
557       for(i=0; i<oldsize; i++){
558         if(to_exclude[i]==0){
559           index = smpi_group_index(group, i);
560           smpi_group_set_mapping(*newgroup, index, j);
561           j++;
562         }
563       }
564
565       xbt_free(to_exclude);
566     }
567     smpi_group_use(*newgroup);
568     retval = MPI_SUCCESS;
569   }
570   smpi_bench_begin();
571   return retval;
572 }
573
574 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
575                          MPI_Group * newgroup)
576 {
577   int retval, i, j, rank, size, index;
578
579   smpi_bench_end();
580   if (group == MPI_GROUP_NULL) {
581     retval = MPI_ERR_GROUP;
582   } else if (newgroup == NULL) {
583     retval = MPI_ERR_ARG;
584   } else {
585     if (n == 0) {
586       *newgroup = MPI_GROUP_EMPTY;
587     } else {
588       size = 0;
589       for (i = 0; i < n; i++) {
590         for (rank = ranges[i][0];       /* First */
591              rank >= 0 && rank <= ranges[i][1]; /* Last */
592              rank += ranges[i][2] /* Stride */ ) {
593           size++;
594         }
595       }
596       if (size == smpi_group_size(group)) {
597         *newgroup = group;
598       } else {
599         *newgroup = smpi_group_new(size);
600         j = 0;
601         for (i = 0; i < n; i++) {
602           for (rank = ranges[i][0];     /* First */
603                rank >= 0 && rank <= ranges[i][1];       /* Last */
604                rank += ranges[i][2] /* Stride */ ) {
605             index = smpi_group_index(group, rank);
606             smpi_group_set_mapping(*newgroup, index, j);
607             j++;
608           }
609         }
610       }
611     }
612     smpi_group_use(*newgroup);
613     retval = MPI_SUCCESS;
614   }
615   smpi_bench_begin();
616   return retval;
617 }
618
619 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
620                          MPI_Group * newgroup)
621 {
622   int retval, i, newrank, rank, size, index, add;
623
624   smpi_bench_end();
625   if (group == MPI_GROUP_NULL) {
626     retval = MPI_ERR_GROUP;
627   } else if (newgroup == NULL) {
628     retval = MPI_ERR_ARG;
629   } else {
630     if (n == 0) {
631       *newgroup = group;
632     } else {
633       size = smpi_group_size(group);
634       for (i = 0; i < n; i++) {
635         for (rank = ranges[i][0];       /* First */
636              rank >= 0 && rank <= ranges[i][1]; /* Last */
637              rank += ranges[i][2] /* Stride */ ) {
638           size--;
639         }
640       }
641       if (size == 0) {
642         *newgroup = MPI_GROUP_EMPTY;
643       } else {
644         *newgroup = smpi_group_new(size);
645         newrank = 0;
646         while (newrank < size) {
647           for (i = 0; i < n; i++) {
648             add = 1;
649             for (rank = ranges[i][0];   /* First */
650                  rank >= 0 && rank <= ranges[i][1];     /* Last */
651                  rank += ranges[i][2] /* Stride */ ) {
652               if (rank == newrank) {
653                 add = 0;
654                 break;
655               }
656             }
657             if (add == 1) {
658               index = smpi_group_index(group, newrank);
659               smpi_group_set_mapping(*newgroup, index, newrank);
660             }
661           }
662           newrank++; //added to avoid looping, need to be checked ..
663         }
664       }
665     }
666     smpi_group_use(*newgroup);
667     retval = MPI_SUCCESS;
668   }
669   smpi_bench_begin();
670   return retval;
671 }
672
673 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
674 {
675   int retval;
676
677   smpi_bench_end();
678   if (comm == MPI_COMM_NULL) {
679     retval = MPI_ERR_COMM;
680   } else if (rank == NULL) {
681     retval = MPI_ERR_ARG;
682   } else {
683     *rank = smpi_comm_rank(comm);
684     retval = MPI_SUCCESS;
685   }
686   smpi_bench_begin();
687   return retval;
688 }
689
690 int PMPI_Comm_size(MPI_Comm comm, int *size)
691 {
692   int retval;
693
694   smpi_bench_end();
695   if (comm == MPI_COMM_NULL) {
696     retval = MPI_ERR_COMM;
697   } else if (size == NULL) {
698     retval = MPI_ERR_ARG;
699   } else {
700     *size = smpi_comm_size(comm);
701     retval = MPI_SUCCESS;
702   }
703   smpi_bench_begin();
704   return retval;
705 }
706
707 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
708 {
709   int retval;
710
711   smpi_bench_end();
712   if (comm == MPI_COMM_NULL)  {
713     retval = MPI_ERR_COMM;
714   } else if (name == NULL || len == NULL)  {
715     retval = MPI_ERR_ARG;
716   } else {
717     smpi_comm_get_name(comm, name, len);
718     retval = MPI_SUCCESS;
719   }
720   smpi_bench_begin();
721   return retval;
722 }
723
724 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
725 {
726   int retval;
727
728   smpi_bench_end();
729   if (comm == MPI_COMM_NULL) {
730     retval = MPI_ERR_COMM;
731   } else if (group == NULL) {
732     retval = MPI_ERR_ARG;
733   } else {
734     *group = smpi_comm_group(comm);
735     retval = MPI_SUCCESS;
736   }
737   smpi_bench_begin();
738   return retval;
739 }
740
741 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
742 {
743   int retval;
744
745   smpi_bench_end();
746   if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
747     retval = MPI_ERR_COMM;
748   } else if (result == NULL) {
749     retval = MPI_ERR_ARG;
750   } else {
751     if (comm1 == comm2) {       /* Same communicators means same groups */
752       *result = MPI_IDENT;
753     } else {
754       *result =
755           smpi_group_compare(smpi_comm_group(comm1),
756                              smpi_comm_group(comm2));
757       if (*result == MPI_IDENT) {
758         *result = MPI_CONGRUENT;
759       }
760     }
761     retval = MPI_SUCCESS;
762   }
763   smpi_bench_begin();
764   return retval;
765 }
766
767 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
768 {
769   int retval;
770
771   smpi_bench_end();
772   if (comm == MPI_COMM_NULL) {
773     retval = MPI_ERR_COMM;
774   } else if (newcomm == NULL) {
775     retval = MPI_ERR_ARG;
776   } else {
777     *newcomm = smpi_comm_new(smpi_comm_group(comm));
778     retval = MPI_SUCCESS;
779   }
780   smpi_bench_begin();
781   return retval;
782 }
783
784 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
785 {
786   int retval;
787
788   smpi_bench_end();
789   if (comm == MPI_COMM_NULL) {
790     retval = MPI_ERR_COMM;
791   } else if (group == MPI_GROUP_NULL) {
792     retval = MPI_ERR_GROUP;
793   } else if (newcomm == NULL) {
794     retval = MPI_ERR_ARG;
795   } else if(smpi_group_rank(group,smpi_process_index())==MPI_UNDEFINED){
796     *newcomm= MPI_COMM_NULL;
797     retval = MPI_SUCCESS;
798   }else{
799
800     *newcomm = smpi_comm_new(group);
801     retval = MPI_SUCCESS;
802   }
803   smpi_bench_begin();
804   return retval;
805 }
806
807 int PMPI_Comm_free(MPI_Comm * comm)
808 {
809   int retval;
810
811   smpi_bench_end();
812   if (comm == NULL) {
813     retval = MPI_ERR_ARG;
814   } else if (*comm == MPI_COMM_NULL) {
815     retval = MPI_ERR_COMM;
816   } else {
817     smpi_comm_destroy(*comm);
818     *comm = MPI_COMM_NULL;
819     retval = MPI_SUCCESS;
820   }
821   smpi_bench_begin();
822   return retval;
823 }
824
825 int PMPI_Comm_disconnect(MPI_Comm * comm)
826 {
827   /* TODO: wait until all communication in comm are done */
828   int retval;
829
830   smpi_bench_end();
831   if (comm == NULL) {
832     retval = MPI_ERR_ARG;
833   } else if (*comm == MPI_COMM_NULL) {
834     retval = MPI_ERR_COMM;
835   } else {
836     smpi_comm_destroy(*comm);
837     *comm = MPI_COMM_NULL;
838     retval = MPI_SUCCESS;
839   }
840   smpi_bench_begin();
841   return retval;
842 }
843
844 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
845 {
846   int retval;
847
848   smpi_bench_end();
849   if (comm_out == NULL) {
850     retval = MPI_ERR_ARG;
851   } else if (comm == MPI_COMM_NULL) {
852     retval = MPI_ERR_COMM;
853   } else {
854     *comm_out = smpi_comm_split(comm, color, key);
855     retval = MPI_SUCCESS;
856   }
857   smpi_bench_begin();
858   return retval;
859 }
860
861 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
862                   int tag, MPI_Comm comm, MPI_Request * request)
863 {
864   int retval;
865
866   smpi_bench_end();
867   if (request == NULL) {
868     retval = MPI_ERR_ARG;
869   } else if (comm == MPI_COMM_NULL) {
870     retval = MPI_ERR_COMM;
871   } else if (dst == MPI_PROC_NULL) {
872     retval = MPI_SUCCESS;
873   } else {
874     *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
875     retval = MPI_SUCCESS;
876   }
877   smpi_bench_begin();
878   return retval;
879 }
880
881 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
882                   int tag, MPI_Comm comm, MPI_Request * request)
883 {
884   int retval;
885
886   smpi_bench_end();
887   if (request == NULL) {
888     retval = MPI_ERR_ARG;
889   } else if (comm == MPI_COMM_NULL) {
890     retval = MPI_ERR_COMM;
891   } else if (src == MPI_PROC_NULL) {
892       retval = MPI_SUCCESS;
893   } else {
894     *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
895     retval = MPI_SUCCESS;
896   }
897   smpi_bench_begin();
898   return retval;
899 }
900
901 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
902   int retval;
903
904     smpi_bench_end();
905     if (request == NULL) {
906       retval = MPI_ERR_ARG;
907     } else if (comm == MPI_COMM_NULL) {
908       retval = MPI_ERR_COMM;
909     } else if (dst == MPI_PROC_NULL) {
910       retval = MPI_SUCCESS;
911     } else {
912       *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
913       retval = MPI_SUCCESS;
914     }
915     smpi_bench_begin();
916     return retval;
917 }
918
919 int PMPI_Start(MPI_Request * request)
920 {
921   int retval;
922
923   smpi_bench_end();
924   if (request == NULL || *request == MPI_REQUEST_NULL) {
925     retval = MPI_ERR_ARG;
926   } else {
927     smpi_mpi_start(*request);
928     retval = MPI_SUCCESS;
929   }
930   smpi_bench_begin();
931   return retval;
932 }
933
934 int PMPI_Startall(int count, MPI_Request * requests)
935 {
936   int retval;
937
938   smpi_bench_end();
939   if (requests == NULL) {
940     retval = MPI_ERR_ARG;
941   } else {
942     smpi_mpi_startall(count, requests);
943     retval = MPI_SUCCESS;
944   }
945   smpi_bench_begin();
946   return retval;
947 }
948
949 int PMPI_Request_free(MPI_Request * request)
950 {
951   int retval;
952
953   smpi_bench_end();
954   if (*request == MPI_REQUEST_NULL) {
955     retval = MPI_ERR_ARG;
956   } else {
957     if((*request)->flags & PERSISTENT)(*request)->refcount--;
958     smpi_mpi_request_free(request);
959     retval = MPI_SUCCESS;
960   }
961   smpi_bench_begin();
962   return retval;
963 }
964
965 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
966               int tag, MPI_Comm comm, MPI_Request * request)
967 {
968   int retval;
969
970   smpi_bench_end();
971
972   if (request == NULL) {
973     retval = MPI_ERR_ARG;
974   } else if (comm == MPI_COMM_NULL) {
975     retval = MPI_ERR_COMM;
976   } else if (src == MPI_PROC_NULL) {
977     *request = MPI_REQUEST_NULL;
978     retval = MPI_SUCCESS;
979   } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
980     retval = MPI_ERR_COMM;
981   } else if (count < 0) {
982     retval = MPI_ERR_COUNT;
983   } else if (buf==NULL && count > 0) {
984     retval = MPI_ERR_COUNT;
985   } else if (datatype == MPI_DATATYPE_NULL){
986     retval = MPI_ERR_TYPE;
987   } else if(tag<0 && tag !=  MPI_ANY_TAG){
988     retval = MPI_ERR_TAG;
989   } else {
990
991 #ifdef HAVE_TRACING
992   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
993   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
994   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
995 #endif
996
997     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
998     retval = MPI_SUCCESS;
999
1000 #ifdef HAVE_TRACING
1001   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1002   (*request)->recv = 1;
1003 #endif
1004   }
1005
1006   smpi_bench_begin();
1007   return retval;
1008 }
1009
1010
1011 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
1012               int tag, MPI_Comm comm, MPI_Request * request)
1013 {
1014   int retval;
1015
1016   smpi_bench_end();
1017   if (request == NULL) {
1018     retval = MPI_ERR_ARG;
1019   } else if (comm == MPI_COMM_NULL) {
1020     retval = MPI_ERR_COMM;
1021   } else if (dst == MPI_PROC_NULL) {
1022     *request = MPI_REQUEST_NULL;
1023     retval = MPI_SUCCESS;
1024   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1025     retval = MPI_ERR_RANK;
1026   } else if (count < 0) {
1027     retval = MPI_ERR_COUNT;
1028   } else if (buf==NULL && count > 0) {
1029     retval = MPI_ERR_COUNT;
1030   } else if (datatype == MPI_DATATYPE_NULL){
1031     retval = MPI_ERR_TYPE;
1032   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1033     retval = MPI_ERR_TAG;
1034   } else {
1035
1036 #ifdef HAVE_TRACING
1037   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1038   TRACE_smpi_computing_out(rank);
1039   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1040   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1041   TRACE_smpi_send(rank, rank, dst_traced);
1042 #endif
1043
1044     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1045     retval = MPI_SUCCESS;
1046
1047 #ifdef HAVE_TRACING
1048   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1049   (*request)->send = 1;
1050   TRACE_smpi_computing_in(rank);
1051 #endif
1052   }
1053
1054   smpi_bench_begin();
1055   return retval;
1056 }
1057
1058 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
1059   int retval;
1060
1061   smpi_bench_end();
1062   if (request == NULL) {
1063     retval = MPI_ERR_ARG;
1064   } else if (comm == MPI_COMM_NULL) {
1065     retval = MPI_ERR_COMM;
1066   } else if (dst == MPI_PROC_NULL) {
1067     *request = MPI_REQUEST_NULL;
1068     retval = MPI_SUCCESS;
1069   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1070     retval = MPI_ERR_RANK;
1071   } else if (count < 0) {
1072     retval = MPI_ERR_COUNT;
1073   } else if (buf==NULL && count > 0) {
1074     retval = MPI_ERR_COUNT;
1075   } else if (datatype == MPI_DATATYPE_NULL){
1076     retval = MPI_ERR_TYPE;
1077   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1078     retval = MPI_ERR_TAG;
1079   } else {
1080
1081 #ifdef HAVE_TRACING
1082   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1083   TRACE_smpi_computing_out(rank);
1084   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1085   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1086   TRACE_smpi_send(rank, rank, dst_traced);
1087 #endif
1088
1089     *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1090     retval = MPI_SUCCESS;
1091
1092 #ifdef HAVE_TRACING
1093   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1094   (*request)->send = 1;
1095   TRACE_smpi_computing_in(rank);
1096 #endif
1097   }
1098
1099   smpi_bench_begin();
1100   return retval;
1101 }
1102
1103 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1104              MPI_Comm comm, MPI_Status * status)
1105 {
1106   int retval;
1107
1108   smpi_bench_end();
1109   if (comm == MPI_COMM_NULL) {
1110     retval = MPI_ERR_COMM;
1111   } else if (src == MPI_PROC_NULL) {
1112     smpi_empty_status(status);
1113     status->MPI_SOURCE = MPI_PROC_NULL;
1114     retval = MPI_SUCCESS;
1115   } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1116     retval = MPI_ERR_RANK;
1117   } else if (count < 0) {
1118     retval = MPI_ERR_COUNT;
1119   } else if (buf==NULL && count > 0) {
1120     retval = MPI_ERR_COUNT;
1121   } else if (datatype == MPI_DATATYPE_NULL){
1122     retval = MPI_ERR_TYPE;
1123   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1124     retval = MPI_ERR_TAG;
1125   } else {
1126 #ifdef HAVE_TRACING
1127   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1128   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1129   TRACE_smpi_computing_out(rank);
1130
1131   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
1132 #endif
1133
1134     smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1135     retval = MPI_SUCCESS;
1136
1137 #ifdef HAVE_TRACING
1138   //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1139   if(status!=MPI_STATUS_IGNORE)src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1140   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1141   TRACE_smpi_recv(rank, src_traced, rank);
1142   TRACE_smpi_computing_in(rank);
1143 #endif
1144   }
1145
1146   smpi_bench_begin();
1147   return retval;
1148 }
1149
1150 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1151              MPI_Comm comm)
1152 {
1153   int retval;
1154
1155   smpi_bench_end();
1156
1157   if (comm == MPI_COMM_NULL) {
1158     retval = MPI_ERR_COMM;
1159   } else if (dst == MPI_PROC_NULL) {
1160     retval = MPI_SUCCESS;
1161   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1162     retval = MPI_ERR_RANK;
1163   } else if (count < 0) {
1164     retval = MPI_ERR_COUNT;
1165   } else if (buf==NULL && count > 0) {
1166     retval = MPI_ERR_COUNT;
1167   } else if (datatype == MPI_DATATYPE_NULL){
1168     retval = MPI_ERR_TYPE;
1169   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1170     retval = MPI_ERR_TAG;
1171   } else {
1172
1173 #ifdef HAVE_TRACING
1174   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1175   TRACE_smpi_computing_out(rank);
1176   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1177   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1178   TRACE_smpi_send(rank, rank, dst_traced);
1179 #endif
1180
1181     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1182     retval = MPI_SUCCESS;
1183
1184 #ifdef HAVE_TRACING
1185   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1186   TRACE_smpi_computing_in(rank);
1187 #endif
1188   }
1189
1190   smpi_bench_begin();
1191   return retval;
1192 }
1193
1194
1195
1196 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1197   int retval;
1198
1199    smpi_bench_end();
1200
1201    if (comm == MPI_COMM_NULL) {
1202      retval = MPI_ERR_COMM;
1203    } else if (dst == MPI_PROC_NULL) {
1204      retval = MPI_SUCCESS;
1205    } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1206      retval = MPI_ERR_RANK;
1207    } else if (count < 0) {
1208      retval = MPI_ERR_COUNT;
1209    } else if (buf==NULL && count > 0) {
1210      retval = MPI_ERR_COUNT;
1211    } else if (datatype == MPI_DATATYPE_NULL){
1212      retval = MPI_ERR_TYPE;
1213    } else if(tag<0 && tag !=  MPI_ANY_TAG){
1214      retval = MPI_ERR_TAG;
1215    } else {
1216
1217  #ifdef HAVE_TRACING
1218    int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1219    TRACE_smpi_computing_out(rank);
1220    int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1221    TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1222    TRACE_smpi_send(rank, rank, dst_traced);
1223  #endif
1224
1225      smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1226      retval = MPI_SUCCESS;
1227
1228  #ifdef HAVE_TRACING
1229    TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1230    TRACE_smpi_computing_in(rank);
1231  #endif
1232    }
1233
1234    smpi_bench_begin();
1235    return retval;}
1236
1237
1238 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1239                  int dst, int sendtag, void *recvbuf, int recvcount,
1240                  MPI_Datatype recvtype, int src, int recvtag,
1241                  MPI_Comm comm, MPI_Status * status)
1242 {
1243   int retval;
1244
1245   smpi_bench_end();
1246
1247   if (comm == MPI_COMM_NULL) {
1248     retval = MPI_ERR_COMM;
1249   } else if (sendtype == MPI_DATATYPE_NULL
1250              || recvtype == MPI_DATATYPE_NULL) {
1251     retval = MPI_ERR_TYPE;
1252   } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1253       smpi_empty_status(status);
1254       status->MPI_SOURCE = MPI_PROC_NULL;
1255       retval = MPI_SUCCESS;
1256   }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1257       (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1258     retval = MPI_ERR_RANK;
1259   } else if (sendcount < 0 || recvcount<0) {
1260       retval = MPI_ERR_COUNT;
1261   } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1262     retval = MPI_ERR_COUNT;
1263   } else if((sendtag<0 && sendtag !=  MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1264     retval = MPI_ERR_TAG;
1265   } else {
1266
1267 #ifdef HAVE_TRACING
1268   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1269   TRACE_smpi_computing_out(rank);
1270   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1271   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1272   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1273   TRACE_smpi_send(rank, rank, dst_traced);
1274 #endif
1275
1276
1277     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1278                       recvcount, recvtype, src, recvtag, comm, status);
1279     retval = MPI_SUCCESS;
1280
1281 #ifdef HAVE_TRACING
1282   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1283   TRACE_smpi_recv(rank, src_traced, rank);
1284   TRACE_smpi_computing_in(rank);
1285 #endif
1286
1287   }
1288
1289   smpi_bench_begin();
1290   return retval;
1291 }
1292
1293 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1294                          int dst, int sendtag, int src, int recvtag,
1295                          MPI_Comm comm, MPI_Status * status)
1296 {
1297   //TODO: suboptimal implementation
1298   void *recvbuf;
1299   int retval;
1300   if (datatype == MPI_DATATYPE_NULL) {
1301       retval = MPI_ERR_TYPE;
1302   } else if (count < 0) {
1303       retval = MPI_ERR_COUNT;
1304   } else {
1305     int size = smpi_datatype_get_extent(datatype) * count;
1306     recvbuf = xbt_new(char, size);
1307     retval =
1308         MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1309                      datatype, src, recvtag, comm, status);
1310     if(retval==MPI_SUCCESS){
1311         smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1312     }
1313     xbt_free(recvbuf);
1314
1315   }
1316   return retval;
1317 }
1318
1319 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1320 {
1321   int retval;
1322
1323   smpi_bench_end();
1324   if (request == MPI_REQUEST_NULL || flag == NULL) {
1325     retval = MPI_ERR_ARG;
1326   } else if (*request == MPI_REQUEST_NULL) {
1327     *flag= TRUE;
1328     retval = MPI_ERR_REQUEST;
1329   } else {
1330     *flag = smpi_mpi_test(request, status);
1331     retval = MPI_SUCCESS;
1332   }
1333   smpi_bench_begin();
1334   return retval;
1335 }
1336
1337 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1338                 MPI_Status * status)
1339 {
1340   int retval;
1341
1342   smpi_bench_end();
1343   if (index == NULL || flag == NULL) {
1344     retval = MPI_ERR_ARG;
1345   } else {
1346     *flag = smpi_mpi_testany(count, requests, index, status);
1347     retval = MPI_SUCCESS;
1348   }
1349   smpi_bench_begin();
1350   return retval;
1351 }
1352
1353 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1354 {
1355   int retval;
1356
1357   smpi_bench_end();
1358   if (flag == NULL) {
1359     retval = MPI_ERR_ARG;
1360   } else {
1361     *flag = smpi_mpi_testall(count, requests, statuses);
1362     retval = MPI_SUCCESS;
1363   }
1364   smpi_bench_begin();
1365   return retval;
1366 }
1367
1368 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1369   int retval;
1370   smpi_bench_end();
1371
1372   if (status == NULL) {
1373     retval = MPI_ERR_ARG;
1374   } else if (comm == MPI_COMM_NULL) {
1375     retval = MPI_ERR_COMM;
1376   } else if (source == MPI_PROC_NULL) {
1377     smpi_empty_status(status);
1378     status->MPI_SOURCE = MPI_PROC_NULL;
1379     retval = MPI_SUCCESS;
1380   } else {
1381     smpi_mpi_probe(source, tag, comm, status);
1382     retval = MPI_SUCCESS;
1383   }
1384   smpi_bench_begin();
1385   return retval;
1386 }
1387
1388
1389 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1390   int retval;
1391   smpi_bench_end();
1392
1393   if (flag == NULL) {
1394     retval = MPI_ERR_ARG;
1395   } else if (status == NULL) {
1396     retval = MPI_ERR_ARG;
1397   } else if (comm == MPI_COMM_NULL) {
1398     retval = MPI_ERR_COMM;
1399   } else if (source == MPI_PROC_NULL) {
1400     smpi_empty_status(status);
1401     status->MPI_SOURCE = MPI_PROC_NULL;
1402     retval = MPI_SUCCESS;
1403   } else {
1404     smpi_mpi_iprobe(source, tag, comm, flag, status);
1405     retval = MPI_SUCCESS;
1406   }
1407   smpi_bench_begin();
1408   return retval;
1409 }
1410
1411 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1412 {
1413   int retval;
1414
1415   smpi_bench_end();
1416
1417   if (request == NULL) {
1418     retval = MPI_ERR_ARG;
1419   } else if (*request == MPI_REQUEST_NULL) {
1420     retval = MPI_ERR_REQUEST;
1421   } else {
1422
1423 #ifdef HAVE_TRACING
1424   int rank = request && (*request)->comm != MPI_COMM_NULL
1425       ? smpi_process_index()
1426       : -1;
1427   TRACE_smpi_computing_out(rank);
1428
1429   int src_traced = (*request)->src;
1430   int dst_traced = (*request)->dst;
1431   int is_wait_for_receive = (*request)->recv;
1432   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1433 #endif
1434
1435     smpi_mpi_wait(request, status);
1436     retval = MPI_SUCCESS;
1437
1438 #ifdef HAVE_TRACING
1439   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1440   if (is_wait_for_receive) {
1441     TRACE_smpi_recv(rank, src_traced, dst_traced);
1442   }
1443   TRACE_smpi_computing_in(rank);
1444 #endif
1445
1446   }
1447
1448   smpi_bench_begin();
1449   return retval;
1450 }
1451
1452 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1453 {
1454   int retval;
1455
1456   smpi_bench_end();
1457 #ifdef HAVE_TRACING
1458   //save requests information for tracing
1459   int i;
1460   int *srcs = xbt_new(int, count);
1461   int *dsts = xbt_new(int, count);
1462   int *recvs = xbt_new(int, count);
1463   for (i = 0; i < count; i++) {
1464     MPI_Request req = requests[i];      //already received requests are no longer valid
1465     if (req) {
1466       srcs[i] = req->src;
1467       dsts[i] = req->dst;
1468       recvs[i] = req->recv;
1469     }
1470   }
1471   int rank_traced = smpi_process_index();
1472   TRACE_smpi_computing_out(rank_traced);
1473
1474   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1475
1476 #endif
1477   if (index == NULL) {
1478     retval = MPI_ERR_ARG;
1479   } else {
1480     *index = smpi_mpi_waitany(count, requests, status);
1481     retval = MPI_SUCCESS;
1482   }
1483 #ifdef HAVE_TRACING
1484   if(*index!=MPI_UNDEFINED){
1485     int src_traced = srcs[*index];
1486     int dst_traced = dsts[*index];
1487     int is_wait_for_receive = recvs[*index];
1488     if (is_wait_for_receive) {
1489       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1490     }
1491     TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1492     xbt_free(srcs);
1493     xbt_free(dsts);
1494     xbt_free(recvs);
1495   }
1496   TRACE_smpi_computing_in(rank_traced);
1497 #endif
1498   smpi_bench_begin();
1499   return retval;
1500 }
1501
1502 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1503 {
1504
1505   smpi_bench_end();
1506 #ifdef HAVE_TRACING
1507   //save information from requests
1508   int i;
1509   int *srcs = xbt_new(int, count);
1510   int *dsts = xbt_new(int, count);
1511   int *recvs = xbt_new(int, count);
1512   int valid_count = 0;
1513   for (i = 0; i < count; i++) {
1514     MPI_Request req = requests[i];
1515     if(req!=MPI_REQUEST_NULL){
1516       srcs[valid_count] = req->src;
1517       dsts[valid_count] = req->dst;
1518       recvs[valid_count] = req->recv;
1519       valid_count++;
1520     }
1521   }
1522   int rank_traced = smpi_process_index();
1523   TRACE_smpi_computing_out(rank_traced);
1524
1525   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1526 #endif
1527   int retval = smpi_mpi_waitall(count, requests, status);
1528 #ifdef HAVE_TRACING
1529   for (i = 0; i < valid_count; i++) {
1530     int src_traced = srcs[i];
1531     int dst_traced = dsts[i];
1532     int is_wait_for_receive = recvs[i];
1533     if (is_wait_for_receive) {
1534       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1535     }
1536   }
1537   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1538   xbt_free(srcs);
1539   xbt_free(dsts);
1540   xbt_free(recvs);
1541   TRACE_smpi_computing_in(rank_traced);
1542 #endif
1543   smpi_bench_begin();
1544   return retval;
1545 }
1546
1547 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1548                  int *indices, MPI_Status status[])
1549 {
1550   int retval;
1551
1552   smpi_bench_end();
1553   if (outcount == NULL || indices == NULL) {
1554     retval = MPI_ERR_ARG;
1555   } else {
1556     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1557     retval = MPI_SUCCESS;
1558   }
1559   smpi_bench_begin();
1560   return retval;
1561 }
1562
1563 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1564                  int* indices, MPI_Status status[])
1565 {
1566   int retval;
1567
1568    smpi_bench_end();
1569    if (outcount == NULL || indices == NULL) {
1570      retval = MPI_ERR_ARG;
1571    } else {
1572      *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1573      retval = MPI_SUCCESS;
1574    }
1575    smpi_bench_begin();
1576    return retval;
1577 }
1578
1579
1580 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1581 {
1582   int retval;
1583
1584   smpi_bench_end();
1585 #ifdef HAVE_TRACING
1586   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1587   TRACE_smpi_computing_out(rank);
1588   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1589   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1590 #endif
1591   if (comm == MPI_COMM_NULL) {
1592     retval = MPI_ERR_COMM;
1593   } else {
1594     mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1595     retval = MPI_SUCCESS;
1596   }
1597 #ifdef HAVE_TRACING
1598   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1599   TRACE_smpi_computing_in(rank);
1600 #endif
1601   smpi_bench_begin();
1602   return retval;
1603 }
1604
1605 int PMPI_Barrier(MPI_Comm comm)
1606 {
1607   int retval;
1608
1609   smpi_bench_end();
1610 #ifdef HAVE_TRACING
1611   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1612   TRACE_smpi_computing_out(rank);
1613   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1614 #endif
1615   if (comm == MPI_COMM_NULL) {
1616     retval = MPI_ERR_COMM;
1617   } else {
1618     mpi_coll_barrier_fun(comm);
1619     retval = MPI_SUCCESS;
1620   }
1621 #ifdef HAVE_TRACING
1622   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1623   TRACE_smpi_computing_in(rank);
1624 #endif
1625   smpi_bench_begin();
1626   return retval;
1627 }
1628
1629 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1630                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1631                int root, MPI_Comm comm)
1632 {
1633   int retval;
1634
1635   smpi_bench_end();
1636 #ifdef HAVE_TRACING
1637   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1638   TRACE_smpi_computing_out(rank);
1639   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1640   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1641 #endif
1642   if (comm == MPI_COMM_NULL) {
1643     retval = MPI_ERR_COMM;
1644   } else if (sendtype == MPI_DATATYPE_NULL
1645              || recvtype == MPI_DATATYPE_NULL) {
1646     retval = MPI_ERR_TYPE;
1647   } else {
1648     mpi_coll_gather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1649                     recvtype, root, comm);
1650     retval = MPI_SUCCESS;
1651   }
1652 #ifdef HAVE_TRACING
1653   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1654   TRACE_smpi_computing_in(rank);
1655 #endif
1656   smpi_bench_begin();
1657   return retval;
1658 }
1659
1660 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1661                 void *recvbuf, int *recvcounts, int *displs,
1662                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1663 {
1664   int retval;
1665
1666   smpi_bench_end();
1667 #ifdef HAVE_TRACING
1668   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1669   TRACE_smpi_computing_out(rank);
1670   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1671   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1672 #endif
1673   if (comm == MPI_COMM_NULL) {
1674     retval = MPI_ERR_COMM;
1675   } else if (sendtype == MPI_DATATYPE_NULL
1676              || recvtype == MPI_DATATYPE_NULL) {
1677     retval = MPI_ERR_TYPE;
1678   } else if (recvcounts == NULL || displs == NULL) {
1679     retval = MPI_ERR_ARG;
1680   } else {
1681     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1682                      displs, recvtype, root, comm);
1683     retval = MPI_SUCCESS;
1684   }
1685 #ifdef HAVE_TRACING
1686   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1687   TRACE_smpi_computing_in(rank);
1688 #endif
1689   smpi_bench_begin();
1690   return retval;
1691 }
1692
1693 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1694                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1695                   MPI_Comm comm)
1696 {
1697   int retval;
1698
1699   smpi_bench_end();
1700 #ifdef HAVE_TRACING
1701   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1702   TRACE_smpi_computing_out(rank);
1703   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1704 #endif
1705   if (comm == MPI_COMM_NULL) {
1706     retval = MPI_ERR_COMM;
1707   } else if (sendtype == MPI_DATATYPE_NULL
1708              || recvtype == MPI_DATATYPE_NULL) {
1709     retval = MPI_ERR_TYPE;
1710   } else {
1711     mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1712                            recvtype, comm);
1713     retval = MPI_SUCCESS;
1714   }
1715 #ifdef HAVE_TRACING
1716   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1717 #endif
1718   smpi_bench_begin();
1719   return retval;
1720 }
1721
1722 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1723                    void *recvbuf, int *recvcounts, int *displs,
1724                    MPI_Datatype recvtype, MPI_Comm comm)
1725 {
1726   int retval;
1727
1728   smpi_bench_end();
1729 #ifdef HAVE_TRACING
1730   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1731   TRACE_smpi_computing_out(rank);
1732   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1733 #endif
1734   if (comm == MPI_COMM_NULL) {
1735     retval = MPI_ERR_COMM;
1736   } else if (sendtype == MPI_DATATYPE_NULL
1737              || recvtype == MPI_DATATYPE_NULL) {
1738     retval = MPI_ERR_TYPE;
1739   } else if (recvcounts == NULL || displs == NULL) {
1740     retval = MPI_ERR_ARG;
1741   } else {
1742     mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1743                         displs, recvtype, comm);
1744     retval = MPI_SUCCESS;
1745   }
1746 #ifdef HAVE_TRACING
1747   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1748   TRACE_smpi_computing_in(rank);
1749 #endif
1750   smpi_bench_begin();
1751   return retval;
1752 }
1753
1754 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1755                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1756                 int root, MPI_Comm comm)
1757 {
1758   int retval;
1759
1760   smpi_bench_end();
1761 #ifdef HAVE_TRACING
1762   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1763   TRACE_smpi_computing_out(rank);
1764   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1765
1766   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1767 #endif
1768   if (comm == MPI_COMM_NULL) {
1769     retval = MPI_ERR_COMM;
1770   } else if (sendtype == MPI_DATATYPE_NULL
1771              || recvtype == MPI_DATATYPE_NULL) {
1772     retval = MPI_ERR_TYPE;
1773   } else {
1774     mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1775                      recvtype, root, comm);
1776     retval = MPI_SUCCESS;
1777   }
1778 #ifdef HAVE_TRACING
1779   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1780   TRACE_smpi_computing_in(rank);
1781 #endif
1782   smpi_bench_begin();
1783   return retval;
1784 }
1785
1786 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1787                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1788                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1789 {
1790   int retval;
1791
1792   smpi_bench_end();
1793 #ifdef HAVE_TRACING
1794   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1795   TRACE_smpi_computing_out(rank);
1796   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1797   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1798 #endif
1799   if (comm == MPI_COMM_NULL) {
1800     retval = MPI_ERR_COMM;
1801   } else if (sendtype == MPI_DATATYPE_NULL
1802              || recvtype == MPI_DATATYPE_NULL) {
1803     retval = MPI_ERR_TYPE;
1804   } else if (sendcounts == NULL || displs == NULL) {
1805     retval = MPI_ERR_ARG;
1806   } else {
1807     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1808                       recvcount, recvtype, root, comm);
1809     retval = MPI_SUCCESS;
1810   }
1811 #ifdef HAVE_TRACING
1812   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1813   TRACE_smpi_computing_in(rank);
1814 #endif
1815   smpi_bench_begin();
1816   return retval;
1817 }
1818
1819 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1820                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1821 {
1822   int retval;
1823
1824   smpi_bench_end();
1825 #ifdef HAVE_TRACING
1826   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1827   TRACE_smpi_computing_out(rank);
1828   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1829   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1830 #endif
1831   if (comm == MPI_COMM_NULL) {
1832     retval = MPI_ERR_COMM;
1833   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1834     retval = MPI_ERR_ARG;
1835   } else {
1836     mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
1837     retval = MPI_SUCCESS;
1838   }
1839 #ifdef HAVE_TRACING
1840   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1841   TRACE_smpi_computing_in(rank);
1842 #endif
1843   smpi_bench_begin();
1844   return retval;
1845 }
1846
1847 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
1848     MPI_Datatype datatype, MPI_Op op){
1849   int retval;
1850
1851     smpi_bench_end();
1852     if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1853       retval = MPI_ERR_ARG;
1854     } else {
1855       smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
1856       retval=MPI_SUCCESS;
1857     }
1858     smpi_bench_begin();
1859     return retval;
1860 }
1861
1862 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1863                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1864 {
1865   int retval;
1866
1867   smpi_bench_end();
1868 #ifdef HAVE_TRACING
1869   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1870   TRACE_smpi_computing_out(rank);
1871   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1872 #endif
1873   if (comm == MPI_COMM_NULL) {
1874     retval = MPI_ERR_COMM;
1875   } else if (datatype == MPI_DATATYPE_NULL) {
1876     retval = MPI_ERR_TYPE;
1877   } else if (op == MPI_OP_NULL) {
1878     retval = MPI_ERR_OP;
1879   } else {
1880       mpi_coll_allreduce_fun(sendbuf, recvbuf, count, datatype, op, comm);
1881     retval = MPI_SUCCESS;
1882   }
1883 #ifdef HAVE_TRACING
1884   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1885   TRACE_smpi_computing_in(rank);
1886 #endif
1887   smpi_bench_begin();
1888   return retval;
1889 }
1890
1891 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
1892              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1893 {
1894   int retval;
1895
1896   smpi_bench_end();
1897 #ifdef HAVE_TRACING
1898   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1899   TRACE_smpi_computing_out(rank);
1900   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1901 #endif
1902   if (comm == MPI_COMM_NULL) {
1903     retval = MPI_ERR_COMM;
1904   } else if (datatype == MPI_DATATYPE_NULL) {
1905     retval = MPI_ERR_TYPE;
1906   } else if (op == MPI_OP_NULL) {
1907     retval = MPI_ERR_OP;
1908   } else {
1909     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1910     retval = MPI_SUCCESS;
1911   }
1912 #ifdef HAVE_TRACING
1913   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1914   TRACE_smpi_computing_in(rank);
1915 #endif
1916   smpi_bench_begin();
1917   return retval;
1918 }
1919
1920 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1921                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1922 {
1923   int retval;
1924   smpi_bench_end();
1925 #ifdef HAVE_TRACING
1926   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1927   TRACE_smpi_computing_out(rank);
1928   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1929 #endif
1930   if (comm == MPI_COMM_NULL) {
1931     retval = MPI_ERR_COMM;
1932   } else if (datatype == MPI_DATATYPE_NULL) {
1933     retval = MPI_ERR_TYPE;
1934   } else if (op == MPI_OP_NULL) {
1935     retval = MPI_ERR_OP;
1936   } else if (recvcounts == NULL) {
1937     retval = MPI_ERR_ARG;
1938   } else {
1939
1940     mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
1941                        datatype,  op, comm);
1942     retval = MPI_SUCCESS;
1943   }
1944 #ifdef HAVE_TRACING
1945   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1946   TRACE_smpi_computing_in(rank);
1947 #endif
1948   smpi_bench_begin();
1949   return retval;
1950 }
1951
1952 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
1953                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1954 {
1955   int retval,i;
1956   smpi_bench_end();
1957 #ifdef HAVE_TRACING
1958   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1959   TRACE_smpi_computing_out(rank);
1960   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1961 #endif
1962   if (comm == MPI_COMM_NULL) {
1963     retval = MPI_ERR_COMM;
1964   } else if (datatype == MPI_DATATYPE_NULL) {
1965     retval = MPI_ERR_TYPE;
1966   } else if (op == MPI_OP_NULL) {
1967     retval = MPI_ERR_OP;
1968   } else if (recvcount < 0) {
1969     retval = MPI_ERR_ARG;
1970   } else {
1971     int count=smpi_comm_size(comm);
1972     int* recvcounts=(int*)xbt_malloc(count);
1973     for (i=0; i<count;i++)recvcounts[i]=recvcount;
1974     mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
1975                        datatype,  op, comm);
1976     xbt_free(recvcounts);
1977     retval = MPI_SUCCESS;
1978   }
1979 #ifdef HAVE_TRACING
1980   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1981   TRACE_smpi_computing_in(rank);
1982 #endif
1983   smpi_bench_begin();
1984   return retval;
1985 }
1986
1987 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1988                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
1989                  MPI_Comm comm)
1990 {
1991   int retval;
1992
1993   smpi_bench_end();
1994 #ifdef HAVE_TRACING
1995   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1996   TRACE_smpi_computing_out(rank);
1997   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1998 #endif
1999   if (comm == MPI_COMM_NULL) {
2000     retval = MPI_ERR_COMM;
2001   } else if (sendtype == MPI_DATATYPE_NULL
2002              || recvtype == MPI_DATATYPE_NULL) {
2003     retval = MPI_ERR_TYPE;
2004   } else {
2005     retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2006   }
2007 #ifdef HAVE_TRACING
2008   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2009   TRACE_smpi_computing_in(rank);
2010 #endif
2011   smpi_bench_begin();
2012   return retval;
2013 }
2014
2015 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2016                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2017                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2018 {
2019   int retval;
2020
2021   smpi_bench_end();
2022 #ifdef HAVE_TRACING
2023   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2024   TRACE_smpi_computing_out(rank);
2025   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2026 #endif
2027   if (comm == MPI_COMM_NULL) {
2028     retval = MPI_ERR_COMM;
2029   } else if (sendtype == MPI_DATATYPE_NULL
2030              || recvtype == MPI_DATATYPE_NULL) {
2031     retval = MPI_ERR_TYPE;
2032   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2033              || recvdisps == NULL) {
2034     retval = MPI_ERR_ARG;
2035   } else {
2036     retval =
2037         mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2038                                   recvbuf, recvcounts, recvdisps, recvtype,
2039                                   comm);
2040   }
2041 #ifdef HAVE_TRACING
2042   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2043   TRACE_smpi_computing_in(rank);
2044 #endif
2045   smpi_bench_begin();
2046   return retval;
2047 }
2048
2049
2050 int PMPI_Get_processor_name(char *name, int *resultlen)
2051 {
2052   int retval = MPI_SUCCESS;
2053
2054   smpi_bench_end();
2055   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2056           strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2057           strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2058           MPI_MAX_PROCESSOR_NAME - 1 );
2059   *resultlen =
2060       strlen(name) >
2061       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2062
2063   smpi_bench_begin();
2064   return retval;
2065 }
2066
2067 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2068 {
2069   int retval = MPI_SUCCESS;
2070   size_t size;
2071
2072   smpi_bench_end();
2073   if (status == NULL || count == NULL) {
2074     retval = MPI_ERR_ARG;
2075   } else if (datatype == MPI_DATATYPE_NULL) {
2076     retval = MPI_ERR_TYPE;
2077   } else {
2078     size = smpi_datatype_size(datatype);
2079     if (size == 0) {
2080       *count = 0;
2081     } else if (status->count % size != 0) {
2082       retval = MPI_UNDEFINED;
2083     } else {
2084       *count = smpi_mpi_get_count(status, datatype);
2085     }
2086   }
2087   smpi_bench_begin();
2088   return retval;
2089 }
2090
2091 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2092   int retval;
2093
2094   smpi_bench_end();
2095   if (old_type == MPI_DATATYPE_NULL) {
2096     retval = MPI_ERR_TYPE;
2097   } else if (count<0){
2098     retval = MPI_ERR_COUNT;
2099   } else {
2100     retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2101   }
2102   smpi_bench_begin();
2103   return retval;
2104 }
2105
2106 int PMPI_Type_commit(MPI_Datatype* datatype) {
2107   int retval;
2108
2109   smpi_bench_end();
2110   if (datatype == MPI_DATATYPE_NULL) {
2111     retval = MPI_ERR_TYPE;
2112   } else {
2113     smpi_datatype_commit(datatype);
2114     retval = MPI_SUCCESS;
2115   }
2116   smpi_bench_begin();
2117   return retval;
2118 }
2119
2120
2121 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2122   int retval;
2123
2124   smpi_bench_end();
2125   if (old_type == MPI_DATATYPE_NULL) {
2126     retval = MPI_ERR_TYPE;
2127   } else if (count<0 || blocklen<0){
2128     retval = MPI_ERR_COUNT;
2129   } else {
2130     retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2131   }
2132   smpi_bench_begin();
2133   return retval;
2134 }
2135
2136 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2137   int retval;
2138
2139   smpi_bench_end();
2140   if (old_type == MPI_DATATYPE_NULL) {
2141     retval = MPI_ERR_TYPE;
2142   } else if (count<0 || blocklen<0){
2143     retval = MPI_ERR_COUNT;
2144   } else {
2145     retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2146   }
2147   smpi_bench_begin();
2148   return retval;
2149 }
2150
2151 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2152   return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2153 }
2154
2155 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2156   int retval;
2157
2158   smpi_bench_end();
2159   if (old_type == MPI_DATATYPE_NULL) {
2160     retval = MPI_ERR_TYPE;
2161   } else if (count<0){
2162     retval = MPI_ERR_COUNT;
2163   } else {
2164     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2165   }
2166   smpi_bench_begin();
2167   return retval;
2168 }
2169
2170 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2171   int retval,i;
2172
2173   smpi_bench_end();
2174   if (old_type == MPI_DATATYPE_NULL) {
2175     retval = MPI_ERR_TYPE;
2176   } else if (count<0){
2177     retval = MPI_ERR_COUNT;
2178   } else {
2179     int* blocklens=(int*)xbt_malloc(blocklength*count);
2180     for (i=0; i<count;i++)blocklens[i]=blocklength;
2181     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2182     xbt_free(blocklens);
2183   }
2184   smpi_bench_begin();
2185   return retval;
2186 }
2187
2188
2189 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2190   int retval;
2191
2192   smpi_bench_end();
2193   if (old_type == MPI_DATATYPE_NULL) {
2194     retval = MPI_ERR_TYPE;
2195   } else if (count<0){
2196     retval = MPI_ERR_COUNT;
2197   } else {
2198     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2199   }
2200   smpi_bench_begin();
2201   return retval;
2202 }
2203
2204 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2205   int retval,i;
2206
2207   smpi_bench_end();
2208   if (old_type == MPI_DATATYPE_NULL) {
2209     retval = MPI_ERR_TYPE;
2210   } else if (count<0){
2211     retval = MPI_ERR_COUNT;
2212   } else {
2213     int* blocklens=(int*)xbt_malloc(blocklength*count);
2214     for (i=0; i<count;i++)blocklens[i]=blocklength;
2215     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2216     xbt_free(blocklens);
2217   }
2218   smpi_bench_begin();
2219   return retval;
2220 }
2221
2222
2223 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2224   int retval;
2225
2226   smpi_bench_end();
2227   if (count<0){
2228     retval = MPI_ERR_COUNT;
2229   } else {
2230     retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2231   }
2232   smpi_bench_begin();
2233   return retval;
2234 }
2235
2236 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2237   return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2238 }
2239
2240
2241 int PMPI_Error_class(int errorcode, int* errorclass) {
2242   // assume smpi uses only standard mpi error codes
2243   *errorclass=errorcode;
2244   return MPI_SUCCESS;
2245 }
2246
2247
2248 int PMPI_Initialized(int* flag) {
2249    *flag=(smpi_process_data()!=NULL);
2250    return MPI_SUCCESS;
2251 }
2252
2253 /* The following calls are not yet implemented and will fail at runtime. */
2254 /* Once implemented, please move them above this notice. */
2255
2256 #define NOT_YET_IMPLEMENTED {\
2257         XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__);\
2258         return MPI_SUCCESS;\
2259         }
2260
2261
2262 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2263   NOT_YET_IMPLEMENTED
2264 }
2265
2266 int PMPI_Type_set_name(MPI_Datatype  datatype, char * name)
2267 {
2268   NOT_YET_IMPLEMENTED
2269 }
2270
2271 int PMPI_Type_get_name(MPI_Datatype  datatype, char * name, int* len)
2272 {
2273   NOT_YET_IMPLEMENTED
2274 }
2275
2276 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2277    NOT_YET_IMPLEMENTED
2278 }
2279
2280 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2281    NOT_YET_IMPLEMENTED
2282 }
2283
2284 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
2285    NOT_YET_IMPLEMENTED
2286 }
2287
2288 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2289    NOT_YET_IMPLEMENTED
2290 }
2291
2292 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2293    NOT_YET_IMPLEMENTED
2294 }
2295
2296 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2297    NOT_YET_IMPLEMENTED
2298 }
2299
2300 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2301    NOT_YET_IMPLEMENTED
2302 }
2303
2304 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2305    NOT_YET_IMPLEMENTED
2306 }
2307
2308 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2309    NOT_YET_IMPLEMENTED
2310 }
2311
2312 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2313    NOT_YET_IMPLEMENTED
2314 }
2315
2316 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2317    NOT_YET_IMPLEMENTED
2318 }
2319
2320 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2321    NOT_YET_IMPLEMENTED
2322 }
2323
2324 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2325    NOT_YET_IMPLEMENTED
2326 }
2327
2328 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2329    NOT_YET_IMPLEMENTED
2330 }
2331
2332 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2333    NOT_YET_IMPLEMENTED
2334 }
2335
2336 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2337    NOT_YET_IMPLEMENTED
2338 }
2339
2340 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2341    NOT_YET_IMPLEMENTED
2342 }
2343
2344 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2345    NOT_YET_IMPLEMENTED
2346 }
2347
2348 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2349    NOT_YET_IMPLEMENTED
2350 }
2351
2352 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2353    NOT_YET_IMPLEMENTED
2354 }
2355
2356 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2357    NOT_YET_IMPLEMENTED
2358 }
2359
2360 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2361    NOT_YET_IMPLEMENTED
2362 }
2363
2364 int PMPI_Cancel(MPI_Request* request) {
2365    NOT_YET_IMPLEMENTED
2366 }
2367
2368 int PMPI_Buffer_attach(void* buffer, int size) {
2369    NOT_YET_IMPLEMENTED
2370 }
2371
2372 int PMPI_Buffer_detach(void* buffer, int* size) {
2373    NOT_YET_IMPLEMENTED
2374 }
2375
2376 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2377    NOT_YET_IMPLEMENTED
2378 }
2379
2380 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2381 {
2382    NOT_YET_IMPLEMENTED
2383 }
2384
2385 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2386 {
2387    NOT_YET_IMPLEMENTED
2388 }
2389
2390 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2391 {
2392    NOT_YET_IMPLEMENTED
2393 }
2394
2395 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2396 {
2397    NOT_YET_IMPLEMENTED
2398 }
2399
2400 int PMPI_Comm_free_keyval(int* keyval) {
2401    NOT_YET_IMPLEMENTED
2402 }
2403
2404 int PMPI_Pcontrol(const int level )
2405 {
2406    NOT_YET_IMPLEMENTED
2407 }
2408
2409 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2410    NOT_YET_IMPLEMENTED
2411 }
2412
2413 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2414 {
2415   NOT_YET_IMPLEMENTED
2416 }
2417
2418 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2419 {
2420   NOT_YET_IMPLEMENTED
2421 }
2422
2423 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2424 {
2425   NOT_YET_IMPLEMENTED
2426 }
2427
2428 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2429 {
2430   NOT_YET_IMPLEMENTED
2431 }
2432
2433 int PMPI_Type_free_keyval(int* keyval) {
2434   NOT_YET_IMPLEMENTED
2435 }
2436
2437 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2438    NOT_YET_IMPLEMENTED
2439 }
2440
2441 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2442    NOT_YET_IMPLEMENTED
2443 }
2444
2445 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2446    NOT_YET_IMPLEMENTED
2447 }
2448
2449 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2450    NOT_YET_IMPLEMENTED
2451 }
2452
2453 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2454    NOT_YET_IMPLEMENTED
2455 }
2456
2457 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2458    NOT_YET_IMPLEMENTED
2459 }
2460
2461 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2462    NOT_YET_IMPLEMENTED
2463 }
2464
2465 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2466    NOT_YET_IMPLEMENTED
2467 }
2468
2469 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2470    NOT_YET_IMPLEMENTED
2471 }
2472
2473 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2474    NOT_YET_IMPLEMENTED
2475 }
2476
2477 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2478    NOT_YET_IMPLEMENTED
2479 }
2480
2481 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2482    NOT_YET_IMPLEMENTED
2483 }
2484
2485 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2486    NOT_YET_IMPLEMENTED
2487 }
2488
2489 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2490    NOT_YET_IMPLEMENTED
2491 }
2492
2493 int PMPI_Keyval_free(int* keyval) {
2494    NOT_YET_IMPLEMENTED
2495 }
2496
2497 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2498    NOT_YET_IMPLEMENTED
2499 }
2500
2501 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2502    NOT_YET_IMPLEMENTED
2503 }
2504
2505 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2506   NOT_YET_IMPLEMENTED
2507 }
2508
2509 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2510   NOT_YET_IMPLEMENTED
2511 }
2512
2513 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2514   NOT_YET_IMPLEMENTED
2515 }
2516
2517 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2518    NOT_YET_IMPLEMENTED
2519 }
2520
2521 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2522    NOT_YET_IMPLEMENTED
2523 }
2524
2525 int PMPI_Win_fence( int assert,  MPI_Win win){
2526    NOT_YET_IMPLEMENTED
2527 }
2528
2529 int PMPI_Win_free( MPI_Win* win){
2530    NOT_YET_IMPLEMENTED
2531 }
2532
2533 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2534   NOT_YET_IMPLEMENTED
2535 }
2536
2537 int PMPI_Info_create( MPI_Info *info){
2538   NOT_YET_IMPLEMENTED
2539 }
2540
2541 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2542   NOT_YET_IMPLEMENTED
2543 }
2544
2545 int PMPI_Info_free( MPI_Info *info){
2546   NOT_YET_IMPLEMENTED
2547 }
2548
2549 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2550     MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2551   NOT_YET_IMPLEMENTED
2552 }
2553
2554 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2555                           int *num_addresses, int *num_datatypes, int *combiner){
2556   NOT_YET_IMPLEMENTED
2557 }
2558
2559 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2560                           int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2561                           MPI_Datatype* array_of_datatypes){
2562   NOT_YET_IMPLEMENTED
2563 }
2564
2565 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2566                             int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2567                             int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2568   NOT_YET_IMPLEMENTED
2569 }
2570
2571 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2572   NOT_YET_IMPLEMENTED
2573 }
2574
2575 int PMPI_Type_create_subarray(int ndims,int *array_of_sizes, int *array_of_subsizes, int *array_of_starts, int order, MPI_Datatype oldtype, MPI_Datatype *newtype){
2576   NOT_YET_IMPLEMENTED
2577 }
2578
2579 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
2580   NOT_YET_IMPLEMENTED
2581 }
2582
2583 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
2584                    void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
2585                    MPI_Comm comm){
2586   NOT_YET_IMPLEMENTED
2587 }
2588
2589 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2590                 MPI_Op op, MPI_Comm comm){
2591   NOT_YET_IMPLEMENTED
2592 }
2593
2594 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
2595   NOT_YET_IMPLEMENTED
2596 }
2597
2598 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
2599   NOT_YET_IMPLEMENTED
2600 }
2601
2602 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
2603   NOT_YET_IMPLEMENTED
2604 }
2605
2606 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
2607   NOT_YET_IMPLEMENTED
2608 }
2609
2610 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
2611   NOT_YET_IMPLEMENTED
2612 }
2613
2614 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
2615   NOT_YET_IMPLEMENTED
2616 }
2617
2618 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
2619   NOT_YET_IMPLEMENTED
2620 }
2621
2622 int PMPI_Add_error_class( int *errorclass){
2623   NOT_YET_IMPLEMENTED
2624 }
2625
2626 int PMPI_Add_error_code(  int errorclass, int *errorcode){
2627   NOT_YET_IMPLEMENTED
2628 }
2629
2630 int PMPI_Add_error_string( int errorcode, char *string){
2631   NOT_YET_IMPLEMENTED
2632 }
2633
2634 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
2635   NOT_YET_IMPLEMENTED
2636 }
2637
2638 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
2639   NOT_YET_IMPLEMENTED
2640 }
2641
2642 int PMPI_Info_delete(MPI_Info info, char *key){
2643   NOT_YET_IMPLEMENTED
2644 }
2645
2646 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
2647   NOT_YET_IMPLEMENTED
2648 }
2649
2650 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
2651   NOT_YET_IMPLEMENTED
2652 }
2653
2654 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
2655   NOT_YET_IMPLEMENTED
2656 }
2657
2658 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2659   NOT_YET_IMPLEMENTED
2660 }
2661
2662 int MPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2663   NOT_YET_IMPLEMENTED
2664 }
2665
2666 int PMPI_Grequest_start( MPI_Grequest_query_function *query_fn, MPI_Grequest_free_function *free_fn, MPI_Grequest_cancel_function *cancel_fn, void *extra_state, MPI_Request *request){
2667   NOT_YET_IMPLEMENTED
2668 }
2669
2670 int PMPI_Grequest_complete( MPI_Request request){
2671   NOT_YET_IMPLEMENTED
2672 }
2673
2674 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
2675   NOT_YET_IMPLEMENTED
2676 }
2677
2678 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
2679   NOT_YET_IMPLEMENTED
2680 }
2681
2682 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2683   NOT_YET_IMPLEMENTED
2684 }
2685
2686 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
2687   NOT_YET_IMPLEMENTED
2688 }
2689
2690 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
2691   NOT_YET_IMPLEMENTED
2692 }
2693
2694 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
2695   NOT_YET_IMPLEMENTED
2696 }
2697
2698 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
2699   NOT_YET_IMPLEMENTED
2700 }
2701
2702 int PMPI_Open_port( MPI_Info info, char *port_name){
2703   NOT_YET_IMPLEMENTED
2704 }
2705
2706 int PMPI_Close_port(char *port_name){
2707   NOT_YET_IMPLEMENTED
2708 }
2709
2710 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2711   NOT_YET_IMPLEMENTED
2712 }
2713
2714 int PMPI_Comm_spawn( char *command, char **argv, int maxprocs, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
2715   NOT_YET_IMPLEMENTED
2716 }
2717
2718 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
2719                              int* array_of_maxprocs, MPI_Info* array_of_info, int root,
2720                              MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
2721   NOT_YET_IMPLEMENTED
2722 }
2723
2724 int PMPI_Comm_get_parent( MPI_Comm *parent){
2725   NOT_YET_IMPLEMENTED
2726 }