Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Encode path so that it is not too long.
[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 == NULL || flag == NULL) {
1378     retval = MPI_ERR_ARG;
1379   } else if (*request == MPI_REQUEST_NULL) {
1380     *flag= TRUE;
1381     smpi_empty_status(status);
1382     retval = MPI_ERR_REQUEST;
1383   } else {
1384     *flag = smpi_mpi_test(request, status);
1385     retval = MPI_SUCCESS;
1386   }
1387   smpi_bench_begin();
1388   return retval;
1389 }
1390
1391 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1392                 MPI_Status * status)
1393 {
1394   int retval = 0;
1395
1396   smpi_bench_end();
1397   if (index == NULL || flag == NULL) {
1398     retval = MPI_ERR_ARG;
1399   } else {
1400     *flag = smpi_mpi_testany(count, requests, index, status);
1401     retval = MPI_SUCCESS;
1402   }
1403   smpi_bench_begin();
1404   return retval;
1405 }
1406
1407 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1408 {
1409   int retval = 0;
1410
1411   smpi_bench_end();
1412   if (flag == NULL) {
1413     retval = MPI_ERR_ARG;
1414   } else {
1415     *flag = smpi_mpi_testall(count, requests, statuses);
1416     retval = MPI_SUCCESS;
1417   }
1418   smpi_bench_begin();
1419   return retval;
1420 }
1421
1422 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1423   int retval = 0;
1424   smpi_bench_end();
1425
1426   if (status == NULL) {
1427     retval = MPI_ERR_ARG;
1428   } else if (comm == MPI_COMM_NULL) {
1429     retval = MPI_ERR_COMM;
1430   } else if (source == MPI_PROC_NULL) {
1431     smpi_empty_status(status);
1432     status->MPI_SOURCE = MPI_PROC_NULL;
1433     retval = MPI_SUCCESS;
1434   } else {
1435     smpi_mpi_probe(source, tag, comm, status);
1436     retval = MPI_SUCCESS;
1437   }
1438   smpi_bench_begin();
1439   return retval;
1440 }
1441
1442
1443 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1444   int retval = 0;
1445   smpi_bench_end();
1446
1447   if (flag == NULL) {
1448     retval = MPI_ERR_ARG;
1449   } else if (status == NULL) {
1450     retval = MPI_ERR_ARG;
1451   } else if (comm == MPI_COMM_NULL) {
1452     retval = MPI_ERR_COMM;
1453   } else if (source == MPI_PROC_NULL) {
1454     *flag=TRUE;
1455     smpi_empty_status(status);
1456     status->MPI_SOURCE = MPI_PROC_NULL;
1457     retval = MPI_SUCCESS;
1458   } else {
1459     smpi_mpi_iprobe(source, tag, comm, flag, status);
1460     retval = MPI_SUCCESS;
1461   }
1462   smpi_bench_begin();
1463   return retval;
1464 }
1465
1466 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1467 {
1468   int retval = 0;
1469
1470   smpi_bench_end();
1471
1472   smpi_empty_status(status);
1473
1474   if (request == NULL) {
1475     retval = MPI_ERR_ARG;
1476   } else if (*request == MPI_REQUEST_NULL) {
1477     retval = MPI_ERR_REQUEST;
1478   } else {
1479
1480 #ifdef HAVE_TRACING
1481     int rank = request && (*request)->comm != MPI_COMM_NULL
1482       ? smpi_process_index()
1483       : -1;
1484     TRACE_smpi_computing_out(rank);
1485
1486     int src_traced = (*request)->src;
1487     int dst_traced = (*request)->dst;
1488     MPI_Comm comm = (*request)->comm;
1489     int is_wait_for_receive = (*request)->recv;
1490     TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__,-1);
1491 #endif
1492
1493     smpi_mpi_wait(request, status);
1494     retval = MPI_SUCCESS;
1495
1496 #ifdef HAVE_TRACING
1497     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1498     TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1499     if (is_wait_for_receive) {
1500       if(src_traced==MPI_ANY_SOURCE)
1501         src_traced = (status!=MPI_STATUS_IGNORE) ?
1502           smpi_group_rank(smpi_comm_group(comm), status->MPI_SOURCE) :
1503           src_traced;
1504       TRACE_smpi_recv(rank, src_traced, dst_traced);
1505     }
1506     TRACE_smpi_computing_in(rank);
1507 #endif
1508
1509   }
1510
1511   smpi_bench_begin();
1512   return retval;
1513 }
1514
1515 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1516 {
1517   int retval = 0;
1518
1519   smpi_bench_end();
1520 #ifdef HAVE_TRACING
1521   //save requests information for tracing
1522   int i;
1523   int *srcs = xbt_new(int, count);
1524   int *dsts = xbt_new(int, count);
1525   int *recvs = xbt_new(int, count);
1526   MPI_Comm *comms = xbt_new(MPI_Comm, count);
1527
1528   for (i = 0; i < count; i++) {
1529     MPI_Request req = requests[i];      //already received requests are no longer valid
1530     if (req) {
1531       srcs[i] = req->src;
1532       dsts[i] = req->dst;
1533       recvs[i] = req->recv;
1534       comms[i] = req->comm;
1535     }
1536   }
1537   int rank_traced = smpi_process_index();
1538   TRACE_smpi_computing_out(rank_traced);
1539
1540   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,count);
1541
1542 #endif
1543   if (index == NULL) {
1544     retval = MPI_ERR_ARG;
1545   } else {
1546     *index = smpi_mpi_waitany(count, requests, status);
1547     retval = MPI_SUCCESS;
1548   }
1549 #ifdef HAVE_TRACING
1550   if(*index!=MPI_UNDEFINED){
1551     int src_traced = srcs[*index];
1552     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1553     int dst_traced = dsts[*index];
1554     int is_wait_for_receive = recvs[*index];
1555     if (is_wait_for_receive) {
1556       if(srcs[*index]==MPI_ANY_SOURCE)
1557         src_traced = (status!=MPI_STATUSES_IGNORE) ?
1558                       smpi_group_rank(smpi_comm_group(comms[*index]), status->MPI_SOURCE) :
1559                       srcs[*index];
1560       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1561     }
1562     TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1563     xbt_free(srcs);
1564     xbt_free(dsts);
1565     xbt_free(recvs);
1566     xbt_free(comms);
1567
1568   }
1569   TRACE_smpi_computing_in(rank_traced);
1570 #endif
1571   smpi_bench_begin();
1572   return retval;
1573 }
1574
1575 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1576 {
1577
1578   smpi_bench_end();
1579 #ifdef HAVE_TRACING
1580   //save information from requests
1581   int i;
1582   int *srcs = xbt_new(int, count);
1583   int *dsts = xbt_new(int, count);
1584   int *recvs = xbt_new(int, count);
1585   int *valid = xbt_new(int, count);
1586   MPI_Comm *comms = xbt_new(MPI_Comm, count);
1587
1588   //int valid_count = 0;
1589   for (i = 0; i < count; i++) {
1590     MPI_Request req = requests[i];
1591     if(req!=MPI_REQUEST_NULL){
1592       srcs[i] = req->src;
1593       dsts[i] = req->dst;
1594       recvs[i] = req->recv;
1595       comms[i] = req->comm;
1596       valid[i]=1;;
1597     }else{
1598       valid[i]=0;
1599     }
1600   }
1601   int rank_traced = smpi_process_index();
1602   TRACE_smpi_computing_out(rank_traced);
1603
1604   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,count);
1605 #endif
1606   int retval = smpi_mpi_waitall(count, requests, status);
1607 #ifdef HAVE_TRACING
1608   for (i = 0; i < count; i++) {
1609     if(valid[i]){
1610     //int src_traced = srcs[*index];
1611     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1612       int src_traced = srcs[i];
1613       int dst_traced = dsts[i];
1614       int is_wait_for_receive = recvs[i];
1615       if (is_wait_for_receive) {
1616         if(src_traced==MPI_ANY_SOURCE)
1617         src_traced = (status!=MPI_STATUSES_IGNORE) ?
1618                           smpi_group_rank(smpi_comm_group(comms[i]), status[i].MPI_SOURCE) :
1619                           srcs[i];
1620         TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1621       }
1622     }
1623   }
1624   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1625   xbt_free(srcs);
1626   xbt_free(dsts);
1627   xbt_free(recvs);
1628   xbt_free(valid);
1629   xbt_free(comms);
1630
1631   TRACE_smpi_computing_in(rank_traced);
1632 #endif
1633   smpi_bench_begin();
1634   return retval;
1635 }
1636
1637 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1638                  int *indices, MPI_Status status[])
1639 {
1640   int retval = 0;
1641
1642   smpi_bench_end();
1643   if (outcount == NULL) {
1644     retval = MPI_ERR_ARG;
1645   } else {
1646     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1647     retval = MPI_SUCCESS;
1648   }
1649   smpi_bench_begin();
1650   return retval;
1651 }
1652
1653 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1654                  int* indices, MPI_Status status[])
1655 {
1656   int retval = 0;
1657
1658    smpi_bench_end();
1659    if (outcount == NULL) {
1660      retval = MPI_ERR_ARG;
1661    } else {
1662      *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1663      retval = MPI_SUCCESS;
1664    }
1665    smpi_bench_begin();
1666    return retval;
1667 }
1668
1669
1670 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1671 {
1672   int retval = 0;
1673
1674   smpi_bench_end();
1675
1676   if (comm == MPI_COMM_NULL) {
1677     retval = MPI_ERR_COMM;
1678   } else {
1679 #ifdef HAVE_TRACING
1680   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1681   TRACE_smpi_computing_out(rank);
1682   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1683   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
1684 #endif
1685     mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1686     retval = MPI_SUCCESS;
1687 #ifdef HAVE_TRACING
1688   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1689   TRACE_smpi_computing_in(rank);
1690 #endif
1691   }
1692
1693   smpi_bench_begin();
1694   return retval;
1695 }
1696
1697 int PMPI_Barrier(MPI_Comm comm)
1698 {
1699   int retval = 0;
1700
1701   smpi_bench_end();
1702
1703   if (comm == MPI_COMM_NULL) {
1704     retval = MPI_ERR_COMM;
1705   } else {
1706 #ifdef HAVE_TRACING
1707   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1708   TRACE_smpi_computing_out(rank);
1709   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, smpi_comm_size(comm));
1710 #endif
1711     mpi_coll_barrier_fun(comm);
1712     retval = MPI_SUCCESS;
1713 #ifdef HAVE_TRACING
1714   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1715   TRACE_smpi_computing_in(rank);
1716 #endif
1717   }
1718
1719   smpi_bench_begin();
1720   return retval;
1721 }
1722
1723 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1724                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1725                int root, MPI_Comm comm)
1726 {
1727   int retval = 0;
1728
1729   smpi_bench_end();
1730
1731   if (comm == MPI_COMM_NULL) {
1732     retval = MPI_ERR_COMM;
1733   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1734             ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1735     retval = MPI_ERR_TYPE;
1736   } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1737             ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1738     retval = MPI_ERR_COUNT;
1739   } else {
1740
1741     char* sendtmpbuf = (char*) sendbuf;
1742     int sendtmpcount = sendcount;
1743     MPI_Datatype sendtmptype = sendtype;
1744     if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1745       sendtmpcount=0;
1746       sendtmptype=recvtype;
1747     }
1748 #ifdef HAVE_TRACING
1749   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1750   TRACE_smpi_computing_out(rank);
1751   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1752   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,sendcount*smpi_datatype_size(sendtmptype));
1753 #endif
1754     mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount,
1755                     recvtype, root, comm);
1756
1757
1758     retval = MPI_SUCCESS;
1759 #ifdef HAVE_TRACING
1760   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1761   TRACE_smpi_computing_in(rank);
1762 #endif
1763   }
1764
1765   smpi_bench_begin();
1766   return retval;
1767 }
1768
1769 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1770                 void *recvbuf, int *recvcounts, int *displs,
1771                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1772 {
1773   int retval = 0;
1774
1775   smpi_bench_end();
1776
1777   if (comm == MPI_COMM_NULL) {
1778     retval = MPI_ERR_COMM;
1779   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1780             ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1781     retval = MPI_ERR_TYPE;
1782   } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1783     retval = MPI_ERR_COUNT;
1784   } else if (recvcounts == NULL || displs == NULL) {
1785     retval = MPI_ERR_ARG;
1786   } else {
1787     char* sendtmpbuf = (char*) sendbuf;
1788     int sendtmpcount = sendcount;
1789     MPI_Datatype sendtmptype = sendtype;
1790     if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1791       sendtmpcount=0;
1792       sendtmptype=recvtype;
1793     }
1794
1795 #ifdef HAVE_TRACING
1796   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1797   TRACE_smpi_computing_out(rank);
1798   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1799   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,sendcount*smpi_datatype_size(sendtmptype));
1800 #endif
1801     smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts,
1802                      displs, recvtype, root, comm);
1803     retval = MPI_SUCCESS;
1804 #ifdef HAVE_TRACING
1805   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1806   TRACE_smpi_computing_in(rank);
1807 #endif
1808   }
1809
1810   smpi_bench_begin();
1811   return retval;
1812 }
1813
1814 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1815                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1816                   MPI_Comm comm)
1817 {
1818   int retval = 0;
1819
1820   smpi_bench_end();
1821
1822   if (comm == MPI_COMM_NULL) {
1823     retval = MPI_ERR_COMM;
1824   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1825             (recvtype == MPI_DATATYPE_NULL)){
1826     retval = MPI_ERR_TYPE;
1827   } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1828             (recvcount <0)){
1829     retval = MPI_ERR_COUNT;
1830   } else {
1831     if(sendbuf == MPI_IN_PLACE) {
1832       sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1833       sendcount=recvcount;
1834       sendtype=recvtype;
1835     }
1836 #ifdef HAVE_TRACING
1837   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1838   TRACE_smpi_computing_out(rank);
1839   TRACE_smpi_collective_in(rank, -1, __FUNCTION__,sendcount*smpi_datatype_size(sendtype));
1840 #endif
1841     mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1842                            recvtype, comm);
1843     retval = MPI_SUCCESS;
1844
1845 #ifdef HAVE_TRACING
1846   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1847 #endif
1848   }
1849   smpi_bench_begin();
1850   return retval;
1851 }
1852
1853 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1854                    void *recvbuf, int *recvcounts, int *displs,
1855                    MPI_Datatype recvtype, MPI_Comm comm)
1856 {
1857   int retval = 0;
1858
1859   smpi_bench_end();
1860
1861   if (comm == MPI_COMM_NULL) {
1862     retval = MPI_ERR_COMM;
1863   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1864             (recvtype == MPI_DATATYPE_NULL)){
1865     retval = MPI_ERR_TYPE;
1866   } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1867     retval = MPI_ERR_COUNT;
1868   } else if (recvcounts == NULL || displs == NULL) {
1869     retval = MPI_ERR_ARG;
1870   } else {
1871
1872     if(sendbuf == MPI_IN_PLACE) {
1873       sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1874       sendcount=recvcounts[smpi_comm_rank(comm)];
1875       sendtype=recvtype;
1876     }
1877 #ifdef HAVE_TRACING
1878   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1879   TRACE_smpi_computing_out(rank);
1880   TRACE_smpi_collective_in(rank, -1, __FUNCTION__,sendcount*smpi_datatype_size(sendtype));
1881 #endif
1882     mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1883                         displs, recvtype, comm);
1884     retval = MPI_SUCCESS;
1885 #ifdef HAVE_TRACING
1886   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1887   TRACE_smpi_computing_in(rank);
1888 #endif
1889   }
1890
1891   smpi_bench_begin();
1892   return retval;
1893 }
1894
1895 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1896                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1897                 int root, MPI_Comm comm)
1898 {
1899   int retval = 0;
1900
1901   smpi_bench_end();
1902
1903   if (comm == MPI_COMM_NULL) {
1904     retval = MPI_ERR_COMM;
1905   } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1906              || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1907     retval = MPI_ERR_TYPE;
1908   } else {
1909
1910     if (recvbuf == MPI_IN_PLACE) {
1911         recvtype=sendtype;
1912         recvcount=sendcount;
1913     }
1914 #ifdef HAVE_TRACING
1915   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1916   TRACE_smpi_computing_out(rank);
1917   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1918
1919   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,sendcount*smpi_datatype_size(recvtype));
1920 #endif
1921     mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1922                      recvtype, root, comm);
1923     retval = MPI_SUCCESS;
1924 #ifdef HAVE_TRACING
1925   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1926   TRACE_smpi_computing_in(rank);
1927 #endif
1928   }
1929
1930   smpi_bench_begin();
1931   return retval;
1932 }
1933
1934 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1935                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1936                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1937 {
1938   int retval = 0;
1939
1940   smpi_bench_end();
1941
1942   if (comm == MPI_COMM_NULL) {
1943     retval = MPI_ERR_COMM;
1944   } else if (sendcounts == NULL || displs == NULL) {
1945     retval = MPI_ERR_ARG;
1946   } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1947              || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1948     retval = MPI_ERR_TYPE;
1949   } else {
1950     if (recvbuf == MPI_IN_PLACE) {
1951         recvtype=sendtype;
1952         recvcount=sendcounts[smpi_comm_rank(comm)];
1953     }
1954 #ifdef HAVE_TRACING
1955   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1956   TRACE_smpi_computing_out(rank);
1957   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1958   int count=0, i;
1959   for(i=0; i<smpi_comm_size(comm);i++)count+=sendcounts[i];
1960   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, count*smpi_datatype_size(sendtype));
1961 #endif
1962     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1963                       recvcount, recvtype, root, comm);
1964     retval = MPI_SUCCESS;
1965 #ifdef HAVE_TRACING
1966   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1967   TRACE_smpi_computing_in(rank);
1968 #endif
1969   }
1970
1971   smpi_bench_begin();
1972   return retval;
1973 }
1974
1975 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1976                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1977 {
1978   int retval = 0;
1979
1980   smpi_bench_end();
1981
1982   if (comm == MPI_COMM_NULL) {
1983     retval = MPI_ERR_COMM;
1984   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1985     retval = MPI_ERR_ARG;
1986   } else {
1987 #ifdef HAVE_TRACING
1988   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1989   TRACE_smpi_computing_out(rank);
1990   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1991   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
1992 #endif
1993     mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
1994
1995     retval = MPI_SUCCESS;
1996 #ifdef HAVE_TRACING
1997   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1998   TRACE_smpi_computing_in(rank);
1999 #endif
2000   }
2001
2002   smpi_bench_begin();
2003   return retval;
2004 }
2005
2006 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
2007     MPI_Datatype datatype, MPI_Op op){
2008   int retval = 0;
2009
2010     smpi_bench_end();
2011     if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
2012       retval = MPI_ERR_ARG;
2013     } else {
2014       smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
2015       retval=MPI_SUCCESS;
2016     }
2017     smpi_bench_begin();
2018     return retval;
2019 }
2020
2021 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
2022                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2023 {
2024   int retval = 0;
2025
2026   smpi_bench_end();
2027
2028   if (comm == MPI_COMM_NULL) {
2029     retval = MPI_ERR_COMM;
2030   } else if (datatype == MPI_DATATYPE_NULL) {
2031     retval = MPI_ERR_TYPE;
2032   } else if (op == MPI_OP_NULL) {
2033     retval = MPI_ERR_OP;
2034   } else {
2035
2036     char* sendtmpbuf = (char*) sendbuf;
2037     if( sendbuf == MPI_IN_PLACE ) {
2038       sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
2039       smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
2040     }
2041 #ifdef HAVE_TRACING
2042   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2043   TRACE_smpi_computing_out(rank);
2044   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, count*smpi_datatype_size(datatype));
2045 #endif
2046     mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
2047
2048     if( sendbuf == MPI_IN_PLACE ) {
2049       xbt_free(sendtmpbuf);
2050     }
2051
2052     retval = MPI_SUCCESS;
2053 #ifdef HAVE_TRACING
2054   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2055   TRACE_smpi_computing_in(rank);
2056 #endif
2057   }
2058
2059   smpi_bench_begin();
2060   return retval;
2061 }
2062
2063 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
2064              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2065 {
2066   int retval = 0;
2067
2068   smpi_bench_end();
2069
2070   if (comm == MPI_COMM_NULL) {
2071     retval = MPI_ERR_COMM;
2072   } else if (datatype == MPI_DATATYPE_NULL) {
2073     retval = MPI_ERR_TYPE;
2074   } else if (op == MPI_OP_NULL) {
2075     retval = MPI_ERR_OP;
2076   } else {
2077 #ifdef HAVE_TRACING
2078   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2079   TRACE_smpi_computing_out(rank);
2080   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, count*smpi_datatype_size(datatype));
2081 #endif
2082     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2083     retval = MPI_SUCCESS;
2084 #ifdef HAVE_TRACING
2085   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2086   TRACE_smpi_computing_in(rank);
2087 #endif
2088   }
2089
2090   smpi_bench_begin();
2091   return retval;
2092 }
2093
2094 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2095                 MPI_Op op, MPI_Comm comm){
2096   int retval = 0;
2097
2098   smpi_bench_end();
2099
2100   if (comm == MPI_COMM_NULL) {
2101     retval = MPI_ERR_COMM;
2102   } else if (datatype == MPI_DATATYPE_NULL) {
2103     retval = MPI_ERR_TYPE;
2104   } else if (op == MPI_OP_NULL) {
2105     retval = MPI_ERR_OP;
2106   } else {
2107 #ifdef HAVE_TRACING
2108   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2109   TRACE_smpi_computing_out(rank);
2110   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, count*smpi_datatype_size(datatype));
2111 #endif
2112     smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2113     retval = MPI_SUCCESS;
2114 #ifdef HAVE_TRACING
2115   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2116   TRACE_smpi_computing_in(rank);
2117 #endif
2118   }
2119
2120   smpi_bench_begin();
2121   return retval;
2122 }
2123
2124 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
2125                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2126 {
2127   int retval = 0;
2128   smpi_bench_end();
2129
2130   if (comm == MPI_COMM_NULL) {
2131     retval = MPI_ERR_COMM;
2132   } else if (datatype == MPI_DATATYPE_NULL) {
2133     retval = MPI_ERR_TYPE;
2134   } else if (op == MPI_OP_NULL) {
2135     retval = MPI_ERR_OP;
2136   } else if (recvcounts == NULL) {
2137     retval = MPI_ERR_ARG;
2138   } else {
2139 #ifdef HAVE_TRACING
2140   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2141   TRACE_smpi_computing_out(rank);
2142   int count=0, i;
2143   for(i=0; i<smpi_comm_size(comm);i++)count+=recvcounts[i];
2144   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, count*smpi_datatype_size(datatype));
2145 #endif
2146     void* sendtmpbuf=sendbuf;
2147     if(sendbuf==MPI_IN_PLACE){
2148       sendtmpbuf=recvbuf;
2149     }
2150
2151     mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts,
2152                        datatype,  op, comm);
2153     retval = MPI_SUCCESS;
2154 #ifdef HAVE_TRACING
2155   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2156   TRACE_smpi_computing_in(rank);
2157 #endif
2158   }
2159
2160   smpi_bench_begin();
2161   return retval;
2162 }
2163
2164 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2165                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2166 {
2167   int retval,i;
2168   smpi_bench_end();
2169
2170   if (comm == MPI_COMM_NULL) {
2171     retval = MPI_ERR_COMM;
2172   } else if (datatype == MPI_DATATYPE_NULL) {
2173     retval = MPI_ERR_TYPE;
2174   } else if (op == MPI_OP_NULL) {
2175     retval = MPI_ERR_OP;
2176   } else if (recvcount < 0) {
2177     retval = MPI_ERR_ARG;
2178   } else {
2179 #ifdef HAVE_TRACING
2180   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2181   TRACE_smpi_computing_out(rank);
2182   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, recvcount*smpi_comm_size(comm)*smpi_datatype_size(datatype));
2183 #endif
2184     int count=smpi_comm_size(comm);
2185     int* recvcounts=(int*)xbt_malloc(count);
2186     for (i=0; i<count;i++)recvcounts[i]=recvcount;
2187     mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
2188                        datatype,  op, comm);
2189     xbt_free(recvcounts);
2190     retval = MPI_SUCCESS;
2191 #ifdef HAVE_TRACING
2192   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2193   TRACE_smpi_computing_in(rank);
2194 #endif
2195   }
2196
2197   smpi_bench_begin();
2198   return retval;
2199 }
2200
2201 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2202                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
2203                  MPI_Comm comm)
2204 {
2205   int retval = 0;
2206
2207   smpi_bench_end();
2208
2209   if (comm == MPI_COMM_NULL) {
2210     retval = MPI_ERR_COMM;
2211   } else if (sendtype == MPI_DATATYPE_NULL
2212              || recvtype == MPI_DATATYPE_NULL) {
2213     retval = MPI_ERR_TYPE;
2214   } else {
2215 #ifdef HAVE_TRACING
2216   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2217   TRACE_smpi_computing_out(rank);
2218   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, sendcount*smpi_datatype_size(sendtype));
2219 #endif
2220     retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2221 #ifdef HAVE_TRACING
2222   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2223   TRACE_smpi_computing_in(rank);
2224 #endif
2225   }
2226
2227   smpi_bench_begin();
2228   return retval;
2229 }
2230
2231 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2232                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2233                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2234 {
2235   int retval = 0;
2236
2237   smpi_bench_end();
2238
2239   if (comm == MPI_COMM_NULL) {
2240     retval = MPI_ERR_COMM;
2241   } else if (sendtype == MPI_DATATYPE_NULL
2242              || recvtype == MPI_DATATYPE_NULL) {
2243     retval = MPI_ERR_TYPE;
2244   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2245              || recvdisps == NULL) {
2246     retval = MPI_ERR_ARG;
2247   } else {
2248 #ifdef HAVE_TRACING
2249   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2250   TRACE_smpi_computing_out(rank);
2251   int i, size=0;
2252   for(i=0; i< smpi_comm_size(comm);i++)size+=sendcounts[i];
2253   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, size*smpi_datatype_size(sendtype));
2254 #endif
2255     retval =
2256         mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2257                                   recvbuf, recvcounts, recvdisps, recvtype,
2258                                   comm);
2259 #ifdef HAVE_TRACING
2260   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2261   TRACE_smpi_computing_in(rank);
2262 #endif
2263   }
2264
2265   smpi_bench_begin();
2266   return retval;
2267 }
2268
2269
2270 int PMPI_Get_processor_name(char *name, int *resultlen)
2271 {
2272   int retval = MPI_SUCCESS;
2273
2274   smpi_bench_end();
2275   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2276           strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2277           strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2278           MPI_MAX_PROCESSOR_NAME - 1 );
2279   *resultlen =
2280       strlen(name) >
2281       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2282
2283   smpi_bench_begin();
2284   return retval;
2285 }
2286
2287 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2288 {
2289   int retval = MPI_SUCCESS;
2290   size_t size;
2291
2292   smpi_bench_end();
2293   if (status == NULL || count == NULL) {
2294     retval = MPI_ERR_ARG;
2295   } else if (datatype == MPI_DATATYPE_NULL) {
2296     retval = MPI_ERR_TYPE;
2297   } else {
2298     size = smpi_datatype_size(datatype);
2299     if (size == 0) {
2300       *count = 0;
2301     } else if (status->count % size != 0) {
2302       retval = MPI_UNDEFINED;
2303     } else {
2304       *count = smpi_mpi_get_count(status, datatype);
2305     }
2306   }
2307   smpi_bench_begin();
2308   return retval;
2309 }
2310
2311 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2312   int retval = 0;
2313
2314   smpi_bench_end();
2315   if (old_type == MPI_DATATYPE_NULL) {
2316     retval = MPI_ERR_TYPE;
2317   } else if (count<0){
2318     retval = MPI_ERR_COUNT;
2319   } else {
2320     retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2321   }
2322   smpi_bench_begin();
2323   return retval;
2324 }
2325
2326 int PMPI_Type_commit(MPI_Datatype* datatype) {
2327   int retval = 0;
2328
2329   smpi_bench_end();
2330   if (datatype == NULL || *datatype == MPI_DATATYPE_NULL) {
2331     retval = MPI_ERR_TYPE;
2332   } else {
2333     smpi_datatype_commit(datatype);
2334     retval = MPI_SUCCESS;
2335   }
2336   smpi_bench_begin();
2337   return retval;
2338 }
2339
2340
2341 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2342   int retval = 0;
2343
2344   smpi_bench_end();
2345   if (old_type == MPI_DATATYPE_NULL) {
2346     retval = MPI_ERR_TYPE;
2347   } else if (count<0 || blocklen<0){
2348     retval = MPI_ERR_COUNT;
2349   } else {
2350     retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2351   }
2352   smpi_bench_begin();
2353   return retval;
2354 }
2355
2356 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2357   int retval = 0;
2358
2359   smpi_bench_end();
2360   if (old_type == MPI_DATATYPE_NULL) {
2361     retval = MPI_ERR_TYPE;
2362   } else if (count<0 || blocklen<0){
2363     retval = MPI_ERR_COUNT;
2364   } else {
2365     retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2366   }
2367   smpi_bench_begin();
2368   return retval;
2369 }
2370
2371 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2372   return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2373 }
2374
2375 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2376   int retval = 0;
2377
2378   smpi_bench_end();
2379   if (old_type == MPI_DATATYPE_NULL) {
2380     retval = MPI_ERR_TYPE;
2381   } else if (count<0){
2382     retval = MPI_ERR_COUNT;
2383   } else {
2384     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2385   }
2386   smpi_bench_begin();
2387   return retval;
2388 }
2389
2390 int PMPI_Type_create_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2391   int retval = 0;
2392
2393   smpi_bench_end();
2394   if (old_type == MPI_DATATYPE_NULL) {
2395     retval = MPI_ERR_TYPE;
2396   } else if (count<0){
2397     retval = MPI_ERR_COUNT;
2398   } else {
2399     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2400   }
2401   smpi_bench_begin();
2402   return retval;
2403 }
2404
2405 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2406   int retval,i;
2407
2408   smpi_bench_end();
2409   if (old_type == MPI_DATATYPE_NULL) {
2410     retval = MPI_ERR_TYPE;
2411   } else if (count<0){
2412     retval = MPI_ERR_COUNT;
2413   } else {
2414     int* blocklens=(int*)xbt_malloc(blocklength*count);
2415     for (i=0; i<count;i++)blocklens[i]=blocklength;
2416     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2417     xbt_free(blocklens);
2418   }
2419   smpi_bench_begin();
2420   return retval;
2421 }
2422
2423
2424 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2425   int retval = 0;
2426
2427   smpi_bench_end();
2428   if (old_type == MPI_DATATYPE_NULL) {
2429     retval = MPI_ERR_TYPE;
2430   } else if (count<0){
2431     retval = MPI_ERR_COUNT;
2432   } else {
2433     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2434   }
2435   smpi_bench_begin();
2436   return retval;
2437 }
2438
2439 int PMPI_Type_create_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2440   return PMPI_Type_hindexed(count, blocklens,indices,old_type,new_type);
2441 }
2442
2443 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2444   int retval,i;
2445
2446   smpi_bench_end();
2447   if (old_type == MPI_DATATYPE_NULL) {
2448     retval = MPI_ERR_TYPE;
2449   } else if (count<0){
2450     retval = MPI_ERR_COUNT;
2451   } else {
2452     int* blocklens=(int*)xbt_malloc(blocklength*count);
2453     for (i=0; i<count;i++)blocklens[i]=blocklength;
2454     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2455     xbt_free(blocklens);
2456   }
2457   smpi_bench_begin();
2458   return retval;
2459 }
2460
2461
2462 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2463   int retval = 0;
2464
2465   smpi_bench_end();
2466   if (count<0){
2467     retval = MPI_ERR_COUNT;
2468   } else {
2469     retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2470   }
2471   smpi_bench_begin();
2472   return retval;
2473 }
2474
2475 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2476   return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2477 }
2478
2479
2480 int PMPI_Error_class(int errorcode, int* errorclass) {
2481   // assume smpi uses only standard mpi error codes
2482   *errorclass=errorcode;
2483   return MPI_SUCCESS;
2484 }
2485
2486
2487 int PMPI_Initialized(int* flag) {
2488    *flag=smpi_process_initialized();
2489    return MPI_SUCCESS;
2490 }
2491
2492 /* The following calls are not yet implemented and will fail at runtime. */
2493 /* Once implemented, please move them above this notice. */
2494
2495 #define NOT_YET_IMPLEMENTED {\
2496         XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__);\
2497         return MPI_SUCCESS;\
2498         }
2499
2500
2501 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2502   NOT_YET_IMPLEMENTED
2503 }
2504
2505 int PMPI_Type_set_name(MPI_Datatype  datatype, char * name)
2506 {
2507   NOT_YET_IMPLEMENTED
2508 }
2509
2510 int PMPI_Type_get_name(MPI_Datatype  datatype, char * name, int* len)
2511 {
2512   NOT_YET_IMPLEMENTED
2513 }
2514
2515 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2516    NOT_YET_IMPLEMENTED
2517 }
2518
2519 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2520    NOT_YET_IMPLEMENTED
2521 }
2522
2523 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
2524    NOT_YET_IMPLEMENTED
2525 }
2526
2527 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2528    NOT_YET_IMPLEMENTED
2529 }
2530
2531 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2532    NOT_YET_IMPLEMENTED
2533 }
2534
2535 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2536    NOT_YET_IMPLEMENTED
2537 }
2538
2539 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2540    NOT_YET_IMPLEMENTED
2541 }
2542
2543 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2544    NOT_YET_IMPLEMENTED
2545 }
2546
2547 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2548    NOT_YET_IMPLEMENTED
2549 }
2550
2551 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2552    NOT_YET_IMPLEMENTED
2553 }
2554
2555 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2556    NOT_YET_IMPLEMENTED
2557 }
2558
2559 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2560    NOT_YET_IMPLEMENTED
2561 }
2562
2563 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2564    NOT_YET_IMPLEMENTED
2565 }
2566
2567 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2568    NOT_YET_IMPLEMENTED
2569 }
2570
2571 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2572    NOT_YET_IMPLEMENTED
2573 }
2574
2575 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2576    NOT_YET_IMPLEMENTED
2577 }
2578
2579 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2580    NOT_YET_IMPLEMENTED
2581 }
2582
2583 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2584    NOT_YET_IMPLEMENTED
2585 }
2586
2587 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2588    NOT_YET_IMPLEMENTED
2589 }
2590
2591 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2592    NOT_YET_IMPLEMENTED
2593 }
2594
2595 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2596    NOT_YET_IMPLEMENTED
2597 }
2598
2599 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2600    NOT_YET_IMPLEMENTED
2601 }
2602
2603 int PMPI_Comm_get_errhandler(MPI_Comm comm, MPI_Errhandler* errhandler) {
2604    NOT_YET_IMPLEMENTED
2605 }
2606
2607 int PMPI_Cancel(MPI_Request* request) {
2608    NOT_YET_IMPLEMENTED
2609 }
2610
2611 int PMPI_Buffer_attach(void* buffer, int size) {
2612    NOT_YET_IMPLEMENTED
2613 }
2614
2615 int PMPI_Buffer_detach(void* buffer, int* size) {
2616    NOT_YET_IMPLEMENTED
2617 }
2618
2619 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2620    NOT_YET_IMPLEMENTED
2621 }
2622
2623 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2624 {
2625    NOT_YET_IMPLEMENTED
2626 }
2627
2628 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2629 {
2630    NOT_YET_IMPLEMENTED
2631 }
2632
2633 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2634 {
2635    NOT_YET_IMPLEMENTED
2636 }
2637
2638 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2639 {
2640    NOT_YET_IMPLEMENTED
2641 }
2642
2643 int PMPI_Comm_free_keyval(int* keyval) {
2644    NOT_YET_IMPLEMENTED
2645 }
2646
2647 int PMPI_Pcontrol(const int level )
2648 {
2649    NOT_YET_IMPLEMENTED
2650 }
2651
2652 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2653    NOT_YET_IMPLEMENTED
2654 }
2655
2656 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2657 {
2658   NOT_YET_IMPLEMENTED
2659 }
2660
2661 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2662 {
2663   NOT_YET_IMPLEMENTED
2664 }
2665
2666 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2667 {
2668   NOT_YET_IMPLEMENTED
2669 }
2670
2671 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2672 {
2673   NOT_YET_IMPLEMENTED
2674 }
2675
2676 int PMPI_Type_free_keyval(int* keyval) {
2677   NOT_YET_IMPLEMENTED
2678 }
2679
2680 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2681    NOT_YET_IMPLEMENTED
2682 }
2683
2684 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2685    NOT_YET_IMPLEMENTED
2686 }
2687
2688 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2689    NOT_YET_IMPLEMENTED
2690 }
2691
2692 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2693    NOT_YET_IMPLEMENTED
2694 }
2695
2696 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2697    NOT_YET_IMPLEMENTED
2698 }
2699
2700 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2701    NOT_YET_IMPLEMENTED
2702 }
2703
2704 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2705    NOT_YET_IMPLEMENTED
2706 }
2707
2708 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2709    NOT_YET_IMPLEMENTED
2710 }
2711
2712 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2713    NOT_YET_IMPLEMENTED
2714 }
2715
2716 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2717    NOT_YET_IMPLEMENTED
2718 }
2719
2720 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2721    NOT_YET_IMPLEMENTED
2722 }
2723
2724 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2725    NOT_YET_IMPLEMENTED
2726 }
2727
2728 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2729    NOT_YET_IMPLEMENTED
2730 }
2731
2732 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2733    NOT_YET_IMPLEMENTED
2734 }
2735
2736 int PMPI_Keyval_free(int* keyval) {
2737    NOT_YET_IMPLEMENTED
2738 }
2739
2740 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2741    NOT_YET_IMPLEMENTED
2742 }
2743
2744 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2745    NOT_YET_IMPLEMENTED
2746 }
2747
2748 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2749   NOT_YET_IMPLEMENTED
2750 }
2751
2752 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2753   NOT_YET_IMPLEMENTED
2754 }
2755
2756 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2757   NOT_YET_IMPLEMENTED
2758 }
2759
2760 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2761    NOT_YET_IMPLEMENTED
2762 }
2763
2764 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2765    NOT_YET_IMPLEMENTED
2766 }
2767
2768 int PMPI_Win_fence( int assert,  MPI_Win win){
2769    NOT_YET_IMPLEMENTED
2770 }
2771
2772 int PMPI_Win_free( MPI_Win* win){
2773    NOT_YET_IMPLEMENTED
2774 }
2775
2776 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2777   NOT_YET_IMPLEMENTED
2778 }
2779
2780 int PMPI_Info_create( MPI_Info *info){
2781   NOT_YET_IMPLEMENTED
2782 }
2783
2784 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2785   NOT_YET_IMPLEMENTED
2786 }
2787
2788 int PMPI_Info_free( MPI_Info *info){
2789   NOT_YET_IMPLEMENTED
2790 }
2791
2792 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2793     MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2794   NOT_YET_IMPLEMENTED
2795 }
2796
2797 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2798                           int *num_addresses, int *num_datatypes, int *combiner){
2799   NOT_YET_IMPLEMENTED
2800 }
2801
2802 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2803                           int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2804                           MPI_Datatype* array_of_datatypes){
2805   NOT_YET_IMPLEMENTED
2806 }
2807
2808 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2809                             int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2810                             int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2811   NOT_YET_IMPLEMENTED
2812 }
2813
2814 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2815   NOT_YET_IMPLEMENTED
2816 }
2817
2818 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){
2819   NOT_YET_IMPLEMENTED
2820 }
2821
2822 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
2823   NOT_YET_IMPLEMENTED
2824 }
2825
2826 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
2827                    void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
2828                    MPI_Comm comm){
2829   NOT_YET_IMPLEMENTED
2830 }
2831
2832 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
2833   NOT_YET_IMPLEMENTED
2834 }
2835
2836 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
2837   NOT_YET_IMPLEMENTED
2838 }
2839
2840 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
2841   NOT_YET_IMPLEMENTED
2842 }
2843
2844 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
2845   NOT_YET_IMPLEMENTED
2846 }
2847
2848 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
2849   NOT_YET_IMPLEMENTED
2850 }
2851
2852 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
2853   NOT_YET_IMPLEMENTED
2854 }
2855
2856 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
2857   NOT_YET_IMPLEMENTED
2858 }
2859
2860 int PMPI_Add_error_class( int *errorclass){
2861   NOT_YET_IMPLEMENTED
2862 }
2863
2864 int PMPI_Add_error_code(  int errorclass, int *errorcode){
2865   NOT_YET_IMPLEMENTED
2866 }
2867
2868 int PMPI_Add_error_string( int errorcode, char *string){
2869   NOT_YET_IMPLEMENTED
2870 }
2871
2872 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
2873   NOT_YET_IMPLEMENTED
2874 }
2875
2876 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
2877   NOT_YET_IMPLEMENTED
2878 }
2879
2880 int PMPI_Info_delete(MPI_Info info, char *key){
2881   NOT_YET_IMPLEMENTED
2882 }
2883
2884 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
2885   NOT_YET_IMPLEMENTED
2886 }
2887
2888 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
2889   NOT_YET_IMPLEMENTED
2890 }
2891
2892 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
2893   NOT_YET_IMPLEMENTED
2894 }
2895
2896 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2897   NOT_YET_IMPLEMENTED
2898 }
2899
2900 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){
2901   NOT_YET_IMPLEMENTED
2902 }
2903
2904 int PMPI_Grequest_complete( MPI_Request request){
2905   NOT_YET_IMPLEMENTED
2906 }
2907
2908 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
2909   NOT_YET_IMPLEMENTED
2910 }
2911
2912 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
2913   NOT_YET_IMPLEMENTED
2914 }
2915
2916 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2917   NOT_YET_IMPLEMENTED
2918 }
2919
2920 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
2921   NOT_YET_IMPLEMENTED
2922 }
2923
2924 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
2925   NOT_YET_IMPLEMENTED
2926 }
2927
2928 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
2929   NOT_YET_IMPLEMENTED
2930 }
2931
2932 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
2933   NOT_YET_IMPLEMENTED
2934 }
2935
2936 int PMPI_Open_port( MPI_Info info, char *port_name){
2937   NOT_YET_IMPLEMENTED
2938 }
2939
2940 int PMPI_Close_port(char *port_name){
2941   NOT_YET_IMPLEMENTED
2942 }
2943
2944 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2945   NOT_YET_IMPLEMENTED
2946 }
2947
2948 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){
2949   NOT_YET_IMPLEMENTED
2950 }
2951
2952 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
2953                              int* array_of_maxprocs, MPI_Info* array_of_info, int root,
2954                              MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
2955   NOT_YET_IMPLEMENTED
2956 }
2957
2958 int PMPI_Comm_get_parent( MPI_Comm *parent){
2959   NOT_YET_IMPLEMENTED
2960 }