From d0baf6c96942da897df348d8b4e881694e980465 Mon Sep 17 00:00:00 2001 From: Mark Abraham Date: Wed, 5 Jul 2017 12:59:35 +0200 Subject: [PATCH] Simplify gmx_gpu_opt_t This temporary data about compatibility can be recomputed simply at need. Change-Id: I70352d75a39f4f777823638c462ea108829e46be --- src/gromacs/gpu_utils/gpu_utils.cu | 1 + src/gromacs/hardware/detecthardware.cpp | 15 ++++++---- src/gromacs/hardware/gpu_hw_info.h | 2 -- src/gromacs/hardware/hardwareassign.cpp | 49 ++++++++++++++------------------- src/gromacs/hardware/hardwareassign.h | 11 ++++++++ 5 files changed, 41 insertions(+), 37 deletions(-) diff --git a/src/gromacs/gpu_utils/gpu_utils.cu b/src/gromacs/gpu_utils/gpu_utils.cu index 9e34ed1e7b..c2f5527efb 100644 --- a/src/gromacs/gpu_utils/gpu_utils.cu +++ b/src/gromacs/gpu_utils/gpu_utils.cu @@ -638,6 +638,7 @@ bool isGpuCompatible(const gmx_gpu_info_t *gpu_info, int index) { assert(gpu_info); + assert(gpu_info->n_dev == 0 || gpu_info->gpu_dev); return (index >= gpu_info->n_dev ? false : diff --git a/src/gromacs/hardware/detecthardware.cpp b/src/gromacs/hardware/detecthardware.cpp index 5a47d0a002..ca4e38ca0a 100644 --- a/src/gromacs/hardware/detecthardware.cpp +++ b/src/gromacs/hardware/detecthardware.cpp @@ -54,6 +54,7 @@ #include "gromacs/gpu_utils/gpu_utils.h" #include "gromacs/hardware/cpuinfo.h" #include "gromacs/hardware/gpu_hw_info.h" +#include "gromacs/hardware/hardwareassign.h" #include "gromacs/hardware/hardwaretopology.h" #include "gromacs/hardware/hw_info.h" #include "gromacs/mdtypes/commrec.h" @@ -167,6 +168,8 @@ std::string sprint_gpus(const gmx_gpu_info_t *gpu_info) return gmx::joinStrings(gpuStrings, "\n"); } +// TODO This function should not live in detectharware.cpp + /*! \brief Helper function for reporting GPU usage information * in the mdrun log file * @@ -204,19 +207,19 @@ makeGpuUsageReport(const gmx_gpu_info_t *gpu_info, std::string output; if (!userSetGpuIds) { - // gpu_opt->dev_compatible is only populated during auto-selection - std::string gpuIdsString = - formatAndJoin(gmx::constArrayRefFromArray(gpu_opt->dev_compatible, - gpu_opt->n_dev_compatible), + auto compatibleGpus = getCompatibleGpus(gpu_info); + int numCompatibleGpus = static_cast(compatibleGpus.size()); + std::string gpuIdsString = + formatAndJoin(compatibleGpus, ",", gmx::StringFormatter("%d")); - bool bPluralGpus = gpu_opt->n_dev_compatible > 1; + bool bPluralGpus = numCompatibleGpus > 1; if (bPrintHostName) { output += gmx::formatString("On host %s ", host); } output += gmx::formatString("%d compatible GPU%s %s present, with ID%s %s\n", - gpu_opt->n_dev_compatible, + numCompatibleGpus, bPluralGpus ? "s" : "", bPluralGpus ? "are" : "is", bPluralGpus ? "s" : "", diff --git a/src/gromacs/hardware/gpu_hw_info.h b/src/gromacs/hardware/gpu_hw_info.h index 04921ef934..f2c2c64059 100644 --- a/src/gromacs/hardware/gpu_hw_info.h +++ b/src/gromacs/hardware/gpu_hw_info.h @@ -74,8 +74,6 @@ typedef struct gmx_gpu_opt_t { char *gpu_id; /* GPU id's to use, each specified as chars */ - int n_dev_compatible; /* number of compatible GPU devices that could be used */ - int *dev_compatible; /* array of compatible GPU device IDs, from which automatic selection occurs */ int n_dev_use; /* number of GPU devices selected to be used, either by the user or automatically */ int *dev_use; /* array mapping from PP rank index to GPU device ID; GPU IDs can be listed multiple times when ranks share them */ } gmx_gpu_opt_t; diff --git a/src/gromacs/hardware/hardwareassign.cpp b/src/gromacs/hardware/hardwareassign.cpp index 0ac9f8f997..6366ed1146 100644 --- a/src/gromacs/hardware/hardwareassign.cpp +++ b/src/gromacs/hardware/hardwareassign.cpp @@ -101,6 +101,7 @@ static void print_gpu_detection_stats(const gmx::MDLogger &mdlog, * This function is responsible for mapping the GPUs to the processes on a single node * (filling the gpu_opt->dev_use array). * + * \param[in] compatibleGpus Vector of GPUs that are compatible * \param[in,out] gpu_opt Input/output GPU assignment data. * \param[in] nrank Number of PP GPU ranks on the node. * \param[in] rank Index of PP GPU rank on the node. @@ -108,14 +109,16 @@ static void print_gpu_detection_stats(const gmx::MDLogger &mdlog, * This selects the GPUs we will use. This is an operation local to each physical node. * If we have less MPI ranks than GPUs, we will waste some GPUs. */ -static void assign_rank_gpu_ids(gmx_gpu_opt_t *gpu_opt, int nrank, int rank) +static void assign_rank_gpu_ids(const std::vector &compatibleGpus, + gmx_gpu_opt_t *gpu_opt, int nrank, int rank) { + int numCompatibleGpus = static_cast(compatibleGpus.size()); GMX_RELEASE_ASSERT(gpu_opt, "Invalid gpu_opt pointer passed"); GMX_RELEASE_ASSERT(nrank >= 1, gmx::formatString("Invalid limit (%d) for the number of GPUs (detected %d compatible GPUs)", - rank, gpu_opt->n_dev_compatible).c_str()); + rank, numCompatibleGpus).c_str()); - if (gpu_opt->n_dev_compatible == 0) + if (numCompatibleGpus == 0) { char host[HOSTNAMELEN]; @@ -126,18 +129,18 @@ static void assign_rank_gpu_ids(gmx_gpu_opt_t *gpu_opt, int nrank, int rank) int nshare; nshare = 1; - if (nrank > gpu_opt->n_dev_compatible) + if (nrank > numCompatibleGpus) { - if (nrank % gpu_opt->n_dev_compatible == 0) + if (nrank % numCompatibleGpus == 0) { - nshare = gmx_gpu_sharing_supported() ? nrank/gpu_opt->n_dev_compatible : 1; + nshare = gmx_gpu_sharing_supported() ? nrank/numCompatibleGpus : 1; } else { if (rank == 0) { gmx_fatal(FARGS, "The number of MPI ranks (%d) in a physical node is not a multiple of the number of GPUs (%d). Select a different number of MPI ranks or use the -gpu_id option to manually specify the GPU to be used.", - nrank, gpu_opt->n_dev_compatible); + nrank, numCompatibleGpus); } #if GMX_MPI @@ -149,8 +152,8 @@ static void assign_rank_gpu_ids(gmx_gpu_opt_t *gpu_opt, int nrank, int rank) } } - /* Here we will waste GPUs when nrank < gpu_opt->n_dev_compatible */ - gpu_opt->n_dev_use = std::min(gpu_opt->n_dev_compatible*nshare, nrank); + /* Here we will waste GPUs when nrank < numCompatibleGpus */ + gpu_opt->n_dev_use = std::min(numCompatibleGpus*nshare, nrank); if (!gmx_multiple_gpu_per_node_supported()) { gpu_opt->n_dev_use = std::min(gpu_opt->n_dev_use, 1); @@ -159,7 +162,7 @@ static void assign_rank_gpu_ids(gmx_gpu_opt_t *gpu_opt, int nrank, int rank) for (int i = 0; i != gpu_opt->n_dev_use; ++i) { /* TODO: improve this implementation: either sort GPUs or remove the weakest here */ - gpu_opt->dev_use[i] = gpu_opt->dev_compatible[i/nshare]; + gpu_opt->dev_use[i] = compatibleGpus[i/nshare]; } } @@ -203,33 +206,21 @@ static bool checkGpuSelection(const gmx_gpu_info_t *gpu_info, return allOK; } -/*! \brief Select the compatible GPUs - * - * This function filters gpu_info->gpu_dev for compatible gpus based - * on the previously run compatibility tests. Sets - * gpu_info->dev_compatible and gpu_info->n_dev_compatible. - * - * \param[in] gpu_info pointer to structure holding GPU information - * \param[out] gpu_opt pointer to structure holding GPU options - */ -static void pickCompatibleGpus(const gmx_gpu_info_t *gpu_info, - gmx_gpu_opt_t *gpu_opt) +std::vector getCompatibleGpus(const gmx_gpu_info_t *gpu_info) { GMX_ASSERT(gpu_info, "Invalid gpu_info"); - GMX_ASSERT(gpu_opt, "Invalid gpu_opt"); // Possible minor over-allocation here, but not important for anything - gpu_opt->n_dev_compatible = 0; - snew(gpu_opt->dev_compatible, gpu_info->n_dev); + std::vector compatibleGpus; + compatibleGpus.reserve(gpu_info->n_dev); for (int i = 0; i < gpu_info->n_dev; i++) { - GMX_ASSERT(gpu_info->gpu_dev, "Invalid gpu_info->gpu_dev"); if (isGpuCompatible(gpu_info, i)) { - gpu_opt->dev_compatible[gpu_opt->n_dev_compatible] = i; - gpu_opt->n_dev_compatible++; + compatibleGpus.push_back(i); } } + return compatibleGpus; } void gmx_select_rank_gpu_ids(const gmx::MDLogger &mdlog, const t_commrec *cr, @@ -261,7 +252,7 @@ void gmx_select_rank_gpu_ids(const gmx::MDLogger &mdlog, const t_commrec *cr, } else { - pickCompatibleGpus(gpu_info, gpu_opt); - assign_rank_gpu_ids(gpu_opt, cr->nrank_pp_intranode, cr->rank_pp_intranode); + auto compatibleGpus = getCompatibleGpus(gpu_info); + assign_rank_gpu_ids(compatibleGpus, gpu_opt, cr->nrank_pp_intranode, cr->rank_pp_intranode); } } diff --git a/src/gromacs/hardware/hardwareassign.h b/src/gromacs/hardware/hardwareassign.h index b61dda0ca4..02495e6136 100644 --- a/src/gromacs/hardware/hardwareassign.h +++ b/src/gromacs/hardware/hardwareassign.h @@ -35,6 +35,8 @@ #ifndef GMX_HARDWARE_HARDWAREASSIGN_H #define GMX_HARDWARE_HARDWAREASSIGN_H +#include + #include "gromacs/utility/basedefinitions.h" struct gmx_gpu_info_t; @@ -46,6 +48,15 @@ namespace gmx class MDLogger; } +/*! \brief Select the compatible GPUs + * + * This function filters gpu_info->gpu_dev for compatible GPUs based + * on the previously run compatibility tests. + * + * \param[in] gpu_info pointer to structure holding GPU information, including compatibility + * \return vector of IDs of GPUs already recorded as compatible */ +std::vector getCompatibleGpus(const gmx_gpu_info_t *gpu_info); + void gmx_select_rank_gpu_ids(const gmx::MDLogger &mdlog, const t_commrec *cr, const gmx_gpu_info_t *gpu_info, bool userSetGpuIds, -- 2.11.4.GIT