Separate management of GPU contexts from modules
[gromacs.git] / src / gromacs / taskassignment / decidegpuusage.cpp
blobd20ca674d6de3be6bbd11d42bf1f5079d9fb45fe
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,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.
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 /*! \internal \file
36 * \brief Defines functionality for deciding whether tasks will run on GPUs.
38 * \author Mark Abraham <mark.j.abraham@gmail.com>
39 * \ingroup module_taskassignment
42 #include "gmxpre.h"
44 #include "decidegpuusage.h"
46 #include "config.h"
48 #include <stdlib.h>
49 #include <string.h>
51 #include <algorithm>
52 #include <string>
54 #include "gromacs/hardware/cpuinfo.h"
55 #include "gromacs/hardware/detecthardware.h"
56 #include "gromacs/hardware/hardwaretopology.h"
57 #include "gromacs/hardware/hw_info.h"
58 #include "gromacs/mdlib/gmx_omp_nthreads.h"
59 #include "gromacs/mdlib/nb_verlet.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/taskassignment/taskassignment.h"
64 #include "gromacs/topology/topology.h"
65 #include "gromacs/utility/baseversion.h"
66 #include "gromacs/utility/exceptions.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/gmxassert.h"
69 #include "gromacs/utility/logger.h"
70 #include "gromacs/utility/stringutil.h"
73 namespace gmx
76 namespace
79 //! Helper variable to localise the text of an often repeated message.
80 const char * g_specifyEverythingFormatString =
81 "When you use mdrun -gputasks, %s must be set to non-default "
82 "values, so that the device IDs can be interpreted correctly."
83 #if GMX_GPU != GMX_GPU_NONE
84 " If you simply want to restrict which GPUs are used, then it is "
85 "better to use mdrun -gpu_id. Otherwise, setting the "
86 # if GMX_GPU == GMX_GPU_CUDA
87 "CUDA_VISIBLE_DEVICES"
88 # elif GMX_GPU == GMX_GPU_OPENCL
89 // Technically there is no portable way to do this offered by the
90 // OpenCL standard, but the only current relevant case for GROMACS
91 // is AMD OpenCL, which offers this variable.
92 "GPU_DEVICE_ORDINAL"
93 # else
94 # error "Unreachable branch"
95 # endif
96 " environment variable in your bash profile or job "
97 "script may be more convenient."
98 #endif
101 } // namespace
103 bool
104 decideWhetherToUseGpusForNonbondedWithThreadMpi(const TaskTarget nonbondedTarget,
105 const std::vector<int> &gpuIdsToUse,
106 const std::vector<int> &userGpuTaskAssignment,
107 const EmulateGpuNonbonded emulateGpuNonbonded,
108 const bool usingVerletScheme,
109 const bool nonbondedOnGpuIsUseful,
110 const int numRanksPerSimulation)
112 // First, exclude all cases where we can't run NB on GPUs.
113 if (nonbondedTarget == TaskTarget::Cpu ||
114 emulateGpuNonbonded == EmulateGpuNonbonded::Yes ||
115 !usingVerletScheme ||
116 !nonbondedOnGpuIsUseful)
118 // If the user required NB on GPUs, we issue an error later.
119 return false;
122 // We now know that NB on GPUs makes sense, if we have any.
124 if (!userGpuTaskAssignment.empty())
126 // Specifying -gputasks requires specifying everything.
127 if (nonbondedTarget == TaskTarget::Auto ||
128 numRanksPerSimulation < 1)
130 GMX_THROW(InconsistentInputError(formatString(g_specifyEverythingFormatString, "-nb and -ntmpi")));
132 return true;
135 if (nonbondedTarget == TaskTarget::Gpu)
137 return true;
140 // Because this is thread-MPI, we already know about the GPUs that
141 // all potential ranks can use, and can use that in a global
142 // decision that will later be consistent.
143 auto haveGpus = !gpuIdsToUse.empty();
145 // If we get here, then the user permitted or required GPUs.
146 return haveGpus;
149 bool
150 decideWhetherToUseGpusForPmeWithThreadMpi(const bool useGpuForNonbonded,
151 const TaskTarget pmeTarget,
152 const std::vector<int> &gpuIdsToUse,
153 const std::vector<int> &userGpuTaskAssignment,
154 const bool canUseGpuForPme,
155 const int numRanksPerSimulation)
157 // First, exclude all cases where we can't run PME on GPUs.
158 if ((pmeTarget == TaskTarget::Cpu) ||
159 !useGpuForNonbonded ||
160 !canUseGpuForPme)
162 // PME can't run on a GPU. If the user required that, we issue
163 // an error later.
164 return false;
167 // We now know that PME on GPUs might make sense, if we have any.
169 if (!userGpuTaskAssignment.empty())
171 // Follow the user's choice of GPU task assignment, if we
172 // can. Checking that their IDs are for compatible GPUs comes
173 // later.
175 // Specifying -gputasks requires specifying everything.
176 if (pmeTarget == TaskTarget::Auto ||
177 numRanksPerSimulation < 1)
179 GMX_THROW(InconsistentInputError(formatString(g_specifyEverythingFormatString, "all of -nb, -pme, and -ntmpi")));
182 // PME on GPUs is only supported in a single case
183 if (pmeTarget == TaskTarget::Gpu)
185 if (numRanksPerSimulation > 1)
187 GMX_THROW(InconsistentInputError
188 ("When you run mdrun -pme gpu -gputasks, you must supply a PME .tpr file and use a single rank."));
190 return true;
193 // pmeTarget == TaskTarget::Auto
194 return numRanksPerSimulation == 1;
197 // Because this is thread-MPI, we already know about the GPUs that
198 // all potential ranks can use, and can use that in a global
199 // decision that will later be consistent.
201 if (pmeTarget == TaskTarget::Gpu)
203 if (numRanksPerSimulation > 1)
205 GMX_THROW(NotImplementedError
206 ("PME tasks were required to run on GPUs, but that is not implemented with "
207 "more than one rank. Use a single rank, or permit PME tasks to be assigned "
208 "to the CPU."));
210 return true;
213 if (numRanksPerSimulation == 1)
215 // PME can run well on a GPU shared with NB, and we permit
216 // mdrun to default to try that.
217 return gpuIdsToUse.size() >= 1;
220 if (numRanksPerSimulation < 1)
222 // Full automated mode for thread-MPI (the default). PME can
223 // run well on a GPU shared with NB, and we permit mdrun to
224 // default to it if there is only one GPU available.
225 return (gpuIdsToUse.size() == 1);
228 // Not enough support for PME on GPUs for anything else
229 return false;
232 bool decideWhetherToUseGpusForNonbonded(const TaskTarget nonbondedTarget,
233 const std::vector<int> &gpuIdsToUse,
234 const std::vector<int> &userGpuTaskAssignment,
235 const EmulateGpuNonbonded emulateGpuNonbonded,
236 const bool usingVerletScheme,
237 const bool nonbondedOnGpuIsUseful)
239 if (nonbondedTarget == TaskTarget::Cpu)
241 if (!userGpuTaskAssignment.empty())
243 GMX_THROW(InconsistentInputError
244 ("A GPU task assignment was specified, but nonbonded interactions were "
245 "assigned to the CPU. Make no more than one of these choices."));
248 return false;
251 // TODO refactor all these TaskTarget::Gpu checks into one place?
252 // e.g. use a subfunction that handles only the cases where
253 // TaskTargets are not Cpu?
254 if (emulateGpuNonbonded == EmulateGpuNonbonded::Yes)
256 if (nonbondedTarget == TaskTarget::Gpu)
258 GMX_THROW(InconsistentInputError
259 ("Nonbonded interactions on the GPU were required, which is inconsistent "
260 "with choosing emulation. Make no more than one of these choices."));
262 if (!gpuIdsToUse.empty() || !userGpuTaskAssignment.empty())
264 GMX_THROW(InconsistentInputError
265 ("GPU ID usage was specified, as was GPU emulation. Make no more than one of these choices."));
268 return false;
271 if (!usingVerletScheme)
273 if (nonbondedTarget == TaskTarget::Gpu)
275 GMX_THROW(InconsistentInputError
276 ("Nonbonded interactions on the GPU were required, which requires using "
277 "the Verlet scheme. Either use the Verlet scheme, or do not require using GPUs."));
280 return false;
283 if (!nonbondedOnGpuIsUseful)
285 if (nonbondedTarget == TaskTarget::Gpu)
287 GMX_THROW(InconsistentInputError
288 ("Nonbonded interactions on the GPU were required, but this would not be "
289 "useful. Probably you should not require using GPUs."));
292 return false;
295 if (!userGpuTaskAssignment.empty())
297 // Specifying -gputasks requires specifying everything.
298 if (nonbondedTarget == TaskTarget::Auto)
300 GMX_THROW(InconsistentInputError(formatString(g_specifyEverythingFormatString, "-nb and -ntmpi")));
303 return true;
306 // We still don't know whether it is an error if no GPUs are found
307 // because we don't know the duty of this rank, yet. For example,
308 // a node with only PME ranks and -pme cpu is OK if there are not
309 // GPUs.
311 // If we get here, then the user permitted or required GPUs.
312 return true;
315 bool decideWhetherToUseGpusForPme(const bool useGpuForNonbonded,
316 const TaskTarget pmeTarget,
317 const std::vector<int> &userGpuTaskAssignment,
318 const bool canUseGpuForPme,
319 const int numRanksPerSimulation)
321 if (pmeTarget == TaskTarget::Cpu)
323 return false;
326 if (!useGpuForNonbonded)
328 if (pmeTarget == TaskTarget::Gpu)
330 GMX_THROW(NotImplementedError
331 ("The PME on the GPU is only supported when nonbonded interactions run on GPUs also."));
333 return false;
336 if (!canUseGpuForPme)
338 if (pmeTarget == TaskTarget::Gpu)
340 // TODO Pass in the inputrec so we can give more help here?
341 GMX_THROW(NotImplementedError
342 ("The input simulation did not use PME in a way that is supported on the GPU."));
344 return false;
347 if (pmeTarget == TaskTarget::Cpu)
349 if (!userGpuTaskAssignment.empty())
351 GMX_THROW(InconsistentInputError
352 ("A GPU task assignment was specified, but PME interactions were "
353 "assigned to the CPU. Make no more than one of these choices."));
356 return false;
359 if (!userGpuTaskAssignment.empty())
361 // Specifying -gputasks requires specifying everything.
362 if (pmeTarget == TaskTarget::Auto)
364 GMX_THROW(InconsistentInputError(formatString(g_specifyEverythingFormatString, "all of -nb, -pme, and -ntmpi")));
367 return true;
370 // We still don't know whether it is an error if no GPUs are found
371 // because we don't know the duty of this rank, yet. For example,
372 // a node with only PME ranks and -pme cpu is OK if there are not
373 // GPUs.
375 if (pmeTarget == TaskTarget::Gpu)
377 if (numRanksPerSimulation > 1)
379 GMX_THROW(NotImplementedError
380 ("PME tasks were required to run on GPUs, but that is not implemented with "
381 "more than one rank. Use a single rank, or permit PME tasks to be assigned "
382 "to the CPU."));
384 return true;
387 if (numRanksPerSimulation == 1)
389 // PME can run well on a single GPU shared with NB when
390 // there is one rank, so we permit mdrun to try that.
391 return true;
394 // Not enough support for PME on GPUs for anything else
395 return false;
398 } // namespace