1 /* Copyright (c) 2007, 2008, 2009, 2010. The SimGrid Team.
2 * All rights reserved. */
4 /* This program is free software; you can redistribute it and/or modify it
5 * under the terms of the license (GNU LGPL) which comes with this package. */
8 #include "smpi_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);
31 int rank = smpi_process_index();
32 TRACE_smpi_init(rank);
34 TRACE_smpi_computing_init(rank);
40 int PMPI_Finalize(void)
42 smpi_process_finalize();
45 int rank = smpi_process_index();
46 TRACE_smpi_computing_out(rank);
47 TRACE_smpi_finalize(smpi_process_index());
49 smpi_process_destroy();
53 int PMPI_Finalized(int* flag)
55 *flag=smpi_process_finalized();
59 int PMPI_Get_version (int *version,int *subversion){
60 *version = MPI_VERSION;
61 *subversion= MPI_SUBVERSION;
65 int PMPI_Get_library_version (char *version,int *len){
66 int retval = MPI_SUCCESS;
68 snprintf(version,MPI_MAX_LIBRARY_VERSION_STRING,"SMPI Version %d.%d. Copyright The Simgrid Team 2007-2013",SIMGRID_VERSION_MAJOR,
69 SIMGRID_VERSION_MINOR);
70 *len = strlen(version) > MPI_MAX_LIBRARY_VERSION_STRING ? MPI_MAX_LIBRARY_VERSION_STRING : strlen(version);
75 int PMPI_Init_thread(int *argc, char ***argv, int required, int *provided)
77 if (provided != NULL) {
78 *provided = MPI_THREAD_MULTIPLE;
80 return MPI_Init(argc, argv);
83 int PMPI_Query_thread(int *provided)
88 if (provided == NULL) {
91 *provided = MPI_THREAD_MULTIPLE;
98 int PMPI_Is_thread_main(int *flag)
104 retval = MPI_ERR_ARG;
106 *flag = smpi_process_index() == 0;
107 retval = MPI_SUCCESS;
113 int PMPI_Abort(MPI_Comm comm, int errorcode)
116 smpi_process_destroy();
118 int rank = smpi_process_index();
119 TRACE_smpi_computing_out(rank);
121 // FIXME: should kill all processes in comm instead
122 simcall_process_kill(SIMIX_process_self());
126 double PMPI_Wtime(void)
131 time = SIMIX_get_clock();
135 extern double sg_maxmin_precision;
136 double PMPI_Wtick(void)
138 return sg_maxmin_precision;
141 int PMPI_Address(void *location, MPI_Aint * address)
147 retval = MPI_ERR_ARG;
149 *address = (MPI_Aint) location;
150 retval = MPI_SUCCESS;
156 int PMPI_Get_address(void *location, MPI_Aint * address)
158 return PMPI_Address(location, address);
161 int PMPI_Type_free(MPI_Datatype * datatype)
167 retval = MPI_ERR_ARG;
169 smpi_datatype_free(datatype);
170 retval = MPI_SUCCESS;
176 int PMPI_Type_size(MPI_Datatype datatype, int *size)
181 if (datatype == MPI_DATATYPE_NULL) {
182 retval = MPI_ERR_TYPE;
183 } else if (size == NULL) {
184 retval = MPI_ERR_ARG;
186 *size = (int) smpi_datatype_size(datatype);
187 retval = MPI_SUCCESS;
193 int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
198 if (datatype == MPI_DATATYPE_NULL) {
199 retval = MPI_ERR_TYPE;
200 } else if (lb == NULL || extent == NULL) {
201 retval = MPI_ERR_ARG;
203 retval = smpi_datatype_extent(datatype, lb, extent);
209 int PMPI_Type_get_true_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
211 return PMPI_Type_get_extent(datatype, lb, extent);
214 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
219 if (datatype == MPI_DATATYPE_NULL) {
220 retval = MPI_ERR_TYPE;
221 } else if (extent == NULL) {
222 retval = MPI_ERR_ARG;
224 *extent = smpi_datatype_get_extent(datatype);
225 retval = MPI_SUCCESS;
231 int PMPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
236 if (datatype == MPI_DATATYPE_NULL) {
237 retval = MPI_ERR_TYPE;
238 } else if (disp == NULL) {
239 retval = MPI_ERR_ARG;
241 *disp = smpi_datatype_lb(datatype);
242 retval = MPI_SUCCESS;
248 int PMPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
253 if (datatype == MPI_DATATYPE_NULL) {
254 retval = MPI_ERR_TYPE;
255 } else if (disp == NULL) {
256 retval = MPI_ERR_ARG;
258 *disp = smpi_datatype_ub(datatype);
259 retval = MPI_SUCCESS;
265 int PMPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
270 if (function == NULL || op == NULL) {
271 retval = MPI_ERR_ARG;
273 *op = smpi_op_new(function, commute);
274 retval = MPI_SUCCESS;
280 int PMPI_Op_free(MPI_Op * op)
286 retval = MPI_ERR_ARG;
287 } else if (*op == MPI_OP_NULL) {
290 smpi_op_destroy(*op);
292 retval = MPI_SUCCESS;
298 int PMPI_Group_free(MPI_Group * group)
304 retval = MPI_ERR_ARG;
306 if(*group!= smpi_comm_group(MPI_COMM_WORLD))// do not free the group of the comm_world
307 smpi_group_destroy(*group);
308 *group = MPI_GROUP_NULL;
309 retval = MPI_SUCCESS;
315 int PMPI_Group_size(MPI_Group group, int *size)
320 if (group == MPI_GROUP_NULL) {
321 retval = MPI_ERR_GROUP;
322 } else if (size == NULL) {
323 retval = MPI_ERR_ARG;
325 *size = smpi_group_size(group);
326 retval = MPI_SUCCESS;
332 int PMPI_Group_rank(MPI_Group group, int *rank)
337 if (group == MPI_GROUP_NULL) {
338 retval = MPI_ERR_GROUP;
339 } else if (rank == NULL) {
340 retval = MPI_ERR_ARG;
342 *rank = smpi_group_rank(group, smpi_process_index());
343 retval = MPI_SUCCESS;
349 int PMPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1,
350 MPI_Group group2, int *ranks2)
352 int retval, i, index;
354 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
355 retval = MPI_ERR_GROUP;
357 for (i = 0; i < n; i++) {
358 if(ranks1[i]==MPI_PROC_NULL){
359 ranks2[i]=MPI_PROC_NULL;
361 index = smpi_group_index(group1, ranks1[i]);
362 ranks2[i] = smpi_group_rank(group2, index);
365 retval = MPI_SUCCESS;
371 int PMPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
376 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
377 retval = MPI_ERR_GROUP;
378 } else if (result == NULL) {
379 retval = MPI_ERR_ARG;
381 *result = smpi_group_compare(group1, group2);
382 retval = MPI_SUCCESS;
388 int PMPI_Group_union(MPI_Group group1, MPI_Group group2,
389 MPI_Group * newgroup)
391 int retval, i, proc1, proc2, size, size2;
394 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
395 retval = MPI_ERR_GROUP;
396 } else if (newgroup == NULL) {
397 retval = MPI_ERR_ARG;
399 size = smpi_group_size(group1);
400 size2 = smpi_group_size(group2);
401 for (i = 0; i < size2; i++) {
402 proc2 = smpi_group_index(group2, i);
403 proc1 = smpi_group_rank(group1, proc2);
404 if (proc1 == MPI_UNDEFINED) {
409 *newgroup = MPI_GROUP_EMPTY;
411 *newgroup = smpi_group_new(size);
412 size2 = smpi_group_size(group1);
413 for (i = 0; i < size2; i++) {
414 proc1 = smpi_group_index(group1, i);
415 smpi_group_set_mapping(*newgroup, proc1, i);
417 for (i = size2; i < size; i++) {
418 proc2 = smpi_group_index(group2, i - size2);
419 smpi_group_set_mapping(*newgroup, proc2, i);
422 smpi_group_use(*newgroup);
423 retval = MPI_SUCCESS;
429 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
430 MPI_Group * newgroup)
432 int retval, i, proc1, proc2, size;
435 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
436 retval = MPI_ERR_GROUP;
437 } else if (newgroup == NULL) {
438 retval = MPI_ERR_ARG;
440 size = smpi_group_size(group2);
441 for (i = 0; i < size; i++) {
442 proc2 = smpi_group_index(group2, i);
443 proc1 = smpi_group_rank(group1, proc2);
444 if (proc1 == MPI_UNDEFINED) {
449 *newgroup = MPI_GROUP_EMPTY;
451 *newgroup = smpi_group_new(size);
453 for (i = 0; i < smpi_group_size(group2); i++) {
454 proc2 = smpi_group_index(group2, i);
455 proc1 = smpi_group_rank(group1, proc2);
456 if (proc1 != MPI_UNDEFINED) {
457 smpi_group_set_mapping(*newgroup, proc2, j);
462 smpi_group_use(*newgroup);
463 retval = MPI_SUCCESS;
469 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
471 int retval, i, proc1, proc2, size, size2;
474 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
475 retval = MPI_ERR_GROUP;
476 } else if (newgroup == NULL) {
477 retval = MPI_ERR_ARG;
479 size = size2 = smpi_group_size(group1);
480 for (i = 0; i < size2; i++) {
481 proc1 = smpi_group_index(group1, i);
482 proc2 = smpi_group_rank(group2, proc1);
483 if (proc2 != MPI_UNDEFINED) {
488 *newgroup = MPI_GROUP_EMPTY;
490 *newgroup = smpi_group_new(size);
491 for (i = 0; i < size2; i++) {
492 proc1 = smpi_group_index(group1, i);
493 proc2 = smpi_group_rank(group2, proc1);
494 if (proc2 == MPI_UNDEFINED) {
495 smpi_group_set_mapping(*newgroup, proc1, i);
499 smpi_group_use(*newgroup);
500 retval = MPI_SUCCESS;
506 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
508 int retval, i, index;
511 if (group == MPI_GROUP_NULL) {
512 retval = MPI_ERR_GROUP;
513 } else if (newgroup == NULL) {
514 retval = MPI_ERR_ARG;
517 *newgroup = MPI_GROUP_EMPTY;
518 } else if (n == smpi_group_size(group)) {
521 *newgroup = smpi_group_new(n);
522 for (i = 0; i < n; i++) {
523 index = smpi_group_index(group, ranks[i]);
524 smpi_group_set_mapping(*newgroup, index, i);
527 smpi_group_use(*newgroup);
528 retval = MPI_SUCCESS;
534 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
536 int retval, i, j, newsize, oldsize, index;
539 if (group == MPI_GROUP_NULL) {
540 retval = MPI_ERR_GROUP;
541 } else if (newgroup == NULL) {
542 retval = MPI_ERR_ARG;
546 } else if (n == smpi_group_size(group)) {
547 *newgroup = MPI_GROUP_EMPTY;
549 oldsize=smpi_group_size(group);
550 newsize = oldsize - n;
551 *newgroup = smpi_group_new(newsize);
553 int* to_exclude=xbt_new(int, smpi_group_size(group));
554 for(i=0; i<oldsize; i++)
557 to_exclude[ranks[i]]=1;
560 for(i=0; i<oldsize; i++){
561 if(to_exclude[i]==0){
562 index = smpi_group_index(group, i);
563 smpi_group_set_mapping(*newgroup, index, j);
568 xbt_free(to_exclude);
570 smpi_group_use(*newgroup);
571 retval = MPI_SUCCESS;
577 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
578 MPI_Group * newgroup)
580 int retval, i, j, rank, size, index;
583 if (group == MPI_GROUP_NULL) {
584 retval = MPI_ERR_GROUP;
585 } else if (newgroup == NULL) {
586 retval = MPI_ERR_ARG;
589 *newgroup = MPI_GROUP_EMPTY;
592 for (i = 0; i < n; i++) {
593 for (rank = ranges[i][0]; /* First */
594 rank >= 0; /* Last */
598 rank += ranges[i][2]; /* Stride */
599 if (ranges[i][0]<ranges[i][1]){
600 if(rank > ranges[i][1])
603 if(rank < ranges[i][1])
609 *newgroup = smpi_group_new(size);
611 for (i = 0; i < n; i++) {
612 for (rank = ranges[i][0]; /* First */
613 rank >= 0; /* Last */
615 index = smpi_group_index(group, rank);
616 smpi_group_set_mapping(*newgroup, index, j);
618 rank += ranges[i][2]; /* Stride */
619 if (ranges[i][0]<ranges[i][1]){
620 if(rank > ranges[i][1])
623 if(rank < ranges[i][1])
630 smpi_group_use(*newgroup);
631 retval = MPI_SUCCESS;
637 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
638 MPI_Group * newgroup)
640 int retval, i, newrank, rank, size, index, add;
643 if (group == MPI_GROUP_NULL) {
644 retval = MPI_ERR_GROUP;
645 } else if (newgroup == NULL) {
646 retval = MPI_ERR_ARG;
651 size = smpi_group_size(group);
652 for (i = 0; i < n; i++) {
653 for (rank = ranges[i][0]; /* First */
654 rank >= 0 && rank <= ranges[i][1]; /* Last */
655 rank += ranges[i][2] /* Stride */ ) {
660 *newgroup = MPI_GROUP_EMPTY;
662 *newgroup = smpi_group_new(size);
664 while (newrank < size) {
665 for (i = 0; i < n; i++) {
667 for (rank = ranges[i][0]; /* First */
668 rank >= 0 && rank <= ranges[i][1]; /* Last */
669 rank += ranges[i][2] /* Stride */ ) {
670 if (rank == newrank) {
676 index = smpi_group_index(group, newrank);
677 smpi_group_set_mapping(*newgroup, index, newrank);
680 newrank++; //added to avoid looping, need to be checked ..
684 smpi_group_use(*newgroup);
685 retval = MPI_SUCCESS;
691 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
696 if (comm == MPI_COMM_NULL) {
697 retval = MPI_ERR_COMM;
698 } else if (rank == NULL) {
699 retval = MPI_ERR_ARG;
701 *rank = smpi_comm_rank(comm);
702 retval = MPI_SUCCESS;
708 int PMPI_Comm_size(MPI_Comm comm, int *size)
713 if (comm == MPI_COMM_NULL) {
714 retval = MPI_ERR_COMM;
715 } else if (size == NULL) {
716 retval = MPI_ERR_ARG;
718 *size = smpi_comm_size(comm);
719 retval = MPI_SUCCESS;
725 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
730 if (comm == MPI_COMM_NULL) {
731 retval = MPI_ERR_COMM;
732 } else if (name == NULL || len == NULL) {
733 retval = MPI_ERR_ARG;
735 smpi_comm_get_name(comm, name, len);
736 retval = MPI_SUCCESS;
742 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
747 if (comm == MPI_COMM_NULL) {
748 retval = MPI_ERR_COMM;
749 } else if (group == NULL) {
750 retval = MPI_ERR_ARG;
752 *group = smpi_comm_group(comm);
753 retval = MPI_SUCCESS;
759 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
764 if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
765 retval = MPI_ERR_COMM;
766 } else if (result == NULL) {
767 retval = MPI_ERR_ARG;
769 if (comm1 == comm2) { /* Same communicators means same groups */
773 smpi_group_compare(smpi_comm_group(comm1),
774 smpi_comm_group(comm2));
775 if (*result == MPI_IDENT) {
776 *result = MPI_CONGRUENT;
779 retval = MPI_SUCCESS;
785 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
790 if (comm == MPI_COMM_NULL) {
791 retval = MPI_ERR_COMM;
792 } else if (newcomm == NULL) {
793 retval = MPI_ERR_ARG;
795 *newcomm = smpi_comm_new(smpi_comm_group(comm));
796 retval = MPI_SUCCESS;
802 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
807 if (comm == MPI_COMM_NULL) {
808 retval = MPI_ERR_COMM;
809 } else if (group == MPI_GROUP_NULL) {
810 retval = MPI_ERR_GROUP;
811 } else if (newcomm == NULL) {
812 retval = MPI_ERR_ARG;
813 } else if(smpi_group_rank(group,smpi_process_index())==MPI_UNDEFINED){
814 *newcomm= MPI_COMM_NULL;
815 retval = MPI_SUCCESS;
818 *newcomm = smpi_comm_new(group);
819 retval = MPI_SUCCESS;
825 int PMPI_Comm_free(MPI_Comm * comm)
831 retval = MPI_ERR_ARG;
832 } else if (*comm == MPI_COMM_NULL) {
833 retval = MPI_ERR_COMM;
835 smpi_comm_destroy(*comm);
836 *comm = MPI_COMM_NULL;
837 retval = MPI_SUCCESS;
843 int PMPI_Comm_disconnect(MPI_Comm * comm)
845 /* TODO: wait until all communication in comm are done */
850 retval = MPI_ERR_ARG;
851 } else if (*comm == MPI_COMM_NULL) {
852 retval = MPI_ERR_COMM;
854 smpi_comm_destroy(*comm);
855 *comm = MPI_COMM_NULL;
856 retval = MPI_SUCCESS;
862 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
867 if (comm_out == NULL) {
868 retval = MPI_ERR_ARG;
869 } else if (comm == MPI_COMM_NULL) {
870 retval = MPI_ERR_COMM;
872 *comm_out = smpi_comm_split(comm, color, key);
873 retval = MPI_SUCCESS;
879 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
880 int tag, MPI_Comm comm, MPI_Request * request)
885 if (request == NULL) {
886 retval = MPI_ERR_ARG;
887 } else if (comm == MPI_COMM_NULL) {
888 retval = MPI_ERR_COMM;
889 } else if (dst == MPI_PROC_NULL) {
890 retval = MPI_SUCCESS;
892 *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
893 retval = MPI_SUCCESS;
899 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
900 int tag, MPI_Comm comm, MPI_Request * request)
905 if (request == NULL) {
906 retval = MPI_ERR_ARG;
907 } else if (comm == MPI_COMM_NULL) {
908 retval = MPI_ERR_COMM;
909 } else if (src == MPI_PROC_NULL) {
910 retval = MPI_SUCCESS;
912 *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
913 retval = MPI_SUCCESS;
919 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
923 if (request == NULL) {
924 retval = MPI_ERR_ARG;
925 } else if (comm == MPI_COMM_NULL) {
926 retval = MPI_ERR_COMM;
927 } else if (dst == MPI_PROC_NULL) {
928 retval = MPI_SUCCESS;
930 *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
931 retval = MPI_SUCCESS;
937 int PMPI_Start(MPI_Request * request)
942 if (request == NULL || *request == MPI_REQUEST_NULL) {
943 retval = MPI_ERR_ARG;
945 smpi_mpi_start(*request);
946 retval = MPI_SUCCESS;
952 int PMPI_Startall(int count, MPI_Request * requests)
957 if (requests == NULL) {
958 retval = MPI_ERR_ARG;
960 smpi_mpi_startall(count, requests);
961 retval = MPI_SUCCESS;
967 int PMPI_Request_free(MPI_Request * request)
972 if (*request == MPI_REQUEST_NULL) {
973 retval = MPI_ERR_ARG;
975 if((*request)->flags & PERSISTENT)(*request)->refcount--;
976 smpi_mpi_request_free(request);
977 retval = MPI_SUCCESS;
983 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
984 int tag, MPI_Comm comm, MPI_Request * request)
990 if (request == NULL) {
991 retval = MPI_ERR_ARG;
992 } else if (comm == MPI_COMM_NULL) {
993 retval = MPI_ERR_COMM;
994 } else if (src == MPI_PROC_NULL) {
995 *request = MPI_REQUEST_NULL;
996 retval = MPI_SUCCESS;
997 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
998 retval = MPI_ERR_COMM;
999 } else if (count < 0) {
1000 retval = MPI_ERR_COUNT;
1001 } else if (buf==NULL && count > 0) {
1002 retval = MPI_ERR_COUNT;
1003 } else if (datatype == MPI_DATATYPE_NULL){
1004 retval = MPI_ERR_TYPE;
1005 } else if(tag<0 && tag != MPI_ANY_TAG){
1006 retval = MPI_ERR_TAG;
1010 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1011 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1012 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
1015 *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
1016 retval = MPI_SUCCESS;
1019 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1020 (*request)->recv = 1;
1029 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
1030 int tag, MPI_Comm comm, MPI_Request * request)
1035 if (request == NULL) {
1036 retval = MPI_ERR_ARG;
1037 } else if (comm == MPI_COMM_NULL) {
1038 retval = MPI_ERR_COMM;
1039 } else if (dst == MPI_PROC_NULL) {
1040 *request = MPI_REQUEST_NULL;
1041 retval = MPI_SUCCESS;
1042 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1043 retval = MPI_ERR_RANK;
1044 } else if (count < 0) {
1045 retval = MPI_ERR_COUNT;
1046 } else if (buf==NULL && count > 0) {
1047 retval = MPI_ERR_COUNT;
1048 } else if (datatype == MPI_DATATYPE_NULL){
1049 retval = MPI_ERR_TYPE;
1050 } else if(tag<0 && tag != MPI_ANY_TAG){
1051 retval = MPI_ERR_TAG;
1055 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1056 TRACE_smpi_computing_out(rank);
1057 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1058 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1059 TRACE_smpi_send(rank, rank, dst_traced);
1062 *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1063 retval = MPI_SUCCESS;
1066 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1067 (*request)->send = 1;
1068 TRACE_smpi_computing_in(rank);
1076 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
1080 if (request == NULL) {
1081 retval = MPI_ERR_ARG;
1082 } else if (comm == MPI_COMM_NULL) {
1083 retval = MPI_ERR_COMM;
1084 } else if (dst == MPI_PROC_NULL) {
1085 *request = MPI_REQUEST_NULL;
1086 retval = MPI_SUCCESS;
1087 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1088 retval = MPI_ERR_RANK;
1089 } else if (count < 0) {
1090 retval = MPI_ERR_COUNT;
1091 } else if (buf==NULL && count > 0) {
1092 retval = MPI_ERR_COUNT;
1093 } else if (datatype == MPI_DATATYPE_NULL){
1094 retval = MPI_ERR_TYPE;
1095 } else if(tag<0 && tag != MPI_ANY_TAG){
1096 retval = MPI_ERR_TAG;
1100 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1101 TRACE_smpi_computing_out(rank);
1102 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1103 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1104 TRACE_smpi_send(rank, rank, dst_traced);
1107 *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1108 retval = MPI_SUCCESS;
1111 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1112 (*request)->send = 1;
1113 TRACE_smpi_computing_in(rank);
1121 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1122 MPI_Comm comm, MPI_Status * status)
1127 if (comm == MPI_COMM_NULL) {
1128 retval = MPI_ERR_COMM;
1129 } else if (src == MPI_PROC_NULL) {
1130 smpi_empty_status(status);
1131 status->MPI_SOURCE = MPI_PROC_NULL;
1132 retval = MPI_SUCCESS;
1133 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1134 retval = MPI_ERR_RANK;
1135 } else if (count < 0) {
1136 retval = MPI_ERR_COUNT;
1137 } else if (buf==NULL && count > 0) {
1138 retval = MPI_ERR_COUNT;
1139 } else if (datatype == MPI_DATATYPE_NULL){
1140 retval = MPI_ERR_TYPE;
1141 } else if(tag<0 && tag != MPI_ANY_TAG){
1142 retval = MPI_ERR_TAG;
1145 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1146 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1147 TRACE_smpi_computing_out(rank);
1149 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
1152 smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1153 retval = MPI_SUCCESS;
1156 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1157 if(status!=MPI_STATUS_IGNORE)src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1158 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1159 TRACE_smpi_recv(rank, src_traced, rank);
1160 TRACE_smpi_computing_in(rank);
1168 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1175 if (comm == MPI_COMM_NULL) {
1176 retval = MPI_ERR_COMM;
1177 } else if (dst == MPI_PROC_NULL) {
1178 retval = MPI_SUCCESS;
1179 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1180 retval = MPI_ERR_RANK;
1181 } else if (count < 0) {
1182 retval = MPI_ERR_COUNT;
1183 } else if (buf==NULL && count > 0) {
1184 retval = MPI_ERR_COUNT;
1185 } else if (datatype == MPI_DATATYPE_NULL){
1186 retval = MPI_ERR_TYPE;
1187 } else if(tag<0 && tag != MPI_ANY_TAG){
1188 retval = MPI_ERR_TAG;
1192 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1193 TRACE_smpi_computing_out(rank);
1194 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1195 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1196 TRACE_smpi_send(rank, rank, dst_traced);
1199 smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1200 retval = MPI_SUCCESS;
1203 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1204 TRACE_smpi_computing_in(rank);
1214 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1219 if (comm == MPI_COMM_NULL) {
1220 retval = MPI_ERR_COMM;
1221 } else if (dst == MPI_PROC_NULL) {
1222 retval = MPI_SUCCESS;
1223 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1224 retval = MPI_ERR_RANK;
1225 } else if (count < 0) {
1226 retval = MPI_ERR_COUNT;
1227 } else if (buf==NULL && count > 0) {
1228 retval = MPI_ERR_COUNT;
1229 } else if (datatype == MPI_DATATYPE_NULL){
1230 retval = MPI_ERR_TYPE;
1231 } else if(tag<0 && tag != MPI_ANY_TAG){
1232 retval = MPI_ERR_TAG;
1236 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1237 TRACE_smpi_computing_out(rank);
1238 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1239 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1240 TRACE_smpi_send(rank, rank, dst_traced);
1243 smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1244 retval = MPI_SUCCESS;
1247 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1248 TRACE_smpi_computing_in(rank);
1256 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1257 int dst, int sendtag, void *recvbuf, int recvcount,
1258 MPI_Datatype recvtype, int src, int recvtag,
1259 MPI_Comm comm, MPI_Status * status)
1265 if (comm == MPI_COMM_NULL) {
1266 retval = MPI_ERR_COMM;
1267 } else if (sendtype == MPI_DATATYPE_NULL
1268 || recvtype == MPI_DATATYPE_NULL) {
1269 retval = MPI_ERR_TYPE;
1270 } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1271 smpi_empty_status(status);
1272 status->MPI_SOURCE = MPI_PROC_NULL;
1273 retval = MPI_SUCCESS;
1274 }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1275 (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1276 retval = MPI_ERR_RANK;
1277 } else if (sendcount < 0 || recvcount<0) {
1278 retval = MPI_ERR_COUNT;
1279 } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1280 retval = MPI_ERR_COUNT;
1281 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1282 retval = MPI_ERR_TAG;
1286 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1287 TRACE_smpi_computing_out(rank);
1288 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1289 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1290 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1291 TRACE_smpi_send(rank, rank, dst_traced);
1295 smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1296 recvcount, recvtype, src, recvtag, comm, status);
1297 retval = MPI_SUCCESS;
1300 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1301 TRACE_smpi_recv(rank, src_traced, rank);
1302 TRACE_smpi_computing_in(rank);
1311 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1312 int dst, int sendtag, int src, int recvtag,
1313 MPI_Comm comm, MPI_Status * status)
1315 //TODO: suboptimal implementation
1318 if (datatype == MPI_DATATYPE_NULL) {
1319 retval = MPI_ERR_TYPE;
1320 } else if (count < 0) {
1321 retval = MPI_ERR_COUNT;
1323 int size = smpi_datatype_get_extent(datatype) * count;
1324 recvbuf = xbt_new(char, size);
1326 MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1327 datatype, src, recvtag, comm, status);
1328 if(retval==MPI_SUCCESS){
1329 smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1337 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1342 if (request == MPI_REQUEST_NULL || flag == NULL) {
1343 retval = MPI_ERR_ARG;
1344 } else if (*request == MPI_REQUEST_NULL) {
1346 retval = MPI_ERR_REQUEST;
1348 *flag = smpi_mpi_test(request, status);
1349 retval = MPI_SUCCESS;
1355 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1356 MPI_Status * status)
1361 if (index == NULL || flag == NULL) {
1362 retval = MPI_ERR_ARG;
1364 *flag = smpi_mpi_testany(count, requests, index, status);
1365 retval = MPI_SUCCESS;
1371 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1377 retval = MPI_ERR_ARG;
1379 *flag = smpi_mpi_testall(count, requests, statuses);
1380 retval = MPI_SUCCESS;
1386 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1390 if (status == NULL) {
1391 retval = MPI_ERR_ARG;
1392 } else if (comm == MPI_COMM_NULL) {
1393 retval = MPI_ERR_COMM;
1394 } else if (source == MPI_PROC_NULL) {
1395 smpi_empty_status(status);
1396 status->MPI_SOURCE = MPI_PROC_NULL;
1397 retval = MPI_SUCCESS;
1399 smpi_mpi_probe(source, tag, comm, status);
1400 retval = MPI_SUCCESS;
1407 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1412 retval = MPI_ERR_ARG;
1413 } else if (status == NULL) {
1414 retval = MPI_ERR_ARG;
1415 } else if (comm == MPI_COMM_NULL) {
1416 retval = MPI_ERR_COMM;
1417 } else if (source == MPI_PROC_NULL) {
1418 smpi_empty_status(status);
1419 status->MPI_SOURCE = MPI_PROC_NULL;
1420 retval = MPI_SUCCESS;
1422 smpi_mpi_iprobe(source, tag, comm, flag, status);
1423 retval = MPI_SUCCESS;
1429 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1435 if (request == NULL) {
1436 retval = MPI_ERR_ARG;
1437 } else if (*request == MPI_REQUEST_NULL) {
1438 retval = MPI_ERR_REQUEST;
1442 int rank = request && (*request)->comm != MPI_COMM_NULL
1443 ? smpi_process_index()
1445 TRACE_smpi_computing_out(rank);
1447 int src_traced = (*request)->src;
1448 int dst_traced = (*request)->dst;
1449 int is_wait_for_receive = (*request)->recv;
1450 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1453 smpi_mpi_wait(request, status);
1454 retval = MPI_SUCCESS;
1457 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1458 if (is_wait_for_receive) {
1459 TRACE_smpi_recv(rank, src_traced, dst_traced);
1461 TRACE_smpi_computing_in(rank);
1470 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1476 //save requests information for tracing
1478 int *srcs = xbt_new(int, count);
1479 int *dsts = xbt_new(int, count);
1480 int *recvs = xbt_new(int, count);
1481 for (i = 0; i < count; i++) {
1482 MPI_Request req = requests[i]; //already received requests are no longer valid
1486 recvs[i] = req->recv;
1489 int rank_traced = smpi_process_index();
1490 TRACE_smpi_computing_out(rank_traced);
1492 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1495 if (index == NULL) {
1496 retval = MPI_ERR_ARG;
1498 *index = smpi_mpi_waitany(count, requests, status);
1499 retval = MPI_SUCCESS;
1502 if(*index!=MPI_UNDEFINED){
1503 int src_traced = srcs[*index];
1504 int dst_traced = dsts[*index];
1505 int is_wait_for_receive = recvs[*index];
1506 if (is_wait_for_receive) {
1507 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1509 TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1514 TRACE_smpi_computing_in(rank_traced);
1520 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1525 //save information from requests
1527 int *srcs = xbt_new(int, count);
1528 int *dsts = xbt_new(int, count);
1529 int *recvs = xbt_new(int, count);
1530 int valid_count = 0;
1531 for (i = 0; i < count; i++) {
1532 MPI_Request req = requests[i];
1533 if(req!=MPI_REQUEST_NULL){
1534 srcs[valid_count] = req->src;
1535 dsts[valid_count] = req->dst;
1536 recvs[valid_count] = req->recv;
1540 int rank_traced = smpi_process_index();
1541 TRACE_smpi_computing_out(rank_traced);
1543 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1545 int retval = smpi_mpi_waitall(count, requests, status);
1547 for (i = 0; i < valid_count; i++) {
1548 int src_traced = srcs[i];
1549 int dst_traced = dsts[i];
1550 int is_wait_for_receive = recvs[i];
1551 if (is_wait_for_receive) {
1552 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1555 TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1559 TRACE_smpi_computing_in(rank_traced);
1565 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1566 int *indices, MPI_Status status[])
1571 if (outcount == NULL) {
1572 retval = MPI_ERR_ARG;
1574 *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1575 retval = MPI_SUCCESS;
1581 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1582 int* indices, MPI_Status status[])
1587 if (outcount == NULL) {
1588 retval = MPI_ERR_ARG;
1590 *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1591 retval = MPI_SUCCESS;
1598 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1604 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1605 TRACE_smpi_computing_out(rank);
1606 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1607 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1609 if (comm == MPI_COMM_NULL) {
1610 retval = MPI_ERR_COMM;
1612 mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1613 retval = MPI_SUCCESS;
1616 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1617 TRACE_smpi_computing_in(rank);
1623 int PMPI_Barrier(MPI_Comm comm)
1629 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1630 TRACE_smpi_computing_out(rank);
1631 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1633 if (comm == MPI_COMM_NULL) {
1634 retval = MPI_ERR_COMM;
1636 mpi_coll_barrier_fun(comm);
1637 retval = MPI_SUCCESS;
1640 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1641 TRACE_smpi_computing_in(rank);
1647 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1648 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1649 int root, MPI_Comm comm)
1655 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1656 TRACE_smpi_computing_out(rank);
1657 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1658 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1660 if (comm == MPI_COMM_NULL) {
1661 retval = MPI_ERR_COMM;
1662 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1663 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1664 retval = MPI_ERR_TYPE;
1665 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1666 ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1667 retval = MPI_ERR_COUNT;
1670 char* sendtmpbuf = (char*) sendbuf;
1671 int sendtmpcount = sendcount;
1672 MPI_Datatype sendtmptype = sendtype;
1673 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1675 sendtmptype=recvtype;
1678 mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount,
1679 recvtype, root, comm);
1682 retval = MPI_SUCCESS;
1685 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1686 TRACE_smpi_computing_in(rank);
1692 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1693 void *recvbuf, int *recvcounts, int *displs,
1694 MPI_Datatype recvtype, int root, MPI_Comm comm)
1700 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1701 TRACE_smpi_computing_out(rank);
1702 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1703 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1705 if (comm == MPI_COMM_NULL) {
1706 retval = MPI_ERR_COMM;
1707 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1708 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1709 retval = MPI_ERR_TYPE;
1710 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1711 retval = MPI_ERR_COUNT;
1712 } else if (recvcounts == NULL || displs == NULL) {
1713 retval = MPI_ERR_ARG;
1716 char* sendtmpbuf = (char*) sendbuf;
1717 int sendtmpcount = sendcount;
1718 MPI_Datatype sendtmptype = sendtype;
1719 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1721 sendtmptype=recvtype;
1724 smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts,
1725 displs, recvtype, root, comm);
1726 retval = MPI_SUCCESS;
1729 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1730 TRACE_smpi_computing_in(rank);
1736 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1737 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1744 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1745 TRACE_smpi_computing_out(rank);
1746 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1748 if (comm == MPI_COMM_NULL) {
1749 retval = MPI_ERR_COMM;
1750 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1751 (recvtype == MPI_DATATYPE_NULL)){
1752 retval = MPI_ERR_TYPE;
1753 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1755 retval = MPI_ERR_COUNT;
1758 if(sendbuf == MPI_IN_PLACE) {
1759 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1760 sendcount=recvcount;
1764 mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1766 retval = MPI_SUCCESS;
1769 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1775 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1776 void *recvbuf, int *recvcounts, int *displs,
1777 MPI_Datatype recvtype, MPI_Comm comm)
1783 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1784 TRACE_smpi_computing_out(rank);
1785 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1787 if (comm == MPI_COMM_NULL) {
1788 retval = MPI_ERR_COMM;
1789 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1790 (recvtype == MPI_DATATYPE_NULL)){
1791 retval = MPI_ERR_TYPE;
1792 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1793 retval = MPI_ERR_COUNT;
1794 } else if (recvcounts == NULL || displs == NULL) {
1795 retval = MPI_ERR_ARG;
1798 if(sendbuf == MPI_IN_PLACE) {
1799 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1800 sendcount=recvcounts[smpi_comm_rank(comm)];
1804 mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1805 displs, recvtype, comm);
1806 retval = MPI_SUCCESS;
1809 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1810 TRACE_smpi_computing_in(rank);
1816 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1817 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1818 int root, MPI_Comm comm)
1824 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1825 TRACE_smpi_computing_out(rank);
1826 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1828 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1830 if (comm == MPI_COMM_NULL) {
1831 retval = MPI_ERR_COMM;
1832 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1833 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1834 retval = MPI_ERR_TYPE;
1837 if(recvbuf==MPI_IN_PLACE){
1840 mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1841 recvtype, root, comm);
1842 retval = MPI_SUCCESS;
1845 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1846 TRACE_smpi_computing_in(rank);
1852 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1853 MPI_Datatype sendtype, void *recvbuf, int recvcount,
1854 MPI_Datatype recvtype, int root, MPI_Comm comm)
1860 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1861 TRACE_smpi_computing_out(rank);
1862 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1863 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1865 if (comm == MPI_COMM_NULL) {
1866 retval = MPI_ERR_COMM;
1867 } else if (sendcounts == NULL || displs == NULL) {
1868 retval = MPI_ERR_ARG;
1869 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1870 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1871 retval = MPI_ERR_TYPE;
1874 if(recvbuf==MPI_IN_PLACE){
1878 smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1879 recvcount, recvtype, root, comm);
1880 retval = MPI_SUCCESS;
1883 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1884 TRACE_smpi_computing_in(rank);
1890 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1891 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1897 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1898 TRACE_smpi_computing_out(rank);
1899 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1900 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1902 if (comm == MPI_COMM_NULL) {
1903 retval = MPI_ERR_COMM;
1904 } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1905 retval = MPI_ERR_ARG;
1908 char* sendtmpbuf = (char*) sendbuf;
1909 if( sendbuf == MPI_IN_PLACE ) {
1910 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
1911 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
1914 mpi_coll_reduce_fun(sendtmpbuf, recvbuf, count, datatype, op, root, comm);
1916 if( sendbuf == MPI_IN_PLACE ) {
1917 xbt_free(sendtmpbuf);
1920 retval = MPI_SUCCESS;
1923 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1924 TRACE_smpi_computing_in(rank);
1930 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
1931 MPI_Datatype datatype, MPI_Op op){
1935 if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1936 retval = MPI_ERR_ARG;
1938 smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
1945 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1946 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1952 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1953 TRACE_smpi_computing_out(rank);
1954 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1956 if (comm == MPI_COMM_NULL) {
1957 retval = MPI_ERR_COMM;
1958 } else if (datatype == MPI_DATATYPE_NULL) {
1959 retval = MPI_ERR_TYPE;
1960 } else if (op == MPI_OP_NULL) {
1961 retval = MPI_ERR_OP;
1964 char* sendtmpbuf = (char*) sendbuf;
1965 if( sendbuf == MPI_IN_PLACE ) {
1966 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
1967 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
1970 mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
1972 if( sendbuf == MPI_IN_PLACE ) {
1973 xbt_free(sendtmpbuf);
1976 retval = MPI_SUCCESS;
1980 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1981 TRACE_smpi_computing_in(rank);
1987 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
1988 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1994 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1995 TRACE_smpi_computing_out(rank);
1996 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1998 if (comm == MPI_COMM_NULL) {
1999 retval = MPI_ERR_COMM;
2000 } else if (datatype == MPI_DATATYPE_NULL) {
2001 retval = MPI_ERR_TYPE;
2002 } else if (op == MPI_OP_NULL) {
2003 retval = MPI_ERR_OP;
2005 smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2006 retval = MPI_SUCCESS;
2009 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2010 TRACE_smpi_computing_in(rank);
2016 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2017 MPI_Op op, MPI_Comm comm){
2022 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2023 TRACE_smpi_computing_out(rank);
2024 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2026 if (comm == MPI_COMM_NULL) {
2027 retval = MPI_ERR_COMM;
2028 } else if (datatype == MPI_DATATYPE_NULL) {
2029 retval = MPI_ERR_TYPE;
2030 } else if (op == MPI_OP_NULL) {
2031 retval = MPI_ERR_OP;
2033 smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2034 retval = MPI_SUCCESS;
2037 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2038 TRACE_smpi_computing_in(rank);
2044 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
2045 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2050 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2051 TRACE_smpi_computing_out(rank);
2052 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2054 if (comm == MPI_COMM_NULL) {
2055 retval = MPI_ERR_COMM;
2056 } else if (datatype == MPI_DATATYPE_NULL) {
2057 retval = MPI_ERR_TYPE;
2058 } else if (op == MPI_OP_NULL) {
2059 retval = MPI_ERR_OP;
2060 } else if (recvcounts == NULL) {
2061 retval = MPI_ERR_ARG;
2063 void* sendtmpbuf=sendbuf;
2064 if(sendbuf==MPI_IN_PLACE){
2068 mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts,
2069 datatype, op, comm);
2070 retval = MPI_SUCCESS;
2073 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2074 TRACE_smpi_computing_in(rank);
2080 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2081 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2086 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2087 TRACE_smpi_computing_out(rank);
2088 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2090 if (comm == MPI_COMM_NULL) {
2091 retval = MPI_ERR_COMM;
2092 } else if (datatype == MPI_DATATYPE_NULL) {
2093 retval = MPI_ERR_TYPE;
2094 } else if (op == MPI_OP_NULL) {
2095 retval = MPI_ERR_OP;
2096 } else if (recvcount < 0) {
2097 retval = MPI_ERR_ARG;
2099 int count=smpi_comm_size(comm);
2100 int* recvcounts=(int*)xbt_malloc(count);
2101 for (i=0; i<count;i++)recvcounts[i]=recvcount;
2102 mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
2103 datatype, op, comm);
2104 xbt_free(recvcounts);
2105 retval = MPI_SUCCESS;
2108 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2109 TRACE_smpi_computing_in(rank);
2115 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2116 void *recvbuf, int recvcount, MPI_Datatype recvtype,
2123 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2124 TRACE_smpi_computing_out(rank);
2125 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2127 if (comm == MPI_COMM_NULL) {
2128 retval = MPI_ERR_COMM;
2129 } else if (sendtype == MPI_DATATYPE_NULL
2130 || recvtype == MPI_DATATYPE_NULL) {
2131 retval = MPI_ERR_TYPE;
2133 retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2136 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2137 TRACE_smpi_computing_in(rank);
2143 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2144 MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2145 int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2151 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2152 TRACE_smpi_computing_out(rank);
2153 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2155 if (comm == MPI_COMM_NULL) {
2156 retval = MPI_ERR_COMM;
2157 } else if (sendtype == MPI_DATATYPE_NULL
2158 || recvtype == MPI_DATATYPE_NULL) {
2159 retval = MPI_ERR_TYPE;
2160 } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2161 || recvdisps == NULL) {
2162 retval = MPI_ERR_ARG;
2165 mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2166 recvbuf, recvcounts, recvdisps, recvtype,
2170 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2171 TRACE_smpi_computing_in(rank);
2178 int PMPI_Get_processor_name(char *name, int *resultlen)
2180 int retval = MPI_SUCCESS;
2183 strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2184 strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2185 strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2186 MPI_MAX_PROCESSOR_NAME - 1 );
2189 MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2195 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2197 int retval = MPI_SUCCESS;
2201 if (status == NULL || count == NULL) {
2202 retval = MPI_ERR_ARG;
2203 } else if (datatype == MPI_DATATYPE_NULL) {
2204 retval = MPI_ERR_TYPE;
2206 size = smpi_datatype_size(datatype);
2209 } else if (status->count % size != 0) {
2210 retval = MPI_UNDEFINED;
2212 *count = smpi_mpi_get_count(status, datatype);
2219 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2223 if (old_type == MPI_DATATYPE_NULL) {
2224 retval = MPI_ERR_TYPE;
2225 } else if (count<0){
2226 retval = MPI_ERR_COUNT;
2228 retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2234 int PMPI_Type_commit(MPI_Datatype* datatype) {
2238 if (datatype == MPI_DATATYPE_NULL) {
2239 retval = MPI_ERR_TYPE;
2241 smpi_datatype_commit(datatype);
2242 retval = MPI_SUCCESS;
2249 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2253 if (old_type == MPI_DATATYPE_NULL) {
2254 retval = MPI_ERR_TYPE;
2255 } else if (count<0 || blocklen<0){
2256 retval = MPI_ERR_COUNT;
2258 retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2264 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2268 if (old_type == MPI_DATATYPE_NULL) {
2269 retval = MPI_ERR_TYPE;
2270 } else if (count<0 || blocklen<0){
2271 retval = MPI_ERR_COUNT;
2273 retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2279 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2280 return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2283 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2287 if (old_type == MPI_DATATYPE_NULL) {
2288 retval = MPI_ERR_TYPE;
2289 } else if (count<0){
2290 retval = MPI_ERR_COUNT;
2292 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2298 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2302 if (old_type == MPI_DATATYPE_NULL) {
2303 retval = MPI_ERR_TYPE;
2304 } else if (count<0){
2305 retval = MPI_ERR_COUNT;
2307 int* blocklens=(int*)xbt_malloc(blocklength*count);
2308 for (i=0; i<count;i++)blocklens[i]=blocklength;
2309 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2310 xbt_free(blocklens);
2317 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2321 if (old_type == MPI_DATATYPE_NULL) {
2322 retval = MPI_ERR_TYPE;
2323 } else if (count<0){
2324 retval = MPI_ERR_COUNT;
2326 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2332 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2336 if (old_type == MPI_DATATYPE_NULL) {
2337 retval = MPI_ERR_TYPE;
2338 } else if (count<0){
2339 retval = MPI_ERR_COUNT;
2341 int* blocklens=(int*)xbt_malloc(blocklength*count);
2342 for (i=0; i<count;i++)blocklens[i]=blocklength;
2343 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2344 xbt_free(blocklens);
2351 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2356 retval = MPI_ERR_COUNT;
2358 retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2364 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2365 return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2369 int PMPI_Error_class(int errorcode, int* errorclass) {
2370 // assume smpi uses only standard mpi error codes
2371 *errorclass=errorcode;
2376 int PMPI_Initialized(int* flag) {
2377 *flag=(smpi_process_data()!=NULL);
2381 /* The following calls are not yet implemented and will fail at runtime. */
2382 /* Once implemented, please move them above this notice. */
2384 #define NOT_YET_IMPLEMENTED {\
2385 XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__);\
2386 return MPI_SUCCESS;\
2390 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2394 int PMPI_Type_set_name(MPI_Datatype datatype, char * name)
2399 int PMPI_Type_get_name(MPI_Datatype datatype, char * name, int* len)
2404 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2408 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2412 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
2416 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2420 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2424 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2428 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2432 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2436 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2440 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2444 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2448 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2452 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2456 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2460 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2464 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2468 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2472 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2476 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2480 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2484 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2488 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2492 int PMPI_Cancel(MPI_Request* request) {
2496 int PMPI_Buffer_attach(void* buffer, int size) {
2500 int PMPI_Buffer_detach(void* buffer, int* size) {
2504 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2508 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2513 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2518 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2523 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2528 int PMPI_Comm_free_keyval(int* keyval) {
2532 int PMPI_Pcontrol(const int level )
2537 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2541 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2546 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2551 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2556 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2561 int PMPI_Type_free_keyval(int* keyval) {
2565 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2569 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2573 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2577 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2581 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2585 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2589 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2593 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2597 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2601 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2605 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2609 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2613 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2617 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2621 int PMPI_Keyval_free(int* keyval) {
2625 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2629 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2633 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2637 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2641 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2645 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2649 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2653 int PMPI_Win_fence( int assert, MPI_Win win){
2657 int PMPI_Win_free( MPI_Win* win){
2661 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2665 int PMPI_Info_create( MPI_Info *info){
2669 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2673 int PMPI_Info_free( MPI_Info *info){
2677 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2678 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2682 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2683 int *num_addresses, int *num_datatypes, int *combiner){
2687 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2688 int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2689 MPI_Datatype* array_of_datatypes){
2693 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2694 int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2695 int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2699 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2703 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){
2707 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
2711 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
2712 void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
2717 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
2721 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
2725 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
2729 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
2733 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
2737 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
2741 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
2745 int PMPI_Add_error_class( int *errorclass){
2749 int PMPI_Add_error_code( int errorclass, int *errorcode){
2753 int PMPI_Add_error_string( int errorcode, char *string){
2757 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
2761 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
2765 int PMPI_Info_delete(MPI_Info info, char *key){
2769 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
2773 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
2777 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
2781 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2785 int MPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2789 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){
2793 int PMPI_Grequest_complete( MPI_Request request){
2797 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
2801 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
2805 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2809 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
2813 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
2817 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
2821 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
2825 int PMPI_Open_port( MPI_Info info, char *port_name){
2829 int PMPI_Close_port(char *port_name){
2833 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2837 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){
2841 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
2842 int* array_of_maxprocs, MPI_Info* array_of_info, int root,
2843 MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
2847 int PMPI_Comm_get_parent( MPI_Comm *parent){