2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2017,2018,2019,2020, 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.
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
50 #ifndef GMX_TASKASSIGNMENT_TASKASSIGNMENT_H
51 #define GMX_TASKASSIGNMENT_TASKASSIGNMENT_H
55 #include "gromacs/utility/basedefinitions.h"
56 #include "gromacs/utility/gmxmpi.h"
58 struct DeviceInformation
;
62 enum class PmeRunMode
;
67 enum class TaskTarget
;
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.
79 //! Long-ranged interactions.
81 //! Number of possible tasks.
86 * \brief Specifies the GPU deviceID_ available for task_ to use. */
89 //! The type of this GPU task.
91 //! Device ID on this node to which this GPU task is mapped.
95 //! 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.
96 using GpuTaskAssignment
= std::vector
<GpuTaskMapping
>;
98 class GpuTaskAssignments
;
101 * \brief Builder for the GpuTaskAssignments for all ranks on this
104 * This will coordinate the final stages of task assignment and
105 * reporting, and build the GpuTaskAssignments object used to
106 * configure the modules that might run tasks on GPUs.
108 * Communicates between ranks on a node to coordinate task assignment
109 * between them onto available hardware, e.g. accelerators.
111 * \todo Later, this might become a loop over all registered modules
112 * relevant to the mdp inputs, to find those that have such tasks.
114 * \todo Later we might need the concept of computeTasksOnThisRank,
115 * from which we construct gpuTasksOnThisRank.
117 * Currently the DD code assigns duty to ranks that can
118 * include PP work that currently can be executed on a single
119 * GPU, if present and compatible. This has to be coordinated
120 * across PP ranks on a node, with possible multiple devices
121 * or sharing devices on a node, either from the user
122 * selection, or automatically. */
123 class GpuTaskAssignmentsBuilder
127 GpuTaskAssignmentsBuilder();
129 /*! \brief Builds a GpuTaskAssignments
131 * This method reconciles
133 * - user mdrun command-line options,
134 * - the results of hardware detection
135 * - the duty assigned by the DD setup,
136 * - the requested simulation modules, and
137 * - the possible existence of multi-simulations
139 * to assign the GPUs on each physical node to the tasks on
140 * the ranks of that node. It throws InconsistentInputError
141 * when a/the useful GPU task assignment is not possible.
143 * \param[in] gpuIdsToUse The compatible GPUs that the user permitted us to use.
144 * \param[in] userGpuTaskAssignment The user-specified assignment of GPU tasks to device IDs.
145 * \param[in] hardwareInfo The detected hardware
146 * \param[in] gromacsWorldComm MPI communicator for all ranks in the current GROMACS run
147 * \param[in] physicalNodeComm Communication object for this physical node.
148 * \param[in] nonbondedTarget The user's choice for mdrun -nb for where to assign
149 * short-ranged nonbonded interaction tasks.
150 * \param[in] pmeTarget The user's choice for mdrun -pme for where to assign
151 * long-ranged PME nonbonded interaction tasks.
152 * \param[in] bondedTarget The user's choice for mdrun -bonded for where to assign tasks.
153 * \param[in] updateTarget The user's choice for mdrun -update for where to assign tasks.
154 * \param[in] useGpuForNonbonded Whether GPUs will be used for nonbonded interactions.
155 * \param[in] useGpuForPme Whether GPUs will be used for PME interactions.
156 * \param[in] rankHasPpTask Whether this rank has a PP task
157 * \param[in] rankHasPmeTask Whether this rank has a PME task
159 * \throws std::bad_alloc If out of memory.
160 * InconsistentInputError If user and/or detected inputs are inconsistent.
162 static GpuTaskAssignments
build(const std::vector
<int>& gpuIdsToUse
,
163 const std::vector
<int>& userGpuTaskAssignment
,
164 const gmx_hw_info_t
& hardwareInfo
,
165 MPI_Comm gromacsWorldComm
,
166 const PhysicalNodeCommunicator
& physicalNodeComm
,
167 TaskTarget nonbondedTarget
,
168 TaskTarget pmeTarget
,
169 TaskTarget bondedTarget
,
170 TaskTarget updateTarget
,
171 bool useGpuForNonbonded
,
174 bool rankHasPmeTask
);
178 * \brief Contains the GPU task assignment for all ranks on this
181 * This can be used to configure the modules that might run tasks on
184 * This assignment is made by a GpuTaskAssignmentsBuilder object. */
185 class GpuTaskAssignments
188 //! Public move constructor to use with the builder
189 GpuTaskAssignments(GpuTaskAssignments
&& source
) noexcept
= default;
192 // Let the builder handle construction
193 friend class GpuTaskAssignmentsBuilder
;
194 //! Private constructor so only the builder can construct
195 GpuTaskAssignments(const gmx_hw_info_t
& hardwareInfo
);
196 /*! \brief Information about hardware on this physical node
198 * The lifetime of the object referred to must exceed that
200 const gmx_hw_info_t
& hardwareInfo_
;
201 //! The GPU task assignment for all ranks on this node
202 std::vector
<GpuTaskAssignment
> assignmentForAllRanksOnThisNode_
;
203 /*! \brief The index of this rank within those on this node.
205 * This is useful for indexing into \c
206 * assignmentForAllRanksOnThisNode_. */
207 index indexOfThisRank_
= -1;
208 //! Number of GPU tasks on this node.
209 size_t numGpuTasksOnThisNode_
= 0;
210 //! Number of ranks on this physical node.
211 size_t numRanksOnThisNode_
= 0;
214 /*! \brief Log a report on how GPUs are being used on
215 * the ranks of the physical node of rank 0 of the simulation.
217 * \todo It could be useful to report also whether any nodes differed,
220 * \param[in] mdlog Logging object.
221 * \param[in] printHostName Print the hostname in the usage information.
222 * \param[in] useGpuForBonded Whether GPU PP tasks will do bonded work on the GPU.
223 * \param[in] pmeRunMode Describes the execution of PME tasks.
224 * \param[in] useGpuForUpdate Whether the update is offloaded on the GPU.
226 * \throws std::bad_alloc if out of memory
228 void reportGpuUsage(const MDLogger
& mdlog
,
230 bool useGpuForBonded
,
231 PmeRunMode pmeRunMode
,
232 bool useGpuForUpdate
);
234 /*! \brief Logs to \c mdlog information that may help a user
235 * learn how to let mdrun make a task assignment that runs
238 * \param[in] mdlog Logging object.
239 * \param[in] numCompatibleGpusOnThisNode The number of compatible GPUs on this node.
241 void logPerformanceHints(const MDLogger
& mdlog
, size_t numCompatibleGpusOnThisNode
);
242 /*! \brief Return handle to the initialized GPU to use in the this rank.
244 * \param[out] deviceId Index of the assigned device.
246 * \returns Device information on the selected devicce. Returns nullptr if no GPU task
247 * is assigned to this rank.
249 DeviceInformation
* initDevice(int* deviceId
) const;
250 //! Return whether this rank has a PME task running on a GPU
251 bool thisRankHasPmeGpuTask() const;
252 //! Return whether this rank has any task running on a GPU
253 bool thisRankHasAnyGpuTask() const;