From 3d42ec45c971024c0b85497e1b9105cc5a6dd77a Mon Sep 17 00:00:00 2001 From: Teemu Murtola Date: Fri, 13 Nov 2015 11:06:07 +0200 Subject: [PATCH] Remove t_commrec dependency from cmdlineinit.cpp This removes remaining legacyheaders dependencies from commandline/. Change-Id: Iff17ee731823fc4979231944a1f53595e8a9f7fc --- src/gromacs/commandline/cmdlineinit.cpp | 50 +++++++++++++++------------------ src/gromacs/utility/basenetwork.cpp | 24 +++++++++++----- src/gromacs/utility/basenetwork.h | 34 +++++----------------- 3 files changed, 47 insertions(+), 61 deletions(-) diff --git a/src/gromacs/commandline/cmdlineinit.cpp b/src/gromacs/commandline/cmdlineinit.cpp index 4b428acbe5..2915f0aacb 100644 --- a/src/gromacs/commandline/cmdlineinit.cpp +++ b/src/gromacs/commandline/cmdlineinit.cpp @@ -43,8 +43,6 @@ #include "cmdlineinit.h" -#include "config.h" - #include #include @@ -52,8 +50,7 @@ #include "gromacs/commandline/cmdlinemodulemanager.h" #include "gromacs/commandline/cmdlineoptionsmodule.h" #include "gromacs/commandline/cmdlineprogramcontext.h" -#include "gromacs/legacyheaders/network.h" -#include "gromacs/legacyheaders/types/commrec.h" +#include "gromacs/utility/basenetwork.h" #include "gromacs/utility/datafilefinder.h" #include "gromacs/utility/exceptions.h" #include "gromacs/utility/futil.h" @@ -77,31 +74,41 @@ std::unique_ptr g_commandLineContext; //! Global library data file finder that respects GMXLIB. std::unique_ptr g_libFileFinder; -#ifdef GMX_LIB_MPI -void broadcastArguments(const t_commrec *cr, int *argc, char ***argv) +/*! \brief + * Broadcasts command-line arguments to all ranks. + * + * MPI does not ensure that command-line arguments would be passed on any + * other rank than zero, but our code wants to parse them on each rank + * separately. + */ +void broadcastArguments(int *argc, char ***argv) { - gmx_bcast(sizeof(*argc), argc, cr); + if (gmx_node_num() <= 1) + { + return; + } + gmx_broadcast_world(sizeof(*argc), argc); - if (!MASTER(cr)) + const bool isMaster = (gmx_node_rank() == 0); + if (!isMaster) { snew(*argv, *argc+1); } for (int i = 0; i < *argc; i++) { int len; - if (MASTER(cr)) + if (isMaster) { len = std::strlen((*argv)[i])+1; } - gmx_bcast(sizeof(len), &len, cr); - if (!MASTER(cr)) + gmx_broadcast_world(sizeof(len), &len); + if (!isMaster) { snew((*argv)[i], len); } - gmx_bcast(len, (*argv)[i], cr); + gmx_broadcast_world(len, (*argv)[i]); } } -#endif //! \} @@ -112,20 +119,9 @@ CommandLineProgramContext &initForCommandLine(int *argc, char ***argv) gmx::init(argc, argv); GMX_RELEASE_ASSERT(!g_commandLineContext, "initForCommandLine() calls cannot be nested"); -#ifdef GMX_LIB_MPI - // TODO: Rewrite this to not use t_commrec once there is clarity on - // the approach for MPI in C++ code. // TODO: Consider whether the argument broadcast would better be done // in CommandLineModuleManager. - t_commrec cr; - std::memset(&cr, 0, sizeof(cr)); - - gmx_fill_commrec_from_mpi(&cr); - if (PAR(&cr)) - { - broadcastArguments(&cr, argc, argv); - } -#endif + broadcastArguments(argc, argv); try { g_commandLineContext.reset(new CommandLineProgramContext(*argc, *argv)); @@ -153,8 +149,8 @@ void finalizeForCommandLine() int processExceptionAtExitForCommandLine(const std::exception &ex) { - int rc = processExceptionAtExit(ex); //Currently this aborts for GMX_LIB_MPI - finalizeForCommandLine(); //thus this MPI_Finalize doesn't matter + int rc = processExceptionAtExit(ex); // Currently this aborts for real MPI + finalizeForCommandLine(); // thus this MPI_Finalize doesn't matter. return rc; } diff --git a/src/gromacs/utility/basenetwork.cpp b/src/gromacs/utility/basenetwork.cpp index 57e7794bff..38151a9a2e 100644 --- a/src/gromacs/utility/basenetwork.cpp +++ b/src/gromacs/utility/basenetwork.cpp @@ -49,7 +49,7 @@ #include "gromacs/utility/fatalerror.h" #include "gromacs/utility/gmxmpi.h" -gmx_bool gmx_mpi_initialized(void) +bool gmx_mpi_initialized() { #ifndef GMX_MPI return 0; @@ -61,7 +61,7 @@ gmx_bool gmx_mpi_initialized(void) #endif } -int gmx_node_num(void) +int gmx_node_num() { #ifndef GMX_MPI return 1; @@ -78,7 +78,7 @@ int gmx_node_num(void) #endif } -int gmx_node_rank(void) +int gmx_node_rank() { #ifndef GMX_MPI return 0; @@ -95,7 +95,7 @@ int gmx_node_rank(void) #endif } -static int mpi_hostname_hash(void) +static int mpi_hostname_hash() { int hash_int; @@ -141,7 +141,7 @@ static uint64_t Kernel_GetJobID(); #endif #include -static int bgq_nodenum(void) +static int bgq_nodenum() { int hostnum; Personality_t personality; @@ -153,7 +153,7 @@ static int bgq_nodenum(void) threads, so 0 <= T <= 63 (but the maximum value of T depends on the confituration of ranks and OpenMP threads per node). However, T is irrelevant for computing a suitable return - value for gmx_hostname_num(). + value for gmx_physicalnode_id_hash(). */ hostnum = personality.Network_Config.Acoord; hostnum *= personality.Network_Config.Bnodes; @@ -191,7 +191,7 @@ static int bgq_nodenum(void) } #endif -int gmx_physicalnode_id_hash(void) +int gmx_physicalnode_id_hash() { int hash; @@ -221,6 +221,16 @@ int gmx_physicalnode_id_hash(void) return hash; } +void gmx_broadcast_world(int size, void *buffer) +{ +#ifdef GMX_MPI + MPI_Bcast(buffer, size, MPI_BYTE, 0, MPI_COMM_WORLD); +#else + GMX_UNUSED_VALUE(size); + GMX_UNUSED_VALUE(buffer); +#endif +} + #ifdef GMX_LIB_MPI void gmx_abort(int errorno) { diff --git a/src/gromacs/utility/basenetwork.h b/src/gromacs/utility/basenetwork.h index a03d6f698f..79125f5e3b 100644 --- a/src/gromacs/utility/basenetwork.h +++ b/src/gromacs/utility/basenetwork.h @@ -3,7 +3,7 @@ * * Copyright (c) 1991-2000, University of Groningen, The Netherlands. * Copyright (c) 2001-2004, The GROMACS development team. - * Copyright (c) 2013,2014, by the GROMACS development team, led by + * Copyright (c) 2013,2014,2015, by the GROMACS development team, led by * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl, * and including many others, as listed in the AUTHORS file in the * top-level source directory and at http://www.gromacs.org. @@ -44,14 +44,6 @@ #ifndef GMX_UTILITY_BASENETWORK_H #define GMX_UTILITY_BASENETWORK_H -#include - -#include "gromacs/utility/basedefinitions.h" - -#ifdef __cplusplus -extern "C" { -#endif - /*! \brief * Returns whether MPI has been initialized. * @@ -63,7 +55,7 @@ extern "C" { * thread start where the return value is different depending on compilation * options. */ -gmx_bool gmx_mpi_initialized(void); +bool gmx_mpi_initialized(); /*! \brief * Returns the number of nodes. @@ -74,7 +66,7 @@ gmx_bool gmx_mpi_initialized(void); * expected: for thread-MPI, at that point they should behave as if the run was * serial. */ -int gmx_node_num(void); +int gmx_node_num(); /*! \brief * Returns the rank of the node. @@ -85,7 +77,7 @@ int gmx_node_num(void); * for thread-MPI, at that point the only thread of execution should behave as * if it the master node. */ -int gmx_node_rank(void); +int gmx_node_rank(); /*! \brief * Return a non-negative hash that is, hopefully, unique for each physical @@ -93,26 +85,14 @@ int gmx_node_rank(void); * * This hash is useful for determining hardware locality. */ -int gmx_physicalnode_id_hash(void); +int gmx_physicalnode_id_hash(); /*! \brief - * Returns an integer characteristic of and, hopefully, unique to each - * physical node in the MPI system. - * - * If the first part of the MPI hostname (up to the first dot) ends with a - * number, returns this number. If the first part of the MPI hostname does not - * ends in a number (0-9 characters), returns 0. - * - * \todo - * This function should be fully replaced by gmx_physicalnode_id_hash(). + * Broadcasts given data from rank zero to all other ranks. */ -int gmx_hostname_num(void); +void gmx_broadcast_world(int size, void *buffer); /** Abort the parallel run */ void gmx_abort(int errorno); -#ifdef __cplusplus -} -#endif - #endif -- 2.11.4.GIT