Logo AND Algorithmique Numérique Distribuée

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