Logo AND Algorithmique Numérique Distribuée

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