Logo AND Algorithmique Numérique Distribuée

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