Logo AND Algorithmique Numérique Distribuée

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