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