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 smpi_group_destroy(*group);
307 *group = MPI_GROUP_NULL;
308 retval = MPI_SUCCESS;
314 int PMPI_Group_size(MPI_Group group, int *size)
319 if (group == MPI_GROUP_NULL) {
320 retval = MPI_ERR_GROUP;
321 } else if (size == NULL) {
322 retval = MPI_ERR_ARG;
324 *size = smpi_group_size(group);
325 retval = MPI_SUCCESS;
331 int PMPI_Group_rank(MPI_Group group, int *rank)
336 if (group == MPI_GROUP_NULL) {
337 retval = MPI_ERR_GROUP;
338 } else if (rank == NULL) {
339 retval = MPI_ERR_ARG;
341 *rank = smpi_group_rank(group, smpi_process_index());
342 retval = MPI_SUCCESS;
348 int PMPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1,
349 MPI_Group group2, int *ranks2)
351 int retval, i, index;
353 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
354 retval = MPI_ERR_GROUP;
356 for (i = 0; i < n; i++) {
357 if(ranks1[i]==MPI_PROC_NULL){
358 ranks2[i]=MPI_PROC_NULL;
360 index = smpi_group_index(group1, ranks1[i]);
361 ranks2[i] = smpi_group_rank(group2, index);
364 retval = MPI_SUCCESS;
370 int PMPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
375 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
376 retval = MPI_ERR_GROUP;
377 } else if (result == NULL) {
378 retval = MPI_ERR_ARG;
380 *result = smpi_group_compare(group1, group2);
381 retval = MPI_SUCCESS;
387 int PMPI_Group_union(MPI_Group group1, MPI_Group group2,
388 MPI_Group * newgroup)
390 int retval, i, proc1, proc2, size, size2;
393 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
394 retval = MPI_ERR_GROUP;
395 } else if (newgroup == NULL) {
396 retval = MPI_ERR_ARG;
398 size = smpi_group_size(group1);
399 size2 = smpi_group_size(group2);
400 for (i = 0; i < size2; i++) {
401 proc2 = smpi_group_index(group2, i);
402 proc1 = smpi_group_rank(group1, proc2);
403 if (proc1 == MPI_UNDEFINED) {
408 *newgroup = MPI_GROUP_EMPTY;
410 *newgroup = smpi_group_new(size);
411 size2 = smpi_group_size(group1);
412 for (i = 0; i < size2; i++) {
413 proc1 = smpi_group_index(group1, i);
414 smpi_group_set_mapping(*newgroup, proc1, i);
416 for (i = size2; i < size; i++) {
417 proc2 = smpi_group_index(group2, i - size2);
418 smpi_group_set_mapping(*newgroup, proc2, i);
421 retval = MPI_SUCCESS;
427 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
428 MPI_Group * newgroup)
430 int retval, i, proc1, proc2, size;
433 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
434 retval = MPI_ERR_GROUP;
435 } else if (newgroup == NULL) {
436 retval = MPI_ERR_ARG;
438 size = smpi_group_size(group2);
439 for (i = 0; i < size; i++) {
440 proc2 = smpi_group_index(group2, i);
441 proc1 = smpi_group_rank(group1, proc2);
442 if (proc1 == MPI_UNDEFINED) {
447 *newgroup = MPI_GROUP_EMPTY;
449 *newgroup = smpi_group_new(size);
451 for (i = 0; i < smpi_group_size(group2); i++) {
452 proc2 = smpi_group_index(group2, i);
453 proc1 = smpi_group_rank(group1, proc2);
454 if (proc1 != MPI_UNDEFINED) {
455 smpi_group_set_mapping(*newgroup, proc2, j);
460 retval = MPI_SUCCESS;
466 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
468 int retval, i, proc1, proc2, size, size2;
471 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
472 retval = MPI_ERR_GROUP;
473 } else if (newgroup == NULL) {
474 retval = MPI_ERR_ARG;
476 size = size2 = smpi_group_size(group1);
477 for (i = 0; i < size2; i++) {
478 proc1 = smpi_group_index(group1, i);
479 proc2 = smpi_group_rank(group2, proc1);
480 if (proc2 != MPI_UNDEFINED) {
485 *newgroup = MPI_GROUP_EMPTY;
487 *newgroup = smpi_group_new(size);
488 for (i = 0; i < size2; i++) {
489 proc1 = smpi_group_index(group1, i);
490 proc2 = smpi_group_rank(group2, proc1);
491 if (proc2 == MPI_UNDEFINED) {
492 smpi_group_set_mapping(*newgroup, proc1, i);
496 retval = MPI_SUCCESS;
502 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
504 int retval, i, index;
507 if (group == MPI_GROUP_NULL) {
508 retval = MPI_ERR_GROUP;
509 } else if (newgroup == NULL) {
510 retval = MPI_ERR_ARG;
513 *newgroup = MPI_GROUP_EMPTY;
514 } else if (n == smpi_group_size(group)) {
516 if(group!= smpi_comm_group(MPI_COMM_WORLD)
517 && group != MPI_GROUP_NULL
518 && group != smpi_comm_group(MPI_COMM_SELF)
519 && group != MPI_GROUP_EMPTY)
520 smpi_group_use(group);
522 *newgroup = smpi_group_new(n);
523 for (i = 0; i < n; i++) {
524 index = smpi_group_index(group, ranks[i]);
525 smpi_group_set_mapping(*newgroup, index, i);
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 if(group!= smpi_comm_group(MPI_COMM_WORLD)
547 && group != MPI_GROUP_NULL
548 && group != smpi_comm_group(MPI_COMM_SELF)
549 && group != MPI_GROUP_EMPTY)
550 smpi_group_use(group);
551 } else if (n == smpi_group_size(group)) {
552 *newgroup = MPI_GROUP_EMPTY;
554 oldsize=smpi_group_size(group);
555 newsize = oldsize - n;
556 *newgroup = smpi_group_new(newsize);
558 int* to_exclude=xbt_new(int, smpi_group_size(group));
559 for(i=0; i<oldsize; i++)
562 to_exclude[ranks[i]]=1;
565 for(i=0; i<oldsize; i++){
566 if(to_exclude[i]==0){
567 index = smpi_group_index(group, i);
568 smpi_group_set_mapping(*newgroup, index, j);
573 xbt_free(to_exclude);
575 retval = MPI_SUCCESS;
581 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
582 MPI_Group * newgroup)
584 int retval, i, j, rank, size, index;
587 if (group == MPI_GROUP_NULL) {
588 retval = MPI_ERR_GROUP;
589 } else if (newgroup == NULL) {
590 retval = MPI_ERR_ARG;
593 *newgroup = MPI_GROUP_EMPTY;
596 for (i = 0; i < n; i++) {
597 for (rank = ranges[i][0]; /* First */
598 rank >= 0; /* Last */
602 rank += ranges[i][2]; /* Stride */
603 if (ranges[i][0]<ranges[i][1]){
604 if(rank > ranges[i][1])
607 if(rank < ranges[i][1])
613 *newgroup = smpi_group_new(size);
615 for (i = 0; i < n; i++) {
616 for (rank = ranges[i][0]; /* First */
617 rank >= 0; /* Last */
619 index = smpi_group_index(group, rank);
620 smpi_group_set_mapping(*newgroup, index, j);
622 rank += ranges[i][2]; /* Stride */
623 if (ranges[i][0]<ranges[i][1]){
624 if(rank > ranges[i][1])
627 if(rank < ranges[i][1])
633 retval = MPI_SUCCESS;
639 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
640 MPI_Group * newgroup)
642 int retval, i, rank, newrank,oldrank, size, index, add;
645 if (group == MPI_GROUP_NULL) {
646 retval = MPI_ERR_GROUP;
647 } else if (newgroup == NULL) {
648 retval = MPI_ERR_ARG;
652 if(group!= smpi_comm_group(MPI_COMM_WORLD)
653 && group != MPI_GROUP_NULL
654 && group != smpi_comm_group(MPI_COMM_SELF)
655 && group != MPI_GROUP_EMPTY)
656 smpi_group_use(group);
658 size = smpi_group_size(group);
659 for (i = 0; i < n; i++) {
660 for (rank = ranges[i][0]; /* First */
661 rank >= 0; /* Last */
665 rank += ranges[i][2]; /* Stride */
666 if (ranges[i][0]<ranges[i][1]){
667 if(rank > ranges[i][1])
670 if(rank < ranges[i][1])
676 *newgroup = MPI_GROUP_EMPTY;
678 *newgroup = smpi_group_new(size);
681 while (newrank < size) {
683 for (i = 0; i < n; i++) {
684 for (rank = ranges[i][0];rank >= 0;){
690 rank += ranges[i][2]; /* Stride */
692 if (ranges[i][0]<ranges[i][1]){
693 if(rank > ranges[i][1])
696 if(rank < ranges[i][1])
702 index = smpi_group_index(group, oldrank);
703 smpi_group_set_mapping(*newgroup, index, newrank);
711 retval = MPI_SUCCESS;
717 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
722 if (comm == MPI_COMM_NULL) {
723 retval = MPI_ERR_COMM;
724 } else if (rank == NULL) {
725 retval = MPI_ERR_ARG;
727 *rank = smpi_comm_rank(comm);
728 retval = MPI_SUCCESS;
734 int PMPI_Comm_size(MPI_Comm comm, int *size)
739 if (comm == MPI_COMM_NULL) {
740 retval = MPI_ERR_COMM;
741 } else if (size == NULL) {
742 retval = MPI_ERR_ARG;
744 *size = smpi_comm_size(comm);
745 retval = MPI_SUCCESS;
751 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
756 if (comm == MPI_COMM_NULL) {
757 retval = MPI_ERR_COMM;
758 } else if (name == NULL || len == NULL) {
759 retval = MPI_ERR_ARG;
761 smpi_comm_get_name(comm, name, len);
762 retval = MPI_SUCCESS;
768 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
773 if (comm == MPI_COMM_NULL) {
774 retval = MPI_ERR_COMM;
775 } else if (group == NULL) {
776 retval = MPI_ERR_ARG;
778 *group = smpi_comm_group(comm);
779 if(*group!= smpi_comm_group(MPI_COMM_WORLD)
780 && *group != MPI_GROUP_NULL
781 && *group != smpi_comm_group(MPI_COMM_SELF)
782 && *group != MPI_GROUP_EMPTY)
783 smpi_group_use(*group);
784 retval = MPI_SUCCESS;
790 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
795 if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
796 retval = MPI_ERR_COMM;
797 } else if (result == NULL) {
798 retval = MPI_ERR_ARG;
800 if (comm1 == comm2) { /* Same communicators means same groups */
804 smpi_group_compare(smpi_comm_group(comm1),
805 smpi_comm_group(comm2));
806 if (*result == MPI_IDENT) {
807 *result = MPI_CONGRUENT;
810 retval = MPI_SUCCESS;
816 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
821 if (comm == MPI_COMM_NULL) {
822 retval = MPI_ERR_COMM;
823 } else if (newcomm == NULL) {
824 retval = MPI_ERR_ARG;
826 *newcomm = smpi_comm_new(smpi_comm_group(comm));
827 retval = MPI_SUCCESS;
833 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
838 if (comm == MPI_COMM_NULL) {
839 retval = MPI_ERR_COMM;
840 } else if (group == MPI_GROUP_NULL) {
841 retval = MPI_ERR_GROUP;
842 } else if (newcomm == NULL) {
843 retval = MPI_ERR_ARG;
844 } else if(smpi_group_rank(group,smpi_process_index())==MPI_UNDEFINED){
845 *newcomm= MPI_COMM_NULL;
846 retval = MPI_SUCCESS;
849 *newcomm = smpi_comm_new(group);
850 retval = MPI_SUCCESS;
856 int PMPI_Comm_free(MPI_Comm * comm)
862 retval = MPI_ERR_ARG;
863 } else if (*comm == MPI_COMM_NULL) {
864 retval = MPI_ERR_COMM;
866 smpi_comm_destroy(*comm);
867 *comm = MPI_COMM_NULL;
868 retval = MPI_SUCCESS;
874 int PMPI_Comm_disconnect(MPI_Comm * comm)
876 /* TODO: wait until all communication in comm are done */
881 retval = MPI_ERR_ARG;
882 } else if (*comm == MPI_COMM_NULL) {
883 retval = MPI_ERR_COMM;
885 smpi_comm_destroy(*comm);
886 *comm = MPI_COMM_NULL;
887 retval = MPI_SUCCESS;
893 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
898 if (comm_out == NULL) {
899 retval = MPI_ERR_ARG;
900 } else if (comm == MPI_COMM_NULL) {
901 retval = MPI_ERR_COMM;
903 *comm_out = smpi_comm_split(comm, color, key);
904 retval = MPI_SUCCESS;
910 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
911 int tag, MPI_Comm comm, MPI_Request * request)
916 if (request == NULL) {
917 retval = MPI_ERR_ARG;
918 } else if (comm == MPI_COMM_NULL) {
919 retval = MPI_ERR_COMM;
920 } else if (dst == MPI_PROC_NULL) {
921 retval = MPI_SUCCESS;
923 *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
924 retval = MPI_SUCCESS;
930 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
931 int tag, MPI_Comm comm, MPI_Request * request)
936 if (request == NULL) {
937 retval = MPI_ERR_ARG;
938 } else if (comm == MPI_COMM_NULL) {
939 retval = MPI_ERR_COMM;
940 } else if (src == MPI_PROC_NULL) {
941 retval = MPI_SUCCESS;
943 *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
944 retval = MPI_SUCCESS;
950 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
954 if (request == NULL) {
955 retval = MPI_ERR_ARG;
956 } else if (comm == MPI_COMM_NULL) {
957 retval = MPI_ERR_COMM;
958 } else if (dst == MPI_PROC_NULL) {
959 retval = MPI_SUCCESS;
961 *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
962 retval = MPI_SUCCESS;
968 int PMPI_Start(MPI_Request * request)
973 if (request == NULL || *request == MPI_REQUEST_NULL) {
974 retval = MPI_ERR_ARG;
976 smpi_mpi_start(*request);
977 retval = MPI_SUCCESS;
983 int PMPI_Startall(int count, MPI_Request * requests)
988 if (requests == NULL) {
989 retval = MPI_ERR_ARG;
991 smpi_mpi_startall(count, requests);
992 retval = MPI_SUCCESS;
998 int PMPI_Request_free(MPI_Request * request)
1003 if (*request == MPI_REQUEST_NULL) {
1004 retval = MPI_ERR_ARG;
1006 if((*request)->flags & PERSISTENT)(*request)->refcount--;
1007 smpi_mpi_request_free(request);
1008 retval = MPI_SUCCESS;
1014 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
1015 int tag, MPI_Comm comm, MPI_Request * request)
1021 if (request == NULL) {
1022 retval = MPI_ERR_ARG;
1023 } else if (comm == MPI_COMM_NULL) {
1024 retval = MPI_ERR_COMM;
1025 } else if (src == MPI_PROC_NULL) {
1026 *request = MPI_REQUEST_NULL;
1027 retval = MPI_SUCCESS;
1028 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1029 retval = MPI_ERR_COMM;
1030 } else if (count < 0) {
1031 retval = MPI_ERR_COUNT;
1032 } else if (buf==NULL && count > 0) {
1033 retval = MPI_ERR_COUNT;
1034 } else if (datatype == MPI_DATATYPE_NULL){
1035 retval = MPI_ERR_TYPE;
1036 } else if(tag<0 && tag != MPI_ANY_TAG){
1037 retval = MPI_ERR_TAG;
1041 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1042 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1043 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
1046 *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
1047 retval = MPI_SUCCESS;
1050 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1051 (*request)->recv = 1;
1060 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
1061 int tag, MPI_Comm comm, MPI_Request * request)
1066 if (request == NULL) {
1067 retval = MPI_ERR_ARG;
1068 } else if (comm == MPI_COMM_NULL) {
1069 retval = MPI_ERR_COMM;
1070 } else if (dst == MPI_PROC_NULL) {
1071 *request = MPI_REQUEST_NULL;
1072 retval = MPI_SUCCESS;
1073 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1074 retval = MPI_ERR_RANK;
1075 } else if (count < 0) {
1076 retval = MPI_ERR_COUNT;
1077 } else if (buf==NULL && count > 0) {
1078 retval = MPI_ERR_COUNT;
1079 } else if (datatype == MPI_DATATYPE_NULL){
1080 retval = MPI_ERR_TYPE;
1081 } else if(tag<0 && tag != MPI_ANY_TAG){
1082 retval = MPI_ERR_TAG;
1086 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1087 TRACE_smpi_computing_out(rank);
1088 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1089 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1090 TRACE_smpi_send(rank, rank, dst_traced);
1093 *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1094 retval = MPI_SUCCESS;
1097 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1098 (*request)->send = 1;
1099 TRACE_smpi_computing_in(rank);
1107 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
1111 if (request == NULL) {
1112 retval = MPI_ERR_ARG;
1113 } else if (comm == MPI_COMM_NULL) {
1114 retval = MPI_ERR_COMM;
1115 } else if (dst == MPI_PROC_NULL) {
1116 *request = MPI_REQUEST_NULL;
1117 retval = MPI_SUCCESS;
1118 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1119 retval = MPI_ERR_RANK;
1120 } else if (count < 0) {
1121 retval = MPI_ERR_COUNT;
1122 } else if (buf==NULL && count > 0) {
1123 retval = MPI_ERR_COUNT;
1124 } else if (datatype == MPI_DATATYPE_NULL){
1125 retval = MPI_ERR_TYPE;
1126 } else if(tag<0 && tag != MPI_ANY_TAG){
1127 retval = MPI_ERR_TAG;
1131 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1132 TRACE_smpi_computing_out(rank);
1133 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1134 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1135 TRACE_smpi_send(rank, rank, dst_traced);
1138 *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1139 retval = MPI_SUCCESS;
1142 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1143 (*request)->send = 1;
1144 TRACE_smpi_computing_in(rank);
1152 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1153 MPI_Comm comm, MPI_Status * status)
1158 if (comm == MPI_COMM_NULL) {
1159 retval = MPI_ERR_COMM;
1160 } else if (src == MPI_PROC_NULL) {
1161 smpi_empty_status(status);
1162 status->MPI_SOURCE = MPI_PROC_NULL;
1163 retval = MPI_SUCCESS;
1164 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1165 retval = MPI_ERR_RANK;
1166 } else if (count < 0) {
1167 retval = MPI_ERR_COUNT;
1168 } else if (buf==NULL && count > 0) {
1169 retval = MPI_ERR_COUNT;
1170 } else if (datatype == MPI_DATATYPE_NULL){
1171 retval = MPI_ERR_TYPE;
1172 } else if(tag<0 && tag != MPI_ANY_TAG){
1173 retval = MPI_ERR_TAG;
1176 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1177 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1178 TRACE_smpi_computing_out(rank);
1180 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
1183 smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1184 retval = MPI_SUCCESS;
1187 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1188 if(status!=MPI_STATUS_IGNORE)src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1189 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1190 TRACE_smpi_recv(rank, src_traced, rank);
1191 TRACE_smpi_computing_in(rank);
1199 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1206 if (comm == MPI_COMM_NULL) {
1207 retval = MPI_ERR_COMM;
1208 } else if (dst == MPI_PROC_NULL) {
1209 retval = MPI_SUCCESS;
1210 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1211 retval = MPI_ERR_RANK;
1212 } else if (count < 0) {
1213 retval = MPI_ERR_COUNT;
1214 } else if (buf==NULL && count > 0) {
1215 retval = MPI_ERR_COUNT;
1216 } else if (datatype == MPI_DATATYPE_NULL){
1217 retval = MPI_ERR_TYPE;
1218 } else if(tag<0 && tag != MPI_ANY_TAG){
1219 retval = MPI_ERR_TAG;
1223 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1224 TRACE_smpi_computing_out(rank);
1225 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1226 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1227 TRACE_smpi_send(rank, rank, dst_traced);
1230 smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1231 retval = MPI_SUCCESS;
1234 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1235 TRACE_smpi_computing_in(rank);
1245 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
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 TRACE_smpi_computing_out(rank);
1269 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1270 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1271 TRACE_smpi_send(rank, rank, dst_traced);
1274 smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1275 retval = MPI_SUCCESS;
1278 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1279 TRACE_smpi_computing_in(rank);
1287 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1288 int dst, int sendtag, void *recvbuf, int recvcount,
1289 MPI_Datatype recvtype, int src, int recvtag,
1290 MPI_Comm comm, MPI_Status * status)
1296 if (comm == MPI_COMM_NULL) {
1297 retval = MPI_ERR_COMM;
1298 } else if (sendtype == MPI_DATATYPE_NULL
1299 || recvtype == MPI_DATATYPE_NULL) {
1300 retval = MPI_ERR_TYPE;
1301 } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1302 smpi_empty_status(status);
1303 status->MPI_SOURCE = MPI_PROC_NULL;
1304 retval = MPI_SUCCESS;
1305 }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1306 (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1307 retval = MPI_ERR_RANK;
1308 } else if (sendcount < 0 || recvcount<0) {
1309 retval = MPI_ERR_COUNT;
1310 } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1311 retval = MPI_ERR_COUNT;
1312 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1313 retval = MPI_ERR_TAG;
1317 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1318 TRACE_smpi_computing_out(rank);
1319 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1320 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1321 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1322 TRACE_smpi_send(rank, rank, dst_traced);
1326 smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1327 recvcount, recvtype, src, recvtag, comm, status);
1328 retval = MPI_SUCCESS;
1331 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1332 TRACE_smpi_recv(rank, src_traced, rank);
1333 TRACE_smpi_computing_in(rank);
1342 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1343 int dst, int sendtag, int src, int recvtag,
1344 MPI_Comm comm, MPI_Status * status)
1346 //TODO: suboptimal implementation
1349 if (datatype == MPI_DATATYPE_NULL) {
1350 retval = MPI_ERR_TYPE;
1351 } else if (count < 0) {
1352 retval = MPI_ERR_COUNT;
1354 int size = smpi_datatype_get_extent(datatype) * count;
1355 recvbuf = xbt_new(char, size);
1357 MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1358 datatype, src, recvtag, comm, status);
1359 if(retval==MPI_SUCCESS){
1360 smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1368 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1373 if (request == MPI_REQUEST_NULL || flag == NULL) {
1374 retval = MPI_ERR_ARG;
1375 } else if (*request == MPI_REQUEST_NULL) {
1377 retval = MPI_ERR_REQUEST;
1379 *flag = smpi_mpi_test(request, status);
1380 retval = MPI_SUCCESS;
1386 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1387 MPI_Status * status)
1392 if (index == NULL || flag == NULL) {
1393 retval = MPI_ERR_ARG;
1395 *flag = smpi_mpi_testany(count, requests, index, status);
1396 retval = MPI_SUCCESS;
1402 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1408 retval = MPI_ERR_ARG;
1410 *flag = smpi_mpi_testall(count, requests, statuses);
1411 retval = MPI_SUCCESS;
1417 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1421 if (status == NULL) {
1422 retval = MPI_ERR_ARG;
1423 } else if (comm == MPI_COMM_NULL) {
1424 retval = MPI_ERR_COMM;
1425 } else if (source == MPI_PROC_NULL) {
1426 smpi_empty_status(status);
1427 status->MPI_SOURCE = MPI_PROC_NULL;
1428 retval = MPI_SUCCESS;
1430 smpi_mpi_probe(source, tag, comm, status);
1431 retval = MPI_SUCCESS;
1438 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1443 retval = MPI_ERR_ARG;
1444 } else if (status == NULL) {
1445 retval = MPI_ERR_ARG;
1446 } else if (comm == MPI_COMM_NULL) {
1447 retval = MPI_ERR_COMM;
1448 } else if (source == MPI_PROC_NULL) {
1450 smpi_empty_status(status);
1451 status->MPI_SOURCE = MPI_PROC_NULL;
1452 retval = MPI_SUCCESS;
1454 smpi_mpi_iprobe(source, tag, comm, flag, status);
1455 retval = MPI_SUCCESS;
1461 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1467 if (request == NULL) {
1468 retval = MPI_ERR_ARG;
1469 } else if (*request == MPI_REQUEST_NULL) {
1470 retval = MPI_ERR_REQUEST;
1474 int rank = request && (*request)->comm != MPI_COMM_NULL
1475 ? smpi_process_index()
1477 TRACE_smpi_computing_out(rank);
1479 int src_traced = (*request)->src;
1480 int dst_traced = (*request)->dst;
1481 int is_wait_for_receive = (*request)->recv;
1482 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1485 smpi_mpi_wait(request, status);
1486 retval = MPI_SUCCESS;
1489 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1490 if (is_wait_for_receive) {
1491 TRACE_smpi_recv(rank, src_traced, dst_traced);
1493 TRACE_smpi_computing_in(rank);
1502 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1508 //save requests information for tracing
1510 int *srcs = xbt_new(int, count);
1511 int *dsts = xbt_new(int, count);
1512 int *recvs = xbt_new(int, count);
1513 for (i = 0; i < count; i++) {
1514 MPI_Request req = requests[i]; //already received requests are no longer valid
1518 recvs[i] = req->recv;
1521 int rank_traced = smpi_process_index();
1522 TRACE_smpi_computing_out(rank_traced);
1524 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1527 if (index == NULL) {
1528 retval = MPI_ERR_ARG;
1530 *index = smpi_mpi_waitany(count, requests, status);
1531 retval = MPI_SUCCESS;
1534 if(*index!=MPI_UNDEFINED){
1535 int src_traced = srcs[*index];
1536 int dst_traced = dsts[*index];
1537 int is_wait_for_receive = recvs[*index];
1538 if (is_wait_for_receive) {
1539 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1541 TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1546 TRACE_smpi_computing_in(rank_traced);
1552 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1557 //save information from requests
1559 int *srcs = xbt_new(int, count);
1560 int *dsts = xbt_new(int, count);
1561 int *recvs = xbt_new(int, count);
1562 int valid_count = 0;
1563 for (i = 0; i < count; i++) {
1564 MPI_Request req = requests[i];
1565 if(req!=MPI_REQUEST_NULL){
1566 srcs[valid_count] = req->src;
1567 dsts[valid_count] = req->dst;
1568 recvs[valid_count] = req->recv;
1572 int rank_traced = smpi_process_index();
1573 TRACE_smpi_computing_out(rank_traced);
1575 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1577 int retval = smpi_mpi_waitall(count, requests, status);
1579 for (i = 0; i < valid_count; i++) {
1580 int src_traced = srcs[i];
1581 int dst_traced = dsts[i];
1582 int is_wait_for_receive = recvs[i];
1583 if (is_wait_for_receive) {
1584 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1587 TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1591 TRACE_smpi_computing_in(rank_traced);
1597 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1598 int *indices, MPI_Status status[])
1603 if (outcount == NULL) {
1604 retval = MPI_ERR_ARG;
1606 *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1607 retval = MPI_SUCCESS;
1613 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1614 int* indices, MPI_Status status[])
1619 if (outcount == NULL) {
1620 retval = MPI_ERR_ARG;
1622 *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1623 retval = MPI_SUCCESS;
1630 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1636 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1637 TRACE_smpi_computing_out(rank);
1638 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1639 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1641 if (comm == MPI_COMM_NULL) {
1642 retval = MPI_ERR_COMM;
1644 mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1645 retval = MPI_SUCCESS;
1648 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1649 TRACE_smpi_computing_in(rank);
1655 int PMPI_Barrier(MPI_Comm comm)
1661 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1662 TRACE_smpi_computing_out(rank);
1663 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1665 if (comm == MPI_COMM_NULL) {
1666 retval = MPI_ERR_COMM;
1668 mpi_coll_barrier_fun(comm);
1669 retval = MPI_SUCCESS;
1672 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1673 TRACE_smpi_computing_in(rank);
1679 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1680 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1681 int root, MPI_Comm comm)
1687 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1688 TRACE_smpi_computing_out(rank);
1689 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1690 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1692 if (comm == MPI_COMM_NULL) {
1693 retval = MPI_ERR_COMM;
1694 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1695 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1696 retval = MPI_ERR_TYPE;
1697 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1698 ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1699 retval = MPI_ERR_COUNT;
1702 char* sendtmpbuf = (char*) sendbuf;
1703 int sendtmpcount = sendcount;
1704 MPI_Datatype sendtmptype = sendtype;
1705 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1707 sendtmptype=recvtype;
1710 mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount,
1711 recvtype, root, comm);
1714 retval = MPI_SUCCESS;
1717 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1718 TRACE_smpi_computing_in(rank);
1724 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1725 void *recvbuf, int *recvcounts, int *displs,
1726 MPI_Datatype recvtype, int root, MPI_Comm comm)
1732 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1733 TRACE_smpi_computing_out(rank);
1734 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1735 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1737 if (comm == MPI_COMM_NULL) {
1738 retval = MPI_ERR_COMM;
1739 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1740 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1741 retval = MPI_ERR_TYPE;
1742 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1743 retval = MPI_ERR_COUNT;
1744 } else if (recvcounts == NULL || displs == NULL) {
1745 retval = MPI_ERR_ARG;
1748 char* sendtmpbuf = (char*) sendbuf;
1749 int sendtmpcount = sendcount;
1750 MPI_Datatype sendtmptype = sendtype;
1751 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1753 sendtmptype=recvtype;
1756 smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts,
1757 displs, recvtype, root, comm);
1758 retval = MPI_SUCCESS;
1761 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1762 TRACE_smpi_computing_in(rank);
1768 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1769 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1776 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1777 TRACE_smpi_computing_out(rank);
1778 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1780 if (comm == MPI_COMM_NULL) {
1781 retval = MPI_ERR_COMM;
1782 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1783 (recvtype == MPI_DATATYPE_NULL)){
1784 retval = MPI_ERR_TYPE;
1785 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1787 retval = MPI_ERR_COUNT;
1790 if(sendbuf == MPI_IN_PLACE) {
1791 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1792 sendcount=recvcount;
1796 mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1798 retval = MPI_SUCCESS;
1801 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1807 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1808 void *recvbuf, int *recvcounts, int *displs,
1809 MPI_Datatype recvtype, MPI_Comm comm)
1815 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1816 TRACE_smpi_computing_out(rank);
1817 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1819 if (comm == MPI_COMM_NULL) {
1820 retval = MPI_ERR_COMM;
1821 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1822 (recvtype == MPI_DATATYPE_NULL)){
1823 retval = MPI_ERR_TYPE;
1824 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1825 retval = MPI_ERR_COUNT;
1826 } else if (recvcounts == NULL || displs == NULL) {
1827 retval = MPI_ERR_ARG;
1830 if(sendbuf == MPI_IN_PLACE) {
1831 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1832 sendcount=recvcounts[smpi_comm_rank(comm)];
1836 mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1837 displs, recvtype, comm);
1838 retval = MPI_SUCCESS;
1841 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1842 TRACE_smpi_computing_in(rank);
1848 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1849 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1850 int root, MPI_Comm comm)
1856 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1857 TRACE_smpi_computing_out(rank);
1858 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1860 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1862 if (comm == MPI_COMM_NULL) {
1863 retval = MPI_ERR_COMM;
1864 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1865 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1866 retval = MPI_ERR_TYPE;
1869 if(recvbuf==MPI_IN_PLACE){
1872 mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1873 recvtype, root, comm);
1874 retval = MPI_SUCCESS;
1877 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1878 TRACE_smpi_computing_in(rank);
1884 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1885 MPI_Datatype sendtype, void *recvbuf, int recvcount,
1886 MPI_Datatype recvtype, int root, MPI_Comm comm)
1892 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1893 TRACE_smpi_computing_out(rank);
1894 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1895 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1897 if (comm == MPI_COMM_NULL) {
1898 retval = MPI_ERR_COMM;
1899 } else if (sendcounts == NULL || displs == NULL) {
1900 retval = MPI_ERR_ARG;
1901 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1902 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1903 retval = MPI_ERR_TYPE;
1906 if(recvbuf==MPI_IN_PLACE){
1910 smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1911 recvcount, recvtype, root, comm);
1912 retval = MPI_SUCCESS;
1915 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1916 TRACE_smpi_computing_in(rank);
1922 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1923 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1929 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1930 TRACE_smpi_computing_out(rank);
1931 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1932 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1934 if (comm == MPI_COMM_NULL) {
1935 retval = MPI_ERR_COMM;
1936 } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1937 retval = MPI_ERR_ARG;
1940 char* sendtmpbuf = (char*) sendbuf;
1941 if( sendbuf == MPI_IN_PLACE ) {
1942 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
1943 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
1946 mpi_coll_reduce_fun(sendtmpbuf, recvbuf, count, datatype, op, root, comm);
1948 if( sendbuf == MPI_IN_PLACE ) {
1949 xbt_free(sendtmpbuf);
1952 retval = MPI_SUCCESS;
1955 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1956 TRACE_smpi_computing_in(rank);
1962 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
1963 MPI_Datatype datatype, MPI_Op op){
1967 if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1968 retval = MPI_ERR_ARG;
1970 smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
1977 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1978 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1984 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1985 TRACE_smpi_computing_out(rank);
1986 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1988 if (comm == MPI_COMM_NULL) {
1989 retval = MPI_ERR_COMM;
1990 } else if (datatype == MPI_DATATYPE_NULL) {
1991 retval = MPI_ERR_TYPE;
1992 } else if (op == MPI_OP_NULL) {
1993 retval = MPI_ERR_OP;
1996 char* sendtmpbuf = (char*) sendbuf;
1997 if( sendbuf == MPI_IN_PLACE ) {
1998 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
1999 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
2002 mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
2004 if( sendbuf == MPI_IN_PLACE ) {
2005 xbt_free(sendtmpbuf);
2008 retval = MPI_SUCCESS;
2012 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2013 TRACE_smpi_computing_in(rank);
2019 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
2020 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2026 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2027 TRACE_smpi_computing_out(rank);
2028 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2030 if (comm == MPI_COMM_NULL) {
2031 retval = MPI_ERR_COMM;
2032 } else if (datatype == MPI_DATATYPE_NULL) {
2033 retval = MPI_ERR_TYPE;
2034 } else if (op == MPI_OP_NULL) {
2035 retval = MPI_ERR_OP;
2037 smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2038 retval = MPI_SUCCESS;
2041 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2042 TRACE_smpi_computing_in(rank);
2048 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2049 MPI_Op op, MPI_Comm comm){
2054 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2055 TRACE_smpi_computing_out(rank);
2056 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2058 if (comm == MPI_COMM_NULL) {
2059 retval = MPI_ERR_COMM;
2060 } else if (datatype == MPI_DATATYPE_NULL) {
2061 retval = MPI_ERR_TYPE;
2062 } else if (op == MPI_OP_NULL) {
2063 retval = MPI_ERR_OP;
2065 smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2066 retval = MPI_SUCCESS;
2069 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2070 TRACE_smpi_computing_in(rank);
2076 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
2077 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2082 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2083 TRACE_smpi_computing_out(rank);
2084 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2086 if (comm == MPI_COMM_NULL) {
2087 retval = MPI_ERR_COMM;
2088 } else if (datatype == MPI_DATATYPE_NULL) {
2089 retval = MPI_ERR_TYPE;
2090 } else if (op == MPI_OP_NULL) {
2091 retval = MPI_ERR_OP;
2092 } else if (recvcounts == NULL) {
2093 retval = MPI_ERR_ARG;
2095 void* sendtmpbuf=sendbuf;
2096 if(sendbuf==MPI_IN_PLACE){
2100 mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts,
2101 datatype, op, comm);
2102 retval = MPI_SUCCESS;
2105 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2106 TRACE_smpi_computing_in(rank);
2112 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2113 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2118 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2119 TRACE_smpi_computing_out(rank);
2120 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2122 if (comm == MPI_COMM_NULL) {
2123 retval = MPI_ERR_COMM;
2124 } else if (datatype == MPI_DATATYPE_NULL) {
2125 retval = MPI_ERR_TYPE;
2126 } else if (op == MPI_OP_NULL) {
2127 retval = MPI_ERR_OP;
2128 } else if (recvcount < 0) {
2129 retval = MPI_ERR_ARG;
2131 int count=smpi_comm_size(comm);
2132 int* recvcounts=(int*)xbt_malloc(count);
2133 for (i=0; i<count;i++)recvcounts[i]=recvcount;
2134 mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
2135 datatype, op, comm);
2136 xbt_free(recvcounts);
2137 retval = MPI_SUCCESS;
2140 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2141 TRACE_smpi_computing_in(rank);
2147 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2148 void *recvbuf, int recvcount, MPI_Datatype recvtype,
2155 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2156 TRACE_smpi_computing_out(rank);
2157 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2159 if (comm == MPI_COMM_NULL) {
2160 retval = MPI_ERR_COMM;
2161 } else if (sendtype == MPI_DATATYPE_NULL
2162 || recvtype == MPI_DATATYPE_NULL) {
2163 retval = MPI_ERR_TYPE;
2165 retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2168 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2169 TRACE_smpi_computing_in(rank);
2175 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2176 MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2177 int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2183 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2184 TRACE_smpi_computing_out(rank);
2185 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2187 if (comm == MPI_COMM_NULL) {
2188 retval = MPI_ERR_COMM;
2189 } else if (sendtype == MPI_DATATYPE_NULL
2190 || recvtype == MPI_DATATYPE_NULL) {
2191 retval = MPI_ERR_TYPE;
2192 } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2193 || recvdisps == NULL) {
2194 retval = MPI_ERR_ARG;
2197 mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2198 recvbuf, recvcounts, recvdisps, recvtype,
2202 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2203 TRACE_smpi_computing_in(rank);
2210 int PMPI_Get_processor_name(char *name, int *resultlen)
2212 int retval = MPI_SUCCESS;
2215 strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2216 strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2217 strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2218 MPI_MAX_PROCESSOR_NAME - 1 );
2221 MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2227 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2229 int retval = MPI_SUCCESS;
2233 if (status == NULL || count == NULL) {
2234 retval = MPI_ERR_ARG;
2235 } else if (datatype == MPI_DATATYPE_NULL) {
2236 retval = MPI_ERR_TYPE;
2238 size = smpi_datatype_size(datatype);
2241 } else if (status->count % size != 0) {
2242 retval = MPI_UNDEFINED;
2244 *count = smpi_mpi_get_count(status, datatype);
2251 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2255 if (old_type == MPI_DATATYPE_NULL) {
2256 retval = MPI_ERR_TYPE;
2257 } else if (count<0){
2258 retval = MPI_ERR_COUNT;
2260 retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2266 int PMPI_Type_commit(MPI_Datatype* datatype) {
2270 if (datatype == MPI_DATATYPE_NULL) {
2271 retval = MPI_ERR_TYPE;
2273 smpi_datatype_commit(datatype);
2274 retval = MPI_SUCCESS;
2281 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2285 if (old_type == MPI_DATATYPE_NULL) {
2286 retval = MPI_ERR_TYPE;
2287 } else if (count<0 || blocklen<0){
2288 retval = MPI_ERR_COUNT;
2290 retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2296 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2300 if (old_type == MPI_DATATYPE_NULL) {
2301 retval = MPI_ERR_TYPE;
2302 } else if (count<0 || blocklen<0){
2303 retval = MPI_ERR_COUNT;
2305 retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2311 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2312 return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2315 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2319 if (old_type == MPI_DATATYPE_NULL) {
2320 retval = MPI_ERR_TYPE;
2321 } else if (count<0){
2322 retval = MPI_ERR_COUNT;
2324 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2330 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2334 if (old_type == MPI_DATATYPE_NULL) {
2335 retval = MPI_ERR_TYPE;
2336 } else if (count<0){
2337 retval = MPI_ERR_COUNT;
2339 int* blocklens=(int*)xbt_malloc(blocklength*count);
2340 for (i=0; i<count;i++)blocklens[i]=blocklength;
2341 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2342 xbt_free(blocklens);
2349 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2353 if (old_type == MPI_DATATYPE_NULL) {
2354 retval = MPI_ERR_TYPE;
2355 } else if (count<0){
2356 retval = MPI_ERR_COUNT;
2358 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2364 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2368 if (old_type == MPI_DATATYPE_NULL) {
2369 retval = MPI_ERR_TYPE;
2370 } else if (count<0){
2371 retval = MPI_ERR_COUNT;
2373 int* blocklens=(int*)xbt_malloc(blocklength*count);
2374 for (i=0; i<count;i++)blocklens[i]=blocklength;
2375 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2376 xbt_free(blocklens);
2383 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2388 retval = MPI_ERR_COUNT;
2390 retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2396 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2397 return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2401 int PMPI_Error_class(int errorcode, int* errorclass) {
2402 // assume smpi uses only standard mpi error codes
2403 *errorclass=errorcode;
2408 int PMPI_Initialized(int* flag) {
2409 *flag=(smpi_process_data()!=NULL);
2413 /* The following calls are not yet implemented and will fail at runtime. */
2414 /* Once implemented, please move them above this notice. */
2416 #define NOT_YET_IMPLEMENTED {\
2417 XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__);\
2418 return MPI_SUCCESS;\
2422 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2426 int PMPI_Type_set_name(MPI_Datatype datatype, char * name)
2431 int PMPI_Type_get_name(MPI_Datatype datatype, char * name, int* len)
2436 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2440 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2444 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
2448 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2452 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2456 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2460 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2464 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2468 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2472 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2476 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2480 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2484 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2488 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2492 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2496 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2500 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2504 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2508 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2512 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2516 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2520 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2524 int PMPI_Cancel(MPI_Request* request) {
2528 int PMPI_Buffer_attach(void* buffer, int size) {
2532 int PMPI_Buffer_detach(void* buffer, int* size) {
2536 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2540 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2545 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2550 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2555 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2560 int PMPI_Comm_free_keyval(int* keyval) {
2564 int PMPI_Pcontrol(const int level )
2569 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2573 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2578 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2583 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2588 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2593 int PMPI_Type_free_keyval(int* keyval) {
2597 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2601 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2605 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2609 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2613 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2617 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2621 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2625 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2629 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2633 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2637 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2641 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2645 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2649 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2653 int PMPI_Keyval_free(int* keyval) {
2657 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2661 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2665 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2669 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2673 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2677 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2681 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2685 int PMPI_Win_fence( int assert, MPI_Win win){
2689 int PMPI_Win_free( MPI_Win* win){
2693 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2697 int PMPI_Info_create( MPI_Info *info){
2701 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2705 int PMPI_Info_free( MPI_Info *info){
2709 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2710 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2714 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2715 int *num_addresses, int *num_datatypes, int *combiner){
2719 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2720 int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2721 MPI_Datatype* array_of_datatypes){
2725 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2726 int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2727 int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2731 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2735 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){
2739 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
2743 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
2744 void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
2749 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
2753 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
2757 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
2761 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
2765 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
2769 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
2773 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
2777 int PMPI_Add_error_class( int *errorclass){
2781 int PMPI_Add_error_code( int errorclass, int *errorcode){
2785 int PMPI_Add_error_string( int errorcode, char *string){
2789 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
2793 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
2797 int PMPI_Info_delete(MPI_Info info, char *key){
2801 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
2805 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
2809 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
2813 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2817 int MPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2821 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){
2825 int PMPI_Grequest_complete( MPI_Request request){
2829 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
2833 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
2837 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2841 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
2845 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
2849 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
2853 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
2857 int PMPI_Open_port( MPI_Info info, char *port_name){
2861 int PMPI_Close_port(char *port_name){
2865 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2869 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){
2873 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
2874 int* array_of_maxprocs, MPI_Info* array_of_info, int root,
2875 MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
2879 int PMPI_Comm_get_parent( MPI_Comm *parent){