Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
model-checker : cleanups for heap comparison from global and local variables
[simgrid.git] / src / smpi / smpi_pmpi.c
1 /* Copyright (c) 2007, 2008, 2009, 2010. The SimGrid Team.
2  * All rights reserved.                                                     */
3
4 /* This program is free software; you can redistribute it and/or modify it
5   * under the terms of the license (GNU LGPL) which comes with this package. */
6
7 #include "private.h"
8 #include "smpi_mpi_dt_private.h"
9
10 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_pmpi, smpi,
11                                 "Logging specific to SMPI (pmpi)");
12
13 #ifdef HAVE_TRACING
14 //this function need to be here because of the calls to smpi_bench
15 void TRACE_smpi_set_category(const char *category)
16 {
17   //need to end bench otherwise categories for execution tasks are wrong
18   smpi_bench_end();
19   TRACE_internal_smpi_set_category (category);
20   //begin bench after changing process's category
21   smpi_bench_begin();
22 }
23 #endif
24
25 /* PMPI User level calls */
26
27 int PMPI_Init(int *argc, char ***argv)
28 {
29   smpi_process_init(argc, argv);
30 #ifdef HAVE_TRACING
31   int rank = smpi_process_index();
32   TRACE_smpi_init(rank);
33
34   TRACE_smpi_computing_init(rank);
35 #endif
36   smpi_bench_begin();
37   return MPI_SUCCESS;
38 }
39
40 int PMPI_Finalize(void)
41 {
42   smpi_process_finalize();
43   smpi_bench_end();
44 #ifdef HAVE_TRACING
45   int rank = smpi_process_index();
46   TRACE_smpi_computing_out(rank);
47   TRACE_smpi_finalize(smpi_process_index());
48 #endif
49   smpi_process_destroy();
50   return MPI_SUCCESS;
51 }
52
53 int PMPI_Init_thread(int *argc, char ***argv, int required, int *provided)
54 {
55   if (provided != NULL) {
56     *provided = MPI_THREAD_MULTIPLE;
57   }
58   return MPI_Init(argc, argv);
59 }
60
61 int PMPI_Query_thread(int *provided)
62 {
63   int retval;
64
65   smpi_bench_end();
66   if (provided == NULL) {
67     retval = MPI_ERR_ARG;
68   } else {
69     *provided = MPI_THREAD_MULTIPLE;
70     retval = MPI_SUCCESS;
71   }
72   smpi_bench_begin();
73   return retval;
74 }
75
76 int PMPI_Is_thread_main(int *flag)
77 {
78   int retval;
79
80   smpi_bench_end();
81   if (flag == NULL) {
82     retval = MPI_ERR_ARG;
83   } else {
84     *flag = smpi_process_index() == 0;
85     retval = MPI_SUCCESS;
86   }
87   smpi_bench_begin();
88   return retval;
89 }
90
91 int PMPI_Abort(MPI_Comm comm, int errorcode)
92 {
93   smpi_bench_end();
94   smpi_process_destroy();
95 #ifdef HAVE_TRACING
96   int rank = smpi_process_index();
97   TRACE_smpi_computing_out(rank);
98 #endif
99   // FIXME: should kill all processes in comm instead
100   simcall_process_kill(SIMIX_process_self());
101   return MPI_SUCCESS;
102 }
103
104 double PMPI_Wtime(void)
105 {
106   double time;
107
108   smpi_bench_end();
109   time = SIMIX_get_clock();
110   smpi_bench_begin();
111   return time;
112 }
113 extern double sg_maxmin_precision;
114 double PMPI_Wtick(void)
115 {
116   return sg_maxmin_precision;
117 }
118
119 int PMPI_Address(void *location, MPI_Aint * address)
120 {
121   int retval;
122
123   smpi_bench_end();
124   if (!address) {
125     retval = MPI_ERR_ARG;
126   } else {
127     *address = (MPI_Aint) location;
128     retval = MPI_SUCCESS;
129   }
130   smpi_bench_begin();
131   return retval;
132 }
133
134 int PMPI_Type_free(MPI_Datatype * datatype)
135 {
136   int retval;
137
138   smpi_bench_end();
139   if (!datatype) {
140     retval = MPI_ERR_ARG;
141   } else {
142     smpi_datatype_free(datatype);
143     retval = MPI_SUCCESS;
144   }
145   smpi_bench_begin();
146   return retval;
147 }
148
149 int PMPI_Type_size(MPI_Datatype datatype, int *size)
150 {
151   int retval;
152
153   smpi_bench_end();
154   if (datatype == MPI_DATATYPE_NULL) {
155     retval = MPI_ERR_TYPE;
156   } else if (size == NULL) {
157     retval = MPI_ERR_ARG;
158   } else {
159     *size = (int) smpi_datatype_size(datatype);
160     retval = MPI_SUCCESS;
161   }
162   smpi_bench_begin();
163   return retval;
164 }
165
166 int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
167 {
168   int retval;
169
170   smpi_bench_end();
171   if (datatype == MPI_DATATYPE_NULL) {
172     retval = MPI_ERR_TYPE;
173   } else if (lb == NULL || extent == NULL) {
174     retval = MPI_ERR_ARG;
175   } else {
176     retval = smpi_datatype_extent(datatype, lb, extent);
177   }
178   smpi_bench_begin();
179   return retval;
180 }
181
182 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
183 {
184   int retval;
185
186   smpi_bench_end();
187   if (datatype == MPI_DATATYPE_NULL) {
188     retval = MPI_ERR_TYPE;
189   } else if (extent == NULL) {
190     retval = MPI_ERR_ARG;
191   } else {
192     *extent = smpi_datatype_get_extent(datatype);
193     retval = MPI_SUCCESS;
194   }
195   smpi_bench_begin();
196   return retval;
197 }
198
199 int PMPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
200 {
201   int retval;
202
203   smpi_bench_end();
204   if (datatype == MPI_DATATYPE_NULL) {
205     retval = MPI_ERR_TYPE;
206   } else if (disp == NULL) {
207     retval = MPI_ERR_ARG;
208   } else {
209     *disp = smpi_datatype_lb(datatype);
210     retval = MPI_SUCCESS;
211   }
212   smpi_bench_begin();
213   return retval;
214 }
215
216 int PMPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
217 {
218   int retval;
219
220   smpi_bench_end();
221   if (datatype == MPI_DATATYPE_NULL) {
222     retval = MPI_ERR_TYPE;
223   } else if (disp == NULL) {
224     retval = MPI_ERR_ARG;
225   } else {
226     *disp = smpi_datatype_ub(datatype);
227     retval = MPI_SUCCESS;
228   }
229   smpi_bench_begin();
230   return retval;
231 }
232
233 int PMPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
234 {
235   int retval;
236
237   smpi_bench_end();
238   if (function == NULL || op == NULL) {
239     retval = MPI_ERR_ARG;
240   } else {
241     *op = smpi_op_new(function, commute);
242     retval = MPI_SUCCESS;
243   }
244   smpi_bench_begin();
245   return retval;
246 }
247
248 int PMPI_Op_free(MPI_Op * op)
249 {
250   int retval;
251
252   smpi_bench_end();
253   if (op == NULL) {
254     retval = MPI_ERR_ARG;
255   } else if (*op == MPI_OP_NULL) {
256     retval = MPI_ERR_OP;
257   } else {
258     smpi_op_destroy(*op);
259     *op = MPI_OP_NULL;
260     retval = MPI_SUCCESS;
261   }
262   smpi_bench_begin();
263   return retval;
264 }
265
266 int PMPI_Group_free(MPI_Group * group)
267 {
268   int retval;
269
270   smpi_bench_end();
271   if (group == NULL) {
272     retval = MPI_ERR_ARG;
273   } else {
274     if(*group!= smpi_comm_group(MPI_COMM_WORLD))// do not free the group of the comm_world
275     smpi_group_destroy(*group);
276     *group = MPI_GROUP_NULL;
277     retval = MPI_SUCCESS;
278   }
279   smpi_bench_begin();
280   return retval;
281 }
282
283 int PMPI_Group_size(MPI_Group group, int *size)
284 {
285   int retval;
286
287   smpi_bench_end();
288   if (group == MPI_GROUP_NULL) {
289     retval = MPI_ERR_GROUP;
290   } else if (size == NULL) {
291     retval = MPI_ERR_ARG;
292   } else {
293     *size = smpi_group_size(group);
294     retval = MPI_SUCCESS;
295   }
296   smpi_bench_begin();
297   return retval;
298 }
299
300 int PMPI_Group_rank(MPI_Group group, int *rank)
301 {
302   int retval;
303
304   smpi_bench_end();
305   if (group == MPI_GROUP_NULL) {
306     retval = MPI_ERR_GROUP;
307   } else if (rank == NULL) {
308     retval = MPI_ERR_ARG;
309   } else {
310     *rank = smpi_group_rank(group, smpi_process_index());
311     retval = MPI_SUCCESS;
312   }
313   smpi_bench_begin();
314   return retval;
315 }
316
317 int PMPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1,
318                               MPI_Group group2, int *ranks2)
319 {
320   int retval, i, index;
321
322   smpi_bench_end();
323   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
324     retval = MPI_ERR_GROUP;
325   } else {
326     for (i = 0; i < n; i++) {
327       index = smpi_group_index(group1, ranks1[i]);
328       ranks2[i] = smpi_group_rank(group2, index);
329     }
330     retval = MPI_SUCCESS;
331   }
332   smpi_bench_begin();
333   return retval;
334 }
335
336 int PMPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
337 {
338   int retval;
339
340   smpi_bench_end();
341   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
342     retval = MPI_ERR_GROUP;
343   } else if (result == NULL) {
344     retval = MPI_ERR_ARG;
345   } else {
346     *result = smpi_group_compare(group1, group2);
347     retval = MPI_SUCCESS;
348   }
349   smpi_bench_begin();
350   return retval;
351 }
352
353 int PMPI_Group_union(MPI_Group group1, MPI_Group group2,
354                     MPI_Group * newgroup)
355 {
356   int retval, i, proc1, proc2, size, size2;
357
358   smpi_bench_end();
359   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
360     retval = MPI_ERR_GROUP;
361   } else if (newgroup == NULL) {
362     retval = MPI_ERR_ARG;
363   } else {
364     size = smpi_group_size(group1);
365     size2 = smpi_group_size(group2);
366     for (i = 0; i < size2; i++) {
367       proc2 = smpi_group_index(group2, i);
368       proc1 = smpi_group_rank(group1, proc2);
369       if (proc1 == MPI_UNDEFINED) {
370         size++;
371       }
372     }
373     if (size == 0) {
374       *newgroup = MPI_GROUP_EMPTY;
375     } else {
376       *newgroup = smpi_group_new(size);
377       size2 = smpi_group_size(group1);
378       for (i = 0; i < size2; i++) {
379         proc1 = smpi_group_index(group1, i);
380         smpi_group_set_mapping(*newgroup, proc1, i);
381       }
382       for (i = size2; i < size; i++) {
383         proc2 = smpi_group_index(group2, i - size2);
384         smpi_group_set_mapping(*newgroup, proc2, i);
385       }
386     }
387     smpi_group_use(*newgroup);
388     retval = MPI_SUCCESS;
389   }
390   smpi_bench_begin();
391   return retval;
392 }
393
394 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
395                            MPI_Group * newgroup)
396 {
397   int retval, i, proc1, proc2, size;
398
399   smpi_bench_end();
400   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
401     retval = MPI_ERR_GROUP;
402   } else if (newgroup == NULL) {
403     retval = MPI_ERR_ARG;
404   } else {
405     size = smpi_group_size(group2);
406     for (i = 0; i < size; i++) {
407       proc2 = smpi_group_index(group2, i);
408       proc1 = smpi_group_rank(group1, proc2);
409       if (proc1 == MPI_UNDEFINED) {
410         size--;
411       }
412     }
413     if (size == 0) {
414       *newgroup = MPI_GROUP_EMPTY;
415     } else {
416       *newgroup = smpi_group_new(size);
417       int j=0;
418       for (i = 0; i < smpi_group_size(group2); i++) {
419         proc2 = smpi_group_index(group2, i);
420         proc1 = smpi_group_rank(group1, proc2);
421         if (proc1 != MPI_UNDEFINED) {
422           smpi_group_set_mapping(*newgroup, proc2, j);
423           j++;
424         }
425       }
426     }
427     smpi_group_use(*newgroup);
428     retval = MPI_SUCCESS;
429   }
430   smpi_bench_begin();
431   return retval;
432 }
433
434 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
435 {
436   int retval, i, proc1, proc2, size, size2;
437
438   smpi_bench_end();
439   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
440     retval = MPI_ERR_GROUP;
441   } else if (newgroup == NULL) {
442     retval = MPI_ERR_ARG;
443   } else {
444     size = size2 = smpi_group_size(group1);
445     for (i = 0; i < size2; i++) {
446       proc1 = smpi_group_index(group1, i);
447       proc2 = smpi_group_rank(group2, proc1);
448       if (proc2 != MPI_UNDEFINED) {
449         size--;
450       }
451     }
452     if (size == 0) {
453       *newgroup = MPI_GROUP_EMPTY;
454     } else {
455       *newgroup = smpi_group_new(size);
456       for (i = 0; i < size2; i++) {
457         proc1 = smpi_group_index(group1, i);
458         proc2 = smpi_group_rank(group2, proc1);
459         if (proc2 == MPI_UNDEFINED) {
460           smpi_group_set_mapping(*newgroup, proc1, i);
461         }
462       }
463     }
464     smpi_group_use(*newgroup);
465     retval = MPI_SUCCESS;
466   }
467   smpi_bench_begin();
468   return retval;
469 }
470
471 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
472 {
473   int retval, i, index;
474
475   smpi_bench_end();
476   if (group == MPI_GROUP_NULL) {
477     retval = MPI_ERR_GROUP;
478   } else if (newgroup == NULL) {
479     retval = MPI_ERR_ARG;
480   } else {
481     if (n == 0) {
482       *newgroup = MPI_GROUP_EMPTY;
483     } else if (n == smpi_group_size(group)) {
484       *newgroup = group;
485     } else {
486       *newgroup = smpi_group_new(n);
487       for (i = 0; i < n; i++) {
488         index = smpi_group_index(group, ranks[i]);
489         smpi_group_set_mapping(*newgroup, index, i);
490       }
491     }
492     smpi_group_use(*newgroup);
493     retval = MPI_SUCCESS;
494   }
495   smpi_bench_begin();
496   return retval;
497 }
498
499 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
500 {
501   int retval, i, j, newsize, oldsize, index;
502
503   smpi_bench_end();
504   if (group == MPI_GROUP_NULL) {
505     retval = MPI_ERR_GROUP;
506   } else if (newgroup == NULL) {
507     retval = MPI_ERR_ARG;
508   } else {
509     if (n == 0) {
510       *newgroup = group;
511     } else if (n == smpi_group_size(group)) {
512       *newgroup = MPI_GROUP_EMPTY;
513     } else {
514       oldsize=smpi_group_size(group);
515       newsize = oldsize - n;
516       *newgroup = smpi_group_new(newsize);
517
518       int* to_exclude=xbt_new(int, smpi_group_size(group));
519       for(i=0; i<oldsize; i++)
520         to_exclude[i]=0;
521       for(i=0; i<n; i++)
522         to_exclude[ranks[i]]=1;
523
524       j=0;
525       for(i=0; i<oldsize; i++){
526         if(to_exclude[i]==0){
527           index = smpi_group_index(group, i);
528           smpi_group_set_mapping(*newgroup, index, j);
529           j++;
530         }
531       }
532
533       xbt_free(to_exclude);
534     }
535     smpi_group_use(*newgroup);
536     retval = MPI_SUCCESS;
537   }
538   smpi_bench_begin();
539   return retval;
540 }
541
542 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
543                          MPI_Group * newgroup)
544 {
545   int retval, i, j, rank, size, index;
546
547   smpi_bench_end();
548   if (group == MPI_GROUP_NULL) {
549     retval = MPI_ERR_GROUP;
550   } else if (newgroup == NULL) {
551     retval = MPI_ERR_ARG;
552   } else {
553     if (n == 0) {
554       *newgroup = MPI_GROUP_EMPTY;
555     } else {
556       size = 0;
557       for (i = 0; i < n; i++) {
558         for (rank = ranges[i][0];       /* First */
559              rank >= 0 && rank <= ranges[i][1]; /* Last */
560              rank += ranges[i][2] /* Stride */ ) {
561           size++;
562         }
563       }
564       if (size == smpi_group_size(group)) {
565         *newgroup = group;
566       } else {
567         *newgroup = smpi_group_new(size);
568         j = 0;
569         for (i = 0; i < n; i++) {
570           for (rank = ranges[i][0];     /* First */
571                rank >= 0 && rank <= ranges[i][1];       /* Last */
572                rank += ranges[i][2] /* Stride */ ) {
573             index = smpi_group_index(group, rank);
574             smpi_group_set_mapping(*newgroup, index, j);
575             j++;
576           }
577         }
578       }
579     }
580     smpi_group_use(*newgroup);
581     retval = MPI_SUCCESS;
582   }
583   smpi_bench_begin();
584   return retval;
585 }
586
587 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
588                          MPI_Group * newgroup)
589 {
590   int retval, i, newrank, rank, size, index, add;
591
592   smpi_bench_end();
593   if (group == MPI_GROUP_NULL) {
594     retval = MPI_ERR_GROUP;
595   } else if (newgroup == NULL) {
596     retval = MPI_ERR_ARG;
597   } else {
598     if (n == 0) {
599       *newgroup = group;
600     } else {
601       size = smpi_group_size(group);
602       for (i = 0; i < n; i++) {
603         for (rank = ranges[i][0];       /* First */
604              rank >= 0 && rank <= ranges[i][1]; /* Last */
605              rank += ranges[i][2] /* Stride */ ) {
606           size--;
607         }
608       }
609       if (size == 0) {
610         *newgroup = MPI_GROUP_EMPTY;
611       } else {
612         *newgroup = smpi_group_new(size);
613         newrank = 0;
614         while (newrank < size) {
615           for (i = 0; i < n; i++) {
616             add = 1;
617             for (rank = ranges[i][0];   /* First */
618                  rank >= 0 && rank <= ranges[i][1];     /* Last */
619                  rank += ranges[i][2] /* Stride */ ) {
620               if (rank == newrank) {
621                 add = 0;
622                 break;
623               }
624             }
625             if (add == 1) {
626               index = smpi_group_index(group, newrank);
627               smpi_group_set_mapping(*newgroup, index, newrank);
628             }
629           }
630           newrank++; //added to avoid looping, need to be checked ..
631         }
632       }
633     }
634     smpi_group_use(*newgroup);
635     retval = MPI_SUCCESS;
636   }
637   smpi_bench_begin();
638   return retval;
639 }
640
641 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
642 {
643   int retval;
644
645   smpi_bench_end();
646   if (comm == MPI_COMM_NULL) {
647     retval = MPI_ERR_COMM;
648   } else if (rank == NULL) {
649     retval = MPI_ERR_ARG;
650   } else {
651     *rank = smpi_comm_rank(comm);
652     retval = MPI_SUCCESS;
653   }
654   smpi_bench_begin();
655   return retval;
656 }
657
658 int PMPI_Comm_size(MPI_Comm comm, int *size)
659 {
660   int retval;
661
662   smpi_bench_end();
663   if (comm == MPI_COMM_NULL) {
664     retval = MPI_ERR_COMM;
665   } else if (size == NULL) {
666     retval = MPI_ERR_ARG;
667   } else {
668     *size = smpi_comm_size(comm);
669     retval = MPI_SUCCESS;
670   }
671   smpi_bench_begin();
672   return retval;
673 }
674
675 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
676 {
677   int retval;
678
679   smpi_bench_end();
680   if (comm == MPI_COMM_NULL)  {
681     retval = MPI_ERR_COMM;
682   } else if (name == NULL || len == NULL)  {
683     retval = MPI_ERR_ARG;
684   } else {
685     smpi_comm_get_name(comm, name, len);
686     retval = MPI_SUCCESS;
687   }
688   smpi_bench_begin();
689   return retval;
690 }
691
692 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
693 {
694   int retval;
695
696   smpi_bench_end();
697   if (comm == MPI_COMM_NULL) {
698     retval = MPI_ERR_COMM;
699   } else if (group == NULL) {
700     retval = MPI_ERR_ARG;
701   } else {
702     *group = smpi_comm_group(comm);
703     retval = MPI_SUCCESS;
704   }
705   smpi_bench_begin();
706   return retval;
707 }
708
709 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
710 {
711   int retval;
712
713   smpi_bench_end();
714   if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
715     retval = MPI_ERR_COMM;
716   } else if (result == NULL) {
717     retval = MPI_ERR_ARG;
718   } else {
719     if (comm1 == comm2) {       /* Same communicators means same groups */
720       *result = MPI_IDENT;
721     } else {
722       *result =
723           smpi_group_compare(smpi_comm_group(comm1),
724                              smpi_comm_group(comm2));
725       if (*result == MPI_IDENT) {
726         *result = MPI_CONGRUENT;
727       }
728     }
729     retval = MPI_SUCCESS;
730   }
731   smpi_bench_begin();
732   return retval;
733 }
734
735 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
736 {
737   int retval;
738
739   smpi_bench_end();
740   if (comm == MPI_COMM_NULL) {
741     retval = MPI_ERR_COMM;
742   } else if (newcomm == NULL) {
743     retval = MPI_ERR_ARG;
744   } else {
745     *newcomm = smpi_comm_new(smpi_comm_group(comm));
746     retval = MPI_SUCCESS;
747   }
748   smpi_bench_begin();
749   return retval;
750 }
751
752 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
753 {
754   int retval;
755
756   smpi_bench_end();
757   if (comm == MPI_COMM_NULL) {
758     retval = MPI_ERR_COMM;
759   } else if (group == MPI_GROUP_NULL) {
760     retval = MPI_ERR_GROUP;
761   } else if (newcomm == NULL) {
762     retval = MPI_ERR_ARG;
763   } else {
764     *newcomm = smpi_comm_new(group);
765     retval = MPI_SUCCESS;
766   }
767   smpi_bench_begin();
768   return retval;
769 }
770
771 int PMPI_Comm_free(MPI_Comm * comm)
772 {
773   int retval;
774
775   smpi_bench_end();
776   if (comm == NULL) {
777     retval = MPI_ERR_ARG;
778   } else if (*comm == MPI_COMM_NULL) {
779     retval = MPI_ERR_COMM;
780   } else {
781     smpi_comm_destroy(*comm);
782     *comm = MPI_COMM_NULL;
783     retval = MPI_SUCCESS;
784   }
785   smpi_bench_begin();
786   return retval;
787 }
788
789 int PMPI_Comm_disconnect(MPI_Comm * comm)
790 {
791   /* TODO: wait until all communication in comm are done */
792   int retval;
793
794   smpi_bench_end();
795   if (comm == NULL) {
796     retval = MPI_ERR_ARG;
797   } else if (*comm == MPI_COMM_NULL) {
798     retval = MPI_ERR_COMM;
799   } else {
800     smpi_comm_destroy(*comm);
801     *comm = MPI_COMM_NULL;
802     retval = MPI_SUCCESS;
803   }
804   smpi_bench_begin();
805   return retval;
806 }
807
808 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
809 {
810   int retval;
811
812   smpi_bench_end();
813   if (comm_out == NULL) {
814     retval = MPI_ERR_ARG;
815   } else if (comm == MPI_COMM_NULL) {
816     retval = MPI_ERR_COMM;
817   } else {
818     *comm_out = smpi_comm_split(comm, color, key);
819     retval = MPI_SUCCESS;
820   }
821   smpi_bench_begin();
822   return retval;
823 }
824
825 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
826                   int tag, MPI_Comm comm, MPI_Request * request)
827 {
828   int retval;
829
830   smpi_bench_end();
831   if (request == NULL) {
832     retval = MPI_ERR_ARG;
833   } else if (comm == MPI_COMM_NULL) {
834     retval = MPI_ERR_COMM;
835   } else if (dst == MPI_PROC_NULL) {
836     retval = MPI_SUCCESS;
837   } else {
838     *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
839     retval = MPI_SUCCESS;
840   }
841   smpi_bench_begin();
842   return retval;
843 }
844
845 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
846                   int tag, MPI_Comm comm, MPI_Request * request)
847 {
848   int retval;
849
850   smpi_bench_end();
851   if (request == NULL) {
852     retval = MPI_ERR_ARG;
853   } else if (comm == MPI_COMM_NULL) {
854     retval = MPI_ERR_COMM;
855   } else if (src == MPI_PROC_NULL) {
856       retval = MPI_SUCCESS;
857   } else {
858     *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
859     retval = MPI_SUCCESS;
860   }
861   smpi_bench_begin();
862   return retval;
863 }
864
865 int PMPI_Start(MPI_Request * request)
866 {
867   int retval;
868
869   smpi_bench_end();
870   if (request == NULL || *request == MPI_REQUEST_NULL) {
871     retval = MPI_ERR_ARG;
872   } else {
873     smpi_mpi_start(*request);
874     retval = MPI_SUCCESS;
875   }
876   smpi_bench_begin();
877   return retval;
878 }
879
880 int PMPI_Startall(int count, MPI_Request * requests)
881 {
882   int retval;
883
884   smpi_bench_end();
885   if (requests == NULL) {
886     retval = MPI_ERR_ARG;
887   } else {
888     smpi_mpi_startall(count, requests);
889     retval = MPI_SUCCESS;
890   }
891   smpi_bench_begin();
892   return retval;
893 }
894
895 int PMPI_Request_free(MPI_Request * request)
896 {
897   int retval;
898
899   smpi_bench_end();
900   if (*request == MPI_REQUEST_NULL) {
901     retval = MPI_ERR_ARG;
902   } else {
903     if((*request)->flags & PERSISTENT)(*request)->refcount--;
904     smpi_mpi_request_free(request);
905     retval = MPI_SUCCESS;
906   }
907   smpi_bench_begin();
908   return retval;
909 }
910
911 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
912               int tag, MPI_Comm comm, MPI_Request * request)
913 {
914   int retval;
915
916   smpi_bench_end();
917
918   if (request == NULL) {
919     retval = MPI_ERR_ARG;
920   } else if (comm == MPI_COMM_NULL) {
921     retval = MPI_ERR_COMM;
922   } else if (src == MPI_PROC_NULL) {
923     *request = MPI_REQUEST_NULL;
924     retval = MPI_SUCCESS;
925   } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
926     retval = MPI_ERR_COMM;
927   } else if (count < 0) {
928     retval = MPI_ERR_COUNT;
929   } else if (buf==NULL && count > 0) {
930     retval = MPI_ERR_COUNT;
931   } else if (datatype == MPI_DATATYPE_NULL){
932     retval = MPI_ERR_TYPE;
933   } else if(tag<0 && tag !=  MPI_ANY_TAG){
934     retval = MPI_ERR_TAG;
935   } else {
936
937 #ifdef HAVE_TRACING
938   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
939   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
940   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
941 #endif
942
943     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
944     retval = MPI_SUCCESS;
945
946 #ifdef HAVE_TRACING
947   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
948   (*request)->recv = 1;
949 #endif
950   }
951
952   smpi_bench_begin();
953   return retval;
954 }
955
956
957 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
958               int tag, MPI_Comm comm, MPI_Request * request)
959 {
960   int retval;
961
962   smpi_bench_end();
963   if (request == NULL) {
964     retval = MPI_ERR_ARG;
965   } else if (comm == MPI_COMM_NULL) {
966     retval = MPI_ERR_COMM;
967   } else if (dst == MPI_PROC_NULL) {
968     *request = MPI_REQUEST_NULL;
969     retval = MPI_SUCCESS;
970   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
971     retval = MPI_ERR_COMM;
972   } else if (count < 0) {
973     retval = MPI_ERR_COUNT;
974   } else if (buf==NULL && count > 0) {
975     retval = MPI_ERR_COUNT;
976   } else if (datatype == MPI_DATATYPE_NULL){
977     retval = MPI_ERR_TYPE;
978   } else if(tag<0 && tag !=  MPI_ANY_TAG){
979     retval = MPI_ERR_TAG;
980   } else {
981
982 #ifdef HAVE_TRACING
983   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
984   TRACE_smpi_computing_out(rank);
985   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
986   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
987   TRACE_smpi_send(rank, rank, dst_traced);
988 #endif
989
990     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
991     retval = MPI_SUCCESS;
992
993 #ifdef HAVE_TRACING
994   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
995   (*request)->send = 1;
996   TRACE_smpi_computing_in(rank);
997 #endif
998   }
999
1000   smpi_bench_begin();
1001   return retval;
1002 }
1003
1004
1005
1006 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1007              MPI_Comm comm, MPI_Status * status)
1008 {
1009   int retval;
1010
1011   smpi_bench_end();
1012   if (comm == MPI_COMM_NULL) {
1013     retval = MPI_ERR_COMM;
1014   } else if (src == MPI_PROC_NULL) {
1015     smpi_empty_status(status);
1016     status->MPI_SOURCE = MPI_PROC_NULL;
1017     retval = MPI_SUCCESS;
1018   } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1019     retval = MPI_ERR_COMM;
1020   } else if (count < 0) {
1021     retval = MPI_ERR_COUNT;
1022   } else if (buf==NULL && count > 0) {
1023     retval = MPI_ERR_COUNT;
1024   } else if (datatype == MPI_DATATYPE_NULL){
1025     retval = MPI_ERR_TYPE;
1026   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1027     retval = MPI_ERR_TAG;
1028   } else {
1029 #ifdef HAVE_TRACING
1030   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1031   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1032   TRACE_smpi_computing_out(rank);
1033
1034   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
1035 #endif
1036
1037     smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1038     retval = MPI_SUCCESS;
1039
1040 #ifdef HAVE_TRACING
1041   //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1042   if(status!=MPI_STATUS_IGNORE)src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1043   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1044   TRACE_smpi_recv(rank, src_traced, rank);
1045   TRACE_smpi_computing_in(rank);
1046 #endif
1047   }
1048
1049   smpi_bench_begin();
1050   return retval;
1051 }
1052
1053 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1054              MPI_Comm comm)
1055 {
1056   int retval;
1057
1058   smpi_bench_end();
1059
1060   if (comm == MPI_COMM_NULL) {
1061     retval = MPI_ERR_COMM;
1062   } else if (dst == MPI_PROC_NULL) {
1063     retval = MPI_SUCCESS;
1064   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1065     retval = MPI_ERR_COMM;
1066   } else if (count < 0) {
1067     retval = MPI_ERR_COUNT;
1068   } else if (buf==NULL && count > 0) {
1069     retval = MPI_ERR_COUNT;
1070   } else if (datatype == MPI_DATATYPE_NULL){
1071     retval = MPI_ERR_TYPE;
1072   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1073     retval = MPI_ERR_TAG;
1074   } else {
1075
1076 #ifdef HAVE_TRACING
1077   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1078   TRACE_smpi_computing_out(rank);
1079   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1080   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1081   TRACE_smpi_send(rank, rank, dst_traced);
1082 #endif
1083
1084     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1085     retval = MPI_SUCCESS;
1086
1087 #ifdef HAVE_TRACING
1088   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1089   TRACE_smpi_computing_in(rank);
1090 #endif
1091   }
1092
1093   smpi_bench_begin();
1094   return retval;
1095 }
1096
1097 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1098                  int dst, int sendtag, void *recvbuf, int recvcount,
1099                  MPI_Datatype recvtype, int src, int recvtag,
1100                  MPI_Comm comm, MPI_Status * status)
1101 {
1102   int retval;
1103
1104   smpi_bench_end();
1105
1106   if (comm == MPI_COMM_NULL) {
1107     retval = MPI_ERR_COMM;
1108   } else if (sendtype == MPI_DATATYPE_NULL
1109              || recvtype == MPI_DATATYPE_NULL) {
1110     retval = MPI_ERR_TYPE;
1111   } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1112       smpi_empty_status(status);
1113       status->MPI_SOURCE = MPI_PROC_NULL;
1114       retval = MPI_SUCCESS;
1115   }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1116       (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1117     retval = MPI_ERR_COMM;
1118   } else if (sendcount < 0 || recvcount<0) {
1119       retval = MPI_ERR_COUNT;
1120   } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1121     retval = MPI_ERR_COUNT;
1122   } else if((sendtag<0 && sendtag !=  MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1123     retval = MPI_ERR_TAG;
1124   } else {
1125
1126 #ifdef HAVE_TRACING
1127   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1128   TRACE_smpi_computing_out(rank);
1129   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1130   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1131   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1132   TRACE_smpi_send(rank, rank, dst_traced);
1133   TRACE_smpi_send(rank, src_traced, rank);
1134 #endif
1135
1136
1137     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1138                       recvcount, recvtype, src, recvtag, comm, status);
1139     retval = MPI_SUCCESS;
1140
1141 #ifdef HAVE_TRACING
1142   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1143   TRACE_smpi_recv(rank, rank, dst_traced);
1144   TRACE_smpi_recv(rank, src_traced, rank);
1145   TRACE_smpi_computing_in(rank);
1146 #endif
1147
1148   }
1149
1150   smpi_bench_begin();
1151   return retval;
1152 }
1153
1154 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1155                          int dst, int sendtag, int src, int recvtag,
1156                          MPI_Comm comm, MPI_Status * status)
1157 {
1158   //TODO: suboptimal implementation
1159   void *recvbuf;
1160   int retval;
1161   if ((datatype == MPI_DATATYPE_NULL)||(datatype->has_subtype==1)) {
1162       retval = MPI_ERR_TYPE;
1163   } else if (count < 0) {
1164       retval = MPI_ERR_COUNT;
1165   } else {
1166     int size = smpi_datatype_size(datatype) * count;
1167     recvbuf = xbt_new(char, size);
1168     retval =
1169         MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1170                      datatype, src, recvtag, comm, status);
1171     if(retval==MPI_SUCCESS){
1172         memcpy(buf, recvbuf, size * sizeof(char));
1173     }
1174     xbt_free(recvbuf);
1175
1176   }
1177   return retval;
1178 }
1179
1180 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1181 {
1182   int retval;
1183
1184   smpi_bench_end();
1185   if (request == MPI_REQUEST_NULL || flag == NULL) {
1186     retval = MPI_ERR_ARG;
1187   } else if (*request == MPI_REQUEST_NULL) {
1188     *flag= TRUE;
1189     retval = MPI_ERR_REQUEST;
1190   } else {
1191     *flag = smpi_mpi_test(request, status);
1192     retval = MPI_SUCCESS;
1193   }
1194   smpi_bench_begin();
1195   return retval;
1196 }
1197
1198 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1199                 MPI_Status * status)
1200 {
1201   int retval;
1202
1203   smpi_bench_end();
1204   if (index == NULL || flag == NULL) {
1205     retval = MPI_ERR_ARG;
1206   } else {
1207     *flag = smpi_mpi_testany(count, requests, index, status);
1208     retval = MPI_SUCCESS;
1209   }
1210   smpi_bench_begin();
1211   return retval;
1212 }
1213
1214 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1215 {
1216   int retval;
1217
1218   smpi_bench_end();
1219   if (flag == NULL) {
1220     retval = MPI_ERR_ARG;
1221   } else {
1222     *flag = smpi_mpi_testall(count, requests, statuses);
1223     retval = MPI_SUCCESS;
1224   }
1225   smpi_bench_begin();
1226   return retval;
1227 }
1228
1229 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1230   int retval;
1231   smpi_bench_end();
1232
1233   if (status == NULL) {
1234     retval = MPI_ERR_ARG;
1235   } else if (comm == MPI_COMM_NULL) {
1236     retval = MPI_ERR_COMM;
1237   } else if (source == MPI_PROC_NULL) {
1238     smpi_empty_status(status);
1239     status->MPI_SOURCE = MPI_PROC_NULL;
1240     retval = MPI_SUCCESS;
1241   } else {
1242     smpi_mpi_probe(source, tag, comm, status);
1243     retval = MPI_SUCCESS;
1244   }
1245   smpi_bench_begin();
1246   return retval;
1247 }
1248
1249
1250 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1251   int retval;
1252   smpi_bench_end();
1253
1254   if (flag == NULL) {
1255     retval = MPI_ERR_ARG;
1256   } else if (status == NULL) {
1257     retval = MPI_ERR_ARG;
1258   } else if (comm == MPI_COMM_NULL) {
1259     retval = MPI_ERR_COMM;
1260   } else if (source == MPI_PROC_NULL) {
1261     smpi_empty_status(status);
1262     status->MPI_SOURCE = MPI_PROC_NULL;
1263     retval = MPI_SUCCESS;
1264   } else {
1265     smpi_mpi_iprobe(source, tag, comm, flag, status);
1266     retval = MPI_SUCCESS;
1267   }
1268   smpi_bench_begin();
1269   return retval;
1270 }
1271
1272 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1273 {
1274   int retval;
1275
1276   smpi_bench_end();
1277
1278   if (request == NULL) {
1279     retval = MPI_ERR_ARG;
1280   } else if (*request == MPI_REQUEST_NULL) {
1281     retval = MPI_ERR_REQUEST;
1282   } else {
1283
1284 #ifdef HAVE_TRACING
1285   int rank = request && (*request)->comm != MPI_COMM_NULL
1286       ? smpi_process_index()
1287       : -1;
1288   TRACE_smpi_computing_out(rank);
1289
1290   MPI_Group group = smpi_comm_group((*request)->comm);
1291   int src_traced = smpi_group_index(group, (*request)->src);
1292   int dst_traced = smpi_group_index(group, (*request)->dst);
1293   int is_wait_for_receive = (*request)->recv;
1294   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1295 #endif
1296
1297     smpi_mpi_wait(request, status);
1298     retval = MPI_SUCCESS;
1299
1300 #ifdef HAVE_TRACING
1301   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1302   if (is_wait_for_receive) {
1303     TRACE_smpi_recv(rank, src_traced, dst_traced);
1304   }
1305   TRACE_smpi_computing_in(rank);
1306 #endif
1307
1308   }
1309
1310   smpi_bench_begin();
1311   return retval;
1312 }
1313
1314 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1315 {
1316   int retval;
1317
1318   smpi_bench_end();
1319 #ifdef HAVE_TRACING
1320   //save requests information for tracing
1321   int i;
1322   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), NULL);
1323   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), NULL);
1324   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), NULL);
1325   for (i = 0; i < count; i++) {
1326     MPI_Request req = requests[i];      //already received requests are no longer valid
1327     if (req) {
1328       int *asrc = xbt_new(int, 1);
1329       int *adst = xbt_new(int, 1);
1330       int *arecv = xbt_new(int, 1);
1331       *asrc = req->src;
1332       *adst = req->dst;
1333       *arecv = req->recv;
1334       xbt_dynar_insert_at(srcs, i, asrc);
1335       xbt_dynar_insert_at(dsts, i, adst);
1336       xbt_dynar_insert_at(recvs, i, arecv);
1337       xbt_free(asrc);
1338       xbt_free(adst);
1339       xbt_free(arecv);
1340     } else {
1341       int *t = xbt_new(int, 1);
1342       xbt_dynar_insert_at(srcs, i, t);
1343       xbt_dynar_insert_at(dsts, i, t);
1344       xbt_dynar_insert_at(recvs, i, t);
1345       xbt_free(t);
1346     }
1347   }
1348   int rank_traced = smpi_process_index();
1349   TRACE_smpi_computing_out(rank_traced);
1350
1351   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1352
1353 #endif
1354   if (index == NULL) {
1355     retval = MPI_ERR_ARG;
1356   } else {
1357     *index = smpi_mpi_waitany(count, requests, status);
1358     retval = MPI_SUCCESS;
1359   }
1360 #ifdef HAVE_TRACING
1361   if(*index!=MPI_UNDEFINED){
1362     int src_traced, dst_traced, is_wait_for_receive;
1363     xbt_dynar_get_cpy(srcs, *index, &src_traced);
1364     xbt_dynar_get_cpy(dsts, *index, &dst_traced);
1365     xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
1366     if (is_wait_for_receive) {
1367       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1368     }
1369     TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1370     //clean-up of dynars
1371     xbt_dynar_free(&srcs);
1372     xbt_dynar_free(&dsts);
1373     xbt_dynar_free(&recvs);
1374   }
1375   TRACE_smpi_computing_in(rank_traced);
1376 #endif
1377   smpi_bench_begin();
1378   return retval;
1379 }
1380
1381 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1382 {
1383
1384   smpi_bench_end();
1385 #ifdef HAVE_TRACING
1386   //save information from requests
1387   int i;
1388   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), NULL);
1389   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), NULL);
1390   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), NULL);
1391   for (i = 0; i < count; i++) {
1392     MPI_Request req = requests[i];
1393     if(req){
1394       int *asrc = xbt_new(int, 1);
1395       int *adst = xbt_new(int, 1);
1396       int *arecv = xbt_new(int, 1);
1397       *asrc = req->src;
1398       *adst = req->dst;
1399       *arecv = req->recv;
1400       xbt_dynar_insert_at(srcs, i, asrc);
1401       xbt_dynar_insert_at(dsts, i, adst);
1402       xbt_dynar_insert_at(recvs, i, arecv);
1403       xbt_free(asrc);
1404       xbt_free(adst);
1405       xbt_free(arecv);
1406     }else {
1407       int *t = xbt_new(int, 1);
1408       xbt_dynar_insert_at(srcs, i, t);
1409       xbt_dynar_insert_at(dsts, i, t);
1410       xbt_dynar_insert_at(recvs, i, t);
1411       xbt_free(t);
1412     }
1413   }
1414   int rank_traced = smpi_process_index();
1415   TRACE_smpi_computing_out(rank_traced);
1416
1417   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1418 #endif
1419   int retval = smpi_mpi_waitall(count, requests, status);
1420 #ifdef HAVE_TRACING
1421   for (i = 0; i < count; i++) {
1422     int src_traced, dst_traced, is_wait_for_receive;
1423     xbt_dynar_get_cpy(srcs, i, &src_traced);
1424     xbt_dynar_get_cpy(dsts, i, &dst_traced);
1425     xbt_dynar_get_cpy(recvs, i, &is_wait_for_receive);
1426     if (is_wait_for_receive) {
1427       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1428     }
1429   }
1430   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1431   //clean-up of dynars
1432   xbt_dynar_free(&srcs);
1433   xbt_dynar_free(&dsts);
1434   xbt_dynar_free(&recvs);
1435   TRACE_smpi_computing_in(rank_traced);
1436 #endif
1437   smpi_bench_begin();
1438   return retval;
1439 }
1440
1441 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1442                  int *indices, MPI_Status status[])
1443 {
1444   int retval;
1445
1446   smpi_bench_end();
1447   if (outcount == NULL || indices == NULL) {
1448     retval = MPI_ERR_ARG;
1449   } else {
1450     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1451     retval = MPI_SUCCESS;
1452   }
1453   smpi_bench_begin();
1454   return retval;
1455 }
1456
1457 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1458                  int* indices, MPI_Status status[])
1459 {
1460   int retval;
1461
1462    smpi_bench_end();
1463    if (outcount == NULL || indices == NULL) {
1464      retval = MPI_ERR_ARG;
1465    } else {
1466      *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1467      retval = MPI_SUCCESS;
1468    }
1469    smpi_bench_begin();
1470    return retval;
1471 }
1472
1473
1474 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1475 {
1476   int retval;
1477
1478   smpi_bench_end();
1479 #ifdef HAVE_TRACING
1480   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1481   TRACE_smpi_computing_out(rank);
1482   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1483   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1484 #endif
1485   if (comm == MPI_COMM_NULL) {
1486     retval = MPI_ERR_COMM;
1487   } else {
1488     smpi_mpi_bcast(buf, count, datatype, root, comm);
1489     retval = MPI_SUCCESS;
1490   }
1491 #ifdef HAVE_TRACING
1492   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1493   TRACE_smpi_computing_in(rank);
1494 #endif
1495   smpi_bench_begin();
1496   return retval;
1497 }
1498
1499 int PMPI_Barrier(MPI_Comm comm)
1500 {
1501   int retval;
1502
1503   smpi_bench_end();
1504 #ifdef HAVE_TRACING
1505   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1506   TRACE_smpi_computing_out(rank);
1507   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1508 #endif
1509   if (comm == MPI_COMM_NULL) {
1510     retval = MPI_ERR_COMM;
1511   } else {
1512     smpi_mpi_barrier(comm);
1513     retval = MPI_SUCCESS;
1514   }
1515 #ifdef HAVE_TRACING
1516   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1517   TRACE_smpi_computing_in(rank);
1518 #endif
1519   smpi_bench_begin();
1520   return retval;
1521 }
1522
1523 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1524                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1525                int root, MPI_Comm comm)
1526 {
1527   int retval;
1528
1529   smpi_bench_end();
1530 #ifdef HAVE_TRACING
1531   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1532   TRACE_smpi_computing_out(rank);
1533   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1534   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1535 #endif
1536   if (comm == MPI_COMM_NULL) {
1537     retval = MPI_ERR_COMM;
1538   } else if (sendtype == MPI_DATATYPE_NULL
1539              || recvtype == MPI_DATATYPE_NULL) {
1540     retval = MPI_ERR_TYPE;
1541   } else {
1542     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1543                     recvtype, root, comm);
1544     retval = MPI_SUCCESS;
1545   }
1546 #ifdef HAVE_TRACING
1547   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1548   TRACE_smpi_computing_in(rank);
1549 #endif
1550   smpi_bench_begin();
1551   return retval;
1552 }
1553
1554 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1555                 void *recvbuf, int *recvcounts, int *displs,
1556                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1557 {
1558   int retval;
1559
1560   smpi_bench_end();
1561 #ifdef HAVE_TRACING
1562   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1563   TRACE_smpi_computing_out(rank);
1564   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1565   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1566 #endif
1567   if (comm == MPI_COMM_NULL) {
1568     retval = MPI_ERR_COMM;
1569   } else if (sendtype == MPI_DATATYPE_NULL
1570              || recvtype == MPI_DATATYPE_NULL) {
1571     retval = MPI_ERR_TYPE;
1572   } else if (recvcounts == NULL || displs == NULL) {
1573     retval = MPI_ERR_ARG;
1574   } else {
1575     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1576                      displs, recvtype, root, comm);
1577     retval = MPI_SUCCESS;
1578   }
1579 #ifdef HAVE_TRACING
1580   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1581   TRACE_smpi_computing_in(rank);
1582 #endif
1583   smpi_bench_begin();
1584   return retval;
1585 }
1586
1587 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1588                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1589                   MPI_Comm comm)
1590 {
1591   int retval;
1592
1593   smpi_bench_end();
1594 #ifdef HAVE_TRACING
1595   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1596   TRACE_smpi_computing_out(rank);
1597   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1598 #endif
1599   if (comm == MPI_COMM_NULL) {
1600     retval = MPI_ERR_COMM;
1601   } else if (sendtype == MPI_DATATYPE_NULL
1602              || recvtype == MPI_DATATYPE_NULL) {
1603     retval = MPI_ERR_TYPE;
1604   } else {
1605     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1606                        recvtype, comm);
1607     retval = MPI_SUCCESS;
1608   }
1609 #ifdef HAVE_TRACING
1610   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1611 #endif
1612   smpi_bench_begin();
1613   return retval;
1614 }
1615
1616 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1617                    void *recvbuf, int *recvcounts, int *displs,
1618                    MPI_Datatype recvtype, MPI_Comm comm)
1619 {
1620   int retval;
1621
1622   smpi_bench_end();
1623 #ifdef HAVE_TRACING
1624   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1625   TRACE_smpi_computing_out(rank);
1626   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1627 #endif
1628   if (comm == MPI_COMM_NULL) {
1629     retval = MPI_ERR_COMM;
1630   } else if (sendtype == MPI_DATATYPE_NULL
1631              || recvtype == MPI_DATATYPE_NULL) {
1632     retval = MPI_ERR_TYPE;
1633   } else if (recvcounts == NULL || displs == NULL) {
1634     retval = MPI_ERR_ARG;
1635   } else {
1636     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1637                         displs, recvtype, comm);
1638     retval = MPI_SUCCESS;
1639   }
1640 #ifdef HAVE_TRACING
1641   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1642   TRACE_smpi_computing_in(rank);
1643 #endif
1644   smpi_bench_begin();
1645   return retval;
1646 }
1647
1648 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1649                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1650                 int root, MPI_Comm comm)
1651 {
1652   int retval;
1653
1654   smpi_bench_end();
1655 #ifdef HAVE_TRACING
1656   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1657   TRACE_smpi_computing_out(rank);
1658   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1659
1660   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1661 #endif
1662   if (comm == MPI_COMM_NULL) {
1663     retval = MPI_ERR_COMM;
1664   } else if (sendtype == MPI_DATATYPE_NULL
1665              || recvtype == MPI_DATATYPE_NULL) {
1666     retval = MPI_ERR_TYPE;
1667   } else {
1668     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1669                      recvtype, root, comm);
1670     retval = MPI_SUCCESS;
1671   }
1672 #ifdef HAVE_TRACING
1673   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1674   TRACE_smpi_computing_in(rank);
1675 #endif
1676   smpi_bench_begin();
1677   return retval;
1678 }
1679
1680 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1681                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1682                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1683 {
1684   int retval;
1685
1686   smpi_bench_end();
1687 #ifdef HAVE_TRACING
1688   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1689   TRACE_smpi_computing_out(rank);
1690   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1691   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1692 #endif
1693   if (comm == MPI_COMM_NULL) {
1694     retval = MPI_ERR_COMM;
1695   } else if (sendtype == MPI_DATATYPE_NULL
1696              || recvtype == MPI_DATATYPE_NULL) {
1697     retval = MPI_ERR_TYPE;
1698   } else if (sendcounts == NULL || displs == NULL) {
1699     retval = MPI_ERR_ARG;
1700   } else {
1701     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1702                       recvcount, recvtype, root, comm);
1703     retval = MPI_SUCCESS;
1704   }
1705 #ifdef HAVE_TRACING
1706   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1707   TRACE_smpi_computing_in(rank);
1708 #endif
1709   smpi_bench_begin();
1710   return retval;
1711 }
1712
1713 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1714                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1715 {
1716   int retval;
1717
1718   smpi_bench_end();
1719 #ifdef HAVE_TRACING
1720   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1721   TRACE_smpi_computing_out(rank);
1722   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1723   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1724 #endif
1725   if (comm == MPI_COMM_NULL) {
1726     retval = MPI_ERR_COMM;
1727   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1728     retval = MPI_ERR_ARG;
1729   } else {
1730     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1731     retval = MPI_SUCCESS;
1732   }
1733 #ifdef HAVE_TRACING
1734   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1735   TRACE_smpi_computing_in(rank);
1736 #endif
1737   smpi_bench_begin();
1738   return retval;
1739 }
1740
1741 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1742                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1743 {
1744   int retval;
1745
1746   smpi_bench_end();
1747 #ifdef HAVE_TRACING
1748   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1749   TRACE_smpi_computing_out(rank);
1750   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1751 #endif
1752   if (comm == MPI_COMM_NULL) {
1753     retval = MPI_ERR_COMM;
1754   } else if (datatype == MPI_DATATYPE_NULL) {
1755     retval = MPI_ERR_TYPE;
1756   } else if (op == MPI_OP_NULL) {
1757     retval = MPI_ERR_OP;
1758   } else {
1759     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
1760     retval = MPI_SUCCESS;
1761   }
1762 #ifdef HAVE_TRACING
1763   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1764   TRACE_smpi_computing_in(rank);
1765 #endif
1766   smpi_bench_begin();
1767   return retval;
1768 }
1769
1770 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
1771              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1772 {
1773   int retval;
1774
1775   smpi_bench_end();
1776 #ifdef HAVE_TRACING
1777   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1778   TRACE_smpi_computing_out(rank);
1779   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1780 #endif
1781   if (comm == MPI_COMM_NULL) {
1782     retval = MPI_ERR_COMM;
1783   } else if (datatype == MPI_DATATYPE_NULL) {
1784     retval = MPI_ERR_TYPE;
1785   } else if (op == MPI_OP_NULL) {
1786     retval = MPI_ERR_OP;
1787   } else {
1788     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1789     retval = MPI_SUCCESS;
1790   }
1791 #ifdef HAVE_TRACING
1792   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1793   TRACE_smpi_computing_in(rank);
1794 #endif
1795   smpi_bench_begin();
1796   return retval;
1797 }
1798
1799 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1800                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1801 {
1802   int retval, i, size, count;
1803   int *displs;
1804   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1805
1806   smpi_bench_end();
1807 #ifdef HAVE_TRACING
1808   TRACE_smpi_computing_out(rank);
1809   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1810 #endif
1811   if (comm == MPI_COMM_NULL) {
1812     retval = MPI_ERR_COMM;
1813   } else if (datatype == MPI_DATATYPE_NULL) {
1814     retval = MPI_ERR_TYPE;
1815   } else if (op == MPI_OP_NULL) {
1816     retval = MPI_ERR_OP;
1817   } else if (recvcounts == NULL) {
1818     retval = MPI_ERR_ARG;
1819   } else {
1820     /* arbitrarily choose root as rank 0 */
1821     /* TODO: faster direct implementation ? */
1822     size = smpi_comm_size(comm);
1823     count = 0;
1824     displs = xbt_new(int, size);
1825     for (i = 0; i < size; i++) {
1826       count += recvcounts[i];
1827       displs[i] = 0;
1828     }
1829     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1830     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf,
1831                       recvcounts[rank], datatype, 0, comm);
1832     xbt_free(displs);
1833     retval = MPI_SUCCESS;
1834   }
1835 #ifdef HAVE_TRACING
1836   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1837   TRACE_smpi_computing_in(rank);
1838 #endif
1839   smpi_bench_begin();
1840   return retval;
1841 }
1842
1843 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1844                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
1845                  MPI_Comm comm)
1846 {
1847   int retval, size, sendsize;
1848
1849   smpi_bench_end();
1850 #ifdef HAVE_TRACING
1851   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1852   TRACE_smpi_computing_out(rank);
1853   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1854 #endif
1855   if (comm == MPI_COMM_NULL) {
1856     retval = MPI_ERR_COMM;
1857   } else if (sendtype == MPI_DATATYPE_NULL
1858              || recvtype == MPI_DATATYPE_NULL) {
1859     retval = MPI_ERR_TYPE;
1860   } else {
1861     size = smpi_comm_size(comm);
1862     sendsize = smpi_datatype_size(sendtype) * sendcount;
1863     if (sendsize < 200 && size > 12) {
1864       retval =
1865           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
1866                                          recvbuf, recvcount, recvtype,
1867                                          comm);
1868     } else if (sendsize < 3000) {
1869       retval =
1870           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
1871                                                 sendtype, recvbuf,
1872                                                 recvcount, recvtype, comm);
1873     } else {
1874       retval =
1875           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
1876                                             recvbuf, recvcount, recvtype,
1877                                             comm);
1878     }
1879   }
1880 #ifdef HAVE_TRACING
1881   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1882   TRACE_smpi_computing_in(rank);
1883 #endif
1884   smpi_bench_begin();
1885   return retval;
1886 }
1887
1888 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
1889                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
1890                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
1891 {
1892   int retval;
1893
1894   smpi_bench_end();
1895 #ifdef HAVE_TRACING
1896   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1897   TRACE_smpi_computing_out(rank);
1898   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1899 #endif
1900   if (comm == MPI_COMM_NULL) {
1901     retval = MPI_ERR_COMM;
1902   } else if (sendtype == MPI_DATATYPE_NULL
1903              || recvtype == MPI_DATATYPE_NULL) {
1904     retval = MPI_ERR_TYPE;
1905   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
1906              || recvdisps == NULL) {
1907     retval = MPI_ERR_ARG;
1908   } else {
1909     retval =
1910         smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
1911                                   recvbuf, recvcounts, recvdisps, recvtype,
1912                                   comm);
1913   }
1914 #ifdef HAVE_TRACING
1915   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1916   TRACE_smpi_computing_in(rank);
1917 #endif
1918   smpi_bench_begin();
1919   return retval;
1920 }
1921
1922
1923 int PMPI_Get_processor_name(char *name, int *resultlen)
1924 {
1925   int retval = MPI_SUCCESS;
1926
1927   smpi_bench_end();
1928   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
1929           strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
1930           strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
1931           MPI_MAX_PROCESSOR_NAME - 1 );
1932   *resultlen =
1933       strlen(name) >
1934       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1935
1936   smpi_bench_begin();
1937   return retval;
1938 }
1939
1940 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
1941 {
1942   int retval = MPI_SUCCESS;
1943   size_t size;
1944
1945   smpi_bench_end();
1946   if (status == NULL || count == NULL) {
1947     retval = MPI_ERR_ARG;
1948   } else if (datatype == MPI_DATATYPE_NULL) {
1949     retval = MPI_ERR_TYPE;
1950   } else {
1951     size = smpi_datatype_size(datatype);
1952     if (size == 0) {
1953       *count = 0;
1954     } else if (status->count % size != 0) {
1955       retval = MPI_UNDEFINED;
1956     } else {
1957       *count = smpi_mpi_get_count(status, datatype);
1958     }
1959   }
1960   smpi_bench_begin();
1961   return retval;
1962 }
1963
1964 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
1965   int retval;
1966
1967   smpi_bench_end();
1968   if (old_type == MPI_DATATYPE_NULL) {
1969     retval = MPI_ERR_TYPE;
1970   } else if (count<0){
1971     retval = MPI_ERR_COUNT;
1972   } else {
1973     retval = smpi_datatype_contiguous(count, old_type, new_type);
1974   }
1975   smpi_bench_begin();
1976   return retval;
1977 }
1978
1979 int PMPI_Type_commit(MPI_Datatype* datatype) {
1980   int retval;
1981
1982   smpi_bench_end();
1983   if (datatype == MPI_DATATYPE_NULL) {
1984     retval = MPI_ERR_TYPE;
1985   } else {
1986     smpi_datatype_commit(datatype);
1987     retval = MPI_SUCCESS;
1988   }
1989   smpi_bench_begin();
1990   return retval;
1991 }
1992
1993
1994 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
1995   int retval;
1996
1997   smpi_bench_end();
1998   if (old_type == MPI_DATATYPE_NULL) {
1999     retval = MPI_ERR_TYPE;
2000   } else if (count<0 || blocklen<0){
2001     retval = MPI_ERR_COUNT;
2002   } else {
2003     retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2004   }
2005   smpi_bench_begin();
2006   return retval;
2007 }
2008
2009 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2010   int retval;
2011
2012   smpi_bench_end();
2013   if (old_type == MPI_DATATYPE_NULL) {
2014     retval = MPI_ERR_TYPE;
2015   } else if (count<0 || blocklen<0){
2016     retval = MPI_ERR_COUNT;
2017   } else {
2018     retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2019   }
2020   smpi_bench_begin();
2021   return retval;
2022 }
2023
2024
2025 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2026   int retval;
2027
2028   smpi_bench_end();
2029   if (old_type == MPI_DATATYPE_NULL) {
2030     retval = MPI_ERR_TYPE;
2031   } else if (count<0){
2032     retval = MPI_ERR_COUNT;
2033   } else {
2034     retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2035   }
2036   smpi_bench_begin();
2037   return retval;
2038 }
2039
2040 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2041   int retval;
2042
2043   smpi_bench_end();
2044   if (old_type == MPI_DATATYPE_NULL) {
2045     retval = MPI_ERR_TYPE;
2046   } else if (count<0){
2047     retval = MPI_ERR_COUNT;
2048   } else {
2049     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2050   }
2051   smpi_bench_begin();
2052   return retval;
2053 }
2054
2055
2056 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2057   int retval;
2058
2059   smpi_bench_end();
2060   if (count<0){
2061     retval = MPI_ERR_COUNT;
2062   } else {
2063     retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2064   }
2065   smpi_bench_begin();
2066   return retval;}
2067
2068 int PMPI_Error_class(int errorcode, int* errorclass) {
2069   // assume smpi uses only standard mpi error codes
2070   *errorclass=errorcode;
2071   return MPI_SUCCESS;
2072 }
2073
2074
2075 int PMPI_Initialized(int* flag) {
2076    *flag=(smpi_process_data()!=NULL);
2077    return MPI_SUCCESS;
2078 }
2079
2080 /* The following calls are not yet implemented and will fail at runtime. */
2081 /* Once implemented, please move them above this notice. */
2082
2083 static int not_yet_implemented(void) {
2084           XBT_WARN("Not yet implemented");
2085    return MPI_SUCCESS;
2086 }
2087
2088 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2089    return not_yet_implemented();
2090 }
2091
2092 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2093    return not_yet_implemented();
2094 }
2095
2096 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
2097    return not_yet_implemented();
2098 }
2099
2100 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2101    return not_yet_implemented();
2102 }
2103
2104 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2105    return not_yet_implemented();
2106 }
2107
2108 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2109    return not_yet_implemented();
2110 }
2111
2112 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2113    return not_yet_implemented();
2114 }
2115
2116 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2117    return not_yet_implemented();
2118 }
2119
2120 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2121    return not_yet_implemented();
2122 }
2123
2124 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2125    return not_yet_implemented();
2126 }
2127
2128 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2129    return not_yet_implemented();
2130 }
2131
2132 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2133    return not_yet_implemented();
2134 }
2135
2136 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2137    return not_yet_implemented();
2138 }
2139
2140 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2141    return not_yet_implemented();
2142 }
2143
2144 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2145    return not_yet_implemented();
2146 }
2147
2148 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2149    return not_yet_implemented();
2150 }
2151
2152 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2153    return not_yet_implemented();
2154 }
2155
2156 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2157    return not_yet_implemented();
2158 }
2159
2160 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2161    return not_yet_implemented();
2162 }
2163
2164 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2165    return not_yet_implemented();
2166 }
2167
2168 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2169    return not_yet_implemented();
2170 }
2171
2172
2173 int PMPI_Cancel(MPI_Request* request) {
2174    return not_yet_implemented();
2175 }
2176
2177 int PMPI_Buffer_attach(void* buffer, int size) {
2178    return not_yet_implemented();
2179 }
2180
2181 int PMPI_Buffer_detach(void* buffer, int* size) {
2182    return not_yet_implemented();
2183 }
2184
2185 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2186    return not_yet_implemented();
2187 }
2188
2189 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2190 {
2191    return not_yet_implemented();
2192 }
2193
2194 int PMPI_Pcontrol(const int level )
2195 {
2196    return not_yet_implemented();
2197 }
2198
2199 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2200    return not_yet_implemented();
2201 }
2202
2203 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2204    return not_yet_implemented();
2205 }
2206
2207 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2208    return not_yet_implemented();
2209 }
2210
2211 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2212    return not_yet_implemented();
2213 }
2214
2215 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2216    return not_yet_implemented();
2217 }
2218
2219 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2220    return not_yet_implemented();
2221 }
2222
2223 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2224    return not_yet_implemented();
2225 }
2226
2227 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2228    return not_yet_implemented();
2229 }
2230
2231 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2232    return not_yet_implemented();
2233 }
2234
2235 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2236    return not_yet_implemented();
2237 }
2238
2239 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2240    return not_yet_implemented();
2241 }
2242
2243
2244 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2245    return not_yet_implemented();
2246 }
2247
2248 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2249    return not_yet_implemented();
2250 }
2251
2252 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2253    return not_yet_implemented();
2254 }
2255
2256 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2257    return not_yet_implemented();
2258 }
2259
2260 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2261    return not_yet_implemented();
2262 }
2263
2264 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2265    return not_yet_implemented();
2266 }
2267
2268 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2269    return not_yet_implemented();
2270 }
2271
2272 int PMPI_Keyval_free(int* keyval) {
2273    return not_yet_implemented();
2274 }
2275
2276 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2277    return not_yet_implemented();
2278 }
2279
2280 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2281    return not_yet_implemented();
2282 }
2283
2284 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2285    return not_yet_implemented();
2286 }
2287
2288 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2289    return not_yet_implemented();
2290 }
2291
2292