2 /* Copyright (c) 2007-2015. The SimGrid Team.
3 * All rights reserved. */
5 /* This program is free software; you can redistribute it and/or modify it
6 * under the terms of the license (GNU LGPL) which comes with this package. */
9 #include "smpi_mpi_dt_private.h"
11 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_pmpi, smpi, "Logging specific to SMPI (pmpi)");
13 //this function need to be here because of the calls to smpi_bench
14 void TRACE_smpi_set_category(const char *category)
16 //need to end bench otherwise categories for execution tasks are wrong
18 TRACE_internal_smpi_set_category (category);
19 //begin bench after changing process's category
23 /* PMPI User level calls */
25 int PMPI_Init(int *argc, char ***argv)
27 // PMPI_Init is call only one time by only by SMPI process
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__);
45 int PMPI_Finalize(void)
48 int rank = smpi_process_index();
49 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
50 extra->type = TRACING_FINALIZE;
51 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
53 smpi_process_finalize();
55 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
56 TRACE_smpi_finalize(smpi_process_index());
57 smpi_process_destroy();
61 int PMPI_Finalized(int* flag)
63 *flag=smpi_process_finalized();
67 int PMPI_Get_version (int *version,int *subversion){
68 *version = MPI_VERSION;
69 *subversion= MPI_SUBVERSION;
73 int PMPI_Get_library_version (char *version,int *len){
74 int retval = MPI_SUCCESS;
76 snprintf(version,MPI_MAX_LIBRARY_VERSION_STRING,"SMPI Version %d.%d. Copyright The Simgrid Team 2007-2015",
77 SIMGRID_VERSION_MAJOR, SIMGRID_VERSION_MINOR);
78 *len = strlen(version) > MPI_MAX_LIBRARY_VERSION_STRING ? MPI_MAX_LIBRARY_VERSION_STRING : strlen(version);
83 int PMPI_Init_thread(int *argc, char ***argv, int required, int *provided)
85 if (provided != NULL) {
86 *provided = MPI_THREAD_SINGLE;
88 return MPI_Init(argc, argv);
91 int PMPI_Query_thread(int *provided)
95 if (provided == NULL) {
98 *provided = MPI_THREAD_SINGLE;
104 int PMPI_Is_thread_main(int *flag)
109 retval = MPI_ERR_ARG;
111 *flag = smpi_process_index() == 0;
112 retval = MPI_SUCCESS;
117 int PMPI_Abort(MPI_Comm comm, int errorcode)
120 smpi_process_destroy();
121 // FIXME: should kill all processes in comm instead
122 simcall_process_kill(SIMIX_process_self());
126 double PMPI_Wtime(void)
128 return smpi_mpi_wtime();
131 extern double sg_maxmin_precision;
132 double PMPI_Wtick(void)
134 return sg_maxmin_precision;
137 int PMPI_Address(void *location, MPI_Aint * address)
142 retval = MPI_ERR_ARG;
144 *address = reinterpret_cast<MPI_Aint>(location);
145 retval = MPI_SUCCESS;
150 int PMPI_Get_address(void *location, MPI_Aint * address)
152 return PMPI_Address(location, address);
155 int PMPI_Type_free(MPI_Datatype * datatype)
158 /* Free a predefined datatype is an error according to the standard, and should be checked for */
159 if (*datatype == MPI_DATATYPE_NULL) {
160 retval = MPI_ERR_ARG;
162 smpi_datatype_unuse(*datatype);
163 retval = MPI_SUCCESS;
168 int PMPI_Type_size(MPI_Datatype datatype, int *size)
172 if (datatype == MPI_DATATYPE_NULL) {
173 retval = MPI_ERR_TYPE;
174 } else if (size == NULL) {
175 retval = MPI_ERR_ARG;
177 *size = static_cast<int>(smpi_datatype_size(datatype));
178 retval = MPI_SUCCESS;
183 int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
187 if (datatype == MPI_DATATYPE_NULL) {
188 retval = MPI_ERR_TYPE;
189 } else if (lb == NULL || extent == NULL) {
190 retval = MPI_ERR_ARG;
192 retval = smpi_datatype_extent(datatype, lb, extent);
197 int PMPI_Type_get_true_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
199 return PMPI_Type_get_extent(datatype, lb, extent);
202 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
206 if (datatype == MPI_DATATYPE_NULL) {
207 retval = MPI_ERR_TYPE;
208 } else if (extent == NULL) {
209 retval = MPI_ERR_ARG;
211 *extent = smpi_datatype_get_extent(datatype);
212 retval = MPI_SUCCESS;
217 int PMPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
221 if (datatype == MPI_DATATYPE_NULL) {
222 retval = MPI_ERR_TYPE;
223 } else if (disp == NULL) {
224 retval = MPI_ERR_ARG;
226 *disp = smpi_datatype_lb(datatype);
227 retval = MPI_SUCCESS;
232 int PMPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
236 if (datatype == MPI_DATATYPE_NULL) {
237 retval = MPI_ERR_TYPE;
238 } else if (disp == NULL) {
239 retval = MPI_ERR_ARG;
241 *disp = smpi_datatype_ub(datatype);
242 retval = MPI_SUCCESS;
247 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
250 if (datatype == MPI_DATATYPE_NULL) {
251 retval = MPI_ERR_TYPE;
253 retval = smpi_datatype_dup(datatype, newtype);
258 int PMPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
262 if (function == NULL || op == NULL) {
263 retval = MPI_ERR_ARG;
265 *op = smpi_op_new(function, commute);
266 retval = MPI_SUCCESS;
271 int PMPI_Op_free(MPI_Op * op)
276 retval = MPI_ERR_ARG;
277 } else if (*op == MPI_OP_NULL) {
280 smpi_op_destroy(*op);
282 retval = MPI_SUCCESS;
287 int PMPI_Group_free(MPI_Group * group)
292 retval = MPI_ERR_ARG;
294 smpi_group_destroy(*group);
295 *group = MPI_GROUP_NULL;
296 retval = MPI_SUCCESS;
301 int PMPI_Group_size(MPI_Group group, int *size)
305 if (group == MPI_GROUP_NULL) {
306 retval = MPI_ERR_GROUP;
307 } else if (size == NULL) {
308 retval = MPI_ERR_ARG;
310 *size = smpi_group_size(group);
311 retval = MPI_SUCCESS;
316 int PMPI_Group_rank(MPI_Group group, int *rank)
320 if (group == MPI_GROUP_NULL) {
321 retval = MPI_ERR_GROUP;
322 } else if (rank == NULL) {
323 retval = MPI_ERR_ARG;
325 *rank = smpi_group_rank(group, smpi_process_index());
326 retval = MPI_SUCCESS;
331 int PMPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1, MPI_Group group2, int *ranks2)
333 int retval, i, index;
334 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
335 retval = MPI_ERR_GROUP;
337 for (i = 0; i < n; i++) {
338 if(ranks1[i]==MPI_PROC_NULL){
339 ranks2[i]=MPI_PROC_NULL;
341 index = smpi_group_index(group1, ranks1[i]);
342 ranks2[i] = smpi_group_rank(group2, index);
345 retval = MPI_SUCCESS;
350 int PMPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
354 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
355 retval = MPI_ERR_GROUP;
356 } else if (result == NULL) {
357 retval = MPI_ERR_ARG;
359 *result = smpi_group_compare(group1, group2);
360 retval = MPI_SUCCESS;
365 int PMPI_Group_union(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
367 int retval, i, proc1, proc2, size, size2;
369 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
370 retval = MPI_ERR_GROUP;
371 } else if (newgroup == NULL) {
372 retval = MPI_ERR_ARG;
374 size = smpi_group_size(group1);
375 size2 = smpi_group_size(group2);
376 for (i = 0; i < size2; i++) {
377 proc2 = smpi_group_index(group2, i);
378 proc1 = smpi_group_rank(group1, proc2);
379 if (proc1 == MPI_UNDEFINED) {
384 *newgroup = MPI_GROUP_EMPTY;
386 *newgroup = smpi_group_new(size);
387 size2 = smpi_group_size(group1);
388 for (i = 0; i < size2; i++) {
389 proc1 = smpi_group_index(group1, i);
390 smpi_group_set_mapping(*newgroup, proc1, i);
392 for (i = size2; i < size; i++) {
393 proc2 = smpi_group_index(group2, i - size2);
394 smpi_group_set_mapping(*newgroup, proc2, i);
397 retval = MPI_SUCCESS;
402 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
404 int retval, i, proc1, proc2, size;
406 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
407 retval = MPI_ERR_GROUP;
408 } else if (newgroup == NULL) {
409 retval = MPI_ERR_ARG;
411 size = smpi_group_size(group2);
412 for (i = 0; i < size; i++) {
413 proc2 = smpi_group_index(group2, i);
414 proc1 = smpi_group_rank(group1, proc2);
415 if (proc1 == MPI_UNDEFINED) {
420 *newgroup = MPI_GROUP_EMPTY;
422 *newgroup = smpi_group_new(size);
424 for (i = 0; i < smpi_group_size(group2); i++) {
425 proc2 = smpi_group_index(group2, i);
426 proc1 = smpi_group_rank(group1, proc2);
427 if (proc1 != MPI_UNDEFINED) {
428 smpi_group_set_mapping(*newgroup, proc2, j);
433 retval = MPI_SUCCESS;
438 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
440 int retval, i, proc1, proc2, size, size2;
442 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
443 retval = MPI_ERR_GROUP;
444 } else if (newgroup == NULL) {
445 retval = MPI_ERR_ARG;
447 size = smpi_group_size(group1);
449 for (i = 0; i < size2; i++) {
450 proc1 = smpi_group_index(group1, i);
451 proc2 = smpi_group_rank(group2, proc1);
452 if (proc2 != MPI_UNDEFINED) {
457 *newgroup = MPI_GROUP_EMPTY;
459 *newgroup = smpi_group_new(size);
460 for (i = 0; i < size2; i++) {
461 proc1 = smpi_group_index(group1, i);
462 proc2 = smpi_group_rank(group2, proc1);
463 if (proc2 == MPI_UNDEFINED) {
464 smpi_group_set_mapping(*newgroup, proc1, i);
468 retval = MPI_SUCCESS;
473 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
477 if (group == MPI_GROUP_NULL) {
478 retval = MPI_ERR_GROUP;
479 } else if (newgroup == NULL) {
480 retval = MPI_ERR_ARG;
482 retval = smpi_group_incl(group, n, ranks, newgroup);
487 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
489 int retval, i, j, newsize, oldsize, index;
491 if (group == MPI_GROUP_NULL) {
492 retval = MPI_ERR_GROUP;
493 } else if (newgroup == NULL) {
494 retval = MPI_ERR_ARG;
498 if(group!= smpi_comm_group(MPI_COMM_WORLD) && group != MPI_GROUP_NULL
499 && group != smpi_comm_group(MPI_COMM_SELF) && group != MPI_GROUP_EMPTY)
500 smpi_group_use(group);
501 } else if (n == smpi_group_size(group)) {
502 *newgroup = MPI_GROUP_EMPTY;
504 oldsize=smpi_group_size(group);
505 newsize = oldsize - n;
506 *newgroup = smpi_group_new(newsize);
508 int* to_exclude=xbt_new0(int, smpi_group_size(group));
509 for(i=0; i<oldsize; i++)
512 to_exclude[ranks[i]]=1;
515 for(i=0; i<oldsize; i++){
516 if(to_exclude[i]==0){
517 index = smpi_group_index(group, i);
518 smpi_group_set_mapping(*newgroup, index, j);
523 xbt_free(to_exclude);
525 retval = MPI_SUCCESS;
530 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3], MPI_Group * newgroup)
532 int retval, i, j, rank, size, index;
534 if (group == MPI_GROUP_NULL) {
535 retval = MPI_ERR_GROUP;
536 } else if (newgroup == NULL) {
537 retval = MPI_ERR_ARG;
540 *newgroup = MPI_GROUP_EMPTY;
543 for (i = 0; i < n; i++) {
544 for (rank = ranges[i][0]; /* First */
545 rank >= 0 && rank < smpi_group_size(group); /* Last */
548 if(rank == ranges[i][1]){/*already last ?*/
551 rank += ranges[i][2]; /* Stride */
552 if (ranges[i][0]<ranges[i][1]){
553 if(rank > ranges[i][1])
556 if(rank < ranges[i][1])
562 *newgroup = smpi_group_new(size);
564 for (i = 0; i < n; i++) {
565 for (rank = ranges[i][0]; /* First */
566 rank >= 0 && rank < smpi_group_size(group); /* Last */
568 index = smpi_group_index(group, rank);
569 smpi_group_set_mapping(*newgroup, index, j);
571 if(rank == ranges[i][1]){/*already last ?*/
574 rank += ranges[i][2]; /* Stride */
575 if (ranges[i][0]<ranges[i][1]){
576 if(rank > ranges[i][1])
579 if(rank < ranges[i][1])
585 retval = MPI_SUCCESS;
590 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3], MPI_Group * newgroup)
592 int retval, i, rank, newrank,oldrank, size, index, add;
594 if (group == MPI_GROUP_NULL) {
595 retval = MPI_ERR_GROUP;
596 } else if (newgroup == NULL) {
597 retval = MPI_ERR_ARG;
601 if(group!= smpi_comm_group(MPI_COMM_WORLD) && group != MPI_GROUP_NULL
602 && group != smpi_comm_group(MPI_COMM_SELF) && group != MPI_GROUP_EMPTY)
603 smpi_group_use(group);
605 size = smpi_group_size(group);
606 for (i = 0; i < n; i++) {
607 for (rank = ranges[i][0]; /* First */
608 rank >= 0 && rank < smpi_group_size(group); /* Last */
611 if(rank == ranges[i][1]){/*already last ?*/
614 rank += ranges[i][2]; /* Stride */
615 if (ranges[i][0]<ranges[i][1]){
616 if(rank > ranges[i][1])
619 if(rank < ranges[i][1])
625 *newgroup = MPI_GROUP_EMPTY;
627 *newgroup = smpi_group_new(size);
630 while (newrank < size) {
632 for (i = 0; i < n; i++) {
633 for (rank = ranges[i][0];
634 rank >= 0 && rank < smpi_group_size(group);
640 if(rank == ranges[i][1]){/*already last ?*/
643 rank += ranges[i][2]; /* Stride */
644 if (ranges[i][0]<ranges[i][1]){
645 if(rank > ranges[i][1])
648 if(rank < ranges[i][1])
654 index = smpi_group_index(group, oldrank);
655 smpi_group_set_mapping(*newgroup, index, newrank);
663 retval = MPI_SUCCESS;
668 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
671 if (comm == MPI_COMM_NULL) {
672 retval = MPI_ERR_COMM;
673 } else if (rank == NULL) {
674 retval = MPI_ERR_ARG;
676 *rank = smpi_comm_rank(comm);
677 retval = MPI_SUCCESS;
682 int PMPI_Comm_size(MPI_Comm comm, int *size)
685 if (comm == MPI_COMM_NULL) {
686 retval = MPI_ERR_COMM;
687 } else if (size == NULL) {
688 retval = MPI_ERR_ARG;
690 *size = smpi_comm_size(comm);
691 retval = MPI_SUCCESS;
696 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
700 if (comm == MPI_COMM_NULL) {
701 retval = MPI_ERR_COMM;
702 } else if (name == NULL || len == NULL) {
703 retval = MPI_ERR_ARG;
705 smpi_comm_get_name(comm, name, len);
706 retval = MPI_SUCCESS;
711 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
715 if (comm == MPI_COMM_NULL) {
716 retval = MPI_ERR_COMM;
717 } else if (group == NULL) {
718 retval = MPI_ERR_ARG;
720 *group = smpi_comm_group(comm);
721 if(*group!= smpi_comm_group(MPI_COMM_WORLD) && *group != MPI_GROUP_NULL
722 && *group != MPI_GROUP_EMPTY)
723 smpi_group_use(*group);
724 retval = MPI_SUCCESS;
729 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
733 if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
734 retval = MPI_ERR_COMM;
735 } else if (result == NULL) {
736 retval = MPI_ERR_ARG;
738 if (comm1 == comm2) { /* Same communicators means same groups */
741 *result = smpi_group_compare(smpi_comm_group(comm1), smpi_comm_group(comm2));
742 if (*result == MPI_IDENT) {
743 *result = MPI_CONGRUENT;
746 retval = MPI_SUCCESS;
751 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
755 if (comm == MPI_COMM_NULL) {
756 retval = MPI_ERR_COMM;
757 } else if (newcomm == NULL) {
758 retval = MPI_ERR_ARG;
760 retval = smpi_comm_dup(comm, newcomm);
765 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
769 if (comm == MPI_COMM_NULL) {
770 retval = MPI_ERR_COMM;
771 } else if (group == MPI_GROUP_NULL) {
772 retval = MPI_ERR_GROUP;
773 } else if (newcomm == NULL) {
774 retval = MPI_ERR_ARG;
775 } else if(smpi_group_rank(group,smpi_process_index())==MPI_UNDEFINED){
776 *newcomm= MPI_COMM_NULL;
777 retval = MPI_SUCCESS;
779 smpi_group_use(group);
780 *newcomm = smpi_comm_new(group, NULL);
781 retval = MPI_SUCCESS;
786 int PMPI_Comm_free(MPI_Comm * comm)
791 retval = MPI_ERR_ARG;
792 } else if (*comm == MPI_COMM_NULL) {
793 retval = MPI_ERR_COMM;
795 smpi_comm_destroy(*comm);
796 *comm = MPI_COMM_NULL;
797 retval = MPI_SUCCESS;
802 int PMPI_Comm_disconnect(MPI_Comm * comm)
804 /* TODO: wait until all communication in comm are done */
808 retval = MPI_ERR_ARG;
809 } else if (*comm == MPI_COMM_NULL) {
810 retval = MPI_ERR_COMM;
812 smpi_comm_destroy(*comm);
813 *comm = MPI_COMM_NULL;
814 retval = MPI_SUCCESS;
819 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
824 if (comm_out == NULL) {
825 retval = MPI_ERR_ARG;
826 } else if (comm == MPI_COMM_NULL) {
827 retval = MPI_ERR_COMM;
829 *comm_out = smpi_comm_split(comm, color, key);
830 retval = MPI_SUCCESS;
837 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
842 if (request == NULL) {
843 retval = MPI_ERR_ARG;
844 } else if (comm == MPI_COMM_NULL) {
845 retval = MPI_ERR_COMM;
846 } else if (!is_datatype_valid(datatype)) {
847 retval = MPI_ERR_TYPE;
848 } else if (dst == MPI_PROC_NULL) {
849 retval = MPI_SUCCESS;
851 *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
852 retval = MPI_SUCCESS;
855 if (retval != MPI_SUCCESS && request != NULL)
856 *request = MPI_REQUEST_NULL;
860 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
865 if (request == NULL) {
866 retval = MPI_ERR_ARG;
867 } else if (comm == MPI_COMM_NULL) {
868 retval = MPI_ERR_COMM;
869 } else if (!is_datatype_valid(datatype)) {
870 retval = MPI_ERR_TYPE;
871 } else if (src == MPI_PROC_NULL) {
872 retval = MPI_SUCCESS;
874 *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
875 retval = MPI_SUCCESS;
878 if (retval != MPI_SUCCESS && request != NULL)
879 *request = MPI_REQUEST_NULL;
883 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
888 if (request == NULL) {
889 retval = MPI_ERR_ARG;
890 } else if (comm == MPI_COMM_NULL) {
891 retval = MPI_ERR_COMM;
892 } else if (!is_datatype_valid(datatype)) {
893 retval = MPI_ERR_TYPE;
894 } else if (dst == MPI_PROC_NULL) {
895 retval = MPI_SUCCESS;
897 *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
898 retval = MPI_SUCCESS;
901 if (retval != MPI_SUCCESS && request != NULL)
902 *request = MPI_REQUEST_NULL;
906 int PMPI_Start(MPI_Request * request)
911 if (request == NULL || *request == MPI_REQUEST_NULL) {
912 retval = MPI_ERR_REQUEST;
914 smpi_mpi_start(*request);
915 retval = MPI_SUCCESS;
921 int PMPI_Startall(int count, MPI_Request * requests)
926 if (requests == NULL) {
927 retval = MPI_ERR_ARG;
929 retval = MPI_SUCCESS;
930 for (i = 0 ; i < count ; i++) {
931 if(requests[i] == MPI_REQUEST_NULL) {
932 retval = MPI_ERR_REQUEST;
935 if(retval != MPI_ERR_REQUEST) {
936 smpi_mpi_startall(count, requests);
943 int PMPI_Request_free(MPI_Request * request)
948 if (*request == MPI_REQUEST_NULL) {
949 retval = MPI_ERR_ARG;
951 smpi_mpi_request_free(request);
952 retval = MPI_SUCCESS;
958 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
964 if (request == NULL) {
965 retval = MPI_ERR_ARG;
966 } else if (comm == MPI_COMM_NULL) {
967 retval = MPI_ERR_COMM;
968 } else if (src == MPI_PROC_NULL) {
969 *request = MPI_REQUEST_NULL;
970 retval = MPI_SUCCESS;
971 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
972 retval = MPI_ERR_RANK;
973 } else if ((count < 0) || (buf==NULL && count > 0)) {
974 retval = MPI_ERR_COUNT;
975 } else if (!is_datatype_valid(datatype)) {
976 retval = MPI_ERR_TYPE;
977 } else if(tag<0 && tag != MPI_ANY_TAG){
978 retval = MPI_ERR_TAG;
981 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
982 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
984 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
985 extra->type = TRACING_IRECV;
986 extra->src = src_traced;
989 extra->datatype1 = encode_datatype(datatype, &known);
990 int dt_size_send = 1;
992 dt_size_send = smpi_datatype_size(datatype);
993 extra->send_size = count*dt_size_send;
994 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
996 *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
997 retval = MPI_SUCCESS;
999 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1000 (*request)->recv = 1;
1004 if (retval != MPI_SUCCESS && request != NULL)
1005 *request = MPI_REQUEST_NULL;
1010 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
1015 if (request == NULL) {
1016 retval = MPI_ERR_ARG;
1017 } else if (comm == MPI_COMM_NULL) {
1018 retval = MPI_ERR_COMM;
1019 } else if (dst == MPI_PROC_NULL) {
1020 *request = MPI_REQUEST_NULL;
1021 retval = MPI_SUCCESS;
1022 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1023 retval = MPI_ERR_RANK;
1024 } else if ((count < 0) || (buf==NULL && count > 0)) {
1025 retval = MPI_ERR_COUNT;
1026 } else if (!is_datatype_valid(datatype)) {
1027 retval = MPI_ERR_TYPE;
1028 } else if(tag<0 && tag != MPI_ANY_TAG){
1029 retval = MPI_ERR_TAG;
1032 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1033 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1034 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1035 extra->type = TRACING_ISEND;
1037 extra->dst = dst_traced;
1039 extra->datatype1 = encode_datatype(datatype, &known);
1040 int dt_size_send = 1;
1042 dt_size_send = smpi_datatype_size(datatype);
1043 extra->send_size = count*dt_size_send;
1044 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1045 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1047 *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1048 retval = MPI_SUCCESS;
1050 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1051 (*request)->send = 1;
1055 if (retval != MPI_SUCCESS && request!=NULL)
1056 *request = MPI_REQUEST_NULL;
1060 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
1065 if (request == NULL) {
1066 retval = MPI_ERR_ARG;
1067 } else if (comm == MPI_COMM_NULL) {
1068 retval = MPI_ERR_COMM;
1069 } else if (dst == MPI_PROC_NULL) {
1070 *request = MPI_REQUEST_NULL;
1071 retval = MPI_SUCCESS;
1072 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1073 retval = MPI_ERR_RANK;
1074 } else if ((count < 0)|| (buf==NULL && count > 0)) {
1075 retval = MPI_ERR_COUNT;
1076 } else if (!is_datatype_valid(datatype)) {
1077 retval = MPI_ERR_TYPE;
1078 } else if(tag<0 && tag != MPI_ANY_TAG){
1079 retval = MPI_ERR_TAG;
1082 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1083 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1084 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1085 extra->type = TRACING_ISSEND;
1087 extra->dst = dst_traced;
1089 extra->datatype1 = encode_datatype(datatype, &known);
1090 int dt_size_send = 1;
1092 dt_size_send = smpi_datatype_size(datatype);
1093 extra->send_size = count*dt_size_send;
1094 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1095 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1097 *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1098 retval = MPI_SUCCESS;
1100 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1101 (*request)->send = 1;
1105 if (retval != MPI_SUCCESS && request!=NULL)
1106 *request = MPI_REQUEST_NULL;
1110 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
1115 if (comm == MPI_COMM_NULL) {
1116 retval = MPI_ERR_COMM;
1117 } else if (src == MPI_PROC_NULL) {
1118 smpi_empty_status(status);
1119 status->MPI_SOURCE = MPI_PROC_NULL;
1120 retval = MPI_SUCCESS;
1121 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1122 retval = MPI_ERR_RANK;
1123 } else if ((count < 0) || (buf==NULL && count > 0)) {
1124 retval = MPI_ERR_COUNT;
1125 } else if (!is_datatype_valid(datatype)) {
1126 retval = MPI_ERR_TYPE;
1127 } else if(tag<0 && tag != MPI_ANY_TAG){
1128 retval = MPI_ERR_TAG;
1130 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1131 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1132 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1133 extra->type = TRACING_RECV;
1134 extra->src = src_traced;
1137 extra->datatype1 = encode_datatype(datatype, &known);
1138 int dt_size_send = 1;
1140 dt_size_send = smpi_datatype_size(datatype);
1141 extra->send_size = count*dt_size_send;
1142 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
1144 smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1145 retval = MPI_SUCCESS;
1147 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1148 if(status!=MPI_STATUS_IGNORE){
1149 src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1150 if (!TRACE_smpi_view_internals()) {
1151 TRACE_smpi_recv(rank, src_traced, rank);
1154 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1161 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
1167 if (comm == MPI_COMM_NULL) {
1168 retval = MPI_ERR_COMM;
1169 } else if (dst == MPI_PROC_NULL) {
1170 retval = MPI_SUCCESS;
1171 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1172 retval = MPI_ERR_RANK;
1173 } else if ((count < 0) || (buf==NULL && count > 0)) {
1174 retval = MPI_ERR_COUNT;
1175 } else if (!is_datatype_valid(datatype)) {
1176 retval = MPI_ERR_TYPE;
1177 } else if(tag<0 && tag != MPI_ANY_TAG){
1178 retval = MPI_ERR_TAG;
1181 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1182 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1183 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1184 extra->type = TRACING_SEND;
1186 extra->dst = dst_traced;
1188 extra->datatype1 = encode_datatype(datatype, &known);
1189 int dt_size_send = 1;
1191 dt_size_send = smpi_datatype_size(datatype);
1192 extra->send_size = count*dt_size_send;
1193 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1194 if (!TRACE_smpi_view_internals()) {
1195 TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1198 smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1199 retval = MPI_SUCCESS;
1201 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1210 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1215 if (comm == MPI_COMM_NULL) {
1216 retval = MPI_ERR_COMM;
1217 } else if (dst == MPI_PROC_NULL) {
1218 retval = MPI_SUCCESS;
1219 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1220 retval = MPI_ERR_RANK;
1221 } else if ((count < 0) || (buf==NULL && count > 0)) {
1222 retval = MPI_ERR_COUNT;
1223 } else if (!is_datatype_valid(datatype)){
1224 retval = MPI_ERR_TYPE;
1225 } else if(tag<0 && tag != MPI_ANY_TAG){
1226 retval = MPI_ERR_TAG;
1229 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1230 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1231 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1232 extra->type = TRACING_SSEND;
1234 extra->dst = dst_traced;
1236 extra->datatype1 = encode_datatype(datatype, &known);
1237 int dt_size_send = 1;
1239 dt_size_send = smpi_datatype_size(datatype);
1240 extra->send_size = count*dt_size_send;
1241 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1242 TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1244 smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1245 retval = MPI_SUCCESS;
1247 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1253 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void *recvbuf,
1254 int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status * status)
1260 if (comm == MPI_COMM_NULL) {
1261 retval = MPI_ERR_COMM;
1262 } else if (!is_datatype_valid(sendtype)
1263 || !is_datatype_valid(recvtype)) {
1264 retval = MPI_ERR_TYPE;
1265 } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1266 smpi_empty_status(status);
1267 status->MPI_SOURCE = MPI_PROC_NULL;
1268 retval = MPI_SUCCESS;
1269 }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1270 (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1271 retval = MPI_ERR_RANK;
1272 } else if ((sendcount < 0 || recvcount<0) ||
1273 (sendbuf==NULL && sendcount > 0) || (recvbuf==NULL && recvcount>0)) {
1274 retval = MPI_ERR_COUNT;
1275 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1276 retval = MPI_ERR_TAG;
1279 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1280 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1281 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1282 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1283 extra->type = TRACING_SENDRECV;
1284 extra->src = src_traced;
1285 extra->dst = dst_traced;
1287 extra->datatype1 = encode_datatype(sendtype, &known);
1288 int dt_size_send = 1;
1290 dt_size_send = smpi_datatype_size(sendtype);
1291 extra->send_size = sendcount*dt_size_send;
1292 extra->datatype2 = encode_datatype(recvtype, &known);
1293 int dt_size_recv = 1;
1295 dt_size_recv = smpi_datatype_size(recvtype);
1296 extra->recv_size = recvcount*dt_size_recv;
1298 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1299 TRACE_smpi_send(rank, rank, dst_traced,sendcount*smpi_datatype_size(sendtype));
1301 smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1302 recvcount, recvtype, src, recvtag, comm, status);
1303 retval = MPI_SUCCESS;
1305 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1306 TRACE_smpi_recv(rank, src_traced, rank);
1313 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
1314 MPI_Comm comm, MPI_Status * status)
1318 if (!is_datatype_valid(datatype)) {
1319 retval = MPI_ERR_TYPE;
1320 } else if (count < 0) {
1321 retval = MPI_ERR_COUNT;
1323 int size = smpi_datatype_get_extent(datatype) * count;
1324 recvbuf = xbt_new0(char, size);
1325 retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
1326 if(retval==MPI_SUCCESS){
1327 smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1335 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1339 if (request == NULL || flag == NULL) {
1340 retval = MPI_ERR_ARG;
1341 } else if (*request == MPI_REQUEST_NULL) {
1343 smpi_empty_status(status);
1344 retval = MPI_SUCCESS;
1346 int rank = (request!=NULL && (*request)->comm != MPI_COMM_NULL) ? smpi_process_index() : -1;
1348 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1349 extra->type = TRACING_TEST;
1350 TRACE_smpi_testing_in(rank, extra);
1352 *flag = smpi_mpi_test(request, status);
1354 TRACE_smpi_testing_out(rank);
1355 retval = MPI_SUCCESS;
1361 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
1366 if (index == NULL || flag == NULL) {
1367 retval = MPI_ERR_ARG;
1369 *flag = smpi_mpi_testany(count, requests, index, status);
1370 retval = MPI_SUCCESS;
1376 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1382 retval = MPI_ERR_ARG;
1384 *flag = smpi_mpi_testall(count, requests, statuses);
1385 retval = MPI_SUCCESS;
1391 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1395 if (status == NULL) {
1396 retval = MPI_ERR_ARG;
1397 } else if (comm == MPI_COMM_NULL) {
1398 retval = MPI_ERR_COMM;
1399 } else if (source == MPI_PROC_NULL) {
1400 smpi_empty_status(status);
1401 status->MPI_SOURCE = MPI_PROC_NULL;
1402 retval = MPI_SUCCESS;
1404 smpi_mpi_probe(source, tag, comm, status);
1405 retval = MPI_SUCCESS;
1411 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1415 if ((flag == NULL) || (status == NULL)) {
1416 retval = MPI_ERR_ARG;
1417 } else if (comm == MPI_COMM_NULL) {
1418 retval = MPI_ERR_COMM;
1419 } else if (source == MPI_PROC_NULL) {
1421 smpi_empty_status(status);
1422 status->MPI_SOURCE = MPI_PROC_NULL;
1423 retval = MPI_SUCCESS;
1425 smpi_mpi_iprobe(source, tag, comm, flag, status);
1426 retval = MPI_SUCCESS;
1432 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1438 smpi_empty_status(status);
1440 if (request == NULL) {
1441 retval = MPI_ERR_ARG;
1442 } else if (*request == MPI_REQUEST_NULL) {
1443 retval = MPI_SUCCESS;
1446 int rank = (request!=NULL && (*request)->comm != MPI_COMM_NULL) ? smpi_process_index() : -1;
1448 int src_traced = (*request)->src;
1449 int dst_traced = (*request)->dst;
1450 MPI_Comm comm = (*request)->comm;
1451 int is_wait_for_receive = (*request)->recv;
1452 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1453 extra->type = TRACING_WAIT;
1454 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1456 smpi_mpi_wait(request, status);
1457 retval = MPI_SUCCESS;
1459 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1460 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1461 if (is_wait_for_receive) {
1462 if(src_traced==MPI_ANY_SOURCE)
1463 src_traced = (status!=MPI_STATUS_IGNORE) ?
1464 smpi_group_rank(smpi_comm_group(comm), status->MPI_SOURCE) :
1466 TRACE_smpi_recv(rank, src_traced, dst_traced);
1474 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1480 //save requests information for tracing
1482 int *srcs = NULL, *dsts = NULL, *recvs = NULL;
1483 MPI_Comm* comms = NULL;
1485 srcs = xbt_new0(int, count);
1486 dsts = xbt_new0(int, count);
1487 recvs = xbt_new0(int, count);
1488 comms = xbt_new0(MPI_Comm, count);
1490 for (i = 0; i < count; i++) {
1491 MPI_Request req = requests[i]; //already received requests are no longer valid
1495 recvs[i] = req->recv;
1496 comms[i] = req->comm;
1499 int rank_traced = smpi_process_index();
1500 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1501 extra->type = TRACING_WAITANY;
1502 extra->send_size=count;
1503 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1505 *index = smpi_mpi_waitany(count, requests, status);
1507 if(*index!=MPI_UNDEFINED){
1508 int src_traced = srcs[*index];
1509 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1510 int dst_traced = dsts[*index];
1511 int is_wait_for_receive = recvs[*index];
1512 if (is_wait_for_receive) {
1513 if(srcs[*index]==MPI_ANY_SOURCE)
1514 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1515 smpi_group_rank(smpi_comm_group(comms[*index]), status->MPI_SOURCE) : srcs[*index];
1516 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1518 TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1529 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1532 //save information from requests
1534 int *srcs = xbt_new0(int, count);
1535 int *dsts = xbt_new0(int, count);
1536 int *recvs = xbt_new0(int, count);
1537 int *valid = xbt_new0(int, count);
1538 MPI_Comm *comms = xbt_new0(MPI_Comm, count);
1540 for (i = 0; i < count; i++) {
1541 MPI_Request req = requests[i];
1542 if(req!=MPI_REQUEST_NULL){
1545 recvs[i] = req->recv;
1546 comms[i] = req->comm;
1552 int rank_traced = smpi_process_index();
1553 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1554 extra->type = TRACING_WAITALL;
1555 extra->send_size=count;
1556 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1558 int retval = smpi_mpi_waitall(count, requests, status);
1560 for (i = 0; i < count; i++) {
1562 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1563 int src_traced = srcs[i];
1564 int dst_traced = dsts[i];
1565 int is_wait_for_receive = recvs[i];
1566 if (is_wait_for_receive) {
1567 if(src_traced==MPI_ANY_SOURCE)
1568 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1569 smpi_group_rank(smpi_comm_group(comms[i]), status[i].MPI_SOURCE) : srcs[i];
1570 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1574 TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1585 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
1590 if (outcount == NULL) {
1591 retval = MPI_ERR_ARG;
1593 *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1594 retval = MPI_SUCCESS;
1600 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
1605 if (outcount == NULL) {
1606 retval = MPI_ERR_ARG;
1608 *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1609 retval = MPI_SUCCESS;
1616 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1622 if (comm == MPI_COMM_NULL) {
1623 retval = MPI_ERR_COMM;
1624 } else if (!is_datatype_valid(datatype)) {
1625 retval = MPI_ERR_ARG;
1627 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1628 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1630 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1631 extra->type = TRACING_BCAST;
1632 extra->root = root_traced;
1634 extra->datatype1 = encode_datatype(datatype, &known);
1635 int dt_size_send = 1;
1637 dt_size_send = smpi_datatype_size(datatype);
1638 extra->send_size = count*dt_size_send;
1639 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1640 if(smpi_comm_size(comm)>1)
1641 mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1642 retval = MPI_SUCCESS;
1644 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1651 int PMPI_Barrier(MPI_Comm comm)
1657 if (comm == MPI_COMM_NULL) {
1658 retval = MPI_ERR_COMM;
1660 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1661 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1662 extra->type = TRACING_BARRIER;
1663 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1665 mpi_coll_barrier_fun(comm);
1666 retval = MPI_SUCCESS;
1668 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1675 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,void *recvbuf, int recvcount, MPI_Datatype recvtype,
1676 int root, MPI_Comm comm)
1682 if (comm == MPI_COMM_NULL) {
1683 retval = MPI_ERR_COMM;
1684 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1685 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1686 retval = MPI_ERR_TYPE;
1687 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) || ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1688 retval = MPI_ERR_COUNT;
1691 char* sendtmpbuf = static_cast<char*>(sendbuf);
1692 int sendtmpcount = sendcount;
1693 MPI_Datatype sendtmptype = sendtype;
1694 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1696 sendtmptype=recvtype;
1698 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1699 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1700 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1701 extra->type = TRACING_GATHER;
1702 extra->root = root_traced;
1704 extra->datatype1 = encode_datatype(sendtmptype, &known);
1705 int dt_size_send = 1;
1707 dt_size_send = smpi_datatype_size(sendtmptype);
1708 extra->send_size = sendtmpcount*dt_size_send;
1709 extra->datatype2 = encode_datatype(recvtype, &known);
1710 int dt_size_recv = 1;
1711 if((smpi_comm_rank(comm)==root) && known==0)
1712 dt_size_recv = smpi_datatype_size(recvtype);
1713 extra->recv_size = recvcount*dt_size_recv;
1715 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1717 mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount, recvtype, root, comm);
1719 retval = MPI_SUCCESS;
1720 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1727 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int *recvcounts, int *displs,
1728 MPI_Datatype recvtype, int root, MPI_Comm comm)
1734 if (comm == MPI_COMM_NULL) {
1735 retval = MPI_ERR_COMM;
1736 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1737 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1738 retval = MPI_ERR_TYPE;
1739 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1740 retval = MPI_ERR_COUNT;
1741 } else if (recvcounts == NULL || displs == NULL) {
1742 retval = MPI_ERR_ARG;
1744 char* sendtmpbuf = static_cast<char*>(sendbuf);
1745 int sendtmpcount = sendcount;
1746 MPI_Datatype sendtmptype = sendtype;
1747 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1749 sendtmptype=recvtype;
1752 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1753 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1755 int size = smpi_comm_size(comm);
1756 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1757 extra->type = TRACING_GATHERV;
1758 extra->num_processes = size;
1759 extra->root = root_traced;
1761 extra->datatype1 = encode_datatype(sendtmptype, &known);
1762 int dt_size_send = 1;
1764 dt_size_send = smpi_datatype_size(sendtype);
1765 extra->send_size = sendtmpcount*dt_size_send;
1766 extra->datatype2 = encode_datatype(recvtype, &known);
1767 int dt_size_recv = 1;
1769 dt_size_recv = smpi_datatype_size(recvtype);
1770 if((smpi_comm_rank(comm)==root)){
1771 extra->recvcounts= xbt_new(int,size);
1772 for(i=0; i< size; i++)//copy data to avoid bad free
1773 extra->recvcounts[i] = recvcounts[i]*dt_size_recv;
1775 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1777 smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts, displs, recvtype, root, comm);
1778 retval = MPI_SUCCESS;
1779 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1786 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1787 void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
1793 if (comm == MPI_COMM_NULL) {
1794 retval = MPI_ERR_COMM;
1795 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1796 (recvtype == MPI_DATATYPE_NULL)){
1797 retval = MPI_ERR_TYPE;
1798 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1800 retval = MPI_ERR_COUNT;
1802 if(sendbuf == MPI_IN_PLACE) {
1803 sendbuf=static_cast<char*>(recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1804 sendcount=recvcount;
1807 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1808 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1809 extra->type = TRACING_ALLGATHER;
1811 extra->datatype1 = encode_datatype(sendtype, &known);
1812 int dt_size_send = 1;
1814 dt_size_send = smpi_datatype_size(sendtype);
1815 extra->send_size = sendcount*dt_size_send;
1816 extra->datatype2 = encode_datatype(recvtype, &known);
1817 int dt_size_recv = 1;
1819 dt_size_recv = smpi_datatype_size(recvtype);
1820 extra->recv_size = recvcount*dt_size_recv;
1822 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1824 mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1825 retval = MPI_SUCCESS;
1826 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1832 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1833 void *recvbuf, int *recvcounts, int *displs, MPI_Datatype recvtype, MPI_Comm comm)
1839 if (comm == MPI_COMM_NULL) {
1840 retval = MPI_ERR_COMM;
1841 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1842 (recvtype == MPI_DATATYPE_NULL)){
1843 retval = MPI_ERR_TYPE;
1844 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1845 retval = MPI_ERR_COUNT;
1846 } else if (recvcounts == NULL || displs == NULL) {
1847 retval = MPI_ERR_ARG;
1850 if(sendbuf == MPI_IN_PLACE) {
1851 sendbuf=static_cast<char*>(recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1852 sendcount=recvcounts[smpi_comm_rank(comm)];
1855 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1857 int size = smpi_comm_size(comm);
1858 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1859 extra->type = TRACING_ALLGATHERV;
1860 extra->num_processes = size;
1862 extra->datatype1 = encode_datatype(sendtype, &known);
1863 int dt_size_send = 1;
1865 dt_size_send = smpi_datatype_size(sendtype);
1866 extra->send_size = sendcount*dt_size_send;
1867 extra->datatype2 = encode_datatype(recvtype, &known);
1868 int dt_size_recv = 1;
1870 dt_size_recv = smpi_datatype_size(recvtype);
1871 extra->recvcounts= xbt_new(int, size);
1872 for(i=0; i< size; i++)//copy data to avoid bad free
1873 extra->recvcounts[i] = recvcounts[i]*dt_size_recv;
1875 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
1877 mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
1878 retval = MPI_SUCCESS;
1879 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1886 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1887 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
1893 if (comm == MPI_COMM_NULL) {
1894 retval = MPI_ERR_COMM;
1895 } else if (((smpi_comm_rank(comm)==root) && (!is_datatype_valid(sendtype)))
1896 || ((recvbuf !=MPI_IN_PLACE) && (!is_datatype_valid(recvtype)))){
1897 retval = MPI_ERR_TYPE;
1898 } else if ((sendbuf == recvbuf) ||
1899 ((smpi_comm_rank(comm)==root) && sendcount>0 && (sendbuf == NULL))){
1900 retval = MPI_ERR_BUFFER;
1903 if (recvbuf == MPI_IN_PLACE) {
1905 recvcount=sendcount;
1907 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1908 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1909 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1910 extra->type = TRACING_SCATTER;
1911 extra->root = root_traced;
1913 extra->datatype1 = encode_datatype(sendtype, &known);
1914 int dt_size_send = 1;
1915 if((smpi_comm_rank(comm)==root) && known==0)
1916 dt_size_send = smpi_datatype_size(sendtype);
1917 extra->send_size = sendcount*dt_size_send;
1918 extra->datatype2 = encode_datatype(recvtype, &known);
1919 int dt_size_recv = 1;
1921 dt_size_recv = smpi_datatype_size(recvtype);
1922 extra->recv_size = recvcount*dt_size_recv;
1923 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1925 mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
1926 retval = MPI_SUCCESS;
1927 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1934 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1935 MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
1941 if (comm == MPI_COMM_NULL) {
1942 retval = MPI_ERR_COMM;
1943 } else if (sendcounts == NULL || displs == NULL) {
1944 retval = MPI_ERR_ARG;
1945 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1946 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1947 retval = MPI_ERR_TYPE;
1949 if (recvbuf == MPI_IN_PLACE) {
1951 recvcount=sendcounts[smpi_comm_rank(comm)];
1953 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1954 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1956 int size = smpi_comm_size(comm);
1957 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1958 extra->type = TRACING_SCATTERV;
1959 extra->num_processes = size;
1960 extra->root = root_traced;
1962 extra->datatype1 = encode_datatype(sendtype, &known);
1963 int dt_size_send = 1;
1965 dt_size_send = smpi_datatype_size(sendtype);
1966 if((smpi_comm_rank(comm)==root)){
1967 extra->sendcounts= xbt_new(int, size);
1968 for(i=0; i< size; i++)//copy data to avoid bad free
1969 extra->sendcounts[i] = sendcounts[i]*dt_size_send;
1971 extra->datatype2 = encode_datatype(recvtype, &known);
1972 int dt_size_recv = 1;
1974 dt_size_recv = smpi_datatype_size(recvtype);
1975 extra->recv_size = recvcount*dt_size_recv;
1976 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1978 smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
1980 retval = MPI_SUCCESS;
1981 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1988 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1994 if (comm == MPI_COMM_NULL) {
1995 retval = MPI_ERR_COMM;
1996 } else if (!is_datatype_valid(datatype) || op == MPI_OP_NULL) {
1997 retval = MPI_ERR_ARG;
1999 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2000 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2001 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2002 extra->type = TRACING_REDUCE;
2004 extra->datatype1 = encode_datatype(datatype, &known);
2005 int dt_size_send = 1;
2007 dt_size_send = smpi_datatype_size(datatype);
2008 extra->send_size = count*dt_size_send;
2009 extra->root = root_traced;
2011 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2013 mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
2015 retval = MPI_SUCCESS;
2016 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2023 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count, MPI_Datatype datatype, MPI_Op op){
2027 if (!is_datatype_valid(datatype) || op == MPI_OP_NULL) {
2028 retval = MPI_ERR_ARG;
2030 smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
2037 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2043 if (comm == MPI_COMM_NULL) {
2044 retval = MPI_ERR_COMM;
2045 } else if (!is_datatype_valid(datatype)) {
2046 retval = MPI_ERR_TYPE;
2047 } else if (op == MPI_OP_NULL) {
2048 retval = MPI_ERR_OP;
2051 char* sendtmpbuf = static_cast<char*>(sendbuf);
2052 if( sendbuf == MPI_IN_PLACE ) {
2053 sendtmpbuf = static_cast<char*>(xbt_malloc(count*smpi_datatype_get_extent(datatype)));
2054 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
2056 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2057 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2058 extra->type = TRACING_ALLREDUCE;
2060 extra->datatype1 = encode_datatype(datatype, &known);
2061 int dt_size_send = 1;
2063 dt_size_send = smpi_datatype_size(datatype);
2064 extra->send_size = count*dt_size_send;
2066 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2068 mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
2070 if( sendbuf == MPI_IN_PLACE )
2071 xbt_free(sendtmpbuf);
2073 retval = MPI_SUCCESS;
2074 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2081 int PMPI_Scan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2087 if (comm == MPI_COMM_NULL) {
2088 retval = MPI_ERR_COMM;
2089 } else if (!is_datatype_valid(datatype)) {
2090 retval = MPI_ERR_TYPE;
2091 } else if (op == MPI_OP_NULL) {
2092 retval = MPI_ERR_OP;
2094 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2095 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2096 extra->type = TRACING_SCAN;
2098 extra->datatype1 = encode_datatype(datatype, &known);
2099 int dt_size_send = 1;
2101 dt_size_send = smpi_datatype_size(datatype);
2102 extra->send_size = count*dt_size_send;
2104 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2106 smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2108 retval = MPI_SUCCESS;
2109 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2116 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm){
2121 if (comm == MPI_COMM_NULL) {
2122 retval = MPI_ERR_COMM;
2123 } else if (!is_datatype_valid(datatype)) {
2124 retval = MPI_ERR_TYPE;
2125 } else if (op == MPI_OP_NULL) {
2126 retval = MPI_ERR_OP;
2128 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2129 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2130 extra->type = TRACING_EXSCAN;
2132 extra->datatype1 = encode_datatype(datatype, &known);
2133 int dt_size_send = 1;
2135 dt_size_send = smpi_datatype_size(datatype);
2136 extra->send_size = count*dt_size_send;
2137 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2139 smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2140 retval = MPI_SUCCESS;
2141 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2148 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2153 if (comm == MPI_COMM_NULL) {
2154 retval = MPI_ERR_COMM;
2155 } else if (!is_datatype_valid(datatype)) {
2156 retval = MPI_ERR_TYPE;
2157 } else if (op == MPI_OP_NULL) {
2158 retval = MPI_ERR_OP;
2159 } else if (recvcounts == NULL) {
2160 retval = MPI_ERR_ARG;
2162 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2164 int size = smpi_comm_size(comm);
2165 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2166 extra->type = TRACING_REDUCE_SCATTER;
2167 extra->num_processes = size;
2169 extra->datatype1 = encode_datatype(datatype, &known);
2170 int dt_size_send = 1;
2172 dt_size_send = smpi_datatype_size(datatype);
2173 extra->send_size = 0;
2174 extra->recvcounts= xbt_new(int, size);
2175 for(i=0; i< size; i++)//copy data to avoid bad free
2176 extra->recvcounts[i] = recvcounts[i]*dt_size_send;
2177 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2179 void* sendtmpbuf=sendbuf;
2180 if(sendbuf==MPI_IN_PLACE)
2183 mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts, datatype, op, comm);
2184 retval = MPI_SUCCESS;
2185 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2192 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2193 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2198 if (comm == MPI_COMM_NULL) {
2199 retval = MPI_ERR_COMM;
2200 } else if (!is_datatype_valid(datatype)) {
2201 retval = MPI_ERR_TYPE;
2202 } else if (op == MPI_OP_NULL) {
2203 retval = MPI_ERR_OP;
2204 } else if (recvcount < 0) {
2205 retval = MPI_ERR_ARG;
2207 int count=smpi_comm_size(comm);
2209 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2210 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2211 extra->type = TRACING_REDUCE_SCATTER;
2212 extra->num_processes = count;
2214 extra->datatype1 = encode_datatype(datatype, &known);
2215 int dt_size_send = 1;
2217 dt_size_send = smpi_datatype_size(datatype);
2218 extra->send_size = 0;
2219 extra->recvcounts= xbt_new(int, count);
2220 for(i=0; i< count; i++)//copy data to avoid bad free
2221 extra->recvcounts[i] = recvcount*dt_size_send;
2223 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2225 int* recvcounts=static_cast<int*>(xbt_malloc(count));
2226 for (i=0; i<count;i++)
2227 recvcounts[i]=recvcount;
2228 mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts, datatype, op, comm);
2229 xbt_free(recvcounts);
2230 retval = MPI_SUCCESS;
2232 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2239 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2240 void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
2246 if (comm == MPI_COMM_NULL) {
2247 retval = MPI_ERR_COMM;
2248 } else if (sendtype == MPI_DATATYPE_NULL
2249 || recvtype == MPI_DATATYPE_NULL) {
2250 retval = MPI_ERR_TYPE;
2252 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2253 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2254 extra->type = TRACING_ALLTOALL;
2256 extra->datatype1 = encode_datatype(sendtype, &known);
2258 extra->send_size = sendcount*smpi_datatype_size(sendtype);
2260 extra->send_size = sendcount;
2261 extra->datatype2 = encode_datatype(recvtype, &known);
2263 extra->recv_size = recvcount*smpi_datatype_size(recvtype);
2265 extra->recv_size = recvcount;
2266 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2268 retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2270 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2277 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2278 int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2284 if (comm == MPI_COMM_NULL) {
2285 retval = MPI_ERR_COMM;
2286 } else if (sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
2287 retval = MPI_ERR_TYPE;
2288 } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL || recvdisps == NULL) {
2289 retval = MPI_ERR_ARG;
2291 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2293 int size = smpi_comm_size(comm);
2294 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2295 extra->type = TRACING_ALLTOALLV;
2296 extra->send_size = 0;
2297 extra->recv_size = 0;
2298 extra->recvcounts= xbt_new(int, size);
2299 extra->sendcounts= xbt_new(int, size);
2301 extra->datatype1 = encode_datatype(sendtype, &known);
2302 int dt_size_send = 1;
2304 dt_size_send = smpi_datatype_size(sendtype);
2305 int dt_size_recv = 1;
2306 extra->datatype2 = encode_datatype(recvtype, &known);
2308 dt_size_recv = smpi_datatype_size(recvtype);
2309 for(i=0; i< size; i++){//copy data to avoid bad free
2310 extra->send_size += sendcounts[i]*dt_size_send;
2311 extra->recv_size += recvcounts[i]*dt_size_recv;
2313 extra->sendcounts[i] = sendcounts[i]*dt_size_send;
2314 extra->recvcounts[i] = recvcounts[i]*dt_size_recv;
2316 extra->num_processes = size;
2317 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2319 retval = mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype, recvbuf, recvcounts, recvdisps, recvtype,
2321 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2329 int PMPI_Get_processor_name(char *name, int *resultlen)
2331 int retval = MPI_SUCCESS;
2333 strncpy(name, sg_host_get_name(SIMIX_host_self()),
2334 strlen(sg_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2335 strlen(sg_host_get_name(SIMIX_host_self())) +1 : MPI_MAX_PROCESSOR_NAME - 1 );
2336 *resultlen = strlen(name) > MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2341 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2343 int retval = MPI_SUCCESS;
2346 if (status == NULL || count == NULL) {
2347 retval = MPI_ERR_ARG;
2348 } else if (!is_datatype_valid(datatype)) {
2349 retval = MPI_ERR_TYPE;
2351 size = smpi_datatype_size(datatype);
2354 } else if (status->count % size != 0) {
2355 retval = MPI_UNDEFINED;
2357 *count = smpi_mpi_get_count(status, datatype);
2363 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2366 if (old_type == MPI_DATATYPE_NULL) {
2367 retval = MPI_ERR_TYPE;
2368 } else if (count<0){
2369 retval = MPI_ERR_COUNT;
2371 retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2376 int PMPI_Type_commit(MPI_Datatype* datatype) {
2379 if (datatype == NULL || *datatype == MPI_DATATYPE_NULL) {
2380 retval = MPI_ERR_TYPE;
2382 smpi_datatype_commit(datatype);
2383 retval = MPI_SUCCESS;
2388 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2391 if (old_type == MPI_DATATYPE_NULL) {
2392 retval = MPI_ERR_TYPE;
2393 } else if (count<0 || blocklen<0){
2394 retval = MPI_ERR_COUNT;
2396 retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2401 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2404 if (old_type == MPI_DATATYPE_NULL) {
2405 retval = MPI_ERR_TYPE;
2406 } else if (count<0 || blocklen<0){
2407 retval = MPI_ERR_COUNT;
2409 retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2414 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2415 return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2418 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2421 if (old_type == MPI_DATATYPE_NULL) {
2422 retval = MPI_ERR_TYPE;
2423 } else if (count<0){
2424 retval = MPI_ERR_COUNT;
2426 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2431 int PMPI_Type_create_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2434 if (old_type == MPI_DATATYPE_NULL) {
2435 retval = MPI_ERR_TYPE;
2436 } else if (count<0){
2437 retval = MPI_ERR_COUNT;
2439 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2444 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type,
2445 MPI_Datatype* new_type) {
2448 if (old_type == MPI_DATATYPE_NULL) {
2449 retval = MPI_ERR_TYPE;
2450 } else if (count<0){
2451 retval = MPI_ERR_COUNT;
2453 int* blocklens=static_cast<int*>(xbt_malloc(blocklength*count));
2454 for (i=0; i<count;i++)
2455 blocklens[i]=blocklength;
2456 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2457 xbt_free(blocklens);
2462 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2465 if (old_type == MPI_DATATYPE_NULL) {
2466 retval = MPI_ERR_TYPE;
2467 } else if (count<0){
2468 retval = MPI_ERR_COUNT;
2470 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2475 int PMPI_Type_create_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type,
2476 MPI_Datatype* new_type) {
2477 return PMPI_Type_hindexed(count, blocklens,indices,old_type,new_type);
2480 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type,
2481 MPI_Datatype* new_type) {
2484 if (old_type == MPI_DATATYPE_NULL) {
2485 retval = MPI_ERR_TYPE;
2486 } else if (count<0){
2487 retval = MPI_ERR_COUNT;
2489 int* blocklens=(int*)xbt_malloc(blocklength*count);
2490 for (i=0; i<count;i++)blocklens[i]=blocklength;
2491 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2492 xbt_free(blocklens);
2497 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2501 retval = MPI_ERR_COUNT;
2503 retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2508 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types,
2509 MPI_Datatype* new_type) {
2510 return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2513 int PMPI_Error_class(int errorcode, int* errorclass) {
2514 // assume smpi uses only standard mpi error codes
2515 *errorclass=errorcode;
2519 int PMPI_Initialized(int* flag) {
2520 *flag=smpi_process_initialized();
2524 /* The topo part of MPI_COMM_WORLD should always be NULL. When other topologies will be implemented, not only should we
2525 * check if the topology is NULL, but we should check if it is the good topology type (so we have to add a
2526 * MPIR_Topo_Type field, and replace the MPI_Topology field by an union)*/
2528 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periodic, int reorder, MPI_Comm* comm_cart) {
2530 if (comm_old == MPI_COMM_NULL){
2531 retval = MPI_ERR_COMM;
2532 } else if (ndims < 0 || (ndims > 0 && (dims == NULL || periodic == NULL)) || comm_cart == NULL) {
2533 retval = MPI_ERR_ARG;
2535 retval = smpi_mpi_cart_create(comm_old, ndims, dims, periodic, reorder, comm_cart);
2540 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2541 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2542 return MPI_ERR_TOPOLOGY;
2544 if (coords == NULL) {
2547 return smpi_mpi_cart_rank(comm, coords, rank);
2550 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2551 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2552 return MPI_ERR_TOPOLOGY;
2554 if (source == NULL || dest == NULL || direction < 0 ) {
2557 return smpi_mpi_cart_shift(comm, direction, displ, source, dest);
2560 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2561 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2562 return MPI_ERR_TOPOLOGY;
2564 if (rank < 0 || rank >= smpi_comm_size(comm)) {
2565 return MPI_ERR_RANK;
2570 if(coords == NULL) {
2573 return smpi_mpi_cart_coords(comm, rank, maxdims, coords);
2576 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2577 if(comm == NULL || smpi_comm_topo(comm) == NULL) {
2578 return MPI_ERR_TOPOLOGY;
2580 if(maxdims <= 0 || dims == NULL || periods == NULL || coords == NULL) {
2583 return smpi_mpi_cart_get(comm, maxdims, dims, periods, coords);
2586 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2587 if (comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2588 return MPI_ERR_TOPOLOGY;
2590 if (ndims == NULL) {
2593 return smpi_mpi_cartdim_get(comm, ndims);
2596 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2600 if (ndims < 1 || nnodes < 1) {
2601 return MPI_ERR_DIMS;
2604 return smpi_mpi_dims_create(nnodes, ndims, dims);
2607 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2608 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2609 return MPI_ERR_TOPOLOGY;
2611 if (comm_new == NULL) {
2614 return smpi_mpi_cart_sub(comm, remain_dims, comm_new);
2617 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2618 if(oldtype == MPI_DATATYPE_NULL) {
2619 return MPI_ERR_TYPE;
2621 int blocks[3] = { 1, 1, 1 };
2622 MPI_Aint disps[3] = { lb, 0, lb+extent };
2623 MPI_Datatype types[3] = { MPI_LB, oldtype, MPI_UB };
2625 s_smpi_mpi_struct_t* subtype = smpi_datatype_struct_create( blocks, disps, 3, types);
2626 smpi_datatype_create(newtype,oldtype->size, lb, lb + extent, sizeof(s_smpi_mpi_struct_t) , subtype, DT_FLAG_VECTOR);
2628 (*newtype)->flags &= ~DT_FLAG_COMMITED;
2632 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2635 if (comm == MPI_COMM_NULL) {
2636 retval= MPI_ERR_COMM;
2637 }else if ((base == NULL && size != 0) || disp_unit <= 0 || size < 0 ){
2638 retval= MPI_ERR_OTHER;
2640 *win = smpi_mpi_win_create( base, size, disp_unit, info, comm);
2641 retval = MPI_SUCCESS;
2647 int PMPI_Win_free( MPI_Win* win){
2650 if (win == NULL || *win == MPI_WIN_NULL) {
2651 retval = MPI_ERR_WIN;
2653 retval=smpi_mpi_win_free(win);
2659 int PMPI_Win_set_name(MPI_Win win, char * name)
2662 if (win == MPI_WIN_NULL) {
2663 retval = MPI_ERR_TYPE;
2664 } else if (name == NULL) {
2665 retval = MPI_ERR_ARG;
2667 smpi_mpi_win_set_name(win, name);
2668 retval = MPI_SUCCESS;
2673 int PMPI_Win_get_name(MPI_Win win, char * name, int* len)
2675 int retval = MPI_SUCCESS;
2677 if (win == MPI_WIN_NULL) {
2678 retval = MPI_ERR_WIN;
2679 } else if (name == NULL) {
2680 retval = MPI_ERR_ARG;
2682 smpi_mpi_win_get_name(win, name, len);
2687 int PMPI_Win_get_group(MPI_Win win, MPI_Group * group){
2688 int retval = MPI_SUCCESS;
2689 if (win == MPI_WIN_NULL) {
2690 retval = MPI_ERR_WIN;
2692 smpi_mpi_win_get_group(win, group);
2693 smpi_group_use(*group);
2699 int PMPI_Win_fence( int assert, MPI_Win win){
2702 if (win == MPI_WIN_NULL) {
2703 retval = MPI_ERR_WIN;
2705 int rank = smpi_process_index();
2706 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, NULL);
2707 retval = smpi_mpi_win_fence(assert, win);
2708 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2714 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2715 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2718 if (win == MPI_WIN_NULL) {
2719 retval = MPI_ERR_WIN;
2720 } else if (target_rank == MPI_PROC_NULL) {
2721 retval = MPI_SUCCESS;
2722 } else if (target_rank <0){
2723 retval = MPI_ERR_RANK;
2724 } else if (target_disp <0){
2725 retval = MPI_ERR_ARG;
2726 } else if ((origin_count < 0 || target_count < 0) ||
2727 (origin_addr==NULL && origin_count > 0)){
2728 retval = MPI_ERR_COUNT;
2729 } else if ((!is_datatype_valid(origin_datatype)) || (!is_datatype_valid(target_datatype))) {
2730 retval = MPI_ERR_TYPE;
2732 int rank = smpi_process_index();
2734 smpi_mpi_win_get_group(win, &group);
2735 int src_traced = smpi_group_index(group, target_rank);
2736 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, NULL);
2738 retval = smpi_mpi_get( origin_addr, origin_count, origin_datatype, target_rank, target_disp, target_count,
2739 target_datatype, win);
2741 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
2747 int PMPI_Put( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2748 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2751 if (win == MPI_WIN_NULL) {
2752 retval = MPI_ERR_WIN;
2753 } else if (target_rank == MPI_PROC_NULL) {
2754 retval = MPI_SUCCESS;
2755 } else if (target_rank <0){
2756 retval = MPI_ERR_RANK;
2757 } else if (target_disp <0){
2758 retval = MPI_ERR_ARG;
2759 } else if ((origin_count < 0 || target_count < 0) ||
2760 (origin_addr==NULL && origin_count > 0)){
2761 retval = MPI_ERR_COUNT;
2762 } else if ((!is_datatype_valid(origin_datatype)) || (!is_datatype_valid(target_datatype))) {
2763 retval = MPI_ERR_TYPE;
2765 int rank = smpi_process_index();
2767 smpi_mpi_win_get_group(win, &group);
2768 int dst_traced = smpi_group_index(group, target_rank);
2769 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, NULL);
2770 TRACE_smpi_send(rank, rank, dst_traced, origin_count*smpi_datatype_size(origin_datatype));
2772 retval = smpi_mpi_put( origin_addr, origin_count, origin_datatype, target_rank, target_disp, target_count,
2773 target_datatype, win);
2775 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
2781 int PMPI_Accumulate( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2782 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Op op, MPI_Win win){
2785 if (win == MPI_WIN_NULL) {
2786 retval = MPI_ERR_WIN;
2787 } else if (target_rank == MPI_PROC_NULL) {
2788 retval = MPI_SUCCESS;
2789 } else if (target_rank <0){
2790 retval = MPI_ERR_RANK;
2791 } else if (target_disp <0){
2792 retval = MPI_ERR_ARG;
2793 } else if ((origin_count < 0 || target_count < 0) ||
2794 (origin_addr==NULL && origin_count > 0)){
2795 retval = MPI_ERR_COUNT;
2796 } else if ((!is_datatype_valid(origin_datatype)) ||
2797 (!is_datatype_valid(target_datatype))) {
2798 retval = MPI_ERR_TYPE;
2799 } else if (op == MPI_OP_NULL) {
2800 retval = MPI_ERR_OP;
2802 int rank = smpi_process_index();
2804 smpi_mpi_win_get_group(win, &group);
2805 int src_traced = smpi_group_index(group, target_rank);
2806 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, NULL);
2808 retval = smpi_mpi_accumulate( origin_addr, origin_count, origin_datatype, target_rank, target_disp, target_count,
2809 target_datatype, op, win);
2811 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
2817 int PMPI_Win_post(MPI_Group group, int assert, MPI_Win win){
2820 if (win == MPI_WIN_NULL) {
2821 retval = MPI_ERR_WIN;
2822 } else if (group==MPI_GROUP_NULL){
2823 retval = MPI_ERR_GROUP;
2826 int rank = smpi_process_index();
2827 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, NULL);
2828 retval = smpi_mpi_win_post(group,assert,win);
2829 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2835 int PMPI_Win_start(MPI_Group group, int assert, MPI_Win win){
2838 if (win == MPI_WIN_NULL) {
2839 retval = MPI_ERR_WIN;
2840 } else if (group==MPI_GROUP_NULL){
2841 retval = MPI_ERR_GROUP;
2844 int rank = smpi_process_index();
2845 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, NULL);
2846 retval = smpi_mpi_win_start(group,assert,win);
2847 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2853 int PMPI_Win_complete(MPI_Win win){
2856 if (win == MPI_WIN_NULL) {
2857 retval = MPI_ERR_WIN;
2860 int rank = smpi_process_index();
2861 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, NULL);
2863 retval = smpi_mpi_win_complete(win);
2865 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2871 int PMPI_Win_wait(MPI_Win win){
2874 if (win == MPI_WIN_NULL) {
2875 retval = MPI_ERR_WIN;
2878 int rank = smpi_process_index();
2879 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, NULL);
2881 retval = smpi_mpi_win_wait(win);
2883 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2889 int PMPI_Alloc_mem(MPI_Aint size, MPI_Info info, void *baseptr){
2890 void *ptr = xbt_malloc(size);
2892 return MPI_ERR_NO_MEM;
2894 *static_cast<void**>(baseptr) = ptr;
2899 int PMPI_Free_mem(void *baseptr){
2904 int PMPI_Type_set_name(MPI_Datatype datatype, char * name)
2907 if (datatype == MPI_DATATYPE_NULL) {
2908 retval = MPI_ERR_TYPE;
2909 } else if (name == NULL) {
2910 retval = MPI_ERR_ARG;
2912 smpi_datatype_set_name(datatype, name);
2913 retval = MPI_SUCCESS;
2918 int PMPI_Type_get_name(MPI_Datatype datatype, char * name, int* len)
2922 if (datatype == MPI_DATATYPE_NULL) {
2923 retval = MPI_ERR_TYPE;
2924 } else if (name == NULL) {
2925 retval = MPI_ERR_ARG;
2927 smpi_datatype_get_name(datatype, name, len);
2928 retval = MPI_SUCCESS;
2933 MPI_Datatype PMPI_Type_f2c(MPI_Fint datatype){
2934 return smpi_type_f2c(datatype);
2937 MPI_Fint PMPI_Type_c2f(MPI_Datatype datatype){
2938 return smpi_type_c2f( datatype);
2941 MPI_Group PMPI_Group_f2c(MPI_Fint group){
2942 return smpi_group_f2c( group);
2945 MPI_Fint PMPI_Group_c2f(MPI_Group group){
2946 return smpi_group_c2f(group);
2949 MPI_Request PMPI_Request_f2c(MPI_Fint request){
2950 return smpi_request_f2c(request);
2953 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
2954 return smpi_request_c2f(request);
2957 MPI_Win PMPI_Win_f2c(MPI_Fint win){
2958 return smpi_win_f2c(win);
2961 MPI_Fint PMPI_Win_c2f(MPI_Win win){
2962 return smpi_win_c2f(win);
2965 MPI_Op PMPI_Op_f2c(MPI_Fint op){
2966 return smpi_op_f2c(op);
2969 MPI_Fint PMPI_Op_c2f(MPI_Op op){
2970 return smpi_op_c2f(op);
2973 MPI_Comm PMPI_Comm_f2c(MPI_Fint comm){
2974 return smpi_comm_f2c(comm);
2977 MPI_Fint PMPI_Comm_c2f(MPI_Comm comm){
2978 return smpi_comm_c2f(comm);
2981 MPI_Info PMPI_Info_f2c(MPI_Fint info){
2982 return smpi_info_f2c(info);
2985 MPI_Fint PMPI_Info_c2f(MPI_Info info){
2986 return smpi_info_c2f(info);
2989 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2990 return smpi_comm_keyval_create(copy_fn, delete_fn, keyval, extra_state);
2993 int PMPI_Keyval_free(int* keyval) {
2994 return smpi_comm_keyval_free(keyval);
2997 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2998 if(keyval == MPI_TAG_UB||keyval == MPI_HOST||keyval == MPI_IO ||keyval == MPI_WTIME_IS_GLOBAL||keyval == MPI_APPNUM
2999 ||keyval == MPI_UNIVERSE_SIZE||keyval == MPI_LASTUSEDCODE)
3001 else if (comm==MPI_COMM_NULL)
3002 return MPI_ERR_COMM;
3004 return smpi_comm_attr_delete(comm, keyval);
3007 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
3009 static int zero = 0;
3010 static int tag_ub = 1000000;
3011 static int last_used_code = MPI_ERR_LASTCODE;
3013 if (comm==MPI_COMM_NULL){
3015 return MPI_ERR_COMM;
3023 *static_cast<int**>(attr_value) = &zero;
3025 case MPI_UNIVERSE_SIZE:
3027 *static_cast<int**>(attr_value) = &smpi_universe_size;
3029 case MPI_LASTUSEDCODE:
3031 *static_cast<int**>(attr_value) = &last_used_code;
3035 *static_cast<int**>(attr_value) = &tag_ub;
3037 case MPI_WTIME_IS_GLOBAL:
3039 *static_cast<int**>(attr_value) = &one;
3042 return smpi_comm_attr_get(comm, keyval, attr_value, flag);
3046 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
3047 if(keyval == MPI_TAG_UB||keyval == MPI_HOST||keyval == MPI_IO ||keyval == MPI_WTIME_IS_GLOBAL||keyval == MPI_APPNUM
3048 ||keyval == MPI_UNIVERSE_SIZE||keyval == MPI_LASTUSEDCODE)
3050 else if (comm==MPI_COMM_NULL)
3051 return MPI_ERR_COMM;
3053 return smpi_comm_attr_put(comm, keyval, attr_value);
3056 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
3058 return PMPI_Attr_get(comm, comm_keyval, attribute_val,flag);
3061 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
3063 return PMPI_Attr_put(comm, comm_keyval, attribute_val);
3066 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
3068 return PMPI_Attr_delete(comm, comm_keyval);
3071 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval,
3074 return PMPI_Keyval_create(copy_fn, delete_fn, keyval, extra_state);
3077 int PMPI_Comm_free_keyval(int* keyval) {
3078 return PMPI_Keyval_free(keyval);
3081 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
3083 if (type==MPI_DATATYPE_NULL)
3084 return MPI_ERR_TYPE;
3086 return smpi_type_attr_get(type, type_keyval, attribute_val, flag);
3089 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
3091 if (type==MPI_DATATYPE_NULL)
3092 return MPI_ERR_TYPE;
3094 return smpi_type_attr_put(type, type_keyval, attribute_val);
3097 int PMPI_Type_delete_attr (MPI_Datatype type, int type_keyval)
3099 if (type==MPI_DATATYPE_NULL)
3100 return MPI_ERR_TYPE;
3102 return smpi_type_attr_delete(type, type_keyval);
3105 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval,
3108 return smpi_type_keyval_create(copy_fn, delete_fn, keyval, extra_state);
3111 int PMPI_Type_free_keyval(int* keyval) {
3112 return smpi_type_keyval_free(keyval);
3115 int PMPI_Info_create( MPI_Info *info){
3118 *info = xbt_new(s_smpi_mpi_info_t, 1);
3119 (*info)->info_dict= xbt_dict_new_homogeneous(NULL);
3120 (*info)->refcount=1;
3124 int PMPI_Info_set( MPI_Info info, char *key, char *value){
3125 if (info == NULL || key == NULL || value == NULL)
3128 xbt_dict_set(info->info_dict, key, (void*)value, NULL);
3132 int PMPI_Info_free( MPI_Info *info){
3133 if (info == NULL || *info==NULL)
3135 (*info)->refcount--;
3136 if((*info)->refcount==0){
3137 xbt_dict_free(&((*info)->info_dict));
3140 *info=MPI_INFO_NULL;
3144 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
3146 if (info == NULL || key == NULL || valuelen <0)
3149 return MPI_ERR_INFO_VALUE;
3150 char* tmpvalue=static_cast<char*>(xbt_dict_get_or_null(info->info_dict, key));
3152 memset(value, 0, valuelen);
3153 memcpy(value,tmpvalue, (strlen(tmpvalue) + 1 < static_cast<size_t>(valuelen)) ? strlen(tmpvalue) + 1 : valuelen);
3159 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
3160 if (info == NULL || newinfo==NULL)
3162 *newinfo = xbt_new(s_smpi_mpi_info_t, 1);
3163 (*newinfo)->info_dict= xbt_dict_new_homogeneous(NULL);
3164 (*newinfo)->refcount=1;
3165 xbt_dict_cursor_t cursor = NULL;
3168 xbt_dict_foreach(info->info_dict,cursor,key,data){
3169 xbt_dict_set((*newinfo)->info_dict, reinterpret_cast<char*>(key), data, NULL);
3174 int PMPI_Info_delete(MPI_Info info, char *key){
3176 if (info == NULL || key==NULL)
3179 xbt_dict_remove(info->info_dict, key);
3182 return MPI_ERR_INFO_NOKEY;
3187 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
3188 if (info == NULL || nkeys==NULL)
3190 *nkeys=xbt_dict_size(info->info_dict);
3194 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
3195 if (info == NULL || key==NULL || n<0 || n> MPI_MAX_INFO_KEY)
3198 xbt_dict_cursor_t cursor = NULL;
3202 xbt_dict_foreach(info->info_dict,cursor,keyn,data){
3204 strncpy(key,keyn,strlen(keyn)+1);
3205 xbt_dict_cursor_free(&cursor);
3213 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
3215 if (info == NULL || key == NULL || valuelen==NULL)
3217 char* tmpvalue=(char*)xbt_dict_get_or_null(info->info_dict, key);
3219 *valuelen=strlen(tmpvalue);
3225 int PMPI_Unpack(void* inbuf, int incount, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
3226 if(incount<0 || outcount < 0 || inbuf==NULL || outbuf==NULL)
3228 if(!is_datatype_valid(type))
3229 return MPI_ERR_TYPE;
3230 if(comm==MPI_COMM_NULL)
3231 return MPI_ERR_COMM;
3232 return smpi_mpi_unpack(inbuf, incount, position, outbuf,outcount,type, comm);
3235 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
3236 if(incount<0 || outcount < 0|| inbuf==NULL || outbuf==NULL)
3238 if(!is_datatype_valid(type))
3239 return MPI_ERR_TYPE;
3240 if(comm==MPI_COMM_NULL)
3241 return MPI_ERR_COMM;
3242 return smpi_mpi_pack(inbuf, incount, type, outbuf,outcount,position, comm);
3245 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
3248 if(!is_datatype_valid(datatype))
3249 return MPI_ERR_TYPE;
3250 if(comm==MPI_COMM_NULL)
3251 return MPI_ERR_COMM;
3253 *size=incount*smpi_datatype_size(datatype);