1
1
openmpi/ompi/mpi/c/scatter.c
Jeff Squyres e08c6e81f5 * Fixes trac:338: Only look at root-significant values at the root (e.g.,
recvbuf in MPI_GATHER).
 * Minor style updates (constants on the left of == and !=)
 * Fix a minor buglet that crept in r11904: had a recvbuf where it
   should have been recvcount.  Thankfully, this would have only
   affected erroneous programs.  ;-)

This commit was SVN r11980.

The following SVN revision numbers were found above:
  r11904 --> open-mpi/ompi@17539dc154

The following Trac tickets were found above:
  Ticket 338 --> https://svn.open-mpi.org/trac/ompi/ticket/338
2006-10-04 22:36:01 +00:00

127 строки
4.3 KiB
C

/*
* 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) 2006 Cisco Systems, Inc. All rights reserved.
* $COPYRIGHT$
*
* Additional copyrights may follow
*
* $HEADER$
*/
#include "ompi_config.h"
#include <stdio.h>
#include "ompi/mpi/c/bindings.h"
/*#include "ompi/mca/coll/coll.h"*/
#include "ompi/datatype/datatype.h"
#if OMPI_HAVE_WEAK_SYMBOLS && OMPI_PROFILING_DEFINES
#pragma weak MPI_Scatter = PMPI_Scatter
#endif
#if OMPI_PROFILING_DEFINES
#include "ompi/mpi/c/profile/defines.h"
#endif
static const char FUNC_NAME[] = "MPI_Scatter";
int MPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
void *recvbuf, int recvcount, MPI_Datatype recvtype,
int root, MPI_Comm comm)
{
int err;
if (MPI_PARAM_CHECK) {
err = MPI_SUCCESS;
OMPI_ERR_INIT_FINALIZE(FUNC_NAME);
if (ompi_comm_invalid(comm)) {
return OMPI_ERRHANDLER_INVOKE(MPI_COMM_WORLD, MPI_ERR_COMM,
FUNC_NAME);
} else if ((ompi_comm_rank(comm) != root && MPI_IN_PLACE == recvbuf) ||
(ompi_comm_rank(comm) == root && MPI_IN_PLACE == sendbuf)) {
return OMPI_ERRHANDLER_INVOKE(comm, MPI_ERR_ARG, FUNC_NAME);
}
/* Errors for intracommunicators */
if (OMPI_COMM_IS_INTRA(comm)) {
/* Errors for all ranks */
if ((root >= ompi_comm_size(comm)) || (root < 0)) {
err = MPI_ERR_ROOT;
} else if (MPI_IN_PLACE != recvbuf) {
if (recvcount < 0) {
err = MPI_ERR_COUNT;
} else if (MPI_DATATYPE_NULL == recvtype) {
err = MPI_ERR_TYPE;
}
}
/* Errors for the root. Some of these could have been
combined into compound if statements above, but since
this whole section can be compiled out (or turned off at
run time) for efficiency, it's more clear to separate
them out into individual tests. */
else if (ompi_comm_rank(comm) == root) {
OMPI_CHECK_DATATYPE_FOR_SEND(err, sendtype, sendcount);
}
OMPI_ERRHANDLER_CHECK(err, comm, err, FUNC_NAME);
}
/* Errors for intercommunicators */
else {
if (! ((root >= 0 && root < ompi_comm_remote_size(comm)) ||
MPI_ROOT == root || MPI_PROC_NULL == root)) {
err = MPI_ERR_ROOT;
}
/* Errors for the receivers */
else if (MPI_ROOT != root && MPI_PROC_NULL != root) {
if (recvcount < 0) {
err = MPI_ERR_COUNT;
} else if (MPI_DATATYPE_NULL == recvtype) {
err = MPI_ERR_TYPE;
}
}
/* Errors for the root. Ditto on the comment above -- these
error checks could have been combined above, but let's
make the code easier to read. */
else if (MPI_ROOT == root) {
OMPI_CHECK_DATATYPE_FOR_SEND(err, sendtype, sendcount);
}
OMPI_ERRHANDLER_CHECK(err, comm, err, FUNC_NAME);
}
}
/* Do we need to do anything? */
if ((0 == recvcount &&
(ompi_comm_rank(comm) != root ||
(ompi_comm_rank(comm) == root && MPI_IN_PLACE != recvbuf))) ||
(ompi_comm_rank(comm) == root && MPI_IN_PLACE == recvbuf &&
0 == sendcount)) {
return MPI_SUCCESS;
}
/* Invoke the coll component to perform the back-end operation */
err = comm->c_coll.coll_scatter(sendbuf, sendcount, sendtype, recvbuf,
recvcount, recvtype, root, comm);
OMPI_ERRHANDLER_RETURN(err, comm, err, FUNC_NAME);
}