1 /* Copyright (c) 2007-2014. The SimGrid Team.
2 * All rights reserved. */
4 /* This program is free software; you can redistribute it and/or modify it
5 * under the terms of the license (GNU LGPL) which comes with this package. */
8 #include "smpi_mpi_dt_private.h"
10 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_pmpi, smpi,
11 "Logging specific to SMPI (pmpi)");
14 //this function need to be here because of the calls to smpi_bench
15 void TRACE_smpi_set_category(const char *category)
17 //need to end bench otherwise categories for execution tasks are wrong
19 TRACE_internal_smpi_set_category (category);
20 //begin bench after changing process's category
25 /* PMPI User level calls */
27 int PMPI_Init(int *argc, char ***argv)
29 smpi_process_init(argc, argv);
30 smpi_process_mark_as_initialized();
32 int rank = smpi_process_index();
33 TRACE_smpi_init(rank);
34 TRACE_smpi_computing_init(rank);
35 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
36 extra->type = TRACING_INIT;
37 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
38 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
44 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());
58 smpi_process_destroy();
62 int PMPI_Finalized(int* flag)
64 *flag=smpi_process_finalized();
68 int PMPI_Get_version (int *version,int *subversion){
69 *version = MPI_VERSION;
70 *subversion= MPI_SUBVERSION;
74 int PMPI_Get_library_version (char *version,int *len){
75 int retval = MPI_SUCCESS;
77 snprintf(version,MPI_MAX_LIBRARY_VERSION_STRING,"SMPI Version %d.%d. Copyright The Simgrid Team 2007-2013",SIMGRID_VERSION_MAJOR,
78 SIMGRID_VERSION_MINOR);
79 *len = strlen(version) > MPI_MAX_LIBRARY_VERSION_STRING ? MPI_MAX_LIBRARY_VERSION_STRING : strlen(version);
84 int PMPI_Init_thread(int *argc, char ***argv, int required, int *provided)
86 if (provided != NULL) {
87 *provided = MPI_THREAD_MULTIPLE;
89 return MPI_Init(argc, argv);
92 int PMPI_Query_thread(int *provided)
97 if (provided == NULL) {
100 *provided = MPI_THREAD_MULTIPLE;
101 retval = MPI_SUCCESS;
107 int PMPI_Is_thread_main(int *flag)
113 retval = MPI_ERR_ARG;
115 *flag = smpi_process_index() == 0;
116 retval = MPI_SUCCESS;
122 int PMPI_Abort(MPI_Comm comm, int errorcode)
125 smpi_process_destroy();
126 // FIXME: should kill all processes in comm instead
127 simcall_process_kill(SIMIX_process_self());
131 double PMPI_Wtime(void)
134 if (smpi_process_initialized() && !smpi_process_finalized() && !smpi_process_get_sampling()) {
136 time = SIMIX_get_clock();
139 time = SIMIX_get_clock();
144 extern double sg_maxmin_precision;
145 double PMPI_Wtick(void)
147 return sg_maxmin_precision;
150 int PMPI_Address(void *location, MPI_Aint * address)
156 retval = MPI_ERR_ARG;
158 *address = (MPI_Aint) location;
159 retval = MPI_SUCCESS;
165 int PMPI_Get_address(void *location, MPI_Aint * address)
167 return PMPI_Address(location, address);
170 int PMPI_Type_free(MPI_Datatype * datatype)
176 retval = MPI_ERR_ARG;
178 smpi_datatype_free(datatype);
179 retval = MPI_SUCCESS;
185 int PMPI_Type_size(MPI_Datatype datatype, int *size)
190 if (datatype == MPI_DATATYPE_NULL) {
191 retval = MPI_ERR_TYPE;
192 } else if (size == NULL) {
193 retval = MPI_ERR_ARG;
195 *size = (int) smpi_datatype_size(datatype);
196 retval = MPI_SUCCESS;
202 int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
207 if (datatype == MPI_DATATYPE_NULL) {
208 retval = MPI_ERR_TYPE;
209 } else if (lb == NULL || extent == NULL) {
210 retval = MPI_ERR_ARG;
212 retval = smpi_datatype_extent(datatype, lb, extent);
218 int PMPI_Type_get_true_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
220 return PMPI_Type_get_extent(datatype, lb, extent);
223 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
228 if (datatype == MPI_DATATYPE_NULL) {
229 retval = MPI_ERR_TYPE;
230 } else if (extent == NULL) {
231 retval = MPI_ERR_ARG;
233 *extent = smpi_datatype_get_extent(datatype);
234 retval = MPI_SUCCESS;
240 int PMPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
245 if (datatype == MPI_DATATYPE_NULL) {
246 retval = MPI_ERR_TYPE;
247 } else if (disp == NULL) {
248 retval = MPI_ERR_ARG;
250 *disp = smpi_datatype_lb(datatype);
251 retval = MPI_SUCCESS;
257 int PMPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
262 if (datatype == MPI_DATATYPE_NULL) {
263 retval = MPI_ERR_TYPE;
264 } else if (disp == NULL) {
265 retval = MPI_ERR_ARG;
267 *disp = smpi_datatype_ub(datatype);
268 retval = MPI_SUCCESS;
274 int PMPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
279 if (function == NULL || op == NULL) {
280 retval = MPI_ERR_ARG;
282 *op = smpi_op_new(function, commute);
283 retval = MPI_SUCCESS;
289 int PMPI_Op_free(MPI_Op * op)
295 retval = MPI_ERR_ARG;
296 } else if (*op == MPI_OP_NULL) {
299 smpi_op_destroy(*op);
301 retval = MPI_SUCCESS;
307 int PMPI_Group_free(MPI_Group * group)
313 retval = MPI_ERR_ARG;
315 smpi_group_destroy(*group);
316 *group = MPI_GROUP_NULL;
317 retval = MPI_SUCCESS;
323 int PMPI_Group_size(MPI_Group group, int *size)
328 if (group == MPI_GROUP_NULL) {
329 retval = MPI_ERR_GROUP;
330 } else if (size == NULL) {
331 retval = MPI_ERR_ARG;
333 *size = smpi_group_size(group);
334 retval = MPI_SUCCESS;
340 int PMPI_Group_rank(MPI_Group group, int *rank)
345 if (group == MPI_GROUP_NULL) {
346 retval = MPI_ERR_GROUP;
347 } else if (rank == NULL) {
348 retval = MPI_ERR_ARG;
350 *rank = smpi_group_rank(group, smpi_process_index());
351 retval = MPI_SUCCESS;
357 int PMPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1,
358 MPI_Group group2, int *ranks2)
360 int retval, i, index;
362 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
363 retval = MPI_ERR_GROUP;
365 for (i = 0; i < n; i++) {
366 if(ranks1[i]==MPI_PROC_NULL){
367 ranks2[i]=MPI_PROC_NULL;
369 index = smpi_group_index(group1, ranks1[i]);
370 ranks2[i] = smpi_group_rank(group2, index);
373 retval = MPI_SUCCESS;
379 int PMPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
384 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
385 retval = MPI_ERR_GROUP;
386 } else if (result == NULL) {
387 retval = MPI_ERR_ARG;
389 *result = smpi_group_compare(group1, group2);
390 retval = MPI_SUCCESS;
396 int PMPI_Group_union(MPI_Group group1, MPI_Group group2,
397 MPI_Group * newgroup)
399 int retval, i, proc1, proc2, size, size2;
402 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
403 retval = MPI_ERR_GROUP;
404 } else if (newgroup == NULL) {
405 retval = MPI_ERR_ARG;
407 size = smpi_group_size(group1);
408 size2 = smpi_group_size(group2);
409 for (i = 0; i < size2; i++) {
410 proc2 = smpi_group_index(group2, i);
411 proc1 = smpi_group_rank(group1, proc2);
412 if (proc1 == MPI_UNDEFINED) {
417 *newgroup = MPI_GROUP_EMPTY;
419 *newgroup = smpi_group_new(size);
420 size2 = smpi_group_size(group1);
421 for (i = 0; i < size2; i++) {
422 proc1 = smpi_group_index(group1, i);
423 smpi_group_set_mapping(*newgroup, proc1, i);
425 for (i = size2; i < size; i++) {
426 proc2 = smpi_group_index(group2, i - size2);
427 smpi_group_set_mapping(*newgroup, proc2, i);
430 retval = MPI_SUCCESS;
436 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
437 MPI_Group * newgroup)
439 int retval, i, proc1, proc2, size;
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(group2);
448 for (i = 0; i < size; i++) {
449 proc2 = smpi_group_index(group2, i);
450 proc1 = smpi_group_rank(group1, proc2);
451 if (proc1 == MPI_UNDEFINED) {
456 *newgroup = MPI_GROUP_EMPTY;
458 *newgroup = smpi_group_new(size);
460 for (i = 0; i < smpi_group_size(group2); i++) {
461 proc2 = smpi_group_index(group2, i);
462 proc1 = smpi_group_rank(group1, proc2);
463 if (proc1 != MPI_UNDEFINED) {
464 smpi_group_set_mapping(*newgroup, proc2, j);
469 retval = MPI_SUCCESS;
475 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
477 int retval, i, proc1, proc2, size, size2;
480 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
481 retval = MPI_ERR_GROUP;
482 } else if (newgroup == NULL) {
483 retval = MPI_ERR_ARG;
485 size = size2 = smpi_group_size(group1);
486 for (i = 0; i < size2; i++) {
487 proc1 = smpi_group_index(group1, i);
488 proc2 = smpi_group_rank(group2, proc1);
489 if (proc2 != MPI_UNDEFINED) {
494 *newgroup = MPI_GROUP_EMPTY;
496 *newgroup = smpi_group_new(size);
497 for (i = 0; i < size2; i++) {
498 proc1 = smpi_group_index(group1, i);
499 proc2 = smpi_group_rank(group2, proc1);
500 if (proc2 == MPI_UNDEFINED) {
501 smpi_group_set_mapping(*newgroup, proc1, i);
505 retval = MPI_SUCCESS;
511 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
513 int retval, i, index;
516 if (group == MPI_GROUP_NULL) {
517 retval = MPI_ERR_GROUP;
518 } else if (newgroup == NULL) {
519 retval = MPI_ERR_ARG;
522 *newgroup = MPI_GROUP_EMPTY;
523 } else if (n == smpi_group_size(group)) {
525 if(group!= smpi_comm_group(MPI_COMM_WORLD)
526 && group != MPI_GROUP_NULL
527 && group != smpi_comm_group(MPI_COMM_SELF)
528 && group != MPI_GROUP_EMPTY)
529 smpi_group_use(group);
531 *newgroup = smpi_group_new(n);
532 for (i = 0; i < n; i++) {
533 index = smpi_group_index(group, ranks[i]);
534 smpi_group_set_mapping(*newgroup, index, i);
537 retval = MPI_SUCCESS;
543 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
545 int retval, i, j, newsize, oldsize, index;
548 if (group == MPI_GROUP_NULL) {
549 retval = MPI_ERR_GROUP;
550 } else if (newgroup == NULL) {
551 retval = MPI_ERR_ARG;
555 if(group!= smpi_comm_group(MPI_COMM_WORLD)
556 && group != MPI_GROUP_NULL
557 && group != smpi_comm_group(MPI_COMM_SELF)
558 && group != MPI_GROUP_EMPTY)
559 smpi_group_use(group);
560 } else if (n == smpi_group_size(group)) {
561 *newgroup = MPI_GROUP_EMPTY;
563 oldsize=smpi_group_size(group);
564 newsize = oldsize - n;
565 *newgroup = smpi_group_new(newsize);
567 int* to_exclude=xbt_new0(int, smpi_group_size(group));
568 for(i=0; i<oldsize; i++)
571 to_exclude[ranks[i]]=1;
574 for(i=0; i<oldsize; i++){
575 if(to_exclude[i]==0){
576 index = smpi_group_index(group, i);
577 smpi_group_set_mapping(*newgroup, index, j);
582 xbt_free(to_exclude);
584 retval = MPI_SUCCESS;
590 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
591 MPI_Group * newgroup)
593 int retval, i, j, rank, size, index;
596 if (group == MPI_GROUP_NULL) {
597 retval = MPI_ERR_GROUP;
598 } else if (newgroup == NULL) {
599 retval = MPI_ERR_ARG;
602 *newgroup = MPI_GROUP_EMPTY;
605 for (i = 0; i < n; i++) {
606 for (rank = ranges[i][0]; /* First */
607 rank >= 0; /* Last */
611 rank += ranges[i][2]; /* Stride */
612 if (ranges[i][0]<ranges[i][1]){
613 if(rank > ranges[i][1])
616 if(rank < ranges[i][1])
622 *newgroup = smpi_group_new(size);
624 for (i = 0; i < n; i++) {
625 for (rank = ranges[i][0]; /* First */
626 rank >= 0; /* Last */
628 index = smpi_group_index(group, rank);
629 smpi_group_set_mapping(*newgroup, index, j);
631 rank += ranges[i][2]; /* Stride */
632 if (ranges[i][0]<ranges[i][1]){
633 if(rank > ranges[i][1])
636 if(rank < ranges[i][1])
642 retval = MPI_SUCCESS;
648 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
649 MPI_Group * newgroup)
651 int retval, i, rank, newrank,oldrank, size, index, add;
654 if (group == MPI_GROUP_NULL) {
655 retval = MPI_ERR_GROUP;
656 } else if (newgroup == NULL) {
657 retval = MPI_ERR_ARG;
661 if(group!= smpi_comm_group(MPI_COMM_WORLD)
662 && group != MPI_GROUP_NULL
663 && group != smpi_comm_group(MPI_COMM_SELF)
664 && group != MPI_GROUP_EMPTY)
665 smpi_group_use(group);
667 size = smpi_group_size(group);
668 for (i = 0; i < n; i++) {
669 for (rank = ranges[i][0]; /* First */
670 rank >= 0; /* Last */
674 rank += ranges[i][2]; /* Stride */
675 if (ranges[i][0]<ranges[i][1]){
676 if(rank > ranges[i][1])
679 if(rank < ranges[i][1])
685 *newgroup = MPI_GROUP_EMPTY;
687 *newgroup = smpi_group_new(size);
690 while (newrank < size) {
692 for (i = 0; i < n; i++) {
693 for (rank = ranges[i][0];rank >= 0;){
699 rank += ranges[i][2]; /* Stride */
701 if (ranges[i][0]<ranges[i][1]){
702 if(rank > ranges[i][1])
705 if(rank < ranges[i][1])
711 index = smpi_group_index(group, oldrank);
712 smpi_group_set_mapping(*newgroup, index, newrank);
720 retval = MPI_SUCCESS;
726 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
731 if (comm == MPI_COMM_NULL) {
732 retval = MPI_ERR_COMM;
733 } else if (rank == NULL) {
734 retval = MPI_ERR_ARG;
736 *rank = smpi_comm_rank(comm);
737 retval = MPI_SUCCESS;
743 int PMPI_Comm_size(MPI_Comm comm, int *size)
748 if (comm == MPI_COMM_NULL) {
749 retval = MPI_ERR_COMM;
750 } else if (size == NULL) {
751 retval = MPI_ERR_ARG;
753 *size = smpi_comm_size(comm);
754 retval = MPI_SUCCESS;
760 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
765 if (comm == MPI_COMM_NULL) {
766 retval = MPI_ERR_COMM;
767 } else if (name == NULL || len == NULL) {
768 retval = MPI_ERR_ARG;
770 smpi_comm_get_name(comm, name, len);
771 retval = MPI_SUCCESS;
777 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
782 if (comm == MPI_COMM_NULL) {
783 retval = MPI_ERR_COMM;
784 } else if (group == NULL) {
785 retval = MPI_ERR_ARG;
787 *group = smpi_comm_group(comm);
788 if(*group!= smpi_comm_group(MPI_COMM_WORLD)
789 && *group != MPI_GROUP_NULL
790 && *group != smpi_comm_group(MPI_COMM_SELF)
791 && *group != MPI_GROUP_EMPTY)
792 smpi_group_use(*group);
793 retval = MPI_SUCCESS;
799 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
804 if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
805 retval = MPI_ERR_COMM;
806 } else if (result == NULL) {
807 retval = MPI_ERR_ARG;
809 if (comm1 == comm2) { /* Same communicators means same groups */
813 smpi_group_compare(smpi_comm_group(comm1),
814 smpi_comm_group(comm2));
815 if (*result == MPI_IDENT) {
816 *result = MPI_CONGRUENT;
819 retval = MPI_SUCCESS;
825 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
830 if (comm == MPI_COMM_NULL) {
831 retval = MPI_ERR_COMM;
832 } else if (newcomm == NULL) {
833 retval = MPI_ERR_ARG;
835 *newcomm = smpi_comm_new(smpi_comm_group(comm));
836 retval = MPI_SUCCESS;
842 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
847 if (comm == MPI_COMM_NULL) {
848 retval = MPI_ERR_COMM;
849 } else if (group == MPI_GROUP_NULL) {
850 retval = MPI_ERR_GROUP;
851 } else if (newcomm == NULL) {
852 retval = MPI_ERR_ARG;
853 } else if(smpi_group_rank(group,smpi_process_index())==MPI_UNDEFINED){
854 *newcomm= MPI_COMM_NULL;
855 retval = MPI_SUCCESS;
858 *newcomm = smpi_comm_new(group);
859 retval = MPI_SUCCESS;
865 int PMPI_Comm_free(MPI_Comm * comm)
871 retval = MPI_ERR_ARG;
872 } else if (*comm == MPI_COMM_NULL) {
873 retval = MPI_ERR_COMM;
875 smpi_comm_destroy(*comm);
876 *comm = MPI_COMM_NULL;
877 retval = MPI_SUCCESS;
883 int PMPI_Comm_disconnect(MPI_Comm * comm)
885 /* TODO: wait until all communication in comm are done */
890 retval = MPI_ERR_ARG;
891 } else if (*comm == MPI_COMM_NULL) {
892 retval = MPI_ERR_COMM;
894 smpi_comm_destroy(*comm);
895 *comm = MPI_COMM_NULL;
896 retval = MPI_SUCCESS;
902 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
907 if (comm_out == NULL) {
908 retval = MPI_ERR_ARG;
909 } else if (comm == MPI_COMM_NULL) {
910 retval = MPI_ERR_COMM;
912 *comm_out = smpi_comm_split(comm, color, key);
913 retval = MPI_SUCCESS;
919 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
920 int tag, MPI_Comm comm, MPI_Request * request)
925 if (request == NULL) {
926 retval = MPI_ERR_ARG;
927 } else if (comm == MPI_COMM_NULL) {
928 retval = MPI_ERR_COMM;
929 } else if (dst == MPI_PROC_NULL) {
930 retval = MPI_SUCCESS;
932 *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
933 retval = MPI_SUCCESS;
936 if (retval != MPI_SUCCESS && request)
937 *request = MPI_REQUEST_NULL;
941 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
942 int tag, MPI_Comm comm, MPI_Request * request)
947 if (request == NULL) {
948 retval = MPI_ERR_ARG;
949 } else if (comm == MPI_COMM_NULL) {
950 retval = MPI_ERR_COMM;
951 } else if (src == MPI_PROC_NULL) {
952 retval = MPI_SUCCESS;
954 *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
955 retval = MPI_SUCCESS;
958 if (retval != MPI_SUCCESS && request)
959 *request = MPI_REQUEST_NULL;
963 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype,
964 int dst, int tag, MPI_Comm comm, MPI_Request* request)
969 if (request == NULL) {
970 retval = MPI_ERR_ARG;
971 } else if (comm == MPI_COMM_NULL) {
972 retval = MPI_ERR_COMM;
973 } else if (dst == MPI_PROC_NULL) {
974 retval = MPI_SUCCESS;
976 *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
977 retval = MPI_SUCCESS;
980 if (retval != MPI_SUCCESS && request)
981 *request = MPI_REQUEST_NULL;
985 int PMPI_Start(MPI_Request * request)
990 if (request == NULL || *request == MPI_REQUEST_NULL) {
991 retval = MPI_ERR_ARG;
993 smpi_mpi_start(*request);
994 retval = MPI_SUCCESS;
1000 int PMPI_Startall(int count, MPI_Request * requests)
1005 if (requests == NULL) {
1006 retval = MPI_ERR_ARG;
1008 smpi_mpi_startall(count, requests);
1009 retval = MPI_SUCCESS;
1015 int PMPI_Request_free(MPI_Request * request)
1020 if (*request == MPI_REQUEST_NULL) {
1021 retval = MPI_ERR_ARG;
1023 smpi_mpi_request_free(request);
1024 retval = MPI_SUCCESS;
1030 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
1031 int tag, MPI_Comm comm, MPI_Request * request)
1037 if (request == NULL) {
1038 retval = MPI_ERR_ARG;
1039 } else if (comm == MPI_COMM_NULL) {
1040 retval = MPI_ERR_COMM;
1041 } else if (src == MPI_PROC_NULL) {
1042 *request = MPI_REQUEST_NULL;
1043 retval = MPI_SUCCESS;
1044 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1045 retval = MPI_ERR_COMM;
1046 } else if (count < 0) {
1047 retval = MPI_ERR_COUNT;
1048 } else if (buf==NULL && count > 0) {
1049 retval = MPI_ERR_COUNT;
1050 } else if (datatype == MPI_DATATYPE_NULL){
1051 retval = MPI_ERR_TYPE;
1052 } else if(tag<0 && tag != MPI_ANY_TAG){
1053 retval = MPI_ERR_TAG;
1057 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1058 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1060 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1061 extra->type = TRACING_IRECV;
1062 extra->send_size = count;
1063 extra->src = src_traced;
1065 extra->datatype1 = encode_datatype(datatype);
1066 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
1069 *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
1070 retval = MPI_SUCCESS;
1073 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1074 (*request)->recv = 1;
1079 if (retval != MPI_SUCCESS && request)
1080 *request = MPI_REQUEST_NULL;
1085 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
1086 int tag, MPI_Comm comm, MPI_Request * request)
1091 if (request == NULL) {
1092 retval = MPI_ERR_ARG;
1093 } else if (comm == MPI_COMM_NULL) {
1094 retval = MPI_ERR_COMM;
1095 } else if (dst == MPI_PROC_NULL) {
1096 *request = MPI_REQUEST_NULL;
1097 retval = MPI_SUCCESS;
1098 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1099 retval = MPI_ERR_RANK;
1100 } else if (count < 0) {
1101 retval = MPI_ERR_COUNT;
1102 } else if (buf==NULL && count > 0) {
1103 retval = MPI_ERR_COUNT;
1104 } else if (datatype == MPI_DATATYPE_NULL){
1105 retval = MPI_ERR_TYPE;
1106 } else if(tag<0 && tag != MPI_ANY_TAG){
1107 retval = MPI_ERR_TAG;
1111 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1112 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1114 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1115 extra->type = TRACING_ISEND;
1116 extra->send_size = count;
1118 extra->dst = dst_traced;
1119 extra->datatype1 = encode_datatype(datatype);
1120 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1121 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1124 *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1125 retval = MPI_SUCCESS;
1128 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1129 (*request)->send = 1;
1134 if (retval != MPI_SUCCESS && request)
1135 *request = MPI_REQUEST_NULL;
1139 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype,
1140 int dst, int tag, MPI_Comm comm, MPI_Request* request)
1145 if (request == NULL) {
1146 retval = MPI_ERR_ARG;
1147 } else if (comm == MPI_COMM_NULL) {
1148 retval = MPI_ERR_COMM;
1149 } else if (dst == MPI_PROC_NULL) {
1150 *request = MPI_REQUEST_NULL;
1151 retval = MPI_SUCCESS;
1152 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1153 retval = MPI_ERR_RANK;
1154 } else if (count < 0) {
1155 retval = MPI_ERR_COUNT;
1156 } else if (buf==NULL && count > 0) {
1157 retval = MPI_ERR_COUNT;
1158 } else if (datatype == MPI_DATATYPE_NULL){
1159 retval = MPI_ERR_TYPE;
1160 } else if(tag<0 && tag != MPI_ANY_TAG){
1161 retval = MPI_ERR_TAG;
1165 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1166 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1167 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1168 extra->type = TRACING_ISSEND;
1169 extra->send_size = count;
1171 extra->dst = dst_traced;
1172 extra->datatype1 = encode_datatype(datatype);
1173 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1174 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1177 *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1178 retval = MPI_SUCCESS;
1181 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1182 (*request)->send = 1;
1187 if (retval != MPI_SUCCESS && request)
1188 *request = MPI_REQUEST_NULL;
1192 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1193 MPI_Comm comm, MPI_Status * status)
1198 if (comm == MPI_COMM_NULL) {
1199 retval = MPI_ERR_COMM;
1200 } else if (src == MPI_PROC_NULL) {
1201 smpi_empty_status(status);
1202 status->MPI_SOURCE = MPI_PROC_NULL;
1203 retval = MPI_SUCCESS;
1204 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1205 retval = MPI_ERR_RANK;
1206 } else if (count < 0) {
1207 retval = MPI_ERR_COUNT;
1208 } else if (buf==NULL && count > 0) {
1209 retval = MPI_ERR_COUNT;
1210 } else if (datatype == MPI_DATATYPE_NULL){
1211 retval = MPI_ERR_TYPE;
1212 } else if(tag<0 && tag != MPI_ANY_TAG){
1213 retval = MPI_ERR_TAG;
1216 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1217 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1218 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1219 extra->type = TRACING_RECV;
1220 extra->send_size = count;
1221 extra->src = src_traced;
1223 extra->datatype1 = encode_datatype(datatype);
1224 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
1227 smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1228 retval = MPI_SUCCESS;
1231 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1232 if(status!=MPI_STATUS_IGNORE){
1233 src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1234 TRACE_smpi_recv(rank, src_traced, rank);
1236 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1244 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1251 if (comm == MPI_COMM_NULL) {
1252 retval = MPI_ERR_COMM;
1253 } else if (dst == MPI_PROC_NULL) {
1254 retval = MPI_SUCCESS;
1255 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1256 retval = MPI_ERR_RANK;
1257 } else if (count < 0) {
1258 retval = MPI_ERR_COUNT;
1259 } else if (buf==NULL && count > 0) {
1260 retval = MPI_ERR_COUNT;
1261 } else if (datatype == MPI_DATATYPE_NULL){
1262 retval = MPI_ERR_TYPE;
1263 } else if(tag<0 && tag != MPI_ANY_TAG){
1264 retval = MPI_ERR_TAG;
1268 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1269 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1270 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1271 extra->type = TRACING_SEND;
1272 extra->send_size = count;
1274 extra->dst = dst_traced;
1275 extra->datatype1 = encode_datatype(datatype);
1276 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1277 TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1280 smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1281 retval = MPI_SUCCESS;
1284 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1294 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1299 if (comm == MPI_COMM_NULL) {
1300 retval = MPI_ERR_COMM;
1301 } else if (dst == MPI_PROC_NULL) {
1302 retval = MPI_SUCCESS;
1303 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1304 retval = MPI_ERR_RANK;
1305 } else if (count < 0) {
1306 retval = MPI_ERR_COUNT;
1307 } else if (buf==NULL && count > 0) {
1308 retval = MPI_ERR_COUNT;
1309 } else if (datatype == MPI_DATATYPE_NULL){
1310 retval = MPI_ERR_TYPE;
1311 } else if(tag<0 && tag != MPI_ANY_TAG){
1312 retval = MPI_ERR_TAG;
1316 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1317 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1318 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1319 extra->type = TRACING_SSEND;
1320 extra->send_size = count;
1322 extra->dst = dst_traced;
1323 extra->datatype1 = encode_datatype(datatype);
1324 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra); TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1327 smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1328 retval = MPI_SUCCESS;
1331 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1339 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1340 int dst, int sendtag, void *recvbuf, int recvcount,
1341 MPI_Datatype recvtype, int src, int recvtag,
1342 MPI_Comm comm, MPI_Status * status)
1348 if (comm == MPI_COMM_NULL) {
1349 retval = MPI_ERR_COMM;
1350 } else if (sendtype == MPI_DATATYPE_NULL
1351 || recvtype == MPI_DATATYPE_NULL) {
1352 retval = MPI_ERR_TYPE;
1353 } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1354 smpi_empty_status(status);
1355 status->MPI_SOURCE = MPI_PROC_NULL;
1356 retval = MPI_SUCCESS;
1357 }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1358 (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1359 retval = MPI_ERR_RANK;
1360 } else if (sendcount < 0 || recvcount<0) {
1361 retval = MPI_ERR_COUNT;
1362 } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1363 retval = MPI_ERR_COUNT;
1364 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1365 retval = MPI_ERR_TAG;
1369 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1370 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1371 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1372 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1373 extra->type = TRACING_SENDRECV;
1374 extra->send_size = sendcount;
1375 extra->recv_size = recvcount;
1376 extra->src = src_traced;
1377 extra->dst = dst_traced;
1378 extra->datatype1 = encode_datatype(sendtype);
1379 extra->datatype2 = encode_datatype(recvtype);
1381 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1382 TRACE_smpi_send(rank, rank, dst_traced,sendcount*smpi_datatype_size(sendtype));
1386 smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1387 recvcount, recvtype, src, recvtag, comm, status);
1388 retval = MPI_SUCCESS;
1391 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1392 TRACE_smpi_recv(rank, src_traced, rank);
1401 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1402 int dst, int sendtag, int src, int recvtag,
1403 MPI_Comm comm, MPI_Status * status)
1405 //TODO: suboptimal implementation
1408 if (datatype == MPI_DATATYPE_NULL) {
1409 retval = MPI_ERR_TYPE;
1410 } else if (count < 0) {
1411 retval = MPI_ERR_COUNT;
1413 int size = smpi_datatype_get_extent(datatype) * count;
1414 recvbuf = xbt_new0(char, size);
1416 MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1417 datatype, src, recvtag, comm, status);
1418 if(retval==MPI_SUCCESS){
1419 smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1427 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1432 if (request == NULL || flag == NULL) {
1433 retval = MPI_ERR_ARG;
1434 } else if (*request == MPI_REQUEST_NULL) {
1436 smpi_empty_status(status);
1437 retval = MPI_ERR_REQUEST;
1439 *flag = smpi_mpi_test(request, status);
1440 retval = MPI_SUCCESS;
1446 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1447 MPI_Status * status)
1452 if (index == NULL || flag == NULL) {
1453 retval = MPI_ERR_ARG;
1455 *flag = smpi_mpi_testany(count, requests, index, status);
1456 retval = MPI_SUCCESS;
1462 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1468 retval = MPI_ERR_ARG;
1470 *flag = smpi_mpi_testall(count, requests, statuses);
1471 retval = MPI_SUCCESS;
1477 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1481 if (status == NULL) {
1482 retval = MPI_ERR_ARG;
1483 } else if (comm == MPI_COMM_NULL) {
1484 retval = MPI_ERR_COMM;
1485 } else if (source == MPI_PROC_NULL) {
1486 smpi_empty_status(status);
1487 status->MPI_SOURCE = MPI_PROC_NULL;
1488 retval = MPI_SUCCESS;
1490 smpi_mpi_probe(source, tag, comm, status);
1491 retval = MPI_SUCCESS;
1498 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1503 retval = MPI_ERR_ARG;
1504 } else if (status == NULL) {
1505 retval = MPI_ERR_ARG;
1506 } else if (comm == MPI_COMM_NULL) {
1507 retval = MPI_ERR_COMM;
1508 } else if (source == MPI_PROC_NULL) {
1510 smpi_empty_status(status);
1511 status->MPI_SOURCE = MPI_PROC_NULL;
1512 retval = MPI_SUCCESS;
1514 smpi_mpi_iprobe(source, tag, comm, flag, status);
1515 retval = MPI_SUCCESS;
1521 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1527 smpi_empty_status(status);
1529 if (request == NULL) {
1530 retval = MPI_ERR_ARG;
1531 } else if (*request == MPI_REQUEST_NULL) {
1532 retval = MPI_ERR_REQUEST;
1536 int rank = request && (*request)->comm != MPI_COMM_NULL
1537 ? smpi_process_index()
1540 int src_traced = (*request)->src;
1541 int dst_traced = (*request)->dst;
1542 MPI_Comm comm = (*request)->comm;
1543 int is_wait_for_receive = (*request)->recv;
1544 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1545 extra->type = TRACING_WAIT;
1546 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1549 smpi_mpi_wait(request, status);
1550 retval = MPI_SUCCESS;
1553 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1554 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1555 if (is_wait_for_receive) {
1556 if(src_traced==MPI_ANY_SOURCE)
1557 src_traced = (status!=MPI_STATUS_IGNORE) ?
1558 smpi_group_rank(smpi_comm_group(comm), status->MPI_SOURCE) :
1560 TRACE_smpi_recv(rank, src_traced, dst_traced);
1570 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1577 //save requests information for tracing
1579 int *srcs = xbt_new0(int, count);
1580 int *dsts = xbt_new0(int, count);
1581 int *recvs = xbt_new0(int, count);
1582 MPI_Comm *comms = xbt_new0(MPI_Comm, count);
1584 for (i = 0; i < count; i++) {
1585 MPI_Request req = requests[i]; //already received requests are no longer valid
1589 recvs[i] = req->recv;
1590 comms[i] = req->comm;
1593 int rank_traced = smpi_process_index();
1594 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1595 extra->type = TRACING_WAITANY;
1596 extra->send_size=count;
1597 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1600 *index = smpi_mpi_waitany(count, requests, status);
1602 if(*index!=MPI_UNDEFINED){
1603 int src_traced = srcs[*index];
1604 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1605 int dst_traced = dsts[*index];
1606 int is_wait_for_receive = recvs[*index];
1607 if (is_wait_for_receive) {
1608 if(srcs[*index]==MPI_ANY_SOURCE)
1609 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1610 smpi_group_rank(smpi_comm_group(comms[*index]), status->MPI_SOURCE) :
1612 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1614 TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1626 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1631 //save information from requests
1633 int *srcs = xbt_new0(int, count);
1634 int *dsts = xbt_new0(int, count);
1635 int *recvs = xbt_new0(int, count);
1636 int *valid = xbt_new0(int, count);
1637 MPI_Comm *comms = xbt_new0(MPI_Comm, count);
1639 //int valid_count = 0;
1640 for (i = 0; i < count; i++) {
1641 MPI_Request req = requests[i];
1642 if(req!=MPI_REQUEST_NULL){
1645 recvs[i] = req->recv;
1646 comms[i] = req->comm;
1652 int rank_traced = smpi_process_index();
1653 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1654 extra->type = TRACING_WAITALL;
1655 extra->send_size=count;
1656 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1658 int retval = smpi_mpi_waitall(count, requests, status);
1660 for (i = 0; i < count; i++) {
1662 //int src_traced = srcs[*index];
1663 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1664 int src_traced = srcs[i];
1665 int dst_traced = dsts[i];
1666 int is_wait_for_receive = recvs[i];
1667 if (is_wait_for_receive) {
1668 if(src_traced==MPI_ANY_SOURCE)
1669 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1670 smpi_group_rank(smpi_comm_group(comms[i]), status[i].MPI_SOURCE) :
1672 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1676 TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1688 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1689 int *indices, MPI_Status status[])
1694 if (outcount == NULL) {
1695 retval = MPI_ERR_ARG;
1697 *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1698 retval = MPI_SUCCESS;
1704 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1705 int* indices, MPI_Status status[])
1710 if (outcount == NULL) {
1711 retval = MPI_ERR_ARG;
1713 *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1714 retval = MPI_SUCCESS;
1721 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1727 if (comm == MPI_COMM_NULL) {
1728 retval = MPI_ERR_COMM;
1731 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1732 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1734 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1735 extra->type = TRACING_BCAST;
1736 extra->send_size = count;
1737 extra->root = root_traced;
1738 extra->datatype1 = encode_datatype(datatype);
1739 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1742 mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1743 retval = MPI_SUCCESS;
1745 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1753 int PMPI_Barrier(MPI_Comm comm)
1759 if (comm == MPI_COMM_NULL) {
1760 retval = MPI_ERR_COMM;
1763 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1764 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1765 extra->type = TRACING_BARRIER;
1766 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1768 mpi_coll_barrier_fun(comm);
1769 retval = MPI_SUCCESS;
1771 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1779 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1780 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1781 int root, MPI_Comm comm)
1787 if (comm == MPI_COMM_NULL) {
1788 retval = MPI_ERR_COMM;
1789 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1790 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1791 retval = MPI_ERR_TYPE;
1792 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1793 ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1794 retval = MPI_ERR_COUNT;
1797 char* sendtmpbuf = (char*) sendbuf;
1798 int sendtmpcount = sendcount;
1799 MPI_Datatype sendtmptype = sendtype;
1800 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1802 sendtmptype=recvtype;
1805 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1806 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1807 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1808 extra->type = TRACING_GATHER;
1809 extra->send_size = sendtmpcount;
1810 extra->recv_size = recvcount;
1811 extra->root = root_traced;
1812 extra->datatype1 = encode_datatype(sendtmptype);
1813 extra->datatype2 = encode_datatype(recvtype);
1815 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1817 mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount,
1818 recvtype, root, comm);
1821 retval = MPI_SUCCESS;
1823 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1831 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1832 void *recvbuf, int *recvcounts, int *displs,
1833 MPI_Datatype recvtype, int root, 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 ((smpi_comm_rank(comm) == root) && (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;
1849 char* sendtmpbuf = (char*) sendbuf;
1850 int sendtmpcount = sendcount;
1851 MPI_Datatype sendtmptype = sendtype;
1852 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1854 sendtmptype=recvtype;
1858 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1859 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1861 int size = smpi_comm_size(comm);
1862 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1863 extra->type = TRACING_GATHERV;
1864 extra->send_size = sendtmpcount;
1865 extra->recvcounts= xbt_malloc(size*sizeof(int));
1866 for(i=0; i< size; i++)//copy data to avoid bad free
1867 extra->recvcounts[i] = recvcounts[i];
1868 extra->num_processes = size;
1869 extra->root = root_traced;
1870 extra->datatype1 = encode_datatype(sendtmptype);
1871 extra->datatype2 = encode_datatype(recvtype);
1873 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1875 smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts,
1876 displs, recvtype, root, comm);
1877 retval = MPI_SUCCESS;
1879 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1887 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1888 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1895 if (comm == MPI_COMM_NULL) {
1896 retval = MPI_ERR_COMM;
1897 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1898 (recvtype == MPI_DATATYPE_NULL)){
1899 retval = MPI_ERR_TYPE;
1900 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1902 retval = MPI_ERR_COUNT;
1904 if(sendbuf == MPI_IN_PLACE) {
1905 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1906 sendcount=recvcount;
1910 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1911 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1912 extra->type = TRACING_ALLGATHER;
1913 extra->send_size = sendcount;
1914 extra->recv_size = recvcount;
1915 extra->datatype1 = encode_datatype(sendtype);
1916 extra->datatype2 = encode_datatype(recvtype);
1918 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1920 mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1922 retval = MPI_SUCCESS;
1925 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1932 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1933 void *recvbuf, int *recvcounts, int *displs,
1934 MPI_Datatype recvtype, MPI_Comm comm)
1940 if (comm == MPI_COMM_NULL) {
1941 retval = MPI_ERR_COMM;
1942 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1943 (recvtype == MPI_DATATYPE_NULL)){
1944 retval = MPI_ERR_TYPE;
1945 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1946 retval = MPI_ERR_COUNT;
1947 } else if (recvcounts == NULL || displs == NULL) {
1948 retval = MPI_ERR_ARG;
1951 if(sendbuf == MPI_IN_PLACE) {
1952 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1953 sendcount=recvcounts[smpi_comm_rank(comm)];
1957 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1959 int size = smpi_comm_size(comm);
1960 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1961 extra->type = TRACING_ALLGATHERV;
1962 extra->send_size = sendcount;
1963 extra->recvcounts= xbt_malloc(size*sizeof(int));
1964 for(i=0; i< size; i++)//copy data to avoid bad free
1965 extra->recvcounts[i] = recvcounts[i];
1966 extra->num_processes = size;
1967 extra->datatype1 = encode_datatype(sendtype);
1968 extra->datatype2 = encode_datatype(recvtype);
1970 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
1972 mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1973 displs, recvtype, comm);
1974 retval = MPI_SUCCESS;
1976 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1984 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1985 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1986 int root, MPI_Comm comm)
1992 if (comm == MPI_COMM_NULL) {
1993 retval = MPI_ERR_COMM;
1994 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1995 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1996 retval = MPI_ERR_TYPE;
1999 if (recvbuf == MPI_IN_PLACE) {
2001 recvcount=sendcount;
2004 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2005 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2006 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2007 extra->type = TRACING_SCATTER;
2008 extra->send_size = sendcount;
2009 extra->recv_size= recvcount;
2010 extra->root = root_traced;
2011 extra->datatype1 = encode_datatype(sendtype);
2012 extra->datatype2 = encode_datatype(recvtype);
2014 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2016 mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
2017 recvtype, root, comm);
2018 retval = MPI_SUCCESS;
2020 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2028 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
2029 MPI_Datatype sendtype, void *recvbuf, int recvcount,
2030 MPI_Datatype recvtype, int root, MPI_Comm comm)
2036 if (comm == MPI_COMM_NULL) {
2037 retval = MPI_ERR_COMM;
2038 } else if (sendcounts == NULL || displs == NULL) {
2039 retval = MPI_ERR_ARG;
2040 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
2041 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
2042 retval = MPI_ERR_TYPE;
2044 if (recvbuf == MPI_IN_PLACE) {
2046 recvcount=sendcounts[smpi_comm_rank(comm)];
2049 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2050 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2052 int size = smpi_comm_size(comm);
2053 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2054 extra->type = TRACING_SCATTERV;
2055 extra->recv_size = recvcount;
2056 extra->sendcounts= xbt_malloc(size*sizeof(int));
2057 for(i=0; i< size; i++)//copy data to avoid bad free
2058 extra->sendcounts[i] = sendcounts[i];
2059 extra->num_processes = size;
2060 extra->root = root_traced;
2061 extra->datatype1 = encode_datatype(sendtype);
2062 extra->datatype2 = encode_datatype(recvtype);
2064 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2067 smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
2068 recvcount, recvtype, root, comm);
2069 retval = MPI_SUCCESS;
2071 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2079 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
2080 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
2086 if (comm == MPI_COMM_NULL) {
2087 retval = MPI_ERR_COMM;
2088 } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
2089 retval = MPI_ERR_ARG;
2092 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2093 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2094 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2095 extra->type = TRACING_REDUCE;
2096 extra->send_size = count;
2097 extra->datatype1 = encode_datatype(datatype);
2098 extra->root = root_traced;
2100 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2102 mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
2104 retval = MPI_SUCCESS;
2106 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2114 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
2115 MPI_Datatype datatype, MPI_Op op){
2119 if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
2120 retval = MPI_ERR_ARG;
2122 smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
2129 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
2130 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2136 if (comm == MPI_COMM_NULL) {
2137 retval = MPI_ERR_COMM;
2138 } else if (datatype == MPI_DATATYPE_NULL) {
2139 retval = MPI_ERR_TYPE;
2140 } else if (op == MPI_OP_NULL) {
2141 retval = MPI_ERR_OP;
2144 char* sendtmpbuf = (char*) sendbuf;
2145 if( sendbuf == MPI_IN_PLACE ) {
2146 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
2147 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
2150 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2151 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2152 extra->type = TRACING_ALLREDUCE;
2153 extra->send_size = count;
2154 extra->datatype1 = encode_datatype(datatype);
2156 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2158 mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
2160 if( sendbuf == MPI_IN_PLACE ) {
2161 xbt_free(sendtmpbuf);
2164 retval = MPI_SUCCESS;
2166 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2174 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
2175 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2181 if (comm == MPI_COMM_NULL) {
2182 retval = MPI_ERR_COMM;
2183 } else if (datatype == MPI_DATATYPE_NULL) {
2184 retval = MPI_ERR_TYPE;
2185 } else if (op == MPI_OP_NULL) {
2186 retval = MPI_ERR_OP;
2189 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2190 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2191 extra->type = TRACING_SCAN;
2192 extra->send_size = count;
2193 extra->datatype1 = encode_datatype(datatype);
2195 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2197 smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2198 retval = MPI_SUCCESS;
2200 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2208 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2209 MPI_Op op, MPI_Comm comm){
2214 if (comm == MPI_COMM_NULL) {
2215 retval = MPI_ERR_COMM;
2216 } else if (datatype == MPI_DATATYPE_NULL) {
2217 retval = MPI_ERR_TYPE;
2218 } else if (op == MPI_OP_NULL) {
2219 retval = MPI_ERR_OP;
2222 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2223 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2224 extra->type = TRACING_EXSCAN;
2225 extra->send_size = count;
2226 extra->datatype1 = encode_datatype(datatype);
2228 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2230 smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2231 retval = MPI_SUCCESS;
2233 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2241 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
2242 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2247 if (comm == MPI_COMM_NULL) {
2248 retval = MPI_ERR_COMM;
2249 } else if (datatype == MPI_DATATYPE_NULL) {
2250 retval = MPI_ERR_TYPE;
2251 } else if (op == MPI_OP_NULL) {
2252 retval = MPI_ERR_OP;
2253 } else if (recvcounts == NULL) {
2254 retval = MPI_ERR_ARG;
2257 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2259 int size = smpi_comm_size(comm);
2260 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2261 extra->type = TRACING_REDUCE_SCATTER;
2262 extra->send_size = 0;
2263 extra->recvcounts= xbt_malloc(size*sizeof(int));
2264 for(i=0; i< size; i++)//copy data to avoid bad free
2265 extra->recvcounts[i] = recvcounts[i];
2266 extra->num_processes = size;
2267 extra->datatype1 = encode_datatype(datatype);
2269 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2271 void* sendtmpbuf=sendbuf;
2272 if(sendbuf==MPI_IN_PLACE){
2276 mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts,
2277 datatype, op, comm);
2278 retval = MPI_SUCCESS;
2280 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2288 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2289 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2294 if (comm == MPI_COMM_NULL) {
2295 retval = MPI_ERR_COMM;
2296 } else if (datatype == MPI_DATATYPE_NULL) {
2297 retval = MPI_ERR_TYPE;
2298 } else if (op == MPI_OP_NULL) {
2299 retval = MPI_ERR_OP;
2300 } else if (recvcount < 0) {
2301 retval = MPI_ERR_ARG;
2303 int count=smpi_comm_size(comm);
2306 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2307 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2308 extra->type = TRACING_REDUCE_SCATTER;
2309 extra->send_size = 0;
2310 extra->recvcounts= xbt_malloc(count*sizeof(int));
2311 for(i=0; i< count; i++)//copy data to avoid bad free
2312 extra->recvcounts[i] = recvcount;
2313 extra->num_processes = count;
2314 extra->datatype1 = encode_datatype(datatype);
2316 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2318 int* recvcounts=(int*)xbt_malloc(count);
2319 for (i=0; i<count;i++)recvcounts[i]=recvcount;
2320 mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
2321 datatype, op, comm);
2322 xbt_free(recvcounts);
2323 retval = MPI_SUCCESS;
2325 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2333 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2334 void *recvbuf, int recvcount, MPI_Datatype recvtype,
2341 if (comm == MPI_COMM_NULL) {
2342 retval = MPI_ERR_COMM;
2343 } else if (sendtype == MPI_DATATYPE_NULL
2344 || recvtype == MPI_DATATYPE_NULL) {
2345 retval = MPI_ERR_TYPE;
2348 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2349 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2350 extra->type = TRACING_ALLTOALL;
2351 extra->send_size = sendcount;
2352 extra->recv_size = recvcount;
2353 extra->datatype1 = encode_datatype(sendtype);
2354 extra->datatype2 = encode_datatype(recvtype);
2356 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2358 retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2360 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2368 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2369 MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2370 int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2376 if (comm == MPI_COMM_NULL) {
2377 retval = MPI_ERR_COMM;
2378 } else if (sendtype == MPI_DATATYPE_NULL
2379 || recvtype == MPI_DATATYPE_NULL) {
2380 retval = MPI_ERR_TYPE;
2381 } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2382 || recvdisps == NULL) {
2383 retval = MPI_ERR_ARG;
2386 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2388 int size = smpi_comm_size(comm);
2389 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2390 extra->type = TRACING_ALLTOALLV;
2391 extra->send_size = 0;
2392 extra->recv_size = 0;
2393 extra->recvcounts= xbt_malloc(size*sizeof(int));
2394 extra->sendcounts= xbt_malloc(size*sizeof(int));
2396 for(i=0; i< size; i++){//copy data to avoid bad free
2397 extra->send_size += sendcounts[i];
2398 extra->recv_size += recvcounts[i];
2400 extra->sendcounts[i] = sendcounts[i];
2401 extra->recvcounts[i] = recvcounts[i];
2403 extra->num_processes = size;
2405 extra->datatype1 = encode_datatype(sendtype);
2406 extra->datatype2 = encode_datatype(recvtype);
2408 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2411 mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2412 recvbuf, recvcounts, recvdisps, recvtype,
2415 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2424 int PMPI_Get_processor_name(char *name, int *resultlen)
2426 int retval = MPI_SUCCESS;
2429 strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2430 strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2431 strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2432 MPI_MAX_PROCESSOR_NAME - 1 );
2435 MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2441 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2443 int retval = MPI_SUCCESS;
2447 if (status == NULL || count == NULL) {
2448 retval = MPI_ERR_ARG;
2449 } else if (datatype == MPI_DATATYPE_NULL) {
2450 retval = MPI_ERR_TYPE;
2452 size = smpi_datatype_size(datatype);
2455 } else if (status->count % size != 0) {
2456 retval = MPI_UNDEFINED;
2458 *count = smpi_mpi_get_count(status, datatype);
2465 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2469 if (old_type == MPI_DATATYPE_NULL) {
2470 retval = MPI_ERR_TYPE;
2471 } else if (count<0){
2472 retval = MPI_ERR_COUNT;
2474 retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2480 int PMPI_Type_commit(MPI_Datatype* datatype) {
2484 if (datatype == NULL || *datatype == MPI_DATATYPE_NULL) {
2485 retval = MPI_ERR_TYPE;
2487 smpi_datatype_commit(datatype);
2488 retval = MPI_SUCCESS;
2495 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2499 if (old_type == MPI_DATATYPE_NULL) {
2500 retval = MPI_ERR_TYPE;
2501 } else if (count<0 || blocklen<0){
2502 retval = MPI_ERR_COUNT;
2504 retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2510 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2514 if (old_type == MPI_DATATYPE_NULL) {
2515 retval = MPI_ERR_TYPE;
2516 } else if (count<0 || blocklen<0){
2517 retval = MPI_ERR_COUNT;
2519 retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2525 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2526 return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2529 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2533 if (old_type == MPI_DATATYPE_NULL) {
2534 retval = MPI_ERR_TYPE;
2535 } else if (count<0){
2536 retval = MPI_ERR_COUNT;
2538 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2544 int PMPI_Type_create_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2548 if (old_type == MPI_DATATYPE_NULL) {
2549 retval = MPI_ERR_TYPE;
2550 } else if (count<0){
2551 retval = MPI_ERR_COUNT;
2553 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2559 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2563 if (old_type == MPI_DATATYPE_NULL) {
2564 retval = MPI_ERR_TYPE;
2565 } else if (count<0){
2566 retval = MPI_ERR_COUNT;
2568 int* blocklens=(int*)xbt_malloc(blocklength*count);
2569 for (i=0; i<count;i++)blocklens[i]=blocklength;
2570 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2571 xbt_free(blocklens);
2578 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2582 if (old_type == MPI_DATATYPE_NULL) {
2583 retval = MPI_ERR_TYPE;
2584 } else if (count<0){
2585 retval = MPI_ERR_COUNT;
2587 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2593 int PMPI_Type_create_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2594 return PMPI_Type_hindexed(count, blocklens,indices,old_type,new_type);
2597 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2601 if (old_type == MPI_DATATYPE_NULL) {
2602 retval = MPI_ERR_TYPE;
2603 } else if (count<0){
2604 retval = MPI_ERR_COUNT;
2606 int* blocklens=(int*)xbt_malloc(blocklength*count);
2607 for (i=0; i<count;i++)blocklens[i]=blocklength;
2608 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2609 xbt_free(blocklens);
2616 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2621 retval = MPI_ERR_COUNT;
2623 retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2629 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2630 return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2634 int PMPI_Error_class(int errorcode, int* errorclass) {
2635 // assume smpi uses only standard mpi error codes
2636 *errorclass=errorcode;
2641 int PMPI_Initialized(int* flag) {
2642 *flag=smpi_process_initialized();
2646 /* The following calls are not yet implemented and will fail at runtime. */
2647 /* Once implemented, please move them above this notice. */
2649 #define NOT_YET_IMPLEMENTED {\
2650 XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__);\
2651 return MPI_SUCCESS;\
2655 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2659 int PMPI_Type_set_name(MPI_Datatype datatype, char * name)
2664 int PMPI_Type_get_name(MPI_Datatype datatype, char * name, int* len)
2669 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2673 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2677 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
2681 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2685 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2689 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2693 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2697 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2701 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2705 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2709 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2713 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2717 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2721 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2725 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2729 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2733 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2737 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2741 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2745 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2749 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2753 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2757 int PMPI_Comm_get_errhandler(MPI_Comm comm, MPI_Errhandler* errhandler) {
2761 int PMPI_Cancel(MPI_Request* request) {
2765 int PMPI_Buffer_attach(void* buffer, int size) {
2769 int PMPI_Buffer_detach(void* buffer, int* size) {
2773 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2777 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2782 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2787 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2792 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2797 int PMPI_Comm_free_keyval(int* keyval) {
2801 int PMPI_Pcontrol(const int level )
2806 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2810 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2815 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2820 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2825 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2830 int PMPI_Type_free_keyval(int* keyval) {
2834 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2838 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2842 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2846 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2850 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2854 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2858 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2862 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2866 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2870 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2874 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2878 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2882 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2886 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2890 int PMPI_Keyval_free(int* keyval) {
2894 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2898 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2902 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2906 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2910 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2914 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2918 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2922 int PMPI_Win_fence( int assert, MPI_Win win){
2926 int PMPI_Win_free( MPI_Win* win){
2930 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2934 int PMPI_Info_create( MPI_Info *info){
2938 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2942 int PMPI_Info_free( MPI_Info *info){
2946 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2947 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2951 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2952 int *num_addresses, int *num_datatypes, int *combiner){
2956 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2957 int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2958 MPI_Datatype* array_of_datatypes){
2962 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2963 int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2964 int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2968 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2972 int PMPI_Type_create_subarray(int ndims,int *array_of_sizes, int *array_of_subsizes, int *array_of_starts, int order, MPI_Datatype oldtype, MPI_Datatype *newtype){
2976 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
2980 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
2981 void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
2986 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
2990 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
2994 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
2998 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
3002 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
3006 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
3010 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
3014 int PMPI_Add_error_class( int *errorclass){
3018 int PMPI_Add_error_code( int errorclass, int *errorcode){
3022 int PMPI_Add_error_string( int errorcode, char *string){
3026 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
3030 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
3034 int PMPI_Info_delete(MPI_Info info, char *key){
3038 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
3042 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
3046 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
3050 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
3054 int PMPI_Grequest_start( MPI_Grequest_query_function *query_fn, MPI_Grequest_free_function *free_fn, MPI_Grequest_cancel_function *cancel_fn, void *extra_state, MPI_Request *request){
3058 int PMPI_Grequest_complete( MPI_Request request){
3062 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
3066 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
3070 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
3074 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
3078 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
3082 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
3086 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
3090 int PMPI_Open_port( MPI_Info info, char *port_name){
3094 int PMPI_Close_port(char *port_name){
3098 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
3102 int PMPI_Comm_spawn( char *command, char **argv, int maxprocs, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
3106 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
3107 int* array_of_maxprocs, MPI_Info* array_of_info, int root,
3108 MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
3112 int PMPI_Comm_get_parent( MPI_Comm *parent){