Extend task assignment code
[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, 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/stringutil.h"
69 #include "gromacs/utility/sysinfo.h"
70
71 #include "findallgputasks.h"
72 #include "reportgpuusage.h"
73
74 namespace gmx
75 {
76
77 namespace
78 {
79
80 /*! \brief Build data structure of types of GPU tasks on a rank,
81  * together with the mapped GPU device IDs, for all GPU tasks on all
82  * the ranks of this node.
83  *
84  * \param[in]   gpuTasksOnRanksOfThisNode  For each rank on this node, the set of tasks
85  *                                         that are eligible to run on GPUs.
86  * \param[in]   gpuIds                     The user-supplied GPU IDs.
87  */
88 static GpuTaskAssignments
89 buildTaskAssignment(const GpuTasksOnRanks  &gpuTasksOnRanksOfThisNode,
90                     ArrayRef<const int>     gpuIds)
91 {
92     GpuTaskAssignments gpuTaskAssignmentOnRanksOfThisNode(gpuTasksOnRanksOfThisNode.size());
93
94     // Loop over the ranks on this node, and the tasks on each
95     // rank. For each task, take the next device ID from those
96     // provided by the user, to build a vector of mappings of task to
97     // ID, for each rank on this node. Note that if there have not
98     // been any GPU tasks identified, then gpuIds can be empty.
99     auto               currentGpuId            = gpuIds.begin();
100     auto               gpuTaskAssignmentOnRank = gpuTaskAssignmentOnRanksOfThisNode.begin();
101     for (const auto &gpuTasksOnRank : gpuTasksOnRanksOfThisNode)
102     {
103         gpuTaskAssignmentOnRank->reserve(gpuTasksOnRank.size());
104         for (const auto &gpuTaskType : gpuTasksOnRank)
105         {
106             GMX_RELEASE_ASSERT(currentGpuId != gpuIds.end(), "Indexing out of range for GPU tasks");
107             gpuTaskAssignmentOnRank->push_back({gpuTaskType, *currentGpuId});
108             ++currentGpuId;
109         }
110         GMX_RELEASE_ASSERT(gpuTaskAssignmentOnRank->size() == gpuTasksOnRank.size(),
111                            "Mismatch in number of GPU tasks on a rank with the number of elements in the resulting task assignment");
112         ++gpuTaskAssignmentOnRank;
113     }
114
115     return gpuTaskAssignmentOnRanksOfThisNode;
116 }
117
118 /*! \brief Return whether a GPU device is shared between any ranks.
119  *
120  * Sharing GPUs among multiple ranks is possible via either user or
121  * automated selection. */
122 static bool isAnyGpuSharedBetweenRanks(const GpuTaskAssignments &gpuTaskAssignments)
123 {
124     // Loop over all ranks i, looking on all higher ranks j whether
125     // any tasks on them share GPU device IDs.
126     //
127     // TODO Should this functionality also consider whether tasks on
128     // the same rank are sharing a device?
129     for (size_t i = 0; i < gpuTaskAssignments.size(); ++i)
130     {
131         for (const auto &taskOnRankI : gpuTaskAssignments[i])
132         {
133             for (size_t j = i+1; j < gpuTaskAssignments.size(); ++j)
134             {
135                 for (const auto &taskOnRankJ : gpuTaskAssignments[j])
136                 {
137                     if (taskOnRankI.deviceId_ == taskOnRankJ.deviceId_)
138                     {
139                         return true;
140                     }
141                 }
142             }
143         }
144     }
145     return false;
146 }
147
148 //! Logs to \c mdlog information that may help a user learn how to let mdrun make a task assignment that runs faster.
149 void logPerformanceHints(const MDLogger           &mdlog,
150                          size_t                    numCompatibleGpus,
151                          size_t                    numGpuTasksOnThisNode,
152                          const GpuTaskAssignments &gpuTaskAssignments)
153 {
154     if (numCompatibleGpus > numGpuTasksOnThisNode)
155     {
156         /* TODO In principle, this warning could be warranted only on
157          * some nodes, but we lack the infrastructure to do a good job
158          * of reporting that. */
159         GMX_LOG(mdlog.warning).asParagraph().
160             appendText("NOTE: You assigned the GPU tasks on a node such that some GPUs "
161                        "available on that node are unused, which might not be optimal.");
162     }
163
164     if (isAnyGpuSharedBetweenRanks(gpuTaskAssignments))
165     {
166         GMX_LOG(mdlog.warning).asParagraph().
167             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.");
168     }
169 }
170
171 //! Counts all the GPU tasks on this node.
172 size_t countGpuTasksOnThisNode(const GpuTasksOnRanks &gpuTasksOnRanksOfThisNode)
173 {
174     size_t numGpuTasksOnThisNode = 0;
175     for (const auto &gpuTasksOnRank : gpuTasksOnRanksOfThisNode)
176     {
177         numGpuTasksOnThisNode += gpuTasksOnRank.size();
178     }
179     return numGpuTasksOnThisNode;
180 }
181
182 }   // namespace
183
184 GpuTaskAssignments::value_type
185 runTaskAssignment(const std::vector<int>     &gpuIdsToUse,
186                   const std::vector<int>     &userGpuTaskAssignment,
187                   const gmx_hw_info_t        &hardwareInfo,
188                   const MDLogger             &mdlog,
189                   const t_commrec            *cr,
190                   const std::vector<GpuTask> &gpuTasksOnThisRank)
191 {
192     /* Communicate among ranks on this node to find each task that can
193      * be executed on a GPU, on each rank. */
194     auto gpuTasksOnRanksOfThisNode = findAllGpuTasksOnThisNode(gpuTasksOnThisRank,
195                                                                cr->nrank_intranode,
196                                                                cr->mpi_comm_physicalnode);
197     auto               numGpuTasksOnThisNode = countGpuTasksOnThisNode(gpuTasksOnRanksOfThisNode);
198
199     GpuTaskAssignments taskAssignmentOnRanksOfThisNode;
200     try
201     {
202         // Use the GPU IDs from the user if they supplied
203         // them. Otherwise, choose from the compatible GPUs.
204         //
205         // GPU ID assignment strings, if provided, cover all the ranks
206         // on a node. If nodes or the process placement on them are
207         // heterogeneous, then the GMX_GPU_ID environment variable
208         // must be set by a user who also wishes to direct GPU ID
209         // assignment.  Thus this implementation of task assignment
210         // can assume it has a GPU ID assignment appropriate for the
211         // node upon which its process is running.
212         //
213         // Valid GPU ID assignments are `an ordered set of digits that
214         // identify GPU device IDs (e.g. as understood by the GPU
215         // runtime, and subject to environment modification such as
216         // with CUDA_VISIBLE_DEVICES) that will be used for the
217         // GPU-suitable tasks on all of the ranks of that node.
218         ArrayRef<const int> gpuIdsForTaskAssignment;
219         std::vector<int>    generatedGpuIds;
220         if (userGpuTaskAssignment.empty())
221         {
222             generatedGpuIds         = makeGpuIds(gpuIdsToUse, numGpuTasksOnThisNode);
223             gpuIdsForTaskAssignment = generatedGpuIds;
224         }
225         else
226         {
227             if (numGpuTasksOnThisNode != userGpuTaskAssignment.size())
228             {
229                 // TODO Decorating the message with hostname should be
230                 // the job of an error-reporting module.
231                 char host[STRLEN];
232                 gmx_gethostname(host, STRLEN);
233
234                 GMX_THROW(InconsistentInputError
235                               (formatString("There were %zu GPU tasks assigned on node %s, but %zu GPU tasks were "
236                                             "identified, and these must match. Reconsider your GPU task assignment, "
237                                             "number of ranks, or your use of the -nb, -pme, and -npme options.", userGpuTaskAssignment.size(),
238                                             host, numGpuTasksOnThisNode)));
239             }
240             // Did the user choose compatible GPUs?
241             checkUserGpuIds(hardwareInfo.gpu_info, gpuIdsToUse, userGpuTaskAssignment);
242
243             gpuIdsForTaskAssignment = userGpuTaskAssignment;
244         }
245         taskAssignmentOnRanksOfThisNode =
246             buildTaskAssignment(gpuTasksOnRanksOfThisNode, gpuIdsForTaskAssignment);
247
248     }
249     catch (const std::exception &ex)
250     {
251         // TODO This implementation is quite similar to that of
252         // processExceptionAsFatalError (which implements
253         // GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR), but it is unclear
254         // how we should involve MPI in the implementation of error
255         // handling.
256         if (cr->rank_intranode == 0)
257         {
258             printFatalErrorMessage(stderr, ex);
259         }
260
261         if (PAR(cr))
262         {
263 #if GMX_MPI
264             MPI_Barrier(cr->mpi_comm_mysim);
265 #endif
266         }
267         if (MULTISIM(cr))
268         {
269 #if GMX_MPI
270             MPI_Barrier(cr->ms->mpi_comm_masters);
271 #endif
272         }
273
274         gmx_exit_on_fatal_error(ExitType_Abort, 1);
275     }
276
277     reportGpuUsage(mdlog, !userGpuTaskAssignment.empty(), taskAssignmentOnRanksOfThisNode,
278                    numGpuTasksOnThisNode, cr->nrank_intranode, cr->nnodes > 1);
279
280     // If the user chose a task assignment, give them some hints where appropriate.
281     if (!userGpuTaskAssignment.empty())
282     {
283         logPerformanceHints(mdlog, gpuIdsToUse.size(),
284                             numGpuTasksOnThisNode,
285                             taskAssignmentOnRanksOfThisNode);
286     }
287
288     return taskAssignmentOnRanksOfThisNode[cr->rank_intranode];
289
290     // TODO There is no check that mdrun -nb gpu or -pme gpu or
291     // -gpu_id is actually being implemented such that nonbonded tasks
292     // are being run on compatible GPUs, on all applicable ranks. That
293     // would require communication.
294 }
295
296 } // namespace