Logo AND Algorithmique Numérique Distribuée

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