Logo AND Algorithmique Numérique Distribuée

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