Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
[SMPI] NULL -> nullptr substitution
[simgrid.git] / src / smpi / smpi_pmpi.cpp
1
2 /* Copyright (c) 2007-2015. The SimGrid Team.
3  * All rights reserved.                                                     */
4
5 /* This program is free software; you can redistribute it and/or modify it
6  * under the terms of the license (GNU LGPL) which comes with this package. */
7
8 #include "private.h"
9 #include "smpi_mpi_dt_private.h"
10
11 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_pmpi, smpi, "Logging specific to SMPI (pmpi)");
12
13 //this function need to be here because of the calls to smpi_bench
14 void TRACE_smpi_set_category(const char *category)
15 {
16   //need to end bench otherwise categories for execution tasks are wrong
17   smpi_bench_end();
18   TRACE_internal_smpi_set_category (category);
19   //begin bench after changing process's category
20   smpi_bench_begin();
21 }
22
23 /* PMPI User level calls */
24
25 int PMPI_Init(int *argc, char ***argv)
26 {
27   // PMPI_Init is call only one time by only by SMPI process
28   int already_init;
29   MPI_Initialized(&already_init);
30   if(already_init == 0){
31     smpi_process_init(argc, argv);
32     smpi_process_mark_as_initialized();
33     int rank = smpi_process_index();
34     TRACE_smpi_init(rank);
35     TRACE_smpi_computing_init(rank);
36     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
37     extra->type = TRACING_INIT;
38     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
39     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
40     smpi_bench_begin();
41   }
42
43   smpi_mpi_init();
44
45   return MPI_SUCCESS;
46 }
47
48 int PMPI_Finalize(void)
49 {
50   smpi_bench_end();
51   int rank = smpi_process_index();
52   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
53   extra->type = TRACING_FINALIZE;
54   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
55
56   smpi_process_finalize();
57
58   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
59   TRACE_smpi_finalize(smpi_process_index());
60   smpi_process_destroy();
61   return MPI_SUCCESS;
62 }
63
64 int PMPI_Finalized(int* flag)
65 {
66   *flag=smpi_process_finalized();
67   return MPI_SUCCESS;
68 }
69
70 int PMPI_Get_version (int *version,int *subversion){
71   *version = MPI_VERSION;
72   *subversion= MPI_SUBVERSION;
73   return MPI_SUCCESS;
74 }
75
76 int PMPI_Get_library_version (char *version,int *len){
77   int retval = MPI_SUCCESS;
78   smpi_bench_end();
79   snprintf(version,MPI_MAX_LIBRARY_VERSION_STRING,"SMPI Version %d.%d. Copyright The Simgrid Team 2007-2015",
80            SIMGRID_VERSION_MAJOR, SIMGRID_VERSION_MINOR);
81   *len = strlen(version) > MPI_MAX_LIBRARY_VERSION_STRING ? MPI_MAX_LIBRARY_VERSION_STRING : strlen(version);
82   smpi_bench_begin();
83   return retval;
84 }
85
86 int PMPI_Init_thread(int *argc, char ***argv, int required, int *provided)
87 {
88   if (provided != nullptr) {
89     *provided = MPI_THREAD_SINGLE;
90   }
91   return MPI_Init(argc, argv);
92 }
93
94 int PMPI_Query_thread(int *provided)
95 {
96   int retval = 0;
97
98   if (provided == nullptr) {
99     retval = MPI_ERR_ARG;
100   } else {
101     *provided = MPI_THREAD_SINGLE;
102     retval = MPI_SUCCESS;
103   }
104   return retval;
105 }
106
107 int PMPI_Is_thread_main(int *flag)
108 {
109   int retval = 0;
110
111   if (flag == nullptr) {
112     retval = MPI_ERR_ARG;
113   } else {
114     *flag = smpi_process_index() == 0;
115     retval = MPI_SUCCESS;
116   }
117   return retval;
118 }
119
120 int PMPI_Abort(MPI_Comm comm, int errorcode)
121 {
122   smpi_bench_end();
123   smpi_process_destroy();
124   // FIXME: should kill all processes in comm instead
125   simcall_process_kill(SIMIX_process_self());
126   return MPI_SUCCESS;
127 }
128
129 double PMPI_Wtime(void)
130 {
131   return smpi_mpi_wtime();
132 }
133
134 extern double sg_maxmin_precision;
135 double PMPI_Wtick(void)
136 {
137   return sg_maxmin_precision;
138 }
139
140 int PMPI_Address(void *location, MPI_Aint * address)
141 {
142   int retval = 0;
143
144   if (address==nullptr) {
145     retval = MPI_ERR_ARG;
146   } else {
147     *address = reinterpret_cast<MPI_Aint>(location);
148     retval = MPI_SUCCESS;
149   }
150   return retval;
151 }
152
153 int PMPI_Get_address(void *location, MPI_Aint * address)
154 {
155   return PMPI_Address(location, address);
156 }
157
158 int PMPI_Type_free(MPI_Datatype * datatype)
159 {
160   int retval = 0;
161   /* Free a predefined datatype is an error according to the standard, and should be checked for */
162   if (*datatype == MPI_DATATYPE_NULL) {
163     retval = MPI_ERR_ARG;
164   } else {
165     smpi_datatype_unuse(*datatype);
166     retval = MPI_SUCCESS;
167   }
168   return retval;
169 }
170
171 int PMPI_Type_size(MPI_Datatype datatype, int *size)
172 {
173   int retval = 0;
174
175   if (datatype == MPI_DATATYPE_NULL) {
176     retval = MPI_ERR_TYPE;
177   } else if (size == nullptr) {
178     retval = MPI_ERR_ARG;
179   } else {
180     *size = static_cast<int>(smpi_datatype_size(datatype));
181     retval = MPI_SUCCESS;
182   }
183   return retval;
184 }
185
186 int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
187 {
188   int retval = 0;
189
190   if (datatype == MPI_DATATYPE_NULL) {
191     retval = MPI_ERR_TYPE;
192   } else if (lb == nullptr || extent == nullptr) {
193     retval = MPI_ERR_ARG;
194   } else {
195     retval = smpi_datatype_extent(datatype, lb, extent);
196   }
197   return retval;
198 }
199
200 int PMPI_Type_get_true_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
201 {
202   return PMPI_Type_get_extent(datatype, lb, extent);
203 }
204
205 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
206 {
207   int retval = 0;
208
209   if (datatype == MPI_DATATYPE_NULL) {
210     retval = MPI_ERR_TYPE;
211   } else if (extent == nullptr) {
212     retval = MPI_ERR_ARG;
213   } else {
214     *extent = smpi_datatype_get_extent(datatype);
215     retval = MPI_SUCCESS;
216   }
217   return retval;
218 }
219
220 int PMPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
221 {
222   int retval = 0;
223
224   if (datatype == MPI_DATATYPE_NULL) {
225     retval = MPI_ERR_TYPE;
226   } else if (disp == nullptr) {
227     retval = MPI_ERR_ARG;
228   } else {
229     *disp = smpi_datatype_lb(datatype);
230     retval = MPI_SUCCESS;
231   }
232   return retval;
233 }
234
235 int PMPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
236 {
237   int retval = 0;
238
239   if (datatype == MPI_DATATYPE_NULL) {
240     retval = MPI_ERR_TYPE;
241   } else if (disp == nullptr) {
242     retval = MPI_ERR_ARG;
243   } else {
244     *disp = smpi_datatype_ub(datatype);
245     retval = MPI_SUCCESS;
246   }
247   return retval;
248 }
249
250 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
251   int retval = 0;
252
253   if (datatype == MPI_DATATYPE_NULL) {
254     retval = MPI_ERR_TYPE;
255   } else {
256     retval = smpi_datatype_dup(datatype, newtype);
257   }
258   return retval;
259 }
260
261 int PMPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
262 {
263   int retval = 0;
264
265   if (function == nullptr || op == nullptr) {
266     retval = MPI_ERR_ARG;
267   } else {
268     *op = smpi_op_new(function, (commute!=0));
269     retval = MPI_SUCCESS;
270   }
271   return retval;
272 }
273
274 int PMPI_Op_free(MPI_Op * op)
275 {
276   int retval = 0;
277
278   if (op == nullptr) {
279     retval = MPI_ERR_ARG;
280   } else if (*op == MPI_OP_NULL) {
281     retval = MPI_ERR_OP;
282   } else {
283     smpi_op_destroy(*op);
284     *op = MPI_OP_NULL;
285     retval = MPI_SUCCESS;
286   }
287   return retval;
288 }
289
290 int PMPI_Group_free(MPI_Group * group)
291 {
292   int retval = 0;
293
294   if (group == nullptr) {
295     retval = MPI_ERR_ARG;
296   } else {
297     smpi_group_destroy(*group);
298     *group = MPI_GROUP_NULL;
299     retval = MPI_SUCCESS;
300   }
301   return retval;
302 }
303
304 int PMPI_Group_size(MPI_Group group, int *size)
305 {
306   int retval = 0;
307
308   if (group == MPI_GROUP_NULL) {
309     retval = MPI_ERR_GROUP;
310   } else if (size == nullptr) {
311     retval = MPI_ERR_ARG;
312   } else {
313     *size = smpi_group_size(group);
314     retval = MPI_SUCCESS;
315   }
316   return retval;
317 }
318
319 int PMPI_Group_rank(MPI_Group group, int *rank)
320 {
321   int retval = 0;
322
323   if (group == MPI_GROUP_NULL) {
324     retval = MPI_ERR_GROUP;
325   } else if (rank == nullptr) {
326     retval = MPI_ERR_ARG;
327   } else {
328     *rank = smpi_group_rank(group, smpi_process_index());
329     retval = MPI_SUCCESS;
330   }
331   return retval;
332 }
333
334 int PMPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1, MPI_Group group2, int *ranks2)
335 {
336   int retval, i, index;
337   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
338     retval = MPI_ERR_GROUP;
339   } else {
340     for (i = 0; i < n; i++) {
341       if(ranks1[i]==MPI_PROC_NULL){
342         ranks2[i]=MPI_PROC_NULL;
343       }else{
344         index = smpi_group_index(group1, ranks1[i]);
345         ranks2[i] = smpi_group_rank(group2, index);
346       }
347     }
348     retval = MPI_SUCCESS;
349   }
350   return retval;
351 }
352
353 int PMPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
354 {
355   int retval = 0;
356
357   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
358     retval = MPI_ERR_GROUP;
359   } else if (result == nullptr) {
360     retval = MPI_ERR_ARG;
361   } else {
362     *result = smpi_group_compare(group1, group2);
363     retval = MPI_SUCCESS;
364   }
365   return retval;
366 }
367
368 int PMPI_Group_union(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
369 {
370   int retval, i, proc1, proc2, size, size2;
371
372   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
373     retval = MPI_ERR_GROUP;
374   } else if (newgroup == nullptr) {
375     retval = MPI_ERR_ARG;
376   } else {
377     size = smpi_group_size(group1);
378     size2 = smpi_group_size(group2);
379     for (i = 0; i < size2; i++) {
380       proc2 = smpi_group_index(group2, i);
381       proc1 = smpi_group_rank(group1, proc2);
382       if (proc1 == MPI_UNDEFINED) {
383         size++;
384       }
385     }
386     if (size == 0) {
387       *newgroup = MPI_GROUP_EMPTY;
388     } else {
389       *newgroup = smpi_group_new(size);
390       size2 = smpi_group_size(group1);
391       for (i = 0; i < size2; i++) {
392         proc1 = smpi_group_index(group1, i);
393         smpi_group_set_mapping(*newgroup, proc1, i);
394       }
395       for (i = size2; i < size; i++) {
396         proc2 = smpi_group_index(group2, i - size2);
397         smpi_group_set_mapping(*newgroup, proc2, i);
398       }
399     }
400     retval = MPI_SUCCESS;
401   }
402   return retval;
403 }
404
405 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
406 {
407   int retval, i, proc1, proc2, size;
408
409   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
410     retval = MPI_ERR_GROUP;
411   } else if (newgroup == nullptr) {
412     retval = MPI_ERR_ARG;
413   } else {
414     size = smpi_group_size(group2);
415     for (i = 0; i < size; i++) {
416       proc2 = smpi_group_index(group2, i);
417       proc1 = smpi_group_rank(group1, proc2);
418       if (proc1 == MPI_UNDEFINED) {
419         size--;
420       }
421     }
422     if (size == 0) {
423       *newgroup = MPI_GROUP_EMPTY;
424     } else {
425       *newgroup = smpi_group_new(size);
426       int j=0;
427       for (i = 0; i < smpi_group_size(group2); i++) {
428         proc2 = smpi_group_index(group2, i);
429         proc1 = smpi_group_rank(group1, proc2);
430         if (proc1 != MPI_UNDEFINED) {
431           smpi_group_set_mapping(*newgroup, proc2, j);
432           j++;
433         }
434       }
435     }
436     retval = MPI_SUCCESS;
437   }
438   return retval;
439 }
440
441 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
442 {
443   int retval, i, proc1, proc2, size, size2;
444
445   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
446     retval = MPI_ERR_GROUP;
447   } else if (newgroup == nullptr) {
448     retval = MPI_ERR_ARG;
449   } else {
450     size = smpi_group_size(group1);
451     size2 = size;
452     for (i = 0; i < size2; i++) {
453       proc1 = smpi_group_index(group1, i);
454       proc2 = smpi_group_rank(group2, proc1);
455       if (proc2 != MPI_UNDEFINED) {
456         size--;
457       }
458     }
459     if (size == 0) {
460       *newgroup = MPI_GROUP_EMPTY;
461     } else {
462       *newgroup = smpi_group_new(size);
463       for (i = 0; i < size2; i++) {
464         proc1 = smpi_group_index(group1, i);
465         proc2 = smpi_group_rank(group2, proc1);
466         if (proc2 == MPI_UNDEFINED) {
467           smpi_group_set_mapping(*newgroup, proc1, i);
468         }
469       }
470     }
471     retval = MPI_SUCCESS;
472   }
473   return retval;
474 }
475
476 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
477 {
478   int retval;
479
480   if (group == MPI_GROUP_NULL) {
481     retval = MPI_ERR_GROUP;
482   } else if (newgroup == nullptr) {
483     retval = MPI_ERR_ARG;
484   } else {
485     retval = smpi_group_incl(group, n, ranks, newgroup);
486   }
487   return retval;
488 }
489
490 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
491 {
492   int retval, i, j, newsize, oldsize, index;
493
494   if (group == MPI_GROUP_NULL) {
495     retval = MPI_ERR_GROUP;
496   } else if (newgroup == nullptr) {
497     retval = MPI_ERR_ARG;
498   } else {
499     if (n == 0) {
500       *newgroup = group;
501       if(group!= smpi_comm_group(MPI_COMM_WORLD) && group != MPI_GROUP_NULL
502                 && group != smpi_comm_group(MPI_COMM_SELF) && group != MPI_GROUP_EMPTY)
503       smpi_group_use(group);
504     } else if (n == smpi_group_size(group)) {
505       *newgroup = MPI_GROUP_EMPTY;
506     } else {
507       oldsize=smpi_group_size(group);
508       newsize = oldsize - n;
509       *newgroup = smpi_group_new(newsize);
510
511       int* to_exclude=xbt_new0(int, smpi_group_size(group));
512       for(i=0; i<oldsize; i++)
513         to_exclude[i]=0;
514       for(i=0; i<n; i++)
515         to_exclude[ranks[i]]=1;
516
517       j=0;
518       for(i=0; i<oldsize; i++){
519         if(to_exclude[i]==0){
520           index = smpi_group_index(group, i);
521           smpi_group_set_mapping(*newgroup, index, j);
522           j++;
523         }
524       }
525
526       xbt_free(to_exclude);
527     }
528     retval = MPI_SUCCESS;
529   }
530   return retval;
531 }
532
533 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3], MPI_Group * newgroup)
534 {
535   int retval, i, j, rank, size, index;
536
537   if (group == MPI_GROUP_NULL) {
538     retval = MPI_ERR_GROUP;
539   } else if (newgroup == nullptr) {
540     retval = MPI_ERR_ARG;
541   } else {
542     if (n == 0) {
543       *newgroup = MPI_GROUP_EMPTY;
544     } else {
545       size = 0;
546       for (i = 0; i < n; i++) {
547         for (rank = ranges[i][0];       /* First */
548              rank >= 0 && rank < smpi_group_size(group); /* Last */
549               ) {
550           size++;
551           if(rank == ranges[i][1]){/*already last ?*/
552             break;
553           }
554           rank += ranges[i][2]; /* Stride */
555     if (ranges[i][0]<ranges[i][1]){
556         if(rank > ranges[i][1])
557           break;
558     }else{
559         if(rank < ranges[i][1])
560           break;
561     }
562         }
563       }
564
565       *newgroup = smpi_group_new(size);
566       j = 0;
567       for (i = 0; i < n; i++) {
568         for (rank = ranges[i][0];     /* First */
569              rank >= 0 && rank < smpi_group_size(group); /* Last */
570              ) {
571           index = smpi_group_index(group, rank);
572           smpi_group_set_mapping(*newgroup, index, j);
573           j++;
574           if(rank == ranges[i][1]){/*already last ?*/
575             break;
576           }
577           rank += ranges[i][2]; /* Stride */
578     if (ranges[i][0]<ranges[i][1]){
579       if(rank > ranges[i][1])
580         break;
581     }else{
582       if(rank < ranges[i][1])
583         break;
584     }
585         }
586       }
587     }
588     retval = MPI_SUCCESS;
589   }
590   return retval;
591 }
592
593 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3], MPI_Group * newgroup)
594 {
595   int retval, i, rank, newrank,oldrank, size, index, add;
596
597   if (group == MPI_GROUP_NULL) {
598     retval = MPI_ERR_GROUP;
599   } else if (newgroup == nullptr) {
600     retval = MPI_ERR_ARG;
601   } else {
602     if (n == 0) {
603       *newgroup = group;
604       if(group!= smpi_comm_group(MPI_COMM_WORLD) && group != MPI_GROUP_NULL
605                 && group != smpi_comm_group(MPI_COMM_SELF) && group != MPI_GROUP_EMPTY)
606       smpi_group_use(group);
607     } else {
608       size = smpi_group_size(group);
609       for (i = 0; i < n; i++) {
610         for (rank = ranges[i][0];       /* First */
611              rank >= 0 && rank < smpi_group_size(group); /* Last */
612               ) {
613           size--;
614           if(rank == ranges[i][1]){/*already last ?*/
615             break;
616           }
617           rank += ranges[i][2]; /* Stride */
618     if (ranges[i][0]<ranges[i][1]){
619         if(rank > ranges[i][1])
620           break;
621     }else{
622         if(rank < ranges[i][1])
623           break;
624     }
625         }
626       }
627       if (size == 0) {
628         *newgroup = MPI_GROUP_EMPTY;
629       } else {
630         *newgroup = smpi_group_new(size);
631         newrank=0;
632         oldrank=0;
633         while (newrank < size) {
634           add=1;
635           for (i = 0; i < n; i++) {
636             for (rank = ranges[i][0];
637                 rank >= 0 && rank < smpi_group_size(group);
638                 ){
639               if(rank==oldrank){
640                   add=0;
641                   break;
642               }
643               if(rank == ranges[i][1]){/*already last ?*/
644                 break;
645               }
646               rank += ranges[i][2]; /* Stride */
647               if (ranges[i][0]<ranges[i][1]){
648                   if(rank > ranges[i][1])
649                     break;
650               }else{
651                   if(rank < ranges[i][1])
652                     break;
653               }
654             }
655           }
656           if(add==1){
657             index = smpi_group_index(group, oldrank);
658             smpi_group_set_mapping(*newgroup, index, newrank);
659             newrank++;
660           }
661           oldrank++;
662         }
663       }
664     }
665
666     retval = MPI_SUCCESS;
667   }
668   return retval;
669 }
670
671 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
672 {
673   int retval = 0;
674   if (comm == MPI_COMM_NULL) {
675     retval = MPI_ERR_COMM;
676   } else if (rank == nullptr) {
677     retval = MPI_ERR_ARG;
678   } else {
679     *rank = smpi_comm_rank(comm);
680     retval = MPI_SUCCESS;
681   }
682   return retval;
683 }
684
685 int PMPI_Comm_size(MPI_Comm comm, int *size)
686 {
687   int retval = 0;
688   if (comm == MPI_COMM_NULL) {
689     retval = MPI_ERR_COMM;
690   } else if (size == nullptr) {
691     retval = MPI_ERR_ARG;
692   } else {
693     *size = smpi_comm_size(comm);
694     retval = MPI_SUCCESS;
695   }
696   return retval;
697 }
698
699 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
700 {
701   int retval = 0;
702
703   if (comm == MPI_COMM_NULL)  {
704     retval = MPI_ERR_COMM;
705   } else if (name == nullptr || len == nullptr)  {
706     retval = MPI_ERR_ARG;
707   } else {
708     smpi_comm_get_name(comm, name, len);
709     retval = MPI_SUCCESS;
710   }
711   return retval;
712 }
713
714 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
715 {
716   int retval = 0;
717
718   if (comm == MPI_COMM_NULL) {
719     retval = MPI_ERR_COMM;
720   } else if (group == nullptr) {
721     retval = MPI_ERR_ARG;
722   } else {
723     *group = smpi_comm_group(comm);
724     if(*group!= smpi_comm_group(MPI_COMM_WORLD) && *group != MPI_GROUP_NULL
725               && *group != MPI_GROUP_EMPTY)
726     smpi_group_use(*group);
727     retval = MPI_SUCCESS;
728   }
729   return retval;
730 }
731
732 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
733 {
734   int retval = 0;
735
736   if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
737     retval = MPI_ERR_COMM;
738   } else if (result == nullptr) {
739     retval = MPI_ERR_ARG;
740   } else {
741     if (comm1 == comm2) {       /* Same communicators means same groups */
742       *result = MPI_IDENT;
743     } else {
744       *result = smpi_group_compare(smpi_comm_group(comm1), smpi_comm_group(comm2));
745       if (*result == MPI_IDENT) {
746         *result = MPI_CONGRUENT;
747       }
748     }
749     retval = MPI_SUCCESS;
750   }
751   return retval;
752 }
753
754 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
755 {
756   int retval = 0;
757
758   if (comm == MPI_COMM_NULL) {
759     retval = MPI_ERR_COMM;
760   } else if (newcomm == nullptr) {
761     retval = MPI_ERR_ARG;
762   } else {
763     retval = smpi_comm_dup(comm, newcomm);
764   }
765   return retval;
766 }
767
768 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
769 {
770   int retval = 0;
771
772   if (comm == MPI_COMM_NULL) {
773     retval = MPI_ERR_COMM;
774   } else if (group == MPI_GROUP_NULL) {
775     retval = MPI_ERR_GROUP;
776   } else if (newcomm == nullptr) {
777     retval = MPI_ERR_ARG;
778   } else if(smpi_group_rank(group,smpi_process_index())==MPI_UNDEFINED){
779     *newcomm= MPI_COMM_NULL;
780     retval = MPI_SUCCESS;
781   }else{
782     smpi_group_use(group);
783     *newcomm = smpi_comm_new(group, nullptr);
784     retval = MPI_SUCCESS;
785   }
786   return retval;
787 }
788
789 int PMPI_Comm_free(MPI_Comm * comm)
790 {
791   int retval = 0;
792
793   if (comm == nullptr) {
794     retval = MPI_ERR_ARG;
795   } else if (*comm == MPI_COMM_NULL) {
796     retval = MPI_ERR_COMM;
797   } else {
798     smpi_comm_destroy(*comm);
799     *comm = MPI_COMM_NULL;
800     retval = MPI_SUCCESS;
801   }
802   return retval;
803 }
804
805 int PMPI_Comm_disconnect(MPI_Comm * comm)
806 {
807   /* TODO: wait until all communication in comm are done */
808   int retval = 0;
809
810   if (comm == nullptr) {
811     retval = MPI_ERR_ARG;
812   } else if (*comm == MPI_COMM_NULL) {
813     retval = MPI_ERR_COMM;
814   } else {
815     smpi_comm_destroy(*comm);
816     *comm = MPI_COMM_NULL;
817     retval = MPI_SUCCESS;
818   }
819   return retval;
820 }
821
822 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
823 {
824   int retval = 0;
825   smpi_bench_end();
826
827   if (comm_out == nullptr) {
828     retval = MPI_ERR_ARG;
829   } else if (comm == MPI_COMM_NULL) {
830     retval = MPI_ERR_COMM;
831   } else {
832     *comm_out = smpi_comm_split(comm, color, key);
833     retval = MPI_SUCCESS;
834   }
835   smpi_bench_begin();
836
837   return retval;
838 }
839
840 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
841 {
842   int retval = 0;
843
844   smpi_bench_end();
845   if (request == nullptr) {
846       retval = MPI_ERR_ARG;
847   } else if (comm == MPI_COMM_NULL) {
848       retval = MPI_ERR_COMM;
849   } else if (!is_datatype_valid(datatype)) {
850       retval = MPI_ERR_TYPE;
851   } else if (dst == MPI_PROC_NULL) {
852       retval = MPI_SUCCESS;
853   } else {
854       *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
855       retval = MPI_SUCCESS;
856   }
857   smpi_bench_begin();
858   if (retval != MPI_SUCCESS && request != nullptr)
859     *request = MPI_REQUEST_NULL;
860   return retval;
861 }
862
863 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
864 {
865   int retval = 0;
866
867   smpi_bench_end();
868   if (request == nullptr) {
869     retval = MPI_ERR_ARG;
870   } else if (comm == MPI_COMM_NULL) {
871     retval = MPI_ERR_COMM;
872   } else if (!is_datatype_valid(datatype)) {
873       retval = MPI_ERR_TYPE;
874   } else if (src == MPI_PROC_NULL) {
875     retval = MPI_SUCCESS;
876   } else {
877     *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
878     retval = MPI_SUCCESS;
879   }
880   smpi_bench_begin();
881   if (retval != MPI_SUCCESS && request != nullptr)
882     *request = MPI_REQUEST_NULL;
883   return retval;
884 }
885
886 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
887 {
888   int retval = 0;
889
890   smpi_bench_end();
891   if (request == nullptr) {
892     retval = MPI_ERR_ARG;
893   } else if (comm == MPI_COMM_NULL) {
894     retval = MPI_ERR_COMM;
895   } else if (!is_datatype_valid(datatype)) {
896       retval = MPI_ERR_TYPE;
897   } else if (dst == MPI_PROC_NULL) {
898     retval = MPI_SUCCESS;
899   } else {
900     *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
901     retval = MPI_SUCCESS;
902   }
903   smpi_bench_begin();
904   if (retval != MPI_SUCCESS && request != nullptr)
905     *request = MPI_REQUEST_NULL;
906   return retval;
907 }
908
909 int PMPI_Start(MPI_Request * request)
910 {
911   int retval = 0;
912
913   smpi_bench_end();
914   if (request == nullptr || *request == MPI_REQUEST_NULL) {
915     retval = MPI_ERR_REQUEST;
916   } else {
917     smpi_mpi_start(*request);
918     retval = MPI_SUCCESS;
919   }
920   smpi_bench_begin();
921   return retval;
922 }
923
924 int PMPI_Startall(int count, MPI_Request * requests)
925 {
926   int retval;
927   int i = 0;
928   smpi_bench_end();
929   if (requests == nullptr) {
930     retval = MPI_ERR_ARG;
931   } else {
932     retval = MPI_SUCCESS;
933     for (i = 0 ;  i < count ; i++) {
934       if(requests[i] == MPI_REQUEST_NULL) {
935         retval = MPI_ERR_REQUEST;
936       }
937     }
938     if(retval != MPI_ERR_REQUEST) {
939       smpi_mpi_startall(count, requests);
940     }
941   }
942   smpi_bench_begin();
943   return retval;
944 }
945
946 int PMPI_Request_free(MPI_Request * request)
947 {
948   int retval = 0;
949
950   smpi_bench_end();
951   if (*request == MPI_REQUEST_NULL) {
952     retval = MPI_ERR_ARG;
953   } else {
954     smpi_mpi_request_free(request);
955     retval = MPI_SUCCESS;
956   }
957   smpi_bench_begin();
958   return retval;
959 }
960
961 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
962 {
963   int retval = 0;
964
965   smpi_bench_end();
966
967   if (request == nullptr) {
968     retval = MPI_ERR_ARG;
969   } else if (comm == MPI_COMM_NULL) {
970     retval = MPI_ERR_COMM;
971   } else if (src == MPI_PROC_NULL) {
972     *request = MPI_REQUEST_NULL;
973     retval = MPI_SUCCESS;
974   } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
975     retval = MPI_ERR_RANK;
976   } else if ((count < 0) || (buf==nullptr && count > 0)) {
977     retval = MPI_ERR_COUNT;
978   } else if (!is_datatype_valid(datatype)) {
979       retval = MPI_ERR_TYPE;
980   } else if(tag<0 && tag !=  MPI_ANY_TAG){
981     retval = MPI_ERR_TAG;
982   } else {
983
984     int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
985     int src_traced = smpi_group_index(smpi_comm_group(comm), src);
986
987     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
988     extra->type = TRACING_IRECV;
989     extra->src = src_traced;
990     extra->dst = rank;
991     int known=0;
992     extra->datatype1 = encode_datatype(datatype, &known);
993     int dt_size_send = 1;
994     if(known==0)
995       dt_size_send = smpi_datatype_size(datatype);
996     extra->send_size = count*dt_size_send;
997     TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
998
999     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
1000     retval = MPI_SUCCESS;
1001
1002     TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1003     (*request)->recv = 1;
1004   }
1005
1006   smpi_bench_begin();
1007   if (retval != MPI_SUCCESS && request != nullptr)
1008     *request = MPI_REQUEST_NULL;
1009   return retval;
1010 }
1011
1012
1013 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
1014 {
1015   int retval = 0;
1016
1017   smpi_bench_end();
1018   if (request == nullptr) {
1019     retval = MPI_ERR_ARG;
1020   } else if (comm == MPI_COMM_NULL) {
1021     retval = MPI_ERR_COMM;
1022   } else if (dst == MPI_PROC_NULL) {
1023     *request = MPI_REQUEST_NULL;
1024     retval = MPI_SUCCESS;
1025   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1026     retval = MPI_ERR_RANK;
1027   } else if ((count < 0) || (buf==nullptr && count > 0)) {
1028     retval = MPI_ERR_COUNT;
1029   } else if (!is_datatype_valid(datatype)) {
1030       retval = MPI_ERR_TYPE;
1031   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1032     retval = MPI_ERR_TAG;
1033   } else {
1034
1035     int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1036     int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1037     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1038     extra->type = TRACING_ISEND;
1039     extra->src = rank;
1040     extra->dst = dst_traced;
1041     int known=0;
1042     extra->datatype1 = encode_datatype(datatype, &known);
1043     int dt_size_send = 1;
1044     if(known==0)
1045       dt_size_send = smpi_datatype_size(datatype);
1046     extra->send_size = count*dt_size_send;
1047     TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1048     TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1049
1050     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1051     retval = MPI_SUCCESS;
1052
1053     TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1054     (*request)->send = 1;
1055   }
1056
1057   smpi_bench_begin();
1058   if (retval != MPI_SUCCESS && request!=nullptr)
1059     *request = MPI_REQUEST_NULL;
1060   return retval;
1061 }
1062
1063 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
1064 {
1065   int retval = 0;
1066
1067   smpi_bench_end();
1068   if (request == nullptr) {
1069     retval = MPI_ERR_ARG;
1070   } else if (comm == MPI_COMM_NULL) {
1071     retval = MPI_ERR_COMM;
1072   } else if (dst == MPI_PROC_NULL) {
1073     *request = MPI_REQUEST_NULL;
1074     retval = MPI_SUCCESS;
1075   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1076     retval = MPI_ERR_RANK;
1077   } else if ((count < 0)|| (buf==nullptr && count > 0)) {
1078     retval = MPI_ERR_COUNT;
1079   } else if (!is_datatype_valid(datatype)) {
1080       retval = MPI_ERR_TYPE;
1081   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1082     retval = MPI_ERR_TAG;
1083   } else {
1084
1085     int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1086     int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1087     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1088     extra->type = TRACING_ISSEND;
1089     extra->src = rank;
1090     extra->dst = dst_traced;
1091     int known=0;
1092     extra->datatype1 = encode_datatype(datatype, &known);
1093     int dt_size_send = 1;
1094     if(known==0)
1095       dt_size_send = smpi_datatype_size(datatype);
1096     extra->send_size = count*dt_size_send;
1097     TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1098     TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1099
1100     *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1101     retval = MPI_SUCCESS;
1102
1103     TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1104     (*request)->send = 1;
1105   }
1106
1107   smpi_bench_begin();
1108   if (retval != MPI_SUCCESS && request!=nullptr)
1109     *request = MPI_REQUEST_NULL;
1110   return retval;
1111 }
1112
1113 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
1114 {
1115   int retval = 0;
1116
1117   smpi_bench_end();
1118   if (comm == MPI_COMM_NULL) {
1119     retval = MPI_ERR_COMM;
1120   } else if (src == MPI_PROC_NULL) {
1121     smpi_empty_status(status);
1122     status->MPI_SOURCE = MPI_PROC_NULL;
1123     retval = MPI_SUCCESS;
1124   } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1125     retval = MPI_ERR_RANK;
1126   } else if ((count < 0) || (buf==nullptr && count > 0)) {
1127     retval = MPI_ERR_COUNT;
1128   } else if (!is_datatype_valid(datatype)) {
1129       retval = MPI_ERR_TYPE;
1130   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1131     retval = MPI_ERR_TAG;
1132   } else {
1133   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1134   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1135   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1136   extra->type = TRACING_RECV;
1137   extra->src = src_traced;
1138   extra->dst = rank;
1139   int known=0;
1140   extra->datatype1 = encode_datatype(datatype, &known);
1141   int dt_size_send = 1;
1142   if(known==0)
1143     dt_size_send = smpi_datatype_size(datatype);
1144   extra->send_size = count*dt_size_send;
1145   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
1146
1147     smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1148     retval = MPI_SUCCESS;
1149
1150   //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1151   if(status!=MPI_STATUS_IGNORE){
1152     src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1153     if (!TRACE_smpi_view_internals()) {
1154       TRACE_smpi_recv(rank, src_traced, rank);
1155     }
1156   }
1157   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1158   }
1159
1160   smpi_bench_begin();
1161   return retval;
1162 }
1163
1164 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
1165 {
1166   int retval = 0;
1167
1168   smpi_bench_end();
1169
1170   if (comm == MPI_COMM_NULL) {
1171     retval = MPI_ERR_COMM;
1172   } else if (dst == MPI_PROC_NULL) {
1173     retval = MPI_SUCCESS;
1174   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1175     retval = MPI_ERR_RANK;
1176   } else if ((count < 0) || (buf == nullptr && count > 0)) {
1177     retval = MPI_ERR_COUNT;
1178   } else if (!is_datatype_valid(datatype)) {
1179     retval = MPI_ERR_TYPE;
1180   } else if(tag < 0 && tag !=  MPI_ANY_TAG){
1181     retval = MPI_ERR_TAG;
1182   } else {
1183     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1184     int dst_traced         = smpi_group_index(smpi_comm_group(comm), dst);
1185     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1186     extra->type            = TRACING_SEND;
1187     extra->src             = rank;
1188     extra->dst             = dst_traced;
1189     int known              = 0;
1190     extra->datatype1       = encode_datatype(datatype, &known);
1191     int dt_size_send       = 1;
1192     if (known == 0) {
1193       dt_size_send = smpi_datatype_size(datatype);
1194     }
1195     extra->send_size = count*dt_size_send;
1196     TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1197     if (!TRACE_smpi_view_internals()) {
1198       TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1199     }
1200
1201     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1202     retval = MPI_SUCCESS;
1203
1204     TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1205   }
1206
1207   smpi_bench_begin();
1208   return retval;
1209 }
1210
1211 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1212   int retval = 0;
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) || (buf==nullptr && count > 0)) {
1223     retval = MPI_ERR_COUNT;
1224   } else if (!is_datatype_valid(datatype)){
1225     retval = MPI_ERR_TYPE;
1226   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1227     retval = MPI_ERR_TAG;
1228   } else {
1229     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1230     int dst_traced         = smpi_group_index(smpi_comm_group(comm), dst);
1231     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1232     extra->type            = TRACING_SSEND;
1233     extra->src             = rank;
1234     extra->dst             = dst_traced;
1235     int known              = 0;
1236     extra->datatype1       = encode_datatype(datatype, &known);
1237     int dt_size_send       = 1;
1238     if(known == 0) {
1239       dt_size_send = smpi_datatype_size(datatype);
1240     }
1241     extra->send_size = count*dt_size_send;
1242     TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1243     TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1244   
1245     smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1246     retval = MPI_SUCCESS;
1247   
1248     TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1249   }
1250
1251   smpi_bench_begin();
1252   return retval;
1253 }
1254
1255 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void *recvbuf,
1256                  int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status * status)
1257 {
1258   int retval = 0;
1259
1260   smpi_bench_end();
1261
1262   if (comm == MPI_COMM_NULL) {
1263     retval = MPI_ERR_COMM;
1264   } else if (!is_datatype_valid(sendtype)
1265              || !is_datatype_valid(recvtype)) {
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       (sendbuf==nullptr && sendcount > 0) || (recvbuf==nullptr && recvcount>0)) {
1276     retval = MPI_ERR_COUNT;
1277   } else if((sendtag<0 && sendtag !=  MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1278     retval = MPI_ERR_TAG;
1279   } else {
1280
1281   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1282   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1283   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1284   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1285   extra->type = TRACING_SENDRECV;
1286   extra->src = src_traced;
1287   extra->dst = dst_traced;
1288   int known=0;
1289   extra->datatype1 = encode_datatype(sendtype, &known);
1290   int dt_size_send = 1;
1291   if(known==0)
1292     dt_size_send = smpi_datatype_size(sendtype);
1293   extra->send_size = sendcount*dt_size_send;
1294   extra->datatype2 = encode_datatype(recvtype, &known);
1295   int dt_size_recv = 1;
1296   if(known==0)
1297     dt_size_recv = smpi_datatype_size(recvtype);
1298   extra->recv_size = recvcount*dt_size_recv;
1299
1300   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1301   TRACE_smpi_send(rank, rank, dst_traced,sendcount*smpi_datatype_size(sendtype));
1302
1303     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1304                       recvcount, recvtype, src, recvtag, comm, status);
1305     retval = MPI_SUCCESS;
1306
1307   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1308   TRACE_smpi_recv(rank, src_traced, rank);
1309   }
1310
1311   smpi_bench_begin();
1312   return retval;
1313 }
1314
1315 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
1316                          MPI_Comm comm, MPI_Status * status)
1317 {
1318   void *recvbuf;
1319   int retval = 0;
1320   if (!is_datatype_valid(datatype)) {
1321       retval = MPI_ERR_TYPE;
1322   } else if (count < 0) {
1323       retval = MPI_ERR_COUNT;
1324   } else {
1325     int size = smpi_datatype_get_extent(datatype) * count;
1326     recvbuf = xbt_new0(char, size);
1327     retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
1328     if(retval==MPI_SUCCESS){
1329         smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1330     }
1331     xbt_free(recvbuf);
1332
1333   }
1334   return retval;
1335 }
1336
1337 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1338 {
1339   int retval = 0;
1340   smpi_bench_end();
1341   if (request == nullptr || flag == nullptr) {
1342     retval = MPI_ERR_ARG;
1343   } else if (*request == MPI_REQUEST_NULL) {
1344     *flag= true;
1345     smpi_empty_status(status);
1346     retval = MPI_SUCCESS;
1347   } else {
1348     int rank = (request!=nullptr && (*request)->comm != MPI_COMM_NULL) ? smpi_process_index() : -1;
1349
1350     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1351     extra->type = TRACING_TEST;
1352     TRACE_smpi_testing_in(rank, extra);
1353
1354     *flag = smpi_mpi_test(request, status);
1355
1356     TRACE_smpi_testing_out(rank);
1357     retval = MPI_SUCCESS;
1358   }
1359   smpi_bench_begin();
1360   return retval;
1361 }
1362
1363 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
1364 {
1365   int retval = 0;
1366
1367   smpi_bench_end();
1368   if (index == nullptr || flag == nullptr) {
1369     retval = MPI_ERR_ARG;
1370   } else {
1371     *flag = smpi_mpi_testany(count, requests, index, status);
1372     retval = MPI_SUCCESS;
1373   }
1374   smpi_bench_begin();
1375   return retval;
1376 }
1377
1378 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1379 {
1380   int retval = 0;
1381
1382   smpi_bench_end();
1383   if (flag == nullptr) {
1384     retval = MPI_ERR_ARG;
1385   } else {
1386     *flag = smpi_mpi_testall(count, requests, statuses);
1387     retval = MPI_SUCCESS;
1388   }
1389   smpi_bench_begin();
1390   return retval;
1391 }
1392
1393 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1394   int retval = 0;
1395   smpi_bench_end();
1396
1397   if (status == nullptr) {
1398     retval = MPI_ERR_ARG;
1399   } else if (comm == MPI_COMM_NULL) {
1400     retval = MPI_ERR_COMM;
1401   } else if (source == MPI_PROC_NULL) {
1402     smpi_empty_status(status);
1403     status->MPI_SOURCE = MPI_PROC_NULL;
1404     retval = MPI_SUCCESS;
1405   } else {
1406     smpi_mpi_probe(source, tag, comm, status);
1407     retval = MPI_SUCCESS;
1408   }
1409   smpi_bench_begin();
1410   return retval;
1411 }
1412
1413 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1414   int retval = 0;
1415   smpi_bench_end();
1416
1417   if ((flag == nullptr) || (status == nullptr)) {
1418     retval = MPI_ERR_ARG;
1419   } else if (comm == MPI_COMM_NULL) {
1420     retval = MPI_ERR_COMM;
1421   } else if (source == MPI_PROC_NULL) {
1422     *flag=true;
1423     smpi_empty_status(status);
1424     status->MPI_SOURCE = MPI_PROC_NULL;
1425     retval = MPI_SUCCESS;
1426   } else {
1427     smpi_mpi_iprobe(source, tag, comm, flag, status);
1428     retval = MPI_SUCCESS;
1429   }
1430   smpi_bench_begin();
1431   return retval;
1432 }
1433
1434 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1435 {
1436   int retval = 0;
1437
1438   smpi_bench_end();
1439
1440   smpi_empty_status(status);
1441
1442   if (request == nullptr) {
1443     retval = MPI_ERR_ARG;
1444   } else if (*request == MPI_REQUEST_NULL) {
1445     retval = MPI_SUCCESS;
1446   } else {
1447
1448     int rank = (request!=nullptr && (*request)->comm != MPI_COMM_NULL) ? smpi_process_index() : -1;
1449
1450     int src_traced = (*request)->src;
1451     int dst_traced = (*request)->dst;
1452     MPI_Comm comm = (*request)->comm;
1453     int is_wait_for_receive = (*request)->recv;
1454     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1455     extra->type = TRACING_WAIT;
1456     TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1457
1458     smpi_mpi_wait(request, status);
1459     retval = MPI_SUCCESS;
1460
1461     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1462     TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1463     if (is_wait_for_receive) {
1464       if(src_traced==MPI_ANY_SOURCE)
1465         src_traced = (status!=MPI_STATUS_IGNORE) ?
1466           smpi_group_rank(smpi_comm_group(comm), status->MPI_SOURCE) :
1467           src_traced;
1468       TRACE_smpi_recv(rank, src_traced, dst_traced);
1469     }
1470   }
1471
1472   smpi_bench_begin();
1473   return retval;
1474 }
1475
1476 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1477 {
1478   if (index == nullptr)
1479     return MPI_ERR_ARG;
1480
1481   smpi_bench_end();
1482   //save requests information for tracing
1483   int i;
1484   int *srcs = nullptr, *dsts = nullptr, *recvs = nullptr;
1485   MPI_Comm* comms = nullptr;
1486   if(count>0){
1487     srcs = xbt_new0(int, count);
1488     dsts = xbt_new0(int, count);
1489     recvs = xbt_new0(int, count);
1490     comms = xbt_new0(MPI_Comm, count);
1491   }
1492   for (i = 0; i < count; i++) {
1493     MPI_Request req = requests[i];      //already received requests are no longer valid
1494     if (req) {
1495       srcs[i] = req->src;
1496       dsts[i] = req->dst;
1497       recvs[i] = req->recv;
1498       comms[i] = req->comm;
1499     }
1500   }
1501   int rank_traced = smpi_process_index();
1502   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1503   extra->type = TRACING_WAITANY;
1504   extra->send_size=count;
1505   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1506
1507   *index = smpi_mpi_waitany(count, requests, status);
1508
1509   if(*index!=MPI_UNDEFINED){
1510     int src_traced = srcs[*index];
1511     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1512     int dst_traced = dsts[*index];
1513     int is_wait_for_receive = recvs[*index];
1514     if (is_wait_for_receive) {
1515       if(srcs[*index]==MPI_ANY_SOURCE)
1516         src_traced = (status!=MPI_STATUSES_IGNORE) ?
1517                       smpi_group_rank(smpi_comm_group(comms[*index]), status->MPI_SOURCE) : srcs[*index];
1518       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1519     }
1520     TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1521     xbt_free(srcs);
1522     xbt_free(dsts);
1523     xbt_free(recvs);
1524     xbt_free(comms);
1525
1526   }
1527   smpi_bench_begin();
1528   return MPI_SUCCESS;
1529 }
1530
1531 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1532 {
1533   smpi_bench_end();
1534   //save information from requests
1535   int i;
1536   int *srcs = xbt_new0(int, count);
1537   int *dsts = xbt_new0(int, count);
1538   int *recvs = xbt_new0(int, count);
1539   int *valid = xbt_new0(int, count);
1540   MPI_Comm *comms = xbt_new0(MPI_Comm, count);
1541
1542   for (i = 0; i < count; i++) {
1543     MPI_Request req = requests[i];
1544     if(req!=MPI_REQUEST_NULL){
1545       srcs[i] = req->src;
1546       dsts[i] = req->dst;
1547       recvs[i] = req->recv;
1548       comms[i] = req->comm;
1549       valid[i]=1;;
1550     }else{
1551       valid[i]=0;
1552     }
1553   }
1554   int rank_traced = smpi_process_index();
1555   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1556   extra->type = TRACING_WAITALL;
1557   extra->send_size=count;
1558   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1559
1560   int retval = smpi_mpi_waitall(count, requests, status);
1561
1562   for (i = 0; i < count; i++) {
1563     if(valid[i]){
1564     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1565       int src_traced = srcs[i];
1566       int dst_traced = dsts[i];
1567       int is_wait_for_receive = recvs[i];
1568       if (is_wait_for_receive) {
1569         if(src_traced==MPI_ANY_SOURCE)
1570         src_traced = (status!=MPI_STATUSES_IGNORE) ?
1571                           smpi_group_rank(smpi_comm_group(comms[i]), status[i].MPI_SOURCE) : srcs[i];
1572         TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1573       }
1574     }
1575   }
1576   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1577   xbt_free(srcs);
1578   xbt_free(dsts);
1579   xbt_free(recvs);
1580   xbt_free(valid);
1581   xbt_free(comms);
1582
1583   smpi_bench_begin();
1584   return retval;
1585 }
1586
1587 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
1588 {
1589   int retval = 0;
1590
1591   smpi_bench_end();
1592   if (outcount == nullptr) {
1593     retval = MPI_ERR_ARG;
1594   } else {
1595     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1596     retval = MPI_SUCCESS;
1597   }
1598   smpi_bench_begin();
1599   return retval;
1600 }
1601
1602 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
1603 {
1604   int retval = 0;
1605
1606    smpi_bench_end();
1607    if (outcount == nullptr) {
1608      retval = MPI_ERR_ARG;
1609    } else {
1610      *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1611      retval = MPI_SUCCESS;
1612    }
1613    smpi_bench_begin();
1614    return retval;
1615 }
1616
1617
1618 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1619 {
1620   int retval = 0;
1621
1622   smpi_bench_end();
1623
1624   if (comm == MPI_COMM_NULL) {
1625     retval = MPI_ERR_COMM;
1626   } else if (!is_datatype_valid(datatype)) {
1627       retval = MPI_ERR_ARG;
1628   } else {
1629   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1630   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1631
1632   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1633   extra->type = TRACING_BCAST;
1634   extra->root = root_traced;
1635   int known=0;
1636   extra->datatype1 = encode_datatype(datatype, &known);
1637   int dt_size_send = 1;
1638   if(known==0)
1639     dt_size_send = smpi_datatype_size(datatype);
1640   extra->send_size = count*dt_size_send;
1641   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1642   if(smpi_comm_size(comm)>1)
1643     mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1644   retval = MPI_SUCCESS;
1645
1646   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1647   }
1648
1649   smpi_bench_begin();
1650   return retval;
1651 }
1652
1653 int PMPI_Barrier(MPI_Comm comm)
1654 {
1655   int retval = 0;
1656
1657   smpi_bench_end();
1658
1659   if (comm == MPI_COMM_NULL) {
1660     retval = MPI_ERR_COMM;
1661   } else {
1662   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1663   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1664   extra->type = TRACING_BARRIER;
1665   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1666
1667   mpi_coll_barrier_fun(comm);
1668   retval = MPI_SUCCESS;
1669
1670   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1671   }
1672
1673   smpi_bench_begin();
1674   return retval;
1675 }
1676
1677 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,void *recvbuf, int recvcount, MPI_Datatype recvtype,
1678                 int root, MPI_Comm comm)
1679 {
1680   int retval = 0;
1681
1682   smpi_bench_end();
1683
1684   if (comm == MPI_COMM_NULL) {
1685     retval = MPI_ERR_COMM;
1686   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1687             ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1688     retval = MPI_ERR_TYPE;
1689   } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) || ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1690     retval = MPI_ERR_COUNT;
1691   } else {
1692
1693     char* sendtmpbuf = static_cast<char*>(sendbuf);
1694     int sendtmpcount = sendcount;
1695     MPI_Datatype sendtmptype = sendtype;
1696     if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1697       sendtmpcount=0;
1698       sendtmptype=recvtype;
1699     }
1700   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1701   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1702   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1703   extra->type = TRACING_GATHER;
1704   extra->root = root_traced;
1705   int known=0;
1706   extra->datatype1 = encode_datatype(sendtmptype, &known);
1707   int dt_size_send = 1;
1708   if(known==0)
1709     dt_size_send = smpi_datatype_size(sendtmptype);
1710   extra->send_size = sendtmpcount*dt_size_send;
1711   extra->datatype2 = encode_datatype(recvtype, &known);
1712   int dt_size_recv = 1;
1713   if((smpi_comm_rank(comm)==root) && known==0)
1714     dt_size_recv = smpi_datatype_size(recvtype);
1715   extra->recv_size = recvcount*dt_size_recv;
1716
1717   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1718
1719   mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount, recvtype, root, comm);
1720
1721   retval = MPI_SUCCESS;
1722   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1723   }
1724
1725   smpi_bench_begin();
1726   return retval;
1727 }
1728
1729 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int *recvcounts, int *displs,
1730                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1731 {
1732   int retval = 0;
1733
1734   smpi_bench_end();
1735
1736   if (comm == MPI_COMM_NULL) {
1737     retval = MPI_ERR_COMM;
1738   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1739             ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1740     retval = MPI_ERR_TYPE;
1741   } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1742     retval = MPI_ERR_COUNT;
1743   } else if (recvcounts == nullptr || displs == nullptr) {
1744     retval = MPI_ERR_ARG;
1745   } else {
1746     char* sendtmpbuf = static_cast<char*>(sendbuf);
1747     int sendtmpcount = sendcount;
1748     MPI_Datatype sendtmptype = sendtype;
1749     if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1750       sendtmpcount=0;
1751       sendtmptype=recvtype;
1752     }
1753
1754   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1755   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1756   int i=0;
1757   int size = smpi_comm_size(comm);
1758   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1759   extra->type = TRACING_GATHERV;
1760   extra->num_processes = size;
1761   extra->root = root_traced;
1762   int known=0;
1763   extra->datatype1 = encode_datatype(sendtmptype, &known);
1764   int dt_size_send = 1;
1765   if(known==0)
1766     dt_size_send = smpi_datatype_size(sendtype);
1767   extra->send_size = sendtmpcount*dt_size_send;
1768   extra->datatype2 = encode_datatype(recvtype, &known);
1769   int dt_size_recv = 1;
1770   if(known==0)
1771     dt_size_recv = smpi_datatype_size(recvtype);
1772   if((smpi_comm_rank(comm)==root)){
1773   extra->recvcounts= xbt_new(int,size);
1774   for(i=0; i< size; i++)//copy data to avoid bad free
1775     extra->recvcounts[i] = recvcounts[i]*dt_size_recv;
1776   }
1777   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1778
1779   smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts, displs, recvtype, root, comm);
1780     retval = MPI_SUCCESS;
1781   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1782   }
1783
1784   smpi_bench_begin();
1785   return retval;
1786 }
1787
1788 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1789                    void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
1790 {
1791   int retval = 0;
1792
1793   smpi_bench_end();
1794
1795   if (comm == MPI_COMM_NULL) {
1796     retval = MPI_ERR_COMM;
1797   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1798             (recvtype == MPI_DATATYPE_NULL)){
1799     retval = MPI_ERR_TYPE;
1800   } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1801             (recvcount <0)){
1802     retval = MPI_ERR_COUNT;
1803   } else {
1804     if(sendbuf == MPI_IN_PLACE) {
1805       sendbuf=static_cast<char*>(recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1806       sendcount=recvcount;
1807       sendtype=recvtype;
1808     }
1809   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1810   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1811   extra->type = TRACING_ALLGATHER;
1812   int known=0;
1813   extra->datatype1 = encode_datatype(sendtype, &known);
1814   int dt_size_send = 1;
1815   if(known==0)
1816     dt_size_send = smpi_datatype_size(sendtype);
1817   extra->send_size = sendcount*dt_size_send;
1818   extra->datatype2 = encode_datatype(recvtype, &known);
1819   int dt_size_recv = 1;
1820   if(known==0)
1821     dt_size_recv = smpi_datatype_size(recvtype);
1822   extra->recv_size = recvcount*dt_size_recv;
1823
1824   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1825
1826   mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1827     retval = MPI_SUCCESS;
1828   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1829   }
1830   smpi_bench_begin();
1831   return retval;
1832 }
1833
1834 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1835                    void *recvbuf, int *recvcounts, int *displs, MPI_Datatype recvtype, MPI_Comm comm)
1836 {
1837   int retval = 0;
1838
1839   smpi_bench_end();
1840
1841   if (comm == MPI_COMM_NULL) {
1842     retval = MPI_ERR_COMM;
1843   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1844             (recvtype == MPI_DATATYPE_NULL)){
1845     retval = MPI_ERR_TYPE;
1846   } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1847     retval = MPI_ERR_COUNT;
1848   } else if (recvcounts == nullptr || displs == nullptr) {
1849     retval = MPI_ERR_ARG;
1850   } else {
1851
1852     if(sendbuf == MPI_IN_PLACE) {
1853       sendbuf=static_cast<char*>(recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1854       sendcount=recvcounts[smpi_comm_rank(comm)];
1855       sendtype=recvtype;
1856     }
1857   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1858   int i=0;
1859   int size = smpi_comm_size(comm);
1860   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1861   extra->type = TRACING_ALLGATHERV;
1862   extra->num_processes = size;
1863   int known=0;
1864   extra->datatype1 = encode_datatype(sendtype, &known);
1865   int dt_size_send = 1;
1866   if(known==0)
1867     dt_size_send = smpi_datatype_size(sendtype);
1868   extra->send_size = sendcount*dt_size_send;
1869   extra->datatype2 = encode_datatype(recvtype, &known);
1870   int dt_size_recv = 1;
1871   if(known==0)
1872     dt_size_recv = smpi_datatype_size(recvtype);
1873   extra->recvcounts= xbt_new(int, size);
1874   for(i=0; i< size; i++)//copy data to avoid bad free
1875     extra->recvcounts[i] = recvcounts[i]*dt_size_recv;
1876
1877   TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
1878
1879     mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
1880     retval = MPI_SUCCESS;
1881   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1882   }
1883
1884   smpi_bench_begin();
1885   return retval;
1886 }
1887
1888 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1889                 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
1890 {
1891   int retval = 0;
1892
1893   smpi_bench_end();
1894
1895   if (comm == MPI_COMM_NULL) {
1896     retval = MPI_ERR_COMM;
1897   } else if (((smpi_comm_rank(comm)==root) && (!is_datatype_valid(sendtype)))
1898              || ((recvbuf !=MPI_IN_PLACE) && (!is_datatype_valid(recvtype)))){
1899     retval = MPI_ERR_TYPE;
1900   } else if ((sendbuf == recvbuf) ||
1901       ((smpi_comm_rank(comm)==root) && sendcount>0 && (sendbuf == nullptr))){
1902     retval = MPI_ERR_BUFFER;
1903   }else {
1904
1905     if (recvbuf == MPI_IN_PLACE) {
1906         recvtype=sendtype;
1907         recvcount=sendcount;
1908     }
1909   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1910   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1911   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1912   extra->type = TRACING_SCATTER;
1913   extra->root = root_traced;
1914   int known=0;
1915   extra->datatype1 = encode_datatype(sendtype, &known);
1916   int dt_size_send = 1;
1917   if((smpi_comm_rank(comm)==root) && known==0)
1918     dt_size_send = smpi_datatype_size(sendtype);
1919   extra->send_size = sendcount*dt_size_send;
1920   extra->datatype2 = encode_datatype(recvtype, &known);
1921   int dt_size_recv = 1;
1922   if(known==0)
1923     dt_size_recv = smpi_datatype_size(recvtype);
1924   extra->recv_size = recvcount*dt_size_recv;
1925   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1926
1927   mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
1928     retval = MPI_SUCCESS;
1929   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1930   }
1931
1932   smpi_bench_begin();
1933   return retval;
1934 }
1935
1936 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1937                  MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
1938 {
1939   int retval = 0;
1940
1941   smpi_bench_end();
1942
1943   if (comm == MPI_COMM_NULL) {
1944     retval = MPI_ERR_COMM;
1945   } else if (sendcounts == nullptr || displs == nullptr) {
1946     retval = MPI_ERR_ARG;
1947   } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1948              || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1949     retval = MPI_ERR_TYPE;
1950   } else {
1951     if (recvbuf == MPI_IN_PLACE) {
1952         recvtype=sendtype;
1953         recvcount=sendcounts[smpi_comm_rank(comm)];
1954     }
1955   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1956   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1957   int i=0;
1958   int size = smpi_comm_size(comm);
1959   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1960   extra->type = TRACING_SCATTERV;
1961   extra->num_processes = size;
1962   extra->root = root_traced;
1963   int known=0;
1964   extra->datatype1 = encode_datatype(sendtype, &known);
1965   int dt_size_send = 1;
1966   if(known==0)
1967     dt_size_send = smpi_datatype_size(sendtype);
1968   if((smpi_comm_rank(comm)==root)){
1969   extra->sendcounts= xbt_new(int, size);
1970   for(i=0; i< size; i++)//copy data to avoid bad free
1971     extra->sendcounts[i] = sendcounts[i]*dt_size_send;
1972   }
1973   extra->datatype2 = encode_datatype(recvtype, &known);
1974   int dt_size_recv = 1;
1975   if(known==0)
1976     dt_size_recv = smpi_datatype_size(recvtype);
1977   extra->recv_size = recvcount*dt_size_recv;
1978   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1979
1980     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
1981
1982     retval = MPI_SUCCESS;
1983   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1984   }
1985
1986   smpi_bench_begin();
1987   return retval;
1988 }
1989
1990 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1991 {
1992   int retval = 0;
1993
1994   smpi_bench_end();
1995
1996   if (comm == MPI_COMM_NULL) {
1997     retval = MPI_ERR_COMM;
1998   } else if (!is_datatype_valid(datatype) || op == MPI_OP_NULL) {
1999     retval = MPI_ERR_ARG;
2000   } else {
2001   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2002   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2003   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2004   extra->type = TRACING_REDUCE;
2005   int known=0;
2006   extra->datatype1 = encode_datatype(datatype, &known);
2007   int dt_size_send = 1;
2008   if(known==0)
2009     dt_size_send = smpi_datatype_size(datatype);
2010   extra->send_size = count*dt_size_send;
2011   extra->root = root_traced;
2012
2013   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2014
2015   mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
2016
2017     retval = MPI_SUCCESS;
2018   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2019   }
2020
2021   smpi_bench_begin();
2022   return retval;
2023 }
2024
2025 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count, MPI_Datatype datatype, MPI_Op op){
2026   int retval = 0;
2027
2028     smpi_bench_end();
2029     if (!is_datatype_valid(datatype) || op == MPI_OP_NULL) {
2030       retval = MPI_ERR_ARG;
2031     } else {
2032       smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
2033       retval=MPI_SUCCESS;
2034     }
2035     smpi_bench_begin();
2036     return retval;
2037 }
2038
2039 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2040 {
2041   int retval = 0;
2042
2043   smpi_bench_end();
2044
2045   if (comm == MPI_COMM_NULL) {
2046     retval = MPI_ERR_COMM;
2047   } else if (!is_datatype_valid(datatype)) {
2048     retval = MPI_ERR_TYPE;
2049   } else if (op == MPI_OP_NULL) {
2050     retval = MPI_ERR_OP;
2051   } else {
2052
2053     char* sendtmpbuf = static_cast<char*>(sendbuf);
2054     if( sendbuf == MPI_IN_PLACE ) {
2055       sendtmpbuf = static_cast<char*>(xbt_malloc(count*smpi_datatype_get_extent(datatype)));
2056       smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
2057     }
2058   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2059   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2060   extra->type = TRACING_ALLREDUCE;
2061   int known=0;
2062   extra->datatype1 = encode_datatype(datatype, &known);
2063   int dt_size_send = 1;
2064   if(known==0)
2065     dt_size_send = smpi_datatype_size(datatype);
2066   extra->send_size = count*dt_size_send;
2067
2068   TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2069
2070   mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
2071
2072     if( sendbuf == MPI_IN_PLACE )
2073       xbt_free(sendtmpbuf);
2074
2075     retval = MPI_SUCCESS;
2076   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2077   }
2078
2079   smpi_bench_begin();
2080   return retval;
2081 }
2082
2083 int PMPI_Scan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2084 {
2085   int retval = 0;
2086
2087   smpi_bench_end();
2088
2089   if (comm == MPI_COMM_NULL) {
2090     retval = MPI_ERR_COMM;
2091   } else if (!is_datatype_valid(datatype)) {
2092     retval = MPI_ERR_TYPE;
2093   } else if (op == MPI_OP_NULL) {
2094     retval = MPI_ERR_OP;
2095   } else {
2096   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2097   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2098   extra->type = TRACING_SCAN;
2099   int known=0;
2100   extra->datatype1 = encode_datatype(datatype, &known);
2101   int dt_size_send = 1;
2102   if(known==0)
2103     dt_size_send = smpi_datatype_size(datatype);
2104   extra->send_size = count*dt_size_send;
2105
2106   TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2107
2108   smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2109
2110   retval = MPI_SUCCESS;
2111   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2112   }
2113
2114   smpi_bench_begin();
2115   return retval;
2116 }
2117
2118 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm){
2119   int retval = 0;
2120
2121   smpi_bench_end();
2122
2123   if (comm == MPI_COMM_NULL) {
2124     retval = MPI_ERR_COMM;
2125   } else if (!is_datatype_valid(datatype)) {
2126     retval = MPI_ERR_TYPE;
2127   } else if (op == MPI_OP_NULL) {
2128     retval = MPI_ERR_OP;
2129   } else {
2130   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2131   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2132   extra->type = TRACING_EXSCAN;
2133   int known=0;
2134   extra->datatype1 = encode_datatype(datatype, &known);
2135   int dt_size_send = 1;
2136   if(known==0)
2137     dt_size_send = smpi_datatype_size(datatype);
2138   extra->send_size = count*dt_size_send;
2139   TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2140
2141   smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2142     retval = MPI_SUCCESS;
2143   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2144   }
2145
2146   smpi_bench_begin();
2147   return retval;
2148 }
2149
2150 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2151 {
2152   int retval = 0;
2153   smpi_bench_end();
2154
2155   if (comm == MPI_COMM_NULL) {
2156     retval = MPI_ERR_COMM;
2157   } else if (!is_datatype_valid(datatype)) {
2158     retval = MPI_ERR_TYPE;
2159   } else if (op == MPI_OP_NULL) {
2160     retval = MPI_ERR_OP;
2161   } else if (recvcounts == nullptr) {
2162     retval = MPI_ERR_ARG;
2163   } else {
2164   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2165   int i=0;
2166   int size = smpi_comm_size(comm);
2167   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2168   extra->type = TRACING_REDUCE_SCATTER;
2169   extra->num_processes = size;
2170   int known=0;
2171   extra->datatype1 = encode_datatype(datatype, &known);
2172   int dt_size_send = 1;
2173   if(known==0)
2174     dt_size_send = smpi_datatype_size(datatype);
2175   extra->send_size = 0;
2176   extra->recvcounts= xbt_new(int, size);
2177   for(i=0; i< size; i++)//copy data to avoid bad free
2178     extra->recvcounts[i] = recvcounts[i]*dt_size_send;
2179   TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2180
2181   void* sendtmpbuf=sendbuf;
2182     if(sendbuf==MPI_IN_PLACE)
2183       sendtmpbuf=recvbuf;
2184
2185     mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts, datatype,  op, comm);
2186     retval = MPI_SUCCESS;
2187   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2188   }
2189
2190   smpi_bench_begin();
2191   return retval;
2192 }
2193
2194 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2195                               MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2196 {
2197   int retval,i;
2198   smpi_bench_end();
2199
2200   if (comm == MPI_COMM_NULL) {
2201     retval = MPI_ERR_COMM;
2202   } else if (!is_datatype_valid(datatype)) {
2203     retval = MPI_ERR_TYPE;
2204   } else if (op == MPI_OP_NULL) {
2205     retval = MPI_ERR_OP;
2206   } else if (recvcount < 0) {
2207     retval = MPI_ERR_ARG;
2208   } else {
2209     int count=smpi_comm_size(comm);
2210
2211   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2212   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2213   extra->type = TRACING_REDUCE_SCATTER;
2214   extra->num_processes = count;
2215   int known=0;
2216   extra->datatype1 = encode_datatype(datatype, &known);
2217   int dt_size_send = 1;
2218   if(known==0)
2219     dt_size_send = smpi_datatype_size(datatype);
2220   extra->send_size = 0;
2221   extra->recvcounts= xbt_new(int, count);
2222   for(i=0; i< count; i++)//copy data to avoid bad free
2223     extra->recvcounts[i] = recvcount*dt_size_send;
2224
2225   TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2226
2227   int* recvcounts=static_cast<int*>(xbt_malloc(count));
2228     for (i=0; i<count;i++)
2229       recvcounts[i]=recvcount;
2230     mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts, datatype,  op, comm);
2231     xbt_free(recvcounts);
2232     retval = MPI_SUCCESS;
2233
2234     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2235   }
2236
2237   smpi_bench_begin();
2238   return retval;
2239 }
2240
2241 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2242                  void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
2243 {
2244   int retval = 0;
2245
2246   smpi_bench_end();
2247
2248   if (comm == MPI_COMM_NULL) {
2249     retval = MPI_ERR_COMM;
2250   } else if (sendtype == MPI_DATATYPE_NULL
2251              || recvtype == MPI_DATATYPE_NULL) {
2252     retval = MPI_ERR_TYPE;
2253   } else {
2254   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2255   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2256   extra->type = TRACING_ALLTOALL;
2257   int known=0;
2258   extra->datatype1 = encode_datatype(sendtype, &known);
2259   if(known==0)
2260     extra->send_size = sendcount*smpi_datatype_size(sendtype);
2261   else
2262     extra->send_size = sendcount;
2263   extra->datatype2 = encode_datatype(recvtype, &known);
2264   if(known==0)
2265     extra->recv_size = recvcount*smpi_datatype_size(recvtype);
2266   else
2267     extra->recv_size = recvcount;
2268   TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2269
2270   retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2271
2272   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2273   }
2274
2275   smpi_bench_begin();
2276   return retval;
2277 }
2278
2279 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2280                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2281 {
2282   int retval = 0;
2283
2284   smpi_bench_end();
2285
2286   if (comm == MPI_COMM_NULL) {
2287     retval = MPI_ERR_COMM;
2288   } else if (sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
2289     retval = MPI_ERR_TYPE;
2290   } else if (sendcounts == nullptr || senddisps == nullptr || recvcounts == nullptr || recvdisps == nullptr) {
2291     retval = MPI_ERR_ARG;
2292   } else {
2293   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2294   int i=0;
2295   int size = smpi_comm_size(comm);
2296   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2297   extra->type = TRACING_ALLTOALLV;
2298   extra->send_size = 0;
2299   extra->recv_size = 0;
2300   extra->recvcounts= xbt_new(int, size);
2301   extra->sendcounts= xbt_new(int, size);
2302   int known=0;
2303   extra->datatype1 = encode_datatype(sendtype, &known);
2304   int dt_size_send = 1;
2305   if(known==0)
2306     dt_size_send = smpi_datatype_size(sendtype);
2307   int dt_size_recv = 1;
2308   extra->datatype2 = encode_datatype(recvtype, &known);
2309   if(known==0)
2310     dt_size_recv = smpi_datatype_size(recvtype);
2311   for(i=0; i< size; i++){//copy data to avoid bad free
2312     extra->send_size += sendcounts[i]*dt_size_send;
2313     extra->recv_size += recvcounts[i]*dt_size_recv;
2314
2315     extra->sendcounts[i] = sendcounts[i]*dt_size_send;
2316     extra->recvcounts[i] = recvcounts[i]*dt_size_recv;
2317   }
2318   extra->num_processes = size;
2319   TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2320
2321   retval = mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype, recvbuf, recvcounts, recvdisps, recvtype,
2322                                   comm);
2323   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2324   }
2325
2326   smpi_bench_begin();
2327   return retval;
2328 }
2329
2330
2331 int PMPI_Get_processor_name(char *name, int *resultlen)
2332 {
2333   int retval = MPI_SUCCESS;
2334
2335   strncpy(name, sg_host_get_name(SIMIX_host_self()),
2336           strlen(sg_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2337           strlen(sg_host_get_name(SIMIX_host_self())) +1 : MPI_MAX_PROCESSOR_NAME - 1 );
2338   *resultlen = strlen(name) > MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2339
2340   return retval;
2341 }
2342
2343 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2344 {
2345   int retval = MPI_SUCCESS;
2346   size_t size;
2347
2348   if (status == nullptr || count == nullptr) {
2349     retval = MPI_ERR_ARG;
2350   } else if (!is_datatype_valid(datatype)) {
2351     retval = MPI_ERR_TYPE;
2352   } else {
2353     size = smpi_datatype_size(datatype);
2354     if (size == 0) {
2355       *count = 0;
2356     } else if (status->count % size != 0) {
2357       retval = MPI_UNDEFINED;
2358     } else {
2359       *count = smpi_mpi_get_count(status, datatype);
2360     }
2361   }
2362   return retval;
2363 }
2364
2365 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2366   int retval = 0;
2367
2368   if (old_type == MPI_DATATYPE_NULL) {
2369     retval = MPI_ERR_TYPE;
2370   } else if (count<0){
2371     retval = MPI_ERR_COUNT;
2372   } else {
2373     retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2374   }
2375   return retval;
2376 }
2377
2378 int PMPI_Type_commit(MPI_Datatype* datatype) {
2379   int retval = 0;
2380
2381   if (datatype == nullptr || *datatype == MPI_DATATYPE_NULL) {
2382     retval = MPI_ERR_TYPE;
2383   } else {
2384     smpi_datatype_commit(datatype);
2385     retval = MPI_SUCCESS;
2386   }
2387   return retval;
2388 }
2389
2390 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2391   int retval = 0;
2392
2393   if (old_type == MPI_DATATYPE_NULL) {
2394     retval = MPI_ERR_TYPE;
2395   } else if (count<0 || blocklen<0){
2396     retval = MPI_ERR_COUNT;
2397   } else {
2398     retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2399   }
2400   return retval;
2401 }
2402
2403 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2404   int retval = 0;
2405
2406   if (old_type == MPI_DATATYPE_NULL) {
2407     retval = MPI_ERR_TYPE;
2408   } else if (count<0 || blocklen<0){
2409     retval = MPI_ERR_COUNT;
2410   } else {
2411     retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2412   }
2413   return retval;
2414 }
2415
2416 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2417   return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2418 }
2419
2420 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2421   int retval = 0;
2422
2423   if (old_type == MPI_DATATYPE_NULL) {
2424     retval = MPI_ERR_TYPE;
2425   } else if (count<0){
2426     retval = MPI_ERR_COUNT;
2427   } else {
2428     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2429   }
2430   return retval;
2431 }
2432
2433 int PMPI_Type_create_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2434   int retval = 0;
2435
2436   if (old_type == MPI_DATATYPE_NULL) {
2437     retval = MPI_ERR_TYPE;
2438   } else if (count<0){
2439     retval = MPI_ERR_COUNT;
2440   } else {
2441     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2442   }
2443   return retval;
2444 }
2445
2446 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type,
2447                                    MPI_Datatype* new_type) {
2448   int retval,i;
2449
2450   if (old_type == MPI_DATATYPE_NULL) {
2451     retval = MPI_ERR_TYPE;
2452   } else if (count<0){
2453     retval = MPI_ERR_COUNT;
2454   } else {
2455     int* blocklens=static_cast<int*>(xbt_malloc(blocklength*count));
2456     for (i=0; i<count;i++)
2457       blocklens[i]=blocklength;
2458     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2459     xbt_free(blocklens);
2460   }
2461   return retval;
2462 }
2463
2464 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2465   int retval = 0;
2466
2467   if (old_type == MPI_DATATYPE_NULL) {
2468     retval = MPI_ERR_TYPE;
2469   } else if (count<0){
2470     retval = MPI_ERR_COUNT;
2471   } else {
2472     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2473   }
2474   return retval;
2475 }
2476
2477 int PMPI_Type_create_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type,
2478                               MPI_Datatype* new_type) {
2479   return PMPI_Type_hindexed(count, blocklens,indices,old_type,new_type);
2480 }
2481
2482 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type,
2483                                     MPI_Datatype* new_type) {
2484   int retval,i;
2485
2486   if (old_type == MPI_DATATYPE_NULL) {
2487     retval = MPI_ERR_TYPE;
2488   } else if (count<0){
2489     retval = MPI_ERR_COUNT;
2490   } else {
2491     int* blocklens=(int*)xbt_malloc(blocklength*count);
2492     for (i=0; i<count;i++)blocklens[i]=blocklength;
2493     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2494     xbt_free(blocklens);
2495   }
2496   return retval;
2497 }
2498
2499 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2500   int retval = 0;
2501
2502   if (count<0){
2503     retval = MPI_ERR_COUNT;
2504   } else {
2505     retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2506   }
2507   return retval;
2508 }
2509
2510 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types,
2511                             MPI_Datatype* new_type) {
2512   return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2513 }
2514
2515 int PMPI_Error_class(int errorcode, int* errorclass) {
2516   // assume smpi uses only standard mpi error codes
2517   *errorclass=errorcode;
2518   return MPI_SUCCESS;
2519 }
2520
2521 int PMPI_Initialized(int* flag) {
2522    *flag=smpi_process_initialized();
2523    return MPI_SUCCESS;
2524 }
2525
2526 /* The topo part of MPI_COMM_WORLD should always be nullptr. When other topologies will be implemented, not only should we
2527  * check if the topology is nullptr, but we should check if it is the good topology type (so we have to add a
2528  *  MPIR_Topo_Type field, and replace the MPI_Topology field by an union)*/
2529
2530 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periodic, int reorder, MPI_Comm* comm_cart) {
2531   int retval = 0;
2532   if (comm_old == MPI_COMM_NULL){
2533     retval =  MPI_ERR_COMM;
2534   } else if (ndims < 0 || (ndims > 0 && (dims == nullptr || periodic == nullptr)) || comm_cart == nullptr) {
2535     retval = MPI_ERR_ARG;
2536   } else{
2537     retval = smpi_mpi_cart_create(comm_old, ndims, dims, periodic, reorder, comm_cart);
2538   }
2539   return retval;
2540 }
2541
2542 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2543   if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == nullptr) {
2544     return MPI_ERR_TOPOLOGY;
2545   }
2546   if (coords == nullptr) {
2547     return MPI_ERR_ARG;
2548   }
2549   return smpi_mpi_cart_rank(comm, coords, rank);
2550 }
2551
2552 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2553   if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == nullptr) {
2554     return MPI_ERR_TOPOLOGY;
2555   }
2556   if (source == nullptr || dest == nullptr || direction < 0 ) {
2557     return MPI_ERR_ARG;
2558   }
2559   return smpi_mpi_cart_shift(comm, direction, displ, source, dest);
2560 }
2561
2562 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2563   if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == nullptr) {
2564     return MPI_ERR_TOPOLOGY;
2565   }
2566   if (rank < 0 || rank >= smpi_comm_size(comm)) {
2567     return MPI_ERR_RANK;
2568   }
2569   if (maxdims <= 0) {
2570     return MPI_ERR_ARG;
2571   }
2572   if(coords == nullptr) {
2573     return MPI_ERR_ARG;
2574   }
2575   return smpi_mpi_cart_coords(comm, rank, maxdims, coords);
2576 }
2577
2578 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2579   if(comm == nullptr || smpi_comm_topo(comm) == nullptr) {
2580     return MPI_ERR_TOPOLOGY;
2581   }
2582   if(maxdims <= 0 || dims == nullptr || periods == nullptr || coords == nullptr) {
2583     return MPI_ERR_ARG;
2584   }
2585   return smpi_mpi_cart_get(comm, maxdims, dims, periods, coords);
2586 }
2587
2588 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2589   if (comm == MPI_COMM_NULL || smpi_comm_topo(comm) == nullptr) {
2590     return MPI_ERR_TOPOLOGY;
2591   }
2592   if (ndims == nullptr) {
2593     return MPI_ERR_ARG;
2594   }
2595   return smpi_mpi_cartdim_get(comm, ndims);
2596 }
2597
2598 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2599   if(dims == nullptr) {
2600     return MPI_ERR_ARG;
2601   }
2602   if (ndims < 1 || nnodes < 1) {
2603     return MPI_ERR_DIMS;
2604   }
2605
2606   return smpi_mpi_dims_create(nnodes, ndims, dims);
2607 }
2608
2609 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2610   if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == nullptr) {
2611     return MPI_ERR_TOPOLOGY;
2612   }
2613   if (comm_new == nullptr) {
2614     return MPI_ERR_ARG;
2615   }
2616   return smpi_mpi_cart_sub(comm, remain_dims, comm_new);
2617 }
2618
2619 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2620     if(oldtype == MPI_DATATYPE_NULL) {
2621         return MPI_ERR_TYPE;
2622     }
2623     int blocks[3] = { 1, 1, 1 };
2624     MPI_Aint disps[3] = { lb, 0, lb+extent };
2625     MPI_Datatype types[3] = { MPI_LB, oldtype, MPI_UB };
2626
2627     s_smpi_mpi_struct_t* subtype = smpi_datatype_struct_create( blocks, disps, 3, types);
2628     smpi_datatype_create(newtype,oldtype->size, lb, lb + extent, sizeof(s_smpi_mpi_struct_t) , subtype, DT_FLAG_VECTOR);
2629
2630     (*newtype)->flags &= ~DT_FLAG_COMMITED;
2631     return MPI_SUCCESS;
2632 }
2633
2634 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2635   int retval = 0;
2636   smpi_bench_end();
2637   if (comm == MPI_COMM_NULL) {
2638     retval= MPI_ERR_COMM;
2639   }else if ((base == nullptr && size != 0) || disp_unit <= 0 || size < 0 ){
2640     retval= MPI_ERR_OTHER;
2641   }else{
2642     *win = smpi_mpi_win_create( base, size, disp_unit, info, comm);
2643     retval = MPI_SUCCESS;
2644   }
2645   smpi_bench_begin();
2646   return retval;
2647 }
2648
2649 int PMPI_Win_free( MPI_Win* win){
2650   int retval = 0;
2651   smpi_bench_end();
2652   if (win == nullptr || *win == MPI_WIN_NULL) {
2653     retval = MPI_ERR_WIN;
2654   }else{
2655     retval=smpi_mpi_win_free(win);
2656   }
2657   smpi_bench_begin();
2658   return retval;
2659 }
2660
2661 int PMPI_Win_set_name(MPI_Win  win, char * name)
2662 {
2663   int retval = 0;
2664   if (win == MPI_WIN_NULL)  {
2665     retval = MPI_ERR_TYPE;
2666   } else if (name == nullptr)  {
2667     retval = MPI_ERR_ARG;
2668   } else {
2669     smpi_mpi_win_set_name(win, name);
2670     retval = MPI_SUCCESS;
2671   }
2672   return retval;
2673 }
2674
2675 int PMPI_Win_get_name(MPI_Win  win, char * name, int* len)
2676 {
2677   int retval = MPI_SUCCESS;
2678
2679   if (win == MPI_WIN_NULL)  {
2680     retval = MPI_ERR_WIN;
2681   } else if (name == nullptr)  {
2682     retval = MPI_ERR_ARG;
2683   } else {
2684     smpi_mpi_win_get_name(win, name, len);
2685   }
2686   return retval;
2687 }
2688
2689 int PMPI_Win_get_group(MPI_Win  win, MPI_Group * group){
2690   int retval = MPI_SUCCESS;
2691   if (win == MPI_WIN_NULL)  {
2692     retval = MPI_ERR_WIN;
2693   }else {
2694     smpi_mpi_win_get_group(win, group);
2695     smpi_group_use(*group);
2696   }
2697   return retval;
2698 }
2699
2700
2701 int PMPI_Win_fence( int assert,  MPI_Win win){
2702   int retval = 0;
2703   smpi_bench_end();
2704   if (win == MPI_WIN_NULL) {
2705     retval = MPI_ERR_WIN;
2706   } else {
2707   int rank = smpi_process_index();
2708   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, nullptr);
2709   retval = smpi_mpi_win_fence(assert, win);
2710   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2711   }
2712   smpi_bench_begin();
2713   return retval;
2714 }
2715
2716 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2717               MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2718   int retval = 0;
2719   smpi_bench_end();
2720   if (win == MPI_WIN_NULL) {
2721     retval = MPI_ERR_WIN;
2722   } else if (target_rank == MPI_PROC_NULL) {
2723     retval = MPI_SUCCESS;
2724   } else if (target_rank <0){
2725     retval = MPI_ERR_RANK;
2726   } else if (target_disp <0){
2727       retval = MPI_ERR_ARG;
2728   } else if ((origin_count < 0 || target_count < 0) ||
2729              (origin_addr==nullptr && origin_count > 0)){
2730     retval = MPI_ERR_COUNT;
2731   } else if ((!is_datatype_valid(origin_datatype)) || (!is_datatype_valid(target_datatype))) {
2732     retval = MPI_ERR_TYPE;
2733   } else {
2734     int rank = smpi_process_index();
2735     MPI_Group group;
2736     smpi_mpi_win_get_group(win, &group);
2737     int src_traced = smpi_group_index(group, target_rank);
2738     TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, nullptr);
2739
2740     retval = smpi_mpi_get( origin_addr, origin_count, origin_datatype, target_rank, target_disp, target_count,
2741                            target_datatype, win);
2742
2743     TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
2744   }
2745   smpi_bench_begin();
2746   return retval;
2747 }
2748
2749 int PMPI_Put( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2750               MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2751   int retval = 0;
2752   smpi_bench_end();
2753   if (win == MPI_WIN_NULL) {
2754     retval = MPI_ERR_WIN;
2755   } else if (target_rank == MPI_PROC_NULL) {
2756     retval = MPI_SUCCESS;
2757   } else if (target_rank <0){
2758     retval = MPI_ERR_RANK;
2759   } else if (target_disp <0){
2760     retval = MPI_ERR_ARG;
2761   } else if ((origin_count < 0 || target_count < 0) ||
2762             (origin_addr==nullptr && origin_count > 0)){
2763     retval = MPI_ERR_COUNT;
2764   } else if ((!is_datatype_valid(origin_datatype)) || (!is_datatype_valid(target_datatype))) {
2765     retval = MPI_ERR_TYPE;
2766   } else {
2767     int rank = smpi_process_index();
2768     MPI_Group group;
2769     smpi_mpi_win_get_group(win, &group);
2770     int dst_traced = smpi_group_index(group, target_rank);
2771     TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, nullptr);
2772     TRACE_smpi_send(rank, rank, dst_traced, origin_count*smpi_datatype_size(origin_datatype));
2773
2774     retval = smpi_mpi_put( origin_addr, origin_count, origin_datatype, target_rank, target_disp, target_count,
2775                            target_datatype, win);
2776
2777     TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
2778   }
2779   smpi_bench_begin();
2780   return retval;
2781 }
2782
2783 int PMPI_Accumulate( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2784               MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Op op, MPI_Win win){
2785   int retval = 0;
2786   smpi_bench_end();
2787   if (win == MPI_WIN_NULL) {
2788     retval = MPI_ERR_WIN;
2789   } else if (target_rank == MPI_PROC_NULL) {
2790     retval = MPI_SUCCESS;
2791   } else if (target_rank <0){
2792     retval = MPI_ERR_RANK;
2793   } else if (target_disp <0){
2794     retval = MPI_ERR_ARG;
2795   } else if ((origin_count < 0 || target_count < 0) ||
2796              (origin_addr==nullptr && origin_count > 0)){
2797     retval = MPI_ERR_COUNT;
2798   } else if ((!is_datatype_valid(origin_datatype)) ||
2799             (!is_datatype_valid(target_datatype))) {
2800     retval = MPI_ERR_TYPE;
2801   } else if (op == MPI_OP_NULL) {
2802     retval = MPI_ERR_OP;
2803   } else {
2804     int rank = smpi_process_index();
2805     MPI_Group group;
2806     smpi_mpi_win_get_group(win, &group);
2807     int src_traced = smpi_group_index(group, target_rank);
2808     TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, nullptr);
2809
2810     retval = smpi_mpi_accumulate( origin_addr, origin_count, origin_datatype, target_rank, target_disp, target_count,
2811                                   target_datatype, op, win);
2812
2813     TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
2814   }
2815   smpi_bench_begin();
2816   return retval;
2817 }
2818
2819 int PMPI_Win_post(MPI_Group group, int assert, MPI_Win win){
2820   int retval = 0;
2821   smpi_bench_end();
2822   if (win == MPI_WIN_NULL) {
2823     retval = MPI_ERR_WIN;
2824   } else if (group==MPI_GROUP_NULL){
2825     retval = MPI_ERR_GROUP;
2826   }
2827   else {
2828     int rank = smpi_process_index();
2829     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, nullptr);
2830     retval = smpi_mpi_win_post(group,assert,win);
2831     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2832   }
2833   smpi_bench_begin();
2834   return retval;
2835 }
2836
2837 int PMPI_Win_start(MPI_Group group, int assert, MPI_Win win){
2838   int retval = 0;
2839   smpi_bench_end();
2840   if (win == MPI_WIN_NULL) {
2841     retval = MPI_ERR_WIN;
2842   } else if (group==MPI_GROUP_NULL){
2843     retval = MPI_ERR_GROUP;
2844   }
2845   else {
2846     int rank = smpi_process_index();
2847     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, nullptr);
2848     retval = smpi_mpi_win_start(group,assert,win);
2849     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2850   }
2851   smpi_bench_begin();
2852   return retval;
2853 }
2854
2855 int PMPI_Win_complete(MPI_Win win){
2856   int retval = 0;
2857   smpi_bench_end();
2858   if (win == MPI_WIN_NULL) {
2859     retval = MPI_ERR_WIN;
2860   }
2861   else {
2862     int rank = smpi_process_index();
2863     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, nullptr);
2864
2865     retval = smpi_mpi_win_complete(win);
2866
2867     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2868   }
2869   smpi_bench_begin();
2870   return retval;
2871 }
2872
2873 int PMPI_Win_wait(MPI_Win win){
2874   int retval = 0;
2875   smpi_bench_end();
2876   if (win == MPI_WIN_NULL) {
2877     retval = MPI_ERR_WIN;
2878   }
2879   else {
2880     int rank = smpi_process_index();
2881     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, nullptr);
2882
2883     retval = smpi_mpi_win_wait(win);
2884
2885     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2886   }
2887   smpi_bench_begin();
2888   return retval;
2889 }
2890
2891 int PMPI_Alloc_mem(MPI_Aint size, MPI_Info info, void *baseptr){
2892   void *ptr = xbt_malloc(size);
2893   if(ptr==nullptr)
2894     return MPI_ERR_NO_MEM;
2895   else {
2896     *static_cast<void**>(baseptr) = ptr;
2897     return MPI_SUCCESS;
2898   }
2899 }
2900
2901 int PMPI_Free_mem(void *baseptr){
2902   xbt_free(baseptr);
2903   return MPI_SUCCESS;
2904 }
2905
2906 int PMPI_Type_set_name(MPI_Datatype  datatype, char * name)
2907 {
2908   int retval = 0;
2909   if (datatype == MPI_DATATYPE_NULL)  {
2910     retval = MPI_ERR_TYPE;
2911   } else if (name == nullptr)  {
2912     retval = MPI_ERR_ARG;
2913   } else {
2914     smpi_datatype_set_name(datatype, name);
2915     retval = MPI_SUCCESS;
2916   }
2917   return retval;
2918 }
2919
2920 int PMPI_Type_get_name(MPI_Datatype  datatype, char * name, int* len)
2921 {
2922   int retval = 0;
2923
2924   if (datatype == MPI_DATATYPE_NULL)  {
2925     retval = MPI_ERR_TYPE;
2926   } else if (name == nullptr)  {
2927     retval = MPI_ERR_ARG;
2928   } else {
2929     smpi_datatype_get_name(datatype, name, len);
2930     retval = MPI_SUCCESS;
2931   }
2932   return retval;
2933 }
2934
2935 MPI_Datatype PMPI_Type_f2c(MPI_Fint datatype){
2936   return smpi_type_f2c(datatype);
2937 }
2938
2939 MPI_Fint PMPI_Type_c2f(MPI_Datatype datatype){
2940   return smpi_type_c2f( datatype);
2941 }
2942
2943 MPI_Group PMPI_Group_f2c(MPI_Fint group){
2944   return smpi_group_f2c( group);
2945 }
2946
2947 MPI_Fint PMPI_Group_c2f(MPI_Group group){
2948   return smpi_group_c2f(group);
2949 }
2950
2951 MPI_Request PMPI_Request_f2c(MPI_Fint request){
2952   return smpi_request_f2c(request);
2953 }
2954
2955 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
2956   return smpi_request_c2f(request);
2957 }
2958
2959 MPI_Win PMPI_Win_f2c(MPI_Fint win){
2960   return smpi_win_f2c(win);
2961 }
2962
2963 MPI_Fint PMPI_Win_c2f(MPI_Win win){
2964   return smpi_win_c2f(win);
2965 }
2966
2967 MPI_Op PMPI_Op_f2c(MPI_Fint op){
2968   return smpi_op_f2c(op);
2969 }
2970
2971 MPI_Fint PMPI_Op_c2f(MPI_Op op){
2972   return smpi_op_c2f(op);
2973 }
2974
2975 MPI_Comm PMPI_Comm_f2c(MPI_Fint comm){
2976   return smpi_comm_f2c(comm);
2977 }
2978
2979 MPI_Fint PMPI_Comm_c2f(MPI_Comm comm){
2980   return smpi_comm_c2f(comm);
2981 }
2982
2983 MPI_Info PMPI_Info_f2c(MPI_Fint info){
2984   return smpi_info_f2c(info);
2985 }
2986
2987 MPI_Fint PMPI_Info_c2f(MPI_Info info){
2988   return smpi_info_c2f(info);
2989 }
2990
2991 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2992   return smpi_comm_keyval_create(copy_fn, delete_fn, keyval, extra_state);
2993 }
2994
2995 int PMPI_Keyval_free(int* keyval) {
2996   return smpi_comm_keyval_free(keyval);
2997 }
2998
2999 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
3000   if(keyval == MPI_TAG_UB||keyval == MPI_HOST||keyval == MPI_IO ||keyval == MPI_WTIME_IS_GLOBAL||keyval == MPI_APPNUM
3001        ||keyval == MPI_UNIVERSE_SIZE||keyval == MPI_LASTUSEDCODE)
3002     return MPI_ERR_ARG;
3003   else if (comm==MPI_COMM_NULL)
3004     return MPI_ERR_COMM;
3005   else
3006     return smpi_comm_attr_delete(comm, keyval);
3007 }
3008
3009 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
3010   static int one = 1;
3011   static int zero = 0;
3012   static int tag_ub = 1000000;
3013   static int last_used_code = MPI_ERR_LASTCODE;
3014
3015   if (comm==MPI_COMM_NULL){
3016     *flag = 0;
3017     return MPI_ERR_COMM;
3018   }
3019
3020   switch (keyval) {
3021   case MPI_HOST:
3022   case MPI_IO:
3023   case MPI_APPNUM:
3024     *flag = 1;
3025     *static_cast<int**>(attr_value) = &zero;
3026     return MPI_SUCCESS;
3027   case MPI_UNIVERSE_SIZE:
3028     *flag = 1;
3029     *static_cast<int**>(attr_value) = &smpi_universe_size;
3030     return MPI_SUCCESS;
3031   case MPI_LASTUSEDCODE:
3032     *flag = 1;
3033     *static_cast<int**>(attr_value) = &last_used_code;
3034     return MPI_SUCCESS;
3035   case MPI_TAG_UB:
3036     *flag=1;
3037     *static_cast<int**>(attr_value) = &tag_ub;
3038     return MPI_SUCCESS;
3039   case MPI_WTIME_IS_GLOBAL:
3040     *flag = 1;
3041     *static_cast<int**>(attr_value) = &one;
3042     return MPI_SUCCESS;
3043   default:
3044     return smpi_comm_attr_get(comm, keyval, attr_value, flag);
3045   }
3046 }
3047
3048 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
3049   if(keyval == MPI_TAG_UB||keyval == MPI_HOST||keyval == MPI_IO ||keyval == MPI_WTIME_IS_GLOBAL||keyval == MPI_APPNUM
3050        ||keyval == MPI_UNIVERSE_SIZE||keyval == MPI_LASTUSEDCODE)
3051     return MPI_ERR_ARG;
3052   else if (comm==MPI_COMM_NULL)
3053     return MPI_ERR_COMM;
3054   else
3055   return smpi_comm_attr_put(comm, keyval, attr_value);
3056 }
3057
3058 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
3059 {
3060   return PMPI_Attr_get(comm, comm_keyval, attribute_val,flag);
3061 }
3062
3063 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
3064 {
3065   return PMPI_Attr_put(comm, comm_keyval, attribute_val);
3066 }
3067
3068 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
3069 {
3070   return PMPI_Attr_delete(comm, comm_keyval);
3071 }
3072
3073 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval,
3074                             void* extra_state)
3075 {
3076   return PMPI_Keyval_create(copy_fn, delete_fn, keyval, extra_state);
3077 }
3078
3079 int PMPI_Comm_free_keyval(int* keyval) {
3080   return PMPI_Keyval_free(keyval);
3081 }
3082
3083 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
3084 {
3085   if (type==MPI_DATATYPE_NULL)
3086     return MPI_ERR_TYPE;
3087   else
3088     return smpi_type_attr_get(type, type_keyval, attribute_val, flag);
3089 }
3090
3091 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
3092 {
3093   if (type==MPI_DATATYPE_NULL)
3094     return MPI_ERR_TYPE;
3095   else
3096     return smpi_type_attr_put(type, type_keyval, attribute_val);
3097 }
3098
3099 int PMPI_Type_delete_attr (MPI_Datatype type, int type_keyval)
3100 {
3101   if (type==MPI_DATATYPE_NULL)
3102     return MPI_ERR_TYPE;
3103   else
3104     return smpi_type_attr_delete(type, type_keyval);
3105 }
3106
3107 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval,
3108                             void* extra_state)
3109 {
3110   return smpi_type_keyval_create(copy_fn, delete_fn, keyval, extra_state);
3111 }
3112
3113 int PMPI_Type_free_keyval(int* keyval) {
3114   return smpi_type_keyval_free(keyval);
3115 }
3116
3117 int PMPI_Info_create( MPI_Info *info){
3118   if (info == nullptr)
3119     return MPI_ERR_ARG;
3120   *info = xbt_new(s_smpi_mpi_info_t, 1);
3121   (*info)->info_dict= xbt_dict_new_homogeneous(nullptr);
3122   (*info)->refcount=1;
3123   return MPI_SUCCESS;
3124 }
3125
3126 int PMPI_Info_set( MPI_Info info, char *key, char *value){
3127   if (info == nullptr || key == nullptr || value == nullptr)
3128     return MPI_ERR_ARG;
3129
3130   xbt_dict_set(info->info_dict, key, (void*)value, nullptr);
3131   return MPI_SUCCESS;
3132 }
3133
3134 int PMPI_Info_free( MPI_Info *info){
3135   if (info == nullptr || *info==nullptr)
3136     return MPI_ERR_ARG;
3137   (*info)->refcount--;
3138   if((*info)->refcount==0){
3139     xbt_dict_free(&((*info)->info_dict));
3140     xbt_free(*info);
3141   }
3142   *info=MPI_INFO_NULL;
3143   return MPI_SUCCESS;
3144 }
3145
3146 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
3147   *flag=false;
3148   if (info == nullptr || key == nullptr || valuelen <0)
3149     return MPI_ERR_ARG;
3150   if (value == nullptr)
3151     return MPI_ERR_INFO_VALUE;
3152   char* tmpvalue=static_cast<char*>(xbt_dict_get_or_null(info->info_dict, key));
3153   if(tmpvalue){
3154     memset(value, 0, valuelen);
3155     memcpy(value,tmpvalue, (strlen(tmpvalue) + 1 < static_cast<size_t>(valuelen)) ? strlen(tmpvalue) + 1 : valuelen);
3156     *flag=true;
3157   }
3158   return MPI_SUCCESS;
3159 }
3160
3161 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
3162   if (info == nullptr || newinfo==nullptr)
3163     return MPI_ERR_ARG;
3164   *newinfo = xbt_new(s_smpi_mpi_info_t, 1);
3165   (*newinfo)->info_dict= xbt_dict_new_homogeneous(nullptr);
3166   (*newinfo)->refcount=1;
3167   xbt_dict_cursor_t cursor = nullptr;
3168   int *key;
3169   void* data;
3170   xbt_dict_foreach(info->info_dict,cursor,key,data){
3171     xbt_dict_set((*newinfo)->info_dict, reinterpret_cast<char*>(key), data, nullptr);
3172   }
3173   return MPI_SUCCESS;
3174 }
3175
3176 int PMPI_Info_delete(MPI_Info info, char *key){
3177   xbt_ex_t e;
3178   if (info == nullptr || key==nullptr)
3179     return MPI_ERR_ARG;
3180   TRY {
3181     xbt_dict_remove(info->info_dict, key);
3182   }CATCH(e){
3183     xbt_ex_free(e);
3184     return MPI_ERR_INFO_NOKEY;
3185   }
3186   return MPI_SUCCESS;
3187 }
3188
3189 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
3190   if (info == nullptr || nkeys==nullptr)
3191     return MPI_ERR_ARG;
3192   *nkeys=xbt_dict_size(info->info_dict);
3193   return MPI_SUCCESS;
3194 }
3195
3196 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
3197   if (info == nullptr || key==nullptr || n<0 || n> MPI_MAX_INFO_KEY)
3198     return MPI_ERR_ARG;
3199
3200   xbt_dict_cursor_t cursor = nullptr;
3201   char *keyn;
3202   void* data;
3203   int num=0;
3204   xbt_dict_foreach(info->info_dict,cursor,keyn,data){
3205     if(num==n){
3206       strncpy(key,keyn,strlen(keyn)+1);
3207       xbt_dict_cursor_free(&cursor);
3208       return MPI_SUCCESS;
3209     }
3210     num++;
3211   }
3212   return MPI_ERR_ARG;
3213 }
3214
3215 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
3216   *flag=false;
3217   if (info == nullptr || key == nullptr || valuelen==nullptr)
3218     return MPI_ERR_ARG;
3219   char* tmpvalue=(char*)xbt_dict_get_or_null(info->info_dict, key);
3220   if(tmpvalue){
3221     *valuelen=strlen(tmpvalue);
3222     *flag=true;
3223   }
3224   return MPI_SUCCESS;
3225 }
3226
3227 int PMPI_Unpack(void* inbuf, int incount, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
3228   if(incount<0 || outcount < 0 || inbuf==nullptr || outbuf==nullptr)
3229     return MPI_ERR_ARG;
3230   if(!is_datatype_valid(type))
3231     return MPI_ERR_TYPE;
3232   if(comm==MPI_COMM_NULL)
3233     return MPI_ERR_COMM;
3234   return smpi_mpi_unpack(inbuf, incount, position, outbuf,outcount,type, comm);
3235 }
3236
3237 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
3238   if(incount<0 || outcount < 0|| inbuf==nullptr || outbuf==nullptr)
3239     return MPI_ERR_ARG;
3240   if(!is_datatype_valid(type))
3241     return MPI_ERR_TYPE;
3242   if(comm==MPI_COMM_NULL)
3243     return MPI_ERR_COMM;
3244   return smpi_mpi_pack(inbuf, incount, type, outbuf,outcount,position, comm);
3245 }
3246
3247 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
3248   if(incount<0)
3249     return MPI_ERR_ARG;
3250   if(!is_datatype_valid(datatype))
3251     return MPI_ERR_TYPE;
3252   if(comm==MPI_COMM_NULL)
3253     return MPI_ERR_COMM;
3254
3255   *size=incount*smpi_datatype_size(datatype);
3256
3257   return MPI_SUCCESS;
3258 }
3259