25bf5f9bc142c3f63767bd076d9a0df4b38fbe68
[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.
5  * Copyright (c) 2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \internal \file
37  * \brief Defines functionality for deciding whether tasks will run on GPUs.
38  *
39  * \author Mark Abraham <mark.j.abraham@gmail.com>
40  * \ingroup module_taskassignment
41  */
42
43 #include "gmxpre.h"
44
45 #include "decidegpuusage.h"
46
47 #include "config.h"
48
49 #include <cstdlib>
50 #include <cstring>
51
52 #include <algorithm>
53 #include <string>
54
55 #include "gromacs/ewald/pme.h"
56 #include "gromacs/hardware/cpuinfo.h"
57 #include "gromacs/hardware/detecthardware.h"
58 #include "gromacs/hardware/hardwaretopology.h"
59 #include "gromacs/hardware/hw_info.h"
60 #include "gromacs/listed_forces/listed_forces_gpu.h"
61 #include "gromacs/mdlib/gmx_omp_nthreads.h"
62 #include "gromacs/mdlib/update_constrain_gpu.h"
63 #include "gromacs/mdtypes/commrec.h"
64 #include "gromacs/mdtypes/inputrec.h"
65 #include "gromacs/mdtypes/md_enums.h"
66 #include "gromacs/mdtypes/mdrunoptions.h"
67 #include "gromacs/pulling/pull.h"
68 #include "gromacs/taskassignment/taskassignment.h"
69 #include "gromacs/topology/mtop_util.h"
70 #include "gromacs/topology/topology.h"
71 #include "gromacs/utility/baseversion.h"
72 #include "gromacs/utility/exceptions.h"
73 #include "gromacs/utility/fatalerror.h"
74 #include "gromacs/utility/gmxassert.h"
75 #include "gromacs/utility/logger.h"
76 #include "gromacs/utility/stringutil.h"
77
78
79 namespace gmx
80 {
81
82 namespace
83 {
84
85 //! Helper variable to localise the text of an often repeated message.
86 const char* const g_specifyEverythingFormatString =
87         "When you use mdrun -gputasks, %s must be set to non-default "
88         "values, so that the device IDs can be interpreted correctly."
89 #if GMX_GPU
90         " If you simply want to restrict which GPUs are used, then it is "
91         "better to use mdrun -gpu_id. Otherwise, setting the "
92 #    if GMX_GPU_CUDA
93         "CUDA_VISIBLE_DEVICES"
94 #    elif GMX_GPU_OPENCL
95         // Technically there is no portable way to do this offered by the
96         // OpenCL standard, but the only current relevant case for GROMACS
97         // is AMD OpenCL, which offers this variable.
98         "GPU_DEVICE_ORDINAL"
99 #    elif GMX_GPU_SYCL && GMX_SYCL_DPCPP
100         // https://github.com/intel/llvm/blob/sycl/sycl/doc/EnvironmentVariables.md
101         "SYCL_DEVICE_FILTER"
102 #    elif GMX_GPU_SYCL && GMX_SYCL_HIPSYCL
103         // Not true if we use hipSYCL over CUDA or IntelLLVM, but in that case the user probably
104         // knows what they are doing.
105         // https://rocmdocs.amd.com/en/latest/Other_Solutions/Other-Solutions.html#hip-environment-variables
106         "HIP_VISIBLE_DEVICES"
107 #    else
108 #        error "Unreachable branch"
109 #    endif
110         " environment variable in your bash profile or job "
111         "script may be more convenient."
112 #endif
113         ;
114
115 } // namespace
116
117 bool decideWhetherToUseGpusForNonbondedWithThreadMpi(const TaskTarget        nonbondedTarget,
118                                                      const bool              haveAvailableDevices,
119                                                      const std::vector<int>& userGpuTaskAssignment,
120                                                      const EmulateGpuNonbonded emulateGpuNonbonded,
121                                                      const bool buildSupportsNonbondedOnGpu,
122                                                      const bool nonbondedOnGpuIsUseful,
123                                                      const int  numRanksPerSimulation)
124 {
125     // First, exclude all cases where we can't run NB on GPUs.
126     if (nonbondedTarget == TaskTarget::Cpu || emulateGpuNonbonded == EmulateGpuNonbonded::Yes
127         || !nonbondedOnGpuIsUseful || !buildSupportsNonbondedOnGpu)
128     {
129         // If the user required NB on GPUs, we issue an error later.
130         return false;
131     }
132
133     // We now know that NB on GPUs makes sense, if we have any.
134
135     if (!userGpuTaskAssignment.empty())
136     {
137         // Specifying -gputasks requires specifying everything.
138         if (nonbondedTarget == TaskTarget::Auto || numRanksPerSimulation < 1)
139         {
140             GMX_THROW(InconsistentInputError(
141                     formatString(g_specifyEverythingFormatString, "-nb and -ntmpi")));
142         }
143         return true;
144     }
145
146     if (nonbondedTarget == TaskTarget::Gpu)
147     {
148         return true;
149     }
150
151     // Because this is thread-MPI, we already know about the GPUs that
152     // all potential ranks can use, and can use that in a global
153     // decision that will later be consistent.
154     // If we get here, then the user permitted or required GPUs.
155     return haveAvailableDevices;
156 }
157
158 bool decideWhetherToUseGpusForPmeWithThreadMpi(const bool              useGpuForNonbonded,
159                                                const TaskTarget        pmeTarget,
160                                                const int               numDevicesToUse,
161                                                const std::vector<int>& userGpuTaskAssignment,
162                                                const gmx_hw_info_t&    hardwareInfo,
163                                                const t_inputrec&       inputrec,
164                                                const int               numRanksPerSimulation,
165                                                const int               numPmeRanksPerSimulation)
166 {
167     // First, exclude all cases where we can't run PME on GPUs.
168     if ((pmeTarget == TaskTarget::Cpu) || !useGpuForNonbonded || !pme_gpu_supports_build(nullptr)
169         || !pme_gpu_supports_hardware(hardwareInfo, nullptr) || !pme_gpu_supports_input(inputrec, nullptr))
170     {
171         // PME can't run on a GPU. If the user required that, we issue
172         // an error later.
173         return false;
174     }
175
176     // We now know that PME on GPUs might make sense, if we have any.
177
178     if (!userGpuTaskAssignment.empty())
179     {
180         // Follow the user's choice of GPU task assignment, if we
181         // can. Checking that their IDs are for compatible GPUs comes
182         // later.
183
184         // Specifying -gputasks requires specifying everything.
185         if (pmeTarget == TaskTarget::Auto || numRanksPerSimulation < 1)
186         {
187             GMX_THROW(InconsistentInputError(
188                     formatString(g_specifyEverythingFormatString, "all of -nb, -pme, and -ntmpi")));
189         }
190
191         // PME on GPUs is only supported in a single case
192         if (pmeTarget == TaskTarget::Gpu)
193         {
194             if (((numRanksPerSimulation > 1) && (numPmeRanksPerSimulation == 0))
195                 || (numPmeRanksPerSimulation > 1))
196             {
197                 GMX_THROW(InconsistentInputError(
198                         "When you run mdrun -pme gpu -gputasks, you must supply a PME-enabled .tpr "
199                         "file and use a single PME rank."));
200             }
201             return true;
202         }
203
204         // pmeTarget == TaskTarget::Auto
205         return numRanksPerSimulation == 1;
206     }
207
208     // Because this is thread-MPI, we already know about the GPUs that
209     // all potential ranks can use, and can use that in a global
210     // decision that will later be consistent.
211
212     if (pmeTarget == TaskTarget::Gpu)
213     {
214         if (((numRanksPerSimulation > 1) && (numPmeRanksPerSimulation == 0))
215             || (numPmeRanksPerSimulation > 1))
216         {
217             GMX_THROW(NotImplementedError(
218                     "PME tasks were required to run on GPUs, but that is not implemented with "
219                     "more than one PME rank. Use a single rank simulation, or a separate PME rank, "
220                     "or permit PME tasks to be assigned to the CPU."));
221         }
222         return true;
223     }
224
225     if (numRanksPerSimulation == 1)
226     {
227         // PME can run well on a GPU shared with NB, and we permit
228         // mdrun to default to try that.
229         return numDevicesToUse > 0;
230     }
231
232     if (numRanksPerSimulation < 1)
233     {
234         // Full automated mode for thread-MPI (the default). PME can
235         // run well on a GPU shared with NB, and we permit mdrun to
236         // default to it if there is only one GPU available.
237         return (numDevicesToUse == 1);
238     }
239
240     // Not enough support for PME on GPUs for anything else
241     return false;
242 }
243
244 bool decideWhetherToUseGpusForNonbonded(const TaskTarget          nonbondedTarget,
245                                         const std::vector<int>&   userGpuTaskAssignment,
246                                         const EmulateGpuNonbonded emulateGpuNonbonded,
247                                         const bool                buildSupportsNonbondedOnGpu,
248                                         const bool                nonbondedOnGpuIsUseful,
249                                         const bool                gpusWereDetected)
250 {
251     if (nonbondedTarget == TaskTarget::Cpu)
252     {
253         if (!userGpuTaskAssignment.empty())
254         {
255             GMX_THROW(InconsistentInputError(
256                     "A GPU task assignment was specified, but nonbonded interactions were "
257                     "assigned to the CPU. Make no more than one of these choices."));
258         }
259
260         return false;
261     }
262
263     if (!buildSupportsNonbondedOnGpu && nonbondedTarget == TaskTarget::Gpu)
264     {
265         GMX_THROW(InconsistentInputError(
266                 "Nonbonded interactions on the GPU were requested with -nb gpu, "
267                 "but the GROMACS binary has been built without GPU support. "
268                 "Either run without selecting GPU options, or recompile GROMACS "
269                 "with GPU support enabled"));
270     }
271
272     // TODO refactor all these TaskTarget::Gpu checks into one place?
273     // e.g. use a subfunction that handles only the cases where
274     // TaskTargets are not Cpu?
275     if (emulateGpuNonbonded == EmulateGpuNonbonded::Yes)
276     {
277         if (nonbondedTarget == TaskTarget::Gpu)
278         {
279             GMX_THROW(InconsistentInputError(
280                     "Nonbonded interactions on the GPU were required, which is inconsistent "
281                     "with choosing emulation. Make no more than one of these choices."));
282         }
283         if (!userGpuTaskAssignment.empty())
284         {
285             GMX_THROW(
286                     InconsistentInputError("GPU ID usage was specified, as was GPU emulation. Make "
287                                            "no more than one of these choices."));
288         }
289
290         return false;
291     }
292
293     if (!nonbondedOnGpuIsUseful)
294     {
295         if (nonbondedTarget == TaskTarget::Gpu)
296         {
297             GMX_THROW(InconsistentInputError(
298                     "Nonbonded interactions on the GPU were required, but not supported for these "
299                     "simulation settings. Change your settings, or do not require using GPUs."));
300         }
301
302         return false;
303     }
304
305     if (!userGpuTaskAssignment.empty())
306     {
307         // Specifying -gputasks requires specifying everything.
308         if (nonbondedTarget == TaskTarget::Auto)
309         {
310             GMX_THROW(InconsistentInputError(
311                     formatString(g_specifyEverythingFormatString, "-nb and -ntmpi")));
312         }
313
314         return true;
315     }
316
317     if (nonbondedTarget == TaskTarget::Gpu)
318     {
319         // We still don't know whether it is an error if no GPUs are found
320         // because we don't know the duty of this rank, yet. For example,
321         // a node with only PME ranks and -pme cpu is OK if there are not
322         // GPUs.
323         return true;
324     }
325
326     // If we get here, then the user permitted GPUs, which we should
327     // use for nonbonded interactions.
328     return buildSupportsNonbondedOnGpu && gpusWereDetected;
329 }
330
331 bool decideWhetherToUseGpusForPme(const bool              useGpuForNonbonded,
332                                   const TaskTarget        pmeTarget,
333                                   const std::vector<int>& userGpuTaskAssignment,
334                                   const gmx_hw_info_t&    hardwareInfo,
335                                   const t_inputrec&       inputrec,
336                                   const int               numRanksPerSimulation,
337                                   const int               numPmeRanksPerSimulation,
338                                   const bool              gpusWereDetected)
339 {
340     if (pmeTarget == TaskTarget::Cpu)
341     {
342         return false;
343     }
344
345     if (!useGpuForNonbonded)
346     {
347         if (pmeTarget == TaskTarget::Gpu)
348         {
349             GMX_THROW(NotImplementedError(
350                     "PME on GPUs is only supported when nonbonded interactions run on GPUs also."));
351         }
352         return false;
353     }
354
355     std::string message;
356     if (!pme_gpu_supports_build(&message))
357     {
358         if (pmeTarget == TaskTarget::Gpu)
359         {
360             GMX_THROW(NotImplementedError("Cannot compute PME interactions on a GPU, because " + message));
361         }
362         return false;
363     }
364     if (!pme_gpu_supports_hardware(hardwareInfo, &message))
365     {
366         if (pmeTarget == TaskTarget::Gpu)
367         {
368             GMX_THROW(NotImplementedError("Cannot compute PME interactions on a GPU, because " + message));
369         }
370         return false;
371     }
372     if (!pme_gpu_supports_input(inputrec, &message))
373     {
374         if (pmeTarget == TaskTarget::Gpu)
375         {
376             GMX_THROW(NotImplementedError("Cannot compute PME interactions on a GPU, because " + message));
377         }
378         return false;
379     }
380
381     if (pmeTarget == TaskTarget::Cpu)
382     {
383         if (!userGpuTaskAssignment.empty())
384         {
385             GMX_THROW(InconsistentInputError(
386                     "A GPU task assignment was specified, but PME interactions were "
387                     "assigned to the CPU. Make no more than one of these choices."));
388         }
389
390         return false;
391     }
392
393     if (!userGpuTaskAssignment.empty())
394     {
395         // Specifying -gputasks requires specifying everything.
396         if (pmeTarget == TaskTarget::Auto)
397         {
398             GMX_THROW(InconsistentInputError(formatString(
399                     g_specifyEverythingFormatString, "all of -nb, -pme, and -ntmpi"))); // TODO ntmpi?
400         }
401
402         return true;
403     }
404
405     // We still don't know whether it is an error if no GPUs are found
406     // because we don't know the duty of this rank, yet. For example,
407     // a node with only PME ranks and -pme cpu is OK if there are not
408     // GPUs.
409
410     if (pmeTarget == TaskTarget::Gpu)
411     {
412         if (((numRanksPerSimulation > 1) && (numPmeRanksPerSimulation == 0))
413             || (numPmeRanksPerSimulation > 1))
414         {
415             GMX_THROW(NotImplementedError(
416                     "PME tasks were required to run on GPUs, but that is not implemented with "
417                     "more than one PME rank. Use a single rank simulation, or a separate PME rank, "
418                     "or permit PME tasks to be assigned to the CPU."));
419         }
420         return true;
421     }
422
423     // If we get here, then the user permitted GPUs.
424     if (numRanksPerSimulation == 1)
425     {
426         // PME can run well on a single GPU shared with NB when there
427         // is one rank, so we permit mdrun to try that if we have
428         // detected GPUs.
429         return gpusWereDetected;
430     }
431
432     // Not enough support for PME on GPUs for anything else
433     return false;
434 }
435
436
437 PmeRunMode determinePmeRunMode(const bool useGpuForPme, const TaskTarget& pmeFftTarget, const t_inputrec& inputrec)
438 {
439     if (!EEL_PME(inputrec.coulombtype))
440     {
441         return PmeRunMode::None;
442     }
443
444     if (useGpuForPme)
445     {
446         if (pmeFftTarget == TaskTarget::Cpu)
447         {
448             return PmeRunMode::Mixed;
449         }
450         else
451         {
452             return PmeRunMode::GPU;
453         }
454     }
455     else
456     {
457         if (pmeFftTarget == TaskTarget::Gpu)
458         {
459             gmx_fatal(FARGS,
460                       "Assigning FFTs to GPU requires PME to be assigned to GPU as well. With PME "
461                       "on CPU you should not be using -pmefft.");
462         }
463         return PmeRunMode::CPU;
464     }
465 }
466
467 bool decideWhetherToUseGpusForBonded(bool              useGpuForNonbonded,
468                                      bool              useGpuForPme,
469                                      TaskTarget        bondedTarget,
470                                      const t_inputrec& inputrec,
471                                      const gmx_mtop_t& mtop,
472                                      int               numPmeRanksPerSimulation,
473                                      bool              gpusWereDetected)
474 {
475     if (bondedTarget == TaskTarget::Cpu)
476     {
477         return false;
478     }
479
480     std::string errorMessage;
481
482     if (!buildSupportsListedForcesGpu(&errorMessage))
483     {
484         if (bondedTarget == TaskTarget::Gpu)
485         {
486             GMX_THROW(InconsistentInputError(errorMessage.c_str()));
487         }
488
489         return false;
490     }
491
492     if (!inputSupportsListedForcesGpu(inputrec, mtop, &errorMessage))
493     {
494         if (bondedTarget == TaskTarget::Gpu)
495         {
496             GMX_THROW(InconsistentInputError(errorMessage.c_str()));
497         }
498
499         return false;
500     }
501
502     if (!useGpuForNonbonded)
503     {
504         if (bondedTarget == TaskTarget::Gpu)
505         {
506             GMX_THROW(InconsistentInputError(
507                     "Bonded interactions on the GPU were required, but this requires that "
508                     "short-ranged non-bonded interactions are also run on the GPU. Change "
509                     "your settings, or do not require using GPUs."));
510         }
511
512         return false;
513     }
514
515     // TODO If the bonded kernels do not get fused, then performance
516     // overheads might suggest alternative choices here.
517
518     if (bondedTarget == TaskTarget::Gpu)
519     {
520         // We still don't know whether it is an error if no GPUs are
521         // found.
522         return true;
523     }
524
525     // If we get here, then the user permitted GPUs, which we should
526     // use for bonded interactions if any were detected and the CPU
527     // is busy, for which we currently only check PME or Ewald.
528     // (It would be better to dynamically assign bondeds based on timings)
529     // Note that here we assume that the auto setting of PME ranks will not
530     // choose seperate PME ranks when nonBonded are assigned to the GPU.
531     bool usingOurCpuForPmeOrEwald =
532             (EVDW_PME(inputrec.vdwtype)
533              || (EEL_PME_EWALD(inputrec.coulombtype) && !useGpuForPme && numPmeRanksPerSimulation <= 0));
534
535     return gpusWereDetected && usingOurCpuForPmeOrEwald;
536 }
537
538 bool decideWhetherToUseGpuForUpdate(const bool                     isDomainDecomposition,
539                                     const bool                     useUpdateGroups,
540                                     const PmeRunMode               pmeRunMode,
541                                     const bool                     havePmeOnlyRank,
542                                     const bool                     useGpuForNonbonded,
543                                     const TaskTarget               updateTarget,
544                                     const bool                     gpusWereDetected,
545                                     const t_inputrec&              inputrec,
546                                     const gmx_mtop_t&              mtop,
547                                     const bool                     useEssentialDynamics,
548                                     const bool                     doOrientationRestraints,
549                                     const bool                     haveFrozenAtoms,
550                                     const bool                     doRerun,
551                                     const DevelopmentFeatureFlags& devFlags,
552                                     const gmx::MDLogger&           mdlog)
553 {
554
555     // '-update cpu' overrides the environment variable, '-update auto' does not
556     if (updateTarget == TaskTarget::Cpu
557         || (updateTarget == TaskTarget::Auto && !devFlags.forceGpuUpdateDefault))
558     {
559         return false;
560     }
561
562     const bool hasAnyConstraints = gmx_mtop_interaction_count(mtop, IF_CONSTRAINT) > 0;
563     const bool pmeUsesCpu = (pmeRunMode == PmeRunMode::CPU || pmeRunMode == PmeRunMode::Mixed);
564
565     std::string errorMessage;
566
567     if (isDomainDecomposition)
568     {
569         if (hasAnyConstraints && !useUpdateGroups)
570         {
571             errorMessage +=
572                     "Domain decomposition is only supported with constraints when update "
573                     "groups "
574                     "are used. This means constraining all bonds is not supported, except for "
575                     "small molecules, and box sizes close to half the pair-list cutoff are not "
576                     "supported.\n ";
577         }
578     }
579
580     if (havePmeOnlyRank)
581     {
582         if (pmeUsesCpu)
583         {
584             errorMessage += "With separate PME rank(s), PME must run fully on the GPU.\n";
585         }
586     }
587
588     if (inputrec.useMts)
589     {
590         errorMessage += "Multiple time stepping is not supported.\n";
591     }
592
593     if (inputrec.eConstrAlg == ConstraintAlgorithm::Shake && hasAnyConstraints
594         && gmx_mtop_ftype_count(mtop, F_CONSTR) > 0)
595     {
596         errorMessage += "SHAKE constraints are not supported.\n";
597     }
598     // Using the GPU-version of update if:
599     // 1. PME is on the GPU (there should be a copy of coordinates on GPU for PME spread) or inactive, or
600     // 2. Non-bonded interactions are on the GPU.
601     if ((pmeRunMode == PmeRunMode::CPU || pmeRunMode == PmeRunMode::None) && !useGpuForNonbonded)
602     {
603         errorMessage +=
604                 "Either PME or short-ranged non-bonded interaction tasks must run on the GPU.\n";
605     }
606     if (!gpusWereDetected)
607     {
608         errorMessage += "Compatible GPUs must have been found.\n";
609     }
610     if (!GMX_GPU_CUDA)
611     {
612         errorMessage += "Only a CUDA build is supported.\n";
613     }
614     if (inputrec.eI != IntegrationAlgorithm::MD)
615     {
616         errorMessage += "Only the md integrator is supported.\n";
617     }
618     if (inputrec.etc == TemperatureCoupling::NoseHoover)
619     {
620         errorMessage += "Nose-Hoover temperature coupling is not supported.\n";
621     }
622     if (!(inputrec.epc == PressureCoupling::No || inputrec.epc == PressureCoupling::ParrinelloRahman
623           || inputrec.epc == PressureCoupling::Berendsen || inputrec.epc == PressureCoupling::CRescale))
624     {
625         errorMessage +=
626                 "Only Parrinello-Rahman, Berendsen, and C-rescale pressure coupling are "
627                 "supported.\n";
628     }
629     if (EEL_PME_EWALD(inputrec.coulombtype) && inputrec.epsilon_surface != 0)
630     {
631         // The graph is needed, but not supported
632         errorMessage += "Ewald surface correction is not supported.\n";
633     }
634     if (gmx_mtop_interaction_count(mtop, IF_VSITE) > 0)
635     {
636         errorMessage += "Virtual sites are not supported.\n";
637     }
638     if (useEssentialDynamics)
639     {
640         errorMessage += "Essential dynamics is not supported.\n";
641     }
642     if (inputrec.bPull && pull_have_constraint(*inputrec.pull))
643     {
644         errorMessage += "Constraints pulling is not supported.\n";
645     }
646     if (doOrientationRestraints)
647     {
648         // The graph is needed, but not supported
649         errorMessage += "Orientation restraints are not supported.\n";
650     }
651     if (inputrec.efep != FreeEnergyPerturbationType::No
652         && (haveFepPerturbedMasses(mtop) || havePerturbedConstraints(mtop)))
653     {
654         errorMessage += "Free energy perturbation for mass and constraints are not supported.\n";
655     }
656     const auto particleTypes = gmx_mtop_particletype_count(mtop);
657     if (particleTypes[ParticleType::Shell] > 0)
658     {
659         errorMessage += "Shells are not supported.\n";
660     }
661     if (inputrec.eSwapCoords != SwapType::No)
662     {
663         errorMessage += "Swapping the coordinates is not supported.\n";
664     }
665     if (doRerun)
666     {
667         errorMessage += "Re-run is not supported.\n";
668     }
669
670     // TODO: F_CONSTRNC is only unsupported, because isNumCoupledConstraintsSupported()
671     // does not support it, the actual CUDA LINCS code does support it
672     if (gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0)
673     {
674         errorMessage += "Non-connecting constraints are not supported\n";
675     }
676     if (!UpdateConstrainGpu::isNumCoupledConstraintsSupported(mtop))
677     {
678         errorMessage +=
679                 "The number of coupled constraints is higher than supported in the GPU LINCS "
680                 "code.\n";
681     }
682     if (hasAnyConstraints && !UpdateConstrainGpu::areConstraintsSupported())
683     {
684         errorMessage += "Chosen GPU implementation does not support constraints.\n";
685     }
686     if (haveFrozenAtoms)
687     {
688         // There is a known bug with frozen atoms and GPU update, see Issue #3920.
689         errorMessage += "Frozen atoms not supported.\n";
690     }
691
692     if (!errorMessage.empty())
693     {
694         if (updateTarget == TaskTarget::Auto && devFlags.forceGpuUpdateDefault)
695         {
696             GMX_LOG(mdlog.warning)
697                     .asParagraph()
698                     .appendText(
699                             "Update task on the GPU was required, by the "
700                             "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable, but the following "
701                             "condition(s) were not satisfied:");
702             GMX_LOG(mdlog.warning).asParagraph().appendText(errorMessage.c_str());
703             GMX_LOG(mdlog.warning).asParagraph().appendText("Will use CPU version of update.");
704         }
705         else if (updateTarget == TaskTarget::Gpu)
706         {
707             std::string prefix = gmx::formatString(
708                     "Update task on the GPU was required,\n"
709                     "but the following condition(s) were not satisfied:\n");
710             GMX_THROW(InconsistentInputError((prefix + errorMessage).c_str()));
711         }
712         return false;
713     }
714
715     return (updateTarget == TaskTarget::Gpu
716             || (updateTarget == TaskTarget::Auto && devFlags.forceGpuUpdateDefault));
717 }
718
719 bool decideWhetherToUseGpuForHalo(const DevelopmentFeatureFlags& devFlags,
720                                   bool                           havePPDomainDecomposition,
721                                   bool                           useGpuForNonbonded,
722                                   bool                           useModularSimulator,
723                                   bool                           doRerun,
724                                   bool                           haveEnergyMinimization)
725 {
726     return havePPDomainDecomposition && devFlags.enableGpuHaloExchange && useGpuForNonbonded
727            && !useModularSimulator && !doRerun && !haveEnergyMinimization;
728 }
729
730 } // namespace gmx