Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Some cleanups, some smells
[simgrid.git] / src / smpi / smpi_pmpi.cpp
1 /* Copyright (c) 2007-2017. The SimGrid Team. All rights reserved.          */
2
3 /* This program is free software; you can redistribute it and/or modify it
4  * under the terms of the license (GNU LGPL) which comes with this package. */
5
6 #include <simgrid/s4u/host.hpp>
7
8 #include "private.h"
9
10 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_pmpi, smpi, "Logging specific to SMPI (pmpi)");
11
12 //this function need to be here because of the calls to smpi_bench
13 void TRACE_smpi_set_category(const char *category)
14 {
15   //need to end bench otherwise categories for execution tasks are wrong
16   smpi_bench_end();
17   TRACE_internal_smpi_set_category (category);
18   //begin bench after changing process's category
19   smpi_bench_begin();
20 }
21
22 /* PMPI User level calls */
23 extern "C" { // Obviously, the C MPI interface should use the C linkage
24
25 int PMPI_Init(int *argc, char ***argv)
26 {
27   // PMPI_Init is call only one time by only by SMPI process
28   int already_init;
29   MPI_Initialized(&already_init);
30   if(already_init == 0){
31     smpi_process_init(argc, argv);
32     smpi_process_mark_as_initialized();
33     int rank = smpi_process_index();
34     TRACE_smpi_init(rank);
35     TRACE_smpi_computing_init(rank);
36     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
37     extra->type = TRACING_INIT;
38     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
39     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
40     smpi_bench_begin();
41   }
42
43   smpi_mpi_init();
44
45   return MPI_SUCCESS;
46 }
47
48 int PMPI_Finalize()
49 {
50   smpi_bench_end();
51   int rank = smpi_process_index();
52   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
53   extra->type = TRACING_FINALIZE;
54   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
55
56   smpi_process_finalize();
57
58   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
59   TRACE_smpi_finalize(smpi_process_index());
60   smpi_process_destroy();
61   return MPI_SUCCESS;
62 }
63
64 int PMPI_Finalized(int* flag)
65 {
66   *flag=smpi_process_finalized();
67   return MPI_SUCCESS;
68 }
69
70 int PMPI_Get_version (int *version,int *subversion){
71   *version = MPI_VERSION;
72   *subversion= MPI_SUBVERSION;
73   return MPI_SUCCESS;
74 }
75
76 int PMPI_Get_library_version (char *version,int *len){
77   smpi_bench_end();
78   snprintf(version,MPI_MAX_LIBRARY_VERSION_STRING,"SMPI Version %d.%d. Copyright The Simgrid Team 2007-2015",
79            SIMGRID_VERSION_MAJOR, SIMGRID_VERSION_MINOR);
80   *len = strlen(version) > MPI_MAX_LIBRARY_VERSION_STRING ? MPI_MAX_LIBRARY_VERSION_STRING : strlen(version);
81   smpi_bench_begin();
82   return MPI_SUCCESS;
83 }
84
85 int PMPI_Init_thread(int *argc, char ***argv, int required, int *provided)
86 {
87   if (provided != nullptr) {
88     *provided = MPI_THREAD_SINGLE;
89   }
90   return MPI_Init(argc, argv);
91 }
92
93 int PMPI_Query_thread(int *provided)
94 {
95   if (provided == nullptr) {
96     return MPI_ERR_ARG;
97   } else {
98     *provided = MPI_THREAD_SINGLE;
99     return MPI_SUCCESS;
100   }
101 }
102
103 int PMPI_Is_thread_main(int *flag)
104 {
105   if (flag == nullptr) {
106     return MPI_ERR_ARG;
107   } else {
108     *flag = smpi_process_index() == 0;
109     return MPI_SUCCESS;
110   }
111 }
112
113 int PMPI_Abort(MPI_Comm comm, int errorcode)
114 {
115   smpi_bench_end();
116   smpi_process_destroy();
117   // FIXME: should kill all processes in comm instead
118   simcall_process_kill(SIMIX_process_self());
119   return MPI_SUCCESS;
120 }
121
122 double PMPI_Wtime()
123 {
124   return smpi_mpi_wtime();
125 }
126
127 extern double sg_maxmin_precision;
128 double PMPI_Wtick()
129 {
130   return sg_maxmin_precision;
131 }
132
133 int PMPI_Address(void *location, MPI_Aint * address)
134 {
135   if (address==nullptr) {
136     return MPI_ERR_ARG;
137   } else {
138     *address = reinterpret_cast<MPI_Aint>(location);
139     return MPI_SUCCESS;
140   }
141 }
142
143 int PMPI_Get_address(void *location, MPI_Aint * address)
144 {
145   return PMPI_Address(location, address);
146 }
147
148 int PMPI_Type_free(MPI_Datatype * datatype)
149 {
150   /* Free a predefined datatype is an error according to the standard, and should be checked for */
151   if (*datatype == MPI_DATATYPE_NULL) {
152     return MPI_ERR_ARG;
153   } else {
154     Datatype::unref(*datatype);
155     return MPI_SUCCESS;
156   }
157 }
158
159 int PMPI_Type_size(MPI_Datatype datatype, int *size)
160 {
161   if (datatype == MPI_DATATYPE_NULL) {
162     return MPI_ERR_TYPE;
163   } else if (size == nullptr) {
164     return MPI_ERR_ARG;
165   } else {
166     *size = static_cast<int>(datatype->size());
167     return MPI_SUCCESS;
168   }
169 }
170
171 int PMPI_Type_size_x(MPI_Datatype datatype, MPI_Count *size)
172 {
173   if (datatype == MPI_DATATYPE_NULL) {
174     return MPI_ERR_TYPE;
175   } else if (size == nullptr) {
176     return MPI_ERR_ARG;
177   } else {
178     *size = static_cast<MPI_Count>(datatype->size());
179     return MPI_SUCCESS;
180   }
181 }
182
183 int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
184 {
185   if (datatype == MPI_DATATYPE_NULL) {
186     return MPI_ERR_TYPE;
187   } else if (lb == nullptr || extent == nullptr) {
188     return MPI_ERR_ARG;
189   } else {
190     return datatype->extent(lb, extent);
191   }
192 }
193
194 int PMPI_Type_get_true_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
195 {
196   return PMPI_Type_get_extent(datatype, lb, extent);
197 }
198
199 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
200 {
201   if (datatype == MPI_DATATYPE_NULL) {
202     return MPI_ERR_TYPE;
203   } else if (extent == nullptr) {
204     return MPI_ERR_ARG;
205   } else {
206     *extent = datatype->get_extent();
207     return MPI_SUCCESS;
208   }
209 }
210
211 int PMPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
212 {
213   if (datatype == MPI_DATATYPE_NULL) {
214     return MPI_ERR_TYPE;
215   } else if (disp == nullptr) {
216     return MPI_ERR_ARG;
217   } else {
218     *disp = datatype->lb();
219     return MPI_SUCCESS;
220   }
221 }
222
223 int PMPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
224 {
225   if (datatype == MPI_DATATYPE_NULL) {
226     return MPI_ERR_TYPE;
227   } else if (disp == nullptr) {
228     return MPI_ERR_ARG;
229   } else {
230     *disp = datatype->ub();
231     return MPI_SUCCESS;
232   }
233 }
234
235 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
236   int retval = MPI_SUCCESS;
237   if (datatype == MPI_DATATYPE_NULL) {
238     retval=MPI_ERR_TYPE;
239   } else {
240     *newtype = new Datatype(datatype, &retval);
241     //error when duplicating, free the new datatype
242     if(retval!=MPI_SUCCESS){
243       Datatype::unref(*newtype);
244       *newtype = MPI_DATATYPE_NULL;
245     }
246   }
247   return retval;
248 }
249
250 int PMPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
251 {
252   if (function == nullptr || op == nullptr) {
253     return MPI_ERR_ARG;
254   } else {
255     *op = new Op(function, (commute!=0));
256     return MPI_SUCCESS;
257   }
258 }
259
260 int PMPI_Op_free(MPI_Op * op)
261 {
262   if (op == nullptr) {
263     return MPI_ERR_ARG;
264   } else if (*op == MPI_OP_NULL) {
265     return MPI_ERR_OP;
266   } else {
267     delete (*op);
268     *op = MPI_OP_NULL;
269     return MPI_SUCCESS;
270   }
271 }
272
273 int PMPI_Group_free(MPI_Group * group)
274 {
275   if (group == nullptr) {
276     return MPI_ERR_ARG;
277   } else {
278     if(*group != MPI_COMM_WORLD->group() && *group != MPI_GROUP_EMPTY)
279       Group::unref(*group);
280     *group = MPI_GROUP_NULL;
281     return MPI_SUCCESS;
282   }
283 }
284
285 int PMPI_Group_size(MPI_Group group, int *size)
286 {
287   if (group == MPI_GROUP_NULL) {
288     return MPI_ERR_GROUP;
289   } else if (size == nullptr) {
290     return MPI_ERR_ARG;
291   } else {
292     *size = group->size();
293     return MPI_SUCCESS;
294   }
295 }
296
297 int PMPI_Group_rank(MPI_Group group, int *rank)
298 {
299   if (group == MPI_GROUP_NULL) {
300     return MPI_ERR_GROUP;
301   } else if (rank == nullptr) {
302     return MPI_ERR_ARG;
303   } else {
304     *rank = group->rank(smpi_process_index());
305     return MPI_SUCCESS;
306   }
307 }
308
309 int PMPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1, MPI_Group group2, int *ranks2)
310 {
311   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
312     return MPI_ERR_GROUP;
313   } else {
314     for (int i = 0; i < n; i++) {
315       if(ranks1[i]==MPI_PROC_NULL){
316         ranks2[i]=MPI_PROC_NULL;
317       }else{
318         int index = group1->index(ranks1[i]);
319         ranks2[i] = group2->rank(index);
320       }
321     }
322     return MPI_SUCCESS;
323   }
324 }
325
326 int PMPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
327 {
328   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
329     return MPI_ERR_GROUP;
330   } else if (result == nullptr) {
331     return MPI_ERR_ARG;
332   } else {
333     *result = group1->compare(group2);
334     return MPI_SUCCESS;
335   }
336 }
337
338 int PMPI_Group_union(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
339 {
340
341   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
342     return MPI_ERR_GROUP;
343   } else if (newgroup == nullptr) {
344     return MPI_ERR_ARG;
345   } else {
346     return group1->group_union(group2, newgroup);
347   }
348 }
349
350 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
351 {
352
353   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
354     return MPI_ERR_GROUP;
355   } else if (newgroup == nullptr) {
356     return MPI_ERR_ARG;
357   } else {
358     return group1->intersection(group2,newgroup);
359   }
360 }
361
362 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
363 {
364   if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
365     return MPI_ERR_GROUP;
366   } else if (newgroup == nullptr) {
367     return MPI_ERR_ARG;
368   } else {
369     return group1->difference(group2,newgroup);
370   }
371 }
372
373 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
374 {
375   if (group == MPI_GROUP_NULL) {
376     return MPI_ERR_GROUP;
377   } else if (newgroup == nullptr) {
378     return MPI_ERR_ARG;
379   } else {
380     return group->incl(n, ranks, newgroup);
381   }
382 }
383
384 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
385 {
386   if (group == MPI_GROUP_NULL) {
387     return MPI_ERR_GROUP;
388   } else if (newgroup == nullptr) {
389     return MPI_ERR_ARG;
390   } else {
391     if (n == 0) {
392       *newgroup = group;
393       if (group != MPI_COMM_WORLD->group()
394                 && group != MPI_COMM_SELF->group() && group != MPI_GROUP_EMPTY)
395       group->ref();
396       return MPI_SUCCESS;
397     } else if (n == group->size()) {
398       *newgroup = MPI_GROUP_EMPTY;
399       return MPI_SUCCESS;
400     } else {
401       return group->excl(n,ranks,newgroup);
402     }
403   }
404 }
405
406 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3], MPI_Group * newgroup)
407 {
408   if (group == MPI_GROUP_NULL) {
409     return MPI_ERR_GROUP;
410   } else if (newgroup == nullptr) {
411     return MPI_ERR_ARG;
412   } else {
413     if (n == 0) {
414       *newgroup = MPI_GROUP_EMPTY;
415       return MPI_SUCCESS;
416     } else {
417       return group->range_incl(n,ranges,newgroup);
418     }
419   }
420 }
421
422 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3], MPI_Group * newgroup)
423 {
424   if (group == MPI_GROUP_NULL) {
425     return MPI_ERR_GROUP;
426   } else if (newgroup == nullptr) {
427     return MPI_ERR_ARG;
428   } else {
429     if (n == 0) {
430       *newgroup = group;
431       if (group != MPI_COMM_WORLD->group() && group != MPI_COMM_SELF->group() &&
432           group != MPI_GROUP_EMPTY)
433         group->ref();
434       return MPI_SUCCESS;
435     } else {
436       return group->range_excl(n,ranges,newgroup);
437     }
438   }
439 }
440
441 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
442 {
443   if (comm == MPI_COMM_NULL) {
444     return MPI_ERR_COMM;
445   } else if (rank == nullptr) {
446     return MPI_ERR_ARG;
447   } else {
448     *rank = comm->rank();
449     return MPI_SUCCESS;
450   }
451 }
452
453 int PMPI_Comm_size(MPI_Comm comm, int *size)
454 {
455   if (comm == MPI_COMM_NULL) {
456     return MPI_ERR_COMM;
457   } else if (size == nullptr) {
458     return MPI_ERR_ARG;
459   } else {
460     *size = comm->size();
461     return MPI_SUCCESS;
462   }
463 }
464
465 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
466 {
467   if (comm == MPI_COMM_NULL)  {
468     return MPI_ERR_COMM;
469   } else if (name == nullptr || len == nullptr)  {
470     return MPI_ERR_ARG;
471   } else {
472     comm->get_name(name, len);
473     return MPI_SUCCESS;
474   }
475 }
476
477 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
478 {
479   if (comm == MPI_COMM_NULL) {
480     return MPI_ERR_COMM;
481   } else if (group == nullptr) {
482     return MPI_ERR_ARG;
483   } else {
484     *group = comm->group();
485     if (*group != MPI_COMM_WORLD->group() && *group != MPI_GROUP_NULL && *group != MPI_GROUP_EMPTY)
486       (*group)->ref();
487     return MPI_SUCCESS;
488   }
489 }
490
491 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
492 {
493   if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
494     return MPI_ERR_COMM;
495   } else if (result == nullptr) {
496     return MPI_ERR_ARG;
497   } else {
498     if (comm1 == comm2) {       /* Same communicators means same groups */
499       *result = MPI_IDENT;
500     } else {
501       *result = comm1->group()->compare(comm2->group());
502       if (*result == MPI_IDENT) {
503         *result = MPI_CONGRUENT;
504       }
505     }
506     return MPI_SUCCESS;
507   }
508 }
509
510 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
511 {
512   if (comm == MPI_COMM_NULL) {
513     return MPI_ERR_COMM;
514   } else if (newcomm == nullptr) {
515     return MPI_ERR_ARG;
516   } else {
517     return comm->dup(newcomm);
518   }
519 }
520
521 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
522 {
523   if (comm == MPI_COMM_NULL) {
524     return MPI_ERR_COMM;
525   } else if (group == MPI_GROUP_NULL) {
526     return MPI_ERR_GROUP;
527   } else if (newcomm == nullptr) {
528     return MPI_ERR_ARG;
529   } else if(group->rank(smpi_process_index())==MPI_UNDEFINED){
530     *newcomm= MPI_COMM_NULL;
531     return MPI_SUCCESS;
532   }else{
533     group->ref();
534     *newcomm = new Comm(group, nullptr);
535     return MPI_SUCCESS;
536   }
537 }
538
539 int PMPI_Comm_free(MPI_Comm * comm)
540 {
541   if (comm == nullptr) {
542     return MPI_ERR_ARG;
543   } else if (*comm == MPI_COMM_NULL) {
544     return MPI_ERR_COMM;
545   } else {
546     Comm::destroy(*comm);
547     *comm = MPI_COMM_NULL;
548     return MPI_SUCCESS;
549   }
550 }
551
552 int PMPI_Comm_disconnect(MPI_Comm * comm)
553 {
554   /* TODO: wait until all communication in comm are done */
555   if (comm == nullptr) {
556     return MPI_ERR_ARG;
557   } else if (*comm == MPI_COMM_NULL) {
558     return MPI_ERR_COMM;
559   } else {
560     Comm::destroy(*comm);
561     *comm = MPI_COMM_NULL;
562     return MPI_SUCCESS;
563   }
564 }
565
566 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
567 {
568   int retval = 0;
569   smpi_bench_end();
570
571   if (comm_out == nullptr) {
572     retval = MPI_ERR_ARG;
573   } else if (comm == MPI_COMM_NULL) {
574     retval = MPI_ERR_COMM;
575   } else {
576     *comm_out = comm->split(color, key);
577     retval = MPI_SUCCESS;
578   }
579   smpi_bench_begin();
580
581   return retval;
582 }
583
584 int PMPI_Comm_create_group(MPI_Comm comm, MPI_Group group, int, MPI_Comm* comm_out)
585 {
586   int retval = 0;
587   smpi_bench_end();
588
589   if (comm_out == nullptr) {
590     retval = MPI_ERR_ARG;
591   } else if (comm == MPI_COMM_NULL) {
592     retval = MPI_ERR_COMM;
593   } else {
594     retval = MPI_Comm_create(comm, group, comm_out);
595   }
596   smpi_bench_begin();
597
598   return retval;
599 }
600
601 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
602 {
603   int retval = 0;
604
605   smpi_bench_end();
606   if (request == nullptr) {
607       retval = MPI_ERR_ARG;
608   } else if (comm == MPI_COMM_NULL) {
609       retval = MPI_ERR_COMM;
610   } else if (!datatype->is_valid()) {
611       retval = MPI_ERR_TYPE;
612   } else if (dst == MPI_PROC_NULL) {
613       retval = MPI_SUCCESS;
614   } else {
615       *request = Request::send_init(buf, count, datatype, dst, tag, comm);
616       retval = MPI_SUCCESS;
617   }
618   smpi_bench_begin();
619   if (retval != MPI_SUCCESS && request != nullptr)
620     *request = MPI_REQUEST_NULL;
621   return retval;
622 }
623
624 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
625 {
626   int retval = 0;
627
628   smpi_bench_end();
629   if (request == nullptr) {
630     retval = MPI_ERR_ARG;
631   } else if (comm == MPI_COMM_NULL) {
632     retval = MPI_ERR_COMM;
633   } else if (!datatype->is_valid()) {
634       retval = MPI_ERR_TYPE;
635   } else if (src == MPI_PROC_NULL) {
636     retval = MPI_SUCCESS;
637   } else {
638     *request = Request::recv_init(buf, count, datatype, src, tag, comm);
639     retval = MPI_SUCCESS;
640   }
641   smpi_bench_begin();
642   if (retval != MPI_SUCCESS && request != nullptr)
643     *request = MPI_REQUEST_NULL;
644   return retval;
645 }
646
647 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
648 {
649   int retval = 0;
650
651   smpi_bench_end();
652   if (request == nullptr) {
653     retval = MPI_ERR_ARG;
654   } else if (comm == MPI_COMM_NULL) {
655     retval = MPI_ERR_COMM;
656   } else if (!datatype->is_valid()) {
657       retval = MPI_ERR_TYPE;
658   } else if (dst == MPI_PROC_NULL) {
659     retval = MPI_SUCCESS;
660   } else {
661     *request = Request::ssend_init(buf, count, datatype, dst, tag, comm);
662     retval = MPI_SUCCESS;
663   }
664   smpi_bench_begin();
665   if (retval != MPI_SUCCESS && request != nullptr)
666     *request = MPI_REQUEST_NULL;
667   return retval;
668 }
669
670 int PMPI_Start(MPI_Request * request)
671 {
672   int retval = 0;
673
674   smpi_bench_end();
675   if (request == nullptr || *request == MPI_REQUEST_NULL) {
676     retval = MPI_ERR_REQUEST;
677   } else {
678     (*request)->start();
679     retval = MPI_SUCCESS;
680   }
681   smpi_bench_begin();
682   return retval;
683 }
684
685 int PMPI_Startall(int count, MPI_Request * requests)
686 {
687   int retval;
688   smpi_bench_end();
689   if (requests == nullptr) {
690     retval = MPI_ERR_ARG;
691   } else {
692     retval = MPI_SUCCESS;
693     for (int i = 0; i < count; i++) {
694       if(requests[i] == MPI_REQUEST_NULL) {
695         retval = MPI_ERR_REQUEST;
696       }
697     }
698     if(retval != MPI_ERR_REQUEST) {
699       Request::startall(count, requests);
700     }
701   }
702   smpi_bench_begin();
703   return retval;
704 }
705
706 int PMPI_Request_free(MPI_Request * request)
707 {
708   int retval = 0;
709
710   smpi_bench_end();
711   if (*request == MPI_REQUEST_NULL) {
712     retval = MPI_ERR_ARG;
713   } else {
714     Request::unref(request);
715     retval = MPI_SUCCESS;
716   }
717   smpi_bench_begin();
718   return retval;
719 }
720
721 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
722 {
723   int retval = 0;
724
725   smpi_bench_end();
726
727   if (request == nullptr) {
728     retval = MPI_ERR_ARG;
729   } else if (comm == MPI_COMM_NULL) {
730     retval = MPI_ERR_COMM;
731   } else if (src == MPI_PROC_NULL) {
732     *request = MPI_REQUEST_NULL;
733     retval = MPI_SUCCESS;
734   } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
735     retval = MPI_ERR_RANK;
736   } else if ((count < 0) || (buf==nullptr && count > 0)) {
737     retval = MPI_ERR_COUNT;
738   } else if (!datatype->is_valid()) {
739       retval = MPI_ERR_TYPE;
740   } else if(tag<0 && tag !=  MPI_ANY_TAG){
741     retval = MPI_ERR_TAG;
742   } else {
743
744     int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
745     int src_traced = comm->group()->index(src);
746
747     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
748     extra->type = TRACING_IRECV;
749     extra->src = src_traced;
750     extra->dst = rank;
751     int known=0;
752     extra->datatype1 = encode_datatype(datatype, &known);
753     int dt_size_send = 1;
754     if(known==0)
755       dt_size_send = datatype->size();
756     extra->send_size = count*dt_size_send;
757     TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
758
759     *request = Request::irecv(buf, count, datatype, src, tag, comm);
760     retval = MPI_SUCCESS;
761
762     TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
763   }
764
765   smpi_bench_begin();
766   if (retval != MPI_SUCCESS && request != nullptr)
767     *request = MPI_REQUEST_NULL;
768   return retval;
769 }
770
771
772 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
773 {
774   int retval = 0;
775
776   smpi_bench_end();
777   if (request == nullptr) {
778     retval = MPI_ERR_ARG;
779   } else if (comm == MPI_COMM_NULL) {
780     retval = MPI_ERR_COMM;
781   } else if (dst == MPI_PROC_NULL) {
782     *request = MPI_REQUEST_NULL;
783     retval = MPI_SUCCESS;
784   } else if (dst >= comm->group()->size() || dst <0){
785     retval = MPI_ERR_RANK;
786   } else if ((count < 0) || (buf==nullptr && count > 0)) {
787     retval = MPI_ERR_COUNT;
788   } else if (!datatype->is_valid()) {
789       retval = MPI_ERR_TYPE;
790   } else if(tag<0 && tag !=  MPI_ANY_TAG){
791     retval = MPI_ERR_TAG;
792   } else {
793     int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
794     int dst_traced = comm->group()->index(dst);
795     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
796     extra->type = TRACING_ISEND;
797     extra->src = rank;
798     extra->dst = dst_traced;
799     int known=0;
800     extra->datatype1 = encode_datatype(datatype, &known);
801     int dt_size_send = 1;
802     if(known==0)
803       dt_size_send = datatype->size();
804     extra->send_size = count*dt_size_send;
805     TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
806     TRACE_smpi_send(rank, rank, dst_traced, tag, count*datatype->size());
807
808     *request = Request::isend(buf, count, datatype, dst, tag, comm);
809     retval = MPI_SUCCESS;
810
811     TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
812   }
813
814   smpi_bench_begin();
815   if (retval != MPI_SUCCESS && request!=nullptr)
816     *request = MPI_REQUEST_NULL;
817   return retval;
818 }
819
820 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
821 {
822   int retval = 0;
823
824   smpi_bench_end();
825   if (request == nullptr) {
826     retval = MPI_ERR_ARG;
827   } else if (comm == MPI_COMM_NULL) {
828     retval = MPI_ERR_COMM;
829   } else if (dst == MPI_PROC_NULL) {
830     *request = MPI_REQUEST_NULL;
831     retval = MPI_SUCCESS;
832   } else if (dst >= comm->group()->size() || dst <0){
833     retval = MPI_ERR_RANK;
834   } else if ((count < 0)|| (buf==nullptr && count > 0)) {
835     retval = MPI_ERR_COUNT;
836   } else if (!datatype->is_valid()) {
837       retval = MPI_ERR_TYPE;
838   } else if(tag<0 && tag !=  MPI_ANY_TAG){
839     retval = MPI_ERR_TAG;
840   } else {
841     int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
842     int dst_traced = comm->group()->index(dst);
843     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
844     extra->type = TRACING_ISSEND;
845     extra->src = rank;
846     extra->dst = dst_traced;
847     int known=0;
848     extra->datatype1 = encode_datatype(datatype, &known);
849     int dt_size_send = 1;
850     if(known==0)
851       dt_size_send = datatype->size();
852     extra->send_size = count*dt_size_send;
853     TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
854     TRACE_smpi_send(rank, rank, dst_traced, tag, count*datatype->size());
855
856     *request = Request::issend(buf, count, datatype, dst, tag, comm);
857     retval = MPI_SUCCESS;
858
859     TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
860   }
861
862   smpi_bench_begin();
863   if (retval != MPI_SUCCESS && request!=nullptr)
864     *request = MPI_REQUEST_NULL;
865   return retval;
866 }
867
868 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
869 {
870   int retval = 0;
871
872   smpi_bench_end();
873   if (comm == MPI_COMM_NULL) {
874     retval = MPI_ERR_COMM;
875   } else if (src == MPI_PROC_NULL) {
876     smpi_empty_status(status);
877     status->MPI_SOURCE = MPI_PROC_NULL;
878     retval = MPI_SUCCESS;
879   } else if (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0)){
880     retval = MPI_ERR_RANK;
881   } else if ((count < 0) || (buf==nullptr && count > 0)) {
882     retval = MPI_ERR_COUNT;
883   } else if (!datatype->is_valid()) {
884       retval = MPI_ERR_TYPE;
885   } else if(tag<0 && tag !=  MPI_ANY_TAG){
886     retval = MPI_ERR_TAG;
887   } else {
888     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
889     int src_traced         = comm->group()->index(src);
890     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
891     extra->type            = TRACING_RECV;
892     extra->src             = src_traced;
893     extra->dst             = rank;
894     int known              = 0;
895     extra->datatype1       = encode_datatype(datatype, &known);
896     int dt_size_send       = 1;
897     if (known == 0)
898       dt_size_send   = datatype->size();
899     extra->send_size = count * dt_size_send;
900     TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
901
902     Request::recv(buf, count, datatype, src, tag, comm, status);
903     retval = MPI_SUCCESS;
904
905     // the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
906     if (status != MPI_STATUS_IGNORE) {
907       src_traced = comm->group()->index(status->MPI_SOURCE);
908       if (!TRACE_smpi_view_internals()) {
909         TRACE_smpi_recv(rank, src_traced, rank, tag);
910       }
911     }
912     TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
913   }
914
915   smpi_bench_begin();
916   return retval;
917 }
918
919 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
920 {
921   int retval = 0;
922
923   smpi_bench_end();
924
925   if (comm == MPI_COMM_NULL) {
926     retval = MPI_ERR_COMM;
927   } else if (dst == MPI_PROC_NULL) {
928     retval = MPI_SUCCESS;
929   } else if (dst >= comm->group()->size() || dst <0){
930     retval = MPI_ERR_RANK;
931   } else if ((count < 0) || (buf == nullptr && count > 0)) {
932     retval = MPI_ERR_COUNT;
933   } else if (!datatype->is_valid()) {
934     retval = MPI_ERR_TYPE;
935   } else if(tag < 0 && tag !=  MPI_ANY_TAG){
936     retval = MPI_ERR_TAG;
937   } else {
938     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
939     int dst_traced         = comm->group()->index(dst);
940     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
941     extra->type            = TRACING_SEND;
942     extra->src             = rank;
943     extra->dst             = dst_traced;
944     int known              = 0;
945     extra->datatype1       = encode_datatype(datatype, &known);
946     int dt_size_send       = 1;
947     if (known == 0) {
948       dt_size_send = datatype->size();
949     }
950     extra->send_size = count*dt_size_send;
951     TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
952     if (!TRACE_smpi_view_internals()) {
953       TRACE_smpi_send(rank, rank, dst_traced, tag,count*datatype->size());
954     }
955
956     Request::send(buf, count, datatype, dst, tag, comm);
957     retval = MPI_SUCCESS;
958
959     TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
960   }
961
962   smpi_bench_begin();
963   return retval;
964 }
965
966 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
967   int retval = 0;
968
969   smpi_bench_end();
970
971   if (comm == MPI_COMM_NULL) {
972     retval = MPI_ERR_COMM;
973   } else if (dst == MPI_PROC_NULL) {
974     retval = MPI_SUCCESS;
975   } else if (dst >= comm->group()->size() || dst <0){
976     retval = MPI_ERR_RANK;
977   } else if ((count < 0) || (buf==nullptr && count > 0)) {
978     retval = MPI_ERR_COUNT;
979   } else if (!datatype->is_valid()){
980     retval = MPI_ERR_TYPE;
981   } else if(tag<0 && tag !=  MPI_ANY_TAG){
982     retval = MPI_ERR_TAG;
983   } else {
984     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
985     int dst_traced         = comm->group()->index(dst);
986     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
987     extra->type            = TRACING_SSEND;
988     extra->src             = rank;
989     extra->dst             = dst_traced;
990     int known              = 0;
991     extra->datatype1       = encode_datatype(datatype, &known);
992     int dt_size_send       = 1;
993     if(known == 0) {
994       dt_size_send = datatype->size();
995     }
996     extra->send_size = count*dt_size_send;
997     TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
998     TRACE_smpi_send(rank, rank, dst_traced, tag,count*datatype->size());
999   
1000     Request::ssend(buf, count, datatype, dst, tag, comm);
1001     retval = MPI_SUCCESS;
1002   
1003     TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1004   }
1005
1006   smpi_bench_begin();
1007   return retval;
1008 }
1009
1010 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void *recvbuf,
1011                  int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status * status)
1012 {
1013   int retval = 0;
1014
1015   smpi_bench_end();
1016
1017   if (comm == MPI_COMM_NULL) {
1018     retval = MPI_ERR_COMM;
1019   } else if (!sendtype->is_valid() || !recvtype->is_valid()) {
1020     retval = MPI_ERR_TYPE;
1021   } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1022     smpi_empty_status(status);
1023     status->MPI_SOURCE = MPI_PROC_NULL;
1024     retval             = MPI_SUCCESS;
1025   }else if (dst >= comm->group()->size() || dst <0 ||
1026       (src!=MPI_ANY_SOURCE && (src >= comm->group()->size() || src <0))){
1027     retval = MPI_ERR_RANK;
1028   } else if ((sendcount < 0 || recvcount<0) || 
1029       (sendbuf==nullptr && sendcount > 0) || (recvbuf==nullptr && recvcount>0)) {
1030     retval = MPI_ERR_COUNT;
1031   } else if((sendtag<0 && sendtag !=  MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1032     retval = MPI_ERR_TAG;
1033   } else {
1034
1035   int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1036   int dst_traced = comm->group()->index(dst);
1037   int src_traced = comm->group()->index(src);
1038   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1039   extra->type = TRACING_SENDRECV;
1040   extra->src = src_traced;
1041   extra->dst = dst_traced;
1042   int known=0;
1043   extra->datatype1 = encode_datatype(sendtype, &known);
1044   int dt_size_send = 1;
1045   if(known==0)
1046     dt_size_send = sendtype->size();
1047   extra->send_size = sendcount*dt_size_send;
1048   extra->datatype2 = encode_datatype(recvtype, &known);
1049   int dt_size_recv = 1;
1050   if(known==0)
1051     dt_size_recv = recvtype->size();
1052   extra->recv_size = recvcount*dt_size_recv;
1053
1054   TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1055   TRACE_smpi_send(rank, rank, dst_traced, sendtag,sendcount*sendtype->size());
1056
1057   Request::sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src, recvtag, comm,
1058                     status);
1059   retval = MPI_SUCCESS;
1060
1061   TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1062   TRACE_smpi_recv(rank, src_traced, rank, recvtag);
1063   }
1064
1065   smpi_bench_begin();
1066   return retval;
1067 }
1068
1069 int PMPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
1070                           MPI_Comm comm, MPI_Status* status)
1071 {
1072   int retval = 0;
1073   if (!datatype->is_valid()) {
1074     return MPI_ERR_TYPE;
1075   } else if (count < 0) {
1076     return MPI_ERR_COUNT;
1077   } else {
1078     int size = datatype->get_extent() * count;
1079     void* recvbuf = xbt_new0(char, size);
1080     retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
1081     if(retval==MPI_SUCCESS){
1082         Datatype::copy(recvbuf, count, datatype, buf, count, datatype);
1083     }
1084     xbt_free(recvbuf);
1085
1086   }
1087   return retval;
1088 }
1089
1090 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1091 {
1092   int retval = 0;
1093   smpi_bench_end();
1094   if (request == nullptr || flag == nullptr) {
1095     retval = MPI_ERR_ARG;
1096   } else if (*request == MPI_REQUEST_NULL) {
1097     *flag= true;
1098     smpi_empty_status(status);
1099     retval = MPI_SUCCESS;
1100   } else {
1101     int rank = ((*request)->comm() != MPI_COMM_NULL) ? smpi_process_index() : -1;
1102
1103     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1104     extra->type = TRACING_TEST;
1105     TRACE_smpi_testing_in(rank, extra);
1106
1107     *flag = Request::test(request,status);
1108
1109     TRACE_smpi_testing_out(rank);
1110     retval = MPI_SUCCESS;
1111   }
1112   smpi_bench_begin();
1113   return retval;
1114 }
1115
1116 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
1117 {
1118   int retval = 0;
1119
1120   smpi_bench_end();
1121   if (index == nullptr || flag == nullptr) {
1122     retval = MPI_ERR_ARG;
1123   } else {
1124     *flag = Request::testany(count, requests, index, status);
1125     retval = MPI_SUCCESS;
1126   }
1127   smpi_bench_begin();
1128   return retval;
1129 }
1130
1131 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1132 {
1133   int retval = 0;
1134
1135   smpi_bench_end();
1136   if (flag == nullptr) {
1137     retval = MPI_ERR_ARG;
1138   } else {
1139     *flag = Request::testall(count, requests, statuses);
1140     retval = MPI_SUCCESS;
1141   }
1142   smpi_bench_begin();
1143   return retval;
1144 }
1145
1146 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1147   int retval = 0;
1148   smpi_bench_end();
1149
1150   if (status == nullptr) {
1151     retval = MPI_ERR_ARG;
1152   } else if (comm == MPI_COMM_NULL) {
1153     retval = MPI_ERR_COMM;
1154   } else if (source == MPI_PROC_NULL) {
1155     smpi_empty_status(status);
1156     status->MPI_SOURCE = MPI_PROC_NULL;
1157     retval = MPI_SUCCESS;
1158   } else {
1159     Request::probe(source, tag, comm, status);
1160     retval = MPI_SUCCESS;
1161   }
1162   smpi_bench_begin();
1163   return retval;
1164 }
1165
1166 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1167   int retval = 0;
1168   smpi_bench_end();
1169
1170   if ((flag == nullptr) || (status == nullptr)) {
1171     retval = MPI_ERR_ARG;
1172   } else if (comm == MPI_COMM_NULL) {
1173     retval = MPI_ERR_COMM;
1174   } else if (source == MPI_PROC_NULL) {
1175     *flag=true;
1176     smpi_empty_status(status);
1177     status->MPI_SOURCE = MPI_PROC_NULL;
1178     retval = MPI_SUCCESS;
1179   } else {
1180     Request::iprobe(source, tag, comm, flag, status);
1181     retval = MPI_SUCCESS;
1182   }
1183   smpi_bench_begin();
1184   return retval;
1185 }
1186
1187 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1188 {
1189   int retval = 0;
1190
1191   smpi_bench_end();
1192
1193   smpi_empty_status(status);
1194
1195   if (request == nullptr) {
1196     retval = MPI_ERR_ARG;
1197   } else if (*request == MPI_REQUEST_NULL) {
1198     retval = MPI_SUCCESS;
1199   } else {
1200
1201     int rank = (request!=nullptr && (*request)->comm() != MPI_COMM_NULL) ? smpi_process_index() : -1;
1202
1203     int src_traced = (*request)->src();
1204     int dst_traced = (*request)->dst();
1205     int tag_traced= (*request)->tag();
1206     MPI_Comm comm = (*request)->comm();
1207     int is_wait_for_receive = ((*request)->flags() & RECV);
1208     instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1209     extra->type = TRACING_WAIT;
1210     TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1211
1212     Request::wait(request, status);
1213     retval = MPI_SUCCESS;
1214
1215     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1216     TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1217     if (is_wait_for_receive) {
1218       if(src_traced==MPI_ANY_SOURCE)
1219         src_traced = (status!=MPI_STATUS_IGNORE) ?
1220           comm->group()->rank(status->MPI_SOURCE) :
1221           src_traced;
1222       TRACE_smpi_recv(rank, src_traced, dst_traced, tag_traced);
1223     }
1224   }
1225
1226   smpi_bench_begin();
1227   return retval;
1228 }
1229
1230 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1231 {
1232   if (index == nullptr)
1233     return MPI_ERR_ARG;
1234
1235   smpi_bench_end();
1236   //save requests information for tracing
1237   typedef struct {
1238     int src;
1239     int dst;
1240     int recv;
1241     int tag;
1242     MPI_Comm comm;
1243   } savedvalstype;
1244   savedvalstype* savedvals=nullptr;
1245   if(count>0){
1246     savedvals = xbt_new0(savedvalstype, count);
1247   }
1248   for (int i = 0; i < count; i++) {
1249     MPI_Request req = requests[i];      //already received requests are no longer valid
1250     if (req) {
1251       savedvals[i]=(savedvalstype){req->src(), req->dst(), (req->flags() & RECV), req->tag(), req->comm()};
1252     }
1253   }
1254   int rank_traced = smpi_process_index();
1255   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1256   extra->type = TRACING_WAITANY;
1257   extra->send_size=count;
1258   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1259
1260   *index = Request::waitany(count, requests, status);
1261
1262   if(*index!=MPI_UNDEFINED){
1263     int src_traced = savedvals[*index].src;
1264     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1265     int dst_traced = savedvals[*index].dst;
1266     int is_wait_for_receive = savedvals[*index].recv;
1267     if (is_wait_for_receive) {
1268       if(savedvals[*index].src==MPI_ANY_SOURCE)
1269         src_traced = (status != MPI_STATUSES_IGNORE)
1270                          ? savedvals[*index].comm->group()->rank(status->MPI_SOURCE)
1271                          : savedvals[*index].src;
1272       TRACE_smpi_recv(rank_traced, src_traced, dst_traced, savedvals[*index].tag);
1273     }
1274     TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1275   }
1276   xbt_free(savedvals);
1277
1278   smpi_bench_begin();
1279   return MPI_SUCCESS;
1280 }
1281
1282 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1283 {
1284   smpi_bench_end();
1285   //save information from requests
1286   typedef struct {
1287     int src;
1288     int dst;
1289     int recv;
1290     int tag;
1291     int valid;
1292     MPI_Comm comm;
1293   } savedvalstype;
1294   savedvalstype* savedvals=xbt_new0(savedvalstype, count);
1295
1296   for (int i = 0; i < count; i++) {
1297     MPI_Request req = requests[i];
1298     if(req!=MPI_REQUEST_NULL){
1299       savedvals[i]=(savedvalstype){req->src(), req->dst(), (req->flags() & RECV), req->tag(), 1, req->comm()};
1300     }else{
1301       savedvals[i].valid=0;
1302     }
1303   }
1304   int rank_traced = smpi_process_index();
1305   instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1306   extra->type = TRACING_WAITALL;
1307   extra->send_size=count;
1308   TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1309
1310   int retval =Request::waitall(count, requests, status);
1311
1312   for (int i = 0; i < count; i++) {
1313     if(savedvals[i].valid){
1314     //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1315       int src_traced = savedvals[i].src;
1316       int dst_traced = savedvals[i].dst;
1317       int is_wait_for_receive = savedvals[i].recv;
1318       if (is_wait_for_receive) {
1319         if(src_traced==MPI_ANY_SOURCE)
1320         src_traced = (status!=MPI_STATUSES_IGNORE) ?
1321                           savedvals[i].comm->group()->rank(status[i].MPI_SOURCE) : savedvals[i].src;
1322         TRACE_smpi_recv(rank_traced, src_traced, dst_traced,savedvals[i].tag);
1323       }
1324     }
1325   }
1326   TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1327   xbt_free(savedvals);
1328
1329   smpi_bench_begin();
1330   return retval;
1331 }
1332
1333 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
1334 {
1335   int retval = 0;
1336
1337   smpi_bench_end();
1338   if (outcount == nullptr) {
1339     retval = MPI_ERR_ARG;
1340   } else {
1341     *outcount = Request::waitsome(incount, requests, indices, status);
1342     retval = MPI_SUCCESS;
1343   }
1344   smpi_bench_begin();
1345   return retval;
1346 }
1347
1348 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
1349 {
1350   int retval = 0;
1351
1352   smpi_bench_end();
1353   if (outcount == nullptr) {
1354     retval = MPI_ERR_ARG;
1355   } else {
1356     *outcount = Request::testsome(incount, requests, indices, status);
1357     retval    = MPI_SUCCESS;
1358   }
1359   smpi_bench_begin();
1360   return retval;
1361 }
1362
1363
1364 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1365 {
1366   int retval = 0;
1367
1368   smpi_bench_end();
1369
1370   if (comm == MPI_COMM_NULL) {
1371     retval = MPI_ERR_COMM;
1372   } else if (!datatype->is_valid()) {
1373     retval = MPI_ERR_ARG;
1374   } else {
1375     int rank        = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1376     int root_traced = comm->group()->index(root);
1377
1378     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1379     extra->type            = TRACING_BCAST;
1380     extra->root            = root_traced;
1381     int known              = 0;
1382     extra->datatype1       = encode_datatype(datatype, &known);
1383     int dt_size_send       = 1;
1384     if (known == 0)
1385       dt_size_send   = datatype->size();
1386     extra->send_size = count * dt_size_send;
1387     TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1388     if (comm->size() > 1)
1389       Colls::bcast(buf, count, datatype, root, comm);
1390     retval = MPI_SUCCESS;
1391
1392     TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1393   }
1394   smpi_bench_begin();
1395   return retval;
1396 }
1397
1398 int PMPI_Barrier(MPI_Comm comm)
1399 {
1400   int retval = 0;
1401
1402   smpi_bench_end();
1403
1404   if (comm == MPI_COMM_NULL) {
1405     retval = MPI_ERR_COMM;
1406   } else {
1407     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1408     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1409     extra->type            = TRACING_BARRIER;
1410     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1411
1412     Colls::barrier(comm);
1413     retval = MPI_SUCCESS;
1414
1415     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1416   }
1417
1418   smpi_bench_begin();
1419   return retval;
1420 }
1421
1422 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,void *recvbuf, int recvcount, MPI_Datatype recvtype,
1423                 int root, MPI_Comm comm)
1424 {
1425   int retval = 0;
1426
1427   smpi_bench_end();
1428
1429   if (comm == MPI_COMM_NULL) {
1430     retval = MPI_ERR_COMM;
1431   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1432             ((comm->rank() == root) && (recvtype == MPI_DATATYPE_NULL))){
1433     retval = MPI_ERR_TYPE;
1434   } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) || ((comm->rank() == root) && (recvcount <0))){
1435     retval = MPI_ERR_COUNT;
1436   } else {
1437
1438     char* sendtmpbuf = static_cast<char*>(sendbuf);
1439     int sendtmpcount = sendcount;
1440     MPI_Datatype sendtmptype = sendtype;
1441     if( (comm->rank() == root) && (sendbuf == MPI_IN_PLACE )) {
1442       sendtmpcount=0;
1443       sendtmptype=recvtype;
1444     }
1445     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1446     int root_traced        = comm->group()->index(root);
1447     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1448     extra->type            = TRACING_GATHER;
1449     extra->root            = root_traced;
1450     int known              = 0;
1451     extra->datatype1       = encode_datatype(sendtmptype, &known);
1452     int dt_size_send       = 1;
1453     if (known == 0)
1454       dt_size_send   = sendtmptype->size();
1455     extra->send_size = sendtmpcount * dt_size_send;
1456     extra->datatype2 = encode_datatype(recvtype, &known);
1457     int dt_size_recv = 1;
1458     if ((comm->rank() == root) && known == 0)
1459       dt_size_recv   = recvtype->size();
1460     extra->recv_size = recvcount * dt_size_recv;
1461
1462     TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1463
1464     Colls::gather(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount, recvtype, root, comm);
1465
1466     retval = MPI_SUCCESS;
1467     TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1468   }
1469
1470   smpi_bench_begin();
1471   return retval;
1472 }
1473
1474 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int *recvcounts, int *displs,
1475                 MPI_Datatype recvtype, int root, MPI_Comm comm)
1476 {
1477   int retval = 0;
1478
1479   smpi_bench_end();
1480
1481   if (comm == MPI_COMM_NULL) {
1482     retval = MPI_ERR_COMM;
1483   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1484             ((comm->rank() == root) && (recvtype == MPI_DATATYPE_NULL))){
1485     retval = MPI_ERR_TYPE;
1486   } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1487     retval = MPI_ERR_COUNT;
1488   } else if (recvcounts == nullptr || displs == nullptr) {
1489     retval = MPI_ERR_ARG;
1490   } else {
1491     char* sendtmpbuf = static_cast<char*>(sendbuf);
1492     int sendtmpcount = sendcount;
1493     MPI_Datatype sendtmptype = sendtype;
1494     if( (comm->rank() == root) && (sendbuf == MPI_IN_PLACE )) {
1495       sendtmpcount=0;
1496       sendtmptype=recvtype;
1497     }
1498
1499     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1500     int root_traced        = comm->group()->index(root);
1501     int i                  = 0;
1502     int size               = comm->size();
1503     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1504     extra->type            = TRACING_GATHERV;
1505     extra->num_processes   = size;
1506     extra->root            = root_traced;
1507     int known              = 0;
1508     extra->datatype1       = encode_datatype(sendtmptype, &known);
1509     int dt_size_send       = 1;
1510     if (known == 0)
1511       dt_size_send   = sendtype->size();
1512     extra->send_size = sendtmpcount * dt_size_send;
1513     extra->datatype2 = encode_datatype(recvtype, &known);
1514     int dt_size_recv = 1;
1515     if (known == 0)
1516       dt_size_recv = recvtype->size();
1517     if ((comm->rank() == root)) {
1518       extra->recvcounts = xbt_new(int, size);
1519       for (i                 = 0; i < size; i++) // copy data to avoid bad free
1520         extra->recvcounts[i] = recvcounts[i] * dt_size_recv;
1521     }
1522     TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1523
1524     retval = Colls::gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts, displs, recvtype, root, comm);
1525     TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1526   }
1527
1528   smpi_bench_begin();
1529   return retval;
1530 }
1531
1532 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1533                    void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
1534 {
1535   int retval = 0;
1536
1537   smpi_bench_end();
1538
1539   if (comm == MPI_COMM_NULL) {
1540     retval = MPI_ERR_COMM;
1541   } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1542             (recvtype == MPI_DATATYPE_NULL)){
1543     retval = MPI_ERR_TYPE;
1544   } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1545             (recvcount <0)){
1546     retval = MPI_ERR_COUNT;
1547   } else {
1548     if(sendbuf == MPI_IN_PLACE) {
1549       sendbuf=static_cast<char*>(recvbuf)+recvtype->get_extent()*recvcount*comm->rank();
1550       sendcount=recvcount;
1551       sendtype=recvtype;
1552     }
1553     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1554     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1555     extra->type            = TRACING_ALLGATHER;
1556     int known              = 0;
1557     extra->datatype1       = encode_datatype(sendtype, &known);
1558     int dt_size_send       = 1;
1559     if (known == 0)
1560       dt_size_send   = sendtype->size();
1561     extra->send_size = sendcount * dt_size_send;
1562     extra->datatype2 = encode_datatype(recvtype, &known);
1563     int dt_size_recv = 1;
1564     if (known == 0)
1565       dt_size_recv   = recvtype->size();
1566     extra->recv_size = recvcount * dt_size_recv;
1567
1568     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1569
1570     Colls::allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1571     retval = MPI_SUCCESS;
1572     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1573   }
1574   smpi_bench_begin();
1575   return retval;
1576 }
1577
1578 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1579                    void *recvbuf, int *recvcounts, int *displs, MPI_Datatype recvtype, MPI_Comm comm)
1580 {
1581   int retval = 0;
1582
1583   smpi_bench_end();
1584
1585   if (comm == MPI_COMM_NULL) {
1586     retval = MPI_ERR_COMM;
1587   } else if (((sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) || (recvtype == MPI_DATATYPE_NULL)) {
1588     retval = MPI_ERR_TYPE;
1589   } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1590     retval = MPI_ERR_COUNT;
1591   } else if (recvcounts == nullptr || displs == nullptr) {
1592     retval = MPI_ERR_ARG;
1593   } else {
1594
1595     if(sendbuf == MPI_IN_PLACE) {
1596       sendbuf=static_cast<char*>(recvbuf)+recvtype->get_extent()*displs[comm->rank()];
1597       sendcount=recvcounts[comm->rank()];
1598       sendtype=recvtype;
1599     }
1600     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1601     int i                  = 0;
1602     int size               = comm->size();
1603     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1604     extra->type            = TRACING_ALLGATHERV;
1605     extra->num_processes   = size;
1606     int known              = 0;
1607     extra->datatype1       = encode_datatype(sendtype, &known);
1608     int dt_size_send       = 1;
1609     if (known == 0)
1610       dt_size_send   = sendtype->size();
1611     extra->send_size = sendcount * dt_size_send;
1612     extra->datatype2 = encode_datatype(recvtype, &known);
1613     int dt_size_recv = 1;
1614     if (known == 0)
1615       dt_size_recv    = recvtype->size();
1616     extra->recvcounts = xbt_new(int, size);
1617     for (i                 = 0; i < size; i++) // copy data to avoid bad free
1618       extra->recvcounts[i] = recvcounts[i] * dt_size_recv;
1619
1620     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1621
1622     Colls::allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
1623     retval = MPI_SUCCESS;
1624     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1625   }
1626
1627   smpi_bench_begin();
1628   return retval;
1629 }
1630
1631 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1632                 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
1633 {
1634   int retval = 0;
1635
1636   smpi_bench_end();
1637
1638   if (comm == MPI_COMM_NULL) {
1639     retval = MPI_ERR_COMM;
1640   } else if (((comm->rank() == root) && (!sendtype->is_valid())) ||
1641              ((recvbuf != MPI_IN_PLACE) && (!recvtype->is_valid()))) {
1642     retval = MPI_ERR_TYPE;
1643   } else if ((sendbuf == recvbuf) ||
1644       ((comm->rank()==root) && sendcount>0 && (sendbuf == nullptr))){
1645     retval = MPI_ERR_BUFFER;
1646   }else {
1647
1648     if (recvbuf == MPI_IN_PLACE) {
1649       recvtype  = sendtype;
1650       recvcount = sendcount;
1651     }
1652     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1653     int root_traced        = comm->group()->index(root);
1654     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1655     extra->type            = TRACING_SCATTER;
1656     extra->root            = root_traced;
1657     int known              = 0;
1658     extra->datatype1       = encode_datatype(sendtype, &known);
1659     int dt_size_send       = 1;
1660     if ((comm->rank() == root) && known == 0)
1661       dt_size_send   = sendtype->size();
1662     extra->send_size = sendcount * dt_size_send;
1663     extra->datatype2 = encode_datatype(recvtype, &known);
1664     int dt_size_recv = 1;
1665     if (known == 0)
1666       dt_size_recv   = recvtype->size();
1667     extra->recv_size = recvcount * dt_size_recv;
1668     TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1669
1670     Colls::scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
1671     retval = MPI_SUCCESS;
1672     TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1673   }
1674
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, MPI_Datatype recvtype, int root, MPI_Comm comm)
1681 {
1682   int retval = 0;
1683
1684   smpi_bench_end();
1685
1686   if (comm == MPI_COMM_NULL) {
1687     retval = MPI_ERR_COMM;
1688   } else if (sendcounts == nullptr || displs == nullptr) {
1689     retval = MPI_ERR_ARG;
1690   } else if (((comm->rank() == root) && (sendtype == MPI_DATATYPE_NULL)) ||
1691              ((recvbuf != MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1692     retval = MPI_ERR_TYPE;
1693   } else {
1694     if (recvbuf == MPI_IN_PLACE) {
1695       recvtype  = sendtype;
1696       recvcount = sendcounts[comm->rank()];
1697     }
1698     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1699     int root_traced        = comm->group()->index(root);
1700     int i                  = 0;
1701     int size               = comm->size();
1702     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1703     extra->type            = TRACING_SCATTERV;
1704     extra->num_processes   = size;
1705     extra->root            = root_traced;
1706     int known              = 0;
1707     extra->datatype1       = encode_datatype(sendtype, &known);
1708     int dt_size_send       = 1;
1709     if (known == 0)
1710       dt_size_send = sendtype->size();
1711     if ((comm->rank() == root)) {
1712       extra->sendcounts = xbt_new(int, size);
1713       for (i                 = 0; i < size; i++) // copy data to avoid bad free
1714         extra->sendcounts[i] = sendcounts[i] * dt_size_send;
1715     }
1716     extra->datatype2 = encode_datatype(recvtype, &known);
1717     int dt_size_recv = 1;
1718     if (known == 0)
1719       dt_size_recv   = recvtype->size();
1720     extra->recv_size = recvcount * dt_size_recv;
1721     TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1722
1723     retval = Colls::scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
1724
1725     TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1726   }
1727
1728   smpi_bench_begin();
1729   return retval;
1730 }
1731
1732 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1733 {
1734   int retval = 0;
1735
1736   smpi_bench_end();
1737
1738   if (comm == MPI_COMM_NULL) {
1739     retval = MPI_ERR_COMM;
1740   } else if (!datatype->is_valid() || op == MPI_OP_NULL) {
1741     retval = MPI_ERR_ARG;
1742   } else {
1743     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1744     int root_traced        = comm->group()->index(root);
1745     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1746     extra->type            = TRACING_REDUCE;
1747     int known              = 0;
1748     extra->datatype1       = encode_datatype(datatype, &known);
1749     int dt_size_send       = 1;
1750     if (known == 0)
1751       dt_size_send   = datatype->size();
1752     extra->send_size = count * dt_size_send;
1753     extra->root      = root_traced;
1754
1755     TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1756
1757     Colls::reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1758
1759     retval = MPI_SUCCESS;
1760     TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1761   }
1762
1763   smpi_bench_begin();
1764   return retval;
1765 }
1766
1767 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count, MPI_Datatype datatype, MPI_Op op){
1768   int retval = 0;
1769
1770   smpi_bench_end();
1771   if (!datatype->is_valid() || op == MPI_OP_NULL) {
1772     retval = MPI_ERR_ARG;
1773   } else {
1774     op->apply(inbuf, inoutbuf, &count, datatype);
1775     retval = MPI_SUCCESS;
1776   }
1777   smpi_bench_begin();
1778   return retval;
1779 }
1780
1781 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1782 {
1783   int retval = 0;
1784
1785   smpi_bench_end();
1786
1787   if (comm == MPI_COMM_NULL) {
1788     retval = MPI_ERR_COMM;
1789   } else if (!datatype->is_valid()) {
1790     retval = MPI_ERR_TYPE;
1791   } else if (op == MPI_OP_NULL) {
1792     retval = MPI_ERR_OP;
1793   } else {
1794
1795     char* sendtmpbuf = static_cast<char*>(sendbuf);
1796     if( sendbuf == MPI_IN_PLACE ) {
1797       sendtmpbuf = static_cast<char*>(xbt_malloc(count*datatype->get_extent()));
1798       Datatype::copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
1799     }
1800     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1801     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1802     extra->type            = TRACING_ALLREDUCE;
1803     int known              = 0;
1804     extra->datatype1       = encode_datatype(datatype, &known);
1805     int dt_size_send       = 1;
1806     if (known == 0)
1807       dt_size_send   = datatype->size();
1808     extra->send_size = count * dt_size_send;
1809
1810     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1811
1812     Colls::allreduce(sendtmpbuf, recvbuf, count, datatype, op, comm);
1813
1814     if( sendbuf == MPI_IN_PLACE )
1815       xbt_free(sendtmpbuf);
1816
1817     retval = MPI_SUCCESS;
1818     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1819   }
1820
1821   smpi_bench_begin();
1822   return retval;
1823 }
1824
1825 int PMPI_Scan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1826 {
1827   int retval = 0;
1828
1829   smpi_bench_end();
1830
1831   if (comm == MPI_COMM_NULL) {
1832     retval = MPI_ERR_COMM;
1833   } else if (!datatype->is_valid()) {
1834     retval = MPI_ERR_TYPE;
1835   } else if (op == MPI_OP_NULL) {
1836     retval = MPI_ERR_OP;
1837   } else {
1838     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1839     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1840     extra->type            = TRACING_SCAN;
1841     int known              = 0;
1842     extra->datatype1       = encode_datatype(datatype, &known);
1843     int dt_size_send       = 1;
1844     if (known == 0)
1845       dt_size_send   = datatype->size();
1846     extra->send_size = count * dt_size_send;
1847
1848     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1849
1850     retval = Colls::scan(sendbuf, recvbuf, count, datatype, op, comm);
1851
1852     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1853   }
1854
1855   smpi_bench_begin();
1856   return retval;
1857 }
1858
1859 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm){
1860   int retval = 0;
1861
1862   smpi_bench_end();
1863
1864   if (comm == MPI_COMM_NULL) {
1865     retval = MPI_ERR_COMM;
1866   } else if (!datatype->is_valid()) {
1867     retval = MPI_ERR_TYPE;
1868   } else if (op == MPI_OP_NULL) {
1869     retval = MPI_ERR_OP;
1870   } else {
1871     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1872     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1873     extra->type            = TRACING_EXSCAN;
1874     int known              = 0;
1875     extra->datatype1       = encode_datatype(datatype, &known);
1876     int dt_size_send       = 1;
1877     if (known == 0)
1878       dt_size_send   = datatype->size();
1879     extra->send_size = count * dt_size_send;
1880     void* sendtmpbuf = sendbuf;
1881     if (sendbuf == MPI_IN_PLACE) {
1882       sendtmpbuf = static_cast<void*>(xbt_malloc(count * datatype->size()));
1883       memcpy(sendtmpbuf, recvbuf, count * datatype->size());
1884     }
1885     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1886
1887     retval = Colls::exscan(sendtmpbuf, recvbuf, count, datatype, op, comm);
1888
1889     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1890     if (sendbuf == MPI_IN_PLACE)
1891       xbt_free(sendtmpbuf);
1892   }
1893
1894   smpi_bench_begin();
1895   return retval;
1896 }
1897
1898 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1899 {
1900   int retval = 0;
1901   smpi_bench_end();
1902
1903   if (comm == MPI_COMM_NULL) {
1904     retval = MPI_ERR_COMM;
1905   } else if (!datatype->is_valid()) {
1906     retval = MPI_ERR_TYPE;
1907   } else if (op == MPI_OP_NULL) {
1908     retval = MPI_ERR_OP;
1909   } else if (recvcounts == nullptr) {
1910     retval = MPI_ERR_ARG;
1911   } else {
1912     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1913     int i                  = 0;
1914     int size               = comm->size();
1915     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1916     extra->type            = TRACING_REDUCE_SCATTER;
1917     extra->num_processes   = size;
1918     int known              = 0;
1919     extra->datatype1       = encode_datatype(datatype, &known);
1920     int dt_size_send       = 1;
1921     if (known == 0)
1922       dt_size_send    = datatype->size();
1923     extra->send_size  = 0;
1924     extra->recvcounts = xbt_new(int, size);
1925     int totalcount    = 0;
1926     for (i = 0; i < size; i++) { // copy data to avoid bad free
1927       extra->recvcounts[i] = recvcounts[i] * dt_size_send;
1928       totalcount += recvcounts[i];
1929     }
1930     void* sendtmpbuf = sendbuf;
1931     if (sendbuf == MPI_IN_PLACE) {
1932       sendtmpbuf = static_cast<void*>(xbt_malloc(totalcount * datatype->size()));
1933       memcpy(sendtmpbuf, recvbuf, totalcount * datatype->size());
1934     }
1935
1936     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1937
1938     Colls::reduce_scatter(sendtmpbuf, recvbuf, recvcounts, datatype, op, comm);
1939     retval = MPI_SUCCESS;
1940     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1941
1942     if (sendbuf == MPI_IN_PLACE)
1943       xbt_free(sendtmpbuf);
1944   }
1945
1946   smpi_bench_begin();
1947   return retval;
1948 }
1949
1950 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
1951                               MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1952 {
1953   int retval;
1954   smpi_bench_end();
1955
1956   if (comm == MPI_COMM_NULL) {
1957     retval = MPI_ERR_COMM;
1958   } else if (!datatype->is_valid()) {
1959     retval = MPI_ERR_TYPE;
1960   } else if (op == MPI_OP_NULL) {
1961     retval = MPI_ERR_OP;
1962   } else if (recvcount < 0) {
1963     retval = MPI_ERR_ARG;
1964   } else {
1965     int count = comm->size();
1966
1967     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1968     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
1969     extra->type            = TRACING_REDUCE_SCATTER;
1970     extra->num_processes   = count;
1971     int known              = 0;
1972     extra->datatype1       = encode_datatype(datatype, &known);
1973     int dt_size_send       = 1;
1974     if (known == 0)
1975       dt_size_send    = datatype->size();
1976     extra->send_size  = 0;
1977     extra->recvcounts = xbt_new(int, count);
1978     for (int i             = 0; i < count; i++) // copy data to avoid bad free
1979       extra->recvcounts[i] = recvcount * dt_size_send;
1980     void* sendtmpbuf       = sendbuf;
1981     if (sendbuf == MPI_IN_PLACE) {
1982       sendtmpbuf = static_cast<void*>(xbt_malloc(recvcount * count * datatype->size()));
1983       memcpy(sendtmpbuf, recvbuf, recvcount * count * datatype->size());
1984     }
1985
1986     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1987
1988     int* recvcounts = static_cast<int*>(xbt_malloc(count * sizeof(int)));
1989     for (int i      = 0; i < count; i++)
1990       recvcounts[i] = recvcount;
1991     Colls::reduce_scatter(sendtmpbuf, recvbuf, recvcounts, datatype, op, comm);
1992     xbt_free(recvcounts);
1993     retval = MPI_SUCCESS;
1994
1995     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1996
1997     if (sendbuf == MPI_IN_PLACE)
1998       xbt_free(sendtmpbuf);
1999   }
2000
2001   smpi_bench_begin();
2002   return retval;
2003 }
2004
2005 int PMPI_Alltoall(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount,
2006                   MPI_Datatype recvtype, MPI_Comm comm)
2007 {
2008   int retval = 0;
2009   smpi_bench_end();
2010
2011   if (comm == MPI_COMM_NULL) {
2012     retval = MPI_ERR_COMM;
2013   } else if ((sendbuf != MPI_IN_PLACE && sendtype == MPI_DATATYPE_NULL) || recvtype == MPI_DATATYPE_NULL) {
2014     retval = MPI_ERR_TYPE;
2015   } else {
2016     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2017     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
2018     extra->type            = TRACING_ALLTOALL;
2019
2020     void* sendtmpbuf         = static_cast<char*>(sendbuf);
2021     int sendtmpcount         = sendcount;
2022     MPI_Datatype sendtmptype = sendtype;
2023     if (sendbuf == MPI_IN_PLACE) {
2024       sendtmpbuf = static_cast<void*>(xbt_malloc(recvcount * comm->size() * recvtype->size()));
2025       memcpy(sendtmpbuf, recvbuf, recvcount * comm->size() * recvtype->size());
2026       sendtmpcount = recvcount;
2027       sendtmptype  = recvtype;
2028     }
2029
2030     int known        = 0;
2031     extra->datatype1 = encode_datatype(sendtmptype, &known);
2032     if (known == 0)
2033       extra->send_size = sendtmpcount * sendtmptype->size();
2034     else
2035       extra->send_size = sendtmpcount;
2036     extra->datatype2   = encode_datatype(recvtype, &known);
2037     if (known == 0)
2038       extra->recv_size = recvcount * recvtype->size();
2039     else
2040       extra->recv_size = recvcount;
2041
2042     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
2043
2044     retval = Colls::alltoall(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount, recvtype, comm);
2045
2046     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2047
2048     if (sendbuf == MPI_IN_PLACE)
2049       xbt_free(sendtmpbuf);
2050   }
2051
2052   smpi_bench_begin();
2053   return retval;
2054 }
2055
2056 int PMPI_Alltoallv(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype sendtype, void* recvbuf,
2057                    int* recvcounts, int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2058 {
2059   int retval = 0;
2060
2061   smpi_bench_end();
2062
2063   if (comm == MPI_COMM_NULL) {
2064     retval = MPI_ERR_COMM;
2065   } else if (sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
2066     retval = MPI_ERR_TYPE;
2067   } else if ((sendbuf != MPI_IN_PLACE && (sendcounts == nullptr || senddisps == nullptr)) || recvcounts == nullptr ||
2068              recvdisps == nullptr) {
2069     retval = MPI_ERR_ARG;
2070   } else {
2071     int rank               = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2072     int i                  = 0;
2073     int size               = comm->size();
2074     instr_extra_data extra = xbt_new0(s_instr_extra_data_t, 1);
2075     extra->type            = TRACING_ALLTOALLV;
2076     extra->send_size       = 0;
2077     extra->recv_size       = 0;
2078     extra->recvcounts      = xbt_new(int, size);
2079     extra->sendcounts      = xbt_new(int, size);
2080     int known              = 0;
2081     int dt_size_recv       = 1;
2082     extra->datatype2       = encode_datatype(recvtype, &known);
2083     dt_size_recv           = recvtype->size();
2084
2085     void* sendtmpbuf         = static_cast<char*>(sendbuf);
2086     int* sendtmpcounts       = sendcounts;
2087     int* sendtmpdisps        = senddisps;
2088     MPI_Datatype sendtmptype = sendtype;
2089     int maxsize              = 0;
2090     for (i = 0; i < size; i++) { // copy data to avoid bad free
2091       extra->recv_size += recvcounts[i] * dt_size_recv;
2092       extra->recvcounts[i] = recvcounts[i] * dt_size_recv;
2093       if (((recvdisps[i] + recvcounts[i]) * dt_size_recv) > maxsize)
2094         maxsize = (recvdisps[i] + recvcounts[i]) * dt_size_recv;
2095     }
2096
2097     if (sendbuf == MPI_IN_PLACE) {
2098       sendtmpbuf = static_cast<void*>(xbt_malloc(maxsize));
2099       memcpy(sendtmpbuf, recvbuf, maxsize);
2100       sendtmpcounts = static_cast<int*>(xbt_malloc(size * sizeof(int)));
2101       memcpy(sendtmpcounts, recvcounts, size * sizeof(int));
2102       sendtmpdisps = static_cast<int*>(xbt_malloc(size * sizeof(int)));
2103       memcpy(sendtmpdisps, recvdisps, size * sizeof(int));
2104       sendtmptype = recvtype;
2105     }
2106
2107     extra->datatype1 = encode_datatype(sendtmptype, &known);
2108     int dt_size_send = 1;
2109     dt_size_send     = sendtmptype->size();
2110
2111     for (i = 0; i < size; i++) { // copy data to avoid bad free
2112       extra->send_size += sendtmpcounts[i] * dt_size_send;
2113       extra->sendcounts[i] = sendtmpcounts[i] * dt_size_send;
2114     }
2115     extra->num_processes = size;
2116     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
2117     retval = Colls::alltoallv(sendtmpbuf, sendtmpcounts, sendtmpdisps, sendtmptype, recvbuf, recvcounts,
2118                                     recvdisps, recvtype, comm);
2119     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2120
2121     if (sendbuf == MPI_IN_PLACE) {
2122       xbt_free(sendtmpbuf);
2123       xbt_free(sendtmpcounts);
2124       xbt_free(sendtmpdisps);
2125     }
2126   }
2127
2128   smpi_bench_begin();
2129   return retval;
2130 }
2131
2132
2133 int PMPI_Get_processor_name(char *name, int *resultlen)
2134 {
2135   strncpy(name, SIMIX_host_self()->cname(), strlen(SIMIX_host_self()->cname()) < MPI_MAX_PROCESSOR_NAME - 1
2136                                                 ? strlen(SIMIX_host_self()->cname()) + 1
2137                                                 : MPI_MAX_PROCESSOR_NAME - 1);
2138   *resultlen = strlen(name) > MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2139
2140   return MPI_SUCCESS;
2141 }
2142
2143 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2144 {
2145   if (status == nullptr || count == nullptr) {
2146     return MPI_ERR_ARG;
2147   } else if (!datatype->is_valid()) {
2148     return MPI_ERR_TYPE;
2149   } else {
2150     size_t size = datatype->size();
2151     if (size == 0) {
2152       *count = 0;
2153       return MPI_SUCCESS;
2154     } else if (status->count % size != 0) {
2155       return MPI_UNDEFINED;
2156     } else {
2157       *count = smpi_mpi_get_count(status, datatype);
2158       return MPI_SUCCESS;
2159     }
2160   }
2161 }
2162
2163 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2164   if (old_type == MPI_DATATYPE_NULL) {
2165     return MPI_ERR_TYPE;
2166   } else if (count<0){
2167     return MPI_ERR_COUNT;
2168   } else {
2169     return Datatype::create_contiguous(count, old_type, 0, new_type);
2170   }
2171 }
2172
2173 int PMPI_Type_commit(MPI_Datatype* datatype) {
2174   if (datatype == nullptr || *datatype == MPI_DATATYPE_NULL) {
2175     return MPI_ERR_TYPE;
2176   } else {
2177     (*datatype)->commit();
2178     return MPI_SUCCESS;
2179   }
2180 }
2181
2182 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2183   if (old_type == MPI_DATATYPE_NULL) {
2184     return MPI_ERR_TYPE;
2185   } else if (count<0 || blocklen<0){
2186     return MPI_ERR_COUNT;
2187   } else {
2188     return Datatype::create_vector(count, blocklen, stride, old_type, new_type);
2189   }
2190 }
2191
2192 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2193   if (old_type == MPI_DATATYPE_NULL) {
2194     return MPI_ERR_TYPE;
2195   } else if (count<0 || blocklen<0){
2196     return MPI_ERR_COUNT;
2197   } else {
2198     return Datatype::create_hvector(count, blocklen, stride, old_type, new_type);
2199   }
2200 }
2201
2202 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2203   return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2204 }
2205
2206 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2207   if (old_type == MPI_DATATYPE_NULL) {
2208     return MPI_ERR_TYPE;
2209   } else if (count<0){
2210     return MPI_ERR_COUNT;
2211   } else {
2212     return Datatype::create_indexed(count, blocklens, indices, old_type, new_type);
2213   }
2214 }
2215
2216 int PMPI_Type_create_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2217   if (old_type == MPI_DATATYPE_NULL) {
2218     return MPI_ERR_TYPE;
2219   } else if (count<0){
2220     return MPI_ERR_COUNT;
2221   } else {
2222     return Datatype::create_indexed(count, blocklens, indices, old_type, new_type);
2223   }
2224 }
2225
2226 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type,
2227                                    MPI_Datatype* new_type)
2228 {
2229   if (old_type == MPI_DATATYPE_NULL) {
2230     return MPI_ERR_TYPE;
2231   } else if (count<0){
2232     return MPI_ERR_COUNT;
2233   } else {
2234     int* blocklens=static_cast<int*>(xbt_malloc(blocklength*count*sizeof(int)));
2235     for (int i    = 0; i < count; i++)
2236       blocklens[i]=blocklength;
2237     int retval    = Datatype::create_indexed(count, blocklens, indices, old_type, new_type);
2238     xbt_free(blocklens);
2239     return retval;
2240   }
2241 }
2242
2243 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type)
2244 {
2245   if (old_type == MPI_DATATYPE_NULL) {
2246     return MPI_ERR_TYPE;
2247   } else if (count<0){
2248     return MPI_ERR_COUNT;
2249   } else {
2250     return Datatype::create_hindexed(count, blocklens, indices, old_type, new_type);
2251   }
2252 }
2253
2254 int PMPI_Type_create_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type,
2255                               MPI_Datatype* new_type) {
2256   return PMPI_Type_hindexed(count, blocklens,indices,old_type,new_type);
2257 }
2258
2259 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type,
2260                                     MPI_Datatype* new_type) {
2261   if (old_type == MPI_DATATYPE_NULL) {
2262     return MPI_ERR_TYPE;
2263   } else if (count<0){
2264     return MPI_ERR_COUNT;
2265   } else {
2266     int* blocklens=(int*)xbt_malloc(blocklength*count*sizeof(int));
2267     for (int i     = 0; i < count; i++)
2268       blocklens[i] = blocklength;
2269     int retval     = Datatype::create_hindexed(count, blocklens, indices, old_type, new_type);
2270     xbt_free(blocklens);
2271     return retval;
2272   }
2273 }
2274
2275 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2276   if (count<0){
2277     return MPI_ERR_COUNT;
2278   } else {
2279     return Datatype::create_struct(count, blocklens, indices, old_types, new_type);
2280   }
2281 }
2282
2283 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types,
2284                             MPI_Datatype* new_type) {
2285   return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2286 }
2287
2288 int PMPI_Error_class(int errorcode, int* errorclass) {
2289   // assume smpi uses only standard mpi error codes
2290   *errorclass=errorcode;
2291   return MPI_SUCCESS;
2292 }
2293
2294 int PMPI_Initialized(int* flag) {
2295    *flag=smpi_process_initialized();
2296    return MPI_SUCCESS;
2297 }
2298
2299 /* The topo part of MPI_COMM_WORLD should always be nullptr. When other topologies will be implemented, not only should we
2300  * check if the topology is nullptr, but we should check if it is the good topology type (so we have to add a
2301  *  MPIR_Topo_Type field, and replace the MPI_Topology field by an union)*/
2302
2303 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periodic, int reorder, MPI_Comm* comm_cart) {
2304   if (comm_old == MPI_COMM_NULL){
2305     return MPI_ERR_COMM;
2306   } else if (ndims < 0 || (ndims > 0 && (dims == nullptr || periodic == nullptr)) || comm_cart == nullptr) {
2307     return MPI_ERR_ARG;
2308   } else{
2309     new Topo_Cart(comm_old, ndims, dims, periodic, reorder, comm_cart);
2310     return MPI_SUCCESS;
2311   }
2312 }
2313
2314 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2315   if(comm == MPI_COMM_NULL || comm->topo() == nullptr) {
2316     return MPI_ERR_TOPOLOGY;
2317   }
2318   if (coords == nullptr) {
2319     return MPI_ERR_ARG;
2320   }
2321   MPIR_Cart_Topology topo = static_cast<MPIR_Cart_Topology>(comm->topo());
2322   if (topo==nullptr) {
2323     return MPI_ERR_ARG;
2324   }
2325   return topo->rank(coords, rank);
2326 }
2327
2328 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2329   if(comm == MPI_COMM_NULL || comm->topo() == nullptr) {
2330     return MPI_ERR_TOPOLOGY;
2331   }
2332   if (source == nullptr || dest == nullptr || direction < 0 ) {
2333     return MPI_ERR_ARG;
2334   }
2335   MPIR_Cart_Topology topo = static_cast<MPIR_Cart_Topology>(comm->topo());
2336   if (topo==nullptr) {
2337     return MPI_ERR_ARG;
2338   }
2339   return topo->shift(direction, displ, source, dest);
2340 }
2341
2342 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2343   if(comm == MPI_COMM_NULL || comm->topo() == nullptr) {
2344     return MPI_ERR_TOPOLOGY;
2345   }
2346   if (rank < 0 || rank >= comm->size()) {
2347     return MPI_ERR_RANK;
2348   }
2349   if (maxdims <= 0) {
2350     return MPI_ERR_ARG;
2351   }
2352   if(coords == nullptr) {
2353     return MPI_ERR_ARG;
2354   }
2355   MPIR_Cart_Topology topo = static_cast<MPIR_Cart_Topology>(comm->topo());
2356   if (topo==nullptr) {
2357     return MPI_ERR_ARG;
2358   }
2359   return topo->coords(rank, maxdims, coords);
2360 }
2361
2362 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2363   if(comm == nullptr || comm->topo() == nullptr) {
2364     return MPI_ERR_TOPOLOGY;
2365   }
2366   if(maxdims <= 0 || dims == nullptr || periods == nullptr || coords == nullptr) {
2367     return MPI_ERR_ARG;
2368   }
2369   MPIR_Cart_Topology topo = static_cast<MPIR_Cart_Topology>(comm->topo());
2370   if (topo==nullptr) {
2371     return MPI_ERR_ARG;
2372   }
2373   return topo->get(maxdims, dims, periods, coords);
2374 }
2375
2376 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2377   if (comm == MPI_COMM_NULL || comm->topo() == nullptr) {
2378     return MPI_ERR_TOPOLOGY;
2379   }
2380   if (ndims == nullptr) {
2381     return MPI_ERR_ARG;
2382   }
2383   MPIR_Cart_Topology topo = static_cast<MPIR_Cart_Topology>(comm->topo());
2384   if (topo==nullptr) {
2385     return MPI_ERR_ARG;
2386   }
2387   return topo->dim_get(ndims);
2388 }
2389
2390 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2391   if(dims == nullptr) {
2392     return MPI_ERR_ARG;
2393   }
2394   if (ndims < 1 || nnodes < 1) {
2395     return MPI_ERR_DIMS;
2396   }
2397   return Dims_create(nnodes, ndims, dims);
2398 }
2399
2400 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2401   if(comm == MPI_COMM_NULL || comm->topo() == nullptr) {
2402     return MPI_ERR_TOPOLOGY;
2403   }
2404   if (comm_new == nullptr) {
2405     return MPI_ERR_ARG;
2406   }
2407   MPIR_Cart_Topology topo = static_cast<MPIR_Cart_Topology>(comm->topo());
2408   if (topo==nullptr) {
2409     return MPI_ERR_ARG;
2410   }
2411   MPIR_Cart_Topology cart = topo->sub(remain_dims, comm_new);
2412   if(cart==nullptr)
2413     return  MPI_ERR_ARG;
2414   return MPI_SUCCESS;
2415 }
2416
2417 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2418   if (oldtype == MPI_DATATYPE_NULL) {
2419     return MPI_ERR_TYPE;
2420   }
2421   int blocks[3]         = {1, 1, 1};
2422   MPI_Aint disps[3]     = {lb, 0, lb + extent};
2423   MPI_Datatype types[3] = {MPI_LB, oldtype, MPI_UB};
2424
2425   *newtype = new Type_Struct(oldtype->size(), lb, lb + extent, DT_FLAG_DERIVED, 3, blocks, disps, types);
2426
2427   (*newtype)->addflag(~DT_FLAG_COMMITED);
2428   return MPI_SUCCESS;
2429 }
2430
2431 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2432   int retval = 0;
2433   smpi_bench_end();
2434   if (comm == MPI_COMM_NULL) {
2435     retval= MPI_ERR_COMM;
2436   }else if ((base == nullptr && size != 0) || disp_unit <= 0 || size < 0 ){
2437     retval= MPI_ERR_OTHER;
2438   }else{
2439     *win = new Win( base, size, disp_unit, info, comm);
2440     retval = MPI_SUCCESS;
2441   }
2442   smpi_bench_begin();
2443   return retval;
2444 }
2445
2446 int PMPI_Win_free( MPI_Win* win){
2447   int retval = 0;
2448   smpi_bench_end();
2449   if (win == nullptr || *win == MPI_WIN_NULL) {
2450     retval = MPI_ERR_WIN;
2451   }else{
2452     delete *win;
2453     retval=MPI_SUCCESS;
2454   }
2455   smpi_bench_begin();
2456   return retval;
2457 }
2458
2459 int PMPI_Win_set_name(MPI_Win  win, char * name)
2460 {
2461   if (win == MPI_WIN_NULL)  {
2462     return MPI_ERR_TYPE;
2463   } else if (name == nullptr)  {
2464     return MPI_ERR_ARG;
2465   } else {
2466     win->set_name(name);
2467     return MPI_SUCCESS;
2468   }
2469 }
2470
2471 int PMPI_Win_get_name(MPI_Win  win, char * name, int* len)
2472 {
2473   if (win == MPI_WIN_NULL)  {
2474     return MPI_ERR_WIN;
2475   } else if (name == nullptr)  {
2476     return MPI_ERR_ARG;
2477   } else {
2478     win->get_name(name, len);
2479     return MPI_SUCCESS;
2480   }
2481 }
2482
2483 int PMPI_Win_get_group(MPI_Win  win, MPI_Group * group){
2484   if (win == MPI_WIN_NULL)  {
2485     return MPI_ERR_WIN;
2486   }else {
2487     win->get_group(group);
2488     (*group)->ref();
2489     return MPI_SUCCESS;
2490   }
2491 }
2492
2493 int PMPI_Win_fence( int assert,  MPI_Win win){
2494   int retval = 0;
2495   smpi_bench_end();
2496   if (win == MPI_WIN_NULL) {
2497     retval = MPI_ERR_WIN;
2498   } else {
2499   int rank = smpi_process_index();
2500   TRACE_smpi_collective_in(rank, -1, __FUNCTION__, nullptr);
2501   retval = win->fence(assert);
2502   TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2503   }
2504   smpi_bench_begin();
2505   return retval;
2506 }
2507
2508 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2509               MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2510   int retval = 0;
2511   smpi_bench_end();
2512   if (win == MPI_WIN_NULL) {
2513     retval = MPI_ERR_WIN;
2514   } else if (target_rank == MPI_PROC_NULL) {
2515     retval = MPI_SUCCESS;
2516   } else if (target_rank <0){
2517     retval = MPI_ERR_RANK;
2518   } else if (target_disp <0){
2519     retval = MPI_ERR_ARG;
2520   } else if ((origin_count < 0 || target_count < 0) ||
2521              (origin_addr==nullptr && origin_count > 0)){
2522     retval = MPI_ERR_COUNT;
2523   } else if ((!origin_datatype->is_valid()) || (!target_datatype->is_valid())) {
2524     retval = MPI_ERR_TYPE;
2525   } else {
2526     int rank = smpi_process_index();
2527     MPI_Group group;
2528     win->get_group(&group);
2529     int src_traced = group->index(target_rank);
2530     TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, nullptr);
2531
2532     retval = win->get( origin_addr, origin_count, origin_datatype, target_rank, target_disp, target_count,
2533                            target_datatype);
2534
2535     TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
2536   }
2537   smpi_bench_begin();
2538   return retval;
2539 }
2540
2541 int PMPI_Put( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2542               MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2543   int retval = 0;
2544   smpi_bench_end();
2545   if (win == MPI_WIN_NULL) {
2546     retval = MPI_ERR_WIN;
2547   } else if (target_rank == MPI_PROC_NULL) {
2548     retval = MPI_SUCCESS;
2549   } else if (target_rank <0){
2550     retval = MPI_ERR_RANK;
2551   } else if (target_disp <0){
2552     retval = MPI_ERR_ARG;
2553   } else if ((origin_count < 0 || target_count < 0) ||
2554             (origin_addr==nullptr && origin_count > 0)){
2555     retval = MPI_ERR_COUNT;
2556   } else if ((!origin_datatype->is_valid()) || (!target_datatype->is_valid())) {
2557     retval = MPI_ERR_TYPE;
2558   } else {
2559     int rank = smpi_process_index();
2560     MPI_Group group;
2561     win->get_group(&group);
2562     int dst_traced = group->index(target_rank);
2563     TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, nullptr);
2564     TRACE_smpi_send(rank, rank, dst_traced, SMPI_RMA_TAG, origin_count*origin_datatype->size());
2565
2566     retval = win->put( origin_addr, origin_count, origin_datatype, target_rank, target_disp, target_count,
2567                            target_datatype);
2568
2569     TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
2570   }
2571   smpi_bench_begin();
2572   return retval;
2573 }
2574
2575 int PMPI_Accumulate( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2576               MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Op op, MPI_Win win){
2577   int retval = 0;
2578   smpi_bench_end();
2579   if (win == MPI_WIN_NULL) {
2580     retval = MPI_ERR_WIN;
2581   } else if (target_rank == MPI_PROC_NULL) {
2582     retval = MPI_SUCCESS;
2583   } else if (target_rank <0){
2584     retval = MPI_ERR_RANK;
2585   } else if (target_disp <0){
2586     retval = MPI_ERR_ARG;
2587   } else if ((origin_count < 0 || target_count < 0) ||
2588              (origin_addr==nullptr && origin_count > 0)){
2589     retval = MPI_ERR_COUNT;
2590   } else if ((!origin_datatype->is_valid()) ||
2591             (!target_datatype->is_valid())) {
2592     retval = MPI_ERR_TYPE;
2593   } else if (op == MPI_OP_NULL) {
2594     retval = MPI_ERR_OP;
2595   } else {
2596     int rank = smpi_process_index();
2597     MPI_Group group;
2598     win->get_group(&group);
2599     int src_traced = group->index(target_rank);
2600     TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, nullptr);
2601
2602     retval = win->accumulate( origin_addr, origin_count, origin_datatype, target_rank, target_disp, target_count,
2603                                   target_datatype, op);
2604
2605     TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
2606   }
2607   smpi_bench_begin();
2608   return retval;
2609 }
2610
2611 int PMPI_Win_post(MPI_Group group, int assert, MPI_Win win){
2612   int retval = 0;
2613   smpi_bench_end();
2614   if (win == MPI_WIN_NULL) {
2615     retval = MPI_ERR_WIN;
2616   } else if (group==MPI_GROUP_NULL){
2617     retval = MPI_ERR_GROUP;
2618   } else {
2619     int rank = smpi_process_index();
2620     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, nullptr);
2621     retval = win->post(group,assert);
2622     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2623   }
2624   smpi_bench_begin();
2625   return retval;
2626 }
2627
2628 int PMPI_Win_start(MPI_Group group, int assert, MPI_Win win){
2629   int retval = 0;
2630   smpi_bench_end();
2631   if (win == MPI_WIN_NULL) {
2632     retval = MPI_ERR_WIN;
2633   } else if (group==MPI_GROUP_NULL){
2634     retval = MPI_ERR_GROUP;
2635   } else {
2636     int rank = smpi_process_index();
2637     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, nullptr);
2638     retval = win->start(group,assert);
2639     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2640   }
2641   smpi_bench_begin();
2642   return retval;
2643 }
2644
2645 int PMPI_Win_complete(MPI_Win win){
2646   int retval = 0;
2647   smpi_bench_end();
2648   if (win == MPI_WIN_NULL) {
2649     retval = MPI_ERR_WIN;
2650   } else {
2651     int rank = smpi_process_index();
2652     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, nullptr);
2653
2654     retval = win->complete();
2655
2656     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2657   }
2658   smpi_bench_begin();
2659   return retval;
2660 }
2661
2662 int PMPI_Win_wait(MPI_Win win){
2663   int retval = 0;
2664   smpi_bench_end();
2665   if (win == MPI_WIN_NULL) {
2666     retval = MPI_ERR_WIN;
2667   } else {
2668     int rank = smpi_process_index();
2669     TRACE_smpi_collective_in(rank, -1, __FUNCTION__, nullptr);
2670
2671     retval = win->wait();
2672
2673     TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2674   }
2675   smpi_bench_begin();
2676   return retval;
2677 }
2678
2679 int PMPI_Alloc_mem(MPI_Aint size, MPI_Info info, void *baseptr){
2680   void *ptr = xbt_malloc(size);
2681   if(ptr==nullptr)
2682     return MPI_ERR_NO_MEM;
2683   else {
2684     *static_cast<void**>(baseptr) = ptr;
2685     return MPI_SUCCESS;
2686   }
2687 }
2688
2689 int PMPI_Free_mem(void *baseptr){
2690   xbt_free(baseptr);
2691   return MPI_SUCCESS;
2692 }
2693
2694 int PMPI_Type_set_name(MPI_Datatype  datatype, char * name)
2695 {
2696   if (datatype == MPI_DATATYPE_NULL)  {
2697     return MPI_ERR_TYPE;
2698   } else if (name == nullptr)  {
2699     return MPI_ERR_ARG;
2700   } else {
2701     datatype->set_name(name);
2702     return MPI_SUCCESS;
2703   }
2704 }
2705
2706 int PMPI_Type_get_name(MPI_Datatype  datatype, char * name, int* len)
2707 {
2708   if (datatype == MPI_DATATYPE_NULL)  {
2709     return MPI_ERR_TYPE;
2710   } else if (name == nullptr)  {
2711     return MPI_ERR_ARG;
2712   } else {
2713     datatype->get_name(name, len);
2714     return MPI_SUCCESS;
2715   }
2716 }
2717
2718 MPI_Datatype PMPI_Type_f2c(MPI_Fint datatype){
2719   return static_cast<MPI_Datatype>(F2C::f2c(datatype));
2720 }
2721
2722 MPI_Fint PMPI_Type_c2f(MPI_Datatype datatype){
2723   return datatype->c2f();
2724 }
2725
2726 MPI_Group PMPI_Group_f2c(MPI_Fint group){
2727   return Group::f2c(group);
2728 }
2729
2730 MPI_Fint PMPI_Group_c2f(MPI_Group group){
2731   return group->c2f();
2732 }
2733
2734 MPI_Request PMPI_Request_f2c(MPI_Fint request){
2735   return static_cast<MPI_Request>(Request::f2c(request));
2736 }
2737
2738 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
2739   return request->c2f();
2740 }
2741
2742 MPI_Win PMPI_Win_f2c(MPI_Fint win){
2743   return static_cast<MPI_Win>(Win::f2c(win));
2744 }
2745
2746 MPI_Fint PMPI_Win_c2f(MPI_Win win){
2747   return win->c2f();
2748 }
2749
2750 MPI_Op PMPI_Op_f2c(MPI_Fint op){
2751   return static_cast<MPI_Op>(Op::f2c(op));
2752 }
2753
2754 MPI_Fint PMPI_Op_c2f(MPI_Op op){
2755   return op->c2f();
2756 }
2757
2758 MPI_Comm PMPI_Comm_f2c(MPI_Fint comm){
2759   return static_cast<MPI_Comm>(Comm::f2c(comm));
2760 }
2761
2762 MPI_Fint PMPI_Comm_c2f(MPI_Comm comm){
2763   return comm->c2f();
2764 }
2765
2766 MPI_Info PMPI_Info_f2c(MPI_Fint info){
2767   return static_cast<MPI_Info>(Info::f2c(info));
2768 }
2769
2770 MPI_Fint PMPI_Info_c2f(MPI_Info info){
2771   return info->c2f();
2772 }
2773
2774 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2775   return Comm::keyval_create(copy_fn, delete_fn, keyval, extra_state);
2776 }
2777
2778 int PMPI_Keyval_free(int* keyval) {
2779   return Comm::keyval_free(keyval);
2780 }
2781
2782 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2783   if(keyval == MPI_TAG_UB||keyval == MPI_HOST||keyval == MPI_IO ||keyval == MPI_WTIME_IS_GLOBAL||keyval == MPI_APPNUM
2784        ||keyval == MPI_UNIVERSE_SIZE||keyval == MPI_LASTUSEDCODE)
2785     return MPI_ERR_ARG;
2786   else if (comm==MPI_COMM_NULL)
2787     return MPI_ERR_COMM;
2788   else
2789     return comm->attr_delete(keyval);
2790 }
2791
2792 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2793   static int one = 1;
2794   static int zero = 0;
2795   static int tag_ub = 1000000;
2796   static int last_used_code = MPI_ERR_LASTCODE;
2797
2798   if (comm==MPI_COMM_NULL){
2799     *flag = 0;
2800     return MPI_ERR_COMM;
2801   }
2802
2803   switch (keyval) {
2804   case MPI_HOST:
2805   case MPI_IO:
2806   case MPI_APPNUM:
2807     *flag = 1;
2808     *static_cast<int**>(attr_value) = &zero;
2809     return MPI_SUCCESS;
2810   case MPI_UNIVERSE_SIZE:
2811     *flag = 1;
2812     *static_cast<int**>(attr_value) = &smpi_universe_size;
2813     return MPI_SUCCESS;
2814   case MPI_LASTUSEDCODE:
2815     *flag = 1;
2816     *static_cast<int**>(attr_value) = &last_used_code;
2817     return MPI_SUCCESS;
2818   case MPI_TAG_UB:
2819     *flag=1;
2820     *static_cast<int**>(attr_value) = &tag_ub;
2821     return MPI_SUCCESS;
2822   case MPI_WTIME_IS_GLOBAL:
2823     *flag = 1;
2824     *static_cast<int**>(attr_value) = &one;
2825     return MPI_SUCCESS;
2826   default:
2827     return comm->attr_get(keyval, attr_value, flag);
2828   }
2829 }
2830
2831 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2832   if(keyval == MPI_TAG_UB||keyval == MPI_HOST||keyval == MPI_IO ||keyval == MPI_WTIME_IS_GLOBAL||keyval == MPI_APPNUM
2833        ||keyval == MPI_UNIVERSE_SIZE||keyval == MPI_LASTUSEDCODE)
2834     return MPI_ERR_ARG;
2835   else if (comm==MPI_COMM_NULL)
2836     return MPI_ERR_COMM;
2837   else
2838   return comm->attr_put(keyval, attr_value);
2839 }
2840
2841 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2842 {
2843   return PMPI_Attr_get(comm, comm_keyval, attribute_val,flag);
2844 }
2845
2846 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2847 {
2848   return PMPI_Attr_put(comm, comm_keyval, attribute_val);
2849 }
2850
2851 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2852 {
2853   return PMPI_Attr_delete(comm, comm_keyval);
2854 }
2855
2856 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval,
2857                             void* extra_state)
2858 {
2859   return PMPI_Keyval_create(copy_fn, delete_fn, keyval, extra_state);
2860 }
2861
2862 int PMPI_Comm_free_keyval(int* keyval) {
2863   return PMPI_Keyval_free(keyval);
2864 }
2865
2866 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2867 {
2868   if (type==MPI_DATATYPE_NULL)
2869     return MPI_ERR_TYPE;
2870   else
2871     return type->attr_get(type_keyval, attribute_val, flag);
2872 }
2873
2874 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2875 {
2876   if (type==MPI_DATATYPE_NULL)
2877     return MPI_ERR_TYPE;
2878   else
2879     return type->attr_put(type_keyval, attribute_val);
2880 }
2881
2882 int PMPI_Type_delete_attr (MPI_Datatype type, int type_keyval)
2883 {
2884   if (type==MPI_DATATYPE_NULL)
2885     return MPI_ERR_TYPE;
2886   else
2887     return type->attr_delete(type_keyval);
2888 }
2889
2890 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval,
2891                             void* extra_state)
2892 {
2893   return Datatype::keyval_create(copy_fn, delete_fn, keyval, extra_state);
2894 }
2895
2896 int PMPI_Type_free_keyval(int* keyval) {
2897   return Datatype::keyval_free(keyval);
2898 }
2899
2900 int PMPI_Info_create( MPI_Info *info){
2901   if (info == nullptr)
2902     return MPI_ERR_ARG;
2903   *info = new Info();
2904   return MPI_SUCCESS;
2905 }
2906
2907 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2908   if (info == nullptr || key == nullptr || value == nullptr)
2909     return MPI_ERR_ARG;
2910   info->set(key, value);
2911   return MPI_SUCCESS;
2912 }
2913
2914 int PMPI_Info_free( MPI_Info *info){
2915   if (info == nullptr || *info==nullptr)
2916     return MPI_ERR_ARG;
2917   Info::unref(*info);
2918   *info=MPI_INFO_NULL;
2919   return MPI_SUCCESS;
2920 }
2921
2922 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
2923   *flag=false;
2924   if (info == nullptr || key == nullptr || valuelen <0)
2925     return MPI_ERR_ARG;
2926   if (value == nullptr)
2927     return MPI_ERR_INFO_VALUE;
2928   return info->get(key, valuelen, value, flag);
2929 }
2930
2931 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
2932   if (info == nullptr || newinfo==nullptr)
2933     return MPI_ERR_ARG;
2934   *newinfo = new Info(info);
2935   return MPI_SUCCESS;
2936 }
2937
2938 int PMPI_Info_delete(MPI_Info info, char *key){
2939   if (info == nullptr || key==nullptr)
2940     return MPI_ERR_ARG;
2941   return info->remove(key);
2942 }
2943
2944 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
2945   if (info == nullptr || nkeys==nullptr)
2946     return MPI_ERR_ARG;
2947   return info->get_nkeys(nkeys);
2948 }
2949
2950 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
2951   if (info == nullptr || key==nullptr || n<0 || n> MPI_MAX_INFO_KEY)
2952     return MPI_ERR_ARG;
2953   return info->get_nthkey(n, key);
2954 }
2955
2956 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
2957   *flag=false;
2958   if (info == nullptr || key == nullptr || valuelen==nullptr)
2959     return MPI_ERR_ARG;
2960   return info->get_valuelen(key, valuelen, flag);
2961 }
2962
2963 int PMPI_Unpack(void* inbuf, int incount, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2964   if(incount<0 || outcount < 0 || inbuf==nullptr || outbuf==nullptr)
2965     return MPI_ERR_ARG;
2966   if(!type->is_valid())
2967     return MPI_ERR_TYPE;
2968   if(comm==MPI_COMM_NULL)
2969     return MPI_ERR_COMM;
2970   return type->unpack(inbuf, incount, position, outbuf,outcount, comm);
2971 }
2972
2973 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2974   if(incount<0 || outcount < 0|| inbuf==nullptr || outbuf==nullptr)
2975     return MPI_ERR_ARG;
2976   if(!type->is_valid())
2977     return MPI_ERR_TYPE;
2978   if(comm==MPI_COMM_NULL)
2979     return MPI_ERR_COMM;
2980   return type->pack(inbuf, incount, outbuf,outcount,position, comm);
2981 }
2982
2983 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2984   if(incount<0)
2985     return MPI_ERR_ARG;
2986   if(!datatype->is_valid())
2987     return MPI_ERR_TYPE;
2988   if(comm==MPI_COMM_NULL)
2989     return MPI_ERR_COMM;
2990
2991   *size=incount*datatype->size();
2992
2993   return MPI_SUCCESS;
2994 }
2995
2996 } // extern "C"