Parse user-supplied GPU task assignment only when needed
[gromacs.git] / src / gromacs / hardware / hardwareassign.cpp
blobfda0b96c69837d9d3adbadc2ef57cfe27deaaa43
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 #include "gmxpre.h"
37 #include "hardwareassign.h"
39 #include "config.h"
41 #include <cstring>
43 #include <algorithm>
44 #include <string>
45 #include <vector>
47 #include "gromacs/gmxlib/network.h"
48 #include "gromacs/gpu_utils/gpu_utils.h"
49 #include "gromacs/hardware/gpu_hw_info.h"
50 #include "gromacs/hardware/hw_info.h"
51 #include "gromacs/mdtypes/commrec.h"
52 #include "gromacs/utility/basenetwork.h"
53 #include "gromacs/utility/exceptions.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/gmxassert.h"
56 #include "gromacs/utility/logger.h"
57 #include "gromacs/utility/smalloc.h"
58 #include "gromacs/utility/stringutil.h"
59 #include "gromacs/utility/sysinfo.h"
61 #define HOSTNAMELEN 80
63 namespace gmx
66 std::vector<int> parseGpuTaskAssignment(const std::string &gpuTaskAssignment)
68 std::vector<int> digits;
69 if (gpuTaskAssignment.empty())
71 return digits;
74 /* Parse a "plain" or comma-separated GPU ID string which contains
75 * a sequence of digits corresponding to GPU IDs; the order will
76 * indicate the assignment of GPU tasks on this node to GPU
77 * device IDs on this node. */
78 try
80 digits = parseDigitsFromString(gpuTaskAssignment);
82 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
84 if (digits.empty())
86 gmx_fatal(FARGS, "Empty GPU ID string encountered.\n"
87 "An empty, delimiter-free, or comma-separated sequence of valid numeric IDs of available GPUs is required.\n");
89 return digits;
92 /*! \brief This function is responsible for the automated mapping the
93 * GPUs to the processes on a single node.
95 * This selects the GPUs we will use. This is an operation local to each physical node.
96 * If we have less MPI ranks than GPUs, we will waste some GPUs.
98 * \param[in] compatibleGpus Vector of GPUs that are compatible
99 * \param[in] nrank Number of PP GPU ranks on the node.
100 * \param[in] rank Index of PP GPU rank on the node.
102 * \returns The assignment of GPU tasks on ranks of this node to GPU devices on this node.
104 static std::vector<int> assign_rank_gpu_ids(const std::vector<int> &compatibleGpus,
105 int nrank, int rank)
107 int numCompatibleGpus = static_cast<int>(compatibleGpus.size());
108 GMX_RELEASE_ASSERT(nrank >= 1,
109 gmx::formatString("Invalid limit (%d) for the number of GPUs (detected %d compatible GPUs)",
110 rank, numCompatibleGpus).c_str());
112 if (numCompatibleGpus == 0)
114 char host[HOSTNAMELEN];
116 gmx_gethostname(host, HOSTNAMELEN);
117 gmx_fatal(FARGS, "A GPU was requested on host %s, but no compatible GPUs were detected. All nodes with PP ranks need to have GPUs. If you intended to use GPU acceleration in a parallel run, you can either avoid using the nodes that don't have GPUs or place PME ranks on these nodes.", host);
120 int nshare;
122 nshare = 1;
123 if (nrank > numCompatibleGpus)
125 if (nrank % numCompatibleGpus == 0)
127 nshare = nrank/numCompatibleGpus;
129 else
131 if (rank == 0)
133 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.",
134 nrank, numCompatibleGpus);
137 #if GMX_MPI
138 /* We use a global barrier to prevent ranks from continuing with
139 * an invalid setup.
141 MPI_Barrier(MPI_COMM_WORLD);
142 #endif
146 /* Here we will waste GPUs when nrank < numCompatibleGpus */
147 std::vector<int> taskAssignment;
148 taskAssignment.resize(std::min(numCompatibleGpus*nshare, nrank));
149 for (size_t i = 0; i != taskAssignment.size(); ++i)
151 /* TODO: improve this implementation: either sort GPUs or remove the weakest here */
152 taskAssignment[i] = compatibleGpus[i/nshare];
154 return taskAssignment;
157 /*! \brief Check that all user-selected GPUs are compatible.
159 * Given the \c userGpuTaskAssignment and \c compatibleGPUs, give a fatal
160 * error if any selected GPUs is not compatible
162 * The error is given with a suitable descriptive message, which will
163 * have context if this check is done after the hardware detection
164 * results have been reported to the user. However, note that only the
165 * GPUs detected on the master rank are reported, because of the
166 * existing limitations of that reporting.
168 * \todo Note that the selected GPUs can be different on each rank,
169 * and the IDs of compatible GPUs can be different on each node, so
170 * this routine ought to do communication to determine whether all
171 * ranks are able to proceed. Currently this relies on the MPI runtime
172 * to kill the other processes because GROMACS lacks the appropriate
173 * infrastructure to do a good job of coordinating error messages and
174 * behaviour across MPMD ranks and multiple simulations.
176 * \param[in] gpu_info GPU information including device description.
177 * \param[in] compatibleGpus Vector of compatible GPUs
178 * \param[in] userGpuTaskAssignment The GPU selection from the user.
180 static void exitUnlessUserGpuTaskAssignmentIsValid(const gmx_gpu_info_t &gpu_info,
181 const std::vector<int> &compatibleGpus,
182 const std::vector<int> &userGpuTaskAssignment)
184 int numIncompatibleGpuIds = 0;
185 std::string message
186 = "Some of the requested GPUs do not exist, behave strangely, or are not compatible:\n";
188 for (const auto &gpuId : userGpuTaskAssignment)
190 if (std::find(compatibleGpus.begin(), compatibleGpus.end(), gpuId) == compatibleGpus.end())
192 numIncompatibleGpuIds++;
193 message += gmx::formatString(" GPU #%d: %s\n",
194 gpuId,
195 getGpuCompatibilityDescription(gpu_info, gpuId));
199 if (numIncompatibleGpuIds > 0)
201 gmx_fatal(FARGS, message.c_str());
205 /*! \brief Filter the compatible GPUs
207 * This function filters gpu_info.gpu_dev for compatible GPUs based
208 * on the previously run compatibility tests.
210 * \param[in] gpu_info Information detected about GPUs, including compatibility
211 * \return vector of IDs of GPUs already recorded as compatible */
212 static std::vector<int> getCompatibleGpus(const gmx_gpu_info_t &gpu_info)
214 // Possible minor over-allocation here, but not important for anything
215 std::vector<int> compatibleGpus;
216 compatibleGpus.reserve(gpu_info.n_dev);
217 for (int i = 0; i < gpu_info.n_dev; i++)
219 GMX_ASSERT(gpu_info.gpu_dev, "Invalid gpu_info.gpu_dev");
220 if (isGpuCompatible(gpu_info, i))
222 compatibleGpus.push_back(i);
225 return compatibleGpus;
228 std::vector<int> mapPpRanksToGpus(bool rankCanUseGpu,
229 const t_commrec *cr,
230 const gmx_gpu_info_t &gpu_info,
231 const gmx_hw_opt_t &hw_opt)
233 std::vector<int> taskAssignment;
235 if (!rankCanUseGpu)
237 return taskAssignment;
240 auto compatibleGpus = getCompatibleGpus(gpu_info);
241 if (!hw_opt.gpuIdTaskAssignment.empty())
243 auto userGpuTaskAssignment = parseGpuTaskAssignment(hw_opt.gpuIdTaskAssignment);
244 exitUnlessUserGpuTaskAssignmentIsValid(gpu_info, compatibleGpus, userGpuTaskAssignment);
245 taskAssignment = userGpuTaskAssignment;
247 else
249 taskAssignment = assign_rank_gpu_ids(compatibleGpus, cr->nrank_pp_intranode, cr->rank_pp_intranode);
251 return taskAssignment;
254 } // namespace