Disable PME Mixed mode with FEP
[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/gpubonded.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* 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
100         // As with OpenCL, there are no portable way to do it.
101         // Intel reference: https://github.com/intel/llvm/blob/sycl/sycl/doc/EnvironmentVariables.md
102         // While SYCL_DEVICE_FILTER is a better option, as of 2021.1-beta10 it is not yet supported.
103         "SYCL_DEVICE_ALLOWLIST"
104 #    else
105 #        error "Unreachable branch"
106 #    endif
107         " environment variable in your bash profile or job "
108         "script may be more convenient."
109 #endif
110         ;
111
112 } // namespace
113
114 bool decideWhetherToUseGpusForNonbondedWithThreadMpi(const TaskTarget        nonbondedTarget,
115                                                      const int               numDevicesToUse,
116                                                      const std::vector<int>& userGpuTaskAssignment,
117                                                      const EmulateGpuNonbonded emulateGpuNonbonded,
118                                                      const bool buildSupportsNonbondedOnGpu,
119                                                      const bool nonbondedOnGpuIsUseful,
120                                                      const int  numRanksPerSimulation)
121 {
122     // First, exclude all cases where we can't run NB on GPUs.
123     if (nonbondedTarget == TaskTarget::Cpu || emulateGpuNonbonded == EmulateGpuNonbonded::Yes
124         || !nonbondedOnGpuIsUseful || !buildSupportsNonbondedOnGpu)
125     {
126         // If the user required NB on GPUs, we issue an error later.
127         return false;
128     }
129
130     // We now know that NB on GPUs makes sense, if we have any.
131
132     if (!userGpuTaskAssignment.empty())
133     {
134         // Specifying -gputasks requires specifying everything.
135         if (nonbondedTarget == TaskTarget::Auto || numRanksPerSimulation < 1)
136         {
137             GMX_THROW(InconsistentInputError(
138                     formatString(g_specifyEverythingFormatString, "-nb and -ntmpi")));
139         }
140         return true;
141     }
142
143     if (nonbondedTarget == TaskTarget::Gpu)
144     {
145         return true;
146     }
147
148     // Because this is thread-MPI, we already know about the GPUs that
149     // all potential ranks can use, and can use that in a global
150     // decision that will later be consistent.
151     // If we get here, then the user permitted or required GPUs.
152     return (numDevicesToUse > 0);
153 }
154
155 bool decideWhetherToUseGpusForPmeWithThreadMpi(const bool              useGpuForNonbonded,
156                                                const TaskTarget        pmeTarget,
157                                                const TaskTarget        pmeFftTarget,
158                                                const int               numDevicesToUse,
159                                                const std::vector<int>& userGpuTaskAssignment,
160                                                const gmx_hw_info_t&    hardwareInfo,
161                                                const t_inputrec&       inputrec,
162                                                const int               numRanksPerSimulation,
163                                                const int               numPmeRanksPerSimulation)
164 {
165     // First, exclude all cases where we can't run PME on GPUs.
166     if ((pmeTarget == TaskTarget::Cpu) || !useGpuForNonbonded || !pme_gpu_supports_build(nullptr)
167         || !pme_gpu_supports_hardware(hardwareInfo, nullptr) || !pme_gpu_supports_input(inputrec, nullptr))
168     {
169         // PME can't run on a GPU. If the user required that, we issue an error later.
170         return false;
171     }
172
173     if (pmeFftTarget == TaskTarget::Cpu && !pme_gpu_mixed_mode_supports_input(inputrec, nullptr))
174     {
175         /* User requested PME FFT on CPU, but the current system is not compatible with Mixed mode,
176          * so we don't use GPUs at all.
177          * If the user had -pme gpu, we issue an error later. */
178         return false;
179     }
180
181     // We now know that PME on GPUs might make sense, if we have any.
182
183     if (!userGpuTaskAssignment.empty())
184     {
185         // Follow the user's choice of GPU task assignment, if we
186         // can. Checking that their IDs are for compatible GPUs comes
187         // later.
188
189         // Specifying -gputasks requires specifying everything.
190         if (pmeTarget == TaskTarget::Auto || numRanksPerSimulation < 1)
191         {
192             GMX_THROW(InconsistentInputError(
193                     formatString(g_specifyEverythingFormatString, "all of -nb, -pme, and -ntmpi")));
194         }
195
196         // PME on GPUs is only supported in a single case
197         if (pmeTarget == TaskTarget::Gpu)
198         {
199             if (((numRanksPerSimulation > 1) && (numPmeRanksPerSimulation == 0))
200                 || (numPmeRanksPerSimulation > 1))
201             {
202                 GMX_THROW(InconsistentInputError(
203                         "When you run mdrun -pme gpu -gputasks, you must supply a PME-enabled .tpr "
204                         "file and use a single PME rank."));
205             }
206             return true;
207         }
208
209         // pmeTarget == TaskTarget::Auto
210         return numRanksPerSimulation == 1;
211     }
212
213     // Because this is thread-MPI, we already know about the GPUs that
214     // all potential ranks can use, and can use that in a global
215     // decision that will later be consistent.
216
217     if (pmeTarget == TaskTarget::Gpu)
218     {
219         if (((numRanksPerSimulation > 1) && (numPmeRanksPerSimulation == 0))
220             || (numPmeRanksPerSimulation > 1))
221         {
222             GMX_THROW(NotImplementedError(
223                     "PME tasks were required to run on GPUs, but that is not implemented with "
224                     "more than one PME rank. Use a single rank simulation, or a separate PME rank, "
225                     "or permit PME tasks to be assigned to the CPU."));
226         }
227         return true;
228     }
229
230     if (numRanksPerSimulation == 1)
231     {
232         // PME can run well on a GPU shared with NB, and we permit
233         // mdrun to default to try that.
234         return numDevicesToUse > 0;
235     }
236
237     if (numRanksPerSimulation < 1)
238     {
239         // Full automated mode for thread-MPI (the default). PME can
240         // run well on a GPU shared with NB, and we permit mdrun to
241         // default to it if there is only one GPU available.
242         return (numDevicesToUse == 1);
243     }
244
245     // Not enough support for PME on GPUs for anything else
246     return false;
247 }
248
249 bool decideWhetherToUseGpusForNonbonded(const TaskTarget          nonbondedTarget,
250                                         const std::vector<int>&   userGpuTaskAssignment,
251                                         const EmulateGpuNonbonded emulateGpuNonbonded,
252                                         const bool                buildSupportsNonbondedOnGpu,
253                                         const bool                nonbondedOnGpuIsUseful,
254                                         const bool                gpusWereDetected)
255 {
256     if (nonbondedTarget == TaskTarget::Cpu)
257     {
258         if (!userGpuTaskAssignment.empty())
259         {
260             GMX_THROW(InconsistentInputError(
261                     "A GPU task assignment was specified, but nonbonded interactions were "
262                     "assigned to the CPU. Make no more than one of these choices."));
263         }
264
265         return false;
266     }
267
268     if (!buildSupportsNonbondedOnGpu && nonbondedTarget == TaskTarget::Gpu)
269     {
270         GMX_THROW(InconsistentInputError(
271                 "Nonbonded interactions on the GPU were requested with -nb gpu, "
272                 "but the GROMACS binary has been built without GPU support. "
273                 "Either run without selecting GPU options, or recompile GROMACS "
274                 "with GPU support enabled"));
275     }
276
277     // TODO refactor all these TaskTarget::Gpu checks into one place?
278     // e.g. use a subfunction that handles only the cases where
279     // TaskTargets are not Cpu?
280     if (emulateGpuNonbonded == EmulateGpuNonbonded::Yes)
281     {
282         if (nonbondedTarget == TaskTarget::Gpu)
283         {
284             GMX_THROW(InconsistentInputError(
285                     "Nonbonded interactions on the GPU were required, which is inconsistent "
286                     "with choosing emulation. Make no more than one of these choices."));
287         }
288         if (!userGpuTaskAssignment.empty())
289         {
290             GMX_THROW(
291                     InconsistentInputError("GPU ID usage was specified, as was GPU emulation. Make "
292                                            "no more than one of these choices."));
293         }
294
295         return false;
296     }
297
298     if (!nonbondedOnGpuIsUseful)
299     {
300         if (nonbondedTarget == TaskTarget::Gpu)
301         {
302             GMX_THROW(InconsistentInputError(
303                     "Nonbonded interactions on the GPU were required, but not supported for these "
304                     "simulation settings. Change your settings, or do not require using GPUs."));
305         }
306
307         return false;
308     }
309
310     if (!userGpuTaskAssignment.empty())
311     {
312         // Specifying -gputasks requires specifying everything.
313         if (nonbondedTarget == TaskTarget::Auto)
314         {
315             GMX_THROW(InconsistentInputError(
316                     formatString(g_specifyEverythingFormatString, "-nb and -ntmpi")));
317         }
318
319         return true;
320     }
321
322     if (nonbondedTarget == TaskTarget::Gpu)
323     {
324         // We still don't know whether it is an error if no GPUs are found
325         // because we don't know the duty of this rank, yet. For example,
326         // a node with only PME ranks and -pme cpu is OK if there are not
327         // GPUs.
328         return true;
329     }
330
331     // If we get here, then the user permitted GPUs, which we should
332     // use for nonbonded interactions.
333     return buildSupportsNonbondedOnGpu && gpusWereDetected;
334 }
335
336 bool decideWhetherToUseGpusForPme(const bool              useGpuForNonbonded,
337                                   const TaskTarget        pmeTarget,
338                                   const TaskTarget        pmeFftTarget,
339                                   const std::vector<int>& userGpuTaskAssignment,
340                                   const gmx_hw_info_t&    hardwareInfo,
341                                   const t_inputrec&       inputrec,
342                                   const int               numRanksPerSimulation,
343                                   const int               numPmeRanksPerSimulation,
344                                   const bool              gpusWereDetected)
345 {
346     if (pmeTarget == TaskTarget::Cpu)
347     {
348         return false;
349     }
350
351     if (!useGpuForNonbonded)
352     {
353         if (pmeTarget == TaskTarget::Gpu)
354         {
355             GMX_THROW(NotImplementedError(
356                     "PME on GPUs is only supported when nonbonded interactions run on GPUs also."));
357         }
358         return false;
359     }
360
361     std::string message;
362     if (!pme_gpu_supports_build(&message))
363     {
364         if (pmeTarget == TaskTarget::Gpu)
365         {
366             GMX_THROW(NotImplementedError("Cannot compute PME interactions on a GPU, because " + message));
367         }
368         return false;
369     }
370     if (!pme_gpu_supports_hardware(hardwareInfo, &message))
371     {
372         if (pmeTarget == TaskTarget::Gpu)
373         {
374             GMX_THROW(NotImplementedError("Cannot compute PME interactions on a GPU, because " + message));
375         }
376         return false;
377     }
378     if (!pme_gpu_supports_input(inputrec, &message))
379     {
380         if (pmeTarget == TaskTarget::Gpu)
381         {
382             GMX_THROW(NotImplementedError("Cannot compute PME interactions on a GPU, because " + message));
383         }
384         return false;
385     }
386     if (pmeFftTarget == TaskTarget::Cpu && !pme_gpu_mixed_mode_supports_input(inputrec, &message))
387     {
388         /* User requested PME FFT on CPU, but the current system is not compatible with Mixed mode,
389          * so we don't use GPUs at all. */
390         if (pmeTarget == TaskTarget::Gpu)
391         {
392             GMX_THROW(NotImplementedError("Cannot compute PME interactions in Mixed mode, because " + message));
393         }
394         return false;
395     }
396
397     if (pmeTarget == TaskTarget::Cpu)
398     {
399         if (!userGpuTaskAssignment.empty())
400         {
401             GMX_THROW(InconsistentInputError(
402                     "A GPU task assignment was specified, but PME interactions were "
403                     "assigned to the CPU. Make no more than one of these choices."));
404         }
405
406         return false;
407     }
408
409     if (!userGpuTaskAssignment.empty())
410     {
411         // Specifying -gputasks requires specifying everything.
412         if (pmeTarget == TaskTarget::Auto)
413         {
414             GMX_THROW(InconsistentInputError(formatString(
415                     g_specifyEverythingFormatString, "all of -nb, -pme, and -ntmpi"))); // TODO ntmpi?
416         }
417
418         return true;
419     }
420
421     // We still don't know whether it is an error if no GPUs are found
422     // because we don't know the duty of this rank, yet. For example,
423     // a node with only PME ranks and -pme cpu is OK if there are not
424     // GPUs.
425
426     if (pmeTarget == TaskTarget::Gpu)
427     {
428         if (((numRanksPerSimulation > 1) && (numPmeRanksPerSimulation == 0))
429             || (numPmeRanksPerSimulation > 1))
430         {
431             GMX_THROW(NotImplementedError(
432                     "PME tasks were required to run on GPUs, but that is not implemented with "
433                     "more than one PME rank. Use a single rank simulation, or a separate PME rank, "
434                     "or permit PME tasks to be assigned to the CPU."));
435         }
436         return true;
437     }
438
439     // If we get here, then the user permitted GPUs.
440     if (numRanksPerSimulation == 1)
441     {
442         // PME can run well on a single GPU shared with NB when there
443         // is one rank, so we permit mdrun to try that if we have
444         // detected GPUs.
445         return gpusWereDetected;
446     }
447
448     // Not enough support for PME on GPUs for anything else
449     return false;
450 }
451
452
453 PmeRunMode determinePmeRunMode(const bool useGpuForPme, const TaskTarget& pmeFftTarget, const t_inputrec& inputrec)
454 {
455     if (!EEL_PME(inputrec.coulombtype))
456     {
457         return PmeRunMode::None;
458     }
459
460     if (useGpuForPme)
461     {
462         if (pmeFftTarget == TaskTarget::Cpu)
463         {
464             return PmeRunMode::Mixed;
465         }
466         else
467         {
468             return PmeRunMode::GPU;
469         }
470     }
471     else
472     {
473         if (pmeFftTarget == TaskTarget::Gpu)
474         {
475             gmx_fatal(FARGS,
476                       "Assigning FFTs to GPU requires PME to be assigned to GPU as well. With PME "
477                       "on CPU you should not be using -pmefft.");
478         }
479         return PmeRunMode::CPU;
480     }
481 }
482
483 bool decideWhetherToUseGpusForBonded(bool              useGpuForNonbonded,
484                                      bool              useGpuForPme,
485                                      TaskTarget        bondedTarget,
486                                      const t_inputrec& inputrec,
487                                      const gmx_mtop_t& mtop,
488                                      int               numPmeRanksPerSimulation,
489                                      bool              gpusWereDetected)
490 {
491     if (bondedTarget == TaskTarget::Cpu)
492     {
493         return false;
494     }
495
496     std::string errorMessage;
497
498     if (!buildSupportsGpuBondeds(&errorMessage))
499     {
500         if (bondedTarget == TaskTarget::Gpu)
501         {
502             GMX_THROW(InconsistentInputError(errorMessage.c_str()));
503         }
504
505         return false;
506     }
507
508     if (!inputSupportsGpuBondeds(inputrec, mtop, &errorMessage))
509     {
510         if (bondedTarget == TaskTarget::Gpu)
511         {
512             GMX_THROW(InconsistentInputError(errorMessage.c_str()));
513         }
514
515         return false;
516     }
517
518     if (!useGpuForNonbonded)
519     {
520         if (bondedTarget == TaskTarget::Gpu)
521         {
522             GMX_THROW(InconsistentInputError(
523                     "Bonded interactions on the GPU were required, but this requires that "
524                     "short-ranged non-bonded interactions are also run on the GPU. Change "
525                     "your settings, or do not require using GPUs."));
526         }
527
528         return false;
529     }
530
531     // TODO If the bonded kernels do not get fused, then performance
532     // overheads might suggest alternative choices here.
533
534     if (bondedTarget == TaskTarget::Gpu)
535     {
536         // We still don't know whether it is an error if no GPUs are
537         // found.
538         return true;
539     }
540
541     // If we get here, then the user permitted GPUs, which we should
542     // use for bonded interactions if any were detected and the CPU
543     // is busy, for which we currently only check PME or Ewald.
544     // (It would be better to dynamically assign bondeds based on timings)
545     // Note that here we assume that the auto setting of PME ranks will not
546     // choose seperate PME ranks when nonBonded are assigned to the GPU.
547     bool usingOurCpuForPmeOrEwald =
548             (EVDW_PME(inputrec.vdwtype)
549              || (EEL_PME_EWALD(inputrec.coulombtype) && !useGpuForPme && numPmeRanksPerSimulation <= 0));
550
551     return gpusWereDetected && usingOurCpuForPmeOrEwald;
552 }
553
554 bool decideWhetherToUseGpuForUpdate(const bool                     isDomainDecomposition,
555                                     const bool                     useUpdateGroups,
556                                     const PmeRunMode               pmeRunMode,
557                                     const bool                     havePmeOnlyRank,
558                                     const bool                     useGpuForNonbonded,
559                                     const TaskTarget               updateTarget,
560                                     const bool                     gpusWereDetected,
561                                     const t_inputrec&              inputrec,
562                                     const gmx_mtop_t&              mtop,
563                                     const bool                     useEssentialDynamics,
564                                     const bool                     doOrientationRestraints,
565                                     const bool                     useReplicaExchange,
566                                     const bool                     haveFrozenAtoms,
567                                     const bool                     doRerun,
568                                     const DevelopmentFeatureFlags& devFlags,
569                                     const gmx::MDLogger&           mdlog)
570 {
571
572     // '-update cpu' overrides the environment variable, '-update auto' does not
573     if (updateTarget == TaskTarget::Cpu
574         || (updateTarget == TaskTarget::Auto && !devFlags.forceGpuUpdateDefault))
575     {
576         return false;
577     }
578
579     const bool hasAnyConstraints = gmx_mtop_interaction_count(mtop, IF_CONSTRAINT) > 0;
580     const bool pmeUsesCpu = (pmeRunMode == PmeRunMode::CPU || pmeRunMode == PmeRunMode::Mixed);
581
582     std::string errorMessage;
583
584     if (isDomainDecomposition)
585     {
586         if (hasAnyConstraints && !useUpdateGroups)
587         {
588             errorMessage +=
589                     "Domain decomposition is only supported with constraints when update "
590                     "groups "
591                     "are used. This means constraining all bonds is not supported, except for "
592                     "small molecules, and box sizes close to half the pair-list cutoff are not "
593                     "supported.\n ";
594         }
595     }
596
597     if (havePmeOnlyRank)
598     {
599         if (pmeUsesCpu)
600         {
601             errorMessage += "With separate PME rank(s), PME must run fully on the GPU.\n";
602         }
603     }
604
605     if (inputrec.useMts)
606     {
607         errorMessage += "Multiple time stepping is not supported.\n";
608     }
609
610     if (inputrec.eConstrAlg == econtSHAKE && hasAnyConstraints && gmx_mtop_ftype_count(mtop, F_CONSTR) > 0)
611     {
612         errorMessage += "SHAKE constraints are not supported.\n";
613     }
614     // Using the GPU-version of update if:
615     // 1. PME is on the GPU (there should be a copy of coordinates on GPU for PME spread) or inactive, or
616     // 2. Non-bonded interactions are on the GPU.
617     if ((pmeRunMode == PmeRunMode::CPU || pmeRunMode == PmeRunMode::None) && !useGpuForNonbonded)
618     {
619         errorMessage +=
620                 "Either PME or short-ranged non-bonded interaction tasks must run on the GPU.\n";
621     }
622
623     if (!gpusWereDetected)
624     {
625         errorMessage += "Compatible GPUs must have been found.\n";
626     }
627     if (!GMX_GPU_CUDA)
628     {
629         errorMessage += "Only a CUDA build is supported.\n";
630     }
631     if (inputrec.eI != eiMD)
632     {
633         errorMessage += "Only the md integrator is supported.\n";
634     }
635     if (inputrec.etc == etcNOSEHOOVER)
636     {
637         errorMessage += "Nose-Hoover temperature coupling is not supported.\n";
638     }
639     if (!(inputrec.epc == epcNO || inputrec.epc == epcPARRINELLORAHMAN
640           || inputrec.epc == epcBERENDSEN || inputrec.epc == epcCRESCALE))
641     {
642         errorMessage +=
643                 "Only Parrinello-Rahman, Berendsen, and C-rescale pressure coupling are "
644                 "supported.\n";
645     }
646     if (EEL_PME_EWALD(inputrec.coulombtype) && inputrec.epsilon_surface != 0)
647     {
648         // The graph is needed, but not supported
649         errorMessage += "Ewald surface correction is not supported.\n";
650     }
651     if (gmx_mtop_interaction_count(mtop, IF_VSITE) > 0)
652     {
653         errorMessage += "Virtual sites are not supported.\n";
654     }
655     if (useEssentialDynamics)
656     {
657         errorMessage += "Essential dynamics is not supported.\n";
658     }
659     if (inputrec.bPull && pull_have_constraint(*inputrec.pull))
660     {
661         errorMessage += "Constraints pulling is not supported.\n";
662     }
663     if (doOrientationRestraints)
664     {
665         // The graph is needed, but not supported
666         errorMessage += "Orientation restraints are not supported.\n";
667     }
668     if (inputrec.efep != efepNO && (haveFepPerturbedMasses(mtop) || havePerturbedConstraints(mtop)))
669     {
670         errorMessage += "Free energy perturbation for mass and constraints are not supported.\n";
671     }
672     const auto particleTypes = gmx_mtop_particletype_count(mtop);
673     if (particleTypes[eptShell] > 0)
674     {
675         errorMessage += "Shells are not supported.\n";
676     }
677     if (useReplicaExchange)
678     {
679         errorMessage += "Replica exchange simulations are not supported.\n";
680     }
681     if (inputrec.eSwapCoords != eswapNO)
682     {
683         errorMessage += "Swapping the coordinates is not supported.\n";
684     }
685     if (doRerun)
686     {
687         errorMessage += "Re-run is not supported.\n";
688     }
689
690     // TODO: F_CONSTRNC is only unsupported, because isNumCoupledConstraintsSupported()
691     // does not support it, the actual CUDA LINCS code does support it
692     if (gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0)
693     {
694         errorMessage += "Non-connecting constraints are not supported\n";
695     }
696     if (!UpdateConstrainGpu::isNumCoupledConstraintsSupported(mtop))
697     {
698         errorMessage +=
699                 "The number of coupled constraints is higher than supported in the GPU LINCS "
700                 "code.\n";
701     }
702     if (haveFrozenAtoms)
703     {
704         // There is a known bug with frozen atoms and GPU update, see Issue #3920.
705         errorMessage += "Frozen atoms not supported.\n";
706     }
707
708     if (!errorMessage.empty())
709     {
710         if (updateTarget == TaskTarget::Auto && devFlags.forceGpuUpdateDefault)
711         {
712             GMX_LOG(mdlog.warning)
713                     .asParagraph()
714                     .appendText(
715                             "Update task on the GPU was required, by the "
716                             "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable, but the following "
717                             "condition(s) were not satisfied:");
718             GMX_LOG(mdlog.warning).asParagraph().appendText(errorMessage.c_str());
719             GMX_LOG(mdlog.warning).asParagraph().appendText("Will use CPU version of update.");
720         }
721         else if (updateTarget == TaskTarget::Gpu)
722         {
723             std::string prefix = gmx::formatString(
724                     "Update task on the GPU was required,\n"
725                     "but the following condition(s) were not satisfied:\n");
726             GMX_THROW(InconsistentInputError((prefix + errorMessage).c_str()));
727         }
728         return false;
729     }
730
731     return (updateTarget == TaskTarget::Gpu
732             || (updateTarget == TaskTarget::Auto && devFlags.forceGpuUpdateDefault));
733 }
734
735 bool decideWhetherToUseGpuForHalo(const DevelopmentFeatureFlags& devFlags,
736                                   bool                           havePPDomainDecomposition,
737                                   bool                           useGpuForNonbonded,
738                                   bool                           useModularSimulator,
739                                   bool                           doRerun,
740                                   bool                           haveEnergyMinimization)
741 {
742     return havePPDomainDecomposition && devFlags.enableGpuHaloExchange && useGpuForNonbonded
743            && !useModularSimulator && !doRerun && !haveEnergyMinimization;
744 }
745
746 } // namespace gmx