Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
abd6b90d14f9b21d24f9b36d5f2c26e987bb8120
[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     smpi_mpi_request_free(request);
904     retval = MPI_SUCCESS;
905   }
906   smpi_bench_begin();
907   return retval;
908 }
909
910 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
911               int tag, MPI_Comm comm, MPI_Request * request)
912 {
913   int retval;
914
915   smpi_bench_end();
916
917   if (request == NULL) {
918     retval = MPI_ERR_ARG;
919   } else if (comm == MPI_COMM_NULL) {
920     retval = MPI_ERR_COMM;
921   } else if (src == MPI_PROC_NULL) {
922     *request = MPI_REQUEST_NULL;
923     retval = MPI_SUCCESS;
924   } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
925     retval = MPI_ERR_COMM;
926   } else if (count < 0) {
927     retval = MPI_ERR_COUNT;
928   } else if (buf==NULL && count > 0) {
929     retval = MPI_ERR_COUNT;
930   } else if (datatype == MPI_DATATYPE_NULL){
931     retval = MPI_ERR_TYPE;
932   } else if(tag<0 && tag !=  MPI_ANY_TAG){
933     retval = MPI_ERR_TAG;
934   } else {
935
936 #ifdef HAVE_TRACING
937   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
938   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
939   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
940 #endif
941
942     *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
943     retval = MPI_SUCCESS;
944
945 #ifdef HAVE_TRACING
946   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
947   (*request)->recv = 1;
948 #endif
949   }
950
951   smpi_bench_begin();
952   return retval;
953 }
954
955
956 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
957               int tag, MPI_Comm comm, MPI_Request * request)
958 {
959   int retval;
960
961   smpi_bench_end();
962   if (request == NULL) {
963     retval = MPI_ERR_ARG;
964   } else if (comm == MPI_COMM_NULL) {
965     retval = MPI_ERR_COMM;
966   } else if (dst == MPI_PROC_NULL) {
967     *request = MPI_REQUEST_NULL;
968     retval = MPI_SUCCESS;
969   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
970     retval = MPI_ERR_COMM;
971   } else if (count < 0) {
972     retval = MPI_ERR_COUNT;
973   } else if (buf==NULL && count > 0) {
974     retval = MPI_ERR_COUNT;
975   } else if (datatype == MPI_DATATYPE_NULL){
976     retval = MPI_ERR_TYPE;
977   } else if(tag<0 && tag !=  MPI_ANY_TAG){
978     retval = MPI_ERR_TAG;
979   } else {
980
981 #ifdef HAVE_TRACING
982   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
983   TRACE_smpi_computing_out(rank);
984   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
985   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
986   TRACE_smpi_send(rank, rank, dst_traced);
987 #endif
988
989     *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
990     retval = MPI_SUCCESS;
991
992 #ifdef HAVE_TRACING
993   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
994   (*request)->send = 1;
995   TRACE_smpi_computing_in(rank);
996 #endif
997   }
998
999   smpi_bench_begin();
1000   return retval;
1001 }
1002
1003
1004
1005 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1006              MPI_Comm comm, MPI_Status * status)
1007 {
1008   int retval;
1009
1010   smpi_bench_end();
1011   if (comm == MPI_COMM_NULL) {
1012     retval = MPI_ERR_COMM;
1013   } else if (src == MPI_PROC_NULL) {
1014     smpi_empty_status(status);
1015     status->MPI_SOURCE = MPI_PROC_NULL;
1016     retval = MPI_SUCCESS;
1017   } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1018     retval = MPI_ERR_COMM;
1019   } else if (count < 0) {
1020     retval = MPI_ERR_COUNT;
1021   } else if (buf==NULL && count > 0) {
1022     retval = MPI_ERR_COUNT;
1023   } else if (datatype == MPI_DATATYPE_NULL){
1024     retval = MPI_ERR_TYPE;
1025   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1026     retval = MPI_ERR_TAG;
1027   } else {
1028 #ifdef HAVE_TRACING
1029   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1030   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1031   TRACE_smpi_computing_out(rank);
1032
1033   TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
1034 #endif
1035
1036     smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1037     retval = MPI_SUCCESS;
1038
1039 #ifdef HAVE_TRACING
1040   //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1041   if(status!=MPI_STATUS_IGNORE)src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1042   TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1043   TRACE_smpi_recv(rank, src_traced, rank);
1044   TRACE_smpi_computing_in(rank);
1045 #endif
1046   }
1047
1048   smpi_bench_begin();
1049   return retval;
1050 }
1051
1052 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1053              MPI_Comm comm)
1054 {
1055   int retval;
1056
1057   smpi_bench_end();
1058
1059   if (comm == MPI_COMM_NULL) {
1060     retval = MPI_ERR_COMM;
1061   } else if (dst == MPI_PROC_NULL) {
1062     retval = MPI_SUCCESS;
1063   } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1064     retval = MPI_ERR_COMM;
1065   } else if (count < 0) {
1066     retval = MPI_ERR_COUNT;
1067   } else if (buf==NULL && count > 0) {
1068     retval = MPI_ERR_COUNT;
1069   } else if (datatype == MPI_DATATYPE_NULL){
1070     retval = MPI_ERR_TYPE;
1071   } else if(tag<0 && tag !=  MPI_ANY_TAG){
1072     retval = MPI_ERR_TAG;
1073   } else {
1074
1075 #ifdef HAVE_TRACING
1076   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1077   TRACE_smpi_computing_out(rank);
1078   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1079   TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1080   TRACE_smpi_send(rank, rank, dst_traced);
1081 #endif
1082
1083     smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1084     retval = MPI_SUCCESS;
1085
1086 #ifdef HAVE_TRACING
1087   TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1088   TRACE_smpi_computing_in(rank);
1089 #endif
1090   }
1091
1092   smpi_bench_begin();
1093   return retval;
1094 }
1095
1096 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1097                  int dst, int sendtag, void *recvbuf, int recvcount,
1098                  MPI_Datatype recvtype, int src, int recvtag,
1099                  MPI_Comm comm, MPI_Status * status)
1100 {
1101   int retval;
1102
1103   smpi_bench_end();
1104
1105   if (comm == MPI_COMM_NULL) {
1106     retval = MPI_ERR_COMM;
1107   } else if (sendtype == MPI_DATATYPE_NULL
1108              || recvtype == MPI_DATATYPE_NULL) {
1109     retval = MPI_ERR_TYPE;
1110   } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1111       smpi_empty_status(status);
1112       status->MPI_SOURCE = MPI_PROC_NULL;
1113       retval = MPI_SUCCESS;
1114   }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1115       (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1116     retval = MPI_ERR_COMM;
1117   } else if (sendcount < 0 || recvcount<0) {
1118       retval = MPI_ERR_COUNT;
1119   } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1120     retval = MPI_ERR_COUNT;
1121   } else if((sendtag<0 && sendtag !=  MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1122     retval = MPI_ERR_TAG;
1123   } else {
1124
1125 #ifdef HAVE_TRACING
1126   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1127   TRACE_smpi_computing_out(rank);
1128   int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1129   int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1130   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1131   TRACE_smpi_send(rank, rank, dst_traced);
1132   TRACE_smpi_send(rank, src_traced, rank);
1133 #endif
1134
1135
1136     smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1137                       recvcount, recvtype, src, recvtag, comm, status);
1138     retval = MPI_SUCCESS;
1139
1140 #ifdef HAVE_TRACING
1141   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1142   TRACE_smpi_recv(rank, rank, dst_traced);
1143   TRACE_smpi_recv(rank, src_traced, rank);
1144   TRACE_smpi_computing_in(rank);
1145 #endif
1146
1147   }
1148
1149   smpi_bench_begin();
1150   return retval;
1151 }
1152
1153 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1154                          int dst, int sendtag, int src, int recvtag,
1155                          MPI_Comm comm, MPI_Status * status)
1156 {
1157   //TODO: suboptimal implementation
1158   void *recvbuf;
1159   int retval;
1160   if ((datatype == MPI_DATATYPE_NULL)||(datatype->has_subtype==1)) {
1161       retval = MPI_ERR_TYPE;
1162   } else if (count < 0) {
1163       retval = MPI_ERR_COUNT;
1164   } else {
1165     int size = smpi_datatype_size(datatype) * count;
1166     recvbuf = xbt_new(char, size);
1167     retval =
1168         MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1169                      datatype, src, recvtag, comm, status);
1170     if(retval==MPI_SUCCESS){
1171         memcpy(buf, recvbuf, size * sizeof(char));
1172     }
1173     xbt_free(recvbuf);
1174
1175   }
1176   return retval;
1177 }
1178
1179 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1180 {
1181   int retval;
1182
1183   smpi_bench_end();
1184   if (request == MPI_REQUEST_NULL || flag == NULL) {
1185     retval = MPI_ERR_ARG;
1186   } else if (*request == MPI_REQUEST_NULL) {
1187     *flag= TRUE;
1188     retval = MPI_ERR_REQUEST;
1189   } else {
1190     *flag = smpi_mpi_test(request, status);
1191     retval = MPI_SUCCESS;
1192   }
1193   smpi_bench_begin();
1194   return retval;
1195 }
1196
1197 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1198                 MPI_Status * status)
1199 {
1200   int retval;
1201
1202   smpi_bench_end();
1203   if (index == NULL || flag == NULL) {
1204     retval = MPI_ERR_ARG;
1205   } else {
1206     *flag = smpi_mpi_testany(count, requests, index, status);
1207     retval = MPI_SUCCESS;
1208   }
1209   smpi_bench_begin();
1210   return retval;
1211 }
1212
1213 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1214 {
1215   int retval;
1216
1217   smpi_bench_end();
1218   if (flag == NULL) {
1219     retval = MPI_ERR_ARG;
1220   } else {
1221     *flag = smpi_mpi_testall(count, requests, statuses);
1222     retval = MPI_SUCCESS;
1223   }
1224   smpi_bench_begin();
1225   return retval;
1226 }
1227
1228 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1229   int retval;
1230   smpi_bench_end();
1231
1232   if (status == NULL) {
1233     retval = MPI_ERR_ARG;
1234   } else if (comm == MPI_COMM_NULL) {
1235     retval = MPI_ERR_COMM;
1236   } else if (source == MPI_PROC_NULL) {
1237     smpi_empty_status(status);
1238     status->MPI_SOURCE = MPI_PROC_NULL;
1239     retval = MPI_SUCCESS;
1240   } else {
1241     smpi_mpi_probe(source, tag, comm, status);
1242     retval = MPI_SUCCESS;
1243   }
1244   smpi_bench_begin();
1245   return retval;
1246 }
1247
1248
1249 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1250   int retval;
1251   smpi_bench_end();
1252
1253   if (flag == NULL) {
1254     retval = MPI_ERR_ARG;
1255   } else if (status == NULL) {
1256     retval = MPI_ERR_ARG;
1257   } else if (comm == MPI_COMM_NULL) {
1258     retval = MPI_ERR_COMM;
1259   } else if (source == MPI_PROC_NULL) {
1260     smpi_empty_status(status);
1261     status->MPI_SOURCE = MPI_PROC_NULL;
1262     retval = MPI_SUCCESS;
1263   } else {
1264     smpi_mpi_iprobe(source, tag, comm, flag, status);
1265     retval = MPI_SUCCESS;
1266   }
1267   smpi_bench_begin();
1268   return retval;
1269 }
1270
1271 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1272 {
1273   int retval;
1274
1275   smpi_bench_end();
1276
1277   if (request == NULL) {
1278     retval = MPI_ERR_ARG;
1279   } else if (*request == MPI_REQUEST_NULL) {
1280     retval = MPI_ERR_REQUEST;
1281   } else {
1282
1283 #ifdef HAVE_TRACING
1284   int rank = request && (*request)->comm != MPI_COMM_NULL
1285       ? smpi_process_index()
1286       : -1;
1287   TRACE_smpi_computing_out(rank);
1288
1289   MPI_Group group = smpi_comm_group((*request)->comm);
1290   int src_traced = smpi_group_index(group, (*request)->src);
1291   int dst_traced = smpi_group_index(group, (*request)->dst);
1292   int is_wait_for_receive = (*request)->recv;
1293   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1294 #endif
1295
1296     smpi_mpi_wait(request, status);
1297     retval = MPI_SUCCESS;
1298
1299 #ifdef HAVE_TRACING
1300   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1301   if (is_wait_for_receive) {
1302     TRACE_smpi_recv(rank, src_traced, dst_traced);
1303   }
1304   TRACE_smpi_computing_in(rank);
1305 #endif
1306
1307   }
1308
1309   smpi_bench_begin();
1310   return retval;
1311 }
1312
1313 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1314 {
1315   int retval;
1316
1317   smpi_bench_end();
1318 #ifdef HAVE_TRACING
1319   //save requests information for tracing
1320   int i;
1321   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), NULL);
1322   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), NULL);
1323   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), NULL);
1324   for (i = 0; i < count; i++) {
1325     MPI_Request req = requests[i];      //already received requests are no longer valid
1326     if (req) {
1327       int *asrc = xbt_new(int, 1);
1328       int *adst = xbt_new(int, 1);
1329       int *arecv = xbt_new(int, 1);
1330       *asrc = req->src;
1331       *adst = req->dst;
1332       *arecv = req->recv;
1333       xbt_dynar_insert_at(srcs, i, asrc);
1334       xbt_dynar_insert_at(dsts, i, adst);
1335       xbt_dynar_insert_at(recvs, i, arecv);
1336       xbt_free(asrc);
1337       xbt_free(adst);
1338       xbt_free(arecv);
1339     } else {
1340       int *t = xbt_new(int, 1);
1341       xbt_dynar_insert_at(srcs, i, t);
1342       xbt_dynar_insert_at(dsts, i, t);
1343       xbt_dynar_insert_at(recvs, i, t);
1344       xbt_free(t);
1345     }
1346   }
1347   int rank_traced = smpi_process_index();
1348   TRACE_smpi_computing_out(rank_traced);
1349
1350   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1351
1352 #endif
1353   if (index == NULL) {
1354     retval = MPI_ERR_ARG;
1355   } else {
1356     *index = smpi_mpi_waitany(count, requests, status);
1357     retval = MPI_SUCCESS;
1358   }
1359 #ifdef HAVE_TRACING
1360   if(*index!=MPI_UNDEFINED){
1361     int src_traced, dst_traced, is_wait_for_receive;
1362     xbt_dynar_get_cpy(srcs, *index, &src_traced);
1363     xbt_dynar_get_cpy(dsts, *index, &dst_traced);
1364     xbt_dynar_get_cpy(recvs, *index, &is_wait_for_receive);
1365     if (is_wait_for_receive) {
1366       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1367     }
1368     TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1369     //clean-up of dynars
1370     xbt_dynar_free(&srcs);
1371     xbt_dynar_free(&dsts);
1372     xbt_dynar_free(&recvs);
1373   }
1374   TRACE_smpi_computing_in(rank_traced);
1375 #endif
1376   smpi_bench_begin();
1377   return retval;
1378 }
1379
1380 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1381 {
1382
1383   smpi_bench_end();
1384 #ifdef HAVE_TRACING
1385   //save information from requests
1386   int i;
1387   xbt_dynar_t srcs = xbt_dynar_new(sizeof(int), NULL);
1388   xbt_dynar_t dsts = xbt_dynar_new(sizeof(int), NULL);
1389   xbt_dynar_t recvs = xbt_dynar_new(sizeof(int), NULL);
1390   for (i = 0; i < count; i++) {
1391     MPI_Request req = requests[i];
1392     if(req){
1393       int *asrc = xbt_new(int, 1);
1394       int *adst = xbt_new(int, 1);
1395       int *arecv = xbt_new(int, 1);
1396       *asrc = req->src;
1397       *adst = req->dst;
1398       *arecv = req->recv;
1399       xbt_dynar_insert_at(srcs, i, asrc);
1400       xbt_dynar_insert_at(dsts, i, adst);
1401       xbt_dynar_insert_at(recvs, i, arecv);
1402       xbt_free(asrc);
1403       xbt_free(adst);
1404       xbt_free(arecv);
1405     }else {
1406       int *t = xbt_new(int, 1);
1407       xbt_dynar_insert_at(srcs, i, t);
1408       xbt_dynar_insert_at(dsts, i, t);
1409       xbt_dynar_insert_at(recvs, i, t);
1410       xbt_free(t);
1411     }
1412   }
1413   int rank_traced = smpi_process_index();
1414   TRACE_smpi_computing_out(rank_traced);
1415
1416   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1417 #endif
1418   int retval = smpi_mpi_waitall(count, requests, status);
1419 #ifdef HAVE_TRACING
1420   for (i = 0; i < count; i++) {
1421     int src_traced, dst_traced, is_wait_for_receive;
1422     xbt_dynar_get_cpy(srcs, i, &src_traced);
1423     xbt_dynar_get_cpy(dsts, i, &dst_traced);
1424     xbt_dynar_get_cpy(recvs, i, &is_wait_for_receive);
1425     if (is_wait_for_receive) {
1426       TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1427     }
1428   }
1429   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1430   //clean-up of dynars
1431   xbt_dynar_free(&srcs);
1432   xbt_dynar_free(&dsts);
1433   xbt_dynar_free(&recvs);
1434   TRACE_smpi_computing_in(rank_traced);
1435 #endif
1436   smpi_bench_begin();
1437   return retval;
1438 }
1439
1440 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1441                  int *indices, MPI_Status status[])
1442 {
1443   int retval;
1444
1445   smpi_bench_end();
1446   if (outcount == NULL || indices == NULL) {
1447     retval = MPI_ERR_ARG;
1448   } else {
1449     *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1450     retval = MPI_SUCCESS;
1451   }
1452   smpi_bench_begin();
1453   return retval;
1454 }
1455
1456 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1457                  int* indices, MPI_Status status[])
1458 {
1459   int retval;
1460
1461    smpi_bench_end();
1462    if (outcount == NULL || indices == NULL) {
1463      retval = MPI_ERR_ARG;
1464    } else {
1465      *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1466      retval = MPI_SUCCESS;
1467    }
1468    smpi_bench_begin();
1469    return retval;
1470 }
1471
1472
1473 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1474 {
1475   int retval;
1476
1477   smpi_bench_end();
1478 #ifdef HAVE_TRACING
1479   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1480   TRACE_smpi_computing_out(rank);
1481   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1482   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1483 #endif
1484   if (comm == MPI_COMM_NULL) {
1485     retval = MPI_ERR_COMM;
1486   } else {
1487     smpi_mpi_bcast(buf, count, datatype, root, comm);
1488     retval = MPI_SUCCESS;
1489   }
1490 #ifdef HAVE_TRACING
1491   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1492   TRACE_smpi_computing_in(rank);
1493 #endif
1494   smpi_bench_begin();
1495   return retval;
1496 }
1497
1498 int PMPI_Barrier(MPI_Comm comm)
1499 {
1500   int retval;
1501
1502   smpi_bench_end();
1503 #ifdef HAVE_TRACING
1504   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1505   TRACE_smpi_computing_out(rank);
1506   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1507 #endif
1508   if (comm == MPI_COMM_NULL) {
1509     retval = MPI_ERR_COMM;
1510   } else {
1511     smpi_mpi_barrier(comm);
1512     retval = MPI_SUCCESS;
1513   }
1514 #ifdef HAVE_TRACING
1515   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1516   TRACE_smpi_computing_in(rank);
1517 #endif
1518   smpi_bench_begin();
1519   return retval;
1520 }
1521
1522 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1523                void *recvbuf, int recvcount, MPI_Datatype recvtype,
1524                int root, MPI_Comm comm)
1525 {
1526   int retval;
1527
1528   smpi_bench_end();
1529 #ifdef HAVE_TRACING
1530   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1531   TRACE_smpi_computing_out(rank);
1532   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1533   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1534 #endif
1535   if (comm == MPI_COMM_NULL) {
1536     retval = MPI_ERR_COMM;
1537   } else if (sendtype == MPI_DATATYPE_NULL
1538              || recvtype == MPI_DATATYPE_NULL) {
1539     retval = MPI_ERR_TYPE;
1540   } else {
1541     smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1542                     recvtype, root, comm);
1543     retval = MPI_SUCCESS;
1544   }
1545 #ifdef HAVE_TRACING
1546   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1547   TRACE_smpi_computing_in(rank);
1548 #endif
1549   smpi_bench_begin();
1550   return retval;
1551 }
1552
1553 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1554                 void *recvbuf, int *recvcounts, int *displs,
1555                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1556 {
1557   int retval;
1558
1559   smpi_bench_end();
1560 #ifdef HAVE_TRACING
1561   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1562   TRACE_smpi_computing_out(rank);
1563   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1564   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1565 #endif
1566   if (comm == MPI_COMM_NULL) {
1567     retval = MPI_ERR_COMM;
1568   } else if (sendtype == MPI_DATATYPE_NULL
1569              || recvtype == MPI_DATATYPE_NULL) {
1570     retval = MPI_ERR_TYPE;
1571   } else if (recvcounts == NULL || displs == NULL) {
1572     retval = MPI_ERR_ARG;
1573   } else {
1574     smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1575                      displs, recvtype, root, comm);
1576     retval = MPI_SUCCESS;
1577   }
1578 #ifdef HAVE_TRACING
1579   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1580   TRACE_smpi_computing_in(rank);
1581 #endif
1582   smpi_bench_begin();
1583   return retval;
1584 }
1585
1586 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1587                   void *recvbuf, int recvcount, MPI_Datatype recvtype,
1588                   MPI_Comm comm)
1589 {
1590   int retval;
1591
1592   smpi_bench_end();
1593 #ifdef HAVE_TRACING
1594   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1595   TRACE_smpi_computing_out(rank);
1596   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1597 #endif
1598   if (comm == MPI_COMM_NULL) {
1599     retval = MPI_ERR_COMM;
1600   } else if (sendtype == MPI_DATATYPE_NULL
1601              || recvtype == MPI_DATATYPE_NULL) {
1602     retval = MPI_ERR_TYPE;
1603   } else {
1604     smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1605                        recvtype, comm);
1606     retval = MPI_SUCCESS;
1607   }
1608 #ifdef HAVE_TRACING
1609   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1610 #endif
1611   smpi_bench_begin();
1612   return retval;
1613 }
1614
1615 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1616                    void *recvbuf, int *recvcounts, int *displs,
1617                    MPI_Datatype recvtype, MPI_Comm comm)
1618 {
1619   int retval;
1620
1621   smpi_bench_end();
1622 #ifdef HAVE_TRACING
1623   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1624   TRACE_smpi_computing_out(rank);
1625   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1626 #endif
1627   if (comm == MPI_COMM_NULL) {
1628     retval = MPI_ERR_COMM;
1629   } else if (sendtype == MPI_DATATYPE_NULL
1630              || recvtype == MPI_DATATYPE_NULL) {
1631     retval = MPI_ERR_TYPE;
1632   } else if (recvcounts == NULL || displs == NULL) {
1633     retval = MPI_ERR_ARG;
1634   } else {
1635     smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1636                         displs, recvtype, comm);
1637     retval = MPI_SUCCESS;
1638   }
1639 #ifdef HAVE_TRACING
1640   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1641   TRACE_smpi_computing_in(rank);
1642 #endif
1643   smpi_bench_begin();
1644   return retval;
1645 }
1646
1647 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1648                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1649                 int root, MPI_Comm comm)
1650 {
1651   int retval;
1652
1653   smpi_bench_end();
1654 #ifdef HAVE_TRACING
1655   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1656   TRACE_smpi_computing_out(rank);
1657   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1658
1659   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1660 #endif
1661   if (comm == MPI_COMM_NULL) {
1662     retval = MPI_ERR_COMM;
1663   } else if (sendtype == MPI_DATATYPE_NULL
1664              || recvtype == MPI_DATATYPE_NULL) {
1665     retval = MPI_ERR_TYPE;
1666   } else {
1667     smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1668                      recvtype, root, comm);
1669     retval = MPI_SUCCESS;
1670   }
1671 #ifdef HAVE_TRACING
1672   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1673   TRACE_smpi_computing_in(rank);
1674 #endif
1675   smpi_bench_begin();
1676   return retval;
1677 }
1678
1679 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1680                  MPI_Datatype sendtype, void *recvbuf, int recvcount,
1681                  MPI_Datatype recvtype, int root, MPI_Comm comm)
1682 {
1683   int retval;
1684
1685   smpi_bench_end();
1686 #ifdef HAVE_TRACING
1687   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1688   TRACE_smpi_computing_out(rank);
1689   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1690   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1691 #endif
1692   if (comm == MPI_COMM_NULL) {
1693     retval = MPI_ERR_COMM;
1694   } else if (sendtype == MPI_DATATYPE_NULL
1695              || recvtype == MPI_DATATYPE_NULL) {
1696     retval = MPI_ERR_TYPE;
1697   } else if (sendcounts == NULL || displs == NULL) {
1698     retval = MPI_ERR_ARG;
1699   } else {
1700     smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1701                       recvcount, recvtype, root, comm);
1702     retval = MPI_SUCCESS;
1703   }
1704 #ifdef HAVE_TRACING
1705   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1706   TRACE_smpi_computing_in(rank);
1707 #endif
1708   smpi_bench_begin();
1709   return retval;
1710 }
1711
1712 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1713                MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1714 {
1715   int retval;
1716
1717   smpi_bench_end();
1718 #ifdef HAVE_TRACING
1719   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1720   TRACE_smpi_computing_out(rank);
1721   int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1722   TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1723 #endif
1724   if (comm == MPI_COMM_NULL) {
1725     retval = MPI_ERR_COMM;
1726   } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1727     retval = MPI_ERR_ARG;
1728   } else {
1729     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1730     retval = MPI_SUCCESS;
1731   }
1732 #ifdef HAVE_TRACING
1733   TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1734   TRACE_smpi_computing_in(rank);
1735 #endif
1736   smpi_bench_begin();
1737   return retval;
1738 }
1739
1740 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1741                   MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1742 {
1743   int retval;
1744
1745   smpi_bench_end();
1746 #ifdef HAVE_TRACING
1747   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1748   TRACE_smpi_computing_out(rank);
1749   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1750 #endif
1751   if (comm == MPI_COMM_NULL) {
1752     retval = MPI_ERR_COMM;
1753   } else if (datatype == MPI_DATATYPE_NULL) {
1754     retval = MPI_ERR_TYPE;
1755   } else if (op == MPI_OP_NULL) {
1756     retval = MPI_ERR_OP;
1757   } else {
1758     smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
1759     retval = MPI_SUCCESS;
1760   }
1761 #ifdef HAVE_TRACING
1762   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1763   TRACE_smpi_computing_in(rank);
1764 #endif
1765   smpi_bench_begin();
1766   return retval;
1767 }
1768
1769 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
1770              MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1771 {
1772   int retval;
1773
1774   smpi_bench_end();
1775 #ifdef HAVE_TRACING
1776   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1777   TRACE_smpi_computing_out(rank);
1778   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1779 #endif
1780   if (comm == MPI_COMM_NULL) {
1781     retval = MPI_ERR_COMM;
1782   } else if (datatype == MPI_DATATYPE_NULL) {
1783     retval = MPI_ERR_TYPE;
1784   } else if (op == MPI_OP_NULL) {
1785     retval = MPI_ERR_OP;
1786   } else {
1787     smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1788     retval = MPI_SUCCESS;
1789   }
1790 #ifdef HAVE_TRACING
1791   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1792   TRACE_smpi_computing_in(rank);
1793 #endif
1794   smpi_bench_begin();
1795   return retval;
1796 }
1797
1798 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
1799                        MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1800 {
1801   int retval, i, size, count;
1802   int *displs;
1803   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1804
1805   smpi_bench_end();
1806 #ifdef HAVE_TRACING
1807   TRACE_smpi_computing_out(rank);
1808   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1809 #endif
1810   if (comm == MPI_COMM_NULL) {
1811     retval = MPI_ERR_COMM;
1812   } else if (datatype == MPI_DATATYPE_NULL) {
1813     retval = MPI_ERR_TYPE;
1814   } else if (op == MPI_OP_NULL) {
1815     retval = MPI_ERR_OP;
1816   } else if (recvcounts == NULL) {
1817     retval = MPI_ERR_ARG;
1818   } else {
1819     /* arbitrarily choose root as rank 0 */
1820     /* TODO: faster direct implementation ? */
1821     size = smpi_comm_size(comm);
1822     count = 0;
1823     displs = xbt_new(int, size);
1824     for (i = 0; i < size; i++) {
1825       count += recvcounts[i];
1826       displs[i] = 0;
1827     }
1828     smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1829     smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf,
1830                       recvcounts[rank], datatype, 0, comm);
1831     xbt_free(displs);
1832     retval = MPI_SUCCESS;
1833   }
1834 #ifdef HAVE_TRACING
1835   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1836   TRACE_smpi_computing_in(rank);
1837 #endif
1838   smpi_bench_begin();
1839   return retval;
1840 }
1841
1842 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1843                  void *recvbuf, int recvcount, MPI_Datatype recvtype,
1844                  MPI_Comm comm)
1845 {
1846   int retval, size, sendsize;
1847
1848   smpi_bench_end();
1849 #ifdef HAVE_TRACING
1850   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1851   TRACE_smpi_computing_out(rank);
1852   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1853 #endif
1854   if (comm == MPI_COMM_NULL) {
1855     retval = MPI_ERR_COMM;
1856   } else if (sendtype == MPI_DATATYPE_NULL
1857              || recvtype == MPI_DATATYPE_NULL) {
1858     retval = MPI_ERR_TYPE;
1859   } else {
1860     size = smpi_comm_size(comm);
1861     sendsize = smpi_datatype_size(sendtype) * sendcount;
1862     if (sendsize < 200 && size > 12) {
1863       retval =
1864           smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype,
1865                                          recvbuf, recvcount, recvtype,
1866                                          comm);
1867     } else if (sendsize < 3000) {
1868       retval =
1869           smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount,
1870                                                 sendtype, recvbuf,
1871                                                 recvcount, recvtype, comm);
1872     } else {
1873       retval =
1874           smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype,
1875                                             recvbuf, recvcount, recvtype,
1876                                             comm);
1877     }
1878   }
1879 #ifdef HAVE_TRACING
1880   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1881   TRACE_smpi_computing_in(rank);
1882 #endif
1883   smpi_bench_begin();
1884   return retval;
1885 }
1886
1887 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
1888                   MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
1889                   int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
1890 {
1891   int retval;
1892
1893   smpi_bench_end();
1894 #ifdef HAVE_TRACING
1895   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1896   TRACE_smpi_computing_out(rank);
1897   TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1898 #endif
1899   if (comm == MPI_COMM_NULL) {
1900     retval = MPI_ERR_COMM;
1901   } else if (sendtype == MPI_DATATYPE_NULL
1902              || recvtype == MPI_DATATYPE_NULL) {
1903     retval = MPI_ERR_TYPE;
1904   } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
1905              || recvdisps == NULL) {
1906     retval = MPI_ERR_ARG;
1907   } else {
1908     retval =
1909         smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype,
1910                                   recvbuf, recvcounts, recvdisps, recvtype,
1911                                   comm);
1912   }
1913 #ifdef HAVE_TRACING
1914   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1915   TRACE_smpi_computing_in(rank);
1916 #endif
1917   smpi_bench_begin();
1918   return retval;
1919 }
1920
1921
1922 int PMPI_Get_processor_name(char *name, int *resultlen)
1923 {
1924   int retval = MPI_SUCCESS;
1925
1926   smpi_bench_end();
1927   strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
1928           strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
1929           strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
1930           MPI_MAX_PROCESSOR_NAME - 1 );
1931   *resultlen =
1932       strlen(name) >
1933       MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1934
1935   smpi_bench_begin();
1936   return retval;
1937 }
1938
1939 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
1940 {
1941   int retval = MPI_SUCCESS;
1942   size_t size;
1943
1944   smpi_bench_end();
1945   if (status == NULL || count == NULL) {
1946     retval = MPI_ERR_ARG;
1947   } else if (datatype == MPI_DATATYPE_NULL) {
1948     retval = MPI_ERR_TYPE;
1949   } else {
1950     size = smpi_datatype_size(datatype);
1951     if (size == 0) {
1952       *count = 0;
1953     } else if (status->count % size != 0) {
1954       retval = MPI_UNDEFINED;
1955     } else {
1956       *count = smpi_mpi_get_count(status, datatype);
1957     }
1958   }
1959   smpi_bench_begin();
1960   return retval;
1961 }
1962
1963 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
1964   int retval;
1965
1966   smpi_bench_end();
1967   if (old_type == MPI_DATATYPE_NULL) {
1968     retval = MPI_ERR_TYPE;
1969   } else if (count<0){
1970     retval = MPI_ERR_COUNT;
1971   } else {
1972     retval = smpi_datatype_contiguous(count, old_type, new_type);
1973   }
1974   smpi_bench_begin();
1975   return retval;
1976 }
1977
1978 int PMPI_Type_commit(MPI_Datatype* datatype) {
1979   int retval;
1980
1981   smpi_bench_end();
1982   if (datatype == MPI_DATATYPE_NULL) {
1983     retval = MPI_ERR_TYPE;
1984   } else {
1985     smpi_datatype_commit(datatype);
1986     retval = MPI_SUCCESS;
1987   }
1988   smpi_bench_begin();
1989   return retval;
1990 }
1991
1992
1993 int PMPI_Type_vector(int count, int blocklen, int 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_vector(count, blocklen, stride, old_type, new_type);
2003   }
2004   smpi_bench_begin();
2005   return retval;
2006 }
2007
2008 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2009   int retval;
2010
2011   smpi_bench_end();
2012   if (old_type == MPI_DATATYPE_NULL) {
2013     retval = MPI_ERR_TYPE;
2014   } else if (count<0 || blocklen<0){
2015     retval = MPI_ERR_COUNT;
2016   } else {
2017     retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2018   }
2019   smpi_bench_begin();
2020   return retval;
2021 }
2022
2023
2024 int PMPI_Type_indexed(int count, int* blocklens, int* 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_indexed(count, blocklens, indices, old_type, new_type);
2034   }
2035   smpi_bench_begin();
2036   return retval;
2037 }
2038
2039 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2040   int retval;
2041
2042   smpi_bench_end();
2043   if (old_type == MPI_DATATYPE_NULL) {
2044     retval = MPI_ERR_TYPE;
2045   } else if (count<0){
2046     retval = MPI_ERR_COUNT;
2047   } else {
2048     retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2049   }
2050   smpi_bench_begin();
2051   return retval;
2052 }
2053
2054
2055 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2056   int retval;
2057
2058   smpi_bench_end();
2059   if (count<0){
2060     retval = MPI_ERR_COUNT;
2061   } else {
2062     retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2063   }
2064   smpi_bench_begin();
2065   return retval;}
2066
2067 int PMPI_Error_class(int errorcode, int* errorclass) {
2068   // assume smpi uses only standard mpi error codes
2069   *errorclass=errorcode;
2070   return MPI_SUCCESS;
2071 }
2072
2073 /* The following calls are not yet implemented and will fail at runtime. */
2074 /* Once implemented, please move them above this notice. */
2075
2076 static int not_yet_implemented(void) {
2077           XBT_WARN("Not yet implemented");
2078    return MPI_SUCCESS;
2079 }
2080
2081 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2082    return not_yet_implemented();
2083 }
2084
2085 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2086    return not_yet_implemented();
2087 }
2088
2089 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
2090    return not_yet_implemented();
2091 }
2092
2093 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2094    return not_yet_implemented();
2095 }
2096
2097 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2098    return not_yet_implemented();
2099 }
2100
2101 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2102    return not_yet_implemented();
2103 }
2104
2105 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2106    return not_yet_implemented();
2107 }
2108
2109 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2110    return not_yet_implemented();
2111 }
2112
2113 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2114    return not_yet_implemented();
2115 }
2116
2117 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2118    return not_yet_implemented();
2119 }
2120
2121 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2122    return not_yet_implemented();
2123 }
2124
2125 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2126    return not_yet_implemented();
2127 }
2128
2129 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2130    return not_yet_implemented();
2131 }
2132
2133 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2134    return not_yet_implemented();
2135 }
2136
2137 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2138    return not_yet_implemented();
2139 }
2140
2141 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2142    return not_yet_implemented();
2143 }
2144
2145 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2146    return not_yet_implemented();
2147 }
2148
2149 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2150    return not_yet_implemented();
2151 }
2152
2153 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2154    return not_yet_implemented();
2155 }
2156
2157 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2158    return not_yet_implemented();
2159 }
2160
2161 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2162    return not_yet_implemented();
2163 }
2164
2165
2166 int PMPI_Cancel(MPI_Request* request) {
2167    return not_yet_implemented();
2168 }
2169
2170 int PMPI_Buffer_attach(void* buffer, int size) {
2171    return not_yet_implemented();
2172 }
2173
2174 int PMPI_Buffer_detach(void* buffer, int* size) {
2175    return not_yet_implemented();
2176 }
2177
2178 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2179    return not_yet_implemented();
2180 }
2181
2182 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2183 {
2184    return not_yet_implemented();
2185 }
2186
2187 int PMPI_Pcontrol(const int level )
2188 {
2189    return not_yet_implemented();
2190 }
2191
2192 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2193    return not_yet_implemented();
2194 }
2195
2196 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2197    return not_yet_implemented();
2198 }
2199
2200 int PMPI_Ssend_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_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2205    return not_yet_implemented();
2206 }
2207
2208 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2209    return not_yet_implemented();
2210 }
2211
2212 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2213    return not_yet_implemented();
2214 }
2215
2216 int PMPI_Bsend_init(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 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2221    return not_yet_implemented();
2222 }
2223
2224 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2225    return not_yet_implemented();
2226 }
2227
2228 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2229    return not_yet_implemented();
2230 }
2231
2232 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2233    return not_yet_implemented();
2234 }
2235
2236
2237 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2238    return not_yet_implemented();
2239 }
2240
2241 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2242    return not_yet_implemented();
2243 }
2244
2245 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2246    return not_yet_implemented();
2247 }
2248
2249 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2250    return not_yet_implemented();
2251 }
2252
2253 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2254    return not_yet_implemented();
2255 }
2256
2257 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2258    return not_yet_implemented();
2259 }
2260
2261 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2262    return not_yet_implemented();
2263 }
2264
2265 int PMPI_Keyval_free(int* keyval) {
2266    return not_yet_implemented();
2267 }
2268
2269 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2270    return not_yet_implemented();
2271 }
2272
2273 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2274    return not_yet_implemented();
2275 }
2276
2277 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2278    return not_yet_implemented();
2279 }
2280
2281 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2282    return not_yet_implemented();
2283 }
2284
2285 int PMPI_Initialized(int* flag) {
2286    return not_yet_implemented();
2287 }