Move GPU task assignment
[gromacs.git] / src / gromacs / taskassignment / taskassignment.h
blob689035806e8d454d6ff9f7e9201268d1252ca236
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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 /*! \defgroup module_taskassignment Assigning simulation tasks to hardware (taskassignment)
36 * \ingroup group_mdrun
37 * \brief Provides code that manages assignment of simulation tasks to hardware.
39 /*! \libinternal
40 * \file
41 * \brief Declares high-level functionality for managing assigning
42 * tasks on ranks of a node to hardware on that node, and the factory
43 * function to build the correct flavours of gmx::INodeTaskAssigner
44 * required to implement the user's requirements.
46 * \author Mark Abraham <mark.j.abraham@gmail.com>
47 * \ingroup module_taskassignment
48 * \inlibraryapi
50 #ifndef GMX_TASKASSIGNMENT_TASKASSIGNMENT_H
51 #define GMX_TASKASSIGNMENT_TASKASSIGNMENT_H
53 #include <vector>
55 #include "gromacs/utility/basedefinitions.h"
57 struct gmx_device_info_t;
58 struct gmx_hw_info_t;
59 struct gmx_multisim_t;
60 struct t_commrec;
62 enum class PmeRunMode;
64 namespace gmx
67 enum class TaskTarget;
68 class MDLogger;
69 class PhysicalNodeCommunicator;
71 /*! \brief Types of compute tasks that can be run on a GPU.
73 * These names refer to existing practice in GROMACS, which is not
74 * strictly accurate. */
75 enum class GpuTask : int
77 //! Short-ranged interactions.
78 Nonbonded,
79 //! Long-ranged interactions.
80 Pme
83 /*! \libinternal
84 * \brief Specifies the GPU deviceID_ available for task_ to use. */
85 struct GpuTaskMapping
87 //! The type of this GPU task.
88 GpuTask task_;
89 //! Device ID on this node to which this GPU task is mapped.
90 int deviceId_;
93 //! Container of GPU tasks on a rank, specifying the task type and GPU device ID, e.g. potentially ready for consumption by the modules on that rank.
94 using GpuTaskAssignment = std::vector <GpuTaskMapping>;
96 class GpuTaskAssignments;
98 /*! \libinternal
99 * \brief Builder for the GpuTaskAssignments for all ranks on this
100 * node.
102 * This will coordinate the final stages of task assignment and
103 * reporting, and build the GpuTaskAssignments object used to
104 * configure the modules that might run tasks on GPUs.
106 * Communicates between ranks on a node to coordinate task assignment
107 * between them onto available hardware, e.g. accelerators.
109 * \todo Later, this might become a loop over all registered modules
110 * relevant to the mdp inputs, to find those that have such tasks.
112 * \todo Later we might need the concept of computeTasksOnThisRank,
113 * from which we construct gpuTasksOnThisRank.
115 * Currently the DD code assigns duty to ranks that can
116 * include PP work that currently can be executed on a single
117 * GPU, if present and compatible. This has to be coordinated
118 * across PP ranks on a node, with possible multiple devices
119 * or sharing devices on a node, either from the user
120 * selection, or automatically. */
121 class GpuTaskAssignmentsBuilder
123 public:
124 //! Constructor
125 GpuTaskAssignmentsBuilder();
127 /*! \brief Builds a GpuTaskAssignments
129 * This method reconciles
131 * - user mdrun command-line options,
132 * - the results of hardware detection
133 * - the duty assigned by the DD setup,
134 * - the requested simulation modules, and
135 * - the possible existence of multi-simulations
137 * to assign the GPUs on each physical node to the tasks on
138 * the ranks of that node. It throws InconsistentInputError
139 * when a/the useful GPU task assignment is not possible.
141 * \param[in] gpuIdsToUse The compatible GPUs that the user permitted us to use.
142 * \param[in] userGpuTaskAssignment The user-specified assignment of GPU tasks to device IDs.
143 * \param[in] hardwareInfo The detected hardware
144 * \param[in] cr Communication object.
145 * \param[in] ms Multi-simulation handler.
146 * \param[in] physicalNodeComm Communication object for this physical node.
147 * \param[in] nonbondedTarget The user's choice for mdrun -nb for where to assign
148 * short-ranged nonbonded interaction tasks.
149 * \param[in] pmeTarget The user's choice for mdrun -pme for where to assign
150 * long-ranged PME nonbonded interaction tasks.
151 * \param[in] bondedTarget The user's choice for mdrun -bonded for where to assign tasks.
152 * \param[in] updateTarget The user's choice for mdrun -update for where to assign tasks.
153 * \param[in] useGpuForNonbonded Whether GPUs will be used for nonbonded interactions.
154 * \param[in] useGpuForPme Whether GPUs will be used for PME interactions.
155 * \param[in] rankHasPpTask Whether this rank has a PP task
156 * \param[in] rankHasPmeTask Whether this rank has a PME task
158 * \throws std::bad_alloc If out of memory.
159 * InconsistentInputError If user and/or detected inputs are inconsistent.
161 GpuTaskAssignments build(const std::vector<int> &gpuIdsToUse,
162 const std::vector<int> &userGpuTaskAssignment,
163 const gmx_hw_info_t &hardwareInfo,
164 const t_commrec *cr,
165 const gmx_multisim_t *ms,
166 const PhysicalNodeCommunicator &physicalNodeComm,
167 TaskTarget nonbondedTarget,
168 TaskTarget pmeTarget,
169 TaskTarget bondedTarget,
170 TaskTarget updateTarget,
171 bool useGpuForNonbonded,
172 bool useGpuForPme,
173 bool rankHasPpTask,
174 bool rankHasPmeTask);
177 /*! \libinternal
178 * \brief Contains the GPU task assignment for all ranks on this
179 * physical node.
181 * This can be used to configure the modules that might run tasks on
182 * GPUs.
184 * This assignment is made by a GpuTaskAssignmentsBuilder object. */
185 class GpuTaskAssignments
187 public:
188 //! Public move constructor to use with the builder
189 GpuTaskAssignments(GpuTaskAssignments &&source) noexcept = default;
190 private:
191 // Let the builder handle construction
192 friend class GpuTaskAssignmentsBuilder;
193 //! Private constructor so only the builder can construct
194 GpuTaskAssignments(const gmx_hw_info_t &hardwareInfo);
195 /*! \brief Information about hardware on this physical node
197 * The lifetime of the object referred to must exceed that
198 * of this object. */
199 const gmx_hw_info_t &hardwareInfo_;
200 //! The GPU task assignment for all ranks on this node
201 std::vector<GpuTaskAssignment> assignmentForAllRanksOnThisNode_;
202 /*! \brief The index of this rank within those on this node.
204 * This is useful for indexing into \c
205 * assignmentForAllRanksOnThisNode_. */
206 index indexOfThisRank_ = -1;
207 //! Number of GPU tasks on this node.
208 size_t numGpuTasksOnThisNode_ = 0;
209 //! Number of ranks on this physical node.
210 size_t numRanksOnThisNode_ = 0;
211 public:
212 /*! \brief Log a report on how GPUs are being used on
213 * the ranks of the physical node of rank 0 of the simulation.
215 * \todo It could be useful to report also whether any nodes differed,
216 * and in what way.
218 * \param[in] mdlog Logging object.
219 * \param[in] printHostName Print the hostname in the usage information
220 * \param[in] useGpuForBonded Whether GPU PP tasks will do bonded work on the GPU
221 * \param[in] pmeRunMode Describes the execution of PME tasks
223 * \throws std::bad_alloc if out of memory */
224 void
225 reportGpuUsage(const MDLogger &mdlog,
226 bool printHostName,
227 bool useGpuForBonded,
228 PmeRunMode pmeRunMode);
229 /*! \brief Logs to \c mdlog information that may help a user
230 * learn how to let mdrun make a task assignment that runs
231 * faster.
233 * \param[in] mdlog Logging object.
234 * \param[in] numCompatibleGpusOnThisNode The number of compatible GPUs on this node.
235 * */
236 void logPerformanceHints(const MDLogger &mdlog,
237 size_t numCompatibleGpusOnThisNode);
238 /*! \brief Return handle to the initialized GPU to use for the
239 * nonbonded task on this rank, if any.
241 * Returns nullptr if no such task is assigned to this rank.
243 * \todo This also sets up DLB for device sharing, where
244 * appropriate, but that responsbility should move
245 * elsewhere. */
246 gmx_device_info_t *initNonbondedDevice(const t_commrec *cr) const;
247 /*! \brief Return handle to the initialized GPU to use for the
248 * PME task on this rank, if any.
250 * Returns nullptr if no such task is assigned to this rank. */
251 gmx_device_info_t *initPmeDevice() const;
252 //! Return whether this rank has a PME task running on a GPU
253 bool thisRankHasPmeGpuTask() const;
254 //! Return whether this rank has any task running on a GPU
255 bool thisRankHasAnyGpuTask() const;
258 //! Function for whether the task of \c mapping has value \c TaskType.
259 template<GpuTask TaskType>
260 bool hasTaskType(const GpuTaskMapping &mapping)
262 return mapping.task_ == TaskType;
265 } // namespace gmx
267 #endif