5a62bfb136633bf58ba42ba7a7f340bab44b715c
[alexxy/gromacs.git] / src / gromacs / taskassignment / taskassignment.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2017,2018,2019, 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 /*! \internal \file
36  * \brief
37  * Defines helper and factory functionality for task assignment.
38  *
39  * Note that the GPU ID assignment could potentially match many
40  * different kinds of simulation setups, including ranks from multiple
41  * simulations, ranks from the same simulation, and/or ranks with duty
42  * only for particular tasks (e.g. PME-only ranks). Which GPU ID
43  * assignments are valid will naturally depend on the other run-time
44  * options given to mdrun, and the current capabilities of the
45  * implementation.
46  *
47  * \author Mark Abraham <mark.j.abraham@gmail.com>
48  * \ingroup module_taskassignment
49  */
50 #include "gmxpre.h"
51
52 #include "taskassignment.h"
53
54 #include "config.h"
55
56 #include <string>
57 #include <vector>
58
59 #include "gromacs/hardware/hw_info.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/taskassignment/usergpuids.h"
62 #include "gromacs/utility/cstringutil.h"
63 #include "gromacs/utility/exceptions.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/gmxassert.h"
66 #include "gromacs/utility/gmxmpi.h"
67 #include "gromacs/utility/logger.h"
68 #include "gromacs/utility/physicalnodecommunicator.h"
69 #include "gromacs/utility/stringutil.h"
70 #include "gromacs/utility/sysinfo.h"
71
72 #include "findallgputasks.h"
73 #include "reportgpuusage.h"
74
75 namespace gmx
76 {
77
78 namespace
79 {
80
81 /*! \brief Build data structure of types of GPU tasks on a rank,
82  * together with the mapped GPU device IDs, for all GPU tasks on all
83  * the ranks of this node.
84  *
85  * \param[in]   gpuTasksOnRanksOfThisNode  For each rank on this node, the set of tasks
86  *                                         that are eligible to run on GPUs.
87  * \param[in]   gpuIds                     The user-supplied GPU IDs.
88  */
89 GpuTaskAssignments
90 buildTaskAssignment(const GpuTasksOnRanks  &gpuTasksOnRanksOfThisNode,
91                     ArrayRef<const int>     gpuIds)
92 {
93     GpuTaskAssignments gpuTaskAssignmentOnRanksOfThisNode(gpuTasksOnRanksOfThisNode.size());
94
95     // Loop over the ranks on this node, and the tasks on each
96     // rank. For each task, take the next device ID from those
97     // provided by the user, to build a vector of mappings of task to
98     // ID, for each rank on this node. Note that if there have not
99     // been any GPU tasks identified, then gpuIds can be empty.
100     auto               currentGpuId            = gpuIds.begin();
101     auto               gpuTaskAssignmentOnRank = gpuTaskAssignmentOnRanksOfThisNode.begin();
102     for (const auto &gpuTasksOnRank : gpuTasksOnRanksOfThisNode)
103     {
104         gpuTaskAssignmentOnRank->reserve(gpuTasksOnRank.size());
105         for (const auto &gpuTaskType : gpuTasksOnRank)
106         {
107             GMX_RELEASE_ASSERT(currentGpuId != gpuIds.end(), "Indexing out of range for GPU tasks");
108             gpuTaskAssignmentOnRank->push_back({gpuTaskType, *currentGpuId});
109             ++currentGpuId;
110         }
111         GMX_RELEASE_ASSERT(gpuTaskAssignmentOnRank->size() == gpuTasksOnRank.size(),
112                            "Mismatch in number of GPU tasks on a rank with the number of elements in the resulting task assignment");
113         ++gpuTaskAssignmentOnRank;
114     }
115
116     return gpuTaskAssignmentOnRanksOfThisNode;
117 }
118
119 /*! \brief Return whether a GPU device is shared between any ranks.
120  *
121  * Sharing GPUs among multiple ranks is possible via either user or
122  * automated selection. */
123 bool isAnyGpuSharedBetweenRanks(const GpuTaskAssignments &gpuTaskAssignments)
124 {
125     // Loop over all ranks i, looking on all higher ranks j whether
126     // any tasks on them share GPU device IDs.
127     //
128     // TODO Should this functionality also consider whether tasks on
129     // the same rank are sharing a device?
130     for (size_t i = 0; i < gpuTaskAssignments.size(); ++i)
131     {
132         for (const auto &taskOnRankI : gpuTaskAssignments[i])
133         {
134             for (size_t j = i+1; j < gpuTaskAssignments.size(); ++j)
135             {
136                 for (const auto &taskOnRankJ : gpuTaskAssignments[j])
137                 {
138                     if (taskOnRankI.deviceId_ == taskOnRankJ.deviceId_)
139                     {
140                         return true;
141                     }
142                 }
143             }
144         }
145     }
146     return false;
147 }
148
149 //! Logs to \c mdlog information that may help a user learn how to let mdrun make a task assignment that runs faster.
150 void logPerformanceHints(const MDLogger           &mdlog,
151                          size_t                    numCompatibleGpus,
152                          size_t                    numGpuTasksOnThisNode,
153                          const GpuTaskAssignments &gpuTaskAssignments)
154 {
155     if (numCompatibleGpus > numGpuTasksOnThisNode)
156     {
157         /* TODO In principle, this warning could be warranted only on
158          * some nodes, but we lack the infrastructure to do a good job
159          * of reporting that. */
160         GMX_LOG(mdlog.warning).asParagraph().
161             appendText("NOTE: You assigned the GPU tasks on a node such that some GPUs "
162                        "available on that node are unused, which might not be optimal.");
163     }
164
165     if (isAnyGpuSharedBetweenRanks(gpuTaskAssignments))
166     {
167         GMX_LOG(mdlog.warning).asParagraph().
168             appendText("NOTE: You assigned the same GPU ID(s) to multiple ranks, which is a good idea if you have measured the performance of alternatives.");
169     }
170 }
171
172 //! Counts all the GPU tasks on this node.
173 size_t countGpuTasksOnThisNode(const GpuTasksOnRanks &gpuTasksOnRanksOfThisNode)
174 {
175     size_t numGpuTasksOnThisNode = 0;
176     for (const auto &gpuTasksOnRank : gpuTasksOnRanksOfThisNode)
177     {
178         numGpuTasksOnThisNode += gpuTasksOnRank.size();
179     }
180     return numGpuTasksOnThisNode;
181 }
182
183 }   // namespace
184
185 GpuTaskAssignments::value_type
186 runTaskAssignment(const std::vector<int>         &gpuIdsToUse,
187                   const std::vector<int>         &userGpuTaskAssignment,
188                   const gmx_hw_info_t            &hardwareInfo,
189                   const MDLogger                 &mdlog,
190                   const t_commrec                *cr,
191                   const gmx_multisim_t           *ms,
192                   const PhysicalNodeCommunicator &physicalNodeComm,
193                   const std::vector<GpuTask>     &gpuTasksOnThisRank,
194                   bool                            useGpuForBonded,
195                   PmeRunMode                      pmeRunMode)
196 {
197     /* Communicate among ranks on this node to find each task that can
198      * be executed on a GPU, on each rank. */
199     auto               gpuTasksOnRanksOfThisNode = findAllGpuTasksOnThisNode(gpuTasksOnThisRank,
200                                                                              physicalNodeComm);
201     auto               numGpuTasksOnThisNode = countGpuTasksOnThisNode(gpuTasksOnRanksOfThisNode);
202
203     GpuTaskAssignments taskAssignmentOnRanksOfThisNode;
204     try
205     {
206         // Use the GPU IDs from the user if they supplied
207         // them. Otherwise, choose from the compatible GPUs.
208         //
209         // GPU ID assignment strings, if provided, cover all the ranks
210         // on a node. If nodes or the process placement on them are
211         // heterogeneous, then the GMX_GPU_ID environment variable
212         // must be set by a user who also wishes to direct GPU ID
213         // assignment.  Thus this implementation of task assignment
214         // can assume it has a GPU ID assignment appropriate for the
215         // node upon which its process is running.
216         //
217         // Valid GPU ID assignments are `an ordered set of digits that
218         // identify GPU device IDs (e.g. as understood by the GPU
219         // runtime, and subject to environment modification such as
220         // with CUDA_VISIBLE_DEVICES) that will be used for the
221         // GPU-suitable tasks on all of the ranks of that node.
222         ArrayRef<const int> gpuIdsForTaskAssignment;
223         std::vector<int>    generatedGpuIds;
224         if (userGpuTaskAssignment.empty())
225         {
226             ArrayRef<const int> compatibleGpusToUse = gpuIdsToUse;
227
228             // enforce the single device/rank restriction
229             if (physicalNodeComm.size_ == 1 && !compatibleGpusToUse.empty())
230             {
231                 compatibleGpusToUse = compatibleGpusToUse.subArray(0, 1);
232             }
233
234             // When doing automated assignment of GPU tasks to GPU
235             // IDs, even if we have more than one kind of GPU task, we
236             // do a simple round-robin assignment. That's not ideal,
237             // but we don't have any way to do a better job reliably.
238             generatedGpuIds         = makeGpuIds(compatibleGpusToUse, numGpuTasksOnThisNode);
239
240             if ((numGpuTasksOnThisNode > gpuIdsToUse.size()) &&
241                 (numGpuTasksOnThisNode % gpuIdsToUse.size() != 0))
242             {
243                 // TODO Decorating the message with hostname should be
244                 // the job of an error-reporting module.
245                 char host[STRLEN];
246                 gmx_gethostname(host, STRLEN);
247
248                 GMX_THROW(InconsistentInputError
249                               (formatString("There were %zu GPU tasks found on node %s, but %zu GPUs were "
250                                             "available. If the GPUs are equivalent, then it is usually best "
251                                             "to have a number of tasks that is a multiple of the number of GPUs. "
252                                             "You should reconsider your GPU task assignment, "
253                                             "number of ranks, or your use of the -nb, -pme, and -npme options, "
254                                             "perhaps after measuring the performance you can get.", numGpuTasksOnThisNode,
255                                             host, gpuIdsToUse.size())));
256             }
257             gpuIdsForTaskAssignment = generatedGpuIds;
258         }
259         else
260         {
261             if (numGpuTasksOnThisNode != userGpuTaskAssignment.size())
262             {
263                 // TODO Decorating the message with hostname should be
264                 // the job of an error-reporting module.
265                 char host[STRLEN];
266                 gmx_gethostname(host, STRLEN);
267
268                 GMX_THROW(InconsistentInputError
269                               (formatString("There were %zu GPU tasks assigned on node %s, but %zu GPU tasks were "
270                                             "identified, and these must match. Reconsider your GPU task assignment, "
271                                             "number of ranks, or your use of the -nb, -pme, and -npme options.", userGpuTaskAssignment.size(),
272                                             host, numGpuTasksOnThisNode)));
273             }
274             // Did the user choose compatible GPUs?
275             checkUserGpuIds(hardwareInfo.gpu_info, gpuIdsToUse, userGpuTaskAssignment);
276
277             gpuIdsForTaskAssignment = userGpuTaskAssignment;
278         }
279         taskAssignmentOnRanksOfThisNode =
280             buildTaskAssignment(gpuTasksOnRanksOfThisNode, gpuIdsForTaskAssignment);
281
282     }
283     catch (const std::exception &ex)
284     {
285         // TODO This implementation is quite similar to that of
286         // processExceptionAsFatalError (which implements
287         // GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR), but it is unclear
288         // how we should involve MPI in the implementation of error
289         // handling.
290         if (physicalNodeComm.rank_ == 0)
291         {
292             printFatalErrorMessage(stderr, ex);
293         }
294
295         if (PAR(cr))
296         {
297 #if GMX_MPI
298             MPI_Barrier(cr->mpi_comm_mysim);
299 #endif
300         }
301         if (isMultiSim(ms))
302         {
303 #if GMX_MPI
304             MPI_Barrier(ms->mpi_comm_masters);
305 #endif
306         }
307
308         gmx_exit_on_fatal_error(ExitType_Abort, 1);
309     }
310
311     reportGpuUsage(mdlog, taskAssignmentOnRanksOfThisNode,
312                    numGpuTasksOnThisNode, physicalNodeComm.size_, cr->nnodes > 1,
313                    useGpuForBonded, pmeRunMode);
314
315     // If the user chose a task assignment, give them some hints where appropriate.
316     if (!userGpuTaskAssignment.empty())
317     {
318         logPerformanceHints(mdlog, gpuIdsToUse.size(),
319                             numGpuTasksOnThisNode,
320                             taskAssignmentOnRanksOfThisNode);
321     }
322
323     return taskAssignmentOnRanksOfThisNode[physicalNodeComm.rank_];
324
325     // TODO There is no check that mdrun -nb gpu or -pme gpu or
326     // -gpu_id is actually being implemented such that nonbonded tasks
327     // are being run on compatible GPUs, on all applicable ranks. That
328     // would require communication.
329 }
330
331 }  // namespace gmx