Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
e79d12669ae9eadac594e3d87bb97f05d86d1360
[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; /* Last */
592               ) {
593           size++;
594
595           rank += ranges[i][2]; /* Stride */
596           if (ranges[i][0]<ranges[i][1]){
597               if(rank > ranges[i][1])
598                 break;
599           }else{
600               if(rank < ranges[i][1])
601                 break;
602           }
603         }
604       }
605
606         *newgroup = smpi_group_new(size);
607         j = 0;
608         for (i = 0; i < n; i++) {
609           for (rank = ranges[i][0];     /* First */
610                rank >= 0; /* Last */
611                ) {
612             index = smpi_group_index(group, rank);
613             smpi_group_set_mapping(*newgroup, index, j);
614             j++;
615             rank += ranges[i][2]; /* Stride */
616             if (ranges[i][0]<ranges[i][1]){
617               if(rank > ranges[i][1])
618                 break;
619             }else{
620               if(rank < ranges[i][1])
621                 break;
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_Group_range_excl(MPI_Group group, int n, int ranges[][3],
635                          MPI_Group * newgroup)
636 {
637   int retval, i, newrank, rank, size, index, add;
638
639   smpi_bench_end();
640   if (group == MPI_GROUP_NULL) {
641     retval = MPI_ERR_GROUP;
642   } else if (newgroup == NULL) {
643     retval = MPI_ERR_ARG;
644   } else {
645     if (n == 0) {
646       *newgroup = group;
647     } else {
648       size = smpi_group_size(group);
649       for (i = 0; i < n; i++) {
650         for (rank = ranges[i][0];       /* First */
651              rank >= 0 && rank <= ranges[i][1]; /* Last */
652              rank += ranges[i][2] /* Stride */ ) {
653           size--;
654         }
655       }
656       if (size == 0) {
657         *newgroup = MPI_GROUP_EMPTY;
658       } else {
659         *newgroup = smpi_group_new(size);
660         newrank = 0;
661         while (newrank < size) {
662           for (i = 0; i < n; i++) {
663             add = 1;
664             for (rank = ranges[i][0];   /* First */
665                  rank >= 0 && rank <= ranges[i][1];     /* Last */
666                  rank += ranges[i][2] /* Stride */ ) {
667               if (rank == newrank) {
668                 add = 0;
669                 break;
670               }
671             }
672             if (add == 1) {
673               index = smpi_group_index(group, newrank);
674               smpi_group_set_mapping(*newgroup, index, newrank);
675             }
676           }
677           newrank++; //added to avoid looping, need to be checked ..
678         }
679       }
680     }
681     smpi_group_use(*newgroup);
682     retval = MPI_SUCCESS;
683   }
684   smpi_bench_begin();
685   return retval;
686 }
687
688 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
689 {
690   int retval;
691
692   smpi_bench_end();
693   if (comm == MPI_COMM_NULL) {
694     retval = MPI_ERR_COMM;
695   } else if (rank == NULL) {
696     retval = MPI_ERR_ARG;
697   } else {
698     *rank = smpi_comm_rank(comm);
699     retval = MPI_SUCCESS;
700   }
701   smpi_bench_begin();
702   return retval;
703 }
704
705 int PMPI_Comm_size(MPI_Comm comm, int *size)
706 {
707   int retval;
708
709   smpi_bench_end();
710   if (comm == MPI_COMM_NULL) {
711     retval = MPI_ERR_COMM;
712   } else if (size == NULL) {
713     retval = MPI_ERR_ARG;
714   } else {
715     *size = smpi_comm_size(comm);
716     retval = MPI_SUCCESS;
717   }
718   smpi_bench_begin();
719   return retval;
720 }
721
722 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
723 {
724   int retval;
725
726   smpi_bench_end();
727   if (comm == MPI_COMM_NULL)  {
728     retval = MPI_ERR_COMM;
729   } else if (name == NULL || len == NULL)  {
730     retval = MPI_ERR_ARG;
731   } else {
732     smpi_comm_get_name(comm, name, len);
733     retval = MPI_SUCCESS;
734   }
735   smpi_bench_begin();
736   return retval;
737 }
738
739 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
740 {
741   int retval;
742
743   smpi_bench_end();
744   if (comm == MPI_COMM_NULL) {
745     retval = MPI_ERR_COMM;
746   } else if (group == NULL) {
747     retval = MPI_ERR_ARG;
748   } else {
749     *group = smpi_comm_group(comm);
750     retval = MPI_SUCCESS;
751   }
752   smpi_bench_begin();
753   return retval;
754 }
755
756 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
757 {
758   int retval;
759
760   smpi_bench_end();
761   if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
762     retval = MPI_ERR_COMM;
763   } else if (result == NULL) {
764     retval = MPI_ERR_ARG;
765   } else {
766     if (comm1 == comm2) {       /* Same communicators means same groups */
767       *result = MPI_IDENT;
768     } else {
769       *result =
770           smpi_group_compare(smpi_comm_group(comm1),
771                              smpi_comm_group(comm2));
772       if (*result == MPI_IDENT) {
773         *result = MPI_CONGRUENT;
774       }
775     }
776     retval = MPI_SUCCESS;
777   }
778   smpi_bench_begin();
779   return retval;
780 }
781
782 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
783 {
784   int retval;
785
786   smpi_bench_end();
787   if (comm == MPI_COMM_NULL) {
788     retval = MPI_ERR_COMM;
789   } else if (newcomm == NULL) {
790     retval = MPI_ERR_ARG;
791   } else {
792     *newcomm = smpi_comm_new(smpi_comm_group(comm));
793     retval = MPI_SUCCESS;
794   }
795   smpi_bench_begin();
796   return retval;
797 }
798
799 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
800 {
801   int retval;
802
803   smpi_bench_end();
804   if (comm == MPI_COMM_NULL) {
805     retval = MPI_ERR_COMM;
806   } else if (group == MPI_GROUP_NULL) {
807     retval = MPI_ERR_GROUP;
808   } else if (newcomm == NULL) {
809     retval = MPI_ERR_ARG;
810   } else if(smpi_group_rank(group,smpi_process_index())==MPI_UNDEFINED){
811     *newcomm= MPI_COMM_NULL;
812     retval = MPI_SUCCESS;
813   }else{
814
815     *newcomm = smpi_comm_new(group);
816     retval = MPI_SUCCESS;
817   }
818   smpi_bench_begin();
819   return retval;
820 }
821
822 int PMPI_Comm_free(MPI_Comm * comm)
823 {
824   int retval;
825
826   smpi_bench_end();
827   if (comm == NULL) {
828     retval = MPI_ERR_ARG;
829   } else if (*comm == MPI_COMM_NULL) {
830     retval = MPI_ERR_COMM;
831   } else {
832     smpi_comm_destroy(*comm);
833     *comm = MPI_COMM_NULL;
834     retval = MPI_SUCCESS;
835   }
836   smpi_bench_begin();
837   return retval;
838 }
839
840 int PMPI_Comm_disconnect(MPI_Comm * comm)
841 {
842   /* TODO: wait until all communication in comm are done */
843   int retval;
844
845   smpi_bench_end();
846   if (comm == NULL) {
847     retval = MPI_ERR_ARG;
848   } else if (*comm == MPI_COMM_NULL) {
849     retval = MPI_ERR_COMM;
850   } else {
851     smpi_comm_destroy(*comm);
852     *comm = MPI_COMM_NULL;
853     retval = MPI_SUCCESS;
854   }
855   smpi_bench_begin();
856   return retval;
857 }
858
859 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
860 {
861   int retval;
862
863   smpi_bench_end();
864   if (comm_out == NULL) {
865     retval = MPI_ERR_ARG;
866   } else if (comm == MPI_COMM_NULL) {
867     retval = MPI_ERR_COMM;
868   } else {
869     *comm_out = smpi_comm_split(comm, color, key);
870     retval = MPI_SUCCESS;
871   }
872   smpi_bench_begin();
873   return retval;
874 }
875
876 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
877                   int tag, MPI_Comm comm, MPI_Request * request)
878 {
879   int retval;
880
881   smpi_bench_end();
882   if (request == NULL) {
883     retval = MPI_ERR_ARG;
884   } else if (comm == MPI_COMM_NULL) {
885     retval = MPI_ERR_COMM;
886   } else if (dst == MPI_PROC_NULL) {
887     retval = MPI_SUCCESS;
888   } else {
889     *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
890     retval = MPI_SUCCESS;
891   }
892   smpi_bench_begin();
893   return retval;
894 }
895
896 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
897                   int tag, MPI_Comm comm, MPI_Request * request)
898 {
899   int retval;
900
901   smpi_bench_end();
902   if (request == NULL) {
903     retval = MPI_ERR_ARG;
904   } else if (comm == MPI_COMM_NULL) {
905     retval = MPI_ERR_COMM;
906   } else if (src == MPI_PROC_NULL) {
907       retval = MPI_SUCCESS;
908   } else {
909     *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
910     retval = MPI_SUCCESS;
911   }
912   smpi_bench_begin();
913   return retval;
914 }
915
916 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
917   int retval;
918
919     smpi_bench_end();
920     if (request == NULL) {
921       retval = MPI_ERR_ARG;
922     } else if (comm == MPI_COMM_NULL) {
923       retval = MPI_ERR_COMM;
924     } else if (dst == MPI_PROC_NULL) {
925       retval = MPI_SUCCESS;
926     } else {
927       *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
928       retval = MPI_SUCCESS;
929     }
930     smpi_bench_begin();
931     return retval;
932 }
933
934 int PMPI_Start(MPI_Request * request)
935 {
936   int retval;
937
938   smpi_bench_end();
939   if (request == NULL || *request == MPI_REQUEST_NULL) {
940     retval = MPI_ERR_ARG;
941   } else {
942     smpi_mpi_start(*request);
943     retval = MPI_SUCCESS;
944   }
945   smpi_bench_begin();
946   return retval;
947 }
948
949 int PMPI_Startall(int count, MPI_Request * requests)
950 {
951   int retval;
952
953   smpi_bench_end();
954   if (requests == NULL) {
955     retval = MPI_ERR_ARG;
956   } else {
957     smpi_mpi_startall(count, requests);
958     retval = MPI_SUCCESS;
959   }
960   smpi_bench_begin();
961   return retval;
962 }
963
964 int PMPI_Request_free(MPI_Request * request)
965 {
966   int retval;
967
968   smpi_bench_end();
969   if (*request == MPI_REQUEST_NULL) {
970     retval = MPI_ERR_ARG;
971   } else {
972     if((*request)->flags & PERSISTENT)(*request)->refcount--;
973     smpi_mpi_request_free(request);
974     retval = MPI_SUCCESS;
975   }
976   smpi_bench_begin();
977   return retval;
978 }
979
980 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
981               int tag, MPI_Comm comm, MPI_Request * request)
982 {
983   int retval;
984
985   smpi_bench_end();
986
987   if (request == NULL) {
988     retval = MPI_ERR_ARG;
989   } else if (comm == MPI_COMM_NULL) {
990     retval = MPI_ERR_COMM;
991   } else if (src == MPI_PROC_NULL) {
992     *request = MPI_REQUEST_NULL;
993     retval = MPI_SUCCESS;
994   } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
995     retval = MPI_ERR_COMM;
996   } else if (count < 0) {
997     retval = MPI_ERR_COUNT;
998   } else if (buf==NULL && count > 0) {
999     retval = MPI_ERR_COUNT;
1000   } else if (datatype == MPI_DATATYPE_NULL){
1001     retval = MPI_ERR_TYPE;
1002   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1003     retval = MPI_ERR_TAG;
1004   } else {
1005
1006 #ifdef HAVE_TRACING
1007   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1008   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1009   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
1010 #endif
1011
1012     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
1013     retval = MPI_SUCCESS;
1014
1015 #ifdef HAVE_TRACING
1016   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1017   (*request)->recv = 1;
1018 #endif
1019   }
1020
1021   smpi_bench_begin();
1022   return retval;
1023 }
1024
1025
1026 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
1027               int tag, MPI_Comm comm, MPI_Request * request)
1028 {
1029   int retval;
1030
1031   smpi_bench_end();
1032   if (request == NULL) {
1033     retval = MPI_ERR_ARG;
1034   } else if (comm == MPI_COMM_NULL) {
1035     retval = MPI_ERR_COMM;
1036   } else if (dst == MPI_PROC_NULL) {
1037     *request = MPI_REQUEST_NULL;
1038     retval = MPI_SUCCESS;
1039   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1040     retval = MPI_ERR_RANK;
1041   } else if (count < 0) {
1042     retval = MPI_ERR_COUNT;
1043   } else if (buf==NULL && count > 0) {
1044     retval = MPI_ERR_COUNT;
1045   } else if (datatype == MPI_DATATYPE_NULL){
1046     retval = MPI_ERR_TYPE;
1047   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1048     retval = MPI_ERR_TAG;
1049   } else {
1050
1051 #ifdef HAVE_TRACING
1052   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1053   TRACE_smpi_computing_out(rank);
1054   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1055   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1056   TRACE_smpi_send(rank, rank, dst_traced);
1057 #endif
1058
1059     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1060     retval = MPI_SUCCESS;
1061
1062 #ifdef HAVE_TRACING
1063   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1064   (*request)->send = 1;
1065   TRACE_smpi_computing_in(rank);
1066 #endif
1067   }
1068
1069   smpi_bench_begin();
1070   return retval;
1071 }
1072
1073 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
1074   int retval;
1075
1076   smpi_bench_end();
1077   if (request == NULL) {
1078     retval = MPI_ERR_ARG;
1079   } else if (comm == MPI_COMM_NULL) {
1080     retval = MPI_ERR_COMM;
1081   } else if (dst == MPI_PROC_NULL) {
1082     *request = MPI_REQUEST_NULL;
1083     retval = MPI_SUCCESS;
1084   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1085     retval = MPI_ERR_RANK;
1086   } else if (count < 0) {
1087     retval = MPI_ERR_COUNT;
1088   } else if (buf==NULL && count > 0) {
1089     retval = MPI_ERR_COUNT;
1090   } else if (datatype == MPI_DATATYPE_NULL){
1091     retval = MPI_ERR_TYPE;
1092   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1093     retval = MPI_ERR_TAG;
1094   } else {
1095
1096 #ifdef HAVE_TRACING
1097   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1098   TRACE_smpi_computing_out(rank);
1099   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1100   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1101   TRACE_smpi_send(rank, rank, dst_traced);
1102 #endif
1103
1104     *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1105     retval = MPI_SUCCESS;
1106
1107 #ifdef HAVE_TRACING
1108   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1109   (*request)->send = 1;
1110   TRACE_smpi_computing_in(rank);
1111 #endif
1112   }
1113
1114   smpi_bench_begin();
1115   return retval;
1116 }
1117
1118 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1119              MPI_Comm comm, MPI_Status * status)
1120 {
1121   int retval;
1122
1123   smpi_bench_end();
1124   if (comm == MPI_COMM_NULL) {
1125     retval = MPI_ERR_COMM;
1126   } else if (src == MPI_PROC_NULL) {
1127     smpi_empty_status(status);
1128     status->MPI_SOURCE = MPI_PROC_NULL;
1129     retval = MPI_SUCCESS;
1130   } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1131     retval = MPI_ERR_RANK;
1132   } else if (count < 0) {
1133     retval = MPI_ERR_COUNT;
1134   } else if (buf==NULL && count > 0) {
1135     retval = MPI_ERR_COUNT;
1136   } else if (datatype == MPI_DATATYPE_NULL){
1137     retval = MPI_ERR_TYPE;
1138   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1139     retval = MPI_ERR_TAG;
1140   } else {
1141 #ifdef HAVE_TRACING
1142   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1143   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1144   TRACE_smpi_computing_out(rank);
1145
1146   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
1147 #endif
1148
1149     smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1150     retval = MPI_SUCCESS;
1151
1152 #ifdef HAVE_TRACING
1153   //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1154   if(status!=MPI_STATUS_IGNORE)src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1155   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1156   TRACE_smpi_recv(rank, src_traced, rank);
1157   TRACE_smpi_computing_in(rank);
1158 #endif
1159   }
1160
1161   smpi_bench_begin();
1162   return retval;
1163 }
1164
1165 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1166              MPI_Comm comm)
1167 {
1168   int retval;
1169
1170   smpi_bench_end();
1171
1172   if (comm == MPI_COMM_NULL) {
1173     retval = MPI_ERR_COMM;
1174   } else if (dst == MPI_PROC_NULL) {
1175     retval = MPI_SUCCESS;
1176   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1177     retval = MPI_ERR_RANK;
1178   } else if (count < 0) {
1179     retval = MPI_ERR_COUNT;
1180   } else if (buf==NULL && count > 0) {
1181     retval = MPI_ERR_COUNT;
1182   } else if (datatype == MPI_DATATYPE_NULL){
1183     retval = MPI_ERR_TYPE;
1184   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1185     retval = MPI_ERR_TAG;
1186   } else {
1187
1188 #ifdef HAVE_TRACING
1189   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1190   TRACE_smpi_computing_out(rank);
1191   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1192   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1193   TRACE_smpi_send(rank, rank, dst_traced);
1194 #endif
1195
1196     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1197     retval = MPI_SUCCESS;
1198
1199 #ifdef HAVE_TRACING
1200   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1201   TRACE_smpi_computing_in(rank);
1202 #endif
1203   }
1204
1205   smpi_bench_begin();
1206   return retval;
1207 }
1208
1209
1210
1211 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1212   int retval;
1213
1214    smpi_bench_end();
1215
1216    if (comm == MPI_COMM_NULL) {
1217      retval = MPI_ERR_COMM;
1218    } else if (dst == MPI_PROC_NULL) {
1219      retval = MPI_SUCCESS;
1220    } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1221      retval = MPI_ERR_RANK;
1222    } else if (count < 0) {
1223      retval = MPI_ERR_COUNT;
1224    } else if (buf==NULL && count > 0) {
1225      retval = MPI_ERR_COUNT;
1226    } else if (datatype == MPI_DATATYPE_NULL){
1227      retval = MPI_ERR_TYPE;
1228    } else if(tag<0 && tag !=  MPI_ANY_TAG){
1229      retval = MPI_ERR_TAG;
1230    } else {
1231
1232  #ifdef HAVE_TRACING
1233    int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1234    TRACE_smpi_computing_out(rank);
1235    int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1236    TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1237    TRACE_smpi_send(rank, rank, dst_traced);
1238  #endif
1239
1240      smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1241      retval = MPI_SUCCESS;
1242
1243  #ifdef HAVE_TRACING
1244    TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1245    TRACE_smpi_computing_in(rank);
1246  #endif
1247    }
1248
1249    smpi_bench_begin();
1250    return retval;}
1251
1252
1253 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1254                  int dst, int sendtag, void *recvbuf, int recvcount,
1255                  MPI_Datatype recvtype, int src, int recvtag,
1256                  MPI_Comm comm, MPI_Status * status)
1257 {
1258   int retval;
1259
1260   smpi_bench_end();
1261
1262   if (comm == MPI_COMM_NULL) {
1263     retval = MPI_ERR_COMM;
1264   } else if (sendtype == MPI_DATATYPE_NULL
1265              || recvtype == MPI_DATATYPE_NULL) {
1266     retval = MPI_ERR_TYPE;
1267   } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1268       smpi_empty_status(status);
1269       status->MPI_SOURCE = MPI_PROC_NULL;
1270       retval = MPI_SUCCESS;
1271   }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1272       (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1273     retval = MPI_ERR_RANK;
1274   } else if (sendcount < 0 || recvcount<0) {
1275       retval = MPI_ERR_COUNT;
1276   } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1277     retval = MPI_ERR_COUNT;
1278   } else if((sendtag<0 && sendtag !=  MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1279     retval = MPI_ERR_TAG;
1280   } else {
1281
1282 #ifdef HAVE_TRACING
1283   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1284   TRACE_smpi_computing_out(rank);
1285   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1286   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1287   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1288   TRACE_smpi_send(rank, rank, dst_traced);
1289 #endif
1290
1291
1292     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1293                       recvcount, recvtype, src, recvtag, comm, status);
1294     retval = MPI_SUCCESS;
1295
1296 #ifdef HAVE_TRACING
1297   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1298   TRACE_smpi_recv(rank, src_traced, rank);
1299   TRACE_smpi_computing_in(rank);
1300 #endif
1301
1302   }
1303
1304   smpi_bench_begin();
1305   return retval;
1306 }
1307
1308 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1309                          int dst, int sendtag, int src, int recvtag,
1310                          MPI_Comm comm, MPI_Status * status)
1311 {
1312   //TODO: suboptimal implementation
1313   void *recvbuf;
1314   int retval;
1315   if (datatype == MPI_DATATYPE_NULL) {
1316       retval = MPI_ERR_TYPE;
1317   } else if (count < 0) {
1318       retval = MPI_ERR_COUNT;
1319   } else {
1320     int size = smpi_datatype_get_extent(datatype) * count;
1321     recvbuf = xbt_new(char, size);
1322     retval =
1323         MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1324                      datatype, src, recvtag, comm, status);
1325     if(retval==MPI_SUCCESS){
1326         smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1327     }
1328     xbt_free(recvbuf);
1329
1330   }
1331   return retval;
1332 }
1333
1334 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1335 {
1336   int retval;
1337
1338   smpi_bench_end();
1339   if (request == MPI_REQUEST_NULL || flag == NULL) {
1340     retval = MPI_ERR_ARG;
1341   } else if (*request == MPI_REQUEST_NULL) {
1342     *flag= TRUE;
1343     retval = MPI_ERR_REQUEST;
1344   } else {
1345     *flag = smpi_mpi_test(request, status);
1346     retval = MPI_SUCCESS;
1347   }
1348   smpi_bench_begin();
1349   return retval;
1350 }
1351
1352 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1353                 MPI_Status * status)
1354 {
1355   int retval;
1356
1357   smpi_bench_end();
1358   if (index == NULL || flag == NULL) {
1359     retval = MPI_ERR_ARG;
1360   } else {
1361     *flag = smpi_mpi_testany(count, requests, index, status);
1362     retval = MPI_SUCCESS;
1363   }
1364   smpi_bench_begin();
1365   return retval;
1366 }
1367
1368 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1369 {
1370   int retval;
1371
1372   smpi_bench_end();
1373   if (flag == NULL) {
1374     retval = MPI_ERR_ARG;
1375   } else {
1376     *flag = smpi_mpi_testall(count, requests, statuses);
1377     retval = MPI_SUCCESS;
1378   }
1379   smpi_bench_begin();
1380   return retval;
1381 }
1382
1383 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1384   int retval;
1385   smpi_bench_end();
1386
1387   if (status == NULL) {
1388     retval = MPI_ERR_ARG;
1389   } else if (comm == MPI_COMM_NULL) {
1390     retval = MPI_ERR_COMM;
1391   } else if (source == MPI_PROC_NULL) {
1392     smpi_empty_status(status);
1393     status->MPI_SOURCE = MPI_PROC_NULL;
1394     retval = MPI_SUCCESS;
1395   } else {
1396     smpi_mpi_probe(source, tag, comm, status);
1397     retval = MPI_SUCCESS;
1398   }
1399   smpi_bench_begin();
1400   return retval;
1401 }
1402
1403
1404 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1405   int retval;
1406   smpi_bench_end();
1407
1408   if (flag == NULL) {
1409     retval = MPI_ERR_ARG;
1410   } else if (status == NULL) {
1411     retval = MPI_ERR_ARG;
1412   } else if (comm == MPI_COMM_NULL) {
1413     retval = MPI_ERR_COMM;
1414   } else if (source == MPI_PROC_NULL) {
1415     smpi_empty_status(status);
1416     status->MPI_SOURCE = MPI_PROC_NULL;
1417     retval = MPI_SUCCESS;
1418   } else {
1419     smpi_mpi_iprobe(source, tag, comm, flag, status);
1420     retval = MPI_SUCCESS;
1421   }
1422   smpi_bench_begin();
1423   return retval;
1424 }
1425
1426 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1427 {
1428   int retval;
1429
1430   smpi_bench_end();
1431
1432   if (request == NULL) {
1433     retval = MPI_ERR_ARG;
1434   } else if (*request == MPI_REQUEST_NULL) {
1435     retval = MPI_ERR_REQUEST;
1436   } else {
1437
1438 #ifdef HAVE_TRACING
1439   int rank = request && (*request)->comm != MPI_COMM_NULL
1440       ? smpi_process_index()
1441       : -1;
1442   TRACE_smpi_computing_out(rank);
1443
1444   int src_traced = (*request)->src;
1445   int dst_traced = (*request)->dst;
1446   int is_wait_for_receive = (*request)->recv;
1447   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1448 #endif
1449
1450     smpi_mpi_wait(request, status);
1451     retval = MPI_SUCCESS;
1452
1453 #ifdef HAVE_TRACING
1454   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1455   if (is_wait_for_receive) {
1456     TRACE_smpi_recv(rank, src_traced, dst_traced);
1457   }
1458   TRACE_smpi_computing_in(rank);
1459 #endif
1460
1461   }
1462
1463   smpi_bench_begin();
1464   return retval;
1465 }
1466
1467 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1468 {
1469   int retval;
1470
1471   smpi_bench_end();
1472 #ifdef HAVE_TRACING
1473   //save requests information for tracing
1474   int i;
1475   int *srcs = xbt_new(int, count);
1476   int *dsts = xbt_new(int, count);
1477   int *recvs = xbt_new(int, count);
1478   for (i = 0; i < count; i++) {
1479     MPI_Request req = requests[i];      //already received requests are no longer valid
1480     if (req) {
1481       srcs[i] = req->src;
1482       dsts[i] = req->dst;
1483       recvs[i] = req->recv;
1484     }
1485   }
1486   int rank_traced = smpi_process_index();
1487   TRACE_smpi_computing_out(rank_traced);
1488
1489   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1490
1491 #endif
1492   if (index == NULL) {
1493     retval = MPI_ERR_ARG;
1494   } else {
1495     *index = smpi_mpi_waitany(count, requests, status);
1496     retval = MPI_SUCCESS;
1497   }
1498 #ifdef HAVE_TRACING
1499   if(*index!=MPI_UNDEFINED){
1500     int src_traced = srcs[*index];
1501     int dst_traced = dsts[*index];
1502     int is_wait_for_receive = recvs[*index];
1503     if (is_wait_for_receive) {
1504       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1505     }
1506     TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1507     xbt_free(srcs);
1508     xbt_free(dsts);
1509     xbt_free(recvs);
1510   }
1511   TRACE_smpi_computing_in(rank_traced);
1512 #endif
1513   smpi_bench_begin();
1514   return retval;
1515 }
1516
1517 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1518 {
1519
1520   smpi_bench_end();
1521 #ifdef HAVE_TRACING
1522   //save information from requests
1523   int i;
1524   int *srcs = xbt_new(int, count);
1525   int *dsts = xbt_new(int, count);
1526   int *recvs = xbt_new(int, count);
1527   int valid_count = 0;
1528   for (i = 0; i < count; i++) {
1529     MPI_Request req = requests[i];
1530     if(req!=MPI_REQUEST_NULL){
1531       srcs[valid_count] = req->src;
1532       dsts[valid_count] = req->dst;
1533       recvs[valid_count] = req->recv;
1534       valid_count++;
1535     }
1536   }
1537   int rank_traced = smpi_process_index();
1538   TRACE_smpi_computing_out(rank_traced);
1539
1540   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1541 #endif
1542   int retval = smpi_mpi_waitall(count, requests, status);
1543 #ifdef HAVE_TRACING
1544   for (i = 0; i < valid_count; i++) {
1545     int src_traced = srcs[i];
1546     int dst_traced = dsts[i];
1547     int is_wait_for_receive = recvs[i];
1548     if (is_wait_for_receive) {
1549       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1550     }
1551   }
1552   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1553   xbt_free(srcs);
1554   xbt_free(dsts);
1555   xbt_free(recvs);
1556   TRACE_smpi_computing_in(rank_traced);
1557 #endif
1558   smpi_bench_begin();
1559   return retval;
1560 }
1561
1562 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1563                  int *indices, MPI_Status status[])
1564 {
1565   int retval;
1566
1567   smpi_bench_end();
1568   if (outcount == NULL || indices == NULL) {
1569     retval = MPI_ERR_ARG;
1570   } else {
1571     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1572     retval = MPI_SUCCESS;
1573   }
1574   smpi_bench_begin();
1575   return retval;
1576 }
1577
1578 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1579                  int* indices, MPI_Status status[])
1580 {
1581   int retval;
1582
1583    smpi_bench_end();
1584    if (outcount == NULL || indices == NULL) {
1585      retval = MPI_ERR_ARG;
1586    } else {
1587      *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1588      retval = MPI_SUCCESS;
1589    }
1590    smpi_bench_begin();
1591    return retval;
1592 }
1593
1594
1595 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1596 {
1597   int retval;
1598
1599   smpi_bench_end();
1600 #ifdef HAVE_TRACING
1601   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1602   TRACE_smpi_computing_out(rank);
1603   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1604   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1605 #endif
1606   if (comm == MPI_COMM_NULL) {
1607     retval = MPI_ERR_COMM;
1608   } else {
1609     mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1610     retval = MPI_SUCCESS;
1611   }
1612 #ifdef HAVE_TRACING
1613   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1614   TRACE_smpi_computing_in(rank);
1615 #endif
1616   smpi_bench_begin();
1617   return retval;
1618 }
1619
1620 int PMPI_Barrier(MPI_Comm comm)
1621 {
1622   int retval;
1623
1624   smpi_bench_end();
1625 #ifdef HAVE_TRACING
1626   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1627   TRACE_smpi_computing_out(rank);
1628   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1629 #endif
1630   if (comm == MPI_COMM_NULL) {
1631     retval = MPI_ERR_COMM;
1632   } else {
1633     mpi_coll_barrier_fun(comm);
1634     retval = MPI_SUCCESS;
1635   }
1636 #ifdef HAVE_TRACING
1637   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1638   TRACE_smpi_computing_in(rank);
1639 #endif
1640   smpi_bench_begin();
1641   return retval;
1642 }
1643
1644 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1645                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1646                int root, MPI_Comm comm)
1647 {
1648   int retval;
1649
1650   smpi_bench_end();
1651 #ifdef HAVE_TRACING
1652   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1653   TRACE_smpi_computing_out(rank);
1654   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1655   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1656 #endif
1657   if (comm == MPI_COMM_NULL) {
1658     retval = MPI_ERR_COMM;
1659   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1660             ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1661     retval = MPI_ERR_TYPE;
1662   } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1663             ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1664     retval = MPI_ERR_COUNT;
1665   } else {
1666
1667     char* sendtmpbuf = (char*) sendbuf;
1668     int sendtmpcount = sendcount;
1669     MPI_Datatype sendtmptype = sendtype;
1670     if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1671       sendtmpcount=0;
1672       sendtmptype=recvtype;
1673     }
1674
1675     mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount,
1676                     recvtype, root, comm);
1677
1678
1679     retval = MPI_SUCCESS;
1680   }
1681 #ifdef HAVE_TRACING
1682   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1683   TRACE_smpi_computing_in(rank);
1684 #endif
1685   smpi_bench_begin();
1686   return retval;
1687 }
1688
1689 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1690                 void *recvbuf, int *recvcounts, int *displs,
1691                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1692 {
1693   int retval;
1694
1695   smpi_bench_end();
1696 #ifdef HAVE_TRACING
1697   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1698   TRACE_smpi_computing_out(rank);
1699   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1700   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1701 #endif
1702   if (comm == MPI_COMM_NULL) {
1703     retval = MPI_ERR_COMM;
1704   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1705             ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1706     retval = MPI_ERR_TYPE;
1707   } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1708     retval = MPI_ERR_COUNT;
1709   } else if (recvcounts == NULL || displs == NULL) {
1710     retval = MPI_ERR_ARG;
1711   } else {
1712
1713     char* sendtmpbuf = (char*) sendbuf;
1714     int sendtmpcount = sendcount;
1715     MPI_Datatype sendtmptype = sendtype;
1716     if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1717       sendtmpcount=0;
1718       sendtmptype=recvtype;
1719     }
1720
1721     smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts,
1722                      displs, recvtype, root, comm);
1723     retval = MPI_SUCCESS;
1724   }
1725 #ifdef HAVE_TRACING
1726   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1727   TRACE_smpi_computing_in(rank);
1728 #endif
1729   smpi_bench_begin();
1730   return retval;
1731 }
1732
1733 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1734                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1735                   MPI_Comm comm)
1736 {
1737   int retval;
1738
1739   smpi_bench_end();
1740 #ifdef HAVE_TRACING
1741   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1742   TRACE_smpi_computing_out(rank);
1743   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1744 #endif
1745   if (comm == MPI_COMM_NULL) {
1746     retval = MPI_ERR_COMM;
1747   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1748             (recvtype == MPI_DATATYPE_NULL)){
1749     retval = MPI_ERR_TYPE;
1750   } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1751             (recvcount <0)){
1752     retval = MPI_ERR_COUNT;
1753   } else {
1754
1755     if(sendbuf == MPI_IN_PLACE) {
1756       sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1757       sendcount=recvcount;
1758       sendtype=recvtype;
1759     }
1760
1761     mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1762                            recvtype, comm);
1763     retval = MPI_SUCCESS;
1764   }
1765 #ifdef HAVE_TRACING
1766   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1767 #endif
1768   smpi_bench_begin();
1769   return retval;
1770 }
1771
1772 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1773                    void *recvbuf, int *recvcounts, int *displs,
1774                    MPI_Datatype recvtype, MPI_Comm comm)
1775 {
1776   int retval;
1777
1778   smpi_bench_end();
1779 #ifdef HAVE_TRACING
1780   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1781   TRACE_smpi_computing_out(rank);
1782   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1783 #endif
1784   if (comm == MPI_COMM_NULL) {
1785     retval = MPI_ERR_COMM;
1786   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1787             (recvtype == MPI_DATATYPE_NULL)){
1788     retval = MPI_ERR_TYPE;
1789   } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1790     retval = MPI_ERR_COUNT;
1791   } else if (recvcounts == NULL || displs == NULL) {
1792     retval = MPI_ERR_ARG;
1793   } else {
1794
1795     if(sendbuf == MPI_IN_PLACE) {
1796       sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1797       sendcount=recvcounts[smpi_comm_rank(comm)];
1798       sendtype=recvtype;
1799     }
1800
1801     mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1802                         displs, recvtype, comm);
1803     retval = MPI_SUCCESS;
1804   }
1805 #ifdef HAVE_TRACING
1806   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1807   TRACE_smpi_computing_in(rank);
1808 #endif
1809   smpi_bench_begin();
1810   return retval;
1811 }
1812
1813 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1814                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1815                 int root, MPI_Comm comm)
1816 {
1817   int retval;
1818
1819   smpi_bench_end();
1820 #ifdef HAVE_TRACING
1821   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1822   TRACE_smpi_computing_out(rank);
1823   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1824
1825   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1826 #endif
1827   if (comm == MPI_COMM_NULL) {
1828     retval = MPI_ERR_COMM;
1829   } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1830              || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1831     retval = MPI_ERR_TYPE;
1832   } else {
1833
1834     if(recvbuf==MPI_IN_PLACE){
1835         recvcount=0;
1836     }
1837     mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1838                      recvtype, root, comm);
1839     retval = MPI_SUCCESS;
1840   }
1841 #ifdef HAVE_TRACING
1842   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1843   TRACE_smpi_computing_in(rank);
1844 #endif
1845   smpi_bench_begin();
1846   return retval;
1847 }
1848
1849 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1850                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1851                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1852 {
1853   int retval;
1854
1855   smpi_bench_end();
1856 #ifdef HAVE_TRACING
1857   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1858   TRACE_smpi_computing_out(rank);
1859   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1860   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1861 #endif
1862   if (comm == MPI_COMM_NULL) {
1863     retval = MPI_ERR_COMM;
1864   } else if (sendcounts == NULL || displs == NULL) {
1865     retval = MPI_ERR_ARG;
1866   } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1867              || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1868     retval = MPI_ERR_TYPE;
1869   } else {
1870
1871     if(recvbuf==MPI_IN_PLACE){
1872         recvcount=0;
1873     }
1874
1875     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1876                       recvcount, recvtype, root, comm);
1877     retval = MPI_SUCCESS;
1878   }
1879 #ifdef HAVE_TRACING
1880   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1881   TRACE_smpi_computing_in(rank);
1882 #endif
1883   smpi_bench_begin();
1884   return retval;
1885 }
1886
1887 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1888                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1889 {
1890   int retval;
1891
1892   smpi_bench_end();
1893 #ifdef HAVE_TRACING
1894   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1895   TRACE_smpi_computing_out(rank);
1896   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1897   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1898 #endif
1899   if (comm == MPI_COMM_NULL) {
1900     retval = MPI_ERR_COMM;
1901   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1902     retval = MPI_ERR_ARG;
1903   } else {
1904
1905     char* sendtmpbuf = (char*) sendbuf;
1906     if( sendbuf == MPI_IN_PLACE ) {
1907       sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
1908       smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
1909     }
1910
1911     mpi_coll_reduce_fun(sendtmpbuf, recvbuf, count, datatype, op, root, comm);
1912
1913     if( sendbuf == MPI_IN_PLACE ) {
1914       xbt_free(sendtmpbuf);
1915     }
1916
1917     retval = MPI_SUCCESS;
1918   }
1919 #ifdef HAVE_TRACING
1920   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1921   TRACE_smpi_computing_in(rank);
1922 #endif
1923   smpi_bench_begin();
1924   return retval;
1925 }
1926
1927 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
1928     MPI_Datatype datatype, MPI_Op op){
1929   int retval;
1930
1931     smpi_bench_end();
1932     if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1933       retval = MPI_ERR_ARG;
1934     } else {
1935       smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
1936       retval=MPI_SUCCESS;
1937     }
1938     smpi_bench_begin();
1939     return retval;
1940 }
1941
1942 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1943                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1944 {
1945   int retval;
1946
1947   smpi_bench_end();
1948 #ifdef HAVE_TRACING
1949   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1950   TRACE_smpi_computing_out(rank);
1951   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1952 #endif
1953   if (comm == MPI_COMM_NULL) {
1954     retval = MPI_ERR_COMM;
1955   } else if (datatype == MPI_DATATYPE_NULL) {
1956     retval = MPI_ERR_TYPE;
1957   } else if (op == MPI_OP_NULL) {
1958     retval = MPI_ERR_OP;
1959   } else {
1960
1961     char* sendtmpbuf = (char*) sendbuf;
1962     if( sendbuf == MPI_IN_PLACE ) {
1963       sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
1964       smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
1965     }
1966
1967     mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
1968
1969     if( sendbuf == MPI_IN_PLACE ) {
1970       xbt_free(sendtmpbuf);
1971     }
1972
1973     retval = MPI_SUCCESS;
1974
1975   }
1976 #ifdef HAVE_TRACING
1977   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1978   TRACE_smpi_computing_in(rank);
1979 #endif
1980   smpi_bench_begin();
1981   return retval;
1982 }
1983
1984 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
1985              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1986 {
1987   int retval;
1988
1989   smpi_bench_end();
1990 #ifdef HAVE_TRACING
1991   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1992   TRACE_smpi_computing_out(rank);
1993   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1994 #endif
1995   if (comm == MPI_COMM_NULL) {
1996     retval = MPI_ERR_COMM;
1997   } else if (datatype == MPI_DATATYPE_NULL) {
1998     retval = MPI_ERR_TYPE;
1999   } else if (op == MPI_OP_NULL) {
2000     retval = MPI_ERR_OP;
2001   } else {
2002     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2003     retval = MPI_SUCCESS;
2004   }
2005 #ifdef HAVE_TRACING
2006   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2007   TRACE_smpi_computing_in(rank);
2008 #endif
2009   smpi_bench_begin();
2010   return retval;
2011 }
2012
2013 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
2014                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2015 {
2016   int retval;
2017   smpi_bench_end();
2018 #ifdef HAVE_TRACING
2019   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2020   TRACE_smpi_computing_out(rank);
2021   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2022 #endif
2023   if (comm == MPI_COMM_NULL) {
2024     retval = MPI_ERR_COMM;
2025   } else if (datatype == MPI_DATATYPE_NULL) {
2026     retval = MPI_ERR_TYPE;
2027   } else if (op == MPI_OP_NULL) {
2028     retval = MPI_ERR_OP;
2029   } else if (recvcounts == NULL) {
2030     retval = MPI_ERR_ARG;
2031   } else {
2032     void* sendtmpbuf=sendbuf;
2033     if(sendbuf==MPI_IN_PLACE){
2034       sendtmpbuf=recvbuf;
2035     }
2036
2037     mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts,
2038                        datatype,  op, comm);
2039     retval = MPI_SUCCESS;
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 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2050                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2051 {
2052   int retval,i;
2053   smpi_bench_end();
2054 #ifdef HAVE_TRACING
2055   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2056   TRACE_smpi_computing_out(rank);
2057   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2058 #endif
2059   if (comm == MPI_COMM_NULL) {
2060     retval = MPI_ERR_COMM;
2061   } else if (datatype == MPI_DATATYPE_NULL) {
2062     retval = MPI_ERR_TYPE;
2063   } else if (op == MPI_OP_NULL) {
2064     retval = MPI_ERR_OP;
2065   } else if (recvcount < 0) {
2066     retval = MPI_ERR_ARG;
2067   } else {
2068     int count=smpi_comm_size(comm);
2069     int* recvcounts=(int*)xbt_malloc(count);
2070     for (i=0; i<count;i++)recvcounts[i]=recvcount;
2071     mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
2072                        datatype,  op, comm);
2073     xbt_free(recvcounts);
2074     retval = MPI_SUCCESS;
2075   }
2076 #ifdef HAVE_TRACING
2077   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2078   TRACE_smpi_computing_in(rank);
2079 #endif
2080   smpi_bench_begin();
2081   return retval;
2082 }
2083
2084 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2085                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
2086                  MPI_Comm comm)
2087 {
2088   int retval;
2089
2090   smpi_bench_end();
2091 #ifdef HAVE_TRACING
2092   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2093   TRACE_smpi_computing_out(rank);
2094   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2095 #endif
2096   if (comm == MPI_COMM_NULL) {
2097     retval = MPI_ERR_COMM;
2098   } else if (sendtype == MPI_DATATYPE_NULL
2099              || recvtype == MPI_DATATYPE_NULL) {
2100     retval = MPI_ERR_TYPE;
2101   } else {
2102     retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2103   }
2104 #ifdef HAVE_TRACING
2105   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2106   TRACE_smpi_computing_in(rank);
2107 #endif
2108   smpi_bench_begin();
2109   return retval;
2110 }
2111
2112 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2113                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2114                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2115 {
2116   int retval;
2117
2118   smpi_bench_end();
2119 #ifdef HAVE_TRACING
2120   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2121   TRACE_smpi_computing_out(rank);
2122   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2123 #endif
2124   if (comm == MPI_COMM_NULL) {
2125     retval = MPI_ERR_COMM;
2126   } else if (sendtype == MPI_DATATYPE_NULL
2127              || recvtype == MPI_DATATYPE_NULL) {
2128     retval = MPI_ERR_TYPE;
2129   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2130              || recvdisps == NULL) {
2131     retval = MPI_ERR_ARG;
2132   } else {
2133     retval =
2134         mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2135                                   recvbuf, recvcounts, recvdisps, recvtype,
2136                                   comm);
2137   }
2138 #ifdef HAVE_TRACING
2139   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2140   TRACE_smpi_computing_in(rank);
2141 #endif
2142   smpi_bench_begin();
2143   return retval;
2144 }
2145
2146
2147 int PMPI_Get_processor_name(char *name, int *resultlen)
2148 {
2149   int retval = MPI_SUCCESS;
2150
2151   smpi_bench_end();
2152   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2153           strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2154           strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2155           MPI_MAX_PROCESSOR_NAME - 1 );
2156   *resultlen =
2157       strlen(name) >
2158       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2159
2160   smpi_bench_begin();
2161   return retval;
2162 }
2163
2164 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2165 {
2166   int retval = MPI_SUCCESS;
2167   size_t size;
2168
2169   smpi_bench_end();
2170   if (status == NULL || count == NULL) {
2171     retval = MPI_ERR_ARG;
2172   } else if (datatype == MPI_DATATYPE_NULL) {
2173     retval = MPI_ERR_TYPE;
2174   } else {
2175     size = smpi_datatype_size(datatype);
2176     if (size == 0) {
2177       *count = 0;
2178     } else if (status->count % size != 0) {
2179       retval = MPI_UNDEFINED;
2180     } else {
2181       *count = smpi_mpi_get_count(status, datatype);
2182     }
2183   }
2184   smpi_bench_begin();
2185   return retval;
2186 }
2187
2188 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2189   int retval;
2190
2191   smpi_bench_end();
2192   if (old_type == MPI_DATATYPE_NULL) {
2193     retval = MPI_ERR_TYPE;
2194   } else if (count<0){
2195     retval = MPI_ERR_COUNT;
2196   } else {
2197     retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2198   }
2199   smpi_bench_begin();
2200   return retval;
2201 }
2202
2203 int PMPI_Type_commit(MPI_Datatype* datatype) {
2204   int retval;
2205
2206   smpi_bench_end();
2207   if (datatype == MPI_DATATYPE_NULL) {
2208     retval = MPI_ERR_TYPE;
2209   } else {
2210     smpi_datatype_commit(datatype);
2211     retval = MPI_SUCCESS;
2212   }
2213   smpi_bench_begin();
2214   return retval;
2215 }
2216
2217
2218 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2219   int retval;
2220
2221   smpi_bench_end();
2222   if (old_type == MPI_DATATYPE_NULL) {
2223     retval = MPI_ERR_TYPE;
2224   } else if (count<0 || blocklen<0){
2225     retval = MPI_ERR_COUNT;
2226   } else {
2227     retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2228   }
2229   smpi_bench_begin();
2230   return retval;
2231 }
2232
2233 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2234   int retval;
2235
2236   smpi_bench_end();
2237   if (old_type == MPI_DATATYPE_NULL) {
2238     retval = MPI_ERR_TYPE;
2239   } else if (count<0 || blocklen<0){
2240     retval = MPI_ERR_COUNT;
2241   } else {
2242     retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2243   }
2244   smpi_bench_begin();
2245   return retval;
2246 }
2247
2248 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2249   return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2250 }
2251
2252 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2253   int retval;
2254
2255   smpi_bench_end();
2256   if (old_type == MPI_DATATYPE_NULL) {
2257     retval = MPI_ERR_TYPE;
2258   } else if (count<0){
2259     retval = MPI_ERR_COUNT;
2260   } else {
2261     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2262   }
2263   smpi_bench_begin();
2264   return retval;
2265 }
2266
2267 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2268   int retval,i;
2269
2270   smpi_bench_end();
2271   if (old_type == MPI_DATATYPE_NULL) {
2272     retval = MPI_ERR_TYPE;
2273   } else if (count<0){
2274     retval = MPI_ERR_COUNT;
2275   } else {
2276     int* blocklens=(int*)xbt_malloc(blocklength*count);
2277     for (i=0; i<count;i++)blocklens[i]=blocklength;
2278     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2279     xbt_free(blocklens);
2280   }
2281   smpi_bench_begin();
2282   return retval;
2283 }
2284
2285
2286 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2287   int retval;
2288
2289   smpi_bench_end();
2290   if (old_type == MPI_DATATYPE_NULL) {
2291     retval = MPI_ERR_TYPE;
2292   } else if (count<0){
2293     retval = MPI_ERR_COUNT;
2294   } else {
2295     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2296   }
2297   smpi_bench_begin();
2298   return retval;
2299 }
2300
2301 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2302   int retval,i;
2303
2304   smpi_bench_end();
2305   if (old_type == MPI_DATATYPE_NULL) {
2306     retval = MPI_ERR_TYPE;
2307   } else if (count<0){
2308     retval = MPI_ERR_COUNT;
2309   } else {
2310     int* blocklens=(int*)xbt_malloc(blocklength*count);
2311     for (i=0; i<count;i++)blocklens[i]=blocklength;
2312     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2313     xbt_free(blocklens);
2314   }
2315   smpi_bench_begin();
2316   return retval;
2317 }
2318
2319
2320 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2321   int retval;
2322
2323   smpi_bench_end();
2324   if (count<0){
2325     retval = MPI_ERR_COUNT;
2326   } else {
2327     retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2328   }
2329   smpi_bench_begin();
2330   return retval;
2331 }
2332
2333 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2334   return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2335 }
2336
2337
2338 int PMPI_Error_class(int errorcode, int* errorclass) {
2339   // assume smpi uses only standard mpi error codes
2340   *errorclass=errorcode;
2341   return MPI_SUCCESS;
2342 }
2343
2344
2345 int PMPI_Initialized(int* flag) {
2346    *flag=(smpi_process_data()!=NULL);
2347    return MPI_SUCCESS;
2348 }
2349
2350 /* The following calls are not yet implemented and will fail at runtime. */
2351 /* Once implemented, please move them above this notice. */
2352
2353 #define NOT_YET_IMPLEMENTED {\
2354         XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__);\
2355         return MPI_SUCCESS;\
2356         }
2357
2358
2359 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2360   NOT_YET_IMPLEMENTED
2361 }
2362
2363 int PMPI_Type_set_name(MPI_Datatype  datatype, char * name)
2364 {
2365   NOT_YET_IMPLEMENTED
2366 }
2367
2368 int PMPI_Type_get_name(MPI_Datatype  datatype, char * name, int* len)
2369 {
2370   NOT_YET_IMPLEMENTED
2371 }
2372
2373 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2374    NOT_YET_IMPLEMENTED
2375 }
2376
2377 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2378    NOT_YET_IMPLEMENTED
2379 }
2380
2381 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
2382    NOT_YET_IMPLEMENTED
2383 }
2384
2385 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2386    NOT_YET_IMPLEMENTED
2387 }
2388
2389 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2390    NOT_YET_IMPLEMENTED
2391 }
2392
2393 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2394    NOT_YET_IMPLEMENTED
2395 }
2396
2397 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2398    NOT_YET_IMPLEMENTED
2399 }
2400
2401 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2402    NOT_YET_IMPLEMENTED
2403 }
2404
2405 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2406    NOT_YET_IMPLEMENTED
2407 }
2408
2409 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2410    NOT_YET_IMPLEMENTED
2411 }
2412
2413 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2414    NOT_YET_IMPLEMENTED
2415 }
2416
2417 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2418    NOT_YET_IMPLEMENTED
2419 }
2420
2421 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2422    NOT_YET_IMPLEMENTED
2423 }
2424
2425 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2426    NOT_YET_IMPLEMENTED
2427 }
2428
2429 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2430    NOT_YET_IMPLEMENTED
2431 }
2432
2433 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2434    NOT_YET_IMPLEMENTED
2435 }
2436
2437 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2438    NOT_YET_IMPLEMENTED
2439 }
2440
2441 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2442    NOT_YET_IMPLEMENTED
2443 }
2444
2445 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2446    NOT_YET_IMPLEMENTED
2447 }
2448
2449 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2450    NOT_YET_IMPLEMENTED
2451 }
2452
2453 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2454    NOT_YET_IMPLEMENTED
2455 }
2456
2457 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2458    NOT_YET_IMPLEMENTED
2459 }
2460
2461 int PMPI_Cancel(MPI_Request* request) {
2462    NOT_YET_IMPLEMENTED
2463 }
2464
2465 int PMPI_Buffer_attach(void* buffer, int size) {
2466    NOT_YET_IMPLEMENTED
2467 }
2468
2469 int PMPI_Buffer_detach(void* buffer, int* size) {
2470    NOT_YET_IMPLEMENTED
2471 }
2472
2473 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2474    NOT_YET_IMPLEMENTED
2475 }
2476
2477 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2478 {
2479    NOT_YET_IMPLEMENTED
2480 }
2481
2482 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2483 {
2484    NOT_YET_IMPLEMENTED
2485 }
2486
2487 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2488 {
2489    NOT_YET_IMPLEMENTED
2490 }
2491
2492 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2493 {
2494    NOT_YET_IMPLEMENTED
2495 }
2496
2497 int PMPI_Comm_free_keyval(int* keyval) {
2498    NOT_YET_IMPLEMENTED
2499 }
2500
2501 int PMPI_Pcontrol(const int level )
2502 {
2503    NOT_YET_IMPLEMENTED
2504 }
2505
2506 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2507    NOT_YET_IMPLEMENTED
2508 }
2509
2510 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2511 {
2512   NOT_YET_IMPLEMENTED
2513 }
2514
2515 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2516 {
2517   NOT_YET_IMPLEMENTED
2518 }
2519
2520 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2521 {
2522   NOT_YET_IMPLEMENTED
2523 }
2524
2525 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2526 {
2527   NOT_YET_IMPLEMENTED
2528 }
2529
2530 int PMPI_Type_free_keyval(int* keyval) {
2531   NOT_YET_IMPLEMENTED
2532 }
2533
2534 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2535    NOT_YET_IMPLEMENTED
2536 }
2537
2538 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2539    NOT_YET_IMPLEMENTED
2540 }
2541
2542 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2543    NOT_YET_IMPLEMENTED
2544 }
2545
2546 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2547    NOT_YET_IMPLEMENTED
2548 }
2549
2550 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2551    NOT_YET_IMPLEMENTED
2552 }
2553
2554 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2555    NOT_YET_IMPLEMENTED
2556 }
2557
2558 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2559    NOT_YET_IMPLEMENTED
2560 }
2561
2562 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2563    NOT_YET_IMPLEMENTED
2564 }
2565
2566 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2567    NOT_YET_IMPLEMENTED
2568 }
2569
2570 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2571    NOT_YET_IMPLEMENTED
2572 }
2573
2574 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2575    NOT_YET_IMPLEMENTED
2576 }
2577
2578 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2579    NOT_YET_IMPLEMENTED
2580 }
2581
2582 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2583    NOT_YET_IMPLEMENTED
2584 }
2585
2586 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2587    NOT_YET_IMPLEMENTED
2588 }
2589
2590 int PMPI_Keyval_free(int* keyval) {
2591    NOT_YET_IMPLEMENTED
2592 }
2593
2594 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2595    NOT_YET_IMPLEMENTED
2596 }
2597
2598 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2599    NOT_YET_IMPLEMENTED
2600 }
2601
2602 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2603   NOT_YET_IMPLEMENTED
2604 }
2605
2606 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2607   NOT_YET_IMPLEMENTED
2608 }
2609
2610 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2611   NOT_YET_IMPLEMENTED
2612 }
2613
2614 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2615    NOT_YET_IMPLEMENTED
2616 }
2617
2618 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2619    NOT_YET_IMPLEMENTED
2620 }
2621
2622 int PMPI_Win_fence( int assert,  MPI_Win win){
2623    NOT_YET_IMPLEMENTED
2624 }
2625
2626 int PMPI_Win_free( MPI_Win* win){
2627    NOT_YET_IMPLEMENTED
2628 }
2629
2630 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2631   NOT_YET_IMPLEMENTED
2632 }
2633
2634 int PMPI_Info_create( MPI_Info *info){
2635   NOT_YET_IMPLEMENTED
2636 }
2637
2638 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2639   NOT_YET_IMPLEMENTED
2640 }
2641
2642 int PMPI_Info_free( MPI_Info *info){
2643   NOT_YET_IMPLEMENTED
2644 }
2645
2646 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2647     MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2648   NOT_YET_IMPLEMENTED
2649 }
2650
2651 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2652                           int *num_addresses, int *num_datatypes, int *combiner){
2653   NOT_YET_IMPLEMENTED
2654 }
2655
2656 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2657                           int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2658                           MPI_Datatype* array_of_datatypes){
2659   NOT_YET_IMPLEMENTED
2660 }
2661
2662 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2663                             int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2664                             int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2665   NOT_YET_IMPLEMENTED
2666 }
2667
2668 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2669   NOT_YET_IMPLEMENTED
2670 }
2671
2672 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){
2673   NOT_YET_IMPLEMENTED
2674 }
2675
2676 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
2677   NOT_YET_IMPLEMENTED
2678 }
2679
2680 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
2681                    void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
2682                    MPI_Comm comm){
2683   NOT_YET_IMPLEMENTED
2684 }
2685
2686 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2687                 MPI_Op op, MPI_Comm comm){
2688   NOT_YET_IMPLEMENTED
2689 }
2690
2691 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
2692   NOT_YET_IMPLEMENTED
2693 }
2694
2695 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
2696   NOT_YET_IMPLEMENTED
2697 }
2698
2699 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
2700   NOT_YET_IMPLEMENTED
2701 }
2702
2703 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
2704   NOT_YET_IMPLEMENTED
2705 }
2706
2707 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
2708   NOT_YET_IMPLEMENTED
2709 }
2710
2711 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
2712   NOT_YET_IMPLEMENTED
2713 }
2714
2715 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
2716   NOT_YET_IMPLEMENTED
2717 }
2718
2719 int PMPI_Add_error_class( int *errorclass){
2720   NOT_YET_IMPLEMENTED
2721 }
2722
2723 int PMPI_Add_error_code(  int errorclass, int *errorcode){
2724   NOT_YET_IMPLEMENTED
2725 }
2726
2727 int PMPI_Add_error_string( int errorcode, char *string){
2728   NOT_YET_IMPLEMENTED
2729 }
2730
2731 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
2732   NOT_YET_IMPLEMENTED
2733 }
2734
2735 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
2736   NOT_YET_IMPLEMENTED
2737 }
2738
2739 int PMPI_Info_delete(MPI_Info info, char *key){
2740   NOT_YET_IMPLEMENTED
2741 }
2742
2743 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
2744   NOT_YET_IMPLEMENTED
2745 }
2746
2747 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
2748   NOT_YET_IMPLEMENTED
2749 }
2750
2751 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
2752   NOT_YET_IMPLEMENTED
2753 }
2754
2755 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2756   NOT_YET_IMPLEMENTED
2757 }
2758
2759 int MPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2760   NOT_YET_IMPLEMENTED
2761 }
2762
2763 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){
2764   NOT_YET_IMPLEMENTED
2765 }
2766
2767 int PMPI_Grequest_complete( MPI_Request request){
2768   NOT_YET_IMPLEMENTED
2769 }
2770
2771 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
2772   NOT_YET_IMPLEMENTED
2773 }
2774
2775 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
2776   NOT_YET_IMPLEMENTED
2777 }
2778
2779 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2780   NOT_YET_IMPLEMENTED
2781 }
2782
2783 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
2784   NOT_YET_IMPLEMENTED
2785 }
2786
2787 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
2788   NOT_YET_IMPLEMENTED
2789 }
2790
2791 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
2792   NOT_YET_IMPLEMENTED
2793 }
2794
2795 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
2796   NOT_YET_IMPLEMENTED
2797 }
2798
2799 int PMPI_Open_port( MPI_Info info, char *port_name){
2800   NOT_YET_IMPLEMENTED
2801 }
2802
2803 int PMPI_Close_port(char *port_name){
2804   NOT_YET_IMPLEMENTED
2805 }
2806
2807 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2808   NOT_YET_IMPLEMENTED
2809 }
2810
2811 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){
2812   NOT_YET_IMPLEMENTED
2813 }
2814
2815 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
2816                              int* array_of_maxprocs, MPI_Info* array_of_info, int root,
2817                              MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
2818   NOT_YET_IMPLEMENTED
2819 }
2820
2821 int PMPI_Comm_get_parent( MPI_Comm *parent){
2822   NOT_YET_IMPLEMENTED
2823 }