root/ompi/mpi/c/scatter.c

/* [<][>][^][v][top][bottom][index][help] */

DEFINITIONS

This source file includes following definitions.
  1. MPI_Scatter

   1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil -*- */
   2 /*
   3  * Copyright (c) 2004-2007 The Trustees of Indiana University and Indiana
   4  *                         University Research and Technology
   5  *                         Corporation.  All rights reserved.
   6  * Copyright (c) 2004-2018 The University of Tennessee and The University
   7  *                         of Tennessee Research Foundation.  All rights
   8  *                         reserved.
   9  * Copyright (c) 2004-2008 High Performance Computing Center Stuttgart,
  10  *                         University of Stuttgart.  All rights reserved.
  11  * Copyright (c) 2004-2005 The Regents of the University of California.
  12  *                         All rights reserved.
  13  * Copyright (c) 2006-2007 Cisco Systems, Inc.  All rights reserved.
  14  * Copyright (c) 2008      University of Houston.  All rights reserved.
  15  * Copyright (c) 2008      Sun Microsystems, Inc.  All rights reserved.
  16  * Copyright (c) 2013      Los Alamos National Security, LLC.  All rights
  17  *                         reserved.
  18  * Copyright (c) 2015      Research Organization for Information Science
  19  *                         and Technology (RIST). All rights reserved.
  20  * $COPYRIGHT$
  21  *
  22  * Additional copyrights may follow
  23  *
  24  * $HEADER$
  25  */
  26 #include "ompi_config.h"
  27 #include <stdio.h>
  28 
  29 #include "ompi/mpi/c/bindings.h"
  30 #include "ompi/runtime/params.h"
  31 #include "ompi/communicator/communicator.h"
  32 #include "ompi/errhandler/errhandler.h"
  33 #include "ompi/datatype/ompi_datatype.h"
  34 #include "ompi/memchecker.h"
  35 #include "ompi/runtime/ompi_spc.h"
  36 
  37 #if OMPI_BUILD_MPI_PROFILING
  38 #if OPAL_HAVE_WEAK_SYMBOLS
  39 #pragma weak MPI_Scatter = PMPI_Scatter
  40 #endif
  41 #define MPI_Scatter PMPI_Scatter
  42 #endif
  43 
  44 static const char FUNC_NAME[] = "MPI_Scatter";
  45 
  46 
  47 int MPI_Scatter(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
  48                 void *recvbuf, int recvcount, MPI_Datatype recvtype,
  49                 int root, MPI_Comm comm)
  50 {
  51     int err;
  52 
  53     SPC_RECORD(OMPI_SPC_SCATTER, 1);
  54 
  55     MEMCHECKER(
  56         memchecker_comm(comm);
  57         if(OMPI_COMM_IS_INTRA(comm)) {
  58             if(ompi_comm_rank(comm) == root) {
  59                 memchecker_datatype(sendtype);
  60                 /* check whether root's send buffer is defined. */
  61                 memchecker_call(&opal_memchecker_base_isdefined, sendbuf, sendcount, sendtype);
  62                 if(MPI_IN_PLACE != recvbuf) {
  63                     memchecker_datatype(recvtype);
  64                     /* check whether receive buffer is addressable. */
  65                     memchecker_call(&opal_memchecker_base_isaddressable, recvbuf, recvcount, recvtype);
  66                 }
  67             } else {
  68                 memchecker_datatype(recvtype);
  69                 /* check whether receive buffer is addressable. */
  70                 memchecker_call(&opal_memchecker_base_isaddressable, recvbuf, recvcount, recvtype);
  71             }
  72         } else {
  73             if(MPI_ROOT == root) {
  74                 memchecker_datatype(sendtype);
  75                 /* check whether root's send buffer is defined. */
  76                 memchecker_call(&opal_memchecker_base_isdefined, sendbuf, sendcount, sendtype);
  77             } else if (MPI_PROC_NULL != root) {
  78                 memchecker_datatype(recvtype);
  79                 /* check whether receive buffer is addressable. */
  80                 memchecker_call(&opal_memchecker_base_isaddressable, recvbuf, recvcount, recvtype);
  81             }
  82         }
  83     );
  84 
  85     if (MPI_PARAM_CHECK) {
  86         err = MPI_SUCCESS;
  87         OMPI_ERR_INIT_FINALIZE(FUNC_NAME);
  88         if (ompi_comm_invalid(comm)) {
  89             return OMPI_ERRHANDLER_INVOKE(MPI_COMM_WORLD, MPI_ERR_COMM,
  90                                           FUNC_NAME);
  91         } else if ((ompi_comm_rank(comm) != root && MPI_IN_PLACE == recvbuf) ||
  92                    (ompi_comm_rank(comm) == root && MPI_IN_PLACE == sendbuf)) {
  93             return OMPI_ERRHANDLER_INVOKE(comm, MPI_ERR_ARG, FUNC_NAME);
  94         }
  95 
  96         /* Errors for intracommunicators */
  97 
  98         if (OMPI_COMM_IS_INTRA(comm)) {
  99 
 100           /* Errors for all ranks */
 101 
 102           if ((root >= ompi_comm_size(comm)) || (root < 0)) {
 103               err = MPI_ERR_ROOT;
 104           } else if (MPI_IN_PLACE != recvbuf) {
 105               if (recvcount < 0) {
 106                   err = MPI_ERR_COUNT;
 107               } else if (MPI_DATATYPE_NULL == recvtype || NULL == recvtype) {
 108                   err = MPI_ERR_TYPE;
 109               }
 110           }
 111 
 112           /* Errors for the root.  Some of these could have been
 113              combined into compound if statements above, but since
 114              this whole section can be compiled out (or turned off at
 115              run time) for efficiency, it's more clear to separate
 116              them out into individual tests. */
 117 
 118           else if (ompi_comm_rank(comm) == root) {
 119               OMPI_CHECK_DATATYPE_FOR_SEND(err, sendtype, sendcount);
 120           }
 121           OMPI_ERRHANDLER_CHECK(err, comm, err, FUNC_NAME);
 122         }
 123 
 124         /* Errors for intercommunicators */
 125 
 126         else {
 127           if (! ((root >= 0 && root < ompi_comm_remote_size(comm)) ||
 128                  MPI_ROOT == root || MPI_PROC_NULL == root)) {
 129               err = MPI_ERR_ROOT;
 130           }
 131 
 132           /* Errors for the receivers */
 133 
 134           else if (MPI_ROOT != root && MPI_PROC_NULL != root) {
 135             if (recvcount < 0) {
 136                 err = MPI_ERR_COUNT;
 137             } else if (MPI_DATATYPE_NULL == recvtype) {
 138                 err = MPI_ERR_TYPE;
 139             }
 140           }
 141 
 142           /* Errors for the root.  Ditto on the comment above -- these
 143              error checks could have been combined above, but let's
 144              make the code easier to read. */
 145 
 146           else if (MPI_ROOT == root) {
 147               OMPI_CHECK_DATATYPE_FOR_SEND(err, sendtype, sendcount);
 148           }
 149           OMPI_ERRHANDLER_CHECK(err, comm, err, FUNC_NAME);
 150         }
 151     }
 152 
 153     /* Do we need to do anything? */
 154 
 155     if ((0 == recvcount && MPI_ROOT != root &&
 156          (ompi_comm_rank(comm) != root ||
 157           (ompi_comm_rank(comm) == root && MPI_IN_PLACE != recvbuf))) ||
 158         (ompi_comm_rank(comm) == root && MPI_IN_PLACE == recvbuf &&
 159          0 == sendcount) ||
 160         (0 == sendcount && (MPI_ROOT == root || MPI_PROC_NULL == root))) {
 161         return MPI_SUCCESS;
 162     }
 163 
 164     OPAL_CR_ENTER_LIBRARY();
 165 
 166     /* Invoke the coll component to perform the back-end operation */
 167     err = comm->c_coll->coll_scatter(sendbuf, sendcount, sendtype, recvbuf,
 168                                     recvcount, recvtype, root, comm,
 169                                     comm->c_coll->coll_scatter_module);
 170     OMPI_ERRHANDLER_RETURN(err, comm, err, FUNC_NAME);
 171 }

/* [<][>][^][v][top][bottom][index][help] */