2004-01-12 00:26:55 +03:00
|
|
|
/*
|
2004-11-22 04:38:40 +03:00
|
|
|
* Copyright (c) 2004-2005 The Trustees of Indiana University.
|
|
|
|
* All rights reserved.
|
|
|
|
* Copyright (c) 2004-2005 The Trustees of the University of Tennessee.
|
|
|
|
* All rights reserved.
|
2004-11-28 23:09:25 +03:00
|
|
|
* Copyright (c) 2004-2005 High Performance Computing Center Stuttgart,
|
|
|
|
* University of Stuttgart. All rights reserved.
|
2005-03-24 15:43:37 +03:00
|
|
|
* Copyright (c) 2004-2005 The Regents of the University of California.
|
|
|
|
* All rights reserved.
|
2004-11-22 04:38:40 +03:00
|
|
|
* $COPYRIGHT$
|
|
|
|
*
|
|
|
|
* Additional copyrights may follow
|
|
|
|
*
|
2004-01-12 00:26:55 +03:00
|
|
|
* $HEADER$
|
|
|
|
*/
|
|
|
|
|
2004-06-07 19:33:53 +04:00
|
|
|
#include "ompi_config.h"
|
2004-02-15 00:26:30 +03:00
|
|
|
#include "coll_basic.h"
|
2004-01-12 00:26:55 +03:00
|
|
|
|
|
|
|
#include "mpi.h"
|
2005-08-13 01:42:07 +04:00
|
|
|
#include "ompi/include/constants.h"
|
2005-09-13 00:21:53 +04:00
|
|
|
#include "ompi/datatype/datatype.h"
|
|
|
|
#include "ompi/mca/coll/coll.h"
|
|
|
|
#include "ompi/mca/coll/base/coll_tags.h"
|
|
|
|
#include "ompi/mca/pml/pml.h"
|
2004-01-12 00:26:55 +03:00
|
|
|
#include "coll_basic.h"
|
|
|
|
|
|
|
|
|
|
|
|
/*
|
2004-06-29 04:02:25 +04:00
|
|
|
* scatter_intra
|
2004-01-12 00:26:55 +03:00
|
|
|
*
|
|
|
|
* Function: - scatter operation
|
|
|
|
* Accepts: - same arguments as MPI_Scatter()
|
|
|
|
* Returns: - MPI_SUCCESS or error code
|
|
|
|
*/
|
2005-08-10 14:51:42 +04:00
|
|
|
int
|
|
|
|
mca_coll_basic_scatter_intra(void *sbuf, int scount,
|
2005-08-10 21:53:43 +04:00
|
|
|
struct ompi_datatype_t *sdtype,
|
|
|
|
void *rbuf, int rcount,
|
|
|
|
struct ompi_datatype_t *rdtype,
|
|
|
|
int root, struct ompi_communicator_t *comm)
|
2004-01-12 00:26:55 +03:00
|
|
|
{
|
2005-08-10 14:51:42 +04:00
|
|
|
int i;
|
|
|
|
int rank;
|
|
|
|
int size;
|
|
|
|
int err;
|
|
|
|
char *ptmp;
|
|
|
|
long lb;
|
|
|
|
long incr;
|
|
|
|
|
|
|
|
/* Initialize */
|
|
|
|
|
|
|
|
rank = ompi_comm_rank(comm);
|
|
|
|
size = ompi_comm_size(comm);
|
|
|
|
|
|
|
|
/* If not root, receive data. */
|
|
|
|
|
|
|
|
if (rank != root) {
|
2005-08-10 21:53:43 +04:00
|
|
|
err = MCA_PML_CALL(recv(rbuf, rcount, rdtype, root,
|
|
|
|
MCA_COLL_BASE_TAG_SCATTER,
|
|
|
|
comm, MPI_STATUS_IGNORE));
|
|
|
|
return err;
|
2005-08-10 14:51:42 +04:00
|
|
|
}
|
2004-01-12 00:26:55 +03:00
|
|
|
|
2005-08-10 14:51:42 +04:00
|
|
|
/* I am the root, loop sending data. */
|
2004-01-12 00:26:55 +03:00
|
|
|
|
2005-08-10 14:51:42 +04:00
|
|
|
err = ompi_ddt_get_extent(rdtype, &lb, &incr);
|
|
|
|
if (OMPI_SUCCESS != err) {
|
2005-08-10 21:53:43 +04:00
|
|
|
return OMPI_ERROR;
|
2004-01-12 00:26:55 +03:00
|
|
|
}
|
2005-08-10 14:51:42 +04:00
|
|
|
|
|
|
|
incr *= scount;
|
|
|
|
for (i = 0, ptmp = (char *) sbuf; i < size; ++i, ptmp += incr) {
|
|
|
|
|
2005-08-10 21:53:43 +04:00
|
|
|
/* simple optimization */
|
|
|
|
|
|
|
|
if (i == rank) {
|
2005-09-13 22:12:10 +04:00
|
|
|
if (MPI_IN_PLACE != rbuf) {
|
2005-09-13 22:02:36 +04:00
|
|
|
err =
|
|
|
|
ompi_ddt_sndrcv(ptmp, scount, sdtype, rbuf, rcount,
|
|
|
|
rdtype);
|
|
|
|
}
|
2005-08-10 21:53:43 +04:00
|
|
|
} else {
|
|
|
|
err = MCA_PML_CALL(send(ptmp, scount, sdtype, i,
|
|
|
|
MCA_COLL_BASE_TAG_SCATTER,
|
|
|
|
MCA_PML_BASE_SEND_STANDARD, comm));
|
|
|
|
}
|
|
|
|
if (MPI_SUCCESS != err) {
|
|
|
|
return err;
|
|
|
|
}
|
2004-01-12 00:26:55 +03:00
|
|
|
}
|
|
|
|
|
2005-08-10 14:51:42 +04:00
|
|
|
/* All done */
|
2004-01-12 00:26:55 +03:00
|
|
|
|
2005-08-10 14:51:42 +04:00
|
|
|
return MPI_SUCCESS;
|
2004-06-29 04:02:25 +04:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
/*
|
|
|
|
* scatter_inter
|
|
|
|
*
|
|
|
|
* Function: - scatter operation
|
|
|
|
* Accepts: - same arguments as MPI_Scatter()
|
|
|
|
* Returns: - MPI_SUCCESS or error code
|
|
|
|
*/
|
2005-08-10 14:51:42 +04:00
|
|
|
int
|
|
|
|
mca_coll_basic_scatter_inter(void *sbuf, int scount,
|
2005-08-10 21:53:43 +04:00
|
|
|
struct ompi_datatype_t *sdtype,
|
|
|
|
void *rbuf, int rcount,
|
|
|
|
struct ompi_datatype_t *rdtype,
|
|
|
|
int root, struct ompi_communicator_t *comm)
|
2004-06-29 04:02:25 +04:00
|
|
|
{
|
2005-08-10 14:51:42 +04:00
|
|
|
int i;
|
|
|
|
int rank;
|
|
|
|
int size;
|
|
|
|
int err;
|
|
|
|
char *ptmp;
|
|
|
|
long lb;
|
|
|
|
long incr;
|
|
|
|
ompi_request_t **reqs = comm->c_coll_basic_data->mccb_reqs;
|
|
|
|
|
|
|
|
/* Initialize */
|
|
|
|
|
|
|
|
rank = ompi_comm_rank(comm);
|
|
|
|
size = ompi_comm_remote_size(comm);
|
|
|
|
|
|
|
|
if (MPI_PROC_NULL == root) {
|
2005-08-10 21:53:43 +04:00
|
|
|
/* do nothing */
|
|
|
|
err = OMPI_SUCCESS;
|
2005-08-10 14:51:42 +04:00
|
|
|
} else if (MPI_ROOT != root) {
|
2005-08-10 21:53:43 +04:00
|
|
|
/* If not root, receive data. */
|
|
|
|
err = MCA_PML_CALL(recv(rbuf, rcount, rdtype, root,
|
|
|
|
MCA_COLL_BASE_TAG_SCATTER,
|
|
|
|
comm, MPI_STATUS_IGNORE));
|
2005-08-10 14:51:42 +04:00
|
|
|
} else {
|
2005-08-10 21:53:43 +04:00
|
|
|
/* I am the root, loop sending data. */
|
|
|
|
err = ompi_ddt_get_extent(rdtype, &lb, &incr);
|
|
|
|
if (OMPI_SUCCESS != err) {
|
|
|
|
return OMPI_ERROR;
|
|
|
|
}
|
|
|
|
|
|
|
|
incr *= scount;
|
|
|
|
for (i = 0, ptmp = (char *) sbuf; i < size; ++i, ptmp += incr) {
|
|
|
|
err = MCA_PML_CALL(isend(ptmp, scount, sdtype, i,
|
|
|
|
MCA_COLL_BASE_TAG_SCATTER,
|
|
|
|
MCA_PML_BASE_SEND_STANDARD, comm,
|
|
|
|
reqs++));
|
|
|
|
if (OMPI_SUCCESS != err) {
|
|
|
|
return err;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
err =
|
|
|
|
ompi_request_wait_all(size, comm->c_coll_basic_data->mccb_reqs,
|
|
|
|
MPI_STATUSES_IGNORE);
|
2005-08-10 14:51:42 +04:00
|
|
|
}
|
|
|
|
|
|
|
|
return err;
|
2004-01-12 00:26:55 +03:00
|
|
|
}
|