Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge branch 'dvfs'
[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 (sendtype == MPI_DATATYPE_NULL
1660              || recvtype == MPI_DATATYPE_NULL) {
1661     retval = MPI_ERR_TYPE;
1662   } else {
1663     mpi_coll_gather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1664                     recvtype, root, comm);
1665     retval = MPI_SUCCESS;
1666   }
1667 #ifdef HAVE_TRACING
1668   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1669   TRACE_smpi_computing_in(rank);
1670 #endif
1671   smpi_bench_begin();
1672   return retval;
1673 }
1674
1675 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1676                 void *recvbuf, int *recvcounts, int *displs,
1677                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1678 {
1679   int retval;
1680
1681   smpi_bench_end();
1682 #ifdef HAVE_TRACING
1683   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1684   TRACE_smpi_computing_out(rank);
1685   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1686   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1687 #endif
1688   if (comm == MPI_COMM_NULL) {
1689     retval = MPI_ERR_COMM;
1690   } else if (sendtype == MPI_DATATYPE_NULL
1691              || recvtype == MPI_DATATYPE_NULL) {
1692     retval = MPI_ERR_TYPE;
1693   } else if (recvcounts == NULL || displs == NULL) {
1694     retval = MPI_ERR_ARG;
1695   } else {
1696     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1697                      displs, recvtype, root, comm);
1698     retval = MPI_SUCCESS;
1699   }
1700 #ifdef HAVE_TRACING
1701   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1702   TRACE_smpi_computing_in(rank);
1703 #endif
1704   smpi_bench_begin();
1705   return retval;
1706 }
1707
1708 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1709                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1710                   MPI_Comm comm)
1711 {
1712   int retval;
1713
1714   smpi_bench_end();
1715 #ifdef HAVE_TRACING
1716   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1717   TRACE_smpi_computing_out(rank);
1718   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1719 #endif
1720   if (comm == MPI_COMM_NULL) {
1721     retval = MPI_ERR_COMM;
1722   } else if (sendtype == MPI_DATATYPE_NULL
1723              || recvtype == MPI_DATATYPE_NULL) {
1724     retval = MPI_ERR_TYPE;
1725   } else {
1726     mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1727                            recvtype, comm);
1728     retval = MPI_SUCCESS;
1729   }
1730 #ifdef HAVE_TRACING
1731   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1732 #endif
1733   smpi_bench_begin();
1734   return retval;
1735 }
1736
1737 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1738                    void *recvbuf, int *recvcounts, int *displs,
1739                    MPI_Datatype recvtype, MPI_Comm comm)
1740 {
1741   int retval;
1742
1743   smpi_bench_end();
1744 #ifdef HAVE_TRACING
1745   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1746   TRACE_smpi_computing_out(rank);
1747   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1748 #endif
1749   if (comm == MPI_COMM_NULL) {
1750     retval = MPI_ERR_COMM;
1751   } else if (sendtype == MPI_DATATYPE_NULL
1752              || recvtype == MPI_DATATYPE_NULL) {
1753     retval = MPI_ERR_TYPE;
1754   } else if (recvcounts == NULL || displs == NULL) {
1755     retval = MPI_ERR_ARG;
1756   } else {
1757     mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1758                         displs, recvtype, comm);
1759     retval = MPI_SUCCESS;
1760   }
1761 #ifdef HAVE_TRACING
1762   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1763   TRACE_smpi_computing_in(rank);
1764 #endif
1765   smpi_bench_begin();
1766   return retval;
1767 }
1768
1769 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1770                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1771                 int root, MPI_Comm comm)
1772 {
1773   int retval;
1774
1775   smpi_bench_end();
1776 #ifdef HAVE_TRACING
1777   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1778   TRACE_smpi_computing_out(rank);
1779   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1780
1781   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1782 #endif
1783   if (comm == MPI_COMM_NULL) {
1784     retval = MPI_ERR_COMM;
1785   } else if (sendtype == MPI_DATATYPE_NULL
1786              || recvtype == MPI_DATATYPE_NULL) {
1787     retval = MPI_ERR_TYPE;
1788   } else {
1789     mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1790                      recvtype, root, comm);
1791     retval = MPI_SUCCESS;
1792   }
1793 #ifdef HAVE_TRACING
1794   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1795   TRACE_smpi_computing_in(rank);
1796 #endif
1797   smpi_bench_begin();
1798   return retval;
1799 }
1800
1801 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1802                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1803                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1804 {
1805   int retval;
1806
1807   smpi_bench_end();
1808 #ifdef HAVE_TRACING
1809   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1810   TRACE_smpi_computing_out(rank);
1811   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1812   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1813 #endif
1814   if (comm == MPI_COMM_NULL) {
1815     retval = MPI_ERR_COMM;
1816   } else if (sendtype == MPI_DATATYPE_NULL
1817              || recvtype == MPI_DATATYPE_NULL) {
1818     retval = MPI_ERR_TYPE;
1819   } else if (sendcounts == NULL || displs == NULL) {
1820     retval = MPI_ERR_ARG;
1821   } else {
1822     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1823                       recvcount, recvtype, root, comm);
1824     retval = MPI_SUCCESS;
1825   }
1826 #ifdef HAVE_TRACING
1827   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1828   TRACE_smpi_computing_in(rank);
1829 #endif
1830   smpi_bench_begin();
1831   return retval;
1832 }
1833
1834 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1835                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1836 {
1837   int retval;
1838
1839   smpi_bench_end();
1840 #ifdef HAVE_TRACING
1841   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1842   TRACE_smpi_computing_out(rank);
1843   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1844   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1845 #endif
1846   if (comm == MPI_COMM_NULL) {
1847     retval = MPI_ERR_COMM;
1848   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1849     retval = MPI_ERR_ARG;
1850   } else {
1851     mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
1852     retval = MPI_SUCCESS;
1853   }
1854 #ifdef HAVE_TRACING
1855   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1856   TRACE_smpi_computing_in(rank);
1857 #endif
1858   smpi_bench_begin();
1859   return retval;
1860 }
1861
1862 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
1863     MPI_Datatype datatype, MPI_Op op){
1864   int retval;
1865
1866     smpi_bench_end();
1867     if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1868       retval = MPI_ERR_ARG;
1869     } else {
1870       smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
1871       retval=MPI_SUCCESS;
1872     }
1873     smpi_bench_begin();
1874     return retval;
1875 }
1876
1877 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1878                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1879 {
1880   int retval;
1881
1882   smpi_bench_end();
1883 #ifdef HAVE_TRACING
1884   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1885   TRACE_smpi_computing_out(rank);
1886   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1887 #endif
1888   if (comm == MPI_COMM_NULL) {
1889     retval = MPI_ERR_COMM;
1890   } else if (datatype == MPI_DATATYPE_NULL) {
1891     retval = MPI_ERR_TYPE;
1892   } else if (op == MPI_OP_NULL) {
1893     retval = MPI_ERR_OP;
1894   } else {
1895       mpi_coll_allreduce_fun(sendbuf, recvbuf, count, datatype, op, comm);
1896     retval = MPI_SUCCESS;
1897   }
1898 #ifdef HAVE_TRACING
1899   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1900   TRACE_smpi_computing_in(rank);
1901 #endif
1902   smpi_bench_begin();
1903   return retval;
1904 }
1905
1906 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
1907              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1908 {
1909   int retval;
1910
1911   smpi_bench_end();
1912 #ifdef HAVE_TRACING
1913   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1914   TRACE_smpi_computing_out(rank);
1915   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1916 #endif
1917   if (comm == MPI_COMM_NULL) {
1918     retval = MPI_ERR_COMM;
1919   } else if (datatype == MPI_DATATYPE_NULL) {
1920     retval = MPI_ERR_TYPE;
1921   } else if (op == MPI_OP_NULL) {
1922     retval = MPI_ERR_OP;
1923   } else {
1924     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1925     retval = MPI_SUCCESS;
1926   }
1927 #ifdef HAVE_TRACING
1928   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1929   TRACE_smpi_computing_in(rank);
1930 #endif
1931   smpi_bench_begin();
1932   return retval;
1933 }
1934
1935 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1936                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1937 {
1938   int retval;
1939   smpi_bench_end();
1940 #ifdef HAVE_TRACING
1941   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1942   TRACE_smpi_computing_out(rank);
1943   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1944 #endif
1945   if (comm == MPI_COMM_NULL) {
1946     retval = MPI_ERR_COMM;
1947   } else if (datatype == MPI_DATATYPE_NULL) {
1948     retval = MPI_ERR_TYPE;
1949   } else if (op == MPI_OP_NULL) {
1950     retval = MPI_ERR_OP;
1951   } else if (recvcounts == NULL) {
1952     retval = MPI_ERR_ARG;
1953   } else {
1954
1955     mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
1956                        datatype,  op, comm);
1957     retval = MPI_SUCCESS;
1958   }
1959 #ifdef HAVE_TRACING
1960   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1961   TRACE_smpi_computing_in(rank);
1962 #endif
1963   smpi_bench_begin();
1964   return retval;
1965 }
1966
1967 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
1968                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1969 {
1970   int retval,i;
1971   smpi_bench_end();
1972 #ifdef HAVE_TRACING
1973   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1974   TRACE_smpi_computing_out(rank);
1975   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1976 #endif
1977   if (comm == MPI_COMM_NULL) {
1978     retval = MPI_ERR_COMM;
1979   } else if (datatype == MPI_DATATYPE_NULL) {
1980     retval = MPI_ERR_TYPE;
1981   } else if (op == MPI_OP_NULL) {
1982     retval = MPI_ERR_OP;
1983   } else if (recvcount < 0) {
1984     retval = MPI_ERR_ARG;
1985   } else {
1986     int count=smpi_comm_size(comm);
1987     int* recvcounts=(int*)xbt_malloc(count);
1988     for (i=0; i<count;i++)recvcounts[i]=recvcount;
1989     mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
1990                        datatype,  op, comm);
1991     xbt_free(recvcounts);
1992     retval = MPI_SUCCESS;
1993   }
1994 #ifdef HAVE_TRACING
1995   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1996   TRACE_smpi_computing_in(rank);
1997 #endif
1998   smpi_bench_begin();
1999   return retval;
2000 }
2001
2002 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2003                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
2004                  MPI_Comm comm)
2005 {
2006   int retval;
2007
2008   smpi_bench_end();
2009 #ifdef HAVE_TRACING
2010   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2011   TRACE_smpi_computing_out(rank);
2012   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2013 #endif
2014   if (comm == MPI_COMM_NULL) {
2015     retval = MPI_ERR_COMM;
2016   } else if (sendtype == MPI_DATATYPE_NULL
2017              || recvtype == MPI_DATATYPE_NULL) {
2018     retval = MPI_ERR_TYPE;
2019   } else {
2020     retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2021   }
2022 #ifdef HAVE_TRACING
2023   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2024   TRACE_smpi_computing_in(rank);
2025 #endif
2026   smpi_bench_begin();
2027   return retval;
2028 }
2029
2030 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2031                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2032                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2033 {
2034   int retval;
2035
2036   smpi_bench_end();
2037 #ifdef HAVE_TRACING
2038   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2039   TRACE_smpi_computing_out(rank);
2040   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2041 #endif
2042   if (comm == MPI_COMM_NULL) {
2043     retval = MPI_ERR_COMM;
2044   } else if (sendtype == MPI_DATATYPE_NULL
2045              || recvtype == MPI_DATATYPE_NULL) {
2046     retval = MPI_ERR_TYPE;
2047   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2048              || recvdisps == NULL) {
2049     retval = MPI_ERR_ARG;
2050   } else {
2051     retval =
2052         mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2053                                   recvbuf, recvcounts, recvdisps, recvtype,
2054                                   comm);
2055   }
2056 #ifdef HAVE_TRACING
2057   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2058   TRACE_smpi_computing_in(rank);
2059 #endif
2060   smpi_bench_begin();
2061   return retval;
2062 }
2063
2064
2065 int PMPI_Get_processor_name(char *name, int *resultlen)
2066 {
2067   int retval = MPI_SUCCESS;
2068
2069   smpi_bench_end();
2070   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2071           strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2072           strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2073           MPI_MAX_PROCESSOR_NAME - 1 );
2074   *resultlen =
2075       strlen(name) >
2076       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2077
2078   smpi_bench_begin();
2079   return retval;
2080 }
2081
2082 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2083 {
2084   int retval = MPI_SUCCESS;
2085   size_t size;
2086
2087   smpi_bench_end();
2088   if (status == NULL || count == NULL) {
2089     retval = MPI_ERR_ARG;
2090   } else if (datatype == MPI_DATATYPE_NULL) {
2091     retval = MPI_ERR_TYPE;
2092   } else {
2093     size = smpi_datatype_size(datatype);
2094     if (size == 0) {
2095       *count = 0;
2096     } else if (status->count % size != 0) {
2097       retval = MPI_UNDEFINED;
2098     } else {
2099       *count = smpi_mpi_get_count(status, datatype);
2100     }
2101   }
2102   smpi_bench_begin();
2103   return retval;
2104 }
2105
2106 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2107   int retval;
2108
2109   smpi_bench_end();
2110   if (old_type == MPI_DATATYPE_NULL) {
2111     retval = MPI_ERR_TYPE;
2112   } else if (count<0){
2113     retval = MPI_ERR_COUNT;
2114   } else {
2115     retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2116   }
2117   smpi_bench_begin();
2118   return retval;
2119 }
2120
2121 int PMPI_Type_commit(MPI_Datatype* datatype) {
2122   int retval;
2123
2124   smpi_bench_end();
2125   if (datatype == MPI_DATATYPE_NULL) {
2126     retval = MPI_ERR_TYPE;
2127   } else {
2128     smpi_datatype_commit(datatype);
2129     retval = MPI_SUCCESS;
2130   }
2131   smpi_bench_begin();
2132   return retval;
2133 }
2134
2135
2136 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2137   int retval;
2138
2139   smpi_bench_end();
2140   if (old_type == MPI_DATATYPE_NULL) {
2141     retval = MPI_ERR_TYPE;
2142   } else if (count<0 || blocklen<0){
2143     retval = MPI_ERR_COUNT;
2144   } else {
2145     retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2146   }
2147   smpi_bench_begin();
2148   return retval;
2149 }
2150
2151 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2152   int retval;
2153
2154   smpi_bench_end();
2155   if (old_type == MPI_DATATYPE_NULL) {
2156     retval = MPI_ERR_TYPE;
2157   } else if (count<0 || blocklen<0){
2158     retval = MPI_ERR_COUNT;
2159   } else {
2160     retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2161   }
2162   smpi_bench_begin();
2163   return retval;
2164 }
2165
2166 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2167   return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2168 }
2169
2170 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2171   int retval;
2172
2173   smpi_bench_end();
2174   if (old_type == MPI_DATATYPE_NULL) {
2175     retval = MPI_ERR_TYPE;
2176   } else if (count<0){
2177     retval = MPI_ERR_COUNT;
2178   } else {
2179     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2180   }
2181   smpi_bench_begin();
2182   return retval;
2183 }
2184
2185 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2186   int retval,i;
2187
2188   smpi_bench_end();
2189   if (old_type == MPI_DATATYPE_NULL) {
2190     retval = MPI_ERR_TYPE;
2191   } else if (count<0){
2192     retval = MPI_ERR_COUNT;
2193   } else {
2194     int* blocklens=(int*)xbt_malloc(blocklength*count);
2195     for (i=0; i<count;i++)blocklens[i]=blocklength;
2196     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2197     xbt_free(blocklens);
2198   }
2199   smpi_bench_begin();
2200   return retval;
2201 }
2202
2203
2204 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2205   int retval;
2206
2207   smpi_bench_end();
2208   if (old_type == MPI_DATATYPE_NULL) {
2209     retval = MPI_ERR_TYPE;
2210   } else if (count<0){
2211     retval = MPI_ERR_COUNT;
2212   } else {
2213     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2214   }
2215   smpi_bench_begin();
2216   return retval;
2217 }
2218
2219 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2220   int retval,i;
2221
2222   smpi_bench_end();
2223   if (old_type == MPI_DATATYPE_NULL) {
2224     retval = MPI_ERR_TYPE;
2225   } else if (count<0){
2226     retval = MPI_ERR_COUNT;
2227   } else {
2228     int* blocklens=(int*)xbt_malloc(blocklength*count);
2229     for (i=0; i<count;i++)blocklens[i]=blocklength;
2230     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2231     xbt_free(blocklens);
2232   }
2233   smpi_bench_begin();
2234   return retval;
2235 }
2236
2237
2238 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2239   int retval;
2240
2241   smpi_bench_end();
2242   if (count<0){
2243     retval = MPI_ERR_COUNT;
2244   } else {
2245     retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2246   }
2247   smpi_bench_begin();
2248   return retval;
2249 }
2250
2251 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2252   return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2253 }
2254
2255
2256 int PMPI_Error_class(int errorcode, int* errorclass) {
2257   // assume smpi uses only standard mpi error codes
2258   *errorclass=errorcode;
2259   return MPI_SUCCESS;
2260 }
2261
2262
2263 int PMPI_Initialized(int* flag) {
2264    *flag=(smpi_process_data()!=NULL);
2265    return MPI_SUCCESS;
2266 }
2267
2268 /* The following calls are not yet implemented and will fail at runtime. */
2269 /* Once implemented, please move them above this notice. */
2270
2271 #define NOT_YET_IMPLEMENTED {\
2272         XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__);\
2273         return MPI_SUCCESS;\
2274         }
2275
2276
2277 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2278   NOT_YET_IMPLEMENTED
2279 }
2280
2281 int PMPI_Type_set_name(MPI_Datatype  datatype, char * name)
2282 {
2283   NOT_YET_IMPLEMENTED
2284 }
2285
2286 int PMPI_Type_get_name(MPI_Datatype  datatype, char * name, int* len)
2287 {
2288   NOT_YET_IMPLEMENTED
2289 }
2290
2291 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2292    NOT_YET_IMPLEMENTED
2293 }
2294
2295 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2296    NOT_YET_IMPLEMENTED
2297 }
2298
2299 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
2300    NOT_YET_IMPLEMENTED
2301 }
2302
2303 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2304    NOT_YET_IMPLEMENTED
2305 }
2306
2307 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2308    NOT_YET_IMPLEMENTED
2309 }
2310
2311 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2312    NOT_YET_IMPLEMENTED
2313 }
2314
2315 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2316    NOT_YET_IMPLEMENTED
2317 }
2318
2319 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2320    NOT_YET_IMPLEMENTED
2321 }
2322
2323 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2324    NOT_YET_IMPLEMENTED
2325 }
2326
2327 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2328    NOT_YET_IMPLEMENTED
2329 }
2330
2331 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2332    NOT_YET_IMPLEMENTED
2333 }
2334
2335 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2336    NOT_YET_IMPLEMENTED
2337 }
2338
2339 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2340    NOT_YET_IMPLEMENTED
2341 }
2342
2343 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2344    NOT_YET_IMPLEMENTED
2345 }
2346
2347 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2348    NOT_YET_IMPLEMENTED
2349 }
2350
2351 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2352    NOT_YET_IMPLEMENTED
2353 }
2354
2355 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2356    NOT_YET_IMPLEMENTED
2357 }
2358
2359 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2360    NOT_YET_IMPLEMENTED
2361 }
2362
2363 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2364    NOT_YET_IMPLEMENTED
2365 }
2366
2367 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2368    NOT_YET_IMPLEMENTED
2369 }
2370
2371 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2372    NOT_YET_IMPLEMENTED
2373 }
2374
2375 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2376    NOT_YET_IMPLEMENTED
2377 }
2378
2379 int PMPI_Cancel(MPI_Request* request) {
2380    NOT_YET_IMPLEMENTED
2381 }
2382
2383 int PMPI_Buffer_attach(void* buffer, int size) {
2384    NOT_YET_IMPLEMENTED
2385 }
2386
2387 int PMPI_Buffer_detach(void* buffer, int* size) {
2388    NOT_YET_IMPLEMENTED
2389 }
2390
2391 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2392    NOT_YET_IMPLEMENTED
2393 }
2394
2395 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2396 {
2397    NOT_YET_IMPLEMENTED
2398 }
2399
2400 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2401 {
2402    NOT_YET_IMPLEMENTED
2403 }
2404
2405 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2406 {
2407    NOT_YET_IMPLEMENTED
2408 }
2409
2410 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2411 {
2412    NOT_YET_IMPLEMENTED
2413 }
2414
2415 int PMPI_Comm_free_keyval(int* keyval) {
2416    NOT_YET_IMPLEMENTED
2417 }
2418
2419 int PMPI_Pcontrol(const int level )
2420 {
2421    NOT_YET_IMPLEMENTED
2422 }
2423
2424 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2425    NOT_YET_IMPLEMENTED
2426 }
2427
2428 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2429 {
2430   NOT_YET_IMPLEMENTED
2431 }
2432
2433 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2434 {
2435   NOT_YET_IMPLEMENTED
2436 }
2437
2438 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2439 {
2440   NOT_YET_IMPLEMENTED
2441 }
2442
2443 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2444 {
2445   NOT_YET_IMPLEMENTED
2446 }
2447
2448 int PMPI_Type_free_keyval(int* keyval) {
2449   NOT_YET_IMPLEMENTED
2450 }
2451
2452 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2453    NOT_YET_IMPLEMENTED
2454 }
2455
2456 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2457    NOT_YET_IMPLEMENTED
2458 }
2459
2460 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2461    NOT_YET_IMPLEMENTED
2462 }
2463
2464 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2465    NOT_YET_IMPLEMENTED
2466 }
2467
2468 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2469    NOT_YET_IMPLEMENTED
2470 }
2471
2472 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2473    NOT_YET_IMPLEMENTED
2474 }
2475
2476 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2477    NOT_YET_IMPLEMENTED
2478 }
2479
2480 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2481    NOT_YET_IMPLEMENTED
2482 }
2483
2484 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2485    NOT_YET_IMPLEMENTED
2486 }
2487
2488 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2489    NOT_YET_IMPLEMENTED
2490 }
2491
2492 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2493    NOT_YET_IMPLEMENTED
2494 }
2495
2496 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2497    NOT_YET_IMPLEMENTED
2498 }
2499
2500 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2501    NOT_YET_IMPLEMENTED
2502 }
2503
2504 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2505    NOT_YET_IMPLEMENTED
2506 }
2507
2508 int PMPI_Keyval_free(int* keyval) {
2509    NOT_YET_IMPLEMENTED
2510 }
2511
2512 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2513    NOT_YET_IMPLEMENTED
2514 }
2515
2516 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2517    NOT_YET_IMPLEMENTED
2518 }
2519
2520 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2521   NOT_YET_IMPLEMENTED
2522 }
2523
2524 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2525   NOT_YET_IMPLEMENTED
2526 }
2527
2528 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2529   NOT_YET_IMPLEMENTED
2530 }
2531
2532 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2533    NOT_YET_IMPLEMENTED
2534 }
2535
2536 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2537    NOT_YET_IMPLEMENTED
2538 }
2539
2540 int PMPI_Win_fence( int assert,  MPI_Win win){
2541    NOT_YET_IMPLEMENTED
2542 }
2543
2544 int PMPI_Win_free( MPI_Win* win){
2545    NOT_YET_IMPLEMENTED
2546 }
2547
2548 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2549   NOT_YET_IMPLEMENTED
2550 }
2551
2552 int PMPI_Info_create( MPI_Info *info){
2553   NOT_YET_IMPLEMENTED
2554 }
2555
2556 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2557   NOT_YET_IMPLEMENTED
2558 }
2559
2560 int PMPI_Info_free( MPI_Info *info){
2561   NOT_YET_IMPLEMENTED
2562 }
2563
2564 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2565     MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2566   NOT_YET_IMPLEMENTED
2567 }
2568
2569 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2570                           int *num_addresses, int *num_datatypes, int *combiner){
2571   NOT_YET_IMPLEMENTED
2572 }
2573
2574 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2575                           int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2576                           MPI_Datatype* array_of_datatypes){
2577   NOT_YET_IMPLEMENTED
2578 }
2579
2580 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2581                             int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2582                             int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2583   NOT_YET_IMPLEMENTED
2584 }
2585
2586 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2587   NOT_YET_IMPLEMENTED
2588 }
2589
2590 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){
2591   NOT_YET_IMPLEMENTED
2592 }
2593
2594 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
2595   NOT_YET_IMPLEMENTED
2596 }
2597
2598 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
2599                    void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
2600                    MPI_Comm comm){
2601   NOT_YET_IMPLEMENTED
2602 }
2603
2604 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2605                 MPI_Op op, MPI_Comm comm){
2606   NOT_YET_IMPLEMENTED
2607 }
2608
2609 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
2610   NOT_YET_IMPLEMENTED
2611 }
2612
2613 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
2614   NOT_YET_IMPLEMENTED
2615 }
2616
2617 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
2618   NOT_YET_IMPLEMENTED
2619 }
2620
2621 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
2622   NOT_YET_IMPLEMENTED
2623 }
2624
2625 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
2626   NOT_YET_IMPLEMENTED
2627 }
2628
2629 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
2630   NOT_YET_IMPLEMENTED
2631 }
2632
2633 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
2634   NOT_YET_IMPLEMENTED
2635 }
2636
2637 int PMPI_Add_error_class( int *errorclass){
2638   NOT_YET_IMPLEMENTED
2639 }
2640
2641 int PMPI_Add_error_code(  int errorclass, int *errorcode){
2642   NOT_YET_IMPLEMENTED
2643 }
2644
2645 int PMPI_Add_error_string( int errorcode, char *string){
2646   NOT_YET_IMPLEMENTED
2647 }
2648
2649 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
2650   NOT_YET_IMPLEMENTED
2651 }
2652
2653 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
2654   NOT_YET_IMPLEMENTED
2655 }
2656
2657 int PMPI_Info_delete(MPI_Info info, char *key){
2658   NOT_YET_IMPLEMENTED
2659 }
2660
2661 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
2662   NOT_YET_IMPLEMENTED
2663 }
2664
2665 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
2666   NOT_YET_IMPLEMENTED
2667 }
2668
2669 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
2670   NOT_YET_IMPLEMENTED
2671 }
2672
2673 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2674   NOT_YET_IMPLEMENTED
2675 }
2676
2677 int MPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2678   NOT_YET_IMPLEMENTED
2679 }
2680
2681 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){
2682   NOT_YET_IMPLEMENTED
2683 }
2684
2685 int PMPI_Grequest_complete( MPI_Request request){
2686   NOT_YET_IMPLEMENTED
2687 }
2688
2689 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
2690   NOT_YET_IMPLEMENTED
2691 }
2692
2693 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
2694   NOT_YET_IMPLEMENTED
2695 }
2696
2697 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2698   NOT_YET_IMPLEMENTED
2699 }
2700
2701 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
2702   NOT_YET_IMPLEMENTED
2703 }
2704
2705 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
2706   NOT_YET_IMPLEMENTED
2707 }
2708
2709 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
2710   NOT_YET_IMPLEMENTED
2711 }
2712
2713 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
2714   NOT_YET_IMPLEMENTED
2715 }
2716
2717 int PMPI_Open_port( MPI_Info info, char *port_name){
2718   NOT_YET_IMPLEMENTED
2719 }
2720
2721 int PMPI_Close_port(char *port_name){
2722   NOT_YET_IMPLEMENTED
2723 }
2724
2725 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2726   NOT_YET_IMPLEMENTED
2727 }
2728
2729 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){
2730   NOT_YET_IMPLEMENTED
2731 }
2732
2733 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
2734                              int* array_of_maxprocs, MPI_Info* array_of_info, int root,
2735                              MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
2736   NOT_YET_IMPLEMENTED
2737 }
2738
2739 int PMPI_Comm_get_parent( MPI_Comm *parent){
2740   NOT_YET_IMPLEMENTED
2741 }