From: Augustin Degomme Date: Thu, 4 Apr 2013 13:50:53 +0000 (+0200) Subject: port a reduce algorithms with internal smpi calls X-Git-Tag: v3_9_90~412^2~49 X-Git-Url: http://info.iut-bm.univ-fcomte.fr/pub/gitweb/simgrid.git/commitdiff_plain/9ce9ffe42a1e578ec44ccf1cb804dff809d32c5d?ds=inline port a reduce algorithms with internal smpi calls --- diff --git a/src/smpi/colls/reduce-arrival-pattern-aware.c b/src/smpi/colls/reduce-arrival-pattern-aware.c index 0182a6b516..cf5df9eae6 100644 --- a/src/smpi/colls/reduce-arrival-pattern-aware.c +++ b/src/smpi/colls/reduce-arrival-pattern-aware.c @@ -42,14 +42,14 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, int header_buf[HEADER_SIZE]; char temp_buf[MAX_NODE]; - MPI_Aint extent; - MPI_Type_extent(datatype, &extent); + MPI_Aint extent, lb; + smpi_datatype_extent(datatype, &lb, &extent); /* source and destination */ int to, from; - MPI_Comm_rank(MPI_COMM_WORLD, &rank); - MPI_Comm_size(MPI_COMM_WORLD, &size); + size=smpi_comm_size(comm); + rank=smpi_comm_rank(comm); /* segment is segment size in number of elements (not bytes) */ @@ -74,7 +74,7 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, char *tmp_buf; tmp_buf = (char *) malloc(count * extent); - MPI_Sendrecv(buf, count, datatype, rank, tag, rbuf, count, datatype, rank, + smpi_mpi_sendrecv(buf, count, datatype, rank, tag, rbuf, count, datatype, rank, tag, comm, &status); @@ -89,7 +89,7 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, for (i = 1; i < size; i++) { if (already_received[i] == 0) - MPI_Iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i], + smpi_mpi_iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i], MPI_STATUSES_IGNORE); } @@ -101,7 +101,7 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, /* 1-byte message arrive */ if ((flag_array[i] == 1) && (already_received[i] == 0)) { - MPI_Recv(temp_buf, 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status); + smpi_mpi_recv(temp_buf, 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status); header_buf[header_index] = i; header_index++; sent_count++; @@ -125,8 +125,8 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, to = header_buf[0]; from = header_buf[header_index - 1]; - MPI_Send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm); - MPI_Recv(tmp_buf, count, datatype, from, tag, comm, &status); + smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm); + smpi_mpi_recv(tmp_buf, count, datatype, from, tag, comm, &status); star_reduction(op, tmp_buf, rbuf, &count, &datatype); } } /* while loop */ @@ -137,10 +137,10 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, else { /* send 1-byte message to root */ - MPI_Send(temp_buf, 1, MPI_CHAR, 0, tag, comm); + smpi_mpi_send(temp_buf, 1, MPI_CHAR, 0, tag, comm); /* wait for header and data, forward when required */ - MPI_Recv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm, + smpi_mpi_recv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm, &status); // MPI_Recv(buf,count,datatype,MPI_ANY_SOURCE,tag,comm,&status); @@ -152,7 +152,7 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, /* forward header */ if (header_buf[myordering + 1] != -1) { - MPI_Send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1], + smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1], tag, comm); } //printf("node %d ordering %d\n",rank,myordering); @@ -166,7 +166,7 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, } else { to = header_buf[myordering + 1]; } - MPI_Send(rbuf, count, datatype, to, tag, comm); + smpi_mpi_send(rbuf, count, datatype, to, tag, comm); } /* recv, reduce, send */ @@ -177,10 +177,10 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, to = header_buf[myordering + 1]; } from = header_buf[myordering - 1]; - MPI_Recv(tmp_buf, count, datatype, header_buf[myordering - 1], tag, + smpi_mpi_recv(tmp_buf, count, datatype, header_buf[myordering - 1], tag, comm, &status); star_reduction(op, tmp_buf, rbuf, &count, &datatype); - MPI_Send(rbuf, count, datatype, to, tag, comm); + smpi_mpi_send(rbuf, count, datatype, to, tag, comm); } } /* non-root */ } @@ -212,11 +212,11 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, //if (i == rank) //continue; if ((already_received[i] == 0) && (will_send[i] == 0)) { - MPI_Iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i], + smpi_mpi_iprobe(i, MPI_ANY_TAG, MPI_COMM_WORLD, &flag_array[i], &temp_status_array[i]); if (flag_array[i] == 1) { will_send[i] = 1; - MPI_Recv(&temp_buf[i], 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, + smpi_mpi_recv(&temp_buf[i], 1, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status); //printf("recv from %d\n",i); i = 1; @@ -248,12 +248,12 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, to = header_buf[0]; /* send header */ - MPI_Send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm); + smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, to, tag, comm); /* recv data - pipeline */ from = header_buf[header_index - 1]; for (i = 0; i < pipe_length; i++) { - MPI_Recv(tmp_buf + (i * increment), segment, datatype, from, tag, + smpi_mpi_recv(tmp_buf + (i * increment), segment, datatype, from, tag, comm, &status); star_reduction(op, tmp_buf + (i * increment), (char *)rbuf + (i * increment), &segment, &datatype); @@ -266,13 +266,12 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, /* none root */ else { /* send 1-byte message to root */ - MPI_Send(temp_buf, 1, MPI_CHAR, 0, tag, comm); + smpi_mpi_send(temp_buf, 1, MPI_CHAR, 0, tag, comm); /* wait for header forward when required */ - MPI_Irecv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm, - &request); - MPI_Wait(&request, MPI_STATUS_IGNORE); + request=smpi_mpi_irecv(header_buf, HEADER_SIZE, MPI_INT, MPI_ANY_SOURCE, tag, comm); + smpi_mpi_wait(&request, MPI_STATUS_IGNORE); /* search for where it is */ int myordering = 0; @@ -283,7 +282,7 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, /* send header when required */ if (header_buf[myordering + 1] != -1) { - MPI_Send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1], + smpi_mpi_send(header_buf, HEADER_SIZE, MPI_INT, header_buf[myordering + 1], tag, comm); } @@ -297,27 +296,24 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, /* send only */ if (myordering == 0) { for (i = 0; i < pipe_length; i++) { - MPI_Isend((char *)rbuf + (i * increment), segment, datatype, to, tag, comm, - &send_request_array[i]); + send_request_array[i]= smpi_mpi_isend((char *)rbuf + (i * increment), segment, datatype, to, tag, comm); } - MPI_Waitall((pipe_length), send_request_array, send_status_array); + smpi_mpi_waitall((pipe_length), send_request_array, send_status_array); } /* receive, reduce, and send */ else { from = header_buf[myordering - 1]; for (i = 0; i < pipe_length; i++) { - MPI_Irecv(tmp_buf + (i * increment), segment, datatype, from, tag, - comm, &recv_request_array[i]); + recv_request_array[i]=smpi_mpi_irecv(tmp_buf + (i * increment), segment, datatype, from, tag, comm); } for (i = 0; i < pipe_length; i++) { - MPI_Wait(&recv_request_array[i], MPI_STATUS_IGNORE); + smpi_mpi_wait(&recv_request_array[i], MPI_STATUS_IGNORE); star_reduction(op, tmp_buf + (i * increment), (char *)rbuf + (i * increment), &segment, &datatype); - MPI_Isend((char *)rbuf + (i * increment), segment, datatype, to, tag, comm, - &send_request_array[i]); + send_request_array[i]=smpi_mpi_isend((char *)rbuf + (i * increment), segment, datatype, to, tag, comm); } - MPI_Waitall((pipe_length), send_request_array, send_status_array); + smpi_mpi_waitall((pipe_length), send_request_array, send_status_array); } } /* non-root */ @@ -338,16 +334,16 @@ int smpi_coll_tuned_reduce_arrival_pattern_aware(void *buf, void *rbuf, */ if (root != 0) { if (rank == 0) { - MPI_Send(rbuf, count, datatype, root, tag, comm); + smpi_mpi_send(rbuf, count, datatype, root, tag, comm); } else if (rank == root) { - MPI_Recv(rbuf, count, datatype, 0, tag, comm, &status); + smpi_mpi_recv(rbuf, count, datatype, 0, tag, comm, &status); } } /* when count is not divisible by block size, use default BCAST for the remainder */ if ((remainder != 0) && (count > segment)) { - MPI_Reduce((char *)buf + (pipe_length * increment), + smpi_mpi_reduce((char *)buf + (pipe_length * increment), (char *)rbuf + (pipe_length * increment), remainder, datatype, op, root, comm); }