Logo AND Algorithmique Numérique Distribuée

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