Logo AND Algorithmique Numérique Distribuée

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