1 /* Copyright (c) 2007-2013. The SimGrid Team.
2 * All rights reserved. */
4 /* This program is free software; you can redistribute it and/or modify it
5 * under the terms of the license (GNU LGPL) which comes with this package. */
8 #include "smpi_mpi_dt_private.h"
10 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_pmpi, smpi,
11 "Logging specific to SMPI (pmpi)");
14 //this function need to be here because of the calls to smpi_bench
15 void TRACE_smpi_set_category(const char *category)
17 //need to end bench otherwise categories for execution tasks are wrong
19 TRACE_internal_smpi_set_category (category);
20 //begin bench after changing process's category
25 /* PMPI User level calls */
27 int PMPI_Init(int *argc, char ***argv)
29 smpi_process_init(argc, argv);
30 smpi_process_mark_as_initialized();
32 int rank = smpi_process_index();
33 TRACE_smpi_init(rank);
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)
130 time = SIMIX_get_clock();
134 extern double sg_maxmin_precision;
135 double PMPI_Wtick(void)
137 return sg_maxmin_precision;
140 int PMPI_Address(void *location, MPI_Aint * address)
146 retval = MPI_ERR_ARG;
148 *address = (MPI_Aint) location;
149 retval = MPI_SUCCESS;
155 int PMPI_Get_address(void *location, MPI_Aint * address)
157 return PMPI_Address(location, address);
160 int PMPI_Type_free(MPI_Datatype * datatype)
166 retval = MPI_ERR_ARG;
168 smpi_datatype_free(datatype);
169 retval = MPI_SUCCESS;
175 int PMPI_Type_size(MPI_Datatype datatype, int *size)
180 if (datatype == MPI_DATATYPE_NULL) {
181 retval = MPI_ERR_TYPE;
182 } else if (size == NULL) {
183 retval = MPI_ERR_ARG;
185 *size = (int) smpi_datatype_size(datatype);
186 retval = MPI_SUCCESS;
192 int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
197 if (datatype == MPI_DATATYPE_NULL) {
198 retval = MPI_ERR_TYPE;
199 } else if (lb == NULL || extent == NULL) {
200 retval = MPI_ERR_ARG;
202 retval = smpi_datatype_extent(datatype, lb, extent);
208 int PMPI_Type_get_true_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
210 return PMPI_Type_get_extent(datatype, lb, extent);
213 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
218 if (datatype == MPI_DATATYPE_NULL) {
219 retval = MPI_ERR_TYPE;
220 } else if (extent == NULL) {
221 retval = MPI_ERR_ARG;
223 *extent = smpi_datatype_get_extent(datatype);
224 retval = MPI_SUCCESS;
230 int PMPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
235 if (datatype == MPI_DATATYPE_NULL) {
236 retval = MPI_ERR_TYPE;
237 } else if (disp == NULL) {
238 retval = MPI_ERR_ARG;
240 *disp = smpi_datatype_lb(datatype);
241 retval = MPI_SUCCESS;
247 int PMPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
252 if (datatype == MPI_DATATYPE_NULL) {
253 retval = MPI_ERR_TYPE;
254 } else if (disp == NULL) {
255 retval = MPI_ERR_ARG;
257 *disp = smpi_datatype_ub(datatype);
258 retval = MPI_SUCCESS;
264 int PMPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
269 if (function == NULL || op == NULL) {
270 retval = MPI_ERR_ARG;
272 *op = smpi_op_new(function, commute);
273 retval = MPI_SUCCESS;
279 int PMPI_Op_free(MPI_Op * op)
285 retval = MPI_ERR_ARG;
286 } else if (*op == MPI_OP_NULL) {
289 smpi_op_destroy(*op);
291 retval = MPI_SUCCESS;
297 int PMPI_Group_free(MPI_Group * group)
303 retval = MPI_ERR_ARG;
305 smpi_group_destroy(*group);
306 *group = MPI_GROUP_NULL;
307 retval = MPI_SUCCESS;
313 int PMPI_Group_size(MPI_Group group, int *size)
318 if (group == MPI_GROUP_NULL) {
319 retval = MPI_ERR_GROUP;
320 } else if (size == NULL) {
321 retval = MPI_ERR_ARG;
323 *size = smpi_group_size(group);
324 retval = MPI_SUCCESS;
330 int PMPI_Group_rank(MPI_Group group, int *rank)
335 if (group == MPI_GROUP_NULL) {
336 retval = MPI_ERR_GROUP;
337 } else if (rank == NULL) {
338 retval = MPI_ERR_ARG;
340 *rank = smpi_group_rank(group, smpi_process_index());
341 retval = MPI_SUCCESS;
347 int PMPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1,
348 MPI_Group group2, int *ranks2)
350 int retval, i, index;
352 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
353 retval = MPI_ERR_GROUP;
355 for (i = 0; i < n; i++) {
356 if(ranks1[i]==MPI_PROC_NULL){
357 ranks2[i]=MPI_PROC_NULL;
359 index = smpi_group_index(group1, ranks1[i]);
360 ranks2[i] = smpi_group_rank(group2, index);
363 retval = MPI_SUCCESS;
369 int PMPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
374 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
375 retval = MPI_ERR_GROUP;
376 } else if (result == NULL) {
377 retval = MPI_ERR_ARG;
379 *result = smpi_group_compare(group1, group2);
380 retval = MPI_SUCCESS;
386 int PMPI_Group_union(MPI_Group group1, MPI_Group group2,
387 MPI_Group * newgroup)
389 int retval, i, proc1, proc2, size, size2;
392 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
393 retval = MPI_ERR_GROUP;
394 } else if (newgroup == NULL) {
395 retval = MPI_ERR_ARG;
397 size = smpi_group_size(group1);
398 size2 = smpi_group_size(group2);
399 for (i = 0; i < size2; i++) {
400 proc2 = smpi_group_index(group2, i);
401 proc1 = smpi_group_rank(group1, proc2);
402 if (proc1 == MPI_UNDEFINED) {
407 *newgroup = MPI_GROUP_EMPTY;
409 *newgroup = smpi_group_new(size);
410 size2 = smpi_group_size(group1);
411 for (i = 0; i < size2; i++) {
412 proc1 = smpi_group_index(group1, i);
413 smpi_group_set_mapping(*newgroup, proc1, i);
415 for (i = size2; i < size; i++) {
416 proc2 = smpi_group_index(group2, i - size2);
417 smpi_group_set_mapping(*newgroup, proc2, i);
420 retval = MPI_SUCCESS;
426 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
427 MPI_Group * newgroup)
429 int retval, i, proc1, proc2, size;
432 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
433 retval = MPI_ERR_GROUP;
434 } else if (newgroup == NULL) {
435 retval = MPI_ERR_ARG;
437 size = smpi_group_size(group2);
438 for (i = 0; i < size; i++) {
439 proc2 = smpi_group_index(group2, i);
440 proc1 = smpi_group_rank(group1, proc2);
441 if (proc1 == MPI_UNDEFINED) {
446 *newgroup = MPI_GROUP_EMPTY;
448 *newgroup = smpi_group_new(size);
450 for (i = 0; i < smpi_group_size(group2); i++) {
451 proc2 = smpi_group_index(group2, i);
452 proc1 = smpi_group_rank(group1, proc2);
453 if (proc1 != MPI_UNDEFINED) {
454 smpi_group_set_mapping(*newgroup, proc2, j);
459 retval = MPI_SUCCESS;
465 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
467 int retval, i, proc1, proc2, size, size2;
470 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
471 retval = MPI_ERR_GROUP;
472 } else if (newgroup == NULL) {
473 retval = MPI_ERR_ARG;
475 size = size2 = smpi_group_size(group1);
476 for (i = 0; i < size2; i++) {
477 proc1 = smpi_group_index(group1, i);
478 proc2 = smpi_group_rank(group2, proc1);
479 if (proc2 != MPI_UNDEFINED) {
484 *newgroup = MPI_GROUP_EMPTY;
486 *newgroup = smpi_group_new(size);
487 for (i = 0; i < size2; i++) {
488 proc1 = smpi_group_index(group1, i);
489 proc2 = smpi_group_rank(group2, proc1);
490 if (proc2 == MPI_UNDEFINED) {
491 smpi_group_set_mapping(*newgroup, proc1, i);
495 retval = MPI_SUCCESS;
501 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
503 int retval, i, index;
506 if (group == MPI_GROUP_NULL) {
507 retval = MPI_ERR_GROUP;
508 } else if (newgroup == NULL) {
509 retval = MPI_ERR_ARG;
512 *newgroup = MPI_GROUP_EMPTY;
513 } else if (n == smpi_group_size(group)) {
515 if(group!= smpi_comm_group(MPI_COMM_WORLD)
516 && group != MPI_GROUP_NULL
517 && group != smpi_comm_group(MPI_COMM_SELF)
518 && group != MPI_GROUP_EMPTY)
519 smpi_group_use(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 retval = MPI_SUCCESS;
533 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
535 int retval, i, j, newsize, oldsize, index;
538 if (group == MPI_GROUP_NULL) {
539 retval = MPI_ERR_GROUP;
540 } else if (newgroup == NULL) {
541 retval = MPI_ERR_ARG;
545 if(group!= smpi_comm_group(MPI_COMM_WORLD)
546 && group != MPI_GROUP_NULL
547 && group != smpi_comm_group(MPI_COMM_SELF)
548 && group != MPI_GROUP_EMPTY)
549 smpi_group_use(group);
550 } else if (n == smpi_group_size(group)) {
551 *newgroup = MPI_GROUP_EMPTY;
553 oldsize=smpi_group_size(group);
554 newsize = oldsize - n;
555 *newgroup = smpi_group_new(newsize);
557 int* to_exclude=xbt_new(int, smpi_group_size(group));
558 for(i=0; i<oldsize; i++)
561 to_exclude[ranks[i]]=1;
564 for(i=0; i<oldsize; i++){
565 if(to_exclude[i]==0){
566 index = smpi_group_index(group, i);
567 smpi_group_set_mapping(*newgroup, index, j);
572 xbt_free(to_exclude);
574 retval = MPI_SUCCESS;
580 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
581 MPI_Group * newgroup)
583 int retval, i, j, rank, size, index;
586 if (group == MPI_GROUP_NULL) {
587 retval = MPI_ERR_GROUP;
588 } else if (newgroup == NULL) {
589 retval = MPI_ERR_ARG;
592 *newgroup = MPI_GROUP_EMPTY;
595 for (i = 0; i < n; i++) {
596 for (rank = ranges[i][0]; /* First */
597 rank >= 0; /* Last */
601 rank += ranges[i][2]; /* Stride */
602 if (ranges[i][0]<ranges[i][1]){
603 if(rank > ranges[i][1])
606 if(rank < ranges[i][1])
612 *newgroup = smpi_group_new(size);
614 for (i = 0; i < n; i++) {
615 for (rank = ranges[i][0]; /* First */
616 rank >= 0; /* Last */
618 index = smpi_group_index(group, rank);
619 smpi_group_set_mapping(*newgroup, index, j);
621 rank += ranges[i][2]; /* Stride */
622 if (ranges[i][0]<ranges[i][1]){
623 if(rank > ranges[i][1])
626 if(rank < ranges[i][1])
632 retval = MPI_SUCCESS;
638 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
639 MPI_Group * newgroup)
641 int retval, i, rank, newrank,oldrank, size, index, add;
644 if (group == MPI_GROUP_NULL) {
645 retval = MPI_ERR_GROUP;
646 } else if (newgroup == NULL) {
647 retval = MPI_ERR_ARG;
651 if(group!= smpi_comm_group(MPI_COMM_WORLD)
652 && group != MPI_GROUP_NULL
653 && group != smpi_comm_group(MPI_COMM_SELF)
654 && group != MPI_GROUP_EMPTY)
655 smpi_group_use(group);
657 size = smpi_group_size(group);
658 for (i = 0; i < n; i++) {
659 for (rank = ranges[i][0]; /* First */
660 rank >= 0; /* Last */
664 rank += ranges[i][2]; /* Stride */
665 if (ranges[i][0]<ranges[i][1]){
666 if(rank > ranges[i][1])
669 if(rank < ranges[i][1])
675 *newgroup = MPI_GROUP_EMPTY;
677 *newgroup = smpi_group_new(size);
680 while (newrank < size) {
682 for (i = 0; i < n; i++) {
683 for (rank = ranges[i][0];rank >= 0;){
689 rank += ranges[i][2]; /* Stride */
691 if (ranges[i][0]<ranges[i][1]){
692 if(rank > ranges[i][1])
695 if(rank < ranges[i][1])
701 index = smpi_group_index(group, oldrank);
702 smpi_group_set_mapping(*newgroup, index, newrank);
710 retval = MPI_SUCCESS;
716 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
721 if (comm == MPI_COMM_NULL) {
722 retval = MPI_ERR_COMM;
723 } else if (rank == NULL) {
724 retval = MPI_ERR_ARG;
726 *rank = smpi_comm_rank(comm);
727 retval = MPI_SUCCESS;
733 int PMPI_Comm_size(MPI_Comm comm, int *size)
738 if (comm == MPI_COMM_NULL) {
739 retval = MPI_ERR_COMM;
740 } else if (size == NULL) {
741 retval = MPI_ERR_ARG;
743 *size = smpi_comm_size(comm);
744 retval = MPI_SUCCESS;
750 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
755 if (comm == MPI_COMM_NULL) {
756 retval = MPI_ERR_COMM;
757 } else if (name == NULL || len == NULL) {
758 retval = MPI_ERR_ARG;
760 smpi_comm_get_name(comm, name, len);
761 retval = MPI_SUCCESS;
767 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
772 if (comm == MPI_COMM_NULL) {
773 retval = MPI_ERR_COMM;
774 } else if (group == NULL) {
775 retval = MPI_ERR_ARG;
777 *group = smpi_comm_group(comm);
778 if(*group!= smpi_comm_group(MPI_COMM_WORLD)
779 && *group != MPI_GROUP_NULL
780 && *group != smpi_comm_group(MPI_COMM_SELF)
781 && *group != MPI_GROUP_EMPTY)
782 smpi_group_use(*group);
783 retval = MPI_SUCCESS;
789 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
794 if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
795 retval = MPI_ERR_COMM;
796 } else if (result == NULL) {
797 retval = MPI_ERR_ARG;
799 if (comm1 == comm2) { /* Same communicators means same groups */
803 smpi_group_compare(smpi_comm_group(comm1),
804 smpi_comm_group(comm2));
805 if (*result == MPI_IDENT) {
806 *result = MPI_CONGRUENT;
809 retval = MPI_SUCCESS;
815 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
820 if (comm == MPI_COMM_NULL) {
821 retval = MPI_ERR_COMM;
822 } else if (newcomm == NULL) {
823 retval = MPI_ERR_ARG;
825 *newcomm = smpi_comm_new(smpi_comm_group(comm));
826 retval = MPI_SUCCESS;
832 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
837 if (comm == MPI_COMM_NULL) {
838 retval = MPI_ERR_COMM;
839 } else if (group == MPI_GROUP_NULL) {
840 retval = MPI_ERR_GROUP;
841 } else if (newcomm == NULL) {
842 retval = MPI_ERR_ARG;
843 } else if(smpi_group_rank(group,smpi_process_index())==MPI_UNDEFINED){
844 *newcomm= MPI_COMM_NULL;
845 retval = MPI_SUCCESS;
848 *newcomm = smpi_comm_new(group);
849 retval = MPI_SUCCESS;
855 int PMPI_Comm_free(MPI_Comm * comm)
861 retval = MPI_ERR_ARG;
862 } else if (*comm == MPI_COMM_NULL) {
863 retval = MPI_ERR_COMM;
865 smpi_comm_destroy(*comm);
866 *comm = MPI_COMM_NULL;
867 retval = MPI_SUCCESS;
873 int PMPI_Comm_disconnect(MPI_Comm * comm)
875 /* TODO: wait until all communication in comm are done */
880 retval = MPI_ERR_ARG;
881 } else if (*comm == MPI_COMM_NULL) {
882 retval = MPI_ERR_COMM;
884 smpi_comm_destroy(*comm);
885 *comm = MPI_COMM_NULL;
886 retval = MPI_SUCCESS;
892 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
897 if (comm_out == NULL) {
898 retval = MPI_ERR_ARG;
899 } else if (comm == MPI_COMM_NULL) {
900 retval = MPI_ERR_COMM;
902 *comm_out = smpi_comm_split(comm, color, key);
903 retval = MPI_SUCCESS;
909 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
910 int tag, MPI_Comm comm, MPI_Request * request)
915 if (request == NULL) {
916 retval = MPI_ERR_ARG;
917 } else if (comm == MPI_COMM_NULL) {
918 retval = MPI_ERR_COMM;
919 } else if (dst == MPI_PROC_NULL) {
920 retval = MPI_SUCCESS;
922 *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
923 retval = MPI_SUCCESS;
926 if (retval != MPI_SUCCESS && request)
927 *request = MPI_REQUEST_NULL;
931 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
932 int tag, MPI_Comm comm, MPI_Request * request)
937 if (request == NULL) {
938 retval = MPI_ERR_ARG;
939 } else if (comm == MPI_COMM_NULL) {
940 retval = MPI_ERR_COMM;
941 } else if (src == MPI_PROC_NULL) {
942 retval = MPI_SUCCESS;
944 *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
945 retval = MPI_SUCCESS;
948 if (retval != MPI_SUCCESS && request)
949 *request = MPI_REQUEST_NULL;
953 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype,
954 int dst, int tag, MPI_Comm comm, MPI_Request* request)
959 if (request == NULL) {
960 retval = MPI_ERR_ARG;
961 } else if (comm == MPI_COMM_NULL) {
962 retval = MPI_ERR_COMM;
963 } else if (dst == MPI_PROC_NULL) {
964 retval = MPI_SUCCESS;
966 *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
967 retval = MPI_SUCCESS;
970 if (retval != MPI_SUCCESS && request)
971 *request = MPI_REQUEST_NULL;
975 int PMPI_Start(MPI_Request * request)
980 if (request == NULL || *request == MPI_REQUEST_NULL) {
981 retval = MPI_ERR_ARG;
983 smpi_mpi_start(*request);
984 retval = MPI_SUCCESS;
990 int PMPI_Startall(int count, MPI_Request * requests)
995 if (requests == NULL) {
996 retval = MPI_ERR_ARG;
998 smpi_mpi_startall(count, requests);
999 retval = MPI_SUCCESS;
1005 int PMPI_Request_free(MPI_Request * request)
1010 if (*request == MPI_REQUEST_NULL) {
1011 retval = MPI_ERR_ARG;
1013 if((*request)->flags & PERSISTENT)(*request)->refcount--;
1014 smpi_mpi_request_free(request);
1015 retval = MPI_SUCCESS;
1021 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
1022 int tag, MPI_Comm comm, MPI_Request * request)
1028 if (request == NULL) {
1029 retval = MPI_ERR_ARG;
1030 } else if (comm == MPI_COMM_NULL) {
1031 retval = MPI_ERR_COMM;
1032 } else if (src == MPI_PROC_NULL) {
1033 *request = MPI_REQUEST_NULL;
1034 retval = MPI_SUCCESS;
1035 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1036 retval = MPI_ERR_COMM;
1037 } else if (count < 0) {
1038 retval = MPI_ERR_COUNT;
1039 } else if (buf==NULL && count > 0) {
1040 retval = MPI_ERR_COUNT;
1041 } else if (datatype == MPI_DATATYPE_NULL){
1042 retval = MPI_ERR_TYPE;
1043 } else if(tag<0 && tag != MPI_ANY_TAG){
1044 retval = MPI_ERR_TAG;
1048 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1049 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1050 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, count*smpi_datatype_size(datatype));
1053 *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
1054 retval = MPI_SUCCESS;
1057 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1058 (*request)->recv = 1;
1063 if (retval != MPI_SUCCESS && request)
1064 *request = MPI_REQUEST_NULL;
1069 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
1070 int tag, MPI_Comm comm, MPI_Request * request)
1075 if (request == NULL) {
1076 retval = MPI_ERR_ARG;
1077 } else if (comm == MPI_COMM_NULL) {
1078 retval = MPI_ERR_COMM;
1079 } else if (dst == MPI_PROC_NULL) {
1080 *request = MPI_REQUEST_NULL;
1081 retval = MPI_SUCCESS;
1082 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1083 retval = MPI_ERR_RANK;
1084 } else if (count < 0) {
1085 retval = MPI_ERR_COUNT;
1086 } else if (buf==NULL && count > 0) {
1087 retval = MPI_ERR_COUNT;
1088 } else if (datatype == MPI_DATATYPE_NULL){
1089 retval = MPI_ERR_TYPE;
1090 } else if(tag<0 && tag != MPI_ANY_TAG){
1091 retval = MPI_ERR_TAG;
1095 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1096 TRACE_smpi_computing_out(rank);
1097 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1098 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
1099 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1102 *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1103 retval = MPI_SUCCESS;
1106 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1107 (*request)->send = 1;
1108 TRACE_smpi_computing_in(rank);
1113 if (retval != MPI_SUCCESS && request)
1114 *request = MPI_REQUEST_NULL;
1118 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype,
1119 int dst, int tag, MPI_Comm comm, MPI_Request* request)
1124 if (request == NULL) {
1125 retval = MPI_ERR_ARG;
1126 } else if (comm == MPI_COMM_NULL) {
1127 retval = MPI_ERR_COMM;
1128 } else if (dst == MPI_PROC_NULL) {
1129 *request = MPI_REQUEST_NULL;
1130 retval = MPI_SUCCESS;
1131 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1132 retval = MPI_ERR_RANK;
1133 } else if (count < 0) {
1134 retval = MPI_ERR_COUNT;
1135 } else if (buf==NULL && count > 0) {
1136 retval = MPI_ERR_COUNT;
1137 } else if (datatype == MPI_DATATYPE_NULL){
1138 retval = MPI_ERR_TYPE;
1139 } else if(tag<0 && tag != MPI_ANY_TAG){
1140 retval = MPI_ERR_TAG;
1144 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1145 TRACE_smpi_computing_out(rank);
1146 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1147 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
1148 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1151 *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1152 retval = MPI_SUCCESS;
1155 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1156 (*request)->send = 1;
1157 TRACE_smpi_computing_in(rank);
1162 if (retval != MPI_SUCCESS && request)
1163 *request = MPI_REQUEST_NULL;
1167 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1168 MPI_Comm comm, MPI_Status * status)
1173 if (comm == MPI_COMM_NULL) {
1174 retval = MPI_ERR_COMM;
1175 } else if (src == MPI_PROC_NULL) {
1176 smpi_empty_status(status);
1177 status->MPI_SOURCE = MPI_PROC_NULL;
1178 retval = MPI_SUCCESS;
1179 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1180 retval = MPI_ERR_RANK;
1181 } else if (count < 0) {
1182 retval = MPI_ERR_COUNT;
1183 } else if (buf==NULL && count > 0) {
1184 retval = MPI_ERR_COUNT;
1185 } else if (datatype == MPI_DATATYPE_NULL){
1186 retval = MPI_ERR_TYPE;
1187 } else if(tag<0 && tag != MPI_ANY_TAG){
1188 retval = MPI_ERR_TAG;
1191 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1192 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1193 TRACE_smpi_computing_out(rank);
1194 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, count*smpi_datatype_size(datatype));
1197 smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1198 retval = MPI_SUCCESS;
1201 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1202 if(status!=MPI_STATUS_IGNORE)src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1203 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1204 TRACE_smpi_recv(rank, src_traced, rank);
1205 TRACE_smpi_computing_in(rank);
1213 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1220 if (comm == MPI_COMM_NULL) {
1221 retval = MPI_ERR_COMM;
1222 } else if (dst == MPI_PROC_NULL) {
1223 retval = MPI_SUCCESS;
1224 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1225 retval = MPI_ERR_RANK;
1226 } else if (count < 0) {
1227 retval = MPI_ERR_COUNT;
1228 } else if (buf==NULL && count > 0) {
1229 retval = MPI_ERR_COUNT;
1230 } else if (datatype == MPI_DATATYPE_NULL){
1231 retval = MPI_ERR_TYPE;
1232 } else if(tag<0 && tag != MPI_ANY_TAG){
1233 retval = MPI_ERR_TAG;
1237 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1238 TRACE_smpi_computing_out(rank);
1239 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1240 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
1241 TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1244 smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1245 retval = MPI_SUCCESS;
1248 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1249 TRACE_smpi_computing_in(rank);
1259 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1264 if (comm == MPI_COMM_NULL) {
1265 retval = MPI_ERR_COMM;
1266 } else if (dst == MPI_PROC_NULL) {
1267 retval = MPI_SUCCESS;
1268 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1269 retval = MPI_ERR_RANK;
1270 } else if (count < 0) {
1271 retval = MPI_ERR_COUNT;
1272 } else if (buf==NULL && count > 0) {
1273 retval = MPI_ERR_COUNT;
1274 } else if (datatype == MPI_DATATYPE_NULL){
1275 retval = MPI_ERR_TYPE;
1276 } else if(tag<0 && tag != MPI_ANY_TAG){
1277 retval = MPI_ERR_TAG;
1281 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1282 TRACE_smpi_computing_out(rank);
1283 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1284 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
1285 TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1288 smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1289 retval = MPI_SUCCESS;
1292 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1293 TRACE_smpi_computing_in(rank);
1301 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1302 int dst, int sendtag, void *recvbuf, int recvcount,
1303 MPI_Datatype recvtype, int src, int recvtag,
1304 MPI_Comm comm, MPI_Status * status)
1310 if (comm == MPI_COMM_NULL) {
1311 retval = MPI_ERR_COMM;
1312 } else if (sendtype == MPI_DATATYPE_NULL
1313 || recvtype == MPI_DATATYPE_NULL) {
1314 retval = MPI_ERR_TYPE;
1315 } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1316 smpi_empty_status(status);
1317 status->MPI_SOURCE = MPI_PROC_NULL;
1318 retval = MPI_SUCCESS;
1319 }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1320 (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1321 retval = MPI_ERR_RANK;
1322 } else if (sendcount < 0 || recvcount<0) {
1323 retval = MPI_ERR_COUNT;
1324 } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1325 retval = MPI_ERR_COUNT;
1326 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1327 retval = MPI_ERR_TAG;
1331 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1332 TRACE_smpi_computing_out(rank);
1333 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1334 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1335 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, sendcount*smpi_datatype_size(sendtype));
1336 TRACE_smpi_send(rank, rank, dst_traced,sendcount*smpi_datatype_size(sendtype));
1340 smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1341 recvcount, recvtype, src, recvtag, comm, status);
1342 retval = MPI_SUCCESS;
1345 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1346 TRACE_smpi_recv(rank, src_traced, rank);
1347 TRACE_smpi_computing_in(rank);
1356 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1357 int dst, int sendtag, int src, int recvtag,
1358 MPI_Comm comm, MPI_Status * status)
1360 //TODO: suboptimal implementation
1363 if (datatype == MPI_DATATYPE_NULL) {
1364 retval = MPI_ERR_TYPE;
1365 } else if (count < 0) {
1366 retval = MPI_ERR_COUNT;
1368 int size = smpi_datatype_get_extent(datatype) * count;
1369 recvbuf = xbt_new(char, size);
1371 MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1372 datatype, src, recvtag, comm, status);
1373 if(retval==MPI_SUCCESS){
1374 smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1382 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1387 if (request == NULL || flag == NULL) {
1388 retval = MPI_ERR_ARG;
1389 } else if (*request == MPI_REQUEST_NULL) {
1391 smpi_empty_status(status);
1392 retval = MPI_ERR_REQUEST;
1394 *flag = smpi_mpi_test(request, status);
1395 retval = MPI_SUCCESS;
1401 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1402 MPI_Status * status)
1407 if (index == NULL || flag == NULL) {
1408 retval = MPI_ERR_ARG;
1410 *flag = smpi_mpi_testany(count, requests, index, status);
1411 retval = MPI_SUCCESS;
1417 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1423 retval = MPI_ERR_ARG;
1425 *flag = smpi_mpi_testall(count, requests, statuses);
1426 retval = MPI_SUCCESS;
1432 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1436 if (status == NULL) {
1437 retval = MPI_ERR_ARG;
1438 } else if (comm == MPI_COMM_NULL) {
1439 retval = MPI_ERR_COMM;
1440 } else if (source == MPI_PROC_NULL) {
1441 smpi_empty_status(status);
1442 status->MPI_SOURCE = MPI_PROC_NULL;
1443 retval = MPI_SUCCESS;
1445 smpi_mpi_probe(source, tag, comm, status);
1446 retval = MPI_SUCCESS;
1453 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1458 retval = MPI_ERR_ARG;
1459 } else if (status == NULL) {
1460 retval = MPI_ERR_ARG;
1461 } else if (comm == MPI_COMM_NULL) {
1462 retval = MPI_ERR_COMM;
1463 } else if (source == MPI_PROC_NULL) {
1465 smpi_empty_status(status);
1466 status->MPI_SOURCE = MPI_PROC_NULL;
1467 retval = MPI_SUCCESS;
1469 smpi_mpi_iprobe(source, tag, comm, flag, status);
1470 retval = MPI_SUCCESS;
1476 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1482 smpi_empty_status(status);
1484 if (request == NULL) {
1485 retval = MPI_ERR_ARG;
1486 } else if (*request == MPI_REQUEST_NULL) {
1487 retval = MPI_ERR_REQUEST;
1491 int rank = request && (*request)->comm != MPI_COMM_NULL
1492 ? smpi_process_index()
1494 TRACE_smpi_computing_out(rank);
1496 int src_traced = (*request)->src;
1497 int dst_traced = (*request)->dst;
1498 MPI_Comm comm = (*request)->comm;
1499 int is_wait_for_receive = (*request)->recv;
1500 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__,-1);
1503 smpi_mpi_wait(request, status);
1504 retval = MPI_SUCCESS;
1507 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1508 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1509 if (is_wait_for_receive) {
1510 if(src_traced==MPI_ANY_SOURCE)
1511 src_traced = (status!=MPI_STATUS_IGNORE) ?
1512 smpi_group_rank(smpi_comm_group(comm), status->MPI_SOURCE) :
1514 TRACE_smpi_recv(rank, src_traced, dst_traced);
1516 TRACE_smpi_computing_in(rank);
1525 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1531 //save requests information for tracing
1533 int *srcs = xbt_new(int, count);
1534 int *dsts = xbt_new(int, count);
1535 int *recvs = xbt_new(int, count);
1536 MPI_Comm *comms = xbt_new(MPI_Comm, count);
1538 for (i = 0; i < count; i++) {
1539 MPI_Request req = requests[i]; //already received requests are no longer valid
1543 recvs[i] = req->recv;
1544 comms[i] = req->comm;
1547 int rank_traced = smpi_process_index();
1548 TRACE_smpi_computing_out(rank_traced);
1550 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,count);
1553 if (index == NULL) {
1554 retval = MPI_ERR_ARG;
1556 *index = smpi_mpi_waitany(count, requests, status);
1557 retval = MPI_SUCCESS;
1560 if(*index!=MPI_UNDEFINED){
1561 int src_traced = srcs[*index];
1562 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1563 int dst_traced = dsts[*index];
1564 int is_wait_for_receive = recvs[*index];
1565 if (is_wait_for_receive) {
1566 if(srcs[*index]==MPI_ANY_SOURCE)
1567 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1568 smpi_group_rank(smpi_comm_group(comms[*index]), status->MPI_SOURCE) :
1570 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1572 TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1579 TRACE_smpi_computing_in(rank_traced);
1585 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1590 //save information from requests
1592 int *srcs = xbt_new(int, count);
1593 int *dsts = xbt_new(int, count);
1594 int *recvs = xbt_new(int, count);
1595 int *valid = xbt_new(int, count);
1596 MPI_Comm *comms = xbt_new(MPI_Comm, count);
1598 //int valid_count = 0;
1599 for (i = 0; i < count; i++) {
1600 MPI_Request req = requests[i];
1601 if(req!=MPI_REQUEST_NULL){
1604 recvs[i] = req->recv;
1605 comms[i] = req->comm;
1611 int rank_traced = smpi_process_index();
1612 TRACE_smpi_computing_out(rank_traced);
1614 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,count);
1616 int retval = smpi_mpi_waitall(count, requests, status);
1618 for (i = 0; i < count; i++) {
1620 //int src_traced = srcs[*index];
1621 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1622 int src_traced = srcs[i];
1623 int dst_traced = dsts[i];
1624 int is_wait_for_receive = recvs[i];
1625 if (is_wait_for_receive) {
1626 if(src_traced==MPI_ANY_SOURCE)
1627 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1628 smpi_group_rank(smpi_comm_group(comms[i]), status[i].MPI_SOURCE) :
1630 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1634 TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1641 TRACE_smpi_computing_in(rank_traced);
1647 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1648 int *indices, MPI_Status status[])
1653 if (outcount == NULL) {
1654 retval = MPI_ERR_ARG;
1656 *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1657 retval = MPI_SUCCESS;
1663 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1664 int* indices, MPI_Status status[])
1669 if (outcount == NULL) {
1670 retval = MPI_ERR_ARG;
1672 *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1673 retval = MPI_SUCCESS;
1680 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1686 if (comm == MPI_COMM_NULL) {
1687 retval = MPI_ERR_COMM;
1690 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1691 TRACE_smpi_computing_out(rank);
1692 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1693 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
1695 mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1696 retval = MPI_SUCCESS;
1698 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1699 TRACE_smpi_computing_in(rank);
1707 int PMPI_Barrier(MPI_Comm comm)
1713 if (comm == MPI_COMM_NULL) {
1714 retval = MPI_ERR_COMM;
1717 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1718 TRACE_smpi_computing_out(rank);
1719 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, smpi_comm_size(comm));
1721 mpi_coll_barrier_fun(comm);
1722 retval = MPI_SUCCESS;
1724 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1725 TRACE_smpi_computing_in(rank);
1733 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1734 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1735 int root, MPI_Comm comm)
1741 if (comm == MPI_COMM_NULL) {
1742 retval = MPI_ERR_COMM;
1743 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1744 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1745 retval = MPI_ERR_TYPE;
1746 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1747 ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1748 retval = MPI_ERR_COUNT;
1751 char* sendtmpbuf = (char*) sendbuf;
1752 int sendtmpcount = sendcount;
1753 MPI_Datatype sendtmptype = sendtype;
1754 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1756 sendtmptype=recvtype;
1759 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1760 TRACE_smpi_computing_out(rank);
1761 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1762 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,sendcount*smpi_datatype_size(sendtmptype));
1764 mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount,
1765 recvtype, root, comm);
1768 retval = MPI_SUCCESS;
1770 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1771 TRACE_smpi_computing_in(rank);
1779 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1780 void *recvbuf, int *recvcounts, int *displs,
1781 MPI_Datatype recvtype, int root, MPI_Comm comm)
1787 if (comm == MPI_COMM_NULL) {
1788 retval = MPI_ERR_COMM;
1789 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1790 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1791 retval = MPI_ERR_TYPE;
1792 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1793 retval = MPI_ERR_COUNT;
1794 } else if (recvcounts == NULL || displs == NULL) {
1795 retval = MPI_ERR_ARG;
1797 char* sendtmpbuf = (char*) sendbuf;
1798 int sendtmpcount = sendcount;
1799 MPI_Datatype sendtmptype = sendtype;
1800 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1802 sendtmptype=recvtype;
1806 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1807 TRACE_smpi_computing_out(rank);
1808 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1809 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,sendcount*smpi_datatype_size(sendtmptype));
1811 smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts,
1812 displs, recvtype, root, comm);
1813 retval = MPI_SUCCESS;
1815 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1816 TRACE_smpi_computing_in(rank);
1824 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1825 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1832 if (comm == MPI_COMM_NULL) {
1833 retval = MPI_ERR_COMM;
1834 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1835 (recvtype == MPI_DATATYPE_NULL)){
1836 retval = MPI_ERR_TYPE;
1837 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1839 retval = MPI_ERR_COUNT;
1841 if(sendbuf == MPI_IN_PLACE) {
1842 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1843 sendcount=recvcount;
1847 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1848 TRACE_smpi_computing_out(rank);
1849 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,sendcount*smpi_datatype_size(sendtype));
1851 mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1853 retval = MPI_SUCCESS;
1856 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1863 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1864 void *recvbuf, int *recvcounts, int *displs,
1865 MPI_Datatype recvtype, MPI_Comm comm)
1871 if (comm == MPI_COMM_NULL) {
1872 retval = MPI_ERR_COMM;
1873 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1874 (recvtype == MPI_DATATYPE_NULL)){
1875 retval = MPI_ERR_TYPE;
1876 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1877 retval = MPI_ERR_COUNT;
1878 } else if (recvcounts == NULL || displs == NULL) {
1879 retval = MPI_ERR_ARG;
1882 if(sendbuf == MPI_IN_PLACE) {
1883 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1884 sendcount=recvcounts[smpi_comm_rank(comm)];
1888 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1889 TRACE_smpi_computing_out(rank);
1890 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,sendcount*smpi_datatype_size(sendtype));
1892 mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1893 displs, recvtype, comm);
1894 retval = MPI_SUCCESS;
1896 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1897 TRACE_smpi_computing_in(rank);
1905 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1906 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1907 int root, MPI_Comm comm)
1913 if (comm == MPI_COMM_NULL) {
1914 retval = MPI_ERR_COMM;
1915 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1916 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1917 retval = MPI_ERR_TYPE;
1920 if (recvbuf == MPI_IN_PLACE) {
1922 recvcount=sendcount;
1925 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1926 TRACE_smpi_computing_out(rank);
1927 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1929 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,sendcount*smpi_datatype_size(recvtype));
1931 mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1932 recvtype, root, comm);
1933 retval = MPI_SUCCESS;
1935 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1936 TRACE_smpi_computing_in(rank);
1944 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1945 MPI_Datatype sendtype, void *recvbuf, int recvcount,
1946 MPI_Datatype recvtype, int root, MPI_Comm comm)
1952 if (comm == MPI_COMM_NULL) {
1953 retval = MPI_ERR_COMM;
1954 } else if (sendcounts == NULL || displs == NULL) {
1955 retval = MPI_ERR_ARG;
1956 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1957 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1958 retval = MPI_ERR_TYPE;
1960 if (recvbuf == MPI_IN_PLACE) {
1962 recvcount=sendcounts[smpi_comm_rank(comm)];
1965 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1966 TRACE_smpi_computing_out(rank);
1967 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1969 for(i=0; i<smpi_comm_size(comm);i++)count+=sendcounts[i];
1970 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, count*smpi_datatype_size(sendtype));
1972 smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
1973 recvcount, recvtype, root, comm);
1974 retval = MPI_SUCCESS;
1976 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1977 TRACE_smpi_computing_in(rank);
1985 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
1986 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
1992 if (comm == MPI_COMM_NULL) {
1993 retval = MPI_ERR_COMM;
1994 } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
1995 retval = MPI_ERR_ARG;
1998 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1999 TRACE_smpi_computing_out(rank);
2000 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2001 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, count*smpi_datatype_size(datatype));
2003 mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
2005 retval = MPI_SUCCESS;
2007 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2008 TRACE_smpi_computing_in(rank);
2016 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
2017 MPI_Datatype datatype, MPI_Op op){
2021 if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
2022 retval = MPI_ERR_ARG;
2024 smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
2031 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
2032 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2038 if (comm == MPI_COMM_NULL) {
2039 retval = MPI_ERR_COMM;
2040 } else if (datatype == MPI_DATATYPE_NULL) {
2041 retval = MPI_ERR_TYPE;
2042 } else if (op == MPI_OP_NULL) {
2043 retval = MPI_ERR_OP;
2046 char* sendtmpbuf = (char*) sendbuf;
2047 if( sendbuf == MPI_IN_PLACE ) {
2048 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
2049 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
2052 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2053 TRACE_smpi_computing_out(rank);
2054 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, count*smpi_datatype_size(datatype));
2056 mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
2058 if( sendbuf == MPI_IN_PLACE ) {
2059 xbt_free(sendtmpbuf);
2062 retval = MPI_SUCCESS;
2064 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2065 TRACE_smpi_computing_in(rank);
2073 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
2074 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2080 if (comm == MPI_COMM_NULL) {
2081 retval = MPI_ERR_COMM;
2082 } else if (datatype == MPI_DATATYPE_NULL) {
2083 retval = MPI_ERR_TYPE;
2084 } else if (op == MPI_OP_NULL) {
2085 retval = MPI_ERR_OP;
2088 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2089 TRACE_smpi_computing_out(rank);
2090 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, count*smpi_datatype_size(datatype));
2092 smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2093 retval = MPI_SUCCESS;
2095 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2096 TRACE_smpi_computing_in(rank);
2104 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2105 MPI_Op op, MPI_Comm comm){
2110 if (comm == MPI_COMM_NULL) {
2111 retval = MPI_ERR_COMM;
2112 } else if (datatype == MPI_DATATYPE_NULL) {
2113 retval = MPI_ERR_TYPE;
2114 } else if (op == MPI_OP_NULL) {
2115 retval = MPI_ERR_OP;
2118 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2119 TRACE_smpi_computing_out(rank);
2120 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, count*smpi_datatype_size(datatype));
2122 smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2123 retval = MPI_SUCCESS;
2125 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2126 TRACE_smpi_computing_in(rank);
2134 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
2135 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2140 if (comm == MPI_COMM_NULL) {
2141 retval = MPI_ERR_COMM;
2142 } else if (datatype == MPI_DATATYPE_NULL) {
2143 retval = MPI_ERR_TYPE;
2144 } else if (op == MPI_OP_NULL) {
2145 retval = MPI_ERR_OP;
2146 } else if (recvcounts == NULL) {
2147 retval = MPI_ERR_ARG;
2150 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2151 TRACE_smpi_computing_out(rank);
2153 for(i=0; i<smpi_comm_size(comm);i++)count+=recvcounts[i];
2154 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, count*smpi_datatype_size(datatype));
2156 void* sendtmpbuf=sendbuf;
2157 if(sendbuf==MPI_IN_PLACE){
2161 mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts,
2162 datatype, op, comm);
2163 retval = MPI_SUCCESS;
2165 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2166 TRACE_smpi_computing_in(rank);
2174 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2175 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2180 if (comm == MPI_COMM_NULL) {
2181 retval = MPI_ERR_COMM;
2182 } else if (datatype == MPI_DATATYPE_NULL) {
2183 retval = MPI_ERR_TYPE;
2184 } else if (op == MPI_OP_NULL) {
2185 retval = MPI_ERR_OP;
2186 } else if (recvcount < 0) {
2187 retval = MPI_ERR_ARG;
2190 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2191 TRACE_smpi_computing_out(rank);
2192 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, recvcount*smpi_comm_size(comm)*smpi_datatype_size(datatype));
2194 int count=smpi_comm_size(comm);
2195 int* recvcounts=(int*)xbt_malloc(count);
2196 for (i=0; i<count;i++)recvcounts[i]=recvcount;
2197 mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
2198 datatype, op, comm);
2199 xbt_free(recvcounts);
2200 retval = MPI_SUCCESS;
2202 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2203 TRACE_smpi_computing_in(rank);
2211 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2212 void *recvbuf, int recvcount, MPI_Datatype recvtype,
2219 if (comm == MPI_COMM_NULL) {
2220 retval = MPI_ERR_COMM;
2221 } else if (sendtype == MPI_DATATYPE_NULL
2222 || recvtype == MPI_DATATYPE_NULL) {
2223 retval = MPI_ERR_TYPE;
2226 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2227 TRACE_smpi_computing_out(rank);
2228 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, sendcount*smpi_datatype_size(sendtype));
2230 retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2232 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2233 TRACE_smpi_computing_in(rank);
2241 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2242 MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2243 int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2249 if (comm == MPI_COMM_NULL) {
2250 retval = MPI_ERR_COMM;
2251 } else if (sendtype == MPI_DATATYPE_NULL
2252 || recvtype == MPI_DATATYPE_NULL) {
2253 retval = MPI_ERR_TYPE;
2254 } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2255 || recvdisps == NULL) {
2256 retval = MPI_ERR_ARG;
2259 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2260 TRACE_smpi_computing_out(rank);
2262 for(i=0; i< smpi_comm_size(comm);i++)size+=sendcounts[i];
2263 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, size*smpi_datatype_size(sendtype));
2266 mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2267 recvbuf, recvcounts, recvdisps, recvtype,
2270 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2271 TRACE_smpi_computing_in(rank);
2280 int PMPI_Get_processor_name(char *name, int *resultlen)
2282 int retval = MPI_SUCCESS;
2285 strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2286 strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2287 strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2288 MPI_MAX_PROCESSOR_NAME - 1 );
2291 MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2297 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2299 int retval = MPI_SUCCESS;
2303 if (status == NULL || count == NULL) {
2304 retval = MPI_ERR_ARG;
2305 } else if (datatype == MPI_DATATYPE_NULL) {
2306 retval = MPI_ERR_TYPE;
2308 size = smpi_datatype_size(datatype);
2311 } else if (status->count % size != 0) {
2312 retval = MPI_UNDEFINED;
2314 *count = smpi_mpi_get_count(status, datatype);
2321 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2325 if (old_type == MPI_DATATYPE_NULL) {
2326 retval = MPI_ERR_TYPE;
2327 } else if (count<0){
2328 retval = MPI_ERR_COUNT;
2330 retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2336 int PMPI_Type_commit(MPI_Datatype* datatype) {
2340 if (datatype == NULL || *datatype == MPI_DATATYPE_NULL) {
2341 retval = MPI_ERR_TYPE;
2343 smpi_datatype_commit(datatype);
2344 retval = MPI_SUCCESS;
2351 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2355 if (old_type == MPI_DATATYPE_NULL) {
2356 retval = MPI_ERR_TYPE;
2357 } else if (count<0 || blocklen<0){
2358 retval = MPI_ERR_COUNT;
2360 retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2366 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2370 if (old_type == MPI_DATATYPE_NULL) {
2371 retval = MPI_ERR_TYPE;
2372 } else if (count<0 || blocklen<0){
2373 retval = MPI_ERR_COUNT;
2375 retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2381 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2382 return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2385 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2389 if (old_type == MPI_DATATYPE_NULL) {
2390 retval = MPI_ERR_TYPE;
2391 } else if (count<0){
2392 retval = MPI_ERR_COUNT;
2394 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2400 int PMPI_Type_create_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2404 if (old_type == MPI_DATATYPE_NULL) {
2405 retval = MPI_ERR_TYPE;
2406 } else if (count<0){
2407 retval = MPI_ERR_COUNT;
2409 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2415 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2419 if (old_type == MPI_DATATYPE_NULL) {
2420 retval = MPI_ERR_TYPE;
2421 } else if (count<0){
2422 retval = MPI_ERR_COUNT;
2424 int* blocklens=(int*)xbt_malloc(blocklength*count);
2425 for (i=0; i<count;i++)blocklens[i]=blocklength;
2426 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2427 xbt_free(blocklens);
2434 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2438 if (old_type == MPI_DATATYPE_NULL) {
2439 retval = MPI_ERR_TYPE;
2440 } else if (count<0){
2441 retval = MPI_ERR_COUNT;
2443 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2449 int PMPI_Type_create_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2450 return PMPI_Type_hindexed(count, blocklens,indices,old_type,new_type);
2453 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2457 if (old_type == MPI_DATATYPE_NULL) {
2458 retval = MPI_ERR_TYPE;
2459 } else if (count<0){
2460 retval = MPI_ERR_COUNT;
2462 int* blocklens=(int*)xbt_malloc(blocklength*count);
2463 for (i=0; i<count;i++)blocklens[i]=blocklength;
2464 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2465 xbt_free(blocklens);
2472 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2477 retval = MPI_ERR_COUNT;
2479 retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2485 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2486 return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2490 int PMPI_Error_class(int errorcode, int* errorclass) {
2491 // assume smpi uses only standard mpi error codes
2492 *errorclass=errorcode;
2497 int PMPI_Initialized(int* flag) {
2498 *flag=smpi_process_initialized();
2502 /* The following calls are not yet implemented and will fail at runtime. */
2503 /* Once implemented, please move them above this notice. */
2505 #define NOT_YET_IMPLEMENTED {\
2506 XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__);\
2507 return MPI_SUCCESS;\
2511 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2515 int PMPI_Type_set_name(MPI_Datatype datatype, char * name)
2520 int PMPI_Type_get_name(MPI_Datatype datatype, char * name, int* len)
2525 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2529 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2533 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periods, int reorder, MPI_Comm* comm_cart) {
2537 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2541 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2545 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2549 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2553 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2557 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2561 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2565 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2569 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2573 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2577 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2581 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2585 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2589 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2593 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2597 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2601 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2605 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2609 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2613 int PMPI_Comm_get_errhandler(MPI_Comm comm, MPI_Errhandler* errhandler) {
2617 int PMPI_Cancel(MPI_Request* request) {
2621 int PMPI_Buffer_attach(void* buffer, int size) {
2625 int PMPI_Buffer_detach(void* buffer, int* size) {
2629 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2633 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2638 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2643 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2648 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2653 int PMPI_Comm_free_keyval(int* keyval) {
2657 int PMPI_Pcontrol(const int level )
2662 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2666 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2671 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2676 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2681 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2686 int PMPI_Type_free_keyval(int* keyval) {
2690 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2694 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2698 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2702 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2706 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2710 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2714 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2718 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2722 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2726 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2730 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2734 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2738 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2742 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2746 int PMPI_Keyval_free(int* keyval) {
2750 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2754 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2758 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2762 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2766 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2770 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2774 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2778 int PMPI_Win_fence( int assert, MPI_Win win){
2782 int PMPI_Win_free( MPI_Win* win){
2786 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2790 int PMPI_Info_create( MPI_Info *info){
2794 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2798 int PMPI_Info_free( MPI_Info *info){
2802 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2803 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2807 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2808 int *num_addresses, int *num_datatypes, int *combiner){
2812 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2813 int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2814 MPI_Datatype* array_of_datatypes){
2818 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2819 int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2820 int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2824 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2828 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){
2832 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
2836 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
2837 void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
2842 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
2846 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
2850 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
2854 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
2858 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
2862 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
2866 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
2870 int PMPI_Add_error_class( int *errorclass){
2874 int PMPI_Add_error_code( int errorclass, int *errorcode){
2878 int PMPI_Add_error_string( int errorcode, char *string){
2882 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
2886 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
2890 int PMPI_Info_delete(MPI_Info info, char *key){
2894 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
2898 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
2902 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
2906 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
2910 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){
2914 int PMPI_Grequest_complete( MPI_Request request){
2918 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
2922 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
2926 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2930 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
2934 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
2938 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
2942 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
2946 int PMPI_Open_port( MPI_Info info, char *port_name){
2950 int PMPI_Close_port(char *port_name){
2954 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
2958 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){
2962 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
2963 int* array_of_maxprocs, MPI_Info* array_of_info, int root,
2964 MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
2968 int PMPI_Comm_get_parent( MPI_Comm *parent){