Logo AND Algorithmique Numérique Distribuée

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