Logo AND Algorithmique Numérique Distribuée

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