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])
629 smpi_group_use(*newgroup);
630 retval = MPI_SUCCESS;
636 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
637 MPI_Group * newgroup)
639 int retval, i, newrank, rank, size, index, add;
642 if (group == MPI_GROUP_NULL) {
643 retval = MPI_ERR_GROUP;
644 } else if (newgroup == NULL) {
645 retval = MPI_ERR_ARG;
650 size = smpi_group_size(group);
651 for (i = 0; i < n; i++) {
652 for (rank = ranges[i][0]; /* First */
653 rank >= 0 && rank <= ranges[i][1]; /* Last */
654 rank += ranges[i][2] /* Stride */ ) {
659 *newgroup = MPI_GROUP_EMPTY;
661 *newgroup = smpi_group_new(size);
663 while (newrank < size) {
664 for (i = 0; i < n; i++) {
666 for (rank = ranges[i][0]; /* First */
667 rank >= 0 && rank <= ranges[i][1]; /* Last */
668 rank += ranges[i][2] /* Stride */ ) {
669 if (rank == newrank) {
675 index = smpi_group_index(group, newrank);
676 smpi_group_set_mapping(*newgroup, index, newrank);
679 newrank++; //added to avoid looping, need to be checked ..
683 smpi_group_use(*newgroup);
684 retval = MPI_SUCCESS;
690 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
695 if (comm == MPI_COMM_NULL) {
696 retval = MPI_ERR_COMM;
697 } else if (rank == NULL) {
698 retval = MPI_ERR_ARG;
700 *rank = smpi_comm_rank(comm);
701 retval = MPI_SUCCESS;
707 int PMPI_Comm_size(MPI_Comm comm, int *size)
712 if (comm == MPI_COMM_NULL) {
713 retval = MPI_ERR_COMM;
714 } else if (size == NULL) {
715 retval = MPI_ERR_ARG;
717 *size = smpi_comm_size(comm);
718 retval = MPI_SUCCESS;
724 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
729 if (comm == MPI_COMM_NULL) {
730 retval = MPI_ERR_COMM;
731 } else if (name == NULL || len == NULL) {
732 retval = MPI_ERR_ARG;
734 smpi_comm_get_name(comm, name, len);
735 retval = MPI_SUCCESS;
741 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
746 if (comm == MPI_COMM_NULL) {
747 retval = MPI_ERR_COMM;
748 } else if (group == NULL) {
749 retval = MPI_ERR_ARG;
751 *group = smpi_comm_group(comm);
752 retval = MPI_SUCCESS;
758 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
763 if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
764 retval = MPI_ERR_COMM;
765 } else if (result == NULL) {
766 retval = MPI_ERR_ARG;
768 if (comm1 == comm2) { /* Same communicators means same groups */
772 smpi_group_compare(smpi_comm_group(comm1),
773 smpi_comm_group(comm2));
774 if (*result == MPI_IDENT) {
775 *result = MPI_CONGRUENT;
778 retval = MPI_SUCCESS;
784 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
789 if (comm == MPI_COMM_NULL) {
790 retval = MPI_ERR_COMM;
791 } else if (newcomm == NULL) {
792 retval = MPI_ERR_ARG;
794 *newcomm = smpi_comm_new(smpi_comm_group(comm));
795 retval = MPI_SUCCESS;
801 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
806 if (comm == MPI_COMM_NULL) {
807 retval = MPI_ERR_COMM;
808 } else if (group == MPI_GROUP_NULL) {
809 retval = MPI_ERR_GROUP;
810 } else if (newcomm == NULL) {
811 retval = MPI_ERR_ARG;
812 } else if(smpi_group_rank(group,smpi_process_index())==MPI_UNDEFINED){
813 *newcomm= MPI_COMM_NULL;
814 retval = MPI_SUCCESS;
817 *newcomm = smpi_comm_new(group);
818 retval = MPI_SUCCESS;
824 int PMPI_Comm_free(MPI_Comm * comm)
830 retval = MPI_ERR_ARG;
831 } else if (*comm == MPI_COMM_NULL) {
832 retval = MPI_ERR_COMM;
834 smpi_comm_destroy(*comm);
835 *comm = MPI_COMM_NULL;
836 retval = MPI_SUCCESS;
842 int PMPI_Comm_disconnect(MPI_Comm * comm)
844 /* TODO: wait until all communication in comm are done */
849 retval = MPI_ERR_ARG;
850 } else if (*comm == MPI_COMM_NULL) {
851 retval = MPI_ERR_COMM;
853 smpi_comm_destroy(*comm);
854 *comm = MPI_COMM_NULL;
855 retval = MPI_SUCCESS;
861 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
866 if (comm_out == NULL) {
867 retval = MPI_ERR_ARG;
868 } else if (comm == MPI_COMM_NULL) {
869 retval = MPI_ERR_COMM;
871 *comm_out = smpi_comm_split(comm, color, key);
872 retval = MPI_SUCCESS;
878 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
879 int tag, MPI_Comm comm, MPI_Request * request)
884 if (request == NULL) {
885 retval = MPI_ERR_ARG;
886 } else if (comm == MPI_COMM_NULL) {
887 retval = MPI_ERR_COMM;
888 } else if (dst == MPI_PROC_NULL) {
889 retval = MPI_SUCCESS;
891 *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
892 retval = MPI_SUCCESS;
898 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
899 int tag, MPI_Comm comm, MPI_Request * request)
904 if (request == NULL) {
905 retval = MPI_ERR_ARG;
906 } else if (comm == MPI_COMM_NULL) {
907 retval = MPI_ERR_COMM;
908 } else if (src == MPI_PROC_NULL) {
909 retval = MPI_SUCCESS;
911 *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
912 retval = MPI_SUCCESS;
918 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
922 if (request == NULL) {
923 retval = MPI_ERR_ARG;
924 } else if (comm == MPI_COMM_NULL) {
925 retval = MPI_ERR_COMM;
926 } else if (dst == MPI_PROC_NULL) {
927 retval = MPI_SUCCESS;
929 *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
930 retval = MPI_SUCCESS;
936 int PMPI_Start(MPI_Request * request)
941 if (request == NULL || *request == MPI_REQUEST_NULL) {
942 retval = MPI_ERR_ARG;
944 smpi_mpi_start(*request);
945 retval = MPI_SUCCESS;
951 int PMPI_Startall(int count, MPI_Request * requests)
956 if (requests == NULL) {
957 retval = MPI_ERR_ARG;
959 smpi_mpi_startall(count, requests);
960 retval = MPI_SUCCESS;
966 int PMPI_Request_free(MPI_Request * request)
971 if (*request == MPI_REQUEST_NULL) {
972 retval = MPI_ERR_ARG;
974 if((*request)->flags & PERSISTENT)(*request)->refcount--;
975 smpi_mpi_request_free(request);
976 retval = MPI_SUCCESS;
982 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
983 int tag, MPI_Comm comm, MPI_Request * request)
989 if (request == NULL) {
990 retval = MPI_ERR_ARG;
991 } else if (comm == MPI_COMM_NULL) {
992 retval = MPI_ERR_COMM;
993 } else if (src == MPI_PROC_NULL) {
994 *request = MPI_REQUEST_NULL;
995 retval = MPI_SUCCESS;
996 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
997 retval = MPI_ERR_COMM;
998 } else if (count < 0) {
999 retval = MPI_ERR_COUNT;
1000 } else if (buf==NULL && count > 0) {
1001 retval = MPI_ERR_COUNT;
1002 } else if (datatype == MPI_DATATYPE_NULL){
1003 retval = MPI_ERR_TYPE;
1004 } else if(tag<0 && tag != MPI_ANY_TAG){
1005 retval = MPI_ERR_TAG;
1009 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1010 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1011 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
1014 *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
1015 retval = MPI_SUCCESS;
1018 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1019 (*request)->recv = 1;
1028 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
1029 int tag, MPI_Comm comm, MPI_Request * request)
1034 if (request == NULL) {
1035 retval = MPI_ERR_ARG;
1036 } else if (comm == MPI_COMM_NULL) {
1037 retval = MPI_ERR_COMM;
1038 } else if (dst == MPI_PROC_NULL) {
1039 *request = MPI_REQUEST_NULL;
1040 retval = MPI_SUCCESS;
1041 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1042 retval = MPI_ERR_RANK;
1043 } else if (count < 0) {
1044 retval = MPI_ERR_COUNT;
1045 } else if (buf==NULL && count > 0) {
1046 retval = MPI_ERR_COUNT;
1047 } else if (datatype == MPI_DATATYPE_NULL){
1048 retval = MPI_ERR_TYPE;
1049 } else if(tag<0 && tag != MPI_ANY_TAG){
1050 retval = MPI_ERR_TAG;
1054 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1055 TRACE_smpi_computing_out(rank);
1056 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1057 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1058 TRACE_smpi_send(rank, rank, dst_traced);
1061 *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1062 retval = MPI_SUCCESS;
1065 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1066 (*request)->send = 1;
1067 TRACE_smpi_computing_in(rank);
1075 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
1079 if (request == NULL) {
1080 retval = MPI_ERR_ARG;
1081 } else if (comm == MPI_COMM_NULL) {
1082 retval = MPI_ERR_COMM;
1083 } else if (dst == MPI_PROC_NULL) {
1084 *request = MPI_REQUEST_NULL;
1085 retval = MPI_SUCCESS;
1086 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1087 retval = MPI_ERR_RANK;
1088 } else if (count < 0) {
1089 retval = MPI_ERR_COUNT;
1090 } else if (buf==NULL && count > 0) {
1091 retval = MPI_ERR_COUNT;
1092 } else if (datatype == MPI_DATATYPE_NULL){
1093 retval = MPI_ERR_TYPE;
1094 } else if(tag<0 && tag != MPI_ANY_TAG){
1095 retval = MPI_ERR_TAG;
1099 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1100 TRACE_smpi_computing_out(rank);
1101 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1102 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1103 TRACE_smpi_send(rank, rank, dst_traced);
1106 *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1107 retval = MPI_SUCCESS;
1110 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1111 (*request)->send = 1;
1112 TRACE_smpi_computing_in(rank);
1120 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1121 MPI_Comm comm, MPI_Status * status)
1126 if (comm == MPI_COMM_NULL) {
1127 retval = MPI_ERR_COMM;
1128 } else if (src == MPI_PROC_NULL) {
1129 smpi_empty_status(status);
1130 status->MPI_SOURCE = MPI_PROC_NULL;
1131 retval = MPI_SUCCESS;
1132 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1133 retval = MPI_ERR_RANK;
1134 } else if (count < 0) {
1135 retval = MPI_ERR_COUNT;
1136 } else if (buf==NULL && count > 0) {
1137 retval = MPI_ERR_COUNT;
1138 } else if (datatype == MPI_DATATYPE_NULL){
1139 retval = MPI_ERR_TYPE;
1140 } else if(tag<0 && tag != MPI_ANY_TAG){
1141 retval = MPI_ERR_TAG;
1144 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1145 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1146 TRACE_smpi_computing_out(rank);
1148 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__);
1151 smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1152 retval = MPI_SUCCESS;
1155 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1156 if(status!=MPI_STATUS_IGNORE)src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1157 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1158 TRACE_smpi_recv(rank, src_traced, rank);
1159 TRACE_smpi_computing_in(rank);
1167 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1174 if (comm == MPI_COMM_NULL) {
1175 retval = MPI_ERR_COMM;
1176 } else if (dst == MPI_PROC_NULL) {
1177 retval = MPI_SUCCESS;
1178 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1179 retval = MPI_ERR_RANK;
1180 } else if (count < 0) {
1181 retval = MPI_ERR_COUNT;
1182 } else if (buf==NULL && count > 0) {
1183 retval = MPI_ERR_COUNT;
1184 } else if (datatype == MPI_DATATYPE_NULL){
1185 retval = MPI_ERR_TYPE;
1186 } else if(tag<0 && tag != MPI_ANY_TAG){
1187 retval = MPI_ERR_TAG;
1191 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1192 TRACE_smpi_computing_out(rank);
1193 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1194 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1195 TRACE_smpi_send(rank, rank, dst_traced);
1198 smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1199 retval = MPI_SUCCESS;
1202 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1203 TRACE_smpi_computing_in(rank);
1213 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1218 if (comm == MPI_COMM_NULL) {
1219 retval = MPI_ERR_COMM;
1220 } else if (dst == MPI_PROC_NULL) {
1221 retval = MPI_SUCCESS;
1222 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1223 retval = MPI_ERR_RANK;
1224 } else if (count < 0) {
1225 retval = MPI_ERR_COUNT;
1226 } else if (buf==NULL && count > 0) {
1227 retval = MPI_ERR_COUNT;
1228 } else if (datatype == MPI_DATATYPE_NULL){
1229 retval = MPI_ERR_TYPE;
1230 } else if(tag<0 && tag != MPI_ANY_TAG){
1231 retval = MPI_ERR_TAG;
1235 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1236 TRACE_smpi_computing_out(rank);
1237 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1238 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__);
1239 TRACE_smpi_send(rank, rank, dst_traced);
1242 smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1243 retval = MPI_SUCCESS;
1246 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1247 TRACE_smpi_computing_in(rank);
1255 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1256 int dst, int sendtag, void *recvbuf, int recvcount,
1257 MPI_Datatype recvtype, int src, int recvtag,
1258 MPI_Comm comm, MPI_Status * status)
1264 if (comm == MPI_COMM_NULL) {
1265 retval = MPI_ERR_COMM;
1266 } else if (sendtype == MPI_DATATYPE_NULL
1267 || recvtype == MPI_DATATYPE_NULL) {
1268 retval = MPI_ERR_TYPE;
1269 } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1270 smpi_empty_status(status);
1271 status->MPI_SOURCE = MPI_PROC_NULL;
1272 retval = MPI_SUCCESS;
1273 }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1274 (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1275 retval = MPI_ERR_RANK;
1276 } else if (sendcount < 0 || recvcount<0) {
1277 retval = MPI_ERR_COUNT;
1278 } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1279 retval = MPI_ERR_COUNT;
1280 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1281 retval = MPI_ERR_TAG;
1285 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1286 TRACE_smpi_computing_out(rank);
1287 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1288 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1289 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1290 TRACE_smpi_send(rank, rank, dst_traced);
1294 smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1295 recvcount, recvtype, src, recvtag, comm, status);
1296 retval = MPI_SUCCESS;
1299 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1300 TRACE_smpi_recv(rank, src_traced, rank);
1301 TRACE_smpi_computing_in(rank);
1310 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1311 int dst, int sendtag, int src, int recvtag,
1312 MPI_Comm comm, MPI_Status * status)
1314 //TODO: suboptimal implementation
1317 if (datatype == MPI_DATATYPE_NULL) {
1318 retval = MPI_ERR_TYPE;
1319 } else if (count < 0) {
1320 retval = MPI_ERR_COUNT;
1322 int size = smpi_datatype_get_extent(datatype) * count;
1323 recvbuf = xbt_new(char, size);
1325 MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1326 datatype, src, recvtag, comm, status);
1327 if(retval==MPI_SUCCESS){
1328 smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1336 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1341 if (request == MPI_REQUEST_NULL || flag == NULL) {
1342 retval = MPI_ERR_ARG;
1343 } else if (*request == MPI_REQUEST_NULL) {
1345 retval = MPI_ERR_REQUEST;
1347 *flag = smpi_mpi_test(request, status);
1348 retval = MPI_SUCCESS;
1354 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1355 MPI_Status * status)
1360 if (index == NULL || flag == NULL) {
1361 retval = MPI_ERR_ARG;
1363 *flag = smpi_mpi_testany(count, requests, index, status);
1364 retval = MPI_SUCCESS;
1370 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1376 retval = MPI_ERR_ARG;
1378 *flag = smpi_mpi_testall(count, requests, statuses);
1379 retval = MPI_SUCCESS;
1385 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1389 if (status == NULL) {
1390 retval = MPI_ERR_ARG;
1391 } else if (comm == MPI_COMM_NULL) {
1392 retval = MPI_ERR_COMM;
1393 } else if (source == MPI_PROC_NULL) {
1394 smpi_empty_status(status);
1395 status->MPI_SOURCE = MPI_PROC_NULL;
1396 retval = MPI_SUCCESS;
1398 smpi_mpi_probe(source, tag, comm, status);
1399 retval = MPI_SUCCESS;
1406 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1411 retval = MPI_ERR_ARG;
1412 } else if (status == NULL) {
1413 retval = MPI_ERR_ARG;
1414 } else if (comm == MPI_COMM_NULL) {
1415 retval = MPI_ERR_COMM;
1416 } else if (source == MPI_PROC_NULL) {
1417 smpi_empty_status(status);
1418 status->MPI_SOURCE = MPI_PROC_NULL;
1419 retval = MPI_SUCCESS;
1421 smpi_mpi_iprobe(source, tag, comm, flag, status);
1422 retval = MPI_SUCCESS;
1428 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1434 if (request == NULL) {
1435 retval = MPI_ERR_ARG;
1436 } else if (*request == MPI_REQUEST_NULL) {
1437 retval = MPI_ERR_REQUEST;
1441 int rank = request && (*request)->comm != MPI_COMM_NULL
1442 ? smpi_process_index()
1444 TRACE_smpi_computing_out(rank);
1446 int src_traced = (*request)->src;
1447 int dst_traced = (*request)->dst;
1448 int is_wait_for_receive = (*request)->recv;
1449 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__);
1452 smpi_mpi_wait(request, status);
1453 retval = MPI_SUCCESS;
1456 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1457 if (is_wait_for_receive) {
1458 TRACE_smpi_recv(rank, src_traced, dst_traced);
1460 TRACE_smpi_computing_in(rank);
1469 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1475 //save requests information for tracing
1477 int *srcs = xbt_new(int, count);
1478 int *dsts = xbt_new(int, count);
1479 int *recvs = xbt_new(int, count);
1480 for (i = 0; i < count; i++) {
1481 MPI_Request req = requests[i]; //already received requests are no longer valid
1485 recvs[i] = req->recv;
1488 int rank_traced = smpi_process_index();
1489 TRACE_smpi_computing_out(rank_traced);
1491 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1494 if (index == NULL) {
1495 retval = MPI_ERR_ARG;
1497 *index = smpi_mpi_waitany(count, requests, status);
1498 retval = MPI_SUCCESS;
1501 if(*index!=MPI_UNDEFINED){
1502 int src_traced = srcs[*index];
1503 int dst_traced = dsts[*index];
1504 int is_wait_for_receive = recvs[*index];
1505 if (is_wait_for_receive) {
1506 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1508 TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1513 TRACE_smpi_computing_in(rank_traced);
1519 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1524 //save information from requests
1526 int *srcs = xbt_new(int, count);
1527 int *dsts = xbt_new(int, count);
1528 int *recvs = xbt_new(int, count);
1529 int valid_count = 0;
1530 for (i = 0; i < count; i++) {
1531 MPI_Request req = requests[i];
1532 if(req!=MPI_REQUEST_NULL){
1533 srcs[valid_count] = req->src;
1534 dsts[valid_count] = req->dst;
1535 recvs[valid_count] = req->recv;
1539 int rank_traced = smpi_process_index();
1540 TRACE_smpi_computing_out(rank_traced);
1542 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__);
1544 int retval = smpi_mpi_waitall(count, requests, status);
1546 for (i = 0; i < valid_count; i++) {
1547 int src_traced = srcs[i];
1548 int dst_traced = dsts[i];
1549 int is_wait_for_receive = recvs[i];
1550 if (is_wait_for_receive) {
1551 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1554 TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1558 TRACE_smpi_computing_in(rank_traced);
1564 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1565 int *indices, MPI_Status status[])
1570 if (outcount == NULL) {
1571 retval = MPI_ERR_ARG;
1573 *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1574 retval = MPI_SUCCESS;
1580 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1581 int* indices, MPI_Status status[])
1586 if (outcount == NULL) {
1587 retval = MPI_ERR_ARG;
1589 *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1590 retval = MPI_SUCCESS;
1597 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1603 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1604 TRACE_smpi_computing_out(rank);
1605 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1606 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1608 if (comm == MPI_COMM_NULL) {
1609 retval = MPI_ERR_COMM;
1611 mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1612 retval = MPI_SUCCESS;
1615 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1616 TRACE_smpi_computing_in(rank);
1622 int PMPI_Barrier(MPI_Comm comm)
1628 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1629 TRACE_smpi_computing_out(rank);
1630 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1632 if (comm == MPI_COMM_NULL) {
1633 retval = MPI_ERR_COMM;
1635 mpi_coll_barrier_fun(comm);
1636 retval = MPI_SUCCESS;
1639 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1640 TRACE_smpi_computing_in(rank);
1646 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1647 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1648 int root, MPI_Comm comm)
1654 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1655 TRACE_smpi_computing_out(rank);
1656 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1657 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1659 if (comm == MPI_COMM_NULL) {
1660 retval = MPI_ERR_COMM;
1661 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1662 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1663 retval = MPI_ERR_TYPE;
1664 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1665 ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1666 retval = MPI_ERR_COUNT;
1669 char* sendtmpbuf = (char*) sendbuf;
1670 int sendtmpcount = sendcount;
1671 MPI_Datatype sendtmptype = sendtype;
1672 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1674 sendtmptype=recvtype;
1677 mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount,
1678 recvtype, root, comm);
1681 retval = MPI_SUCCESS;
1684 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1685 TRACE_smpi_computing_in(rank);
1691 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1692 void *recvbuf, int *recvcounts, int *displs,
1693 MPI_Datatype recvtype, int root, MPI_Comm comm)
1699 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1700 TRACE_smpi_computing_out(rank);
1701 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1702 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1704 if (comm == MPI_COMM_NULL) {
1705 retval = MPI_ERR_COMM;
1706 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1707 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1708 retval = MPI_ERR_TYPE;
1709 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1710 retval = MPI_ERR_COUNT;
1711 } else if (recvcounts == NULL || displs == NULL) {
1712 retval = MPI_ERR_ARG;
1715 char* sendtmpbuf = (char*) sendbuf;
1716 int sendtmpcount = sendcount;
1717 MPI_Datatype sendtmptype = sendtype;
1718 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1720 sendtmptype=recvtype;
1723 smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts,
1724 displs, recvtype, root, comm);
1725 retval = MPI_SUCCESS;
1728 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1729 TRACE_smpi_computing_in(rank);
1735 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1736 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1743 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1744 TRACE_smpi_computing_out(rank);
1745 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1747 if (comm == MPI_COMM_NULL) {
1748 retval = MPI_ERR_COMM;
1749 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1750 (recvtype == MPI_DATATYPE_NULL)){
1751 retval = MPI_ERR_TYPE;
1752 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1754 retval = MPI_ERR_COUNT;
1757 if(sendbuf == MPI_IN_PLACE) {
1758 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1759 sendcount=recvcount;
1763 mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1765 retval = MPI_SUCCESS;
1768 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1774 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1775 void *recvbuf, int *recvcounts, int *displs,
1776 MPI_Datatype recvtype, MPI_Comm comm)
1782 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1783 TRACE_smpi_computing_out(rank);
1784 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1786 if (comm == MPI_COMM_NULL) {
1787 retval = MPI_ERR_COMM;
1788 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1789 (recvtype == MPI_DATATYPE_NULL)){
1790 retval = MPI_ERR_TYPE;
1791 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1792 retval = MPI_ERR_COUNT;
1793 } else if (recvcounts == NULL || displs == NULL) {
1794 retval = MPI_ERR_ARG;
1797 if(sendbuf == MPI_IN_PLACE) {
1798 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1799 sendcount=recvcounts[smpi_comm_rank(comm)];
1803 mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1804 displs, recvtype, comm);
1805 retval = MPI_SUCCESS;
1808 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1809 TRACE_smpi_computing_in(rank);
1815 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1816 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1817 int root, MPI_Comm comm)
1823 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1824 TRACE_smpi_computing_out(rank);
1825 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1827 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1829 if (comm == MPI_COMM_NULL) {
1830 retval = MPI_ERR_COMM;
1831 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1832 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1833 retval = MPI_ERR_TYPE;
1836 if(recvbuf==MPI_IN_PLACE){
1839 mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1840 recvtype, root, comm);
1841 retval = MPI_SUCCESS;
1844 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1845 TRACE_smpi_computing_in(rank);
1851 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1852 MPI_Datatype sendtype, void *recvbuf, int recvcount,
1853 MPI_Datatype recvtype, int root, MPI_Comm comm)
1859 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1860 TRACE_smpi_computing_out(rank);
1861 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1862 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1864 if (comm == MPI_COMM_NULL) {
1865 retval = MPI_ERR_COMM;
1866 } else if (sendcounts == NULL || displs == NULL) {
1867 retval = MPI_ERR_ARG;
1868 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1869 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1870 retval = MPI_ERR_TYPE;
1873 if(recvbuf==MPI_IN_PLACE){
1877 smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1878 recvcount, recvtype, root, comm);
1879 retval = MPI_SUCCESS;
1882 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1883 TRACE_smpi_computing_in(rank);
1889 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1890 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1896 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1897 TRACE_smpi_computing_out(rank);
1898 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1899 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__);
1901 if (comm == MPI_COMM_NULL) {
1902 retval = MPI_ERR_COMM;
1903 } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1904 retval = MPI_ERR_ARG;
1907 char* sendtmpbuf = (char*) sendbuf;
1908 if( sendbuf == MPI_IN_PLACE ) {
1909 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
1910 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
1913 mpi_coll_reduce_fun(sendtmpbuf, recvbuf, count, datatype, op, root, comm);
1915 if( sendbuf == MPI_IN_PLACE ) {
1916 xbt_free(sendtmpbuf);
1919 retval = MPI_SUCCESS;
1922 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1923 TRACE_smpi_computing_in(rank);
1929 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
1930 MPI_Datatype datatype, MPI_Op op){
1934 if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1935 retval = MPI_ERR_ARG;
1937 smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
1944 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
1945 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1951 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1952 TRACE_smpi_computing_out(rank);
1953 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1955 if (comm == MPI_COMM_NULL) {
1956 retval = MPI_ERR_COMM;
1957 } else if (datatype == MPI_DATATYPE_NULL) {
1958 retval = MPI_ERR_TYPE;
1959 } else if (op == MPI_OP_NULL) {
1960 retval = MPI_ERR_OP;
1963 char* sendtmpbuf = (char*) sendbuf;
1964 if( sendbuf == MPI_IN_PLACE ) {
1965 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
1966 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
1969 mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
1971 if( sendbuf == MPI_IN_PLACE ) {
1972 xbt_free(sendtmpbuf);
1975 retval = MPI_SUCCESS;
1979 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1980 TRACE_smpi_computing_in(rank);
1986 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
1987 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
1993 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1994 TRACE_smpi_computing_out(rank);
1995 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
1997 if (comm == MPI_COMM_NULL) {
1998 retval = MPI_ERR_COMM;
1999 } else if (datatype == MPI_DATATYPE_NULL) {
2000 retval = MPI_ERR_TYPE;
2001 } else if (op == MPI_OP_NULL) {
2002 retval = MPI_ERR_OP;
2004 smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2005 retval = MPI_SUCCESS;
2008 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2009 TRACE_smpi_computing_in(rank);
2015 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2016 MPI_Op op, MPI_Comm comm){
2021 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2022 TRACE_smpi_computing_out(rank);
2023 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2025 if (comm == MPI_COMM_NULL) {
2026 retval = MPI_ERR_COMM;
2027 } else if (datatype == MPI_DATATYPE_NULL) {
2028 retval = MPI_ERR_TYPE;
2029 } else if (op == MPI_OP_NULL) {
2030 retval = MPI_ERR_OP;
2032 smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2033 retval = MPI_SUCCESS;
2036 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2037 TRACE_smpi_computing_in(rank);
2043 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
2044 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2049 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2050 TRACE_smpi_computing_out(rank);
2051 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2053 if (comm == MPI_COMM_NULL) {
2054 retval = MPI_ERR_COMM;
2055 } else if (datatype == MPI_DATATYPE_NULL) {
2056 retval = MPI_ERR_TYPE;
2057 } else if (op == MPI_OP_NULL) {
2058 retval = MPI_ERR_OP;
2059 } else if (recvcounts == NULL) {
2060 retval = MPI_ERR_ARG;
2062 void* sendtmpbuf=sendbuf;
2063 if(sendbuf==MPI_IN_PLACE){
2067 mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts,
2068 datatype, op, comm);
2069 retval = MPI_SUCCESS;
2072 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2073 TRACE_smpi_computing_in(rank);
2079 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2080 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2085 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2086 TRACE_smpi_computing_out(rank);
2087 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2089 if (comm == MPI_COMM_NULL) {
2090 retval = MPI_ERR_COMM;
2091 } else if (datatype == MPI_DATATYPE_NULL) {
2092 retval = MPI_ERR_TYPE;
2093 } else if (op == MPI_OP_NULL) {
2094 retval = MPI_ERR_OP;
2095 } else if (recvcount < 0) {
2096 retval = MPI_ERR_ARG;
2098 int count=smpi_comm_size(comm);
2099 int* recvcounts=(int*)xbt_malloc(count);
2100 for (i=0; i<count;i++)recvcounts[i]=recvcount;
2101 mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
2102 datatype, op, comm);
2103 xbt_free(recvcounts);
2104 retval = MPI_SUCCESS;
2107 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2108 TRACE_smpi_computing_in(rank);
2114 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2115 void *recvbuf, int recvcount, MPI_Datatype recvtype,
2122 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2123 TRACE_smpi_computing_out(rank);
2124 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2126 if (comm == MPI_COMM_NULL) {
2127 retval = MPI_ERR_COMM;
2128 } else if (sendtype == MPI_DATATYPE_NULL
2129 || recvtype == MPI_DATATYPE_NULL) {
2130 retval = MPI_ERR_TYPE;
2132 retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2135 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2136 TRACE_smpi_computing_in(rank);
2142 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2143 MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2144 int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2150 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2151 TRACE_smpi_computing_out(rank);
2152 TRACE_smpi_collective_in(rank, -1, __FUNCTION__);
2154 if (comm == MPI_COMM_NULL) {
2155 retval = MPI_ERR_COMM;
2156 } else if (sendtype == MPI_DATATYPE_NULL
2157 || recvtype == MPI_DATATYPE_NULL) {
2158 retval = MPI_ERR_TYPE;
2159 } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2160 || recvdisps == NULL) {
2161 retval = MPI_ERR_ARG;
2164 mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2165 recvbuf, recvcounts, recvdisps, recvtype,
2169 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2170 TRACE_smpi_computing_in(rank);
2177 int PMPI_Get_processor_name(char *name, int *resultlen)
2179 int retval = MPI_SUCCESS;
2182 strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2183 strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2184 strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2185 MPI_MAX_PROCESSOR_NAME - 1 );
2188 MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2194 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2196 int retval = MPI_SUCCESS;
2200 if (status == NULL || count == NULL) {
2201 retval = MPI_ERR_ARG;
2202 } else if (datatype == MPI_DATATYPE_NULL) {
2203 retval = MPI_ERR_TYPE;
2205 size = smpi_datatype_size(datatype);
2208 } else if (status->count % size != 0) {
2209 retval = MPI_UNDEFINED;
2211 *count = smpi_mpi_get_count(status, datatype);
2218 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2222 if (old_type == MPI_DATATYPE_NULL) {
2223 retval = MPI_ERR_TYPE;
2224 } else if (count<0){
2225 retval = MPI_ERR_COUNT;
2227 retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2233 int PMPI_Type_commit(MPI_Datatype* datatype) {
2237 if (datatype == MPI_DATATYPE_NULL) {
2238 retval = MPI_ERR_TYPE;
2240 smpi_datatype_commit(datatype);
2241 retval = MPI_SUCCESS;
2248 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2252 if (old_type == MPI_DATATYPE_NULL) {
2253 retval = MPI_ERR_TYPE;
2254 } else if (count<0 || blocklen<0){
2255 retval = MPI_ERR_COUNT;
2257 retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2263 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2267 if (old_type == MPI_DATATYPE_NULL) {
2268 retval = MPI_ERR_TYPE;
2269 } else if (count<0 || blocklen<0){
2270 retval = MPI_ERR_COUNT;
2272 retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2278 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2279 return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2282 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2286 if (old_type == MPI_DATATYPE_NULL) {
2287 retval = MPI_ERR_TYPE;
2288 } else if (count<0){
2289 retval = MPI_ERR_COUNT;
2291 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2297 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2301 if (old_type == MPI_DATATYPE_NULL) {
2302 retval = MPI_ERR_TYPE;
2303 } else if (count<0){
2304 retval = MPI_ERR_COUNT;
2306 int* blocklens=(int*)xbt_malloc(blocklength*count);
2307 for (i=0; i<count;i++)blocklens[i]=blocklength;
2308 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2309 xbt_free(blocklens);
2316 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2320 if (old_type == MPI_DATATYPE_NULL) {
2321 retval = MPI_ERR_TYPE;
2322 } else if (count<0){
2323 retval = MPI_ERR_COUNT;
2325 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2331 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2335 if (old_type == MPI_DATATYPE_NULL) {
2336 retval = MPI_ERR_TYPE;
2337 } else if (count<0){
2338 retval = MPI_ERR_COUNT;
2340 int* blocklens=(int*)xbt_malloc(blocklength*count);
2341 for (i=0; i<count;i++)blocklens[i]=blocklength;
2342 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2343 xbt_free(blocklens);
2350 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2355 retval = MPI_ERR_COUNT;
2357 retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2363 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2364 return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2368 int PMPI_Error_class(int errorcode, int* errorclass) {
2369 // assume smpi uses only standard mpi error codes
2370 *errorclass=errorcode;
2375 int PMPI_Initialized(int* flag) {
2376 *flag=(smpi_process_data()!=NULL);
2380 /* The following calls are not yet implemented and will fail at runtime. */
2381 /* Once implemented, please move them above this notice. */
2383 #define NOT_YET_IMPLEMENTED {\
2384 XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__);\
2385 return MPI_SUCCESS;\
2389 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2393 int PMPI_Type_set_name(MPI_Datatype datatype, char * name)
2398 int PMPI_Type_get_name(MPI_Datatype datatype, char * name, int* len)
2403 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2407 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2411 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
2415 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2419 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2423 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2427 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2431 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2435 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2439 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2443 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2447 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2451 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2455 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2459 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2463 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2467 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2471 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2475 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2479 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2483 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2487 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2491 int PMPI_Cancel(MPI_Request* request) {
2495 int PMPI_Buffer_attach(void* buffer, int size) {
2499 int PMPI_Buffer_detach(void* buffer, int* size) {
2503 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2507 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2512 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2517 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2522 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2527 int PMPI_Comm_free_keyval(int* keyval) {
2531 int PMPI_Pcontrol(const int level )
2536 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2540 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2545 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2550 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2555 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2560 int PMPI_Type_free_keyval(int* keyval) {
2564 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2568 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2572 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2576 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2580 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2584 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2588 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2592 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2596 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2600 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2604 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2608 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2612 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2616 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2620 int PMPI_Keyval_free(int* keyval) {
2624 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2628 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2632 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2636 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2640 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2644 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2648 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2652 int PMPI_Win_fence( int assert, MPI_Win win){
2656 int PMPI_Win_free( MPI_Win* win){
2660 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2664 int PMPI_Info_create( MPI_Info *info){
2668 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2672 int PMPI_Info_free( MPI_Info *info){
2676 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2677 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2681 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2682 int *num_addresses, int *num_datatypes, int *combiner){
2686 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2687 int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2688 MPI_Datatype* array_of_datatypes){
2692 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2693 int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2694 int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2698 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2702 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){
2706 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
2710 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
2711 void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
2716 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
2720 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
2724 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
2728 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
2732 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
2736 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
2740 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
2744 int PMPI_Add_error_class( int *errorclass){
2748 int PMPI_Add_error_code( int errorclass, int *errorcode){
2752 int PMPI_Add_error_string( int errorcode, char *string){
2756 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
2760 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
2764 int PMPI_Info_delete(MPI_Info info, char *key){
2768 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
2772 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
2776 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
2780 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2784 int MPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2788 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){
2792 int PMPI_Grequest_complete( MPI_Request request){
2796 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
2800 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
2804 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2808 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
2812 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
2816 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
2820 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
2824 int PMPI_Open_port( MPI_Info info, char *port_name){
2828 int PMPI_Close_port(char *port_name){
2832 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2836 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){
2840 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
2841 int* array_of_maxprocs, MPI_Info* array_of_info, int root,
2842 MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
2846 int PMPI_Comm_get_parent( MPI_Comm *parent){