1
1
openmpi/ompi/mpi/c/scatter.c
Edgar Gabriel 5989fa570c Sorry, previous commit was in the wrong directory. This is the real fix (have
to undo 1822).

The verification of recvcount==0 and rank = root was braking
inter-communicator scatter, since the root (root==MPI_ROOT) might very well
have recvcount=0. The same fix has been applied to gather.c just the other way
round. 
 
Fixes the bug reported on the mainling list by Martin Audet. If there is a
1.2.7 this fix might be worthwhile porting it over.

Please note, that while the test works now for basic and for inter, we get a
0byte malloc warning from the inter module, which we still have to fix in a
separate patch.

This commit was SVN r18123.
2008-04-10 15:03:14 +00:00

139 строки
4.7 KiB
C

/*
* Copyright (c) 2004-2007 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-2008 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-2007 Cisco Systems, Inc. All rights reserved.
* Copyright (c) 2008 University of Houston. 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"
#include "ompi/memchecker.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;
MEMCHECKER(
memchecker_datatype(sendtype);
memchecker_datatype(recvtype);
memchecker_call(&opal_memchecker_base_isdefined, sendbuf, sendcount, sendtype);
memchecker_comm(comm);
);
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 || 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 && MPI_ROOT != root &&
(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;
}
OPAL_CR_ENTER_LIBRARY();
/* Invoke the coll component to perform the back-end operation */
err = comm->c_coll.coll_scatter(sendbuf, sendcount, sendtype, recvbuf,
recvcount, recvtype, root, comm,
comm->c_coll.coll_scatter_module);
OMPI_ERRHANDLER_RETURN(err, comm, err, FUNC_NAME);
}