Split simulationWork.useGpuBufferOps into separate x and f flags
[alexxy/gromacs.git] / src / gromacs / taskassignment / taskassignment.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2017,2018,2019,2020,2021, 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 /*! \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.
38  */
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.
45  *
46  * \author Mark Abraham <mark.j.abraham@gmail.com>
47  * \ingroup module_taskassignment
48  * \inlibraryapi
49  */
50 #ifndef GMX_TASKASSIGNMENT_TASKASSIGNMENT_H
51 #define GMX_TASKASSIGNMENT_TASKASSIGNMENT_H
52
53 #include <vector>
54
55 #include "gromacs/utility/arrayref.h"
56 #include "gromacs/utility/basedefinitions.h"
57 #include "gromacs/utility/gmxmpi.h"
58
59 struct DeviceInformation;
60 struct gmx_hw_info_t;
61 struct t_commrec;
62
63 enum class PmeRunMode;
64
65 namespace gmx
66 {
67
68 enum class TaskTarget;
69 class MDLogger;
70 class PhysicalNodeCommunicator;
71
72 /*! \brief Types of compute tasks that can be run on a GPU.
73  *
74  * These names refer to existing practice in GROMACS, which is not
75  * strictly accurate. */
76 enum class GpuTask : int
77 {
78     //! Short-ranged interactions.
79     Nonbonded,
80     //! Long-ranged interactions.
81     Pme,
82     //! Number of possible tasks.
83     Count
84 };
85
86 /*! \libinternal
87  * \brief Specifies the GPU deviceID_ available for task_ to use. */
88 struct GpuTaskMapping
89 {
90     //! The type of this GPU task.
91     GpuTask task_;
92     //! Device ID on this node to which this GPU task is mapped.
93     int deviceId_;
94 };
95
96 //! 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.
97 using GpuTaskAssignment = std::vector<GpuTaskMapping>;
98
99 class GpuTaskAssignments;
100
101 /*! \libinternal
102  * \brief Builder for the GpuTaskAssignments for all ranks on this
103  * node.
104  *
105  * This will coordinate the final stages of task assignment and
106  * reporting, and build the GpuTaskAssignments object used to
107  * configure the modules that might run tasks on GPUs.
108  *
109  * Communicates between ranks on a node to coordinate task assignment
110  * between them onto available hardware, e.g. accelerators.
111  *
112  * \todo Later, this might become a loop over all registered modules
113  * relevant to the mdp inputs, to find those that have such tasks.
114  *
115  * \todo Later we might need the concept of computeTasksOnThisRank,
116  * from which we construct gpuTasksOnThisRank.
117  *
118  * Currently the DD code assigns duty to ranks that can
119  * include PP work that currently can be executed on a single
120  * GPU, if present and compatible.  This has to be coordinated
121  * across PP ranks on a node, with possible multiple devices
122  * or sharing devices on a node, either from the user
123  * selection, or automatically. */
124 class GpuTaskAssignmentsBuilder
125 {
126 public:
127     //! Constructor
128     GpuTaskAssignmentsBuilder();
129
130     /*! \brief Builds a GpuTaskAssignments
131      *
132      * This method reconciles
133      *
134      *   - user mdrun command-line options,
135      *   - the results of hardware detection
136      *   - the duty assigned by the DD setup,
137      *   - the requested simulation modules, and
138      *   - the possible existence of multi-simulations
139      *
140      * to assign the GPUs on each physical node to the tasks on
141      * the ranks of that node. It throws InconsistentInputError
142      * when a/the useful GPU task assignment is not possible.
143      *
144      * \param[in]  availableDevices       The compatible devices that the user permitted us to use.
145      * \param[in]  userGpuTaskAssignment  The user-specified assignment of GPU tasks to device IDs.
146      * \param[in]  hardwareInfo           The detected hardware
147      * \param[in]  gromacsWorldComm       MPI communicator for all ranks in the current GROMACS run
148      * \param[in]  physicalNodeComm       Communication object for this physical node.
149      * \param[in]  nonbondedTarget        The user's choice for mdrun -nb for where to assign
150      *                                    short-ranged nonbonded interaction tasks.
151      * \param[in]  pmeTarget              The user's choice for mdrun -pme for where to assign
152      *                                    long-ranged PME nonbonded interaction tasks.
153      * \param[in]  bondedTarget           The user's choice for mdrun -bonded for where to assign tasks.
154      * \param[in]  updateTarget           The user's choice for mdrun -update for where to assign tasks.
155      * \param[in]  useGpuForNonbonded     Whether GPUs will be used for nonbonded interactions.
156      * \param[in]  useGpuForPme           Whether GPUs will be used for PME interactions.
157      * \param[in]  rankHasPpTask          Whether this rank has a PP task
158      * \param[in]  rankHasPmeTask         Whether this rank has a PME task
159      *
160      * \throws   std::bad_alloc          If out of memory.
161      *           InconsistentInputError  If user and/or detected inputs are inconsistent.
162      */
163     static GpuTaskAssignments build(ArrayRef<const int>             availableDevices,
164                                     ArrayRef<const int>             userGpuTaskAssignment,
165                                     const gmx_hw_info_t&            hardwareInfo,
166                                     MPI_Comm                        gromacsWorldComm,
167                                     const PhysicalNodeCommunicator& physicalNodeComm,
168                                     TaskTarget                      nonbondedTarget,
169                                     TaskTarget                      pmeTarget,
170                                     TaskTarget                      bondedTarget,
171                                     TaskTarget                      updateTarget,
172                                     bool                            useGpuForNonbonded,
173                                     bool                            useGpuForPme,
174                                     bool                            rankHasPpTask,
175                                     bool                            rankHasPmeTask);
176 };
177
178 /*! \libinternal
179  * \brief Contains the GPU task assignment for all ranks on this
180  * physical node.
181  *
182  * This can be used to configure the modules that might run tasks on
183  * GPUs.
184  *
185  * This assignment is made by a GpuTaskAssignmentsBuilder object. */
186 class GpuTaskAssignments
187 {
188 public:
189     //! Public move constructor to use with the builder
190     GpuTaskAssignments(GpuTaskAssignments&& source) noexcept = default;
191
192 private:
193     // Let the builder handle construction
194     friend class GpuTaskAssignmentsBuilder;
195     //! Private constructor so only the builder can construct
196     GpuTaskAssignments(const gmx_hw_info_t& hardwareInfo);
197     /*! \brief Information about hardware on this physical node
198      *
199      * The lifetime of the object referred to must exceed that
200      * of this object. */
201     const gmx_hw_info_t& hardwareInfo_;
202     //! The GPU task assignment for all ranks on this node
203     std::vector<GpuTaskAssignment> assignmentForAllRanksOnThisNode_;
204     /*! \brief The index of this rank within those on this node.
205      *
206      * This is useful for indexing into \c
207      * assignmentForAllRanksOnThisNode_. */
208     index indexOfThisRank_ = -1;
209     //! Number of GPU tasks on this node.
210     size_t numGpuTasksOnThisNode_ = 0;
211     //! Number of ranks on this physical node.
212     size_t numRanksOnThisNode_ = 0;
213
214     //! Vector of device IDs assigned to this node
215     std::vector<int> deviceIdsAssigned_;
216
217 public:
218     /*! \brief Log a report on how GPUs are being used on
219      * the ranks of the physical node of rank 0 of the simulation.
220      *
221      * \todo It could be useful to report also whether any nodes differed,
222      * and in what way.
223      *
224      * \param[in]  mdlog           Logging object.
225      * \param[in]  printHostName   Print the hostname in the usage information.
226      * \param[in]  useGpuForBonded Whether GPU PP tasks will do bonded work on the GPU.
227      * \param[in]  pmeRunMode      Describes the execution of PME tasks.
228      * \param[in]  useGpuForUpdate Whether the update is offloaded on the GPU.
229      *
230      * \throws     std::bad_alloc if out of memory
231      */
232     void reportGpuUsage(const MDLogger& mdlog,
233                         bool            printHostName,
234                         bool            useGpuForBonded,
235                         PmeRunMode      pmeRunMode,
236                         bool            useGpuForUpdate);
237
238     /*! \brief Logs to \c mdlog information that may help a user
239      * learn how to let mdrun make a task assignment that runs
240      * faster.
241      *
242      * \param[in]  mdlog                         Logging object.
243      * \param[in]  numAvailableDevicesOnThisNode The number of compatible devices on this node
244      *                                           that the user permitted us to use.
245      * */
246     void logPerformanceHints(const MDLogger& mdlog, size_t numAvailableDevicesOnThisNode);
247     /*! \brief Return handle to the initialized GPU to use in the this rank.
248      *
249      * \param[out] deviceId Index of the assigned device.
250      *
251      * \returns Device information on the selected devicce. Returns nullptr if no GPU task
252      *          is assigned to this rank.
253      */
254     DeviceInformation* initDevice(int* deviceId) const;
255     //! Return whether this rank has a PME task running on a GPU
256     bool thisRankHasPmeGpuTask() const;
257     //! Return whether this rank has any task running on a GPU
258     bool thisRankHasAnyGpuTask() const;
259     //! Get the list of unique devices that have been assigned tasks on this physical node
260     std::vector<int> deviceIdsAssigned() { return deviceIdsAssigned_; }
261 };
262
263 } // namespace gmx
264
265 #endif