Update GPU update restrictions
[alexxy/gromacs.git] / src / gromacs / taskassignment / decidegpuusage.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2015,2016,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 /*! \internal \file
36  * \brief Defines functionality for deciding whether tasks will run on GPUs.
37  *
38  * \author Mark Abraham <mark.j.abraham@gmail.com>
39  * \ingroup module_taskassignment
40  */
41
42 #include "gmxpre.h"
43
44 #include "decidegpuusage.h"
45
46 #include "config.h"
47
48 #include <cstdlib>
49 #include <cstring>
50
51 #include <algorithm>
52 #include <string>
53
54 #include "gromacs/ewald/pme.h"
55 #include "gromacs/hardware/cpuinfo.h"
56 #include "gromacs/hardware/detecthardware.h"
57 #include "gromacs/hardware/hardwaretopology.h"
58 #include "gromacs/hardware/hw_info.h"
59 #include "gromacs/mdlib/gmx_omp_nthreads.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/mdtypes/mdrunoptions.h"
64 #include "gromacs/taskassignment/taskassignment.h"
65 #include "gromacs/topology/topology.h"
66 #include "gromacs/utility/baseversion.h"
67 #include "gromacs/utility/exceptions.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/gmxassert.h"
70 #include "gromacs/utility/logger.h"
71 #include "gromacs/utility/stringutil.h"
72
73
74 namespace gmx
75 {
76
77 namespace
78 {
79
80 //! Helper variable to localise the text of an often repeated message.
81 const char * g_specifyEverythingFormatString =
82     "When you use mdrun -gputasks, %s must be set to non-default "
83     "values, so that the device IDs can be interpreted correctly."
84 #if GMX_GPU != GMX_GPU_NONE
85     " If you simply want to restrict which GPUs are used, then it is "
86     "better to use mdrun -gpu_id. Otherwise, setting the "
87 #  if GMX_GPU == GMX_GPU_CUDA
88     "CUDA_VISIBLE_DEVICES"
89 #  elif GMX_GPU == GMX_GPU_OPENCL
90     // Technically there is no portable way to do this offered by the
91     // OpenCL standard, but the only current relevant case for GROMACS
92     // is AMD OpenCL, which offers this variable.
93     "GPU_DEVICE_ORDINAL"
94 #  else
95 #  error "Unreachable branch"
96 #  endif
97     " environment variable in your bash profile or job "
98     "script may be more convenient."
99 #endif
100 ;
101
102 }   // namespace
103
104 bool
105 decideWhetherToUseGpusForNonbondedWithThreadMpi(const TaskTarget          nonbondedTarget,
106                                                 const std::vector<int>   &gpuIdsToUse,
107                                                 const std::vector<int>   &userGpuTaskAssignment,
108                                                 const EmulateGpuNonbonded emulateGpuNonbonded,
109                                                 const bool                buildSupportsNonbondedOnGpu,
110                                                 const bool                nonbondedOnGpuIsUseful,
111                                                 const int                 numRanksPerSimulation)
112 {
113     // First, exclude all cases where we can't run NB on GPUs.
114     if (nonbondedTarget == TaskTarget::Cpu ||
115         emulateGpuNonbonded == EmulateGpuNonbonded::Yes ||
116         !nonbondedOnGpuIsUseful ||
117         !buildSupportsNonbondedOnGpu)
118     {
119         // If the user required NB on GPUs, we issue an error later.
120         return false;
121     }
122
123     // We now know that NB on GPUs makes sense, if we have any.
124
125     if (!userGpuTaskAssignment.empty())
126     {
127         // Specifying -gputasks requires specifying everything.
128         if (nonbondedTarget == TaskTarget::Auto ||
129             numRanksPerSimulation < 1)
130         {
131             GMX_THROW(InconsistentInputError(formatString(g_specifyEverythingFormatString, "-nb and -ntmpi")));
132         }
133         return true;
134     }
135
136     if (nonbondedTarget == TaskTarget::Gpu)
137     {
138         return true;
139     }
140
141     // Because this is thread-MPI, we already know about the GPUs that
142     // all potential ranks can use, and can use that in a global
143     // decision that will later be consistent.
144     auto haveGpus = !gpuIdsToUse.empty();
145
146     // If we get here, then the user permitted or required GPUs.
147     return haveGpus;
148 }
149
150 bool
151 decideWhetherToUseGpusForPmeWithThreadMpi(const bool              useGpuForNonbonded,
152                                           const TaskTarget        pmeTarget,
153                                           const std::vector<int> &gpuIdsToUse,
154                                           const std::vector<int> &userGpuTaskAssignment,
155                                           const gmx_hw_info_t    &hardwareInfo,
156                                           const t_inputrec       &inputrec,
157                                           const gmx_mtop_t       &mtop,
158                                           const int               numRanksPerSimulation,
159                                           const int               numPmeRanksPerSimulation)
160 {
161     // First, exclude all cases where we can't run PME on GPUs.
162     if ((pmeTarget == TaskTarget::Cpu) ||
163         !useGpuForNonbonded ||
164         !pme_gpu_supports_build(nullptr) ||
165         !pme_gpu_supports_hardware(hardwareInfo, nullptr) ||
166         !pme_gpu_supports_input(inputrec, mtop, nullptr))
167     {
168         // PME can't run on a GPU. If the user required that, we issue
169         // an error later.
170         return false;
171     }
172
173     // We now know that PME on GPUs might make sense, if we have any.
174
175     if (!userGpuTaskAssignment.empty())
176     {
177         // Follow the user's choice of GPU task assignment, if we
178         // can. Checking that their IDs are for compatible GPUs comes
179         // later.
180
181         // Specifying -gputasks requires specifying everything.
182         if (pmeTarget == TaskTarget::Auto ||
183             numRanksPerSimulation < 1)
184         {
185             GMX_THROW(InconsistentInputError(formatString(g_specifyEverythingFormatString, "all of -nb, -pme, and -ntmpi")));
186         }
187
188         // PME on GPUs is only supported in a single case
189         if (pmeTarget == TaskTarget::Gpu)
190         {
191             if (((numRanksPerSimulation > 1) && (numPmeRanksPerSimulation == 0)) ||
192                 (numPmeRanksPerSimulation > 1))
193             {
194                 GMX_THROW(InconsistentInputError
195                               ("When you run mdrun -pme gpu -gputasks, you must supply a PME-enabled .tpr file and use a single PME rank."));
196             }
197             return true;
198         }
199
200         // pmeTarget == TaskTarget::Auto
201         return numRanksPerSimulation == 1;
202     }
203
204     // Because this is thread-MPI, we already know about the GPUs that
205     // all potential ranks can use, and can use that in a global
206     // decision that will later be consistent.
207
208     if (pmeTarget == TaskTarget::Gpu)
209     {
210         if (((numRanksPerSimulation > 1) && (numPmeRanksPerSimulation == 0)) ||
211             (numPmeRanksPerSimulation > 1))
212         {
213             GMX_THROW(NotImplementedError
214                           ("PME tasks were required to run on GPUs, but that is not implemented with "
215                           "more than one PME rank. Use a single rank simulation, or a separate PME rank, "
216                           "or permit PME tasks to be assigned to the CPU."));
217         }
218         return true;
219     }
220
221     if (numRanksPerSimulation == 1)
222     {
223         // PME can run well on a GPU shared with NB, and we permit
224         // mdrun to default to try that.
225         return !gpuIdsToUse.empty();
226     }
227
228     if (numRanksPerSimulation < 1)
229     {
230         // Full automated mode for thread-MPI (the default). PME can
231         // run well on a GPU shared with NB, and we permit mdrun to
232         // default to it if there is only one GPU available.
233         return (gpuIdsToUse.size() == 1);
234     }
235
236     // Not enough support for PME on GPUs for anything else
237     return false;
238 }
239
240 bool decideWhetherToUseGpusForNonbonded(const TaskTarget           nonbondedTarget,
241                                         const std::vector<int>    &userGpuTaskAssignment,
242                                         const EmulateGpuNonbonded  emulateGpuNonbonded,
243                                         const bool                 buildSupportsNonbondedOnGpu,
244                                         const bool                 nonbondedOnGpuIsUseful,
245                                         const bool                 gpusWereDetected)
246 {
247     if (nonbondedTarget == TaskTarget::Cpu)
248     {
249         if (!userGpuTaskAssignment.empty())
250         {
251             GMX_THROW(InconsistentInputError
252                           ("A GPU task assignment was specified, but nonbonded interactions were "
253                           "assigned to the CPU. Make no more than one of these choices."));
254         }
255
256         return false;
257     }
258
259     if (!buildSupportsNonbondedOnGpu && nonbondedTarget == TaskTarget::Gpu)
260     {
261         GMX_THROW(InconsistentInputError
262                       ("Nonbonded interactions on the GPU were requested with -nb gpu, "
263                       "but the GROMACS binary has been built without GPU support. "
264                       "Either run without selecting GPU options, or recompile GROMACS "
265                       "with GPU support enabled"));
266     }
267
268     // TODO refactor all these TaskTarget::Gpu checks into one place?
269     // e.g. use a subfunction that handles only the cases where
270     // TaskTargets are not Cpu?
271     if (emulateGpuNonbonded == EmulateGpuNonbonded::Yes)
272     {
273         if (nonbondedTarget == TaskTarget::Gpu)
274         {
275             GMX_THROW(InconsistentInputError
276                           ("Nonbonded interactions on the GPU were required, which is inconsistent "
277                           "with choosing emulation. Make no more than one of these choices."));
278         }
279         if (!userGpuTaskAssignment.empty())
280         {
281             GMX_THROW(InconsistentInputError
282                           ("GPU ID usage was specified, as was GPU emulation. Make no more than one of these choices."));
283         }
284
285         return false;
286     }
287
288     if (!nonbondedOnGpuIsUseful)
289     {
290         if (nonbondedTarget == TaskTarget::Gpu)
291         {
292             GMX_THROW(InconsistentInputError
293                           ("Nonbonded interactions on the GPU were required, but not supported for these "
294                           "simulation settings. Change your settings, or do not require using GPUs."));
295         }
296
297         return false;
298     }
299
300     if (!userGpuTaskAssignment.empty())
301     {
302         // Specifying -gputasks requires specifying everything.
303         if (nonbondedTarget == TaskTarget::Auto)
304         {
305             GMX_THROW(InconsistentInputError(formatString(g_specifyEverythingFormatString, "-nb and -ntmpi")));
306         }
307
308         return true;
309     }
310
311     if (nonbondedTarget == TaskTarget::Gpu)
312     {
313         // We still don't know whether it is an error if no GPUs are found
314         // because we don't know the duty of this rank, yet. For example,
315         // a node with only PME ranks and -pme cpu is OK if there are not
316         // GPUs.
317         return true;
318     }
319
320     // If we get here, then the user permitted GPUs, which we should
321     // use for nonbonded interactions.
322     return gpusWereDetected;
323 }
324
325 bool decideWhetherToUseGpusForPme(const bool              useGpuForNonbonded,
326                                   const TaskTarget        pmeTarget,
327                                   const std::vector<int> &userGpuTaskAssignment,
328                                   const gmx_hw_info_t    &hardwareInfo,
329                                   const t_inputrec       &inputrec,
330                                   const gmx_mtop_t       &mtop,
331                                   const int               numRanksPerSimulation,
332                                   const int               numPmeRanksPerSimulation,
333                                   const bool              gpusWereDetected)
334 {
335     if (pmeTarget == TaskTarget::Cpu)
336     {
337         return false;
338     }
339
340     if (!useGpuForNonbonded)
341     {
342         if (pmeTarget == TaskTarget::Gpu)
343         {
344             GMX_THROW(NotImplementedError
345                           ("PME on GPUs is only supported when nonbonded interactions run on GPUs also."));
346         }
347         return false;
348     }
349
350     std::string message;
351     if (!pme_gpu_supports_build(&message))
352     {
353         if (pmeTarget == TaskTarget::Gpu)
354         {
355             GMX_THROW(NotImplementedError
356                           ("Cannot compute PME interactions on a GPU, because " + message));
357         }
358         return false;
359     }
360     if (!pme_gpu_supports_hardware(hardwareInfo, &message))
361     {
362         if (pmeTarget == TaskTarget::Gpu)
363         {
364             GMX_THROW(NotImplementedError
365                           ("Cannot compute PME interactions on a GPU, because " + message));
366         }
367         return false;
368     }
369     if (!pme_gpu_supports_input(inputrec, mtop, &message))
370     {
371         if (pmeTarget == TaskTarget::Gpu)
372         {
373             GMX_THROW(NotImplementedError
374                           ("Cannot compute PME interactions on a GPU, because " + message));
375         }
376         return false;
377     }
378
379     if (pmeTarget == TaskTarget::Cpu)
380     {
381         if (!userGpuTaskAssignment.empty())
382         {
383             GMX_THROW(InconsistentInputError
384                           ("A GPU task assignment was specified, but PME interactions were "
385                           "assigned to the CPU. Make no more than one of these choices."));
386         }
387
388         return false;
389     }
390
391     if (!userGpuTaskAssignment.empty())
392     {
393         // Specifying -gputasks requires specifying everything.
394         if (pmeTarget == TaskTarget::Auto)
395         {
396             GMX_THROW(InconsistentInputError(formatString(g_specifyEverythingFormatString, "all of -nb, -pme, and -ntmpi"))); // TODO ntmpi?
397         }
398
399         return true;
400     }
401
402     // We still don't know whether it is an error if no GPUs are found
403     // because we don't know the duty of this rank, yet. For example,
404     // a node with only PME ranks and -pme cpu is OK if there are not
405     // GPUs.
406
407     if (pmeTarget == TaskTarget::Gpu)
408     {
409         if (((numRanksPerSimulation > 1) && (numPmeRanksPerSimulation == 0)) ||
410             (numPmeRanksPerSimulation > 1))
411         {
412             GMX_THROW(NotImplementedError
413                           ("PME tasks were required to run on GPUs, but that is not implemented with "
414                           "more than one PME rank. Use a single rank simulation, or a separate PME rank, "
415                           "or permit PME tasks to be assigned to the CPU."));
416         }
417         return true;
418     }
419
420     // If we get here, then the user permitted GPUs.
421     if (numRanksPerSimulation == 1)
422     {
423         // PME can run well on a single GPU shared with NB when there
424         // is one rank, so we permit mdrun to try that if we have
425         // detected GPUs.
426         return gpusWereDetected;
427     }
428
429     // Not enough support for PME on GPUs for anything else
430     return false;
431 }
432
433 bool decideWhetherToUseGpusForBonded(const bool       useGpuForNonbonded,
434                                      const bool       useGpuForPme,
435                                      const TaskTarget bondedTarget,
436                                      const bool       canUseGpuForBonded,
437                                      const bool       usingLJPme,
438                                      const bool       usingElecPmeOrEwald,
439                                      const int        numPmeRanksPerSimulation,
440                                      const bool       gpusWereDetected)
441 {
442     if (bondedTarget == TaskTarget::Cpu)
443     {
444         return false;
445     }
446
447     if (!canUseGpuForBonded)
448     {
449         if (bondedTarget == TaskTarget::Gpu)
450         {
451             GMX_THROW(InconsistentInputError
452                           ("Bonded interactions on the GPU were required, but not supported for these "
453                           "simulation settings. Change your settings, or do not require using GPUs."));
454         }
455
456         return false;
457     }
458
459     if (!useGpuForNonbonded)
460     {
461         if (bondedTarget == TaskTarget::Gpu)
462         {
463             GMX_THROW(InconsistentInputError
464                           ("Bonded interactions on the GPU were required, but this requires that "
465                           "short-ranged non-bonded interactions are also run on the GPU. Change "
466                           "your settings, or do not require using GPUs."));
467         }
468
469         return false;
470     }
471
472     // TODO If the bonded kernels do not get fused, then performance
473     // overheads might suggest alternative choices here.
474
475     if (bondedTarget == TaskTarget::Gpu)
476     {
477         // We still don't know whether it is an error if no GPUs are
478         // found.
479         return true;
480     }
481
482     // If we get here, then the user permitted GPUs, which we should
483     // use for bonded interactions if any were detected and the CPU
484     // is busy, for which we currently only check PME or Ewald.
485     // (It would be better to dynamically assign bondeds based on timings)
486     // Note that here we assume that the auto setting of PME ranks will not
487     // choose seperate PME ranks when nonBonded are assigned to the GPU.
488     bool usingOurCpuForPmeOrEwald = (usingLJPme || (usingElecPmeOrEwald && !useGpuForPme && numPmeRanksPerSimulation <= 0));
489
490     return gpusWereDetected && usingOurCpuForPmeOrEwald;
491 }
492
493 bool decideWhetherToUseGpuForUpdate(const bool        forceGpuUpdateDefaultOn,
494                                     const bool        isDomainDecomposition,
495                                     const bool        useGpuForPme,
496                                     const bool        useGpuForNonbonded,
497                                     const TaskTarget  updateTarget,
498                                     const bool        gpusWereDetected,
499                                     const t_inputrec &inputrec,
500                                     const bool        haveVSites,
501                                     const bool        useEssentialDynamics,
502                                     const bool        doOrientationRestraints,
503                                     const bool        useReplicaExchange)
504 {
505
506     if (updateTarget == TaskTarget::Cpu)
507     {
508         return false;
509     }
510
511     std::string errorMessage;
512
513     if (isDomainDecomposition)
514     {
515         errorMessage += "Domain decomposition is not supported.\n";
516     }
517     // Using the GPU-version of update if:
518     // 1. PME is on the GPU (there should be a copy of coordinates on GPU for PME spread), or
519     // 2. Non-bonded interactions are on the GPU.
520     if (!(useGpuForPme || useGpuForNonbonded))
521     {
522         errorMessage += "Either PME or short-ranged non-bonded interaction tasks must run on the GPU.\n";
523     }
524     if (!gpusWereDetected)
525     {
526         errorMessage += "Compatible GPUs must have been found.\n";
527     }
528     if (GMX_GPU != GMX_GPU_CUDA)
529     {
530         errorMessage += "Only a CUDA build is supported.\n";
531     }
532     if (inputrec.eI != eiMD)
533     {
534         errorMessage += "Only the md integrator is supported.\n";
535     }
536     if (inputrec.etc == etcNOSEHOOVER)
537     {
538         errorMessage += "Nose-Hoover temperature coupling is not supported.\n";
539     }
540     if (inputrec.epc != epcNO)
541     {
542         // Coordinate D2H and H2d are missing as well as PBC reinitialization
543         errorMessage += "Pressure coupling is not supported.\n";
544     }
545     if (EEL_PME_EWALD(inputrec.coulombtype) && inputrec.epsilon_surface != 0)
546     {
547         // The graph is needed, but not supported
548         errorMessage += "Ewald surface correction is not supported.\n";
549     }
550     if (haveVSites)
551     {
552         errorMessage += "Virtual sites are not supported.\n";
553     }
554     if (useEssentialDynamics)
555     {
556         errorMessage += "Essential dynamics is not supported.\n";
557     }
558     if (inputrec.bPull || inputrec.pull)
559     {
560         // Pull potentials are actually supported, but constraint pulling is not
561         errorMessage += "Pulling is not supported.\n";
562     }
563     if (doOrientationRestraints)
564     {
565         // The graph is needed, but not supported
566         errorMessage += "Orientation restraints are not supported.\n";
567     }
568     if (inputrec.efep != efepNO)
569     {
570         // Actually all free-energy options except for mass and constraint perturbation are supported
571         errorMessage += "Free energy perturbations are not supported.\n";
572     }
573     if (useReplicaExchange)
574     {
575         errorMessage += "Replica exchange simulations are not supported.\n";
576     }
577     if (inputrec.eSwapCoords != eswapNO)
578     {
579         errorMessage += "Swapping the coordinates is not supported.\n";
580     }
581     if (!errorMessage.empty())
582     {
583         if (updateTarget == TaskTarget::Gpu)
584         {
585             std::string prefix = gmx::formatString("Update task on the GPU was required,\n"
586                                                    "but the following condition(s) were not satisfied:\n");
587             GMX_THROW(InconsistentInputError((prefix + errorMessage).c_str()));
588         }
589         return false;
590     }
591
592     return ((forceGpuUpdateDefaultOn && updateTarget == TaskTarget::Auto) || (updateTarget == TaskTarget::Gpu));
593 }
594
595 }  // namespace gmx