5f686a90d0
of MPI_Alltoall. - add support for MPI_IN_PLACE in the self collective component. - fix the extent usage in the tuned collective component. - correctly use the peer counts instead of local - add support for MPI_IN_PLACE in the self collective component. - fix the extent usage in the tuned collective component. - correctly use the peer counts instead of local. Thanks to Fujitsu for the patch. This commit was SVN r29187.
328 строки
11 KiB
C
328 строки
11 KiB
C
/* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil -*- */
|
|
/*
|
|
* Copyright (c) 2004-2005 The Trustees of Indiana University and Indiana
|
|
* University Research and Technology
|
|
* Corporation. All rights reserved.
|
|
* Copyright (c) 2004-2005 The University of Tennessee and The University
|
|
* of Tennessee Research Foundation. All rights
|
|
* reserved.
|
|
* Copyright (c) 2004-2005 High Performance Computing Center Stuttgart,
|
|
* University of Stuttgart. All rights reserved.
|
|
* Copyright (c) 2004-2005 The Regents of the University of California.
|
|
* All rights reserved.
|
|
* Copyright (c) 2013 Los Alamos National Security, LLC. All rights
|
|
* reserved.
|
|
* Copyright (c) 2013 FUJITSU LIMITED. All rights reserved.
|
|
* $COPYRIGHT$
|
|
*
|
|
* Additional copyrights may follow
|
|
*
|
|
* $HEADER$
|
|
*/
|
|
|
|
#include "ompi_config.h"
|
|
#include "coll_basic.h"
|
|
|
|
#include "mpi.h"
|
|
#include "ompi/constants.h"
|
|
#include "ompi/datatype/ompi_datatype.h"
|
|
#include "ompi/mca/coll/coll.h"
|
|
#include "ompi/mca/coll/base/coll_tags.h"
|
|
#include "ompi/mca/pml/pml.h"
|
|
|
|
|
|
static int
|
|
mca_coll_basic_alltoallv_intra_inplace(void *rbuf, const int *rcounts, const int *rdisps,
|
|
struct ompi_datatype_t *rdtype,
|
|
struct ompi_communicator_t *comm,
|
|
mca_coll_base_module_t *module)
|
|
{
|
|
mca_coll_basic_module_t *basic_module = (mca_coll_basic_module_t*) module;
|
|
int i, j, size, rank, err=MPI_SUCCESS;
|
|
MPI_Request *preq;
|
|
char *tmp_buffer;
|
|
size_t max_size;
|
|
ptrdiff_t ext;
|
|
|
|
/* Initialize. */
|
|
|
|
size = ompi_comm_size(comm);
|
|
rank = ompi_comm_rank(comm);
|
|
|
|
/* If only one process, we're done. */
|
|
if (1 == size) {
|
|
return MPI_SUCCESS;
|
|
}
|
|
|
|
/* Find the largest receive amount */
|
|
ompi_datatype_type_extent (rdtype, &ext);
|
|
for (i = 0, max_size = 0 ; i < size ; ++i) {
|
|
size_t size = ext * rcounts[i];
|
|
|
|
max_size = size > max_size ? size : max_size;
|
|
}
|
|
|
|
/* Allocate a temporary buffer */
|
|
tmp_buffer = calloc (max_size, 1);
|
|
if (NULL == tmp_buffer) {
|
|
return OMPI_ERR_OUT_OF_RESOURCE;
|
|
}
|
|
|
|
/* in-place alltoallv slow algorithm (but works) */
|
|
for (i = 0 ; i < size ; ++i) {
|
|
for (j = i+1 ; j < size ; ++j) {
|
|
/* Initiate all send/recv to/from others. */
|
|
preq = basic_module->mccb_reqs;
|
|
|
|
if (i == rank && rcounts[j]) {
|
|
/* Copy the data into the temporary buffer */
|
|
err = ompi_datatype_copy_content_same_ddt (rdtype, rcounts[j],
|
|
tmp_buffer, (char *) rbuf + rdisps[j] * ext);
|
|
if (MPI_SUCCESS != err) { goto error_hndl; }
|
|
|
|
/* Exchange data with the peer */
|
|
err = MCA_PML_CALL(irecv ((char *) rbuf + rdisps[j] * ext, rcounts[j], rdtype,
|
|
j, MCA_COLL_BASE_TAG_ALLTOALLV, comm, preq++));
|
|
if (MPI_SUCCESS != err) { goto error_hndl; }
|
|
|
|
err = MCA_PML_CALL(isend ((void *) tmp_buffer, rcounts[j], rdtype,
|
|
j, MCA_COLL_BASE_TAG_ALLTOALLV, MCA_PML_BASE_SEND_STANDARD,
|
|
comm, preq++));
|
|
if (MPI_SUCCESS != err) { goto error_hndl; }
|
|
} else if (j == rank && rcounts[i]) {
|
|
/* Copy the data into the temporary buffer */
|
|
err = ompi_datatype_copy_content_same_ddt (rdtype, rcounts[i],
|
|
tmp_buffer, (char *) rbuf + rdisps[i] * ext);
|
|
if (MPI_SUCCESS != err) { goto error_hndl; }
|
|
|
|
/* Exchange data with the peer */
|
|
err = MCA_PML_CALL(irecv ((char *) rbuf + rdisps[i] * ext, rcounts[i], rdtype,
|
|
i, MCA_COLL_BASE_TAG_ALLTOALLV, comm, preq++));
|
|
if (MPI_SUCCESS != err) { goto error_hndl; }
|
|
|
|
err = MCA_PML_CALL(isend ((void *) tmp_buffer, rcounts[i], rdtype,
|
|
i, MCA_COLL_BASE_TAG_ALLTOALLV, MCA_PML_BASE_SEND_STANDARD,
|
|
comm, preq++));
|
|
if (MPI_SUCCESS != err) { goto error_hndl; }
|
|
} else {
|
|
continue;
|
|
}
|
|
|
|
/* Wait for the requests to complete */
|
|
err = ompi_request_wait_all (2, basic_module->mccb_reqs, MPI_STATUS_IGNORE);
|
|
if (MPI_SUCCESS != err) { goto error_hndl; }
|
|
|
|
/* Free the requests. */
|
|
mca_coll_basic_free_reqs(basic_module->mccb_reqs, 2);
|
|
}
|
|
}
|
|
|
|
error_hndl:
|
|
/* Free the temporary buffer */
|
|
free (tmp_buffer);
|
|
|
|
/* All done */
|
|
|
|
return err;
|
|
}
|
|
|
|
/*
|
|
* alltoallv_intra
|
|
*
|
|
* Function: - MPI_Alltoallv
|
|
* Accepts: - same as MPI_Alltoallv()
|
|
* Returns: - MPI_SUCCESS or an MPI error code
|
|
*/
|
|
int
|
|
mca_coll_basic_alltoallv_intra(void *sbuf, int *scounts, int *sdisps,
|
|
struct ompi_datatype_t *sdtype,
|
|
void *rbuf, int *rcounts, int *rdisps,
|
|
struct ompi_datatype_t *rdtype,
|
|
struct ompi_communicator_t *comm,
|
|
mca_coll_base_module_t *module)
|
|
{
|
|
int i;
|
|
int size;
|
|
int rank;
|
|
int err;
|
|
char *psnd;
|
|
char *prcv;
|
|
int nreqs;
|
|
MPI_Aint sndextent;
|
|
MPI_Aint rcvextent;
|
|
MPI_Request *preq;
|
|
|
|
mca_coll_basic_module_t *basic_module = (mca_coll_basic_module_t*) module;
|
|
|
|
/* Initialize. */
|
|
if (MPI_IN_PLACE == sbuf) {
|
|
return mca_coll_basic_alltoallv_intra_inplace (rbuf, rcounts, rdisps,
|
|
rdtype, comm, module);
|
|
}
|
|
|
|
size = ompi_comm_size(comm);
|
|
rank = ompi_comm_rank(comm);
|
|
|
|
ompi_datatype_type_extent(sdtype, &sndextent);
|
|
ompi_datatype_type_extent(rdtype, &rcvextent);
|
|
|
|
/* simple optimization */
|
|
|
|
psnd = ((char *) sbuf) + (sdisps[rank] * sndextent);
|
|
prcv = ((char *) rbuf) + (rdisps[rank] * rcvextent);
|
|
|
|
if (0 != scounts[rank]) {
|
|
err = ompi_datatype_sndrcv(psnd, scounts[rank], sdtype,
|
|
prcv, rcounts[rank], rdtype);
|
|
if (MPI_SUCCESS != err) {
|
|
return err;
|
|
}
|
|
}
|
|
|
|
/* If only one process, we're done. */
|
|
|
|
if (1 == size) {
|
|
return MPI_SUCCESS;
|
|
}
|
|
|
|
/* Initiate all send/recv to/from others. */
|
|
|
|
nreqs = 0;
|
|
preq = basic_module->mccb_reqs;
|
|
|
|
/* Post all receives first -- a simple optimization */
|
|
|
|
for (i = 0; i < size; ++i) {
|
|
if (i == rank || 0 == rcounts[i]) {
|
|
continue;
|
|
}
|
|
|
|
prcv = ((char *) rbuf) + (rdisps[i] * rcvextent);
|
|
err = MCA_PML_CALL(irecv_init(prcv, rcounts[i], rdtype,
|
|
i, MCA_COLL_BASE_TAG_ALLTOALLV, comm,
|
|
preq++));
|
|
++nreqs;
|
|
if (MPI_SUCCESS != err) {
|
|
mca_coll_basic_free_reqs(basic_module->mccb_reqs, nreqs);
|
|
return err;
|
|
}
|
|
}
|
|
|
|
/* Now post all sends */
|
|
|
|
for (i = 0; i < size; ++i) {
|
|
if (i == rank || 0 == scounts[i]) {
|
|
continue;
|
|
}
|
|
|
|
psnd = ((char *) sbuf) + (sdisps[i] * sndextent);
|
|
err = MCA_PML_CALL(isend_init(psnd, scounts[i], sdtype,
|
|
i, MCA_COLL_BASE_TAG_ALLTOALLV,
|
|
MCA_PML_BASE_SEND_STANDARD, comm,
|
|
preq++));
|
|
++nreqs;
|
|
if (MPI_SUCCESS != err) {
|
|
mca_coll_basic_free_reqs(basic_module->mccb_reqs, nreqs);
|
|
return err;
|
|
}
|
|
}
|
|
|
|
/* Start your engines. This will never return an error. */
|
|
|
|
MCA_PML_CALL(start(nreqs, basic_module->mccb_reqs));
|
|
|
|
/* Wait for them all. If there's an error, note that we don't care
|
|
* what the error was -- just that there *was* an error. The PML
|
|
* will finish all requests, even if one or more of them fail.
|
|
* i.e., by the end of this call, all the requests are free-able.
|
|
* So free them anyway -- even if there was an error, and return the
|
|
* error after we free everything. */
|
|
|
|
err = ompi_request_wait_all(nreqs, basic_module->mccb_reqs,
|
|
MPI_STATUSES_IGNORE);
|
|
|
|
/* Free the requests. */
|
|
|
|
mca_coll_basic_free_reqs(basic_module->mccb_reqs, nreqs);
|
|
|
|
/* All done */
|
|
|
|
return err;
|
|
}
|
|
|
|
|
|
/*
|
|
* alltoallv_inter
|
|
*
|
|
* Function: - MPI_Alltoallv
|
|
* Accepts: - same as MPI_Alltoallv()
|
|
* Returns: - MPI_SUCCESS or an MPI error code
|
|
*/
|
|
int
|
|
mca_coll_basic_alltoallv_inter(void *sbuf, int *scounts, int *sdisps,
|
|
struct ompi_datatype_t *sdtype, void *rbuf,
|
|
int *rcounts, int *rdisps,
|
|
struct ompi_datatype_t *rdtype,
|
|
struct ompi_communicator_t *comm,
|
|
mca_coll_base_module_t *module)
|
|
{
|
|
int i;
|
|
int rsize;
|
|
int err;
|
|
char *psnd;
|
|
char *prcv;
|
|
size_t nreqs;
|
|
MPI_Aint sndextent;
|
|
MPI_Aint rcvextent;
|
|
|
|
mca_coll_basic_module_t *basic_module = (mca_coll_basic_module_t*) module;
|
|
ompi_request_t **preq = basic_module->mccb_reqs;
|
|
|
|
/* Initialize. */
|
|
|
|
rsize = ompi_comm_remote_size(comm);
|
|
|
|
ompi_datatype_type_extent(sdtype, &sndextent);
|
|
ompi_datatype_type_extent(rdtype, &rcvextent);
|
|
|
|
/* Initiate all send/recv to/from others. */
|
|
nreqs = rsize * 2;
|
|
|
|
/* Post all receives first */
|
|
/* A simple optimization: do not send and recv msgs of length zero */
|
|
for (i = 0; i < rsize; ++i) {
|
|
prcv = ((char *) rbuf) + (rdisps[i] * rcvextent);
|
|
if (rcounts[i] > 0) {
|
|
err = MCA_PML_CALL(irecv(prcv, rcounts[i], rdtype,
|
|
i, MCA_COLL_BASE_TAG_ALLTOALLV, comm,
|
|
&preq[i]));
|
|
if (MPI_SUCCESS != err) {
|
|
return err;
|
|
}
|
|
} else {
|
|
preq[i] = MPI_REQUEST_NULL;
|
|
}
|
|
}
|
|
|
|
/* Now post all sends */
|
|
for (i = 0; i < rsize; ++i) {
|
|
psnd = ((char *) sbuf) + (sdisps[i] * sndextent);
|
|
if (scounts[i] > 0) {
|
|
err = MCA_PML_CALL(isend(psnd, scounts[i], sdtype,
|
|
i, MCA_COLL_BASE_TAG_ALLTOALLV,
|
|
MCA_PML_BASE_SEND_STANDARD, comm,
|
|
&preq[rsize + i]));
|
|
if (MPI_SUCCESS != err) {
|
|
return err;
|
|
}
|
|
} else {
|
|
preq[rsize + i] = MPI_REQUEST_NULL;
|
|
}
|
|
}
|
|
|
|
err = ompi_request_wait_all(nreqs, preq, MPI_STATUSES_IGNORE);
|
|
|
|
/* All done */
|
|
return err;
|
|
}
|