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);
30 smpi_process_mark_as_initialized();
32 int rank = smpi_process_index();
33 TRACE_smpi_init(rank);
35 TRACE_smpi_computing_init(rank);
41 int PMPI_Finalize(void)
43 smpi_process_finalize();
46 int rank = smpi_process_index();
47 TRACE_smpi_computing_out(rank);
48 TRACE_smpi_finalize(smpi_process_index());
50 smpi_process_destroy();
54 int PMPI_Finalized(int* flag)
56 *flag=smpi_process_finalized();
60 int PMPI_Get_version (int *version,int *subversion){
61 *version = MPI_VERSION;
62 *subversion= MPI_SUBVERSION;
66 int PMPI_Get_library_version (char *version,int *len){
67 int retval = MPI_SUCCESS;
69 snprintf(version,MPI_MAX_LIBRARY_VERSION_STRING,"SMPI Version %d.%d. Copyright The Simgrid Team 2007-2013",SIMGRID_VERSION_MAJOR,
70 SIMGRID_VERSION_MINOR);
71 *len = strlen(version) > MPI_MAX_LIBRARY_VERSION_STRING ? MPI_MAX_LIBRARY_VERSION_STRING : strlen(version);
76 int PMPI_Init_thread(int *argc, char ***argv, int required, int *provided)
78 if (provided != NULL) {
79 *provided = MPI_THREAD_MULTIPLE;
81 return MPI_Init(argc, argv);
84 int PMPI_Query_thread(int *provided)
89 if (provided == NULL) {
92 *provided = MPI_THREAD_MULTIPLE;
99 int PMPI_Is_thread_main(int *flag)
105 retval = MPI_ERR_ARG;
107 *flag = smpi_process_index() == 0;
108 retval = MPI_SUCCESS;
114 int PMPI_Abort(MPI_Comm comm, int errorcode)
117 smpi_process_destroy();
119 int rank = smpi_process_index();
120 TRACE_smpi_computing_out(rank);
122 // FIXME: should kill all processes in comm instead
123 simcall_process_kill(SIMIX_process_self());
127 double PMPI_Wtime(void)
132 time = SIMIX_get_clock();
136 extern double sg_maxmin_precision;
137 double PMPI_Wtick(void)
139 return sg_maxmin_precision;
142 int PMPI_Address(void *location, MPI_Aint * address)
148 retval = MPI_ERR_ARG;
150 *address = (MPI_Aint) location;
151 retval = MPI_SUCCESS;
157 int PMPI_Get_address(void *location, MPI_Aint * address)
159 return PMPI_Address(location, address);
162 int PMPI_Type_free(MPI_Datatype * datatype)
168 retval = MPI_ERR_ARG;
170 smpi_datatype_free(datatype);
171 retval = MPI_SUCCESS;
177 int PMPI_Type_size(MPI_Datatype datatype, int *size)
182 if (datatype == MPI_DATATYPE_NULL) {
183 retval = MPI_ERR_TYPE;
184 } else if (size == NULL) {
185 retval = MPI_ERR_ARG;
187 *size = (int) smpi_datatype_size(datatype);
188 retval = MPI_SUCCESS;
194 int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
199 if (datatype == MPI_DATATYPE_NULL) {
200 retval = MPI_ERR_TYPE;
201 } else if (lb == NULL || extent == NULL) {
202 retval = MPI_ERR_ARG;
204 retval = smpi_datatype_extent(datatype, lb, extent);
210 int PMPI_Type_get_true_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
212 return PMPI_Type_get_extent(datatype, lb, extent);
215 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
220 if (datatype == MPI_DATATYPE_NULL) {
221 retval = MPI_ERR_TYPE;
222 } else if (extent == NULL) {
223 retval = MPI_ERR_ARG;
225 *extent = smpi_datatype_get_extent(datatype);
226 retval = MPI_SUCCESS;
232 int PMPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
237 if (datatype == MPI_DATATYPE_NULL) {
238 retval = MPI_ERR_TYPE;
239 } else if (disp == NULL) {
240 retval = MPI_ERR_ARG;
242 *disp = smpi_datatype_lb(datatype);
243 retval = MPI_SUCCESS;
249 int PMPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
254 if (datatype == MPI_DATATYPE_NULL) {
255 retval = MPI_ERR_TYPE;
256 } else if (disp == NULL) {
257 retval = MPI_ERR_ARG;
259 *disp = smpi_datatype_ub(datatype);
260 retval = MPI_SUCCESS;
266 int PMPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
271 if (function == NULL || op == NULL) {
272 retval = MPI_ERR_ARG;
274 *op = smpi_op_new(function, commute);
275 retval = MPI_SUCCESS;
281 int PMPI_Op_free(MPI_Op * op)
287 retval = MPI_ERR_ARG;
288 } else if (*op == MPI_OP_NULL) {
291 smpi_op_destroy(*op);
293 retval = MPI_SUCCESS;
299 int PMPI_Group_free(MPI_Group * group)
305 retval = MPI_ERR_ARG;
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 retval = MPI_SUCCESS;
428 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
429 MPI_Group * newgroup)
431 int retval, i, proc1, proc2, size;
434 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
435 retval = MPI_ERR_GROUP;
436 } else if (newgroup == NULL) {
437 retval = MPI_ERR_ARG;
439 size = smpi_group_size(group2);
440 for (i = 0; i < size; i++) {
441 proc2 = smpi_group_index(group2, i);
442 proc1 = smpi_group_rank(group1, proc2);
443 if (proc1 == MPI_UNDEFINED) {
448 *newgroup = MPI_GROUP_EMPTY;
450 *newgroup = smpi_group_new(size);
452 for (i = 0; i < smpi_group_size(group2); i++) {
453 proc2 = smpi_group_index(group2, i);
454 proc1 = smpi_group_rank(group1, proc2);
455 if (proc1 != MPI_UNDEFINED) {
456 smpi_group_set_mapping(*newgroup, proc2, j);
461 retval = MPI_SUCCESS;
467 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
469 int retval, i, proc1, proc2, size, size2;
472 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
473 retval = MPI_ERR_GROUP;
474 } else if (newgroup == NULL) {
475 retval = MPI_ERR_ARG;
477 size = size2 = smpi_group_size(group1);
478 for (i = 0; i < size2; i++) {
479 proc1 = smpi_group_index(group1, i);
480 proc2 = smpi_group_rank(group2, proc1);
481 if (proc2 != MPI_UNDEFINED) {
486 *newgroup = MPI_GROUP_EMPTY;
488 *newgroup = smpi_group_new(size);
489 for (i = 0; i < size2; i++) {
490 proc1 = smpi_group_index(group1, i);
491 proc2 = smpi_group_rank(group2, proc1);
492 if (proc2 == MPI_UNDEFINED) {
493 smpi_group_set_mapping(*newgroup, proc1, i);
497 retval = MPI_SUCCESS;
503 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
505 int retval, i, index;
508 if (group == MPI_GROUP_NULL) {
509 retval = MPI_ERR_GROUP;
510 } else if (newgroup == NULL) {
511 retval = MPI_ERR_ARG;
514 *newgroup = MPI_GROUP_EMPTY;
515 } else if (n == smpi_group_size(group)) {
517 if(group!= smpi_comm_group(MPI_COMM_WORLD)
518 && group != MPI_GROUP_NULL
519 && group != smpi_comm_group(MPI_COMM_SELF)
520 && group != MPI_GROUP_EMPTY)
521 smpi_group_use(group);
523 *newgroup = smpi_group_new(n);
524 for (i = 0; i < n; i++) {
525 index = smpi_group_index(group, ranks[i]);
526 smpi_group_set_mapping(*newgroup, index, i);
529 retval = MPI_SUCCESS;
535 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
537 int retval, i, j, newsize, oldsize, index;
540 if (group == MPI_GROUP_NULL) {
541 retval = MPI_ERR_GROUP;
542 } else if (newgroup == NULL) {
543 retval = MPI_ERR_ARG;
547 if(group!= smpi_comm_group(MPI_COMM_WORLD)
548 && group != MPI_GROUP_NULL
549 && group != smpi_comm_group(MPI_COMM_SELF)
550 && group != MPI_GROUP_EMPTY)
551 smpi_group_use(group);
552 } else if (n == smpi_group_size(group)) {
553 *newgroup = MPI_GROUP_EMPTY;
555 oldsize=smpi_group_size(group);
556 newsize = oldsize - n;
557 *newgroup = smpi_group_new(newsize);
559 int* to_exclude=xbt_new(int, smpi_group_size(group));
560 for(i=0; i<oldsize; i++)
563 to_exclude[ranks[i]]=1;
566 for(i=0; i<oldsize; i++){
567 if(to_exclude[i]==0){
568 index = smpi_group_index(group, i);
569 smpi_group_set_mapping(*newgroup, index, j);
574 xbt_free(to_exclude);
576 retval = MPI_SUCCESS;
582 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
583 MPI_Group * newgroup)
585 int retval, i, j, rank, size, index;
588 if (group == MPI_GROUP_NULL) {
589 retval = MPI_ERR_GROUP;
590 } else if (newgroup == NULL) {
591 retval = MPI_ERR_ARG;
594 *newgroup = MPI_GROUP_EMPTY;
597 for (i = 0; i < n; i++) {
598 for (rank = ranges[i][0]; /* First */
599 rank >= 0; /* Last */
603 rank += ranges[i][2]; /* Stride */
604 if (ranges[i][0]<ranges[i][1]){
605 if(rank > ranges[i][1])
608 if(rank < ranges[i][1])
614 *newgroup = smpi_group_new(size);
616 for (i = 0; i < n; i++) {
617 for (rank = ranges[i][0]; /* First */
618 rank >= 0; /* Last */
620 index = smpi_group_index(group, rank);
621 smpi_group_set_mapping(*newgroup, index, j);
623 rank += ranges[i][2]; /* Stride */
624 if (ranges[i][0]<ranges[i][1]){
625 if(rank > ranges[i][1])
628 if(rank < ranges[i][1])
634 retval = MPI_SUCCESS;
640 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
641 MPI_Group * newgroup)
643 int retval, i, rank, newrank,oldrank, size, index, add;
646 if (group == MPI_GROUP_NULL) {
647 retval = MPI_ERR_GROUP;
648 } else if (newgroup == NULL) {
649 retval = MPI_ERR_ARG;
653 if(group!= smpi_comm_group(MPI_COMM_WORLD)
654 && group != MPI_GROUP_NULL
655 && group != smpi_comm_group(MPI_COMM_SELF)
656 && group != MPI_GROUP_EMPTY)
657 smpi_group_use(group);
659 size = smpi_group_size(group);
660 for (i = 0; i < n; i++) {
661 for (rank = ranges[i][0]; /* First */
662 rank >= 0; /* Last */
666 rank += ranges[i][2]; /* Stride */
667 if (ranges[i][0]<ranges[i][1]){
668 if(rank > ranges[i][1])
671 if(rank < ranges[i][1])
677 *newgroup = MPI_GROUP_EMPTY;
679 *newgroup = smpi_group_new(size);
682 while (newrank < size) {
684 for (i = 0; i < n; i++) {
685 for (rank = ranges[i][0];rank >= 0;){
691 rank += ranges[i][2]; /* Stride */
693 if (ranges[i][0]<ranges[i][1]){
694 if(rank > ranges[i][1])
697 if(rank < ranges[i][1])
703 index = smpi_group_index(group, oldrank);
704 smpi_group_set_mapping(*newgroup, index, newrank);
712 retval = MPI_SUCCESS;
718 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
723 if (comm == MPI_COMM_NULL) {
724 retval = MPI_ERR_COMM;
725 } else if (rank == NULL) {
726 retval = MPI_ERR_ARG;
728 *rank = smpi_comm_rank(comm);
729 retval = MPI_SUCCESS;
735 int PMPI_Comm_size(MPI_Comm comm, int *size)
740 if (comm == MPI_COMM_NULL) {
741 retval = MPI_ERR_COMM;
742 } else if (size == NULL) {
743 retval = MPI_ERR_ARG;
745 *size = smpi_comm_size(comm);
746 retval = MPI_SUCCESS;
752 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
757 if (comm == MPI_COMM_NULL) {
758 retval = MPI_ERR_COMM;
759 } else if (name == NULL || len == NULL) {
760 retval = MPI_ERR_ARG;
762 smpi_comm_get_name(comm, name, len);
763 retval = MPI_SUCCESS;
769 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
774 if (comm == MPI_COMM_NULL) {
775 retval = MPI_ERR_COMM;
776 } else if (group == NULL) {
777 retval = MPI_ERR_ARG;
779 *group = smpi_comm_group(comm);
780 if(*group!= smpi_comm_group(MPI_COMM_WORLD)
781 && *group != MPI_GROUP_NULL
782 && *group != smpi_comm_group(MPI_COMM_SELF)
783 && *group != MPI_GROUP_EMPTY)
784 smpi_group_use(*group);
785 retval = MPI_SUCCESS;
791 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
796 if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
797 retval = MPI_ERR_COMM;
798 } else if (result == NULL) {
799 retval = MPI_ERR_ARG;
801 if (comm1 == comm2) { /* Same communicators means same groups */
805 smpi_group_compare(smpi_comm_group(comm1),
806 smpi_comm_group(comm2));
807 if (*result == MPI_IDENT) {
808 *result = MPI_CONGRUENT;
811 retval = MPI_SUCCESS;
817 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
822 if (comm == MPI_COMM_NULL) {
823 retval = MPI_ERR_COMM;
824 } else if (newcomm == NULL) {
825 retval = MPI_ERR_ARG;
827 *newcomm = smpi_comm_new(smpi_comm_group(comm));
828 retval = MPI_SUCCESS;
834 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
839 if (comm == MPI_COMM_NULL) {
840 retval = MPI_ERR_COMM;
841 } else if (group == MPI_GROUP_NULL) {
842 retval = MPI_ERR_GROUP;
843 } else if (newcomm == NULL) {
844 retval = MPI_ERR_ARG;
845 } else if(smpi_group_rank(group,smpi_process_index())==MPI_UNDEFINED){
846 *newcomm= MPI_COMM_NULL;
847 retval = MPI_SUCCESS;
850 *newcomm = smpi_comm_new(group);
851 retval = MPI_SUCCESS;
857 int PMPI_Comm_free(MPI_Comm * comm)
863 retval = MPI_ERR_ARG;
864 } else if (*comm == MPI_COMM_NULL) {
865 retval = MPI_ERR_COMM;
867 smpi_comm_destroy(*comm);
868 *comm = MPI_COMM_NULL;
869 retval = MPI_SUCCESS;
875 int PMPI_Comm_disconnect(MPI_Comm * comm)
877 /* TODO: wait until all communication in comm are done */
882 retval = MPI_ERR_ARG;
883 } else if (*comm == MPI_COMM_NULL) {
884 retval = MPI_ERR_COMM;
886 smpi_comm_destroy(*comm);
887 *comm = MPI_COMM_NULL;
888 retval = MPI_SUCCESS;
894 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
899 if (comm_out == NULL) {
900 retval = MPI_ERR_ARG;
901 } else if (comm == MPI_COMM_NULL) {
902 retval = MPI_ERR_COMM;
904 *comm_out = smpi_comm_split(comm, color, key);
905 retval = MPI_SUCCESS;
911 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
912 int tag, MPI_Comm comm, MPI_Request * request)
917 if (request == NULL) {
918 retval = MPI_ERR_ARG;
919 } else if (comm == MPI_COMM_NULL) {
920 retval = MPI_ERR_COMM;
921 } else if (dst == MPI_PROC_NULL) {
922 retval = MPI_SUCCESS;
924 *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
925 retval = MPI_SUCCESS;
928 if(retval!=MPI_SUCCESS)*request=MPI_REQUEST_NULL;
932 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
933 int tag, MPI_Comm comm, MPI_Request * request)
938 if (request == NULL) {
939 retval = MPI_ERR_ARG;
940 } else if (comm == MPI_COMM_NULL) {
941 retval = MPI_ERR_COMM;
942 } else if (src == MPI_PROC_NULL) {
943 retval = MPI_SUCCESS;
945 *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
946 retval = MPI_SUCCESS;
949 if(retval!=MPI_SUCCESS)*request=MPI_REQUEST_NULL;
953 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
957 if (request == NULL) {
958 retval = MPI_ERR_ARG;
959 } else if (comm == MPI_COMM_NULL) {
960 retval = MPI_ERR_COMM;
961 } else if (dst == MPI_PROC_NULL) {
962 retval = MPI_SUCCESS;
964 *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
965 retval = MPI_SUCCESS;
968 if(retval!=MPI_SUCCESS)*request=MPI_REQUEST_NULL;
972 int PMPI_Start(MPI_Request * request)
977 if (request == NULL || *request == MPI_REQUEST_NULL) {
978 retval = MPI_ERR_ARG;
980 smpi_mpi_start(*request);
981 retval = MPI_SUCCESS;
987 int PMPI_Startall(int count, MPI_Request * requests)
992 if (requests == NULL) {
993 retval = MPI_ERR_ARG;
995 smpi_mpi_startall(count, requests);
996 retval = MPI_SUCCESS;
1002 int PMPI_Request_free(MPI_Request * request)
1007 if (*request == MPI_REQUEST_NULL) {
1008 retval = MPI_ERR_ARG;
1010 if((*request)->flags & PERSISTENT)(*request)->refcount--;
1011 smpi_mpi_request_free(request);
1012 retval = MPI_SUCCESS;
1018 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
1019 int tag, MPI_Comm comm, MPI_Request * request)
1025 if (request == NULL) {
1026 retval = MPI_ERR_ARG;
1027 } else if (comm == MPI_COMM_NULL) {
1028 retval = MPI_ERR_COMM;
1029 } else if (src == MPI_PROC_NULL) {
1030 *request = MPI_REQUEST_NULL;
1031 retval = MPI_SUCCESS;
1032 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1033 retval = MPI_ERR_COMM;
1034 } else if (count < 0) {
1035 retval = MPI_ERR_COUNT;
1036 } else if (buf==NULL && count > 0) {
1037 retval = MPI_ERR_COUNT;
1038 } else if (datatype == MPI_DATATYPE_NULL){
1039 retval = MPI_ERR_TYPE;
1040 } else if(tag<0 && tag != MPI_ANY_TAG){
1041 retval = MPI_ERR_TAG;
1045 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1046 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1047 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, count*smpi_datatype_size(datatype));
1050 *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
1051 retval = MPI_SUCCESS;
1054 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1055 (*request)->recv = 1;
1060 if(retval!=MPI_SUCCESS)*request=MPI_REQUEST_NULL;
1065 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
1066 int tag, MPI_Comm comm, MPI_Request * request)
1071 if (request == NULL) {
1072 retval = MPI_ERR_ARG;
1073 } else if (comm == MPI_COMM_NULL) {
1074 retval = MPI_ERR_COMM;
1075 } else if (dst == MPI_PROC_NULL) {
1076 *request = MPI_REQUEST_NULL;
1077 retval = MPI_SUCCESS;
1078 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1079 retval = MPI_ERR_RANK;
1080 } else if (count < 0) {
1081 retval = MPI_ERR_COUNT;
1082 } else if (buf==NULL && count > 0) {
1083 retval = MPI_ERR_COUNT;
1084 } else if (datatype == MPI_DATATYPE_NULL){
1085 retval = MPI_ERR_TYPE;
1086 } else if(tag<0 && tag != MPI_ANY_TAG){
1087 retval = MPI_ERR_TAG;
1091 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1092 TRACE_smpi_computing_out(rank);
1093 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1094 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
1095 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1098 *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1099 retval = MPI_SUCCESS;
1102 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1103 (*request)->send = 1;
1104 TRACE_smpi_computing_in(rank);
1109 if(retval!=MPI_SUCCESS)*request=MPI_REQUEST_NULL;
1113 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request* request) {
1117 if (request == NULL) {
1118 retval = MPI_ERR_ARG;
1119 } else if (comm == MPI_COMM_NULL) {
1120 retval = MPI_ERR_COMM;
1121 } else if (dst == MPI_PROC_NULL) {
1122 *request = MPI_REQUEST_NULL;
1123 retval = MPI_SUCCESS;
1124 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1125 retval = MPI_ERR_RANK;
1126 } else if (count < 0) {
1127 retval = MPI_ERR_COUNT;
1128 } else if (buf==NULL && count > 0) {
1129 retval = MPI_ERR_COUNT;
1130 } else if (datatype == MPI_DATATYPE_NULL){
1131 retval = MPI_ERR_TYPE;
1132 } else if(tag<0 && tag != MPI_ANY_TAG){
1133 retval = MPI_ERR_TAG;
1137 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1138 TRACE_smpi_computing_out(rank);
1139 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1140 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
1141 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1144 *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1145 retval = MPI_SUCCESS;
1148 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1149 (*request)->send = 1;
1150 TRACE_smpi_computing_in(rank);
1155 if(retval!=MPI_SUCCESS)*request=MPI_REQUEST_NULL;
1159 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1160 MPI_Comm comm, MPI_Status * status)
1165 if (comm == MPI_COMM_NULL) {
1166 retval = MPI_ERR_COMM;
1167 } else if (src == MPI_PROC_NULL) {
1168 smpi_empty_status(status);
1169 status->MPI_SOURCE = MPI_PROC_NULL;
1170 retval = MPI_SUCCESS;
1171 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1172 retval = MPI_ERR_RANK;
1173 } else if (count < 0) {
1174 retval = MPI_ERR_COUNT;
1175 } else if (buf==NULL && count > 0) {
1176 retval = MPI_ERR_COUNT;
1177 } else if (datatype == MPI_DATATYPE_NULL){
1178 retval = MPI_ERR_TYPE;
1179 } else if(tag<0 && tag != MPI_ANY_TAG){
1180 retval = MPI_ERR_TAG;
1183 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1184 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1185 TRACE_smpi_computing_out(rank);
1186 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, count*smpi_datatype_size(datatype));
1189 smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1190 retval = MPI_SUCCESS;
1193 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1194 if(status!=MPI_STATUS_IGNORE)src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1195 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1196 TRACE_smpi_recv(rank, src_traced, rank);
1197 TRACE_smpi_computing_in(rank);
1205 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1212 if (comm == MPI_COMM_NULL) {
1213 retval = MPI_ERR_COMM;
1214 } else if (dst == MPI_PROC_NULL) {
1215 retval = MPI_SUCCESS;
1216 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1217 retval = MPI_ERR_RANK;
1218 } else if (count < 0) {
1219 retval = MPI_ERR_COUNT;
1220 } else if (buf==NULL && count > 0) {
1221 retval = MPI_ERR_COUNT;
1222 } else if (datatype == MPI_DATATYPE_NULL){
1223 retval = MPI_ERR_TYPE;
1224 } else if(tag<0 && tag != MPI_ANY_TAG){
1225 retval = MPI_ERR_TAG;
1229 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1230 TRACE_smpi_computing_out(rank);
1231 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1232 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
1233 TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1236 smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1237 retval = MPI_SUCCESS;
1240 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1241 TRACE_smpi_computing_in(rank);
1251 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1256 if (comm == MPI_COMM_NULL) {
1257 retval = MPI_ERR_COMM;
1258 } else if (dst == MPI_PROC_NULL) {
1259 retval = MPI_SUCCESS;
1260 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1261 retval = MPI_ERR_RANK;
1262 } else if (count < 0) {
1263 retval = MPI_ERR_COUNT;
1264 } else if (buf==NULL && count > 0) {
1265 retval = MPI_ERR_COUNT;
1266 } else if (datatype == MPI_DATATYPE_NULL){
1267 retval = MPI_ERR_TYPE;
1268 } else if(tag<0 && tag != MPI_ANY_TAG){
1269 retval = MPI_ERR_TAG;
1273 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1274 TRACE_smpi_computing_out(rank);
1275 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1276 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
1277 TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1280 smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1281 retval = MPI_SUCCESS;
1284 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1285 TRACE_smpi_computing_in(rank);
1293 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1294 int dst, int sendtag, void *recvbuf, int recvcount,
1295 MPI_Datatype recvtype, int src, int recvtag,
1296 MPI_Comm comm, MPI_Status * status)
1302 if (comm == MPI_COMM_NULL) {
1303 retval = MPI_ERR_COMM;
1304 } else if (sendtype == MPI_DATATYPE_NULL
1305 || recvtype == MPI_DATATYPE_NULL) {
1306 retval = MPI_ERR_TYPE;
1307 } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1308 smpi_empty_status(status);
1309 status->MPI_SOURCE = MPI_PROC_NULL;
1310 retval = MPI_SUCCESS;
1311 }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1312 (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1313 retval = MPI_ERR_RANK;
1314 } else if (sendcount < 0 || recvcount<0) {
1315 retval = MPI_ERR_COUNT;
1316 } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1317 retval = MPI_ERR_COUNT;
1318 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1319 retval = MPI_ERR_TAG;
1323 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1324 TRACE_smpi_computing_out(rank);
1325 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1326 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1327 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, sendcount*smpi_datatype_size(sendtype));
1328 TRACE_smpi_send(rank, rank, dst_traced,sendcount*smpi_datatype_size(sendtype));
1332 smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1333 recvcount, recvtype, src, recvtag, comm, status);
1334 retval = MPI_SUCCESS;
1337 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1338 TRACE_smpi_recv(rank, src_traced, rank);
1339 TRACE_smpi_computing_in(rank);
1348 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1349 int dst, int sendtag, int src, int recvtag,
1350 MPI_Comm comm, MPI_Status * status)
1352 //TODO: suboptimal implementation
1355 if (datatype == MPI_DATATYPE_NULL) {
1356 retval = MPI_ERR_TYPE;
1357 } else if (count < 0) {
1358 retval = MPI_ERR_COUNT;
1360 int size = smpi_datatype_get_extent(datatype) * count;
1361 recvbuf = xbt_new(char, size);
1363 MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1364 datatype, src, recvtag, comm, status);
1365 if(retval==MPI_SUCCESS){
1366 smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1374 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1379 if (request == MPI_REQUEST_NULL || flag == NULL) {
1380 retval = MPI_ERR_ARG;
1381 } else if (*request == MPI_REQUEST_NULL) {
1383 retval = MPI_ERR_REQUEST;
1385 *flag = smpi_mpi_test(request, status);
1386 retval = MPI_SUCCESS;
1392 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1393 MPI_Status * status)
1398 if (index == NULL || flag == NULL) {
1399 retval = MPI_ERR_ARG;
1401 *flag = smpi_mpi_testany(count, requests, index, status);
1402 retval = MPI_SUCCESS;
1408 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1414 retval = MPI_ERR_ARG;
1416 *flag = smpi_mpi_testall(count, requests, statuses);
1417 retval = MPI_SUCCESS;
1423 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1427 if (status == NULL) {
1428 retval = MPI_ERR_ARG;
1429 } else if (comm == MPI_COMM_NULL) {
1430 retval = MPI_ERR_COMM;
1431 } else if (source == MPI_PROC_NULL) {
1432 smpi_empty_status(status);
1433 status->MPI_SOURCE = MPI_PROC_NULL;
1434 retval = MPI_SUCCESS;
1436 smpi_mpi_probe(source, tag, comm, status);
1437 retval = MPI_SUCCESS;
1444 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1449 retval = MPI_ERR_ARG;
1450 } else if (status == NULL) {
1451 retval = MPI_ERR_ARG;
1452 } else if (comm == MPI_COMM_NULL) {
1453 retval = MPI_ERR_COMM;
1454 } else if (source == MPI_PROC_NULL) {
1456 smpi_empty_status(status);
1457 status->MPI_SOURCE = MPI_PROC_NULL;
1458 retval = MPI_SUCCESS;
1460 smpi_mpi_iprobe(source, tag, comm, flag, status);
1461 retval = MPI_SUCCESS;
1467 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1473 if (request == NULL) {
1474 retval = MPI_ERR_ARG;
1475 } else if (*request == MPI_REQUEST_NULL) {
1476 retval = MPI_ERR_REQUEST;
1480 int rank = request && (*request)->comm != MPI_COMM_NULL
1481 ? smpi_process_index()
1483 TRACE_smpi_computing_out(rank);
1485 int src_traced = (*request)->src;
1486 int dst_traced = (*request)->dst;
1487 MPI_Comm comm = (*request)->comm;
1488 int is_wait_for_receive = (*request)->recv;
1489 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__,-1);
1492 smpi_mpi_wait(request, status);
1493 retval = MPI_SUCCESS;
1496 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1497 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1498 if (is_wait_for_receive) {
1499 if(src_traced==MPI_ANY_SOURCE)
1500 src_traced = (status!=MPI_STATUS_IGNORE) ?
1501 smpi_group_rank(smpi_comm_group(comm), status->MPI_SOURCE) :
1503 TRACE_smpi_recv(rank, src_traced, dst_traced);
1505 TRACE_smpi_computing_in(rank);
1514 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1520 //save requests information for tracing
1522 int *srcs = xbt_new(int, count);
1523 int *dsts = xbt_new(int, count);
1524 int *recvs = xbt_new(int, count);
1525 MPI_Comm *comms = xbt_new(MPI_Comm, count);
1527 for (i = 0; i < count; i++) {
1528 MPI_Request req = requests[i]; //already received requests are no longer valid
1532 recvs[i] = req->recv;
1533 comms[i] = req->comm;
1536 int rank_traced = smpi_process_index();
1537 TRACE_smpi_computing_out(rank_traced);
1539 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,count);
1542 if (index == NULL) {
1543 retval = MPI_ERR_ARG;
1545 *index = smpi_mpi_waitany(count, requests, status);
1546 retval = MPI_SUCCESS;
1549 if(*index!=MPI_UNDEFINED){
1550 int src_traced = srcs[*index];
1551 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1552 int dst_traced = dsts[*index];
1553 int is_wait_for_receive = recvs[*index];
1554 if (is_wait_for_receive) {
1555 if(srcs[*index]==MPI_ANY_SOURCE)
1556 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1557 smpi_group_rank(smpi_comm_group(comms[*index]), status->MPI_SOURCE) :
1559 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1561 TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1568 TRACE_smpi_computing_in(rank_traced);
1574 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1579 //save information from requests
1581 int *srcs = xbt_new(int, count);
1582 int *dsts = xbt_new(int, count);
1583 int *recvs = xbt_new(int, count);
1584 int *valid = xbt_new(int, count);
1585 MPI_Comm *comms = xbt_new(MPI_Comm, count);
1587 //int valid_count = 0;
1588 for (i = 0; i < count; i++) {
1589 MPI_Request req = requests[i];
1590 if(req!=MPI_REQUEST_NULL){
1593 recvs[i] = req->recv;
1594 comms[i] = req->comm;
1600 int rank_traced = smpi_process_index();
1601 TRACE_smpi_computing_out(rank_traced);
1603 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,count);
1605 int retval = smpi_mpi_waitall(count, requests, status);
1607 for (i = 0; i < count; i++) {
1609 //int src_traced = srcs[*index];
1610 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1611 int src_traced = srcs[i];
1612 int dst_traced = dsts[i];
1613 int is_wait_for_receive = recvs[i];
1614 if (is_wait_for_receive) {
1615 if(src_traced==MPI_ANY_SOURCE)
1616 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1617 smpi_group_rank(smpi_comm_group(comms[i]), status[i].MPI_SOURCE) :
1619 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1623 TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1630 TRACE_smpi_computing_in(rank_traced);
1636 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1637 int *indices, MPI_Status status[])
1642 if (outcount == NULL) {
1643 retval = MPI_ERR_ARG;
1645 *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1646 retval = MPI_SUCCESS;
1652 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1653 int* indices, MPI_Status status[])
1658 if (outcount == NULL) {
1659 retval = MPI_ERR_ARG;
1661 *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1662 retval = MPI_SUCCESS;
1669 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1675 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1676 TRACE_smpi_computing_out(rank);
1677 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1678 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
1680 if (comm == MPI_COMM_NULL) {
1681 retval = MPI_ERR_COMM;
1683 mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1684 retval = MPI_SUCCESS;
1687 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1688 TRACE_smpi_computing_in(rank);
1694 int PMPI_Barrier(MPI_Comm comm)
1700 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1701 TRACE_smpi_computing_out(rank);
1702 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, smpi_comm_size(comm));
1704 if (comm == MPI_COMM_NULL) {
1705 retval = MPI_ERR_COMM;
1707 mpi_coll_barrier_fun(comm);
1708 retval = MPI_SUCCESS;
1711 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1712 TRACE_smpi_computing_in(rank);
1718 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1719 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1720 int root, MPI_Comm comm)
1726 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1727 TRACE_smpi_computing_out(rank);
1728 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1729 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,sendcount*smpi_datatype_size(sendtype));
1731 if (comm == MPI_COMM_NULL) {
1732 retval = MPI_ERR_COMM;
1733 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1734 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1735 retval = MPI_ERR_TYPE;
1736 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1737 ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1738 retval = MPI_ERR_COUNT;
1741 char* sendtmpbuf = (char*) sendbuf;
1742 int sendtmpcount = sendcount;
1743 MPI_Datatype sendtmptype = sendtype;
1744 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1746 sendtmptype=recvtype;
1749 mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount,
1750 recvtype, root, comm);
1753 retval = MPI_SUCCESS;
1756 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1757 TRACE_smpi_computing_in(rank);
1763 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1764 void *recvbuf, int *recvcounts, int *displs,
1765 MPI_Datatype recvtype, int root, MPI_Comm comm)
1771 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1772 TRACE_smpi_computing_out(rank);
1773 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1774 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,sendcount*smpi_datatype_size(sendtype));
1776 if (comm == MPI_COMM_NULL) {
1777 retval = MPI_ERR_COMM;
1778 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1779 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1780 retval = MPI_ERR_TYPE;
1781 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1782 retval = MPI_ERR_COUNT;
1783 } else if (recvcounts == NULL || displs == NULL) {
1784 retval = MPI_ERR_ARG;
1787 char* sendtmpbuf = (char*) sendbuf;
1788 int sendtmpcount = sendcount;
1789 MPI_Datatype sendtmptype = sendtype;
1790 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1792 sendtmptype=recvtype;
1795 smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts,
1796 displs, recvtype, root, comm);
1797 retval = MPI_SUCCESS;
1800 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1801 TRACE_smpi_computing_in(rank);
1807 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1808 void *recvbuf, int recvcount, MPI_Datatype recvtype,
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__,sendcount*smpi_datatype_size(sendtype));
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)) ||
1826 retval = MPI_ERR_COUNT;
1829 if(sendbuf == MPI_IN_PLACE) {
1830 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1831 sendcount=recvcount;
1835 mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1837 retval = MPI_SUCCESS;
1840 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1846 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1847 void *recvbuf, int *recvcounts, int *displs,
1848 MPI_Datatype recvtype, MPI_Comm comm)
1854 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1855 TRACE_smpi_computing_out(rank);
1856 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,sendcount*smpi_datatype_size(sendtype));
1858 if (comm == MPI_COMM_NULL) {
1859 retval = MPI_ERR_COMM;
1860 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1861 (recvtype == MPI_DATATYPE_NULL)){
1862 retval = MPI_ERR_TYPE;
1863 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1864 retval = MPI_ERR_COUNT;
1865 } else if (recvcounts == NULL || displs == NULL) {
1866 retval = MPI_ERR_ARG;
1869 if(sendbuf == MPI_IN_PLACE) {
1870 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1871 sendcount=recvcounts[smpi_comm_rank(comm)];
1875 mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1876 displs, recvtype, comm);
1877 retval = MPI_SUCCESS;
1880 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1881 TRACE_smpi_computing_in(rank);
1887 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1888 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1889 int root, MPI_Comm comm)
1895 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1896 TRACE_smpi_computing_out(rank);
1897 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1899 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,sendcount*smpi_datatype_size(sendtype));
1901 if (comm == MPI_COMM_NULL) {
1902 retval = MPI_ERR_COMM;
1903 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1904 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1905 retval = MPI_ERR_TYPE;
1907 if (recvbuf == MPI_IN_PLACE) {
1909 recvcount=sendcount;
1911 mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1912 recvtype, root, comm);
1913 retval = MPI_SUCCESS;
1916 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1917 TRACE_smpi_computing_in(rank);
1923 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1924 MPI_Datatype sendtype, void *recvbuf, int recvcount,
1925 MPI_Datatype recvtype, int root, MPI_Comm comm)
1931 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1932 TRACE_smpi_computing_out(rank);
1933 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1935 for(i=0; i<smpi_comm_size(comm);i++)count+=sendcounts[i];
1936 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, count*smpi_datatype_size(sendtype));
1938 if (comm == MPI_COMM_NULL) {
1939 retval = MPI_ERR_COMM;
1940 } else if (sendcounts == NULL || displs == NULL) {
1941 retval = MPI_ERR_ARG;
1942 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1943 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1944 retval = MPI_ERR_TYPE;
1946 if (recvbuf == MPI_IN_PLACE) {
1948 recvcount=sendcounts[smpi_comm_rank(comm)];
1950 smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1951 recvcount, recvtype, root, comm);
1952 retval = MPI_SUCCESS;
1955 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1956 TRACE_smpi_computing_in(rank);
1962 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1963 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1969 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1970 TRACE_smpi_computing_out(rank);
1971 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1972 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
1974 if (comm == MPI_COMM_NULL) {
1975 retval = MPI_ERR_COMM;
1976 } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1977 retval = MPI_ERR_ARG;
1980 mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
1982 retval = MPI_SUCCESS;
1985 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1986 TRACE_smpi_computing_in(rank);
1992 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
1993 MPI_Datatype datatype, MPI_Op op){
1997 if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1998 retval = MPI_ERR_ARG;
2000 smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
2007 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
2008 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2014 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2015 TRACE_smpi_computing_out(rank);
2016 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, count*smpi_datatype_size(datatype));
2018 if (comm == MPI_COMM_NULL) {
2019 retval = MPI_ERR_COMM;
2020 } else if (datatype == MPI_DATATYPE_NULL) {
2021 retval = MPI_ERR_TYPE;
2022 } else if (op == MPI_OP_NULL) {
2023 retval = MPI_ERR_OP;
2026 char* sendtmpbuf = (char*) sendbuf;
2027 if( sendbuf == MPI_IN_PLACE ) {
2028 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
2029 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
2032 mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
2034 if( sendbuf == MPI_IN_PLACE ) {
2035 xbt_free(sendtmpbuf);
2038 retval = MPI_SUCCESS;
2042 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2043 TRACE_smpi_computing_in(rank);
2049 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
2050 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2056 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2057 TRACE_smpi_computing_out(rank);
2058 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, count*smpi_datatype_size(datatype));
2060 if (comm == MPI_COMM_NULL) {
2061 retval = MPI_ERR_COMM;
2062 } else if (datatype == MPI_DATATYPE_NULL) {
2063 retval = MPI_ERR_TYPE;
2064 } else if (op == MPI_OP_NULL) {
2065 retval = MPI_ERR_OP;
2067 smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2068 retval = MPI_SUCCESS;
2071 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2072 TRACE_smpi_computing_in(rank);
2078 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2079 MPI_Op op, MPI_Comm comm){
2084 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2085 TRACE_smpi_computing_out(rank);
2086 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, count*smpi_datatype_size(datatype));
2088 if (comm == MPI_COMM_NULL) {
2089 retval = MPI_ERR_COMM;
2090 } else if (datatype == MPI_DATATYPE_NULL) {
2091 retval = MPI_ERR_TYPE;
2092 } else if (op == MPI_OP_NULL) {
2093 retval = MPI_ERR_OP;
2095 smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2096 retval = MPI_SUCCESS;
2099 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2100 TRACE_smpi_computing_in(rank);
2106 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
2107 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2112 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2113 TRACE_smpi_computing_out(rank);
2115 for(i=0; i<smpi_comm_size(comm);i++)count+=recvcounts[i];
2116 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, count*smpi_datatype_size(datatype));
2118 if (comm == MPI_COMM_NULL) {
2119 retval = MPI_ERR_COMM;
2120 } else if (datatype == MPI_DATATYPE_NULL) {
2121 retval = MPI_ERR_TYPE;
2122 } else if (op == MPI_OP_NULL) {
2123 retval = MPI_ERR_OP;
2124 } else if (recvcounts == NULL) {
2125 retval = MPI_ERR_ARG;
2127 void* sendtmpbuf=sendbuf;
2128 if(sendbuf==MPI_IN_PLACE){
2132 mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts,
2133 datatype, op, comm);
2134 retval = MPI_SUCCESS;
2137 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2138 TRACE_smpi_computing_in(rank);
2144 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2145 MPI_Datatype datatype, MPI_Op op, 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__, recvcount*smpi_comm_size(comm)*smpi_datatype_size(datatype));
2154 if (comm == MPI_COMM_NULL) {
2155 retval = MPI_ERR_COMM;
2156 } else if (datatype == MPI_DATATYPE_NULL) {
2157 retval = MPI_ERR_TYPE;
2158 } else if (op == MPI_OP_NULL) {
2159 retval = MPI_ERR_OP;
2160 } else if (recvcount < 0) {
2161 retval = MPI_ERR_ARG;
2163 int count=smpi_comm_size(comm);
2164 int* recvcounts=(int*)xbt_malloc(count);
2165 for (i=0; i<count;i++)recvcounts[i]=recvcount;
2166 mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
2167 datatype, op, comm);
2168 xbt_free(recvcounts);
2169 retval = MPI_SUCCESS;
2172 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2173 TRACE_smpi_computing_in(rank);
2179 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2180 void *recvbuf, int recvcount, MPI_Datatype recvtype,
2187 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2188 TRACE_smpi_computing_out(rank);
2189 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, sendcount*smpi_datatype_size(sendtype));
2191 if (comm == MPI_COMM_NULL) {
2192 retval = MPI_ERR_COMM;
2193 } else if (sendtype == MPI_DATATYPE_NULL
2194 || recvtype == MPI_DATATYPE_NULL) {
2195 retval = MPI_ERR_TYPE;
2197 retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2200 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2201 TRACE_smpi_computing_in(rank);
2207 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2208 MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2209 int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2215 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2216 TRACE_smpi_computing_out(rank);
2218 for(i=0; i< smpi_comm_size(comm);i++)size+=sendcounts[i];
2219 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, size*smpi_datatype_size(sendtype));
2221 if (comm == MPI_COMM_NULL) {
2222 retval = MPI_ERR_COMM;
2223 } else if (sendtype == MPI_DATATYPE_NULL
2224 || recvtype == MPI_DATATYPE_NULL) {
2225 retval = MPI_ERR_TYPE;
2226 } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2227 || recvdisps == NULL) {
2228 retval = MPI_ERR_ARG;
2231 mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2232 recvbuf, recvcounts, recvdisps, recvtype,
2236 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2237 TRACE_smpi_computing_in(rank);
2244 int PMPI_Get_processor_name(char *name, int *resultlen)
2246 int retval = MPI_SUCCESS;
2249 strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2250 strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2251 strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2252 MPI_MAX_PROCESSOR_NAME - 1 );
2255 MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2261 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2263 int retval = MPI_SUCCESS;
2267 if (status == NULL || count == NULL) {
2268 retval = MPI_ERR_ARG;
2269 } else if (datatype == MPI_DATATYPE_NULL) {
2270 retval = MPI_ERR_TYPE;
2272 size = smpi_datatype_size(datatype);
2275 } else if (status->count % size != 0) {
2276 retval = MPI_UNDEFINED;
2278 *count = smpi_mpi_get_count(status, datatype);
2285 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2289 if (old_type == MPI_DATATYPE_NULL) {
2290 retval = MPI_ERR_TYPE;
2291 } else if (count<0){
2292 retval = MPI_ERR_COUNT;
2294 retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2300 int PMPI_Type_commit(MPI_Datatype* datatype) {
2304 if (datatype == MPI_DATATYPE_NULL) {
2305 retval = MPI_ERR_TYPE;
2307 smpi_datatype_commit(datatype);
2308 retval = MPI_SUCCESS;
2315 int PMPI_Type_vector(int count, int blocklen, int stride, 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 || blocklen<0){
2322 retval = MPI_ERR_COUNT;
2324 retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2330 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, 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 || blocklen<0){
2337 retval = MPI_ERR_COUNT;
2339 retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2345 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2346 return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2349 int PMPI_Type_indexed(int count, int* blocklens, int* 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_indexed(count, blocklens, indices, old_type, new_type);
2364 int PMPI_Type_create_indexed(int count, int* blocklens, int* 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 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2379 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2383 if (old_type == MPI_DATATYPE_NULL) {
2384 retval = MPI_ERR_TYPE;
2385 } else if (count<0){
2386 retval = MPI_ERR_COUNT;
2388 int* blocklens=(int*)xbt_malloc(blocklength*count);
2389 for (i=0; i<count;i++)blocklens[i]=blocklength;
2390 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2391 xbt_free(blocklens);
2398 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2402 if (old_type == MPI_DATATYPE_NULL) {
2403 retval = MPI_ERR_TYPE;
2404 } else if (count<0){
2405 retval = MPI_ERR_COUNT;
2407 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2413 int PMPI_Type_create_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2414 return PMPI_Type_hindexed(count, blocklens,indices,old_type,new_type);
2417 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2421 if (old_type == MPI_DATATYPE_NULL) {
2422 retval = MPI_ERR_TYPE;
2423 } else if (count<0){
2424 retval = MPI_ERR_COUNT;
2426 int* blocklens=(int*)xbt_malloc(blocklength*count);
2427 for (i=0; i<count;i++)blocklens[i]=blocklength;
2428 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2429 xbt_free(blocklens);
2436 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2441 retval = MPI_ERR_COUNT;
2443 retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2449 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2450 return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2454 int PMPI_Error_class(int errorcode, int* errorclass) {
2455 // assume smpi uses only standard mpi error codes
2456 *errorclass=errorcode;
2461 int PMPI_Initialized(int* flag) {
2462 *flag=smpi_process_initialized();
2466 /* The following calls are not yet implemented and will fail at runtime. */
2467 /* Once implemented, please move them above this notice. */
2469 #define NOT_YET_IMPLEMENTED {\
2470 XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__);\
2471 return MPI_SUCCESS;\
2475 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2479 int PMPI_Type_set_name(MPI_Datatype datatype, char * name)
2484 int PMPI_Type_get_name(MPI_Datatype datatype, char * name, int* len)
2489 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2493 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2497 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
2501 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2505 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2509 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2513 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2517 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2521 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2525 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2529 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2533 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2537 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2541 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2545 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2549 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2553 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2557 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2561 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2565 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2569 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2573 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2577 int PMPI_Comm_get_errhandler(MPI_Comm comm, MPI_Errhandler* errhandler) {
2581 int PMPI_Cancel(MPI_Request* request) {
2585 int PMPI_Buffer_attach(void* buffer, int size) {
2589 int PMPI_Buffer_detach(void* buffer, int* size) {
2593 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2597 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2602 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2607 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2612 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2617 int PMPI_Comm_free_keyval(int* keyval) {
2621 int PMPI_Pcontrol(const int level )
2626 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2630 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2635 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2640 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2645 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2650 int PMPI_Type_free_keyval(int* keyval) {
2654 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2658 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2662 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2666 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2670 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2674 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2678 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2682 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2686 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2690 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2694 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2698 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2702 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2706 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2710 int PMPI_Keyval_free(int* keyval) {
2714 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2718 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2722 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2726 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2730 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2734 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2738 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2742 int PMPI_Win_fence( int assert, MPI_Win win){
2746 int PMPI_Win_free( MPI_Win* win){
2750 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2754 int PMPI_Info_create( MPI_Info *info){
2758 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2762 int PMPI_Info_free( MPI_Info *info){
2766 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2767 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2771 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2772 int *num_addresses, int *num_datatypes, int *combiner){
2776 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2777 int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2778 MPI_Datatype* array_of_datatypes){
2782 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2783 int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2784 int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2788 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2792 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){
2796 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
2800 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
2801 void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
2806 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
2810 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
2814 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
2818 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
2822 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
2826 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
2830 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
2834 int PMPI_Add_error_class( int *errorclass){
2838 int PMPI_Add_error_code( int errorclass, int *errorcode){
2842 int PMPI_Add_error_string( int errorcode, char *string){
2846 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
2850 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
2854 int PMPI_Info_delete(MPI_Info info, char *key){
2858 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
2862 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
2866 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
2870 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2874 int MPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2878 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){
2882 int PMPI_Grequest_complete( MPI_Request request){
2886 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
2890 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
2894 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2898 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
2902 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
2906 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
2910 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
2914 int PMPI_Open_port( MPI_Info info, char *port_name){
2918 int PMPI_Close_port(char *port_name){
2922 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2926 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){
2930 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
2931 int* array_of_maxprocs, MPI_Info* array_of_info, int root,
2932 MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
2936 int PMPI_Comm_get_parent( MPI_Comm *parent){