ef9cb744bbbabd541c7b61a81d150e8b879a32cf
[alexxy/gromacs.git] / src / gromacs / taskassignment / hardwareassign.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #include "gmxpre.h"
36
37 #include "hardwareassign.h"
38
39 #include "config.h"
40
41 #include <cstring>
42
43 #include <algorithm>
44 #include <set>
45 #include <string>
46 #include <vector>
47
48 #include "gromacs/gmxlib/network.h"
49 #include "gromacs/gpu_utils/gpu_utils.h"
50 #include "gromacs/hardware/gpu_hw_info.h"
51 #include "gromacs/hardware/hw_info.h"
52 #include "gromacs/mdtypes/commrec.h"
53 #include "gromacs/utility/basenetwork.h"
54 #include "gromacs/utility/cstringutil.h"
55 #include "gromacs/utility/exceptions.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/gmxassert.h"
58 #include "gromacs/utility/logger.h"
59 #include "gromacs/utility/smalloc.h"
60 #include "gromacs/utility/stringutil.h"
61 #include "gromacs/utility/sysinfo.h"
62
63 #define HOSTNAMELEN 80
64
65 namespace gmx
66 {
67
68 std::vector<int> parseGpuTaskAssignment(const std::string &gpuTaskAssignment)
69 {
70     std::vector<int> digits;
71     if (gpuTaskAssignment.empty())
72     {
73         return digits;
74     }
75
76     /* Parse a "plain" or comma-separated GPU ID string which contains
77      * a sequence of digits corresponding to GPU IDs; the order will
78      * indicate the assignment of GPU tasks on this node to GPU
79      * device IDs on this node. */
80     try
81     {
82         digits = parseDigitsFromString(gpuTaskAssignment);
83     }
84     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
85
86     if (digits.empty())
87     {
88         gmx_fatal(FARGS, "Empty GPU ID string encountered.\n"
89                   "An empty, delimiter-free, or comma-separated sequence of valid numeric IDs of available GPUs is required.\n");
90     }
91     return digits;
92 }
93
94 /*! \brief This function is responsible for the automated mapping the
95  * GPUs to the processes on a single node.
96  *
97  * This selects the GPUs we will use. This is an operation local to each physical node.
98  * If we have less MPI ranks than GPUs, we will waste some GPUs.
99  *
100  * \param[in]        compatibleGpus       Vector of GPUs that are compatible
101  * \param[in]        nrank                Number of PP GPU ranks on the node.
102  * \param[in]        rank                 Index of PP GPU rank on the node.
103  *
104  * \returns The assignment of GPU tasks on ranks of this node to GPU devices on this node.
105  */
106 static std::vector<int> assign_rank_gpu_ids(const std::vector<int> &compatibleGpus,
107                                             int nrank, int rank)
108 {
109     int numCompatibleGpus = static_cast<int>(compatibleGpus.size());
110     GMX_RELEASE_ASSERT(nrank >= 1,
111                        gmx::formatString("Invalid limit (%d) for the number of GPUs (detected %d compatible GPUs)",
112                                          rank, numCompatibleGpus).c_str());
113
114     if (numCompatibleGpus == 0)
115     {
116         char host[HOSTNAMELEN];
117
118         gmx_gethostname(host, HOSTNAMELEN);
119         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     }
121
122     int nshare;
123
124     nshare = 1;
125     if (nrank > numCompatibleGpus)
126     {
127         if (nrank % numCompatibleGpus == 0)
128         {
129             nshare = nrank/numCompatibleGpus;
130         }
131         else
132         {
133             if (rank == 0)
134             {
135                 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.",
136                           nrank, numCompatibleGpus);
137             }
138
139 #if GMX_MPI
140             /* We use a global barrier to prevent ranks from continuing with
141              * an invalid setup.
142              */
143             MPI_Barrier(MPI_COMM_WORLD);
144 #endif
145         }
146     }
147
148     /* Here we will waste GPUs when nrank < numCompatibleGpus */
149     std::vector<int> taskAssignment;
150     taskAssignment.resize(std::min(numCompatibleGpus*nshare, nrank));
151     for (size_t i = 0; i != taskAssignment.size(); ++i)
152     {
153         /* TODO: improve this implementation: either sort GPUs or remove the weakest here */
154         taskAssignment[i] = compatibleGpus[i/nshare];
155     }
156     return taskAssignment;
157 }
158
159 /*! \brief Check that all user-selected GPUs are compatible.
160  *
161  * Given the \c userGpuTaskAssignment and \c compatibleGPUs, give a fatal
162  * error if any selected GPUs is not compatible
163  *
164  * The error is given with a suitable descriptive message, which will
165  * have context if this check is done after the hardware detection
166  * results have been reported to the user. However, note that only the
167  * GPUs detected on the master rank are reported, because of the
168  * existing limitations of that reporting.
169  *
170  * \todo Note that the selected GPUs can be different on each rank,
171  * and the IDs of compatible GPUs can be different on each node, so
172  * this routine ought to do communication to determine whether all
173  * ranks are able to proceed. Currently this relies on the MPI runtime
174  * to kill the other processes because GROMACS lacks the appropriate
175  * infrastructure to do a good job of coordinating error messages and
176  * behaviour across MPMD ranks and multiple simulations.
177  *
178  * \param[in]   gpu_info               GPU information including device description.
179  * \param[in]   compatibleGpus         Vector of compatible GPUs
180  * \param[in]   userGpuTaskAssignment  The GPU selection from the user.
181  */
182 static void exitUnlessUserGpuTaskAssignmentIsValid(const gmx_gpu_info_t   &gpu_info,
183                                                    const std::vector<int> &compatibleGpus,
184                                                    const std::vector<int> &userGpuTaskAssignment)
185 {
186     int         numIncompatibleGpuIds = 0;
187     std::string message
188         = "Some of the requested GPUs do not exist, behave strangely, or are not compatible:\n";
189
190     for (const auto &gpuId : userGpuTaskAssignment)
191     {
192         if (std::find(compatibleGpus.begin(), compatibleGpus.end(), gpuId) == compatibleGpus.end())
193         {
194             numIncompatibleGpuIds++;
195             message += gmx::formatString("    GPU #%d: %s\n",
196                                          gpuId,
197                                          getGpuCompatibilityDescription(gpu_info, gpuId));
198         }
199     }
200
201     if (numIncompatibleGpuIds > 0)
202     {
203         gmx_fatal(FARGS, message.c_str());
204     }
205 }
206
207 std::vector<int> getCompatibleGpus(const gmx_gpu_info_t &gpu_info)
208 {
209     // Possible minor over-allocation here, but not important for anything
210     std::vector<int> compatibleGpus;
211     compatibleGpus.reserve(gpu_info.n_dev);
212     for (int i = 0; i < gpu_info.n_dev; i++)
213     {
214         GMX_ASSERT(gpu_info.gpu_dev, "Invalid gpu_info.gpu_dev");
215         if (isGpuCompatible(gpu_info, i))
216         {
217             compatibleGpus.push_back(i);
218         }
219     }
220     return compatibleGpus;
221 }
222
223 std::vector<int> mapPpRanksToGpus(bool                    rankCanUseGpu,
224                                   const t_commrec        *cr,
225                                   const gmx_gpu_info_t   &gpu_info,
226                                   const gmx_hw_opt_t     &hw_opt)
227 {
228     std::vector<int> taskAssignment;
229
230     if (!rankCanUseGpu)
231     {
232         return taskAssignment;
233     }
234
235     auto compatibleGpus = getCompatibleGpus(gpu_info);
236     if (!hw_opt.gpuIdTaskAssignment.empty())
237     {
238         auto userGpuTaskAssignment = parseGpuTaskAssignment(hw_opt.gpuIdTaskAssignment);
239         exitUnlessUserGpuTaskAssignmentIsValid(gpu_info, compatibleGpus, userGpuTaskAssignment);
240         taskAssignment = userGpuTaskAssignment;
241     }
242     else
243     {
244         taskAssignment = assign_rank_gpu_ids(compatibleGpus, cr->nrank_pp_intranode, cr->rank_pp_intranode);
245     }
246     return taskAssignment;
247 }
248
249 } // namespace
250
251 /*! \brief Return the number of PP rank pairs that share a GPU device between them.
252  *
253  * Sharing GPUs among multiple PP ranks is possible via either user or
254  * automated selection. */
255 static int gmx_count_gpu_dev_shared(const std::vector<int> &gpuTaskAssignment,
256                                     bool                    userSetGpuIds)
257 {
258     int      same_count    = 0;
259
260     if (userSetGpuIds)
261     {
262         GMX_RELEASE_ASSERT(!gpuTaskAssignment.empty(),
263                            "The user cannot choose an empty set of GPU IDs, code is wrong somewhere");
264         size_t ngpu = gpuTaskAssignment.size();
265
266         for (size_t i = 0; i < ngpu - 1; i++)
267         {
268             for (size_t j = i + 1; j < ngpu; j++)
269             {
270                 same_count      += (gpuTaskAssignment[i] ==
271                                     gpuTaskAssignment[j]);
272             }
273         }
274     }
275
276     return same_count;
277 }
278
279 /* Count and return the number of unique GPUs (per node) selected.
280  *
281  * As sharing GPUs among multiple PP ranks is possible, the number of
282  * GPUs used (per node) can be different from the number of GPU IDs
283  * used.
284  */
285 static size_t gmx_count_gpu_dev_unique(const std::vector<int> &gpuTaskAssignment)
286 {
287     std::set<int> uniqIds;
288     for (const auto &deviceId : gpuTaskAssignment)
289     {
290         uniqIds.insert(deviceId);
291     }
292     return uniqIds.size();
293 }
294
295 void reportGpuUsage(const gmx::MDLogger    &mdlog,
296                     const gmx_gpu_info_t   &gpu_info,
297                     bool                    userSetGpuIds,
298                     const std::vector<int> &gpuTaskAssignment,
299                     size_t                  numPpRanks,
300                     bool                    bPrintHostName)
301 {
302     if (gpuTaskAssignment.empty())
303     {
304         return;
305     }
306
307     std::string output;
308     {
309         std::string gpuIdsString =
310             formatAndJoin(gpuTaskAssignment, ",", gmx::StringFormatter("%d"));
311         size_t      numGpusInUse = gmx_count_gpu_dev_unique(gpuTaskAssignment);
312         bool        bPluralGpus  = numGpusInUse > 1;
313
314         if (bPrintHostName)
315         {
316             char host[STRLEN];
317             gmx_gethostname(host, STRLEN);
318             output += gmx::formatString("On host %s ", host);
319         }
320         output += gmx::formatString("%zu GPU%s %sselected for this run.\n"
321                                     "Mapping of GPU ID%s to the %d PP rank%s in this node: %s\n",
322                                     numGpusInUse, bPluralGpus ? "s" : "",
323                                     userSetGpuIds ? "user-" : "auto-",
324                                     bPluralGpus ? "s" : "",
325                                     numPpRanks,
326                                     (numPpRanks > 1) ? "s" : "",
327                                     gpuIdsString.c_str());
328     }
329
330     int same_count = gmx_count_gpu_dev_shared(gpuTaskAssignment, userSetGpuIds);
331
332     if (same_count > 0)
333     {
334         output += gmx::formatString("NOTE: You assigned %s to multiple ranks.\n",
335                                     same_count > 1 ? "GPU IDs" : "a GPU ID");
336     }
337
338     if (static_cast<size_t>(gpu_info.n_dev_compatible) > numPpRanks)
339     {
340         /* TODO In principle, this warning could be warranted only on
341          * ranks on some nodes, but we lack the infrastructure to do a
342          * good job of reporting that. */
343         output += gmx::formatString("NOTE: potentially sub-optimal launch configuration using fewer\n"
344                                     "      PP ranks on a node than GPUs available on that node.\n");
345     }
346
347     /* NOTE: this print is only for and on one physical node */
348     GMX_LOG(mdlog.warning).appendText(output);
349 }