Logo AND Algorithmique Numérique Distribuée

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