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);
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 = (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 = (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 = size2 = smpi_group_size(group1);
448 for (i = 0; i < size2; i++) {
449 proc1 = smpi_group_index(group1, i);
450 proc2 = smpi_group_rank(group2, proc1);
451 if (proc2 != MPI_UNDEFINED) {
456 *newgroup = MPI_GROUP_EMPTY;
458 *newgroup = smpi_group_new(size);
459 for (i = 0; i < size2; i++) {
460 proc1 = smpi_group_index(group1, i);
461 proc2 = smpi_group_rank(group2, proc1);
462 if (proc2 == MPI_UNDEFINED) {
463 smpi_group_set_mapping(*newgroup, proc1, i);
467 retval = MPI_SUCCESS;
472 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
476 if (group == MPI_GROUP_NULL) {
477 retval = MPI_ERR_GROUP;
478 } else if (newgroup == NULL) {
479 retval = MPI_ERR_ARG;
481 retval = smpi_group_incl(group, n, ranks, newgroup);
486 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
488 int retval, i, j, newsize, oldsize, index;
490 if (group == MPI_GROUP_NULL) {
491 retval = MPI_ERR_GROUP;
492 } else if (newgroup == NULL) {
493 retval = MPI_ERR_ARG;
497 if(group!= smpi_comm_group(MPI_COMM_WORLD) && group != MPI_GROUP_NULL
498 && group != smpi_comm_group(MPI_COMM_SELF) && group != MPI_GROUP_EMPTY)
499 smpi_group_use(group);
500 } else if (n == smpi_group_size(group)) {
501 *newgroup = MPI_GROUP_EMPTY;
503 oldsize=smpi_group_size(group);
504 newsize = oldsize - n;
505 *newgroup = smpi_group_new(newsize);
507 int* to_exclude=xbt_new0(int, smpi_group_size(group));
508 for(i=0; i<oldsize; i++)
511 to_exclude[ranks[i]]=1;
514 for(i=0; i<oldsize; i++){
515 if(to_exclude[i]==0){
516 index = smpi_group_index(group, i);
517 smpi_group_set_mapping(*newgroup, index, j);
522 xbt_free(to_exclude);
524 retval = MPI_SUCCESS;
529 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3], MPI_Group * newgroup)
531 int retval, i, j, rank, size, index;
533 if (group == MPI_GROUP_NULL) {
534 retval = MPI_ERR_GROUP;
535 } else if (newgroup == NULL) {
536 retval = MPI_ERR_ARG;
539 *newgroup = MPI_GROUP_EMPTY;
542 for (i = 0; i < n; i++) {
543 for (rank = ranges[i][0]; /* First */
544 rank >= 0 && rank < smpi_group_size(group); /* Last */
547 if(rank == ranges[i][1]){/*already last ?*/
550 rank += ranges[i][2]; /* Stride */
551 if (ranges[i][0]<ranges[i][1]){
552 if(rank > ranges[i][1])
555 if(rank < ranges[i][1])
561 *newgroup = smpi_group_new(size);
563 for (i = 0; i < n; i++) {
564 for (rank = ranges[i][0]; /* First */
565 rank >= 0 && rank < smpi_group_size(group); /* Last */
567 index = smpi_group_index(group, rank);
568 smpi_group_set_mapping(*newgroup, index, j);
570 if(rank == ranges[i][1]){/*already last ?*/
573 rank += ranges[i][2]; /* Stride */
574 if (ranges[i][0]<ranges[i][1]){
575 if(rank > ranges[i][1])
578 if(rank < ranges[i][1])
584 retval = MPI_SUCCESS;
589 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3], MPI_Group * newgroup)
591 int retval, i, rank, newrank,oldrank, size, index, add;
593 if (group == MPI_GROUP_NULL) {
594 retval = MPI_ERR_GROUP;
595 } else if (newgroup == NULL) {
596 retval = MPI_ERR_ARG;
600 if(group!= smpi_comm_group(MPI_COMM_WORLD) && group != MPI_GROUP_NULL
601 && group != smpi_comm_group(MPI_COMM_SELF) && group != MPI_GROUP_EMPTY)
602 smpi_group_use(group);
604 size = smpi_group_size(group);
605 for (i = 0; i < n; i++) {
606 for (rank = ranges[i][0]; /* First */
607 rank >= 0 && rank < smpi_group_size(group); /* Last */
610 if(rank == ranges[i][1]){/*already last ?*/
613 rank += ranges[i][2]; /* Stride */
614 if (ranges[i][0]<ranges[i][1]){
615 if(rank > ranges[i][1])
618 if(rank < ranges[i][1])
624 *newgroup = MPI_GROUP_EMPTY;
626 *newgroup = smpi_group_new(size);
629 while (newrank < size) {
631 for (i = 0; i < n; i++) {
632 for (rank = ranges[i][0];
633 rank >= 0 && rank < smpi_group_size(group);
639 if(rank == ranges[i][1]){/*already last ?*/
642 rank += ranges[i][2]; /* Stride */
643 if (ranges[i][0]<ranges[i][1]){
644 if(rank > ranges[i][1])
647 if(rank < ranges[i][1])
653 index = smpi_group_index(group, oldrank);
654 smpi_group_set_mapping(*newgroup, index, newrank);
662 retval = MPI_SUCCESS;
667 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
670 if (comm == MPI_COMM_NULL) {
671 retval = MPI_ERR_COMM;
672 } else if (rank == NULL) {
673 retval = MPI_ERR_ARG;
675 *rank = smpi_comm_rank(comm);
676 retval = MPI_SUCCESS;
681 int PMPI_Comm_size(MPI_Comm comm, int *size)
684 if (comm == MPI_COMM_NULL) {
685 retval = MPI_ERR_COMM;
686 } else if (size == NULL) {
687 retval = MPI_ERR_ARG;
689 *size = smpi_comm_size(comm);
690 retval = MPI_SUCCESS;
695 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
699 if (comm == MPI_COMM_NULL) {
700 retval = MPI_ERR_COMM;
701 } else if (name == NULL || len == NULL) {
702 retval = MPI_ERR_ARG;
704 smpi_comm_get_name(comm, name, len);
705 retval = MPI_SUCCESS;
710 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
714 if (comm == MPI_COMM_NULL) {
715 retval = MPI_ERR_COMM;
716 } else if (group == NULL) {
717 retval = MPI_ERR_ARG;
719 *group = smpi_comm_group(comm);
720 if(*group!= smpi_comm_group(MPI_COMM_WORLD) && *group != MPI_GROUP_NULL
721 && *group != MPI_GROUP_EMPTY)
722 smpi_group_use(*group);
723 retval = MPI_SUCCESS;
728 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
732 if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
733 retval = MPI_ERR_COMM;
734 } else if (result == NULL) {
735 retval = MPI_ERR_ARG;
737 if (comm1 == comm2) { /* Same communicators means same groups */
740 *result = smpi_group_compare(smpi_comm_group(comm1), smpi_comm_group(comm2));
741 if (*result == MPI_IDENT) {
742 *result = MPI_CONGRUENT;
745 retval = MPI_SUCCESS;
750 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
754 if (comm == MPI_COMM_NULL) {
755 retval = MPI_ERR_COMM;
756 } else if (newcomm == NULL) {
757 retval = MPI_ERR_ARG;
759 retval = smpi_comm_dup(comm, newcomm);
764 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
768 if (comm == MPI_COMM_NULL) {
769 retval = MPI_ERR_COMM;
770 } else if (group == MPI_GROUP_NULL) {
771 retval = MPI_ERR_GROUP;
772 } else if (newcomm == NULL) {
773 retval = MPI_ERR_ARG;
774 } else if(smpi_group_rank(group,smpi_process_index())==MPI_UNDEFINED){
775 *newcomm= MPI_COMM_NULL;
776 retval = MPI_SUCCESS;
778 smpi_group_use(group);
779 *newcomm = smpi_comm_new(group, NULL);
780 retval = MPI_SUCCESS;
785 int PMPI_Comm_free(MPI_Comm * comm)
790 retval = MPI_ERR_ARG;
791 } else if (*comm == MPI_COMM_NULL) {
792 retval = MPI_ERR_COMM;
794 smpi_comm_destroy(*comm);
795 *comm = MPI_COMM_NULL;
796 retval = MPI_SUCCESS;
801 int PMPI_Comm_disconnect(MPI_Comm * comm)
803 /* TODO: wait until all communication in comm are done */
807 retval = MPI_ERR_ARG;
808 } else if (*comm == MPI_COMM_NULL) {
809 retval = MPI_ERR_COMM;
811 smpi_comm_destroy(*comm);
812 *comm = MPI_COMM_NULL;
813 retval = MPI_SUCCESS;
818 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
823 if (comm_out == NULL) {
824 retval = MPI_ERR_ARG;
825 } else if (comm == MPI_COMM_NULL) {
826 retval = MPI_ERR_COMM;
828 *comm_out = smpi_comm_split(comm, color, key);
829 retval = MPI_SUCCESS;
836 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
841 if (request == NULL) {
842 retval = MPI_ERR_ARG;
843 } else if (comm == MPI_COMM_NULL) {
844 retval = MPI_ERR_COMM;
845 } else if (!is_datatype_valid(datatype)) {
846 retval = MPI_ERR_TYPE;
847 } else if (dst == MPI_PROC_NULL) {
848 retval = MPI_SUCCESS;
850 *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
851 retval = MPI_SUCCESS;
854 if (retval != MPI_SUCCESS && request)
855 *request = MPI_REQUEST_NULL;
859 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
864 if (request == NULL) {
865 retval = MPI_ERR_ARG;
866 } else if (comm == MPI_COMM_NULL) {
867 retval = MPI_ERR_COMM;
868 } else if (!is_datatype_valid(datatype)) {
869 retval = MPI_ERR_TYPE;
870 } else if (src == MPI_PROC_NULL) {
871 retval = MPI_SUCCESS;
873 *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
874 retval = MPI_SUCCESS;
877 if (retval != MPI_SUCCESS && request)
878 *request = MPI_REQUEST_NULL;
882 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
887 if (request == NULL) {
888 retval = MPI_ERR_ARG;
889 } else if (comm == MPI_COMM_NULL) {
890 retval = MPI_ERR_COMM;
891 } else if (!is_datatype_valid(datatype)) {
892 retval = MPI_ERR_TYPE;
893 } else if (dst == MPI_PROC_NULL) {
894 retval = MPI_SUCCESS;
896 *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
897 retval = MPI_SUCCESS;
900 if (retval != MPI_SUCCESS && request)
901 *request = MPI_REQUEST_NULL;
905 int PMPI_Start(MPI_Request * request)
910 if (request == NULL || *request == MPI_REQUEST_NULL) {
911 retval = MPI_ERR_REQUEST;
913 smpi_mpi_start(*request);
914 retval = MPI_SUCCESS;
920 int PMPI_Startall(int count, MPI_Request * requests)
925 if (requests == NULL) {
926 retval = MPI_ERR_ARG;
928 retval = MPI_SUCCESS;
929 for (i = 0 ; i < count ; i++) {
930 if(requests[i] == MPI_REQUEST_NULL) {
931 retval = MPI_ERR_REQUEST;
934 if(retval != MPI_ERR_REQUEST) {
935 smpi_mpi_startall(count, requests);
942 int PMPI_Request_free(MPI_Request * request)
947 if (*request == MPI_REQUEST_NULL) {
948 retval = MPI_ERR_ARG;
950 smpi_mpi_request_free(request);
951 retval = MPI_SUCCESS;
957 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request * request)
963 if (request == NULL) {
964 retval = MPI_ERR_ARG;
965 } else if (comm == MPI_COMM_NULL) {
966 retval = MPI_ERR_COMM;
967 } else if (src == MPI_PROC_NULL) {
968 *request = MPI_REQUEST_NULL;
969 retval = MPI_SUCCESS;
970 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
971 retval = MPI_ERR_RANK;
972 } else if (count < 0) {
973 retval = MPI_ERR_COUNT;
974 } else if (buf==NULL && count > 0) {
975 retval = MPI_ERR_COUNT;
976 } else if (!is_datatype_valid(datatype)) {
977 retval = MPI_ERR_TYPE;
978 } else if(tag<0 && tag != MPI_ANY_TAG){
979 retval = MPI_ERR_TAG;
982 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
983 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
985 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
986 extra->type = TRACING_IRECV;
987 extra->src = src_traced;
990 extra->datatype1 = encode_datatype(datatype, &known);
991 int dt_size_send = 1;
993 dt_size_send = smpi_datatype_size(datatype);
994 extra->send_size = count*dt_size_send;
995 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
997 *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
998 retval = MPI_SUCCESS;
1000 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1001 (*request)->recv = 1;
1005 if (retval != MPI_SUCCESS && request)
1006 *request = MPI_REQUEST_NULL;
1011 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request * request)
1016 if (request == NULL) {
1017 retval = MPI_ERR_ARG;
1018 } else if (comm == MPI_COMM_NULL) {
1019 retval = MPI_ERR_COMM;
1020 } else if (dst == MPI_PROC_NULL) {
1021 *request = MPI_REQUEST_NULL;
1022 retval = MPI_SUCCESS;
1023 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1024 retval = MPI_ERR_RANK;
1025 } else if (count < 0) {
1026 retval = MPI_ERR_COUNT;
1027 } else if (buf==NULL && count > 0) {
1028 retval = MPI_ERR_COUNT;
1029 } else if (!is_datatype_valid(datatype)) {
1030 retval = MPI_ERR_TYPE;
1031 } else if(tag<0 && tag != MPI_ANY_TAG){
1032 retval = MPI_ERR_TAG;
1035 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1036 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1037 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1038 extra->type = TRACING_ISEND;
1040 extra->dst = dst_traced;
1042 extra->datatype1 = encode_datatype(datatype, &known);
1043 int dt_size_send = 1;
1045 dt_size_send = smpi_datatype_size(datatype);
1046 extra->send_size = count*dt_size_send;
1047 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1048 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1050 *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1051 retval = MPI_SUCCESS;
1053 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1054 (*request)->send = 1;
1058 if (retval != MPI_SUCCESS && request)
1059 *request = MPI_REQUEST_NULL;
1063 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request)
1068 if (request == NULL) {
1069 retval = MPI_ERR_ARG;
1070 } else if (comm == MPI_COMM_NULL) {
1071 retval = MPI_ERR_COMM;
1072 } else if (dst == MPI_PROC_NULL) {
1073 *request = MPI_REQUEST_NULL;
1074 retval = MPI_SUCCESS;
1075 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1076 retval = MPI_ERR_RANK;
1077 } else if (count < 0) {
1078 retval = MPI_ERR_COUNT;
1079 } else if (buf==NULL && count > 0) {
1080 retval = MPI_ERR_COUNT;
1081 } else if (!is_datatype_valid(datatype)) {
1082 retval = MPI_ERR_TYPE;
1083 } else if(tag<0 && tag != MPI_ANY_TAG){
1084 retval = MPI_ERR_TAG;
1087 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1088 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1089 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1090 extra->type = TRACING_ISSEND;
1092 extra->dst = dst_traced;
1094 extra->datatype1 = encode_datatype(datatype, &known);
1095 int dt_size_send = 1;
1097 dt_size_send = smpi_datatype_size(datatype);
1098 extra->send_size = count*dt_size_send;
1099 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1100 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1102 *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1103 retval = MPI_SUCCESS;
1105 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1106 (*request)->send = 1;
1110 if (retval != MPI_SUCCESS && request)
1111 *request = MPI_REQUEST_NULL;
1115 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status * status)
1120 if (comm == MPI_COMM_NULL) {
1121 retval = MPI_ERR_COMM;
1122 } else if (src == MPI_PROC_NULL) {
1123 smpi_empty_status(status);
1124 status->MPI_SOURCE = MPI_PROC_NULL;
1125 retval = MPI_SUCCESS;
1126 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1127 retval = MPI_ERR_RANK;
1128 } else if (count < 0) {
1129 retval = MPI_ERR_COUNT;
1130 } else if (buf==NULL && count > 0) {
1131 retval = MPI_ERR_COUNT;
1132 } else if (!is_datatype_valid(datatype)) {
1133 retval = MPI_ERR_TYPE;
1134 } else if(tag<0 && tag != MPI_ANY_TAG){
1135 retval = MPI_ERR_TAG;
1137 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1138 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1139 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1140 extra->type = TRACING_RECV;
1141 extra->src = src_traced;
1144 extra->datatype1 = encode_datatype(datatype, &known);
1145 int dt_size_send = 1;
1147 dt_size_send = smpi_datatype_size(datatype);
1148 extra->send_size = count*dt_size_send;
1149 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
1151 smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1152 retval = MPI_SUCCESS;
1154 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1155 if(status!=MPI_STATUS_IGNORE){
1156 src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1157 if (!TRACE_smpi_view_internals()) {
1158 TRACE_smpi_recv(rank, src_traced, rank);
1161 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1168 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
1174 if (comm == MPI_COMM_NULL) {
1175 retval = MPI_ERR_COMM;
1176 } else if (dst == MPI_PROC_NULL) {
1177 retval = MPI_SUCCESS;
1178 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1179 retval = MPI_ERR_RANK;
1180 } else if (count < 0) {
1181 retval = MPI_ERR_COUNT;
1182 } else if (buf==NULL && count > 0) {
1183 retval = MPI_ERR_COUNT;
1184 } else if (!is_datatype_valid(datatype)) {
1185 retval = MPI_ERR_TYPE;
1186 } else if(tag<0 && tag != MPI_ANY_TAG){
1187 retval = MPI_ERR_TAG;
1190 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1191 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1192 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1193 extra->type = TRACING_SEND;
1195 extra->dst = dst_traced;
1197 extra->datatype1 = encode_datatype(datatype, &known);
1198 int dt_size_send = 1;
1200 dt_size_send = smpi_datatype_size(datatype);
1201 extra->send_size = count*dt_size_send;
1202 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1203 if (!TRACE_smpi_view_internals()) {
1204 TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1207 smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1208 retval = MPI_SUCCESS;
1210 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1219 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1224 if (comm == MPI_COMM_NULL) {
1225 retval = MPI_ERR_COMM;
1226 } else if (dst == MPI_PROC_NULL) {
1227 retval = MPI_SUCCESS;
1228 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1229 retval = MPI_ERR_RANK;
1230 } else if (count < 0) {
1231 retval = MPI_ERR_COUNT;
1232 } else if (buf==NULL && count > 0) {
1233 retval = MPI_ERR_COUNT;
1234 } else if (!is_datatype_valid(datatype)){
1235 retval = MPI_ERR_TYPE;
1236 } else if(tag<0 && tag != MPI_ANY_TAG){
1237 retval = MPI_ERR_TAG;
1240 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1241 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1242 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1243 extra->type = TRACING_SSEND;
1245 extra->dst = dst_traced;
1247 extra->datatype1 = encode_datatype(datatype, &known);
1248 int dt_size_send = 1;
1250 dt_size_send = smpi_datatype_size(datatype);
1251 extra->send_size = count*dt_size_send;
1252 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1253 TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1255 smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1256 retval = MPI_SUCCESS;
1258 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1264 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype, int dst, int sendtag, void *recvbuf,
1265 int recvcount, MPI_Datatype recvtype, int src, int recvtag, MPI_Comm comm, MPI_Status * status)
1271 if (comm == MPI_COMM_NULL) {
1272 retval = MPI_ERR_COMM;
1273 } else if (!is_datatype_valid(sendtype)
1274 || !is_datatype_valid(recvtype)) {
1275 retval = MPI_ERR_TYPE;
1276 } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1277 smpi_empty_status(status);
1278 status->MPI_SOURCE = MPI_PROC_NULL;
1279 retval = MPI_SUCCESS;
1280 }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1281 (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1282 retval = MPI_ERR_RANK;
1283 } else if (sendcount < 0 || recvcount<0) {
1284 retval = MPI_ERR_COUNT;
1285 } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1286 retval = MPI_ERR_COUNT;
1287 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1288 retval = MPI_ERR_TAG;
1291 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1292 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1293 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1294 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1295 extra->type = TRACING_SENDRECV;
1296 extra->src = src_traced;
1297 extra->dst = dst_traced;
1299 extra->datatype1 = encode_datatype(sendtype, &known);
1300 int dt_size_send = 1;
1302 dt_size_send = smpi_datatype_size(sendtype);
1303 extra->send_size = sendcount*dt_size_send;
1304 extra->datatype2 = encode_datatype(recvtype, &known);
1305 int dt_size_recv = 1;
1307 dt_size_recv = smpi_datatype_size(recvtype);
1308 extra->recv_size = recvcount*dt_size_recv;
1310 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1311 TRACE_smpi_send(rank, rank, dst_traced,sendcount*smpi_datatype_size(sendtype));
1313 smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1314 recvcount, recvtype, src, recvtag, comm, status);
1315 retval = MPI_SUCCESS;
1317 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1318 TRACE_smpi_recv(rank, src_traced, rank);
1325 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype, int dst, int sendtag, int src, int recvtag,
1326 MPI_Comm comm, MPI_Status * status)
1328 //TODO: suboptimal implementation
1331 if (!is_datatype_valid(datatype)) {
1332 retval = MPI_ERR_TYPE;
1333 } else if (count < 0) {
1334 retval = MPI_ERR_COUNT;
1336 int size = smpi_datatype_get_extent(datatype) * count;
1337 recvbuf = xbt_new0(char, size);
1338 retval = MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count, datatype, src, recvtag, comm, status);
1339 if(retval==MPI_SUCCESS){
1340 smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1348 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1352 if (request == NULL || flag == NULL) {
1353 retval = MPI_ERR_ARG;
1354 } else if (*request == MPI_REQUEST_NULL) {
1356 smpi_empty_status(status);
1357 retval = MPI_SUCCESS;
1359 int rank = request && (*request)->comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1361 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1362 extra->type = TRACING_TEST;
1363 TRACE_smpi_testing_in(rank, extra);
1365 *flag = smpi_mpi_test(request, status);
1367 TRACE_smpi_testing_out(rank);
1368 retval = MPI_SUCCESS;
1374 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag, MPI_Status * status)
1379 if (index == NULL || flag == NULL) {
1380 retval = MPI_ERR_ARG;
1382 *flag = smpi_mpi_testany(count, requests, index, status);
1383 retval = MPI_SUCCESS;
1389 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1395 retval = MPI_ERR_ARG;
1397 *flag = smpi_mpi_testall(count, requests, statuses);
1398 retval = MPI_SUCCESS;
1404 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1408 if (status == NULL) {
1409 retval = MPI_ERR_ARG;
1410 } else if (comm == MPI_COMM_NULL) {
1411 retval = MPI_ERR_COMM;
1412 } else if (source == MPI_PROC_NULL) {
1413 smpi_empty_status(status);
1414 status->MPI_SOURCE = MPI_PROC_NULL;
1415 retval = MPI_SUCCESS;
1417 smpi_mpi_probe(source, tag, comm, status);
1418 retval = MPI_SUCCESS;
1424 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1429 retval = MPI_ERR_ARG;
1430 } else if (status == NULL) {
1431 retval = MPI_ERR_ARG;
1432 } else if (comm == MPI_COMM_NULL) {
1433 retval = MPI_ERR_COMM;
1434 } else if (source == MPI_PROC_NULL) {
1436 smpi_empty_status(status);
1437 status->MPI_SOURCE = MPI_PROC_NULL;
1438 retval = MPI_SUCCESS;
1440 smpi_mpi_iprobe(source, tag, comm, flag, status);
1441 retval = MPI_SUCCESS;
1447 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1453 smpi_empty_status(status);
1455 if (request == NULL) {
1456 retval = MPI_ERR_ARG;
1457 } else if (*request == MPI_REQUEST_NULL) {
1458 retval = MPI_SUCCESS;
1461 int rank = request && (*request)->comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1463 int src_traced = (*request)->src;
1464 int dst_traced = (*request)->dst;
1465 MPI_Comm comm = (*request)->comm;
1466 int is_wait_for_receive = (*request)->recv;
1467 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1468 extra->type = TRACING_WAIT;
1469 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1471 smpi_mpi_wait(request, status);
1472 retval = MPI_SUCCESS;
1474 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1475 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1476 if (is_wait_for_receive) {
1477 if(src_traced==MPI_ANY_SOURCE)
1478 src_traced = (status!=MPI_STATUS_IGNORE) ?
1479 smpi_group_rank(smpi_comm_group(comm), status->MPI_SOURCE) :
1481 TRACE_smpi_recv(rank, src_traced, dst_traced);
1489 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1495 //save requests information for tracing
1497 int *srcs = NULL, *dsts = NULL, *recvs = NULL;
1498 MPI_Comm* comms = NULL;
1500 srcs = xbt_new0(int, count);
1501 dsts = xbt_new0(int, count);
1502 recvs = xbt_new0(int, count);
1503 comms = xbt_new0(MPI_Comm, count);
1505 for (i = 0; i < count; i++) {
1506 MPI_Request req = requests[i]; //already received requests are no longer valid
1510 recvs[i] = req->recv;
1511 comms[i] = req->comm;
1514 int rank_traced = smpi_process_index();
1515 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1516 extra->type = TRACING_WAITANY;
1517 extra->send_size=count;
1518 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1520 *index = smpi_mpi_waitany(count, requests, status);
1522 if(*index!=MPI_UNDEFINED){
1523 int src_traced = srcs[*index];
1524 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1525 int dst_traced = dsts[*index];
1526 int is_wait_for_receive = recvs[*index];
1527 if (is_wait_for_receive) {
1528 if(srcs[*index]==MPI_ANY_SOURCE)
1529 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1530 smpi_group_rank(smpi_comm_group(comms[*index]), status->MPI_SOURCE) : srcs[*index];
1531 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1533 TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1544 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1547 //save information from requests
1549 int *srcs = xbt_new0(int, count);
1550 int *dsts = xbt_new0(int, count);
1551 int *recvs = xbt_new0(int, count);
1552 int *valid = xbt_new0(int, count);
1553 MPI_Comm *comms = xbt_new0(MPI_Comm, count);
1555 //int valid_count = 0;
1556 for (i = 0; i < count; i++) {
1557 MPI_Request req = requests[i];
1558 if(req!=MPI_REQUEST_NULL){
1561 recvs[i] = req->recv;
1562 comms[i] = req->comm;
1568 int rank_traced = smpi_process_index();
1569 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1570 extra->type = TRACING_WAITALL;
1571 extra->send_size=count;
1572 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1574 int retval = smpi_mpi_waitall(count, requests, status);
1576 for (i = 0; i < count; i++) {
1578 //int src_traced = srcs[*index];
1579 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1580 int src_traced = srcs[i];
1581 int dst_traced = dsts[i];
1582 int is_wait_for_receive = recvs[i];
1583 if (is_wait_for_receive) {
1584 if(src_traced==MPI_ANY_SOURCE)
1585 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1586 smpi_group_rank(smpi_comm_group(comms[i]), status[i].MPI_SOURCE) : srcs[i];
1587 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1591 TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1602 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount, int *indices, MPI_Status status[])
1607 if (outcount == NULL) {
1608 retval = MPI_ERR_ARG;
1610 *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1611 retval = MPI_SUCCESS;
1617 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount, int* indices, MPI_Status status[])
1622 if (outcount == NULL) {
1623 retval = MPI_ERR_ARG;
1625 *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1626 retval = MPI_SUCCESS;
1633 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1639 if (comm == MPI_COMM_NULL) {
1640 retval = MPI_ERR_COMM;
1641 } else if (!is_datatype_valid(datatype)) {
1642 retval = MPI_ERR_ARG;
1644 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1645 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1647 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1648 extra->type = TRACING_BCAST;
1649 extra->root = root_traced;
1651 extra->datatype1 = encode_datatype(datatype, &known);
1652 int dt_size_send = 1;
1654 dt_size_send = smpi_datatype_size(datatype);
1655 extra->send_size = count*dt_size_send;
1656 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1658 mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1659 retval = MPI_SUCCESS;
1661 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1668 int PMPI_Barrier(MPI_Comm comm)
1674 if (comm == MPI_COMM_NULL) {
1675 retval = MPI_ERR_COMM;
1677 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1678 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1679 extra->type = TRACING_BARRIER;
1680 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1682 mpi_coll_barrier_fun(comm);
1683 retval = MPI_SUCCESS;
1685 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1692 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,void *recvbuf, int recvcount, MPI_Datatype recvtype,
1693 int root, MPI_Comm comm)
1699 if (comm == MPI_COMM_NULL) {
1700 retval = MPI_ERR_COMM;
1701 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1702 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1703 retval = MPI_ERR_TYPE;
1704 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) || ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1705 retval = MPI_ERR_COUNT;
1708 char* sendtmpbuf = (char*) sendbuf;
1709 int sendtmpcount = sendcount;
1710 MPI_Datatype sendtmptype = sendtype;
1711 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1713 sendtmptype=recvtype;
1715 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1716 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1717 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1718 extra->type = TRACING_GATHER;
1719 extra->root = root_traced;
1721 extra->datatype1 = encode_datatype(sendtmptype, &known);
1722 int dt_size_send = 1;
1724 dt_size_send = smpi_datatype_size(sendtmptype);
1725 extra->send_size = sendtmpcount*dt_size_send;
1726 extra->datatype2 = encode_datatype(recvtype, &known);
1727 int dt_size_recv = 1;
1728 if((smpi_comm_rank(comm)==root) && !known)
1729 dt_size_recv = smpi_datatype_size(recvtype);
1730 extra->recv_size = recvcount*dt_size_recv;
1732 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1734 mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount, recvtype, root, comm);
1736 retval = MPI_SUCCESS;
1737 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1744 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int *recvcounts, int *displs,
1745 MPI_Datatype recvtype, int root, MPI_Comm comm)
1751 if (comm == MPI_COMM_NULL) {
1752 retval = MPI_ERR_COMM;
1753 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1754 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1755 retval = MPI_ERR_TYPE;
1756 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1757 retval = MPI_ERR_COUNT;
1758 } else if (recvcounts == NULL || displs == NULL) {
1759 retval = MPI_ERR_ARG;
1761 char* sendtmpbuf = (char*) sendbuf;
1762 int sendtmpcount = sendcount;
1763 MPI_Datatype sendtmptype = sendtype;
1764 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1766 sendtmptype=recvtype;
1769 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1770 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1772 int size = smpi_comm_size(comm);
1773 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1774 extra->type = TRACING_GATHERV;
1775 extra->num_processes = size;
1776 extra->root = root_traced;
1778 extra->datatype1 = encode_datatype(sendtmptype, &known);
1779 int dt_size_send = 1;
1781 dt_size_send = smpi_datatype_size(sendtype);
1782 extra->send_size = sendtmpcount*dt_size_send;
1783 extra->datatype2 = encode_datatype(recvtype, &known);
1784 int dt_size_recv = 1;
1786 dt_size_recv = smpi_datatype_size(recvtype);
1787 if((smpi_comm_rank(comm)==root)){
1788 extra->recvcounts= xbt_new(int,size);
1789 for(i=0; i< size; i++)//copy data to avoid bad free
1790 extra->recvcounts[i] = recvcounts[i]*dt_size_recv;
1792 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1794 smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts, displs, recvtype, root, comm);
1795 retval = MPI_SUCCESS;
1796 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1803 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1804 void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
1810 if (comm == MPI_COMM_NULL) {
1811 retval = MPI_ERR_COMM;
1812 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1813 (recvtype == MPI_DATATYPE_NULL)){
1814 retval = MPI_ERR_TYPE;
1815 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1817 retval = MPI_ERR_COUNT;
1819 if(sendbuf == MPI_IN_PLACE) {
1820 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1821 sendcount=recvcount;
1824 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1825 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1826 extra->type = TRACING_ALLGATHER;
1828 extra->datatype1 = encode_datatype(sendtype, &known);
1829 int dt_size_send = 1;
1831 dt_size_send = smpi_datatype_size(sendtype);
1832 extra->send_size = sendcount*dt_size_send;
1833 extra->datatype2 = encode_datatype(recvtype, &known);
1834 int dt_size_recv = 1;
1836 dt_size_recv = smpi_datatype_size(recvtype);
1837 extra->recv_size = recvcount*dt_size_recv;
1839 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1841 mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
1842 retval = MPI_SUCCESS;
1843 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1849 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1850 void *recvbuf, int *recvcounts, int *displs, MPI_Datatype recvtype, MPI_Comm comm)
1856 if (comm == MPI_COMM_NULL) {
1857 retval = MPI_ERR_COMM;
1858 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1859 (recvtype == MPI_DATATYPE_NULL)){
1860 retval = MPI_ERR_TYPE;
1861 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1862 retval = MPI_ERR_COUNT;
1863 } else if (recvcounts == NULL || displs == NULL) {
1864 retval = MPI_ERR_ARG;
1867 if(sendbuf == MPI_IN_PLACE) {
1868 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1869 sendcount=recvcounts[smpi_comm_rank(comm)];
1872 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1874 int size = smpi_comm_size(comm);
1875 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1876 extra->type = TRACING_ALLGATHERV;
1877 extra->num_processes = size;
1879 extra->datatype1 = encode_datatype(sendtype, &known);
1880 int dt_size_send = 1;
1882 dt_size_send = smpi_datatype_size(sendtype);
1883 extra->send_size = sendcount*dt_size_send;
1884 extra->datatype2 = encode_datatype(recvtype, &known);
1885 int dt_size_recv = 1;
1887 dt_size_recv = smpi_datatype_size(recvtype);
1888 extra->recvcounts= xbt_new(int, size);
1889 for(i=0; i< size; i++)//copy data to avoid bad free
1890 extra->recvcounts[i] = recvcounts[i]*dt_size_recv;
1892 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
1894 mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm);
1895 retval = MPI_SUCCESS;
1896 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1903 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1904 void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
1910 if (comm == MPI_COMM_NULL) {
1911 retval = MPI_ERR_COMM;
1912 } else if (((smpi_comm_rank(comm)==root) && (!is_datatype_valid(sendtype)))
1913 || ((recvbuf !=MPI_IN_PLACE) && (!is_datatype_valid(recvtype)))){
1914 retval = MPI_ERR_TYPE;
1915 } else if ((sendbuf == recvbuf) ||
1916 ((smpi_comm_rank(comm)==root) && sendcount>0 && (sendbuf == NULL))){
1917 retval = MPI_ERR_BUFFER;
1920 if (recvbuf == MPI_IN_PLACE) {
1922 recvcount=sendcount;
1924 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1925 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1926 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1927 extra->type = TRACING_SCATTER;
1928 extra->root = root_traced;
1930 extra->datatype1 = encode_datatype(sendtype, &known);
1931 int dt_size_send = 1;
1932 if((smpi_comm_rank(comm)==root) && !known)
1933 dt_size_send = smpi_datatype_size(sendtype);
1934 extra->send_size = sendcount*dt_size_send;
1935 extra->datatype2 = encode_datatype(recvtype, &known);
1936 int dt_size_recv = 1;
1938 dt_size_recv = smpi_datatype_size(recvtype);
1939 extra->recv_size = recvcount*dt_size_recv;
1940 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1942 mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, root, comm);
1943 retval = MPI_SUCCESS;
1944 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1951 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1952 MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm)
1958 if (comm == MPI_COMM_NULL) {
1959 retval = MPI_ERR_COMM;
1960 } else if (sendcounts == NULL || displs == NULL) {
1961 retval = MPI_ERR_ARG;
1962 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1963 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1964 retval = MPI_ERR_TYPE;
1966 if (recvbuf == MPI_IN_PLACE) {
1968 recvcount=sendcounts[smpi_comm_rank(comm)];
1970 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1971 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1973 int size = smpi_comm_size(comm);
1974 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1975 extra->type = TRACING_SCATTERV;
1976 extra->num_processes = size;
1977 extra->root = root_traced;
1979 extra->datatype1 = encode_datatype(sendtype, &known);
1980 int dt_size_send = 1;
1982 dt_size_send = smpi_datatype_size(sendtype);
1983 if((smpi_comm_rank(comm)==root)){
1984 extra->sendcounts= xbt_new(int, size);
1985 for(i=0; i< size; i++)//copy data to avoid bad free
1986 extra->sendcounts[i] = sendcounts[i]*dt_size_send;
1988 extra->datatype2 = encode_datatype(recvtype, &known);
1989 int dt_size_recv = 1;
1991 dt_size_recv = smpi_datatype_size(recvtype);
1992 extra->recv_size = recvcount*dt_size_recv;
1993 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1995 smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf, recvcount, recvtype, root, comm);
1997 retval = MPI_SUCCESS;
1998 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2005 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
2011 if (comm == MPI_COMM_NULL) {
2012 retval = MPI_ERR_COMM;
2013 } else if (!is_datatype_valid(datatype) || op == MPI_OP_NULL) {
2014 retval = MPI_ERR_ARG;
2016 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2017 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2018 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2019 extra->type = TRACING_REDUCE;
2021 extra->datatype1 = encode_datatype(datatype, &known);
2022 int dt_size_send = 1;
2024 dt_size_send = smpi_datatype_size(datatype);
2025 extra->send_size = count*dt_size_send;
2026 extra->root = root_traced;
2028 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2030 mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
2032 retval = MPI_SUCCESS;
2033 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2040 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count, MPI_Datatype datatype, MPI_Op op){
2044 if (!is_datatype_valid(datatype) || op == MPI_OP_NULL) {
2045 retval = MPI_ERR_ARG;
2047 smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
2054 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2060 if (comm == MPI_COMM_NULL) {
2061 retval = MPI_ERR_COMM;
2062 } else if (!is_datatype_valid(datatype)) {
2063 retval = MPI_ERR_TYPE;
2064 } else if (op == MPI_OP_NULL) {
2065 retval = MPI_ERR_OP;
2068 char* sendtmpbuf = (char*) sendbuf;
2069 if( sendbuf == MPI_IN_PLACE ) {
2070 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
2071 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
2073 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2074 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2075 extra->type = TRACING_ALLREDUCE;
2077 extra->datatype1 = encode_datatype(datatype, &known);
2078 int dt_size_send = 1;
2080 dt_size_send = smpi_datatype_size(datatype);
2081 extra->send_size = count*dt_size_send;
2083 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2085 mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
2087 if( sendbuf == MPI_IN_PLACE )
2088 xbt_free(sendtmpbuf);
2090 retval = MPI_SUCCESS;
2091 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2098 int PMPI_Scan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2104 if (comm == MPI_COMM_NULL) {
2105 retval = MPI_ERR_COMM;
2106 } else if (!is_datatype_valid(datatype)) {
2107 retval = MPI_ERR_TYPE;
2108 } else if (op == MPI_OP_NULL) {
2109 retval = MPI_ERR_OP;
2111 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2112 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2113 extra->type = TRACING_SCAN;
2115 extra->datatype1 = encode_datatype(datatype, &known);
2116 int dt_size_send = 1;
2118 dt_size_send = smpi_datatype_size(datatype);
2119 extra->send_size = count*dt_size_send;
2121 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2123 smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2125 retval = MPI_SUCCESS;
2126 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2133 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm){
2138 if (comm == MPI_COMM_NULL) {
2139 retval = MPI_ERR_COMM;
2140 } else if (!is_datatype_valid(datatype)) {
2141 retval = MPI_ERR_TYPE;
2142 } else if (op == MPI_OP_NULL) {
2143 retval = MPI_ERR_OP;
2145 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2146 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2147 extra->type = TRACING_EXSCAN;
2149 extra->datatype1 = encode_datatype(datatype, &known);
2150 int dt_size_send = 1;
2152 dt_size_send = smpi_datatype_size(datatype);
2153 extra->send_size = count*dt_size_send;
2154 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2156 smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2157 retval = MPI_SUCCESS;
2158 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2165 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2170 if (comm == MPI_COMM_NULL) {
2171 retval = MPI_ERR_COMM;
2172 } else if (!is_datatype_valid(datatype)) {
2173 retval = MPI_ERR_TYPE;
2174 } else if (op == MPI_OP_NULL) {
2175 retval = MPI_ERR_OP;
2176 } else if (recvcounts == NULL) {
2177 retval = MPI_ERR_ARG;
2179 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2181 int size = smpi_comm_size(comm);
2182 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2183 extra->type = TRACING_REDUCE_SCATTER;
2184 extra->num_processes = size;
2186 extra->datatype1 = encode_datatype(datatype, &known);
2187 int dt_size_send = 1;
2189 dt_size_send = smpi_datatype_size(datatype);
2190 extra->send_size = 0;
2191 extra->recvcounts= xbt_new(int, size);
2192 for(i=0; i< size; i++)//copy data to avoid bad free
2193 extra->recvcounts[i] = recvcounts[i]*dt_size_send;
2194 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2196 void* sendtmpbuf=sendbuf;
2197 if(sendbuf==MPI_IN_PLACE)
2200 mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts, datatype, op, comm);
2201 retval = MPI_SUCCESS;
2202 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2209 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2210 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2215 if (comm == MPI_COMM_NULL) {
2216 retval = MPI_ERR_COMM;
2217 } else if (!is_datatype_valid(datatype)) {
2218 retval = MPI_ERR_TYPE;
2219 } else if (op == MPI_OP_NULL) {
2220 retval = MPI_ERR_OP;
2221 } else if (recvcount < 0) {
2222 retval = MPI_ERR_ARG;
2224 int count=smpi_comm_size(comm);
2226 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2227 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2228 extra->type = TRACING_REDUCE_SCATTER;
2229 extra->num_processes = count;
2231 extra->datatype1 = encode_datatype(datatype, &known);
2232 int dt_size_send = 1;
2234 dt_size_send = smpi_datatype_size(datatype);
2235 extra->send_size = 0;
2236 extra->recvcounts= xbt_new(int, count);
2237 for(i=0; i< count; i++)//copy data to avoid bad free
2238 extra->recvcounts[i] = recvcount*dt_size_send;
2240 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2242 int* recvcounts=(int*)xbt_malloc(count);
2243 for (i=0; i<count;i++)recvcounts[i]=recvcount;
2244 mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts, datatype, op, comm);
2245 xbt_free(recvcounts);
2246 retval = MPI_SUCCESS;
2248 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2255 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2256 void *recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm)
2262 if (comm == MPI_COMM_NULL) {
2263 retval = MPI_ERR_COMM;
2264 } else if (sendtype == MPI_DATATYPE_NULL
2265 || recvtype == MPI_DATATYPE_NULL) {
2266 retval = MPI_ERR_TYPE;
2268 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2269 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2270 extra->type = TRACING_ALLTOALL;
2272 extra->datatype1 = encode_datatype(sendtype, &known);
2274 extra->send_size = sendcount*smpi_datatype_size(sendtype);
2276 extra->send_size = sendcount;
2277 extra->datatype2 = encode_datatype(recvtype, &known);
2279 extra->recv_size = recvcount*smpi_datatype_size(recvtype);
2281 extra->recv_size = recvcount;
2282 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2284 retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2286 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2293 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2294 int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2300 if (comm == MPI_COMM_NULL) {
2301 retval = MPI_ERR_COMM;
2302 } else if (sendtype == MPI_DATATYPE_NULL || recvtype == MPI_DATATYPE_NULL) {
2303 retval = MPI_ERR_TYPE;
2304 } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL || recvdisps == NULL) {
2305 retval = MPI_ERR_ARG;
2307 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2309 int size = smpi_comm_size(comm);
2310 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2311 extra->type = TRACING_ALLTOALLV;
2312 extra->send_size = 0;
2313 extra->recv_size = 0;
2314 extra->recvcounts= xbt_new(int, size);
2315 extra->sendcounts= xbt_new(int, size);
2317 extra->datatype1 = encode_datatype(sendtype, &known);
2318 int dt_size_send = 1;
2320 dt_size_send = smpi_datatype_size(sendtype);
2321 int dt_size_recv = 1;
2322 extra->datatype2 = encode_datatype(recvtype, &known);
2324 dt_size_recv = smpi_datatype_size(recvtype);
2325 for(i=0; i< size; i++){//copy data to avoid bad free
2326 extra->send_size += sendcounts[i]*dt_size_send;
2327 extra->recv_size += recvcounts[i]*dt_size_recv;
2329 extra->sendcounts[i] = sendcounts[i]*dt_size_send;
2330 extra->recvcounts[i] = recvcounts[i]*dt_size_recv;
2332 extra->num_processes = size;
2333 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2335 retval = mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype, recvbuf, recvcounts, recvdisps, recvtype,
2337 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2345 int PMPI_Get_processor_name(char *name, int *resultlen)
2347 int retval = MPI_SUCCESS;
2349 strncpy(name, sg_host_get_name(SIMIX_host_self()),
2350 strlen(sg_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2351 strlen(sg_host_get_name(SIMIX_host_self())) +1 : MPI_MAX_PROCESSOR_NAME - 1 );
2352 *resultlen = strlen(name) > MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2357 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2359 int retval = MPI_SUCCESS;
2362 if (status == NULL || count == NULL) {
2363 retval = MPI_ERR_ARG;
2364 } else if (!is_datatype_valid(datatype)) {
2365 retval = MPI_ERR_TYPE;
2367 size = smpi_datatype_size(datatype);
2370 } else if (status->count % size != 0) {
2371 retval = MPI_UNDEFINED;
2373 *count = smpi_mpi_get_count(status, datatype);
2379 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2382 if (old_type == MPI_DATATYPE_NULL) {
2383 retval = MPI_ERR_TYPE;
2384 } else if (count<0){
2385 retval = MPI_ERR_COUNT;
2387 retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2392 int PMPI_Type_commit(MPI_Datatype* datatype) {
2395 if (datatype == NULL || *datatype == MPI_DATATYPE_NULL) {
2396 retval = MPI_ERR_TYPE;
2398 smpi_datatype_commit(datatype);
2399 retval = MPI_SUCCESS;
2404 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2407 if (old_type == MPI_DATATYPE_NULL) {
2408 retval = MPI_ERR_TYPE;
2409 } else if (count<0 || blocklen<0){
2410 retval = MPI_ERR_COUNT;
2412 retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2417 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2420 if (old_type == MPI_DATATYPE_NULL) {
2421 retval = MPI_ERR_TYPE;
2422 } else if (count<0 || blocklen<0){
2423 retval = MPI_ERR_COUNT;
2425 retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2430 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2431 return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2434 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2437 if (old_type == MPI_DATATYPE_NULL) {
2438 retval = MPI_ERR_TYPE;
2439 } else if (count<0){
2440 retval = MPI_ERR_COUNT;
2442 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2447 int PMPI_Type_create_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2450 if (old_type == MPI_DATATYPE_NULL) {
2451 retval = MPI_ERR_TYPE;
2452 } else if (count<0){
2453 retval = MPI_ERR_COUNT;
2455 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2460 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type,
2461 MPI_Datatype* new_type) {
2464 if (old_type == MPI_DATATYPE_NULL) {
2465 retval = MPI_ERR_TYPE;
2466 } else if (count<0){
2467 retval = MPI_ERR_COUNT;
2469 int* blocklens=(int*)xbt_malloc(blocklength*count);
2470 for (i=0; i<count;i++)blocklens[i]=blocklength;
2471 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2472 xbt_free(blocklens);
2477 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2480 if (old_type == MPI_DATATYPE_NULL) {
2481 retval = MPI_ERR_TYPE;
2482 } else if (count<0){
2483 retval = MPI_ERR_COUNT;
2485 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2490 int PMPI_Type_create_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type,
2491 MPI_Datatype* new_type) {
2492 return PMPI_Type_hindexed(count, blocklens,indices,old_type,new_type);
2495 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type,
2496 MPI_Datatype* new_type) {
2499 if (old_type == MPI_DATATYPE_NULL) {
2500 retval = MPI_ERR_TYPE;
2501 } else if (count<0){
2502 retval = MPI_ERR_COUNT;
2504 int* blocklens=(int*)xbt_malloc(blocklength*count);
2505 for (i=0; i<count;i++)blocklens[i]=blocklength;
2506 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2507 xbt_free(blocklens);
2512 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2516 retval = MPI_ERR_COUNT;
2518 retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2523 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types,
2524 MPI_Datatype* new_type) {
2525 return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2528 int PMPI_Error_class(int errorcode, int* errorclass) {
2529 // assume smpi uses only standard mpi error codes
2530 *errorclass=errorcode;
2534 int PMPI_Initialized(int* flag) {
2535 *flag=smpi_process_initialized();
2539 /* The topo part of MPI_COMM_WORLD should always be NULL. When other topologies will be implemented, not only should we
2540 * check if the topology is NULL, but we should check if it is the good topology type (so we have to add a
2541 * MPIR_Topo_Type field, and replace the MPI_Topology field by an union)*/
2543 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periodic, int reorder, MPI_Comm* comm_cart) {
2545 if (comm_old == MPI_COMM_NULL){
2546 retval = MPI_ERR_COMM;
2547 } else if (ndims < 0 || (ndims > 0 && (dims == NULL || periodic == NULL)) || comm_cart == NULL) {
2548 retval = MPI_ERR_ARG;
2550 retval = smpi_mpi_cart_create(comm_old, ndims, dims, periodic, reorder, comm_cart);
2555 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2556 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2557 return MPI_ERR_TOPOLOGY;
2559 if (coords == NULL) {
2562 return smpi_mpi_cart_rank(comm, coords, rank);
2565 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2566 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2567 return MPI_ERR_TOPOLOGY;
2569 if (source == NULL || dest == NULL || direction < 0 ) {
2572 return smpi_mpi_cart_shift(comm, direction, displ, source, dest);
2575 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2576 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2577 return MPI_ERR_TOPOLOGY;
2579 if (rank < 0 || rank >= smpi_comm_size(comm)) {
2580 return MPI_ERR_RANK;
2585 if(coords == NULL) {
2588 return smpi_mpi_cart_coords(comm, rank, maxdims, coords);
2591 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2592 if(comm == NULL || smpi_comm_topo(comm) == NULL) {
2593 return MPI_ERR_TOPOLOGY;
2595 if(maxdims <= 0 || dims == NULL || periods == NULL || coords == NULL) {
2598 return smpi_mpi_cart_get(comm, maxdims, dims, periods, coords);
2601 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2602 if (comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2603 return MPI_ERR_TOPOLOGY;
2605 if (ndims == NULL) {
2608 return smpi_mpi_cartdim_get(comm, ndims);
2611 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2615 if (ndims < 1 || nnodes < 1) {
2616 return MPI_ERR_DIMS;
2619 return smpi_mpi_dims_create(nnodes, ndims, dims);
2622 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2623 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2624 return MPI_ERR_TOPOLOGY;
2626 if (comm_new == NULL) {
2629 return smpi_mpi_cart_sub(comm, remain_dims, comm_new);
2632 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2633 if(oldtype == MPI_DATATYPE_NULL) {
2634 return MPI_ERR_TYPE;
2636 int blocks[3] = { 1, 1, 1 };
2637 MPI_Aint disps[3] = { lb, 0, lb+extent };
2638 MPI_Datatype types[3] = { MPI_LB, oldtype, MPI_UB };
2640 s_smpi_mpi_struct_t* subtype = smpi_datatype_struct_create( blocks, disps, 3, types);
2641 smpi_datatype_create(newtype,oldtype->size, lb, lb + extent, sizeof(s_smpi_mpi_struct_t) , subtype, DT_FLAG_VECTOR);
2643 (*newtype)->flags &= ~DT_FLAG_COMMITED;
2647 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2650 if (comm == MPI_COMM_NULL) {
2651 retval= MPI_ERR_COMM;
2652 }else if ((base == NULL && size != 0) || disp_unit <= 0 || size < 0 ){
2653 retval= MPI_ERR_OTHER;
2655 *win = smpi_mpi_win_create( base, size, disp_unit, info, comm);
2656 retval = MPI_SUCCESS;
2662 int PMPI_Win_free( MPI_Win* win){
2665 if (win == NULL || *win == MPI_WIN_NULL) {
2666 retval = MPI_ERR_WIN;
2668 retval=smpi_mpi_win_free(win);
2674 int PMPI_Win_set_name(MPI_Win win, char * name)
2677 if (win == MPI_WIN_NULL) {
2678 retval = MPI_ERR_TYPE;
2679 } else if (name == NULL) {
2680 retval = MPI_ERR_ARG;
2682 smpi_mpi_win_set_name(win, name);
2683 retval = MPI_SUCCESS;
2688 int PMPI_Win_get_name(MPI_Win win, char * name, int* len)
2690 int retval = MPI_SUCCESS;
2692 if (win == MPI_WIN_NULL) {
2693 retval = MPI_ERR_WIN;
2694 } else if (name == NULL) {
2695 retval = MPI_ERR_ARG;
2697 smpi_mpi_win_get_name(win, name, len);
2702 int PMPI_Win_get_group(MPI_Win win, MPI_Group * group){
2703 int retval = MPI_SUCCESS;
2704 if (win == MPI_WIN_NULL) {
2705 retval = MPI_ERR_WIN;
2707 smpi_mpi_win_get_group(win, group);
2713 int PMPI_Win_fence( int assert, MPI_Win win){
2716 if (win == MPI_WIN_NULL) {
2717 retval = MPI_ERR_WIN;
2719 int rank = smpi_process_index();
2720 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, NULL);
2721 retval = smpi_mpi_win_fence(assert, win);
2722 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2728 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2729 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2732 if (win == MPI_WIN_NULL) {
2733 retval = MPI_ERR_WIN;
2734 } else if (target_rank == MPI_PROC_NULL) {
2735 retval = MPI_SUCCESS;
2736 } else if (target_rank <0){
2737 retval = MPI_ERR_RANK;
2738 } else if (target_disp <0){
2739 retval = MPI_ERR_ARG;
2740 } else if (origin_count < 0 || target_count < 0) {
2741 retval = MPI_ERR_COUNT;
2742 } else if (origin_addr==NULL && origin_count > 0){
2743 retval = MPI_ERR_COUNT;
2744 } else if ((!is_datatype_valid(origin_datatype)) || (!is_datatype_valid(target_datatype))) {
2745 retval = MPI_ERR_TYPE;
2747 int rank = smpi_process_index();
2749 smpi_mpi_win_get_group(win, &group);
2750 int src_traced = smpi_group_index(group, target_rank);
2751 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, NULL);
2753 retval = smpi_mpi_get( origin_addr, origin_count, origin_datatype, target_rank, target_disp, target_count,
2754 target_datatype, win);
2756 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
2762 int PMPI_Put( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2763 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2766 if (win == MPI_WIN_NULL) {
2767 retval = MPI_ERR_WIN;
2768 } else if (target_rank == MPI_PROC_NULL) {
2769 retval = MPI_SUCCESS;
2770 } else if (target_rank <0){
2771 retval = MPI_ERR_RANK;
2772 } else if (target_disp <0){
2773 retval = MPI_ERR_ARG;
2774 } else if (origin_count < 0 || target_count < 0) {
2775 retval = MPI_ERR_COUNT;
2776 } else if (origin_addr==NULL && origin_count > 0){
2777 retval = MPI_ERR_COUNT;
2778 } else if ((!is_datatype_valid(origin_datatype)) || (!is_datatype_valid(target_datatype))) {
2779 retval = MPI_ERR_TYPE;
2781 int rank = smpi_process_index();
2783 smpi_mpi_win_get_group(win, &group);
2784 int dst_traced = smpi_group_index(group, target_rank);
2785 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, NULL);
2786 TRACE_smpi_send(rank, rank, dst_traced, origin_count*smpi_datatype_size(origin_datatype));
2788 retval = smpi_mpi_put( origin_addr, origin_count, origin_datatype, target_rank, target_disp, target_count,
2789 target_datatype, win);
2791 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
2797 int PMPI_Accumulate( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2798 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Op op, MPI_Win win){
2801 if (win == MPI_WIN_NULL) {
2802 retval = MPI_ERR_WIN;
2803 } else if (target_rank == MPI_PROC_NULL) {
2804 retval = MPI_SUCCESS;
2805 } else if (target_rank <0){
2806 retval = MPI_ERR_RANK;
2807 } else if (target_disp <0){
2808 retval = MPI_ERR_ARG;
2809 } else if (origin_count < 0 || target_count < 0) {
2810 retval = MPI_ERR_COUNT;
2811 } else if (origin_addr==NULL && origin_count > 0){
2812 retval = MPI_ERR_COUNT;
2813 } else if ((!is_datatype_valid(origin_datatype)) ||
2814 (!is_datatype_valid(target_datatype))) {
2815 retval = MPI_ERR_TYPE;
2816 } else if (op == MPI_OP_NULL) {
2817 retval = MPI_ERR_OP;
2819 int rank = smpi_process_index();
2821 smpi_mpi_win_get_group(win, &group);
2822 int src_traced = smpi_group_index(group, target_rank);
2823 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, NULL);
2825 retval = smpi_mpi_accumulate( origin_addr, origin_count, origin_datatype, target_rank, target_disp, target_count,
2826 target_datatype, op, win);
2828 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
2834 int PMPI_Win_post(MPI_Group group, int assert, MPI_Win win){
2837 if (win == MPI_WIN_NULL) {
2838 retval = MPI_ERR_WIN;
2839 } else if (group==MPI_GROUP_NULL){
2840 retval = MPI_ERR_GROUP;
2843 int rank = smpi_process_index();
2844 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, NULL);
2845 retval = smpi_mpi_win_post(group,assert,win);
2846 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2852 int PMPI_Win_start(MPI_Group group, int assert, MPI_Win win){
2855 if (win == MPI_WIN_NULL) {
2856 retval = MPI_ERR_WIN;
2857 } else if (group==MPI_GROUP_NULL){
2858 retval = MPI_ERR_GROUP;
2861 int rank = smpi_process_index();
2862 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, NULL);
2863 retval = smpi_mpi_win_start(group,assert,win);
2864 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2870 int PMPI_Win_complete(MPI_Win win){
2873 if (win == MPI_WIN_NULL) {
2874 retval = MPI_ERR_WIN;
2877 int rank = smpi_process_index();
2878 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, NULL);
2880 retval = smpi_mpi_win_complete(win);
2882 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2888 int PMPI_Win_wait(MPI_Win win){
2891 if (win == MPI_WIN_NULL) {
2892 retval = MPI_ERR_WIN;
2895 int rank = smpi_process_index();
2896 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, NULL);
2898 retval = smpi_mpi_win_wait(win);
2900 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2906 int PMPI_Alloc_mem(MPI_Aint size, MPI_Info info, void *baseptr){
2907 void *ptr = xbt_malloc(size);
2909 return MPI_ERR_NO_MEM;
2911 *(void **)baseptr = ptr;
2916 int PMPI_Free_mem(void *baseptr){
2921 int PMPI_Type_set_name(MPI_Datatype datatype, char * name)
2924 if (datatype == MPI_DATATYPE_NULL) {
2925 retval = MPI_ERR_TYPE;
2926 } else if (name == NULL) {
2927 retval = MPI_ERR_ARG;
2929 smpi_datatype_set_name(datatype, name);
2930 retval = MPI_SUCCESS;
2935 int PMPI_Type_get_name(MPI_Datatype datatype, char * name, int* len)
2939 if (datatype == MPI_DATATYPE_NULL) {
2940 retval = MPI_ERR_TYPE;
2941 } else if (name == NULL) {
2942 retval = MPI_ERR_ARG;
2944 smpi_datatype_get_name(datatype, name, len);
2945 retval = MPI_SUCCESS;
2950 MPI_Datatype PMPI_Type_f2c(MPI_Fint datatype){
2951 return smpi_type_f2c(datatype);
2954 MPI_Fint PMPI_Type_c2f(MPI_Datatype datatype){
2955 return smpi_type_c2f( datatype);
2958 MPI_Group PMPI_Group_f2c(MPI_Fint group){
2959 return smpi_group_f2c( group);
2962 MPI_Fint PMPI_Group_c2f(MPI_Group group){
2963 return smpi_group_c2f(group);
2966 MPI_Request PMPI_Request_f2c(MPI_Fint request){
2967 return smpi_request_f2c(request);
2970 MPI_Fint PMPI_Request_c2f(MPI_Request request) {
2971 return smpi_request_c2f(request);
2974 MPI_Win PMPI_Win_f2c(MPI_Fint win){
2975 return smpi_win_f2c(win);
2978 MPI_Fint PMPI_Win_c2f(MPI_Win win){
2979 return smpi_win_c2f(win);
2982 MPI_Op PMPI_Op_f2c(MPI_Fint op){
2983 return smpi_op_f2c(op);
2986 MPI_Fint PMPI_Op_c2f(MPI_Op op){
2987 return smpi_op_c2f(op);
2990 MPI_Comm PMPI_Comm_f2c(MPI_Fint comm){
2991 return smpi_comm_f2c(comm);
2994 MPI_Fint PMPI_Comm_c2f(MPI_Comm comm){
2995 return smpi_comm_c2f(comm);
2998 MPI_Info PMPI_Info_f2c(MPI_Fint info){
2999 return smpi_info_f2c(info);
3002 MPI_Fint PMPI_Info_c2f(MPI_Info info){
3003 return smpi_info_c2f(info);
3006 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
3007 return smpi_comm_keyval_create(copy_fn, delete_fn, keyval, extra_state);
3010 int PMPI_Keyval_free(int* keyval) {
3011 return smpi_comm_keyval_free(keyval);
3014 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
3015 if(keyval == MPI_TAG_UB||keyval == MPI_HOST||keyval == MPI_IO ||keyval == MPI_WTIME_IS_GLOBAL||keyval == MPI_APPNUM
3016 ||keyval == MPI_UNIVERSE_SIZE||keyval == MPI_LASTUSEDCODE)
3018 else if (comm==MPI_COMM_NULL)
3019 return MPI_ERR_COMM;
3021 return smpi_comm_attr_delete(comm, keyval);
3024 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
3026 static int zero = 0;
3027 static int tag_ub = 1000000;
3028 static int last_used_code = MPI_ERR_LASTCODE;
3030 if (comm==MPI_COMM_NULL){
3032 return MPI_ERR_COMM;
3040 *(int**)attr_value = &zero;
3042 case MPI_UNIVERSE_SIZE:
3044 *(int**)attr_value = &smpi_universe_size;
3046 case MPI_LASTUSEDCODE:
3048 *(int**)attr_value = &last_used_code;
3052 *(int**)attr_value = &tag_ub;
3054 case MPI_WTIME_IS_GLOBAL:
3056 *(int**)attr_value = &one;
3059 return smpi_comm_attr_get(comm, keyval, attr_value, flag);
3063 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
3064 if(keyval == MPI_TAG_UB||keyval == MPI_HOST||keyval == MPI_IO ||keyval == MPI_WTIME_IS_GLOBAL||keyval == MPI_APPNUM
3065 ||keyval == MPI_UNIVERSE_SIZE||keyval == MPI_LASTUSEDCODE)
3067 else if (comm==MPI_COMM_NULL)
3068 return MPI_ERR_COMM;
3070 return smpi_comm_attr_put(comm, keyval, attr_value);
3073 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
3075 return PMPI_Attr_get(comm, comm_keyval, attribute_val,flag);
3078 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
3080 return PMPI_Attr_put(comm, comm_keyval, attribute_val);
3083 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
3085 return PMPI_Attr_delete(comm, comm_keyval);
3088 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval,
3091 return PMPI_Keyval_create(copy_fn, delete_fn, keyval, extra_state);
3094 int PMPI_Comm_free_keyval(int* keyval) {
3095 return PMPI_Keyval_free(keyval);
3098 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
3100 if (type==MPI_DATATYPE_NULL)
3101 return MPI_ERR_TYPE;
3103 return smpi_type_attr_get(type, type_keyval, attribute_val, flag);
3106 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
3108 if (type==MPI_DATATYPE_NULL)
3109 return MPI_ERR_TYPE;
3111 return smpi_type_attr_put(type, type_keyval, attribute_val);
3114 int PMPI_Type_delete_attr (MPI_Datatype type, int type_keyval)
3116 if (type==MPI_DATATYPE_NULL)
3117 return MPI_ERR_TYPE;
3119 return smpi_type_attr_delete(type, type_keyval);
3122 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval,
3125 return smpi_type_keyval_create(copy_fn, delete_fn, keyval, extra_state);
3128 int PMPI_Type_free_keyval(int* keyval) {
3129 return smpi_type_keyval_free(keyval);
3132 int PMPI_Info_create( MPI_Info *info){
3135 *info = xbt_new(s_smpi_mpi_info_t, 1);
3136 (*info)->info_dict= xbt_dict_new_homogeneous(NULL);
3137 (*info)->refcount=1;
3141 int PMPI_Info_set( MPI_Info info, char *key, char *value){
3142 if (info == NULL || key == NULL || value == NULL)
3145 xbt_dict_set(info->info_dict, key, (void*)value, NULL);
3149 int PMPI_Info_free( MPI_Info *info){
3150 if (info == NULL || *info==NULL)
3152 (*info)->refcount--;
3153 if((*info)->refcount==0){
3154 xbt_dict_free(&((*info)->info_dict));
3157 *info=MPI_INFO_NULL;
3161 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
3163 if (info == NULL || key == NULL || valuelen <0)
3166 return MPI_ERR_INFO_VALUE;
3167 char* tmpvalue=(char*)xbt_dict_get_or_null(info->info_dict, key);
3169 memcpy(value,tmpvalue, (strlen(tmpvalue) + 1 < static_cast<size_t>(valuelen)) ? strlen(tmpvalue) + 1 : valuelen);
3175 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
3176 if (info == NULL || newinfo==NULL)
3178 *newinfo = xbt_new(s_smpi_mpi_info_t, 1);
3179 (*newinfo)->info_dict= xbt_dict_new_homogeneous(NULL);
3180 xbt_dict_cursor_t cursor = NULL;
3183 xbt_dict_foreach(info->info_dict,cursor,key,data){
3184 xbt_dict_set((*newinfo)->info_dict, (char*)key, data, NULL);
3189 int PMPI_Info_delete(MPI_Info info, char *key){
3191 if (info == NULL || key==NULL)
3194 xbt_dict_remove(info->info_dict, key);
3197 return MPI_ERR_INFO_NOKEY;
3202 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
3203 if (info == NULL || nkeys==NULL)
3205 *nkeys=xbt_dict_size(info->info_dict);
3209 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
3210 if (info == NULL || key==NULL || n<0 || n> MPI_MAX_INFO_KEY)
3213 xbt_dict_cursor_t cursor = NULL;
3217 xbt_dict_foreach(info->info_dict,cursor,keyn,data){
3227 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
3229 if (info == NULL || key == NULL || valuelen==NULL)
3231 char* tmpvalue=(char*)xbt_dict_get_or_null(info->info_dict, key);
3233 *valuelen=strlen(tmpvalue);
3239 int PMPI_Unpack(void* inbuf, int incount, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
3240 if(incount<0 || outcount < 0 || inbuf==NULL || outbuf==NULL)
3242 if(!is_datatype_valid(type))
3243 return MPI_ERR_TYPE;
3244 if(comm==MPI_COMM_NULL)
3245 return MPI_ERR_COMM;
3246 return smpi_mpi_unpack(inbuf, incount, position, outbuf,outcount,type, comm);
3249 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
3250 if(incount<0 || outcount < 0|| inbuf==NULL || outbuf==NULL)
3252 if(!is_datatype_valid(type))
3253 return MPI_ERR_TYPE;
3254 if(comm==MPI_COMM_NULL)
3255 return MPI_ERR_COMM;
3256 return smpi_mpi_pack(inbuf, incount, type, outbuf,outcount,position, comm);
3259 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
3262 if(!is_datatype_valid(datatype))
3263 return MPI_ERR_TYPE;
3264 if(comm==MPI_COMM_NULL)
3265 return MPI_ERR_COMM;
3267 *size=incount*smpi_datatype_size(datatype);