1 /* Copyright (c) 2007-2013. 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_sample_is_running) {
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 if((*request)->flags & PERSISTENT)(*request)->refcount--;
1024 smpi_mpi_request_free(request);
1025 retval = MPI_SUCCESS;
1031 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
1032 int tag, MPI_Comm comm, MPI_Request * request)
1038 if (request == NULL) {
1039 retval = MPI_ERR_ARG;
1040 } else if (comm == MPI_COMM_NULL) {
1041 retval = MPI_ERR_COMM;
1042 } else if (src == MPI_PROC_NULL) {
1043 *request = MPI_REQUEST_NULL;
1044 retval = MPI_SUCCESS;
1045 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1046 retval = MPI_ERR_COMM;
1047 } else if (count < 0) {
1048 retval = MPI_ERR_COUNT;
1049 } else if (buf==NULL && count > 0) {
1050 retval = MPI_ERR_COUNT;
1051 } else if (datatype == MPI_DATATYPE_NULL){
1052 retval = MPI_ERR_TYPE;
1053 } else if(tag<0 && tag != MPI_ANY_TAG){
1054 retval = MPI_ERR_TAG;
1058 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1059 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1061 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1062 extra->type = TRACING_IRECV;
1063 extra->send_size = count;
1064 extra->src = src_traced;
1066 extra->datatype1 = encode_datatype(datatype);
1067 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
1070 *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
1071 retval = MPI_SUCCESS;
1074 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1075 (*request)->recv = 1;
1080 if (retval != MPI_SUCCESS && request)
1081 *request = MPI_REQUEST_NULL;
1086 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
1087 int tag, MPI_Comm comm, MPI_Request * request)
1092 if (request == NULL) {
1093 retval = MPI_ERR_ARG;
1094 } else if (comm == MPI_COMM_NULL) {
1095 retval = MPI_ERR_COMM;
1096 } else if (dst == MPI_PROC_NULL) {
1097 *request = MPI_REQUEST_NULL;
1098 retval = MPI_SUCCESS;
1099 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1100 retval = MPI_ERR_RANK;
1101 } else if (count < 0) {
1102 retval = MPI_ERR_COUNT;
1103 } else if (buf==NULL && count > 0) {
1104 retval = MPI_ERR_COUNT;
1105 } else if (datatype == MPI_DATATYPE_NULL){
1106 retval = MPI_ERR_TYPE;
1107 } else if(tag<0 && tag != MPI_ANY_TAG){
1108 retval = MPI_ERR_TAG;
1112 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1113 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1115 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1116 extra->type = TRACING_ISEND;
1117 extra->send_size = count;
1119 extra->dst = dst_traced;
1120 extra->datatype1 = encode_datatype(datatype);
1121 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1122 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1125 *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1126 retval = MPI_SUCCESS;
1129 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1130 (*request)->send = 1;
1135 if (retval != MPI_SUCCESS && request)
1136 *request = MPI_REQUEST_NULL;
1140 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype,
1141 int dst, int tag, MPI_Comm comm, MPI_Request* request)
1146 if (request == NULL) {
1147 retval = MPI_ERR_ARG;
1148 } else if (comm == MPI_COMM_NULL) {
1149 retval = MPI_ERR_COMM;
1150 } else if (dst == MPI_PROC_NULL) {
1151 *request = MPI_REQUEST_NULL;
1152 retval = MPI_SUCCESS;
1153 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1154 retval = MPI_ERR_RANK;
1155 } else if (count < 0) {
1156 retval = MPI_ERR_COUNT;
1157 } else if (buf==NULL && count > 0) {
1158 retval = MPI_ERR_COUNT;
1159 } else if (datatype == MPI_DATATYPE_NULL){
1160 retval = MPI_ERR_TYPE;
1161 } else if(tag<0 && tag != MPI_ANY_TAG){
1162 retval = MPI_ERR_TAG;
1166 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1167 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1168 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1169 extra->type = TRACING_ISSEND;
1170 extra->send_size = count;
1172 extra->dst = dst_traced;
1173 extra->datatype1 = encode_datatype(datatype);
1174 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1175 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1178 *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1179 retval = MPI_SUCCESS;
1182 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1183 (*request)->send = 1;
1188 if (retval != MPI_SUCCESS && request)
1189 *request = MPI_REQUEST_NULL;
1193 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1194 MPI_Comm comm, MPI_Status * status)
1199 if (comm == MPI_COMM_NULL) {
1200 retval = MPI_ERR_COMM;
1201 } else if (src == MPI_PROC_NULL) {
1202 smpi_empty_status(status);
1203 status->MPI_SOURCE = MPI_PROC_NULL;
1204 retval = MPI_SUCCESS;
1205 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1206 retval = MPI_ERR_RANK;
1207 } else if (count < 0) {
1208 retval = MPI_ERR_COUNT;
1209 } else if (buf==NULL && count > 0) {
1210 retval = MPI_ERR_COUNT;
1211 } else if (datatype == MPI_DATATYPE_NULL){
1212 retval = MPI_ERR_TYPE;
1213 } else if(tag<0 && tag != MPI_ANY_TAG){
1214 retval = MPI_ERR_TAG;
1217 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1218 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1219 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1220 extra->type = TRACING_RECV;
1221 extra->send_size = count;
1222 extra->src = src_traced;
1224 extra->datatype1 = encode_datatype(datatype);
1225 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
1228 smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1229 retval = MPI_SUCCESS;
1232 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1233 if(status!=MPI_STATUS_IGNORE)src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1234 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1235 TRACE_smpi_recv(rank, src_traced, rank);
1243 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1250 if (comm == MPI_COMM_NULL) {
1251 retval = MPI_ERR_COMM;
1252 } else if (dst == MPI_PROC_NULL) {
1253 retval = MPI_SUCCESS;
1254 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1255 retval = MPI_ERR_RANK;
1256 } else if (count < 0) {
1257 retval = MPI_ERR_COUNT;
1258 } else if (buf==NULL && count > 0) {
1259 retval = MPI_ERR_COUNT;
1260 } else if (datatype == MPI_DATATYPE_NULL){
1261 retval = MPI_ERR_TYPE;
1262 } else if(tag<0 && tag != MPI_ANY_TAG){
1263 retval = MPI_ERR_TAG;
1267 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1268 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1269 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1270 extra->type = TRACING_SEND;
1271 extra->send_size = count;
1273 extra->dst = dst_traced;
1274 extra->datatype1 = encode_datatype(datatype);
1275 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1276 TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1279 smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1280 retval = MPI_SUCCESS;
1283 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1293 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1298 if (comm == MPI_COMM_NULL) {
1299 retval = MPI_ERR_COMM;
1300 } else if (dst == MPI_PROC_NULL) {
1301 retval = MPI_SUCCESS;
1302 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1303 retval = MPI_ERR_RANK;
1304 } else if (count < 0) {
1305 retval = MPI_ERR_COUNT;
1306 } else if (buf==NULL && count > 0) {
1307 retval = MPI_ERR_COUNT;
1308 } else if (datatype == MPI_DATATYPE_NULL){
1309 retval = MPI_ERR_TYPE;
1310 } else if(tag<0 && tag != MPI_ANY_TAG){
1311 retval = MPI_ERR_TAG;
1315 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1316 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1317 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1318 extra->type = TRACING_SSEND;
1319 extra->send_size = count;
1321 extra->dst = dst_traced;
1322 extra->datatype1 = encode_datatype(datatype);
1323 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra); TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1326 smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1327 retval = MPI_SUCCESS;
1330 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1338 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1339 int dst, int sendtag, void *recvbuf, int recvcount,
1340 MPI_Datatype recvtype, int src, int recvtag,
1341 MPI_Comm comm, MPI_Status * status)
1347 if (comm == MPI_COMM_NULL) {
1348 retval = MPI_ERR_COMM;
1349 } else if (sendtype == MPI_DATATYPE_NULL
1350 || recvtype == MPI_DATATYPE_NULL) {
1351 retval = MPI_ERR_TYPE;
1352 } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1353 smpi_empty_status(status);
1354 status->MPI_SOURCE = MPI_PROC_NULL;
1355 retval = MPI_SUCCESS;
1356 }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1357 (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1358 retval = MPI_ERR_RANK;
1359 } else if (sendcount < 0 || recvcount<0) {
1360 retval = MPI_ERR_COUNT;
1361 } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1362 retval = MPI_ERR_COUNT;
1363 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1364 retval = MPI_ERR_TAG;
1368 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1369 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1370 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1371 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1372 extra->type = TRACING_SENDRECV;
1373 extra->send_size = sendcount;
1374 extra->recv_size = recvcount;
1375 extra->src = src_traced;
1376 extra->dst = dst_traced;
1377 extra->datatype1 = encode_datatype(sendtype);
1378 extra->datatype2 = encode_datatype(recvtype);
1380 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1381 TRACE_smpi_send(rank, rank, dst_traced,sendcount*smpi_datatype_size(sendtype));
1385 smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1386 recvcount, recvtype, src, recvtag, comm, status);
1387 retval = MPI_SUCCESS;
1390 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1391 TRACE_smpi_recv(rank, src_traced, rank);
1400 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1401 int dst, int sendtag, int src, int recvtag,
1402 MPI_Comm comm, MPI_Status * status)
1404 //TODO: suboptimal implementation
1407 if (datatype == MPI_DATATYPE_NULL) {
1408 retval = MPI_ERR_TYPE;
1409 } else if (count < 0) {
1410 retval = MPI_ERR_COUNT;
1412 int size = smpi_datatype_get_extent(datatype) * count;
1413 recvbuf = xbt_new0(char, size);
1415 MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1416 datatype, src, recvtag, comm, status);
1417 if(retval==MPI_SUCCESS){
1418 smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1426 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1431 if (request == NULL || flag == NULL) {
1432 retval = MPI_ERR_ARG;
1433 } else if (*request == MPI_REQUEST_NULL) {
1435 smpi_empty_status(status);
1436 retval = MPI_ERR_REQUEST;
1438 *flag = smpi_mpi_test(request, status);
1439 retval = MPI_SUCCESS;
1445 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1446 MPI_Status * status)
1451 if (index == NULL || flag == NULL) {
1452 retval = MPI_ERR_ARG;
1454 *flag = smpi_mpi_testany(count, requests, index, status);
1455 retval = MPI_SUCCESS;
1461 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1467 retval = MPI_ERR_ARG;
1469 *flag = smpi_mpi_testall(count, requests, statuses);
1470 retval = MPI_SUCCESS;
1476 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1480 if (status == NULL) {
1481 retval = MPI_ERR_ARG;
1482 } else if (comm == MPI_COMM_NULL) {
1483 retval = MPI_ERR_COMM;
1484 } else if (source == MPI_PROC_NULL) {
1485 smpi_empty_status(status);
1486 status->MPI_SOURCE = MPI_PROC_NULL;
1487 retval = MPI_SUCCESS;
1489 smpi_mpi_probe(source, tag, comm, status);
1490 retval = MPI_SUCCESS;
1497 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1502 retval = MPI_ERR_ARG;
1503 } else if (status == NULL) {
1504 retval = MPI_ERR_ARG;
1505 } else if (comm == MPI_COMM_NULL) {
1506 retval = MPI_ERR_COMM;
1507 } else if (source == MPI_PROC_NULL) {
1509 smpi_empty_status(status);
1510 status->MPI_SOURCE = MPI_PROC_NULL;
1511 retval = MPI_SUCCESS;
1513 smpi_mpi_iprobe(source, tag, comm, flag, status);
1514 retval = MPI_SUCCESS;
1520 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1526 smpi_empty_status(status);
1528 if (request == NULL) {
1529 retval = MPI_ERR_ARG;
1530 } else if (*request == MPI_REQUEST_NULL) {
1531 retval = MPI_ERR_REQUEST;
1535 int rank = request && (*request)->comm != MPI_COMM_NULL
1536 ? smpi_process_index()
1539 int src_traced = (*request)->src;
1540 int dst_traced = (*request)->dst;
1541 MPI_Comm comm = (*request)->comm;
1542 int is_wait_for_receive = (*request)->recv;
1543 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1544 extra->type = TRACING_WAIT;
1545 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1548 smpi_mpi_wait(request, status);
1549 retval = MPI_SUCCESS;
1552 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1553 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1554 if (is_wait_for_receive) {
1555 if(src_traced==MPI_ANY_SOURCE)
1556 src_traced = (status!=MPI_STATUS_IGNORE) ?
1557 smpi_group_rank(smpi_comm_group(comm), status->MPI_SOURCE) :
1559 TRACE_smpi_recv(rank, src_traced, dst_traced);
1569 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1575 //save requests information for tracing
1577 int *srcs = xbt_new0(int, count);
1578 int *dsts = xbt_new0(int, count);
1579 int *recvs = xbt_new0(int, count);
1580 MPI_Comm *comms = xbt_new0(MPI_Comm, count);
1582 for (i = 0; i < count; i++) {
1583 MPI_Request req = requests[i]; //already received requests are no longer valid
1587 recvs[i] = req->recv;
1588 comms[i] = req->comm;
1591 int rank_traced = smpi_process_index();
1592 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1593 extra->type = TRACING_WAITANY;
1594 extra->send_size=count;
1595 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1598 if (index == NULL) {
1599 retval = MPI_ERR_ARG;
1601 *index = smpi_mpi_waitany(count, requests, status);
1602 retval = MPI_SUCCESS;
1605 if(*index!=MPI_UNDEFINED){
1606 int src_traced = srcs[*index];
1607 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1608 int dst_traced = dsts[*index];
1609 int is_wait_for_receive = recvs[*index];
1610 if (is_wait_for_receive) {
1611 if(srcs[*index]==MPI_ANY_SOURCE)
1612 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1613 smpi_group_rank(smpi_comm_group(comms[*index]), status->MPI_SOURCE) :
1615 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1617 TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1629 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1634 //save information from requests
1636 int *srcs = xbt_new0(int, count);
1637 int *dsts = xbt_new0(int, count);
1638 int *recvs = xbt_new0(int, count);
1639 int *valid = xbt_new0(int, count);
1640 MPI_Comm *comms = xbt_new0(MPI_Comm, count);
1642 //int valid_count = 0;
1643 for (i = 0; i < count; i++) {
1644 MPI_Request req = requests[i];
1645 if(req!=MPI_REQUEST_NULL){
1648 recvs[i] = req->recv;
1649 comms[i] = req->comm;
1655 int rank_traced = smpi_process_index();
1656 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1657 extra->type = TRACING_WAITALL;
1658 extra->send_size=count;
1659 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1661 int retval = smpi_mpi_waitall(count, requests, status);
1663 for (i = 0; i < count; i++) {
1665 //int src_traced = srcs[*index];
1666 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1667 int src_traced = srcs[i];
1668 int dst_traced = dsts[i];
1669 int is_wait_for_receive = recvs[i];
1670 if (is_wait_for_receive) {
1671 if(src_traced==MPI_ANY_SOURCE)
1672 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1673 smpi_group_rank(smpi_comm_group(comms[i]), status[i].MPI_SOURCE) :
1675 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1679 TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1691 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1692 int *indices, MPI_Status status[])
1697 if (outcount == NULL) {
1698 retval = MPI_ERR_ARG;
1700 *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1701 retval = MPI_SUCCESS;
1707 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1708 int* indices, MPI_Status status[])
1713 if (outcount == NULL) {
1714 retval = MPI_ERR_ARG;
1716 *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1717 retval = MPI_SUCCESS;
1724 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1730 if (comm == MPI_COMM_NULL) {
1731 retval = MPI_ERR_COMM;
1734 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1735 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1737 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1738 extra->type = TRACING_BCAST;
1739 extra->send_size = count;
1740 extra->root = root_traced;
1741 extra->datatype1 = encode_datatype(datatype);
1742 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1745 mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1746 retval = MPI_SUCCESS;
1748 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1756 int PMPI_Barrier(MPI_Comm comm)
1762 if (comm == MPI_COMM_NULL) {
1763 retval = MPI_ERR_COMM;
1766 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1767 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1768 extra->type = TRACING_BARRIER;
1769 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1771 mpi_coll_barrier_fun(comm);
1772 retval = MPI_SUCCESS;
1774 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1782 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1783 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1784 int root, MPI_Comm comm)
1790 if (comm == MPI_COMM_NULL) {
1791 retval = MPI_ERR_COMM;
1792 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1793 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1794 retval = MPI_ERR_TYPE;
1795 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1796 ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1797 retval = MPI_ERR_COUNT;
1800 char* sendtmpbuf = (char*) sendbuf;
1801 int sendtmpcount = sendcount;
1802 MPI_Datatype sendtmptype = sendtype;
1803 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1805 sendtmptype=recvtype;
1808 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1809 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1810 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1811 extra->type = TRACING_GATHER;
1812 extra->send_size = sendtmpcount;
1813 extra->recv_size = recvcount;
1814 extra->root = root_traced;
1815 extra->datatype1 = encode_datatype(sendtmptype);
1816 extra->datatype2 = encode_datatype(recvtype);
1818 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1820 mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount,
1821 recvtype, root, comm);
1824 retval = MPI_SUCCESS;
1826 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1834 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1835 void *recvbuf, int *recvcounts, int *displs,
1836 MPI_Datatype recvtype, int root, MPI_Comm comm)
1842 if (comm == MPI_COMM_NULL) {
1843 retval = MPI_ERR_COMM;
1844 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1845 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1846 retval = MPI_ERR_TYPE;
1847 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1848 retval = MPI_ERR_COUNT;
1849 } else if (recvcounts == NULL || displs == NULL) {
1850 retval = MPI_ERR_ARG;
1852 char* sendtmpbuf = (char*) sendbuf;
1853 int sendtmpcount = sendcount;
1854 MPI_Datatype sendtmptype = sendtype;
1855 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1857 sendtmptype=recvtype;
1861 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1862 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1864 int size = smpi_comm_size(comm);
1865 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1866 extra->type = TRACING_GATHERV;
1867 extra->send_size = sendtmpcount;
1868 extra->recvcounts= xbt_malloc(size*sizeof(int));
1869 for(i=0; i< size; i++)//copy data to avoid bad free
1870 extra->recvcounts[i] = recvcounts[i];
1871 extra->num_processes = size;
1872 extra->root = root_traced;
1873 extra->datatype1 = encode_datatype(sendtmptype);
1874 extra->datatype2 = encode_datatype(recvtype);
1876 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1878 smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts,
1879 displs, recvtype, root, comm);
1880 retval = MPI_SUCCESS;
1882 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1890 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1891 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1898 if (comm == MPI_COMM_NULL) {
1899 retval = MPI_ERR_COMM;
1900 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1901 (recvtype == MPI_DATATYPE_NULL)){
1902 retval = MPI_ERR_TYPE;
1903 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1905 retval = MPI_ERR_COUNT;
1907 if(sendbuf == MPI_IN_PLACE) {
1908 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1909 sendcount=recvcount;
1913 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1914 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1915 extra->type = TRACING_ALLGATHER;
1916 extra->send_size = sendcount;
1917 extra->recv_size = recvcount;
1918 extra->datatype1 = encode_datatype(sendtype);
1919 extra->datatype2 = encode_datatype(recvtype);
1921 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1923 mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1925 retval = MPI_SUCCESS;
1928 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1935 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1936 void *recvbuf, int *recvcounts, int *displs,
1937 MPI_Datatype recvtype, MPI_Comm comm)
1943 if (comm == MPI_COMM_NULL) {
1944 retval = MPI_ERR_COMM;
1945 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1946 (recvtype == MPI_DATATYPE_NULL)){
1947 retval = MPI_ERR_TYPE;
1948 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1949 retval = MPI_ERR_COUNT;
1950 } else if (recvcounts == NULL || displs == NULL) {
1951 retval = MPI_ERR_ARG;
1954 if(sendbuf == MPI_IN_PLACE) {
1955 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1956 sendcount=recvcounts[smpi_comm_rank(comm)];
1960 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1962 int size = smpi_comm_size(comm);
1963 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1964 extra->type = TRACING_ALLGATHERV;
1965 extra->send_size = sendcount;
1966 extra->recvcounts= xbt_malloc(size*sizeof(int));
1967 for(i=0; i< size; i++)//copy data to avoid bad free
1968 extra->recvcounts[i] = recvcounts[i];
1969 extra->num_processes = size;
1970 extra->datatype1 = encode_datatype(sendtype);
1971 extra->datatype2 = encode_datatype(recvtype);
1973 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
1975 mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1976 displs, recvtype, comm);
1977 retval = MPI_SUCCESS;
1979 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1987 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1988 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1989 int root, MPI_Comm comm)
1995 if (comm == MPI_COMM_NULL) {
1996 retval = MPI_ERR_COMM;
1997 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1998 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1999 retval = MPI_ERR_TYPE;
2002 if (recvbuf == MPI_IN_PLACE) {
2004 recvcount=sendcount;
2007 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2008 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2009 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2010 extra->type = TRACING_SCATTER;
2011 extra->send_size = sendcount;
2012 extra->recv_size= recvcount;
2013 extra->root = root_traced;
2014 extra->datatype1 = encode_datatype(sendtype);
2015 extra->datatype2 = encode_datatype(recvtype);
2017 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2019 mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
2020 recvtype, root, comm);
2021 retval = MPI_SUCCESS;
2023 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2031 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
2032 MPI_Datatype sendtype, void *recvbuf, int recvcount,
2033 MPI_Datatype recvtype, int root, MPI_Comm comm)
2039 if (comm == MPI_COMM_NULL) {
2040 retval = MPI_ERR_COMM;
2041 } else if (sendcounts == NULL || displs == NULL) {
2042 retval = MPI_ERR_ARG;
2043 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
2044 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
2045 retval = MPI_ERR_TYPE;
2047 if (recvbuf == MPI_IN_PLACE) {
2049 recvcount=sendcounts[smpi_comm_rank(comm)];
2052 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2053 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2055 int size = smpi_comm_size(comm);
2056 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2057 extra->type = TRACING_SCATTERV;
2058 extra->recv_size = recvcount;
2059 extra->sendcounts= xbt_malloc(size*sizeof(int));
2060 for(i=0; i< size; i++)//copy data to avoid bad free
2061 extra->sendcounts[i] = sendcounts[i];
2062 extra->num_processes = size;
2063 extra->root = root_traced;
2064 extra->datatype1 = encode_datatype(sendtype);
2065 extra->datatype2 = encode_datatype(recvtype);
2067 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2070 smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
2071 recvcount, recvtype, root, comm);
2072 retval = MPI_SUCCESS;
2074 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2082 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
2083 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
2089 if (comm == MPI_COMM_NULL) {
2090 retval = MPI_ERR_COMM;
2091 } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
2092 retval = MPI_ERR_ARG;
2095 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2096 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2097 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2098 extra->type = TRACING_REDUCE;
2099 extra->send_size = count;
2100 extra->datatype1 = encode_datatype(datatype);
2101 extra->root = root_traced;
2103 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2105 mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
2107 retval = MPI_SUCCESS;
2109 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2117 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
2118 MPI_Datatype datatype, MPI_Op op){
2122 if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
2123 retval = MPI_ERR_ARG;
2125 smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
2132 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
2133 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2139 if (comm == MPI_COMM_NULL) {
2140 retval = MPI_ERR_COMM;
2141 } else if (datatype == MPI_DATATYPE_NULL) {
2142 retval = MPI_ERR_TYPE;
2143 } else if (op == MPI_OP_NULL) {
2144 retval = MPI_ERR_OP;
2147 char* sendtmpbuf = (char*) sendbuf;
2148 if( sendbuf == MPI_IN_PLACE ) {
2149 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
2150 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
2153 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2154 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2155 extra->type = TRACING_ALLREDUCE;
2156 extra->send_size = count;
2157 extra->datatype1 = encode_datatype(datatype);
2159 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2161 mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
2163 if( sendbuf == MPI_IN_PLACE ) {
2164 xbt_free(sendtmpbuf);
2167 retval = MPI_SUCCESS;
2169 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2177 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
2178 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2184 if (comm == MPI_COMM_NULL) {
2185 retval = MPI_ERR_COMM;
2186 } else if (datatype == MPI_DATATYPE_NULL) {
2187 retval = MPI_ERR_TYPE;
2188 } else if (op == MPI_OP_NULL) {
2189 retval = MPI_ERR_OP;
2192 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2193 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2194 extra->type = TRACING_SCAN;
2195 extra->send_size = count;
2196 extra->datatype1 = encode_datatype(datatype);
2198 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2200 smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2201 retval = MPI_SUCCESS;
2203 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2211 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2212 MPI_Op op, MPI_Comm comm){
2217 if (comm == MPI_COMM_NULL) {
2218 retval = MPI_ERR_COMM;
2219 } else if (datatype == MPI_DATATYPE_NULL) {
2220 retval = MPI_ERR_TYPE;
2221 } else if (op == MPI_OP_NULL) {
2222 retval = MPI_ERR_OP;
2225 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2226 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2227 extra->type = TRACING_EXSCAN;
2228 extra->send_size = count;
2229 extra->datatype1 = encode_datatype(datatype);
2231 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2233 smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2234 retval = MPI_SUCCESS;
2236 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2244 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
2245 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2250 if (comm == MPI_COMM_NULL) {
2251 retval = MPI_ERR_COMM;
2252 } else if (datatype == MPI_DATATYPE_NULL) {
2253 retval = MPI_ERR_TYPE;
2254 } else if (op == MPI_OP_NULL) {
2255 retval = MPI_ERR_OP;
2256 } else if (recvcounts == NULL) {
2257 retval = MPI_ERR_ARG;
2260 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2262 int size = smpi_comm_size(comm);
2263 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2264 extra->type = TRACING_REDUCE_SCATTER;
2265 extra->send_size = 0;
2266 extra->recvcounts= xbt_malloc(size*sizeof(int));
2267 for(i=0; i< size; i++)//copy data to avoid bad free
2268 extra->recvcounts[i] = recvcounts[i];
2269 extra->num_processes = size;
2270 extra->datatype1 = encode_datatype(datatype);
2272 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2274 void* sendtmpbuf=sendbuf;
2275 if(sendbuf==MPI_IN_PLACE){
2279 mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts,
2280 datatype, op, comm);
2281 retval = MPI_SUCCESS;
2283 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2291 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2292 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2297 if (comm == MPI_COMM_NULL) {
2298 retval = MPI_ERR_COMM;
2299 } else if (datatype == MPI_DATATYPE_NULL) {
2300 retval = MPI_ERR_TYPE;
2301 } else if (op == MPI_OP_NULL) {
2302 retval = MPI_ERR_OP;
2303 } else if (recvcount < 0) {
2304 retval = MPI_ERR_ARG;
2306 int count=smpi_comm_size(comm);
2309 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2310 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2311 extra->type = TRACING_REDUCE_SCATTER;
2312 extra->send_size = 0;
2313 extra->recvcounts= xbt_malloc(count*sizeof(int));
2314 for(i=0; i< count; i++)//copy data to avoid bad free
2315 extra->recvcounts[i] = recvcount;
2316 extra->num_processes = count;
2317 extra->datatype1 = encode_datatype(datatype);
2319 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2321 int* recvcounts=(int*)xbt_malloc(count);
2322 for (i=0; i<count;i++)recvcounts[i]=recvcount;
2323 mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
2324 datatype, op, comm);
2325 xbt_free(recvcounts);
2326 retval = MPI_SUCCESS;
2328 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2336 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2337 void *recvbuf, int recvcount, MPI_Datatype recvtype,
2344 if (comm == MPI_COMM_NULL) {
2345 retval = MPI_ERR_COMM;
2346 } else if (sendtype == MPI_DATATYPE_NULL
2347 || recvtype == MPI_DATATYPE_NULL) {
2348 retval = MPI_ERR_TYPE;
2351 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2352 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2353 extra->type = TRACING_ALLTOALL;
2354 extra->send_size = sendcount;
2355 extra->recv_size = recvcount;
2356 extra->datatype1 = encode_datatype(sendtype);
2357 extra->datatype2 = encode_datatype(recvtype);
2359 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2361 retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2363 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2371 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2372 MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2373 int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2379 if (comm == MPI_COMM_NULL) {
2380 retval = MPI_ERR_COMM;
2381 } else if (sendtype == MPI_DATATYPE_NULL
2382 || recvtype == MPI_DATATYPE_NULL) {
2383 retval = MPI_ERR_TYPE;
2384 } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2385 || recvdisps == NULL) {
2386 retval = MPI_ERR_ARG;
2389 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2391 int size = smpi_comm_size(comm);
2392 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2393 extra->type = TRACING_ALLTOALLV;
2394 extra->send_size = 0;
2395 extra->recv_size = 0;
2396 extra->recvcounts= xbt_malloc(size*sizeof(int));
2397 extra->sendcounts= xbt_malloc(size*sizeof(int));
2399 for(i=0; i< size; i++){//copy data to avoid bad free
2400 extra->send_size += sendcounts[i];
2401 extra->recv_size += recvcounts[i];
2403 extra->sendcounts[i] = sendcounts[i];
2404 extra->recvcounts[i] = recvcounts[i];
2406 extra->num_processes = size;
2408 extra->datatype1 = encode_datatype(sendtype);
2409 extra->datatype2 = encode_datatype(recvtype);
2411 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2414 mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2415 recvbuf, recvcounts, recvdisps, recvtype,
2418 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2427 int PMPI_Get_processor_name(char *name, int *resultlen)
2429 int retval = MPI_SUCCESS;
2432 strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2433 strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2434 strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2435 MPI_MAX_PROCESSOR_NAME - 1 );
2438 MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2444 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2446 int retval = MPI_SUCCESS;
2450 if (status == NULL || count == NULL) {
2451 retval = MPI_ERR_ARG;
2452 } else if (datatype == MPI_DATATYPE_NULL) {
2453 retval = MPI_ERR_TYPE;
2455 size = smpi_datatype_size(datatype);
2458 } else if (status->count % size != 0) {
2459 retval = MPI_UNDEFINED;
2461 *count = smpi_mpi_get_count(status, datatype);
2468 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2472 if (old_type == MPI_DATATYPE_NULL) {
2473 retval = MPI_ERR_TYPE;
2474 } else if (count<0){
2475 retval = MPI_ERR_COUNT;
2477 retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2483 int PMPI_Type_commit(MPI_Datatype* datatype) {
2487 if (datatype == NULL || *datatype == MPI_DATATYPE_NULL) {
2488 retval = MPI_ERR_TYPE;
2490 smpi_datatype_commit(datatype);
2491 retval = MPI_SUCCESS;
2498 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2502 if (old_type == MPI_DATATYPE_NULL) {
2503 retval = MPI_ERR_TYPE;
2504 } else if (count<0 || blocklen<0){
2505 retval = MPI_ERR_COUNT;
2507 retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2513 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2517 if (old_type == MPI_DATATYPE_NULL) {
2518 retval = MPI_ERR_TYPE;
2519 } else if (count<0 || blocklen<0){
2520 retval = MPI_ERR_COUNT;
2522 retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2528 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2529 return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2532 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2536 if (old_type == MPI_DATATYPE_NULL) {
2537 retval = MPI_ERR_TYPE;
2538 } else if (count<0){
2539 retval = MPI_ERR_COUNT;
2541 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2547 int PMPI_Type_create_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2551 if (old_type == MPI_DATATYPE_NULL) {
2552 retval = MPI_ERR_TYPE;
2553 } else if (count<0){
2554 retval = MPI_ERR_COUNT;
2556 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2562 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2566 if (old_type == MPI_DATATYPE_NULL) {
2567 retval = MPI_ERR_TYPE;
2568 } else if (count<0){
2569 retval = MPI_ERR_COUNT;
2571 int* blocklens=(int*)xbt_malloc(blocklength*count);
2572 for (i=0; i<count;i++)blocklens[i]=blocklength;
2573 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2574 xbt_free(blocklens);
2581 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2585 if (old_type == MPI_DATATYPE_NULL) {
2586 retval = MPI_ERR_TYPE;
2587 } else if (count<0){
2588 retval = MPI_ERR_COUNT;
2590 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2596 int PMPI_Type_create_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2597 return PMPI_Type_hindexed(count, blocklens,indices,old_type,new_type);
2600 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2604 if (old_type == MPI_DATATYPE_NULL) {
2605 retval = MPI_ERR_TYPE;
2606 } else if (count<0){
2607 retval = MPI_ERR_COUNT;
2609 int* blocklens=(int*)xbt_malloc(blocklength*count);
2610 for (i=0; i<count;i++)blocklens[i]=blocklength;
2611 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2612 xbt_free(blocklens);
2619 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2624 retval = MPI_ERR_COUNT;
2626 retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2632 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2633 return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2637 int PMPI_Error_class(int errorcode, int* errorclass) {
2638 // assume smpi uses only standard mpi error codes
2639 *errorclass=errorcode;
2644 int PMPI_Initialized(int* flag) {
2645 *flag=smpi_process_initialized();
2649 /* The following calls are not yet implemented and will fail at runtime. */
2650 /* Once implemented, please move them above this notice. */
2652 #define NOT_YET_IMPLEMENTED {\
2653 XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__);\
2654 return MPI_SUCCESS;\
2658 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2662 int PMPI_Type_set_name(MPI_Datatype datatype, char * name)
2667 int PMPI_Type_get_name(MPI_Datatype datatype, char * name, int* len)
2672 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2676 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2680 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
2684 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2688 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2692 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2696 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2700 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2704 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2708 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2712 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2716 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2720 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2724 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2728 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2732 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2736 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2740 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2744 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2748 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2752 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2756 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2760 int PMPI_Comm_get_errhandler(MPI_Comm comm, MPI_Errhandler* errhandler) {
2764 int PMPI_Cancel(MPI_Request* request) {
2768 int PMPI_Buffer_attach(void* buffer, int size) {
2772 int PMPI_Buffer_detach(void* buffer, int* size) {
2776 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2780 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2785 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2790 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2795 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2800 int PMPI_Comm_free_keyval(int* keyval) {
2804 int PMPI_Pcontrol(const int level )
2809 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2813 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2818 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2823 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2828 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2833 int PMPI_Type_free_keyval(int* keyval) {
2837 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2841 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2845 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2849 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2853 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2857 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2861 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2865 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2869 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2873 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2877 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2881 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2885 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2889 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2893 int PMPI_Keyval_free(int* keyval) {
2897 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2901 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2905 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2909 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2913 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2917 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2921 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2925 int PMPI_Win_fence( int assert, MPI_Win win){
2929 int PMPI_Win_free( MPI_Win* win){
2933 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2937 int PMPI_Info_create( MPI_Info *info){
2941 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2945 int PMPI_Info_free( MPI_Info *info){
2949 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2950 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2954 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2955 int *num_addresses, int *num_datatypes, int *combiner){
2959 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2960 int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2961 MPI_Datatype* array_of_datatypes){
2965 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2966 int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2967 int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2971 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2975 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){
2979 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
2983 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
2984 void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
2989 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
2993 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
2997 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
3001 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
3005 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
3009 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
3013 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
3017 int PMPI_Add_error_class( int *errorclass){
3021 int PMPI_Add_error_code( int errorclass, int *errorcode){
3025 int PMPI_Add_error_string( int errorcode, char *string){
3029 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
3033 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
3037 int PMPI_Info_delete(MPI_Info info, char *key){
3041 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
3045 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
3049 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
3053 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
3057 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){
3061 int PMPI_Grequest_complete( MPI_Request request){
3065 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
3069 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
3073 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
3077 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
3081 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
3085 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
3089 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
3093 int PMPI_Open_port( MPI_Info info, char *port_name){
3097 int PMPI_Close_port(char *port_name){
3101 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
3105 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){
3109 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
3110 int* array_of_maxprocs, MPI_Info* array_of_info, int root,
3111 MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
3115 int PMPI_Comm_get_parent( MPI_Comm *parent){