Logo AND Algorithmique Numérique Distribuée

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