Merge "Merge branch 'release-2019' into master"
[alexxy/gromacs.git] / src / gromacs / taskassignment / decidegpuusage.h
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 /*! \libinternal \file
36  * \brief Declares functionality for deciding whether tasks will run on GPUs.
37  *
38  * \author Mark Abraham <mark.j.abraham@gmail.com>
39  * \ingroup module_taskassignment
40  * \inlibraryapi
41  */
42
43 #ifndef GMX_TASKASSIGNMENT_DECIDEGPUUSAGE_H
44 #define GMX_TASKASSIGNMENT_DECIDEGPUUSAGE_H
45
46 #include <vector>
47
48 struct gmx_hw_info_t;
49 struct gmx_mtop_t;
50 struct t_inputrec;
51
52 enum class EmulateGpuNonbonded : bool;
53
54 namespace gmx
55 {
56
57 //! Record where a compute task is targetted.
58 enum class TaskTarget : int
59 {
60     Auto,
61     Cpu,
62     Gpu
63 };
64
65 /*! \brief Decide whether this thread-MPI simulation will run
66  * nonbonded tasks on GPUs.
67  *
68  * The number of GPU tasks and devices influences both the choice of
69  * the number of ranks, and checks upon any such choice made by the
70  * user. So we need to consider this before any automated choice of
71  * the number of thread-MPI ranks.
72  *
73  * \param[in]  nonbondedTarget             The user's choice for mdrun -nb for where to assign short-ranged nonbonded interaction tasks.
74  * \param[in]  gpuIdsToUse                 The compatible GPUs that the user permitted us to use.
75  * \param[in]  userGpuTaskAssignment       The user-specified assignment of GPU tasks to device IDs.
76  * \param[in]  emulateGpuNonbonded         Whether we will emulate GPU calculation of nonbonded interactions.
77  * \param[in]  buildSupportsNonbondedOnGpu Whether GROMACS was built with GPU support.
78  * \param[in]  usingVerletScheme           Whether the nonbondeds are using the Verlet scheme.
79  * \param[in]  nonbondedOnGpuIsUseful    Whether computing nonbonded interactions on a GPU is useful for this calculation.
80  * \param[in]  numRanksPerSimulation     The number of ranks in each simulation.
81  *
82  * \returns    Whether the simulation will run nonbonded tasks on GPUs.
83  *
84  * \throws     std::bad_alloc          If out of memory
85  *             InconsistentInputError  If the user requirements are inconsistent. */
86 bool decideWhetherToUseGpusForNonbondedWithThreadMpi(TaskTarget              nonbondedTarget,
87                                                      const std::vector<int> &gpuIdsToUse,
88                                                      const std::vector<int> &userGpuTaskAssignment,
89                                                      EmulateGpuNonbonded     emulateGpuNonbonded,
90                                                      bool                    buildSupportsNonbondedOnGpu,
91                                                      bool                    usingVerletScheme,
92                                                      bool                    nonbondedOnGpuIsUseful,
93                                                      int                     numRanksPerSimulation);
94
95 /*! \brief Decide whether this thread-MPI simulation will run
96  * PME tasks on GPUs.
97  *
98  * The number of GPU tasks and devices influences both the choice of
99  * the number of ranks, and checks upon any such choice made by the
100  * user. So we need to consider this before any automated choice of
101  * the number of thread-MPI ranks.
102  *
103  * \param[in]  useGpuForNonbonded        Whether GPUs will be used for nonbonded interactions.
104  * \param[in]  pmeTarget                 The user's choice for mdrun -pme for where to assign long-ranged PME nonbonded interaction tasks.
105  * \param[in]  gpuIdsToUse               The compatible GPUs that the user permitted us to use.
106  * \param[in]  userGpuTaskAssignment     The user-specified assignment of GPU tasks to device IDs.
107  * \param[in]  hardwareInfo              Hardware information
108  * \param[in]  inputrec                  The user input
109  * \param[in]  mtop                      Global system topology
110  * \param[in]  numRanksPerSimulation     The number of ranks in each simulation.
111  * \param[in]  numPmeRanksPerSimulation  The number of PME ranks in each simulation.
112  *
113  * \returns    Whether the simulation will run PME tasks on GPUs.
114  *
115  * \throws     std::bad_alloc          If out of memory
116  *             InconsistentInputError  If the user requirements are inconsistent. */
117 bool decideWhetherToUseGpusForPmeWithThreadMpi(bool                    useGpuForNonbonded,
118                                                TaskTarget              pmeTarget,
119                                                const std::vector<int> &gpuIdsToUse,
120                                                const std::vector<int> &userGpuTaskAssignment,
121                                                const gmx_hw_info_t    &hardwareInfo,
122                                                const t_inputrec       &inputrec,
123                                                const gmx_mtop_t       &mtop,
124                                                int                     numRanksPerSimulation,
125                                                int                     numPmeRanksPerSimulation);
126
127 /*! \brief Decide whether the simulation will try to run nonbonded
128  * tasks on GPUs.
129  *
130  * The final decision cannot be made until after the duty of the rank
131  * is known. But we need to know if nonbonded will run on GPUs for
132  * setting up DD (particularly rlist) and determining duty. If the
133  * user requires GPUs for the tasks of that duty, then it will be an
134  * error when none are found.
135  *
136  * With thread-MPI, calls have been made to
137  * decideWhetherToUseGpusForNonbondedWithThreadMpi() and
138  * decideWhetherToUseGpusForPmeWithThreadMpi() to help determine
139  * the number of ranks and run some checks, but the final
140  * decision is made in this routine, along with many more
141  * consistency checks.
142  *
143  * \param[in]  nonbondedTarget             The user's choice for mdrun -nb for where to assign short-ranged nonbonded interaction tasks.
144  * \param[in]  userGpuTaskAssignment       The user-specified assignment of GPU tasks to device IDs.
145  * \param[in]  emulateGpuNonbonded         Whether we will emulate GPU calculation of nonbonded interactions.
146  * \param[in]  buildSupportsNonbondedOnGpu Whether GROMACS was build with GPU support.
147  * \param[in]  usingVerletScheme           Whether the nonbondeds are using the Verlet scheme.
148  * \param[in]  nonbondedOnGpuIsUseful      Whether computing nonbonded interactions on a GPU is useful for this calculation.
149  * \param[in]  gpusWereDetected            Whether compatible GPUs were detected on any node.
150  *
151  * \returns    Whether the simulation will run nonbonded and PME tasks, respectively, on GPUs.
152  *
153  * \throws     std::bad_alloc          If out of memory
154  *             InconsistentInputError  If the user requirements are inconsistent. */
155 bool decideWhetherToUseGpusForNonbonded(TaskTarget              nonbondedTarget,
156                                         const std::vector<int> &userGpuTaskAssignment,
157                                         EmulateGpuNonbonded     emulateGpuNonbonded,
158                                         bool                    buildSupportsNonbondedOnGpu,
159                                         bool                    usingVerletScheme,
160                                         bool                    nonbondedOnGpuIsUseful,
161                                         bool                    gpusWereDetected);
162
163 /*! \brief Decide whether the simulation will try to run tasks of
164  * different types on GPUs.
165  *
166  * The final decision cannot be made until after the duty of the rank
167  * is known. But we need to know if nonbonded will run on GPUs for
168  * setting up DD (particularly rlist) and determining duty. If the
169  * user requires GPUs for the tasks of that duty, then it will be an
170  * error when none are found.
171  *
172  * With thread-MPI, calls have been made to
173  * decideWhetherToUseGpusForNonbondedWithThreadMpi() and
174  * decideWhetherToUseGpusForPmeWithThreadMpi() to help determine
175  * the number of ranks and run some checks, but the final
176  * decision is made in this routine, along with many more
177  * consistency checks.
178  *
179  * \param[in]  useGpuForNonbonded        Whether GPUs will be used for nonbonded interactions.
180  * \param[in]  pmeTarget                 The user's choice for mdrun -pme for where to assign long-ranged PME nonbonded interaction tasks.
181  * \param[in]  userGpuTaskAssignment     The user-specified assignment of GPU tasks to device IDs.
182  * \param[in]  hardwareInfo              Hardware information
183  * \param[in]  inputrec                  The user input
184  * \param[in]  mtop                      Global system topology
185  * \param[in]  numRanksPerSimulation     The number of ranks in each simulation.
186  * \param[in]  numPmeRanksPerSimulation  The number of PME ranks in each simulation.
187  * \param[in]  gpusWereDetected          Whether compatible GPUs were detected on any node.
188  *
189  * \returns    Whether the simulation will run nonbonded and PME tasks, respectively, on GPUs.
190  *
191  * \throws     std::bad_alloc          If out of memory
192  *             InconsistentInputError  If the user requirements are inconsistent. */
193 bool decideWhetherToUseGpusForPme(bool                    useGpuForNonbonded,
194                                   TaskTarget              pmeTarget,
195                                   const std::vector<int> &userGpuTaskAssignment,
196                                   const gmx_hw_info_t    &hardwareInfo,
197                                   const t_inputrec       &inputrec,
198                                   const gmx_mtop_t       &mtop,
199                                   int                     numRanksPerSimulation,
200                                   int                     numPmeRanksPerSimulation,
201                                   bool                    gpusWereDetected);
202
203 /*! \brief Decide whether the simulation will try to run bonded tasks on GPUs.
204  *
205  * \param[in]  useGpuForNonbonded        Whether GPUs will be used for nonbonded interactions.
206  * \param[in]  useGpuForPme              Whether GPUs will be used for PME interactions.
207  * \param[in]  usingVerletScheme         Whether the nonbondeds are using the Verlet scheme.
208  * \param[in]  bondedTarget              The user's choice for mdrun -bonded for where to assign tasks.
209  * \param[in]  canUseGpuForBonded        Whether the bonded interactions can run on a GPU
210  * \param[in]  usingLJPme                Whether Vdw interactions use LJ-PME.
211  * \param[in]  usingElecPmeOrEwald       Whether a PME or Ewald type method is used for electrostatics.
212  * \param[in]  numPmeRanksPerSimulation  The number of PME ranks in each simulation, can be -1 for auto.
213  * \param[in]  gpusWereDetected          Whether compatible GPUs were detected on any node.
214  *
215  * \returns    Whether the simulation will run bondeded tasks on GPUs.
216  *
217  * \throws     std::bad_alloc          If out of memory
218  *             InconsistentInputError  If the user requirements are inconsistent. */
219 bool decideWhetherToUseGpusForBonded(bool       useGpuForNonbonded,
220                                      bool       useGpuForPme,
221                                      bool       usingVerletScheme,
222                                      TaskTarget bondedTarget,
223                                      bool       canUseGpuForBonded,
224                                      bool       usingLJPme,
225                                      bool       usingElecPmeOrEwald,
226                                      int        numPmeRanksPerSimulation,
227                                      bool       gpusWereDetected);
228
229 }  // namespace gmx
230
231 #endif