1 /* Copyright (c) 2007, 2008, 2009, 2010. The SimGrid Team.
2 * All rights reserved. */
4 /* This program is free software; you can redistribute it and/or modify it
5 * under the terms of the license (GNU LGPL) which comes with this package. */
8 #include "smpi_coll_private.h"
9 #include "smpi_mpi_dt_private.h"
11 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_mpi, smpi,
12 "Logging specific to SMPI (mpi)");
14 /* MPI User level calls */
16 int MPI_Init(int* argc, char*** argv) {
17 smpi_process_init(argc, argv);
19 TRACE_smpi_init(smpi_process_index());
21 smpi_bench_begin(-1, NULL);
25 int MPI_Finalize(void) {
26 smpi_bench_end(-1, NULL);
28 TRACE_smpi_finalize(smpi_process_index());
30 smpi_process_destroy();
34 int MPI_Init_thread(int* argc, char*** argv, int required, int* provided) {
35 if(provided != NULL) {
36 *provided = MPI_THREAD_MULTIPLE;
38 return MPI_Init(argc, argv);
41 int MPI_Query_thread(int* provided) {
44 smpi_bench_end(-1, NULL);
45 if(provided == NULL) {
48 *provided = MPI_THREAD_MULTIPLE;
51 smpi_bench_begin(-1, NULL);
55 int MPI_Is_thread_main(int* flag) {
58 smpi_bench_end(-1, NULL);
62 *flag = smpi_process_index() == 0;
65 smpi_bench_begin(-1, NULL);
69 int MPI_Abort(MPI_Comm comm, int errorcode) {
70 smpi_bench_end(-1, NULL);
71 smpi_process_destroy();
72 // FIXME: should kill all processes in comm instead
73 SIMIX_process_kill(SIMIX_process_self());
77 double MPI_Wtime(void) {
80 smpi_bench_end(-1, NULL);
81 time = SIMIX_get_clock();
82 smpi_bench_begin(-1, NULL);
86 int MPI_Address(void *location, MPI_Aint *address) {
89 smpi_bench_end(-1, NULL);
93 *address = (MPI_Aint)location;
95 smpi_bench_begin(-1, NULL);
99 int MPI_Type_free(MPI_Datatype* datatype) {
102 smpi_bench_end(-1, NULL);
104 retval = MPI_ERR_ARG;
106 // FIXME: always fail for now
107 retval = MPI_ERR_TYPE;
109 smpi_bench_begin(-1, NULL);
113 int MPI_Type_size(MPI_Datatype datatype, int* size) {
116 smpi_bench_end(-1, NULL);
117 if(datatype == MPI_DATATYPE_NULL) {
118 retval = MPI_ERR_TYPE;
119 } else if(size == NULL) {
120 retval = MPI_ERR_ARG;
122 *size = (int)smpi_datatype_size(datatype);
123 retval = MPI_SUCCESS;
125 smpi_bench_begin(-1, NULL);
129 int MPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint* lb, MPI_Aint* extent) {
132 smpi_bench_end(-1, NULL);
133 if(datatype == MPI_DATATYPE_NULL) {
134 retval = MPI_ERR_TYPE;
135 } else if(lb == NULL || extent == NULL) {
136 retval = MPI_ERR_ARG;
138 retval = smpi_datatype_extent(datatype, lb, extent);
140 smpi_bench_begin(-1, NULL);
144 int MPI_Type_extent(MPI_Datatype datatype, MPI_Aint* extent) {
148 smpi_bench_end(-1, NULL);
149 if(datatype == MPI_DATATYPE_NULL) {
150 retval = MPI_ERR_TYPE;
151 } else if(extent == NULL) {
152 retval = MPI_ERR_ARG;
154 retval = smpi_datatype_extent(datatype, &dummy, extent);
156 smpi_bench_begin(-1, NULL);
160 int MPI_Type_lb(MPI_Datatype datatype, MPI_Aint* disp) {
163 smpi_bench_end(-1, NULL);
164 if(datatype == MPI_DATATYPE_NULL) {
165 retval = MPI_ERR_TYPE;
166 } else if(disp == NULL) {
167 retval = MPI_ERR_ARG;
169 *disp = smpi_datatype_lb(datatype);
170 retval = MPI_SUCCESS;
172 smpi_bench_begin(-1, NULL);
176 int MPI_Type_ub(MPI_Datatype datatype, MPI_Aint* disp) {
179 smpi_bench_end(-1, NULL);
180 if(datatype == MPI_DATATYPE_NULL) {
181 retval = MPI_ERR_TYPE;
182 } else if(disp == NULL) {
183 retval = MPI_ERR_ARG;
185 *disp = smpi_datatype_ub(datatype);
186 retval = MPI_SUCCESS;
188 smpi_bench_begin(-1, NULL);
192 int MPI_Op_create(MPI_User_function* function, int commute, MPI_Op* op) {
195 smpi_bench_end(-1, NULL);
196 if(function == NULL || op == NULL) {
197 retval = MPI_ERR_ARG;
199 *op = smpi_op_new(function, commute);
200 retval = MPI_SUCCESS;
202 smpi_bench_begin(-1, NULL);
206 int MPI_Op_free(MPI_Op* op) {
209 smpi_bench_end(-1, NULL);
211 retval = MPI_ERR_ARG;
212 } else if(*op == MPI_OP_NULL) {
215 smpi_op_destroy(*op);
217 retval = MPI_SUCCESS;
219 smpi_bench_begin(-1, NULL);
223 int MPI_Group_free(MPI_Group *group) {
226 smpi_bench_end(-1, NULL);
228 retval = MPI_ERR_ARG;
230 smpi_group_destroy(*group);
231 *group = MPI_GROUP_NULL;
232 retval = MPI_SUCCESS;
234 smpi_bench_begin(-1, NULL);
238 int MPI_Group_size(MPI_Group group, int* size) {
241 smpi_bench_end(-1, NULL);
242 if(group == MPI_GROUP_NULL) {
243 retval = MPI_ERR_GROUP;
244 } else if(size == NULL) {
245 retval = MPI_ERR_ARG;
247 *size = smpi_group_size(group);
248 retval = MPI_SUCCESS;
250 smpi_bench_begin(-1, NULL);
254 int MPI_Group_rank(MPI_Group group, int* rank) {
257 smpi_bench_end(-1, NULL);
258 if(group == MPI_GROUP_NULL) {
259 retval = MPI_ERR_GROUP;
260 } else if(rank == NULL) {
261 retval = MPI_ERR_ARG;
263 *rank = smpi_group_rank(group, smpi_process_index());
264 retval = MPI_SUCCESS;
266 smpi_bench_begin(-1, NULL);
270 int MPI_Group_translate_ranks (MPI_Group group1, int n, int* ranks1, MPI_Group group2, int* ranks2) {
271 int retval, i, index;
273 smpi_bench_end(-1, NULL);
274 if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
275 retval = MPI_ERR_GROUP;
277 for(i = 0; i < n; i++) {
278 index = smpi_group_index(group1, ranks1[i]);
279 ranks2[i] = smpi_group_rank(group2, index);
281 retval = MPI_SUCCESS;
283 smpi_bench_begin(-1, NULL);
287 int MPI_Group_compare(MPI_Group group1, MPI_Group group2, int* result) {
290 smpi_bench_end(-1, NULL);
291 if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
292 retval = MPI_ERR_GROUP;
293 } else if(result == NULL) {
294 retval = MPI_ERR_ARG;
296 *result = smpi_group_compare(group1, group2);
297 retval = MPI_SUCCESS;
299 smpi_bench_begin(-1, NULL);
303 int MPI_Group_union(MPI_Group group1, MPI_Group group2, MPI_Group* newgroup) {
304 int retval, i, proc1, proc2, size, size2;
306 smpi_bench_end(-1, NULL);
307 if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
308 retval = MPI_ERR_GROUP;
309 } else if(newgroup == NULL) {
310 retval = MPI_ERR_ARG;
312 size = smpi_group_size(group1);
313 size2 = smpi_group_size(group2);
314 for(i = 0; i < size2; i++) {
315 proc2 = smpi_group_index(group2, i);
316 proc1 = smpi_group_rank(group1, proc2);
317 if(proc1 == MPI_UNDEFINED) {
322 *newgroup = MPI_GROUP_EMPTY;
324 *newgroup = smpi_group_new(size);
325 size2 = smpi_group_size(group1);
326 for(i = 0; i < size2; i++) {
327 proc1 = smpi_group_index(group1, i);
328 smpi_group_set_mapping(*newgroup, proc1, i);
330 for(i = size2; i < size; i++) {
331 proc2 = smpi_group_index(group2, i - size2);
332 smpi_group_set_mapping(*newgroup, proc2, i);
335 smpi_group_use(*newgroup);
336 retval = MPI_SUCCESS;
338 smpi_bench_begin(-1, NULL);
342 int MPI_Group_intersection(MPI_Group group1, MPI_Group group2, MPI_Group* newgroup) {
343 int retval, i, proc1, proc2, size, size2;
345 smpi_bench_end(-1, NULL);
346 if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
347 retval = MPI_ERR_GROUP;
348 } else if(newgroup == NULL) {
349 retval = MPI_ERR_ARG;
351 size = smpi_group_size(group1);
352 size2 = smpi_group_size(group2);
353 for(i = 0; i < size2; i++) {
354 proc2 = smpi_group_index(group2, i);
355 proc1 = smpi_group_rank(group1, proc2);
356 if(proc1 == MPI_UNDEFINED) {
361 *newgroup = MPI_GROUP_EMPTY;
363 *newgroup = smpi_group_new(size);
364 size2 = smpi_group_size(group1);
365 for(i = 0; i < size2; i++) {
366 proc1 = smpi_group_index(group1, i);
367 proc2 = smpi_group_rank(group2, proc1);
368 if(proc2 != MPI_UNDEFINED) {
369 smpi_group_set_mapping(*newgroup, proc1, i);
373 smpi_group_use(*newgroup);
374 retval = MPI_SUCCESS;
376 smpi_bench_begin(-1, NULL);
380 int MPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group* newgroup) {
381 int retval, i, proc1, proc2, size, size2;
383 smpi_bench_end(-1, NULL);
384 if(group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
385 retval = MPI_ERR_GROUP;
386 } else if(newgroup == NULL) {
387 retval = MPI_ERR_ARG;
389 size = size2 = smpi_group_size(group1);
390 for(i = 0; i < size2; i++) {
391 proc1 = smpi_group_index(group1, i);
392 proc2 = smpi_group_rank(group2, proc1);
393 if(proc2 != MPI_UNDEFINED) {
398 *newgroup = MPI_GROUP_EMPTY;
400 *newgroup = smpi_group_new(size);
401 for(i = 0; i < size2; i++) {
402 proc1 = smpi_group_index(group1, i);
403 proc2 = smpi_group_rank(group2, proc1);
404 if(proc2 == MPI_UNDEFINED) {
405 smpi_group_set_mapping(*newgroup, proc1, i);
409 smpi_group_use(*newgroup);
410 retval = MPI_SUCCESS;
412 smpi_bench_begin(-1, NULL);
416 int MPI_Group_incl(MPI_Group group, int n, int* ranks, MPI_Group* newgroup) {
417 int retval, i, index;
419 smpi_bench_end(-1, NULL);
420 if(group == MPI_GROUP_NULL) {
421 retval = MPI_ERR_GROUP;
422 } else if(newgroup == NULL) {
423 retval = MPI_ERR_ARG;
426 *newgroup = MPI_GROUP_EMPTY;
427 } else if(n == smpi_group_size(group)) {
430 *newgroup = smpi_group_new(n);
431 for(i = 0; i < n; i++) {
432 index = smpi_group_index(group, ranks[i]);
433 smpi_group_set_mapping(*newgroup, index, i);
436 smpi_group_use(*newgroup);
437 retval = MPI_SUCCESS;
439 smpi_bench_begin(-1, NULL);
443 int MPI_Group_excl(MPI_Group group, int n, int* ranks, MPI_Group* newgroup) {
444 int retval, i, size, rank, index;
446 smpi_bench_end(-1, NULL);
447 if(group == MPI_GROUP_NULL) {
448 retval = MPI_ERR_GROUP;
449 } else if(newgroup == NULL) {
450 retval = MPI_ERR_ARG;
454 } else if(n == smpi_group_size(group)) {
455 *newgroup = MPI_GROUP_EMPTY;
457 size = smpi_group_size(group) - n;
458 *newgroup = smpi_group_new(size);
461 for(i = 0; i < n; i++) {
462 if(ranks[i] == rank) {
467 index = smpi_group_index(group, rank);
468 smpi_group_set_mapping(*newgroup, index, rank);
473 smpi_group_use(*newgroup);
474 retval = MPI_SUCCESS;
476 smpi_bench_begin(-1, NULL);
480 int MPI_Group_range_incl(MPI_Group group, int n, int ranges[][3], MPI_Group* newgroup) {
481 int retval, i, j, rank, size, index;
483 smpi_bench_end(-1, NULL);
484 if(group == MPI_GROUP_NULL) {
485 retval = MPI_ERR_GROUP;
486 } else if(newgroup == NULL) {
487 retval = MPI_ERR_ARG;
490 *newgroup = MPI_GROUP_EMPTY;
493 for(i = 0; i < n; i++) {
494 for(rank = ranges[i][0]; /* First */
495 rank >= 0 && rank <= ranges[i][1]; /* Last */
496 rank += ranges[i][2] /* Stride */) {
500 if(size == smpi_group_size(group)) {
503 *newgroup = smpi_group_new(size);
505 for(i = 0; i < n; i++) {
506 for(rank = ranges[i][0]; /* First */
507 rank >= 0 && rank <= ranges[i][1]; /* Last */
508 rank += ranges[i][2] /* Stride */) {
509 index = smpi_group_index(group, rank);
510 smpi_group_set_mapping(*newgroup, index, j);
516 smpi_group_use(*newgroup);
517 retval = MPI_SUCCESS;
519 smpi_bench_begin(-1, NULL);
523 int MPI_Group_range_excl(MPI_Group group, int n, int ranges[][3], MPI_Group* newgroup) {
524 int retval, i, newrank, rank, size, index, add;
526 smpi_bench_end(-1, NULL);
527 if(group == MPI_GROUP_NULL) {
528 retval = MPI_ERR_GROUP;
529 } else if(newgroup == NULL) {
530 retval = MPI_ERR_ARG;
535 size = smpi_group_size(group);
536 for(i = 0; i < n; i++) {
537 for(rank = ranges[i][0]; /* First */
538 rank >= 0 && rank <= ranges[i][1]; /* Last */
539 rank += ranges[i][2] /* Stride */) {
544 *newgroup = MPI_GROUP_EMPTY;
546 *newgroup = smpi_group_new(size);
548 while(newrank < size) {
549 for(i = 0; i < n; i++) {
551 for(rank = ranges[i][0]; /* First */
552 rank >= 0 && rank <= ranges[i][1]; /* Last */
553 rank += ranges[i][2] /* Stride */) {
554 if(rank == newrank) {
560 index = smpi_group_index(group, newrank);
561 smpi_group_set_mapping(*newgroup, index, newrank);
567 smpi_group_use(*newgroup);
568 retval = MPI_SUCCESS;
570 smpi_bench_begin(-1, NULL);
574 int MPI_Comm_rank(MPI_Comm comm, int* rank) {
577 smpi_bench_end(-1, NULL);
578 if(comm == MPI_COMM_NULL) {
579 retval = MPI_ERR_COMM;
581 *rank = smpi_comm_rank(comm);
582 retval = MPI_SUCCESS;
584 smpi_bench_begin(-1, NULL);
588 int MPI_Comm_size(MPI_Comm comm, int* size) {
591 smpi_bench_end(-1, NULL);
592 if(comm == MPI_COMM_NULL) {
593 retval = MPI_ERR_COMM;
594 } else if(size == NULL) {
595 retval = MPI_ERR_ARG;
597 *size = smpi_comm_size(comm);
598 retval = MPI_SUCCESS;
600 smpi_bench_begin(-1, NULL);
604 int MPI_Comm_group(MPI_Comm comm, MPI_Group* group) {
607 smpi_bench_end(-1, NULL);
608 if(comm == MPI_COMM_NULL) {
609 retval = MPI_ERR_COMM;
610 } else if(group == NULL) {
611 retval = MPI_ERR_ARG;
613 *group = smpi_comm_group(comm);
614 retval = MPI_SUCCESS;
616 smpi_bench_begin(-1, NULL);
620 int MPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int* result) {
623 smpi_bench_end(-1, NULL);
624 if(comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
625 retval = MPI_ERR_COMM;
626 } else if(result == NULL) {
627 retval = MPI_ERR_ARG;
629 if(comm1 == comm2) { /* Same communicators means same groups */
632 *result = smpi_group_compare(smpi_comm_group(comm1), smpi_comm_group(comm2));
633 if(*result == MPI_IDENT) {
634 *result = MPI_CONGRUENT;
637 retval = MPI_SUCCESS;
639 smpi_bench_begin(-1, NULL);
643 int MPI_Comm_dup(MPI_Comm comm, MPI_Comm* newcomm) {
646 smpi_bench_end(-1, NULL);
647 if(comm == MPI_COMM_NULL) {
648 retval = MPI_ERR_COMM;
649 } else if(newcomm == NULL) {
650 retval = MPI_ERR_ARG;
652 *newcomm = smpi_comm_new(smpi_comm_group(comm));
653 retval = MPI_SUCCESS;
655 smpi_bench_begin(-1, NULL);
659 int MPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm* newcomm) {
662 smpi_bench_end(-1, NULL);
663 if(comm == MPI_COMM_NULL) {
664 retval = MPI_ERR_COMM;
665 } else if(group == MPI_GROUP_NULL) {
666 retval = MPI_ERR_GROUP;
667 } else if(newcomm == NULL) {
668 retval = MPI_ERR_ARG;
670 *newcomm = smpi_comm_new(group);
671 retval = MPI_SUCCESS;
673 smpi_bench_begin(-1, NULL);
677 int MPI_Comm_free(MPI_Comm* comm) {
680 smpi_bench_end(-1, NULL);
682 retval = MPI_ERR_ARG;
683 } else if(*comm == MPI_COMM_NULL) {
684 retval = MPI_ERR_COMM;
686 smpi_comm_destroy(*comm);
687 *comm = MPI_COMM_NULL;
688 retval = MPI_SUCCESS;
690 smpi_bench_begin(-1, NULL);
694 int MPI_Send_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
696 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
698 smpi_bench_end(rank, "Send_init");
699 if(request == NULL) {
700 retval = MPI_ERR_ARG;
701 } else if (comm == MPI_COMM_NULL) {
702 retval = MPI_ERR_COMM;
704 *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
705 retval = MPI_SUCCESS;
707 smpi_bench_begin(rank, "Send_init");
711 int MPI_Recv_init(void* buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request* request) {
713 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
715 smpi_bench_end(rank, "Recv_init");
716 if(request == NULL) {
717 retval = MPI_ERR_ARG;
718 } else if (comm == MPI_COMM_NULL) {
719 retval = MPI_ERR_COMM;
721 *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
722 retval = MPI_SUCCESS;
724 smpi_bench_begin(rank, "Recv_init");
728 int MPI_Start(MPI_Request* request) {
730 MPI_Comm comm = (*request)->comm;
731 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
733 smpi_bench_end(rank, "Start");
734 if(request == NULL) {
735 retval = MPI_ERR_ARG;
737 smpi_mpi_start(*request);
738 retval = MPI_SUCCESS;
740 smpi_bench_begin(rank, "Start");
744 int MPI_Startall(int count, MPI_Request* requests) {
746 MPI_Comm comm = count > 0 && requests ? requests[0]->comm : MPI_COMM_NULL;
747 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
749 smpi_bench_end(rank, "Startall");
750 if(requests == NULL) {
751 retval = MPI_ERR_ARG;
753 smpi_mpi_startall(count, requests);
754 retval = MPI_SUCCESS;
756 smpi_bench_begin(rank, "Startall");
760 int MPI_Request_free(MPI_Request* request) {
762 MPI_Comm comm = (*request)->comm;
763 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
765 smpi_bench_end(rank, "Request_free");
766 if(request == NULL) {
767 retval = MPI_ERR_ARG;
769 smpi_mpi_request_free(request);
770 retval = MPI_SUCCESS;
772 smpi_bench_begin(rank, "Request_free");
776 int MPI_Irecv(void* buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request* request) {
778 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
780 smpi_bench_end(rank, "Irecv");
782 int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
783 TRACE_smpi_ptp_in (rank, src_traced, rank, __FUNCTION__);
785 if(request == NULL) {
786 retval = MPI_ERR_ARG;
787 } else if (comm == MPI_COMM_NULL) {
788 retval = MPI_ERR_COMM;
790 *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
791 retval = MPI_SUCCESS;
794 TRACE_smpi_ptp_out (rank, src_traced, rank, __FUNCTION__);
795 (*request)->recv = 1;
797 smpi_bench_begin(rank, "Irecv");
801 int MPI_Isend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
803 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
805 smpi_bench_end(rank, "Isend");
807 int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
808 TRACE_smpi_ptp_in (rank, rank, dst_traced, __FUNCTION__);
809 TRACE_smpi_send (rank, rank, dst_traced);
811 if(request == NULL) {
812 retval = MPI_ERR_ARG;
813 } else if (comm == MPI_COMM_NULL) {
814 retval = MPI_ERR_COMM;
816 *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
817 retval = MPI_SUCCESS;
820 TRACE_smpi_ptp_out (rank, rank, dst_traced, __FUNCTION__);
821 (*request)->send = 1;
823 smpi_bench_begin(rank, "Isend");
827 int MPI_Recv(void* buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status* status) {
829 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
831 smpi_bench_end(rank, "Recv");
833 int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
834 TRACE_smpi_ptp_in (rank, src_traced, rank, __FUNCTION__);
836 if (comm == MPI_COMM_NULL) {
837 retval = MPI_ERR_COMM;
839 smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
840 retval = MPI_SUCCESS;
843 TRACE_smpi_ptp_out (rank, src_traced, rank, __FUNCTION__);
844 TRACE_smpi_recv (rank, src_traced, rank);
846 smpi_bench_begin(rank, "Recv");
850 int MPI_Send(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
852 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
854 smpi_bench_end(rank, "Send");
856 int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
857 TRACE_smpi_ptp_in (rank, rank, dst_traced, __FUNCTION__);
858 TRACE_smpi_send (rank, rank, dst_traced);
860 if (comm == MPI_COMM_NULL) {
861 retval = MPI_ERR_COMM;
863 smpi_mpi_send(buf, count, datatype, dst, tag, comm);
864 retval = MPI_SUCCESS;
867 TRACE_smpi_ptp_out (rank, rank, dst_traced, __FUNCTION__);
869 smpi_bench_begin(rank, "Send");
873 int MPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void* recvbuf, int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status* status) {
875 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
877 smpi_bench_end(rank, "Sendrecv");
879 int dst_traced = smpi_group_rank(smpi_comm_group(comm), dst);
880 int src_traced = smpi_group_rank(smpi_comm_group(comm), src);
881 TRACE_smpi_ptp_in (rank, src_traced, dst_traced, __FUNCTION__);
882 TRACE_smpi_send (rank, rank, dst_traced);
883 TRACE_smpi_send (rank, src_traced, rank);
885 if (comm == MPI_COMM_NULL) {
886 retval = MPI_ERR_COMM;
887 } else if (sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
888 retval = MPI_ERR_TYPE;
890 smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf, recvcount, recvtype, src, recvtag, comm, status);
891 retval = MPI_SUCCESS;
894 TRACE_smpi_ptp_out (rank, src_traced, dst_traced, __FUNCTION__);
895 TRACE_smpi_recv (rank, rank, dst_traced);
896 TRACE_smpi_recv (rank, src_traced, rank);
898 smpi_bench_begin(rank, "Sendrecv");
902 int MPI_Sendrecv_replace(void* buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag, MPI_Comm comm, MPI_Status* status) {
903 //TODO: suboptimal implementation
907 size = smpi_datatype_size(datatype) * count;
908 recvbuf = xbt_new(char, size);
909 retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
910 memcpy(buf, recvbuf, size * sizeof(char));
915 int MPI_Test(MPI_Request* request, int* flag, MPI_Status* status) {
917 int rank = request && (*request)->comm != MPI_COMM_NULL
918 ? smpi_comm_rank((*request)->comm)
921 smpi_bench_end(rank, "Test");
922 if(request == NULL || flag == NULL) {
923 retval = MPI_ERR_ARG;
924 } else if(*request == MPI_REQUEST_NULL) {
925 retval = MPI_ERR_REQUEST;
927 *flag = smpi_mpi_test(request, status);
928 retval = MPI_SUCCESS;
930 smpi_bench_begin(rank, "Test");
934 int MPI_Testany(int count, MPI_Request requests[], int* index, int* flag, MPI_Status* status) {
937 smpi_bench_end(-1, NULL); //FIXME
938 if(index == NULL || flag == NULL) {
939 retval = MPI_ERR_ARG;
941 *flag = smpi_mpi_testany(count, requests, index, status);
942 retval = MPI_SUCCESS;
944 smpi_bench_begin(-1, NULL);
948 int MPI_Wait(MPI_Request* request, MPI_Status* status) {
950 int rank = request && (*request)->comm != MPI_COMM_NULL
951 ? smpi_comm_rank((*request)->comm)
954 smpi_bench_end(rank, "Wait");
956 MPI_Group group = smpi_comm_group((*request)->comm);
957 int src_traced = smpi_group_rank (group , (*request)->src);
958 int dst_traced = smpi_group_rank (group , (*request)->dst);
959 int is_wait_for_receive = (*request)->recv;
960 TRACE_smpi_ptp_in (rank, src_traced, dst_traced, __FUNCTION__);
962 if(request == NULL) {
963 retval = MPI_ERR_ARG;
964 } else if(*request == MPI_REQUEST_NULL) {
965 retval = MPI_ERR_REQUEST;
967 smpi_mpi_wait(request, status);
968 retval = MPI_SUCCESS;
971 TRACE_smpi_ptp_out (rank, src_traced, dst_traced, __FUNCTION__);
972 if (is_wait_for_receive){
973 TRACE_smpi_recv (rank, src_traced, dst_traced);
976 smpi_bench_begin(rank, "Wait");
980 int MPI_Waitany(int count, MPI_Request requests[], int* index, MPI_Status* status) {
983 smpi_bench_end(-1, NULL); //FIXME
985 //save requests information for tracing
987 xbt_dynar_t srcs = xbt_dynar_new (sizeof(int), xbt_free);
988 xbt_dynar_t dsts = xbt_dynar_new (sizeof(int), xbt_free);
989 xbt_dynar_t recvs = xbt_dynar_new (sizeof(int), xbt_free);
990 for (i = 0; i < count; i++){
991 MPI_Request req = requests[i]; //already received requests are no longer valid
993 int *asrc = xbt_new(int, 1);
994 int *adst = xbt_new(int, 1);
995 int *arecv = xbt_new(int, 1);
999 xbt_dynar_insert_at (srcs, i, asrc);
1000 xbt_dynar_insert_at (dsts, i, adst);
1001 xbt_dynar_insert_at (recvs, i, arecv);
1003 int *t = xbt_new(int, 1);
1004 xbt_dynar_insert_at (srcs, i, t);
1005 xbt_dynar_insert_at (dsts, i, t);
1006 xbt_dynar_insert_at (recvs, i, t);
1010 //search for a suitable request to give the rank of current mpi proc
1011 MPI_Request req = NULL;
1012 for (i = 0; i < count && req == NULL; i++) {
1015 MPI_Comm comm = (req)->comm;
1016 int rank_traced = smpi_comm_rank(comm);
1017 TRACE_smpi_ptp_in (rank_traced, -1, -1, __FUNCTION__);
1020 retval = MPI_ERR_ARG;
1022 *index = smpi_mpi_waitany(count, requests, status);
1023 retval = MPI_SUCCESS;
1026 int src_traced, dst_traced, is_wait_for_receive;
1027 xbt_dynar_get_cpy (srcs, *index, &src_traced);
1028 xbt_dynar_get_cpy (dsts, *index, &dst_traced);
1029 xbt_dynar_get_cpy (recvs, *index, &is_wait_for_receive);
1030 if (is_wait_for_receive){
1031 TRACE_smpi_recv (rank_traced, src_traced, dst_traced);
1033 TRACE_smpi_ptp_out (rank_traced, src_traced, dst_traced, __FUNCTION__);
1034 //clean-up of dynars
1039 smpi_bench_begin(-1, NULL);
1043 int MPI_Waitall(int count, MPI_Request requests[], MPI_Status status[]) {
1045 smpi_bench_end(-1, NULL); //FIXME
1047 //save information from requests
1049 xbt_dynar_t srcs = xbt_dynar_new (sizeof(int), xbt_free);
1050 xbt_dynar_t dsts = xbt_dynar_new (sizeof(int), xbt_free);
1051 xbt_dynar_t recvs = xbt_dynar_new (sizeof(int), xbt_free);
1052 for (i = 0; i < count; i++){
1053 MPI_Request req = requests[i]; //all req should be valid in Waitall
1054 int *asrc = xbt_new(int, 1);
1055 int *adst = xbt_new(int, 1);
1056 int *arecv = xbt_new(int, 1);
1060 xbt_dynar_insert_at (srcs, i, asrc);
1061 xbt_dynar_insert_at (dsts, i, adst);
1062 xbt_dynar_insert_at (recvs, i, arecv);
1065 // find my rank inside one of MPI_Comm's of the requests
1066 MPI_Request req = NULL;
1067 for (i = 0; i < count && req == NULL; i++) {
1070 MPI_Comm comm = (req)->comm;
1071 int rank_traced = smpi_comm_rank(comm);
1072 TRACE_smpi_ptp_in (rank_traced, -1, -1, __FUNCTION__);
1074 smpi_mpi_waitall(count, requests, status);
1076 for (i = 0; i < count; i++){
1077 int src_traced, dst_traced, is_wait_for_receive;
1078 xbt_dynar_get_cpy (srcs, i, &src_traced);
1079 xbt_dynar_get_cpy (dsts, i, &dst_traced);
1080 xbt_dynar_get_cpy (recvs, i, &is_wait_for_receive);
1081 if (is_wait_for_receive){
1082 TRACE_smpi_recv (rank_traced, src_traced, dst_traced);
1085 TRACE_smpi_ptp_out (rank_traced, -1, -1, __FUNCTION__);
1086 //clean-up of dynars
1091 smpi_bench_begin(-1, NULL);
1095 int MPI_Waitsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[]) {
1098 smpi_bench_end(-1, NULL); //FIXME
1099 if(outcount == NULL || indices == NULL) {
1100 retval = MPI_ERR_ARG;
1102 *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1103 retval = MPI_SUCCESS;
1105 smpi_bench_begin(-1, NULL);
1109 int MPI_Bcast(void* buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm) {
1111 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1113 smpi_bench_end(rank, "Bcast");
1115 int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1116 TRACE_smpi_collective_in (rank, root_traced, __FUNCTION__);
1118 if(comm == MPI_COMM_NULL) {
1119 retval = MPI_ERR_COMM;
1121 smpi_mpi_bcast(buf, count, datatype, root, comm);
1122 retval = MPI_SUCCESS;
1125 TRACE_smpi_collective_out (rank, root_traced, __FUNCTION__);
1127 smpi_bench_begin(rank, "Bcast");
1131 int MPI_Barrier(MPI_Comm comm) {
1133 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1135 smpi_bench_end(rank, "Barrier");
1137 TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
1139 if(comm == MPI_COMM_NULL) {
1140 retval = MPI_ERR_COMM;
1142 smpi_mpi_barrier(comm);
1143 retval = MPI_SUCCESS;
1146 TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
1148 smpi_bench_begin(rank, "Barrier");
1152 int MPI_Gather(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm) {
1154 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1156 smpi_bench_end(rank, "Gather");
1158 int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1159 TRACE_smpi_collective_in (rank, root_traced, __FUNCTION__);
1161 if(comm == MPI_COMM_NULL) {
1162 retval = MPI_ERR_COMM;
1163 } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
1164 retval = MPI_ERR_TYPE;
1166 smpi_mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
1167 retval = MPI_SUCCESS;
1170 TRACE_smpi_collective_out (rank, root_traced, __FUNCTION__);
1172 smpi_bench_begin(rank, "Gather");
1176 int MPI_Gatherv(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int* recvcounts, int* displs, MPI_Datatype recvtype, int root, MPI_Comm comm) {
1178 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1180 smpi_bench_end(rank, "Gatherv");
1182 int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1183 TRACE_smpi_collective_in (rank, root_traced, __FUNCTION__);
1185 if(comm == MPI_COMM_NULL) {
1186 retval = MPI_ERR_COMM;
1187 } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
1188 retval = MPI_ERR_TYPE;
1189 } else if(recvcounts == NULL || displs == NULL) {
1190 retval = MPI_ERR_ARG;
1192 smpi_mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, root, comm);
1193 retval = MPI_SUCCESS;
1196 TRACE_smpi_collective_out (rank, root_traced, __FUNCTION__);
1198 smpi_bench_begin(rank, "Gatherv");
1202 int MPI_Allgather(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm) {
1204 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1206 smpi_bench_end(rank, "Allgather");
1208 TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
1210 if(comm == MPI_COMM_NULL) {
1211 retval = MPI_ERR_COMM;
1212 } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
1213 retval = MPI_ERR_TYPE;
1215 smpi_mpi_allgather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1216 retval = MPI_SUCCESS;
1219 TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
1221 smpi_bench_begin(rank, "Allgather");
1225 int MPI_Allgatherv(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int* recvcounts, int* displs, MPI_Datatype recvtype, MPI_Comm comm) {
1227 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1229 smpi_bench_end(rank, "Allgatherv");
1231 TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
1233 if(comm == MPI_COMM_NULL) {
1234 retval = MPI_ERR_COMM;
1235 } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
1236 retval = MPI_ERR_TYPE;
1237 } else if(recvcounts == NULL || displs == NULL) {
1238 retval = MPI_ERR_ARG;
1240 smpi_mpi_allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
1241 retval = MPI_SUCCESS;
1244 TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
1246 smpi_bench_begin(rank, "Allgatherv");
1250 int MPI_Scatter(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm) {
1252 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1254 smpi_bench_end(rank, "Scatter");
1256 int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1257 TRACE_smpi_collective_in (rank, root_traced, __FUNCTION__);
1259 if(comm == MPI_COMM_NULL) {
1260 retval = MPI_ERR_COMM;
1261 } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
1262 retval = MPI_ERR_TYPE;
1264 smpi_mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
1265 retval = MPI_SUCCESS;
1268 TRACE_smpi_collective_out (rank, root_traced, __FUNCTION__);
1270 smpi_bench_begin(rank, "Scatter");
1274 int MPI_Scatterv(void* sendbuf, int* sendcounts, int* displs, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm) {
1276 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1278 smpi_bench_end(rank, "Scatterv");
1280 int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1281 TRACE_smpi_collective_in (rank, root_traced, __FUNCTION__);
1283 if(comm == MPI_COMM_NULL) {
1284 retval = MPI_ERR_COMM;
1285 } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
1286 retval = MPI_ERR_TYPE;
1287 } else if(sendcounts == NULL || displs == NULL) {
1288 retval = MPI_ERR_ARG;
1290 smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
1291 retval = MPI_SUCCESS;
1294 TRACE_smpi_collective_out (rank, root_traced, __FUNCTION__);
1296 smpi_bench_begin(rank, "Scatterv");
1300 int MPI_Reduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm) {
1302 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1304 smpi_bench_end(rank, "Reduce");
1306 int root_traced = smpi_group_rank(smpi_comm_group(comm), root);
1307 TRACE_smpi_collective_in (rank, root_traced, __FUNCTION__);
1309 if(comm == MPI_COMM_NULL) {
1310 retval = MPI_ERR_COMM;
1311 } else if(datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1312 retval = MPI_ERR_ARG;
1314 smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
1315 retval = MPI_SUCCESS;
1318 TRACE_smpi_collective_out (rank, root_traced, __FUNCTION__);
1320 smpi_bench_begin(rank, "Reduce");
1324 int MPI_Allreduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) {
1326 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1328 smpi_bench_end(rank, "Allreduce");
1330 TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
1332 if(comm == MPI_COMM_NULL) {
1333 retval = MPI_ERR_COMM;
1334 } else if(datatype == MPI_DATATYPE_NULL) {
1335 retval = MPI_ERR_TYPE;
1336 } else if(op == MPI_OP_NULL) {
1337 retval = MPI_ERR_OP;
1339 smpi_mpi_allreduce(sendbuf, recvbuf, count, datatype, op, comm);
1340 retval = MPI_SUCCESS;
1343 TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
1345 smpi_bench_begin(rank, "Allreduce");
1349 int MPI_Scan(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) {
1351 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1353 smpi_bench_end(rank, "Scan");
1355 TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
1357 if(comm == MPI_COMM_NULL) {
1358 retval = MPI_ERR_COMM;
1359 } else if(datatype == MPI_DATATYPE_NULL) {
1360 retval = MPI_ERR_TYPE;
1361 } else if(op == MPI_OP_NULL) {
1362 retval = MPI_ERR_OP;
1364 smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
1365 retval = MPI_SUCCESS;
1368 TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
1370 smpi_bench_begin(rank, "Scan");
1374 int MPI_Reduce_scatter(void* sendbuf, void* recvbuf, int* recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) {
1375 int retval, i, size, count;
1377 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1379 smpi_bench_end(rank, "Reduce_scatter");
1381 TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
1383 if(comm == MPI_COMM_NULL) {
1384 retval = MPI_ERR_COMM;
1385 } else if(datatype == MPI_DATATYPE_NULL) {
1386 retval = MPI_ERR_TYPE;
1387 } else if(op == MPI_OP_NULL) {
1388 retval = MPI_ERR_OP;
1389 } else if(recvcounts == NULL) {
1390 retval = MPI_ERR_ARG;
1392 /* arbitrarily choose root as rank 0 */
1393 /* TODO: faster direct implementation ? */
1394 size = smpi_comm_size(comm);
1396 displs = xbt_new(int, size);
1397 for(i = 0; i < size; i++) {
1398 count += recvcounts[i];
1401 smpi_mpi_reduce(sendbuf, recvbuf, count, datatype, op, 0, comm);
1402 smpi_mpi_scatterv(recvbuf, recvcounts, displs, datatype, recvbuf, recvcounts[rank], datatype, 0, comm);
1404 retval = MPI_SUCCESS;
1407 TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
1409 smpi_bench_begin(rank, "Reduce_scatter");
1413 int MPI_Alltoall(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm) {
1414 int retval, size, sendsize;
1415 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1417 smpi_bench_end(rank, "Alltoall");
1419 TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
1421 if(comm == MPI_COMM_NULL) {
1422 retval = MPI_ERR_COMM;
1423 } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
1424 retval = MPI_ERR_TYPE;
1426 size = smpi_comm_size(comm);
1427 sendsize = smpi_datatype_size(sendtype) * sendcount;
1428 if(sendsize < 200 && size > 12) {
1429 retval = smpi_coll_tuned_alltoall_bruck(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1430 } else if(sendsize < 3000) {
1431 retval = smpi_coll_tuned_alltoall_basic_linear(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1433 retval = smpi_coll_tuned_alltoall_pairwise(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1437 TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
1439 smpi_bench_begin(rank, "Alltoall");
1443 int MPI_Alltoallv(void* sendbuf, int* sendcounts, int* senddisps, MPI_Datatype sendtype, void* recvbuf, int *recvcounts, int* recvdisps, MPI_Datatype recvtype, MPI_Comm comm) {
1445 int rank = comm != MPI_COMM_NULL ? smpi_comm_rank(comm) : -1;
1447 smpi_bench_end(rank, "Alltoallv");
1449 TRACE_smpi_collective_in (rank, -1, __FUNCTION__);
1451 if(comm == MPI_COMM_NULL) {
1452 retval = MPI_ERR_COMM;
1453 } else if(sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
1454 retval = MPI_ERR_TYPE;
1455 } else if(sendcounts == NULL || senddisps == NULL || recvcounts == NULL || recvdisps == NULL) {
1456 retval = MPI_ERR_ARG;
1458 retval = smpi_coll_basic_alltoallv(sendbuf, sendcounts, senddisps, sendtype, recvbuf, recvcounts, recvdisps, recvtype, comm);
1461 TRACE_smpi_collective_out (rank, -1, __FUNCTION__);
1463 smpi_bench_begin(rank, "Alltoallv");
1468 int MPI_Get_processor_name(char* name, int* resultlen) {
1469 int retval = MPI_SUCCESS;
1471 smpi_bench_end(-1, NULL);
1472 strncpy( name , SIMIX_host_get_name(SIMIX_host_self()), MPI_MAX_PROCESSOR_NAME-1);
1473 *resultlen= strlen(name) > MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
1475 smpi_bench_begin(-1, NULL);
1479 int MPI_Get_count(MPI_Status* status, MPI_Datatype datatype, int* count) {
1480 int retval = MPI_SUCCESS;
1483 smpi_bench_end(-1, NULL);
1484 if (status == NULL || count == NULL) {
1485 retval = MPI_ERR_ARG;
1486 } else if (datatype == MPI_DATATYPE_NULL) {
1487 retval = MPI_ERR_TYPE;
1489 size = smpi_datatype_size(datatype);
1492 } else if (status->count % size != 0) {
1493 retval = MPI_UNDEFINED;
1495 *count = smpi_mpi_get_count(status, datatype);
1498 smpi_bench_begin(-1, NULL);