Logo AND Algorithmique Numérique Distribuée

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