X-Git-Url: http://info.iut-bm.univ-fcomte.fr/pub/gitweb/simgrid.git/blobdiff_plain/299c52751d8bf7ed4f9fc40f668025a5cb425901..71af21fb892c4e4d6676b8ac14b836e59a1af4bc:/src/smpi/colls/smpi_nbc_impl.cpp diff --git a/src/smpi/colls/smpi_nbc_impl.cpp b/src/smpi/colls/smpi_nbc_impl.cpp index f1db5eac28..51af4b1705 100644 --- a/src/smpi/colls/smpi_nbc_impl.cpp +++ b/src/smpi/colls/smpi_nbc_impl.cpp @@ -14,14 +14,12 @@ namespace smpi{ int Colls::ibarrier(MPI_Comm comm, MPI_Request* request) { - int i; int size = comm->size(); int rank = comm->rank(); - MPI_Request* requests; (*request) = new Request( nullptr, 0, MPI_BYTE, rank,rank, COLL_TAG_BARRIER, comm, MPI_REQ_PERSISTENT); if (rank > 0) { - requests = new MPI_Request[2]; + MPI_Request* requests = new MPI_Request[2]; requests[0] = Request::isend (nullptr, 0, MPI_BYTE, 0, COLL_TAG_BARRIER, comm); @@ -31,15 +29,10 @@ int Colls::ibarrier(MPI_Comm comm, MPI_Request* request) (*request)->set_nbc_requests(requests, 2); } else { - requests = new MPI_Request[(size-1)*2]; - for (i = 1; i < 2*size-1; i+=2) { - requests[i-1] = Request::irecv(nullptr, 0, MPI_BYTE, MPI_ANY_SOURCE, - COLL_TAG_BARRIER, comm - ); - requests[i] = Request::isend(nullptr, 0, MPI_BYTE, (i+1)/2, - COLL_TAG_BARRIER, - comm - ); + MPI_Request* requests = new MPI_Request[(size - 1) * 2]; + for (int i = 1; i < 2 * size - 1; i += 2) { + requests[i - 1] = Request::irecv(nullptr, 0, MPI_BYTE, MPI_ANY_SOURCE, COLL_TAG_BARRIER, comm); + requests[i] = Request::isend(nullptr, 0, MPI_BYTE, (i + 1) / 2, COLL_TAG_BARRIER, comm); } (*request)->set_nbc_requests(requests, 2*(size-1)); } @@ -48,23 +41,21 @@ int Colls::ibarrier(MPI_Comm comm, MPI_Request* request) int Colls::ibcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm, MPI_Request* request) { - int i; int size = comm->size(); int rank = comm->rank(); - MPI_Request* requests; (*request) = new Request( nullptr, 0, MPI_BYTE, - rank,rank, COLL_TAG_BARRIER, comm, MPI_REQ_PERSISTENT); + rank,rank, COLL_TAG_BCAST, comm, MPI_REQ_PERSISTENT); if (rank != root) { - requests = new MPI_Request[1]; + MPI_Request* requests = new MPI_Request[1]; requests[0] = Request::irecv (buf, count, datatype, root, COLL_TAG_BCAST, comm); (*request)->set_nbc_requests(requests, 1); } else { - requests = new MPI_Request[size-1]; + MPI_Request* requests = new MPI_Request[size - 1]; int n = 0; - for (i = 0; i < size; i++) { + for (int i = 0; i < size; i++) { if(i!=root){ requests[n] = Request::isend(buf, count, datatype, i, COLL_TAG_BCAST, @@ -78,26 +69,25 @@ int Colls::ibcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Com return MPI_SUCCESS; } -int Colls::iallgather(void *sendbuf, int sendcount, MPI_Datatype sendtype, +int Colls::iallgather(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf,int recvcount, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request) { const int system_tag = COLL_TAG_ALLGATHER; MPI_Aint lb = 0; MPI_Aint recvext = 0; - MPI_Request *requests; int rank = comm->rank(); int size = comm->size(); (*request) = new Request( nullptr, 0, MPI_BYTE, - rank,rank, COLL_TAG_BARRIER, comm, MPI_REQ_PERSISTENT); + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT); // FIXME: check for errors recvtype->extent(&lb, &recvext); // Local copy from self Datatype::copy(sendbuf, sendcount, sendtype, static_cast(recvbuf) + rank * recvcount * recvext, recvcount, recvtype); // Send/Recv buffers to/from others; - requests = new MPI_Request[2 * (size - 1)]; + MPI_Request* requests = new MPI_Request[2 * (size - 1)]; int index = 0; for (int other = 0; other < size; other++) { if(other != rank) { @@ -113,20 +103,19 @@ int Colls::iallgather(void *sendbuf, int sendcount, MPI_Datatype sendtype, return MPI_SUCCESS; } -int Colls::iscatter(void *sendbuf, int sendcount, MPI_Datatype sendtype, +int Colls::iscatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request* request) { const int system_tag = COLL_TAG_SCATTER; MPI_Aint lb = 0; MPI_Aint sendext = 0; - MPI_Request *requests; int rank = comm->rank(); int size = comm->size(); (*request) = new Request( nullptr, 0, MPI_BYTE, - rank,rank, COLL_TAG_BARRIER, comm, MPI_REQ_PERSISTENT); + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT); if(rank != root) { - requests = new MPI_Request[1]; + MPI_Request* requests = new MPI_Request[1]; // Recv buffer from root requests[0] = Request::irecv(recvbuf, recvcount, recvtype, root, system_tag, comm); (*request)->set_nbc_requests(requests, 1); @@ -134,15 +123,15 @@ int Colls::iscatter(void *sendbuf, int sendcount, MPI_Datatype sendtype, sendtype->extent(&lb, &sendext); // Local copy from root if(recvbuf!=MPI_IN_PLACE){ - Datatype::copy(static_cast(sendbuf) + root * sendcount * sendext, + Datatype::copy(static_cast(sendbuf) + root * sendcount * sendext, sendcount, sendtype, recvbuf, recvcount, recvtype); } // Send buffers to receivers - requests = new MPI_Request[size - 1]; + MPI_Request* requests = new MPI_Request[size - 1]; int index = 0; for(int dst = 0; dst < size; dst++) { if(dst != root) { - requests[index] = Request::isend_init(static_cast(sendbuf) + dst * sendcount * sendext, sendcount, sendtype, + requests[index] = Request::isend_init(static_cast(sendbuf) + dst * sendcount * sendext, sendcount, sendtype, dst, system_tag, comm); index++; } @@ -154,8 +143,8 @@ int Colls::iscatter(void *sendbuf, int sendcount, MPI_Datatype sendtype, return MPI_SUCCESS; } -int Colls::iallgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, - int *recvcounts, int *displs, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request) +int Colls::iallgatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, + const int *recvcounts, const int *displs, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request) { const int system_tag = COLL_TAG_ALLGATHERV; MPI_Aint lb = 0; @@ -164,7 +153,7 @@ int Colls::iallgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void int rank = comm->rank(); int size = comm->size(); (*request) = new Request( nullptr, 0, MPI_BYTE, - rank,rank, COLL_TAG_BARRIER, comm, MPI_REQ_PERSISTENT); + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT); recvtype->extent(&lb, &recvext); // Local copy from self Datatype::copy(sendbuf, sendcount, sendtype, @@ -188,29 +177,28 @@ int Colls::iallgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void return MPI_SUCCESS; } -int Colls::ialltoall( void *sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request){ -int system_tag = COLL_TAG_ALLTOALL; - int i; - int count; - MPI_Aint lb = 0, sendext = 0, recvext = 0; - MPI_Request *requests; +int Colls::ialltoall( const void *sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request* request){ + int system_tag = COLL_TAG_ALLTOALL; + MPI_Aint lb = 0; + MPI_Aint sendext = 0; + MPI_Aint recvext = 0; /* Initialize. */ int rank = comm->rank(); int size = comm->size(); (*request) = new Request( nullptr, 0, MPI_BYTE, - rank,rank, COLL_TAG_ALLTOALL, comm, MPI_REQ_PERSISTENT); + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT); sendtype->extent(&lb, &sendext); recvtype->extent(&lb, &recvext); /* simple optimization */ - int err = Datatype::copy(static_cast(sendbuf) + rank * sendcount * sendext, sendcount, sendtype, + int err = Datatype::copy(static_cast(sendbuf) + rank * sendcount * sendext, sendcount, sendtype, static_cast(recvbuf) + rank * recvcount * recvext, recvcount, recvtype); if (err == MPI_SUCCESS && size > 1) { /* Initiate all send/recv to/from others. */ - requests = new MPI_Request[2 * (size - 1)]; + MPI_Request* requests = new MPI_Request[2 * (size - 1)]; /* Post all receives first -- a simple optimization */ - count = 0; - for (i = (rank + 1) % size; i != rank; i = (i + 1) % size) { + int count = 0; + for (int i = (rank + 1) % size; i != rank; i = (i + 1) % size) { requests[count] = Request::irecv_init(static_cast(recvbuf) + i * recvcount * recvext, recvcount, recvtype, i, system_tag, comm); count++; @@ -220,8 +208,8 @@ int system_tag = COLL_TAG_ALLTOALL; * when messages actually arrive in the order in which they were posted. * TODO: check the previous assertion */ - for (i = (rank + size - 1) % size; i != rank; i = (i + size - 1) % size) { - requests[count] = Request::isend_init(static_cast(sendbuf) + i * sendcount * sendext, sendcount, + for (int i = (rank + size - 1) % size; i != rank; i = (i + size - 1) % size) { + requests[count] = Request::isend_init(static_cast(sendbuf) + i * sendcount * sendext, sendcount, sendtype, i, system_tag, comm); count++; } @@ -232,31 +220,30 @@ int system_tag = COLL_TAG_ALLTOALL; return MPI_SUCCESS; } -int Colls::ialltoallv(void *sendbuf, int *sendcounts, int *senddisps, MPI_Datatype sendtype, - void *recvbuf, int *recvcounts, int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request *request){ +int Colls::ialltoallv(const void *sendbuf, const int *sendcounts, const int *senddisps, MPI_Datatype sendtype, + void *recvbuf, const int *recvcounts, const int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm, MPI_Request *request){ const int system_tag = COLL_TAG_ALLTOALLV; MPI_Aint lb = 0; MPI_Aint sendext = 0; MPI_Aint recvext = 0; - MPI_Request *requests; /* Initialize. */ int rank = comm->rank(); int size = comm->size(); (*request) = new Request( nullptr, 0, MPI_BYTE, - rank,rank, COLL_TAG_ALLTOALLV, comm, MPI_REQ_PERSISTENT); + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT); sendtype->extent(&lb, &sendext); recvtype->extent(&lb, &recvext); /* Local copy from self */ - int err = Datatype::copy(static_cast(sendbuf) + senddisps[rank] * sendext, sendcounts[rank], sendtype, + int err = Datatype::copy(static_cast(sendbuf) + senddisps[rank] * sendext, sendcounts[rank], sendtype, static_cast(recvbuf) + recvdisps[rank] * recvext, recvcounts[rank], recvtype); if (err == MPI_SUCCESS && size > 1) { /* Initiate all send/recv to/from others. */ - requests = new MPI_Request[2 * (size - 1)]; + MPI_Request* requests = new MPI_Request[2 * (size - 1)]; int count = 0; /* Create all receives that will be posted first */ for (int i = 0; i < size; ++i) { - if (i != rank && recvcounts[i] != 0) { + if (i != rank) { requests[count] = Request::irecv_init(static_cast(recvbuf) + recvdisps[i] * recvext, recvcounts[i], recvtype, i, system_tag, comm); count++; @@ -266,8 +253,8 @@ int Colls::ialltoallv(void *sendbuf, int *sendcounts, int *senddisps, MPI_Dataty } /* Now create all sends */ for (int i = 0; i < size; ++i) { - if (i != rank && sendcounts[i] != 0) { - requests[count] = Request::isend_init(static_cast(sendbuf) + senddisps[i] * sendext, + if (i != rank) { + requests[count] = Request::isend_init(static_cast(sendbuf) + senddisps[i] * sendext, sendcounts[i], sendtype, i, system_tag, comm); count++; }else{ @@ -281,21 +268,63 @@ int Colls::ialltoallv(void *sendbuf, int *sendcounts, int *senddisps, MPI_Dataty return err; } -int Colls::igather(void *sendbuf, int sendcount, MPI_Datatype sendtype, +int Colls::ialltoallw(const void *sendbuf, const int *sendcounts, const int *senddisps, const MPI_Datatype* sendtypes, + void *recvbuf, const int *recvcounts, const int *recvdisps, const MPI_Datatype* recvtypes, MPI_Comm comm, MPI_Request *request){ + const int system_tag = COLL_TAG_ALLTOALLV; + + /* Initialize. */ + int rank = comm->rank(); + int size = comm->size(); + (*request) = new Request( nullptr, 0, MPI_BYTE, + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT); + /* Local copy from self */ + int err = (sendcounts[rank]>0 && recvcounts[rank]) ? Datatype::copy(static_cast(sendbuf) + senddisps[rank], sendcounts[rank], sendtypes[rank], + static_cast(recvbuf) + recvdisps[rank], recvcounts[rank], recvtypes[rank]): MPI_SUCCESS; + if (err == MPI_SUCCESS && size > 1) { + /* Initiate all send/recv to/from others. */ + MPI_Request* requests = new MPI_Request[2 * (size - 1)]; + int count = 0; + /* Create all receives that will be posted first */ + for (int i = 0; i < size; ++i) { + if (i != rank) { + requests[count] = Request::irecv_init(static_cast(recvbuf) + recvdisps[i], + recvcounts[i], recvtypes[i], i, system_tag, comm); + count++; + }else{ + XBT_DEBUG("<%d> skip request creation [src = %d, recvcounts[src] = %d]", rank, i, recvcounts[i]); + } + } + /* Now create all sends */ + for (int i = 0; i < size; ++i) { + if (i != rank) { + requests[count] = Request::isend_init(static_cast(sendbuf) + senddisps[i] , + sendcounts[i], sendtypes[i], i, system_tag, comm); + count++; + }else{ + XBT_DEBUG("<%d> skip request creation [dst = %d, sendcounts[dst] = %d]", rank, i, sendcounts[i]); + } + } + /* Wait for them all. */ + Request::startall(count, requests); + (*request)->set_nbc_requests(requests, count); + } + return err; +} + +int Colls::igather(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request *request) { const int system_tag = COLL_TAG_GATHER; MPI_Aint lb = 0; MPI_Aint recvext = 0; - MPI_Request *requests; int rank = comm->rank(); int size = comm->size(); (*request) = new Request( nullptr, 0, MPI_BYTE, - rank,rank, COLL_TAG_GATHER, comm, MPI_REQ_PERSISTENT); + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT); if(rank != root) { // Send buffer to root - requests = new MPI_Request[1]; + MPI_Request* requests = new MPI_Request[1]; requests[0]=Request::isend(sendbuf, sendcount, sendtype, root, system_tag, comm); (*request)->set_nbc_requests(requests, 1); } else { @@ -304,7 +333,7 @@ int Colls::igather(void *sendbuf, int sendcount, MPI_Datatype sendtype, Datatype::copy(sendbuf, sendcount, sendtype, static_cast(recvbuf) + root * recvcount * recvext, recvcount, recvtype); // Receive buffers from senders - requests = new MPI_Request[size - 1]; + MPI_Request* requests = new MPI_Request[size - 1]; int index = 0; for (int src = 0; src < size; src++) { if(src != root) { @@ -320,21 +349,20 @@ int Colls::igather(void *sendbuf, int sendcount, MPI_Datatype sendtype, return MPI_SUCCESS; } -int Colls::igatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int *recvcounts, int *displs, +int Colls::igatherv(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int *recvcounts, const int *displs, MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request *request) { int system_tag = COLL_TAG_GATHERV; MPI_Aint lb = 0; MPI_Aint recvext = 0; - MPI_Request *requests; int rank = comm->rank(); int size = comm->size(); (*request) = new Request( nullptr, 0, MPI_BYTE, - rank,rank, COLL_TAG_GATHERV, comm, MPI_REQ_PERSISTENT); + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT); if (rank != root) { // Send buffer to root - requests = new MPI_Request[1]; + MPI_Request* requests = new MPI_Request[1]; requests[0]=Request::isend(sendbuf, sendcount, sendtype, root, system_tag, comm); (*request)->set_nbc_requests(requests, 1); } else { @@ -343,7 +371,7 @@ int Colls::igatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *r Datatype::copy(sendbuf, sendcount, sendtype, static_cast(recvbuf) + displs[root] * recvext, recvcounts[root], recvtype); // Receive buffers from senders - requests = new MPI_Request[size - 1]; + MPI_Request* requests = new MPI_Request[size - 1]; int index = 0; for (int src = 0; src < size; src++) { if(src != root) { @@ -358,28 +386,27 @@ int Colls::igatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype, void *r } return MPI_SUCCESS; } -int Colls::iscatterv(void *sendbuf, int *sendcounts, int *displs, MPI_Datatype sendtype, void *recvbuf, int recvcount, +int Colls::iscatterv(const void *sendbuf, const int *sendcounts, const int *displs, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm, MPI_Request *request) { int system_tag = COLL_TAG_SCATTERV; MPI_Aint lb = 0; MPI_Aint sendext = 0; - MPI_Request* requests; int rank = comm->rank(); int size = comm->size(); (*request) = new Request( nullptr, 0, MPI_BYTE, - rank,rank, COLL_TAG_SCATTERV, comm, MPI_REQ_PERSISTENT); + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT); if(rank != root) { // Recv buffer from root - requests = new MPI_Request[1]; + MPI_Request* requests = new MPI_Request[1]; requests[0]=Request::irecv(recvbuf, recvcount, recvtype, root, system_tag, comm); (*request)->set_nbc_requests(requests, 1); } else { sendtype->extent(&lb, &sendext); // Local copy from root if(recvbuf!=MPI_IN_PLACE){ - Datatype::copy(static_cast(sendbuf) + displs[root] * sendext, sendcounts[root], + Datatype::copy(static_cast(sendbuf) + displs[root] * sendext, sendcounts[root], sendtype, recvbuf, recvcount, recvtype); } // Send buffers to receivers @@ -387,7 +414,7 @@ int Colls::iscatterv(void *sendbuf, int *sendcounts, int *displs, MPI_Datatype s int index = 0; for (int dst = 0; dst < size; dst++) { if (dst != root) { - requests[index] = Request::isend_init(static_cast(sendbuf) + displs[dst] * sendext, sendcounts[dst], + requests[index] = Request::isend_init(static_cast(sendbuf) + displs[dst] * sendext, sendcounts[dst], sendtype, dst, system_tag, comm); index++; } @@ -398,5 +425,197 @@ int Colls::iscatterv(void *sendbuf, int *sendcounts, int *displs, MPI_Datatype s } return MPI_SUCCESS; } + +int Colls::ireduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, + MPI_Comm comm, MPI_Request* request) +{ + const int system_tag = COLL_TAG_REDUCE; + MPI_Aint lb = 0; + MPI_Aint dataext = 0; + + const void* real_sendbuf = sendbuf; + + int rank = comm->rank(); + int size = comm->size(); + + if (size <= 0) + return MPI_ERR_COMM; + + unsigned char* tmp_sendbuf = nullptr; + if( sendbuf == MPI_IN_PLACE ) { + tmp_sendbuf = smpi_get_tmp_sendbuffer(count * datatype->get_extent()); + Datatype::copy(recvbuf, count, datatype, tmp_sendbuf, count, datatype); + real_sendbuf = tmp_sendbuf; + } + + if(rank == root){ + (*request) = new Request( recvbuf, count, datatype, + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT, op); + } + else + (*request) = new Request( nullptr, count, datatype, + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT); + + if(rank != root) { + // Send buffer to root + MPI_Request* requests = new MPI_Request[1]; + requests[0] = Request::isend(real_sendbuf, count, datatype, root, system_tag, comm); + (*request)->set_nbc_requests(requests, 1); + } else { + datatype->extent(&lb, &dataext); + // Local copy from root + if (real_sendbuf != nullptr && recvbuf != nullptr) + Datatype::copy(real_sendbuf, count, datatype, recvbuf, count, datatype); + // Receive buffers from senders + MPI_Request *requests = new MPI_Request[size - 1]; + int index = 0; + for (int src = 0; src < size; src++) { + if (src != root) { + requests[index] = + Request::irecv_init(smpi_get_tmp_sendbuffer(count * dataext), count, datatype, src, system_tag, comm); + index++; + } + } + // Wait for completion of irecv's. + Request::startall(size - 1, requests); + (*request)->set_nbc_requests(requests, size - 1); + } + if( sendbuf == MPI_IN_PLACE ) { + smpi_free_tmp_buffer(tmp_sendbuf); + } + return MPI_SUCCESS; +} + +int Colls::iallreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, + MPI_Op op, MPI_Comm comm, MPI_Request* request) +{ + + const int system_tag = COLL_TAG_ALLREDUCE; + MPI_Aint lb = 0; + MPI_Aint dataext = 0; + + int rank = comm->rank(); + int size = comm->size(); + (*request) = new Request( recvbuf, count, datatype, + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT, op); + // FIXME: check for errors + datatype->extent(&lb, &dataext); + // Local copy from self + Datatype::copy(sendbuf, count, datatype, recvbuf, count, datatype); + // Send/Recv buffers to/from others; + MPI_Request* requests = new MPI_Request[2 * (size - 1)]; + int index = 0; + for (int other = 0; other < size; other++) { + if(other != rank) { + requests[index] = Request::isend_init(sendbuf, count, datatype, other, system_tag,comm); + index++; + requests[index] = Request::irecv_init(smpi_get_tmp_sendbuffer(count * dataext), count, datatype, + other, system_tag, comm); + index++; + } + } + Request::startall(2 * (size - 1), requests); + (*request)->set_nbc_requests(requests, 2 * (size - 1)); + return MPI_SUCCESS; +} + +int Colls::iscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request) +{ + int system_tag = -888; + MPI_Aint lb = 0; + MPI_Aint dataext = 0; + + int rank = comm->rank(); + int size = comm->size(); + (*request) = new Request( recvbuf, count, datatype, + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT, op); + datatype->extent(&lb, &dataext); + + // Local copy from self + Datatype::copy(sendbuf, count, datatype, recvbuf, count, datatype); + + // Send/Recv buffers to/from others + MPI_Request *requests = new MPI_Request[size - 1]; + int index = 0; + for (int other = 0; other < rank; other++) { + requests[index] = Request::irecv_init(smpi_get_tmp_sendbuffer(count * dataext), count, datatype, other, system_tag, comm); + index++; + } + for (int other = rank + 1; other < size; other++) { + requests[index] = Request::isend_init(sendbuf, count, datatype, other, system_tag, comm); + index++; + } + // Wait for completion of all comms. + Request::startall(size - 1, requests); + (*request)->set_nbc_requests(requests, size - 1); + return MPI_SUCCESS; +} + +int Colls::iexscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, MPI_Request* request) +{ + int system_tag = -888; + MPI_Aint lb = 0; + MPI_Aint dataext = 0; + int rank = comm->rank(); + int size = comm->size(); + (*request) = new Request( recvbuf, count, datatype, + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT, op); + datatype->extent(&lb, &dataext); + if(rank != 0) + memset(recvbuf, 0, count*dataext); + + // Send/Recv buffers to/from others + MPI_Request *requests = new MPI_Request[size - 1]; + int index = 0; + for (int other = 0; other < rank; other++) { + requests[index] = Request::irecv_init(smpi_get_tmp_sendbuffer(count * dataext), count, datatype, other, system_tag, comm); + index++; + } + for (int other = rank + 1; other < size; other++) { + requests[index] = Request::isend_init(sendbuf, count, datatype, other, system_tag, comm); + index++; + } + // Wait for completion of all comms. + Request::startall(size - 1, requests); + (*request)->set_nbc_requests(requests, size - 1); + return MPI_SUCCESS; +} + +int Colls::ireduce_scatter(const void *sendbuf, void *recvbuf, const int *recvcounts, MPI_Datatype datatype, MPI_Op op, + MPI_Comm comm, MPI_Request* request){ +//Version where each process performs the reduce for its own part. Alltoall pattern for comms. + const int system_tag = COLL_TAG_REDUCE_SCATTER; + MPI_Aint lb = 0; + MPI_Aint dataext = 0; + + int rank = comm->rank(); + int size = comm->size(); + int count=recvcounts[rank]; + (*request) = new Request( recvbuf, count, datatype, + rank,rank, system_tag, comm, MPI_REQ_PERSISTENT, op); + datatype->extent(&lb, &dataext); + + // Send/Recv buffers to/from others; + MPI_Request* requests = new MPI_Request[2 * (size - 1)]; + int index = 0; + int recvdisp=0; + for (int other = 0; other < size; other++) { + if(other != rank) { + requests[index] = Request::isend_init(static_cast(sendbuf) + recvdisp * dataext, recvcounts[other], datatype, other, system_tag,comm); + XBT_VERB("sending with recvdisp %d", recvdisp); + index++; + requests[index] = Request::irecv_init(smpi_get_tmp_sendbuffer(count * dataext), count, datatype, + other, system_tag, comm); + index++; + }else{ + Datatype::copy(static_cast(sendbuf) + recvdisp * dataext, count, datatype, recvbuf, count, datatype); + } + recvdisp+=recvcounts[other]; + } + Request::startall(2 * (size - 1), requests); + (*request)->set_nbc_requests(requests, 2 * (size - 1)); + return MPI_SUCCESS; +} + } }