Split lines with many copyright years
[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_cuda.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 != GMX_GPU_NONE
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 == GMX_GPU_CUDA
92         "CUDA_VISIBLE_DEVICES"
93 #    elif GMX_GPU == 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 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) || !useGpuForNonbonded || !pme_gpu_supports_build(nullptr)
163         || !pme_gpu_supports_hardware(hardwareInfo, nullptr)
164         || !pme_gpu_supports_input(inputrec, mtop, nullptr))
165     {
166         // PME can't run on a GPU. If the user required that, we issue
167         // an error later.
168         return false;
169     }
170
171     // We now know that PME on GPUs might make sense, if we have any.
172
173     if (!userGpuTaskAssignment.empty())
174     {
175         // Follow the user's choice of GPU task assignment, if we
176         // can. Checking that their IDs are for compatible GPUs comes
177         // later.
178
179         // Specifying -gputasks requires specifying everything.
180         if (pmeTarget == TaskTarget::Auto || numRanksPerSimulation < 1)
181         {
182             GMX_THROW(InconsistentInputError(
183                     formatString(g_specifyEverythingFormatString, "all of -nb, -pme, and -ntmpi")));
184         }
185
186         // PME on GPUs is only supported in a single case
187         if (pmeTarget == TaskTarget::Gpu)
188         {
189             if (((numRanksPerSimulation > 1) && (numPmeRanksPerSimulation == 0))
190                 || (numPmeRanksPerSimulation > 1))
191             {
192                 GMX_THROW(InconsistentInputError(
193                         "When you run mdrun -pme gpu -gputasks, you must supply a PME-enabled .tpr "
194                         "file and use a single PME rank."));
195             }
196             return true;
197         }
198
199         // pmeTarget == TaskTarget::Auto
200         return numRanksPerSimulation == 1;
201     }
202
203     // Because this is thread-MPI, we already know about the GPUs that
204     // all potential ranks can use, and can use that in a global
205     // decision that will later be consistent.
206
207     if (pmeTarget == TaskTarget::Gpu)
208     {
209         if (((numRanksPerSimulation > 1) && (numPmeRanksPerSimulation == 0))
210             || (numPmeRanksPerSimulation > 1))
211         {
212             GMX_THROW(NotImplementedError(
213                     "PME tasks were required to run on GPUs, but that is not implemented with "
214                     "more than one PME rank. Use a single rank simulation, or a separate PME rank, "
215                     "or permit PME tasks to be assigned to the CPU."));
216         }
217         return true;
218     }
219
220     if (numRanksPerSimulation == 1)
221     {
222         // PME can run well on a GPU shared with NB, and we permit
223         // mdrun to default to try that.
224         return !gpuIdsToUse.empty();
225     }
226
227     if (numRanksPerSimulation < 1)
228     {
229         // Full automated mode for thread-MPI (the default). PME can
230         // run well on a GPU shared with NB, and we permit mdrun to
231         // default to it if there is only one GPU available.
232         return (gpuIdsToUse.size() == 1);
233     }
234
235     // Not enough support for PME on GPUs for anything else
236     return false;
237 }
238
239 bool decideWhetherToUseGpusForNonbonded(const TaskTarget          nonbondedTarget,
240                                         const std::vector<int>&   userGpuTaskAssignment,
241                                         const EmulateGpuNonbonded emulateGpuNonbonded,
242                                         const bool                buildSupportsNonbondedOnGpu,
243                                         const bool                nonbondedOnGpuIsUseful,
244                                         const bool                gpusWereDetected)
245 {
246     if (nonbondedTarget == TaskTarget::Cpu)
247     {
248         if (!userGpuTaskAssignment.empty())
249         {
250             GMX_THROW(InconsistentInputError(
251                     "A GPU task assignment was specified, but nonbonded interactions were "
252                     "assigned to the CPU. Make no more than one of these choices."));
253         }
254
255         return false;
256     }
257
258     if (!buildSupportsNonbondedOnGpu && nonbondedTarget == TaskTarget::Gpu)
259     {
260         GMX_THROW(InconsistentInputError(
261                 "Nonbonded interactions on the GPU were requested with -nb gpu, "
262                 "but the GROMACS binary has been built without GPU support. "
263                 "Either run without selecting GPU options, or recompile GROMACS "
264                 "with GPU support enabled"));
265     }
266
267     // TODO refactor all these TaskTarget::Gpu checks into one place?
268     // e.g. use a subfunction that handles only the cases where
269     // TaskTargets are not Cpu?
270     if (emulateGpuNonbonded == EmulateGpuNonbonded::Yes)
271     {
272         if (nonbondedTarget == TaskTarget::Gpu)
273         {
274             GMX_THROW(InconsistentInputError(
275                     "Nonbonded interactions on the GPU were required, which is inconsistent "
276                     "with choosing emulation. Make no more than one of these choices."));
277         }
278         if (!userGpuTaskAssignment.empty())
279         {
280             GMX_THROW(
281                     InconsistentInputError("GPU ID usage was specified, as was GPU emulation. Make "
282                                            "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(
306                     formatString(g_specifyEverythingFormatString, "-nb and -ntmpi")));
307         }
308
309         return true;
310     }
311
312     if (nonbondedTarget == TaskTarget::Gpu)
313     {
314         // We still don't know whether it is an error if no GPUs are found
315         // because we don't know the duty of this rank, yet. For example,
316         // a node with only PME ranks and -pme cpu is OK if there are not
317         // GPUs.
318         return true;
319     }
320
321     // If we get here, then the user permitted GPUs, which we should
322     // use for nonbonded interactions.
323     return gpusWereDetected;
324 }
325
326 bool decideWhetherToUseGpusForPme(const bool              useGpuForNonbonded,
327                                   const TaskTarget        pmeTarget,
328                                   const std::vector<int>& userGpuTaskAssignment,
329                                   const gmx_hw_info_t&    hardwareInfo,
330                                   const t_inputrec&       inputrec,
331                                   const gmx_mtop_t&       mtop,
332                                   const int               numRanksPerSimulation,
333                                   const int               numPmeRanksPerSimulation,
334                                   const bool              gpusWereDetected)
335 {
336     if (pmeTarget == TaskTarget::Cpu)
337     {
338         return false;
339     }
340
341     if (!useGpuForNonbonded)
342     {
343         if (pmeTarget == TaskTarget::Gpu)
344         {
345             GMX_THROW(NotImplementedError(
346                     "PME on GPUs is only supported when nonbonded interactions run on GPUs also."));
347         }
348         return false;
349     }
350
351     std::string message;
352     if (!pme_gpu_supports_build(&message))
353     {
354         if (pmeTarget == TaskTarget::Gpu)
355         {
356             GMX_THROW(NotImplementedError("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("Cannot compute PME interactions on a GPU, because " + message));
365         }
366         return false;
367     }
368     if (!pme_gpu_supports_input(inputrec, mtop, &message))
369     {
370         if (pmeTarget == TaskTarget::Gpu)
371         {
372             GMX_THROW(NotImplementedError("Cannot compute PME interactions on a GPU, because " + message));
373         }
374         return false;
375     }
376
377     if (pmeTarget == TaskTarget::Cpu)
378     {
379         if (!userGpuTaskAssignment.empty())
380         {
381             GMX_THROW(InconsistentInputError(
382                     "A GPU task assignment was specified, but PME interactions were "
383                     "assigned to the CPU. Make no more than one of these choices."));
384         }
385
386         return false;
387     }
388
389     if (!userGpuTaskAssignment.empty())
390     {
391         // Specifying -gputasks requires specifying everything.
392         if (pmeTarget == TaskTarget::Auto)
393         {
394             GMX_THROW(InconsistentInputError(formatString(
395                     g_specifyEverythingFormatString, "all of -nb, -pme, and -ntmpi"))); // TODO ntmpi?
396         }
397
398         return true;
399     }
400
401     // We still don't know whether it is an error if no GPUs are found
402     // because we don't know the duty of this rank, yet. For example,
403     // a node with only PME ranks and -pme cpu is OK if there are not
404     // GPUs.
405
406     if (pmeTarget == TaskTarget::Gpu)
407     {
408         if (((numRanksPerSimulation > 1) && (numPmeRanksPerSimulation == 0))
409             || (numPmeRanksPerSimulation > 1))
410         {
411             GMX_THROW(NotImplementedError(
412                     "PME tasks were required to run on GPUs, but that is not implemented with "
413                     "more than one PME rank. Use a single rank simulation, or a separate PME rank, "
414                     "or permit PME tasks to be assigned to the CPU."));
415         }
416         return true;
417     }
418
419     // If we get here, then the user permitted GPUs.
420     if (numRanksPerSimulation == 1)
421     {
422         // PME can run well on a single GPU shared with NB when there
423         // is one rank, so we permit mdrun to try that if we have
424         // detected GPUs.
425         return gpusWereDetected;
426     }
427
428     // Not enough support for PME on GPUs for anything else
429     return false;
430 }
431
432 bool decideWhetherToUseGpusForBonded(const bool       useGpuForNonbonded,
433                                      const bool       useGpuForPme,
434                                      const TaskTarget bondedTarget,
435                                      const bool       canUseGpuForBonded,
436                                      const bool       usingLJPme,
437                                      const bool       usingElecPmeOrEwald,
438                                      const int        numPmeRanksPerSimulation,
439                                      const bool       gpusWereDetected)
440 {
441     if (bondedTarget == TaskTarget::Cpu)
442     {
443         return false;
444     }
445
446     if (!canUseGpuForBonded)
447     {
448         if (bondedTarget == TaskTarget::Gpu)
449         {
450             GMX_THROW(InconsistentInputError(
451                     "Bonded interactions on the GPU were required, but not supported for these "
452                     "simulation settings. Change your settings, or do not require using GPUs."));
453         }
454
455         return false;
456     }
457
458     if (!useGpuForNonbonded)
459     {
460         if (bondedTarget == TaskTarget::Gpu)
461         {
462             GMX_THROW(InconsistentInputError(
463                     "Bonded interactions on the GPU were required, but this requires that "
464                     "short-ranged non-bonded interactions are also run on the GPU. Change "
465                     "your settings, or do not require using GPUs."));
466         }
467
468         return false;
469     }
470
471     // TODO If the bonded kernels do not get fused, then performance
472     // overheads might suggest alternative choices here.
473
474     if (bondedTarget == TaskTarget::Gpu)
475     {
476         // We still don't know whether it is an error if no GPUs are
477         // found.
478         return true;
479     }
480
481     // If we get here, then the user permitted GPUs, which we should
482     // use for bonded interactions if any were detected and the CPU
483     // is busy, for which we currently only check PME or Ewald.
484     // (It would be better to dynamically assign bondeds based on timings)
485     // Note that here we assume that the auto setting of PME ranks will not
486     // choose seperate PME ranks when nonBonded are assigned to the GPU.
487     bool usingOurCpuForPmeOrEwald =
488             (usingLJPme || (usingElecPmeOrEwald && !useGpuForPme && numPmeRanksPerSimulation <= 0));
489
490     return gpusWereDetected && usingOurCpuForPmeOrEwald;
491 }
492
493 bool decideWhetherToUseGpuForUpdate(const bool           forceGpuUpdateDefault,
494                                     const bool           isDomainDecomposition,
495                                     const bool           useUpdateGroups,
496                                     const PmeRunMode     pmeRunMode,
497                                     const bool           havePmeOnlyRank,
498                                     const bool           useGpuForNonbonded,
499                                     const TaskTarget     updateTarget,
500                                     const bool           gpusWereDetected,
501                                     const t_inputrec&    inputrec,
502                                     const gmx_mtop_t&    mtop,
503                                     const bool           useEssentialDynamics,
504                                     const bool           doOrientationRestraints,
505                                     const bool           useReplicaExchange,
506                                     const bool           doRerun,
507                                     const gmx::MDLogger& mdlog)
508 {
509
510     // '-update cpu' overrides the environment variable, '-update auto' does not
511     if (updateTarget == TaskTarget::Cpu || (updateTarget == TaskTarget::Auto && !forceGpuUpdateDefault))
512     {
513         return false;
514     }
515
516     const bool hasAnyConstraints = gmx_mtop_interaction_count(mtop, IF_CONSTRAINT) > 0;
517
518     std::string errorMessage;
519
520     if (isDomainDecomposition)
521     {
522         if (!forceGpuUpdateDefault)
523         {
524             errorMessage += "Domain decomposition is not supported.\n ";
525         }
526         else if (hasAnyConstraints && !useUpdateGroups)
527         {
528             errorMessage +=
529                     "Domain decomposition is only supported with constraints when update groups "
530                     "are used. This means constraining all bonds is not supported, except for "
531                     "small molecules, and box sizes close to half the pair-list cutoff are not "
532                     "supported.\n ";
533         }
534     }
535     if (inputrec.eConstrAlg == econtSHAKE && hasAnyConstraints && gmx_mtop_ftype_count(mtop, F_CONSTR) > 0)
536     {
537         errorMessage += "SHAKE constraints are not supported.\n";
538     }
539     // Using the GPU-version of update if:
540     // 1. PME is on the GPU (there should be a copy of coordinates on GPU for PME spread), or
541     // 2. Non-bonded interactions are on the GPU.
542     if (pmeRunMode == PmeRunMode::CPU && !useGpuForNonbonded)
543     {
544         errorMessage +=
545                 "Either PME or short-ranged non-bonded interaction tasks must run on the GPU.\n";
546     }
547     // Since only direct GPU communications are supported with GPU update, PME should be fully offloaded in DD and PME only cases.
548     if (pmeRunMode != PmeRunMode::GPU && (isDomainDecomposition || havePmeOnlyRank))
549     {
550         errorMessage += "PME should run on GPU.\n";
551     }
552     if (!gpusWereDetected)
553     {
554         errorMessage += "Compatible GPUs must have been found.\n";
555     }
556     if (GMX_GPU != GMX_GPU_CUDA)
557     {
558         errorMessage += "Only a CUDA build is supported.\n";
559     }
560     if (inputrec.eI != eiMD)
561     {
562         errorMessage += "Only the md integrator is supported.\n";
563     }
564     if (inputrec.etc == etcNOSEHOOVER)
565     {
566         errorMessage += "Nose-Hoover temperature coupling is not supported.\n";
567     }
568     if (!(inputrec.epc == epcNO || inputrec.epc == epcPARRINELLORAHMAN || inputrec.epc == epcBERENDSEN))
569     {
570         errorMessage += "Only Parrinello-Rahman and Berendsen pressure coupling are supported.\n";
571     }
572     if (EEL_PME_EWALD(inputrec.coulombtype) && inputrec.epsilon_surface != 0)
573     {
574         // The graph is needed, but not supported
575         errorMessage += "Ewald surface correction is not supported.\n";
576     }
577     if (gmx_mtop_interaction_count(mtop, IF_VSITE) > 0)
578     {
579         errorMessage += "Virtual sites are not supported.\n";
580     }
581     if (useEssentialDynamics)
582     {
583         errorMessage += "Essential dynamics is not supported.\n";
584     }
585     if (inputrec.bPull && pull_have_constraint(inputrec.pull))
586     {
587         errorMessage += "Constraints pulling is not supported.\n";
588     }
589     if (doOrientationRestraints)
590     {
591         // The graph is needed, but not supported
592         errorMessage += "Orientation restraints are not supported.\n";
593     }
594     if (inputrec.efep != efepNO)
595     {
596         // Actually all free-energy options except for mass and constraint perturbation are supported
597         errorMessage += "Free energy perturbations are not supported.\n";
598     }
599     if (useReplicaExchange)
600     {
601         errorMessage += "Replica exchange simulations are not supported.\n";
602     }
603     if (inputrec.eSwapCoords != eswapNO)
604     {
605         errorMessage += "Swapping the coordinates is not supported.\n";
606     }
607     if (doRerun)
608     {
609         errorMessage += "Re-run is not supported.\n";
610     }
611
612     // TODO: F_CONSTRNC is only unsupported, because isNumCoupledConstraintsSupported()
613     // does not support it, the actual CUDA LINCS code does support it
614     if (gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0)
615     {
616         errorMessage += "Non-connecting constraints are not supported";
617     }
618     if (!UpdateConstrainCuda::isNumCoupledConstraintsSupported(mtop))
619     {
620         errorMessage +=
621                 "The number of coupled constraints is higher than supported in the CUDA LINCS "
622                 "code.\n";
623     }
624
625     if (!errorMessage.empty())
626     {
627         if (updateTarget != TaskTarget::Gpu && forceGpuUpdateDefault)
628         {
629             GMX_LOG(mdlog.warning)
630                     .asParagraph()
631                     .appendText(
632                             "Update task on the GPU was required, by the "
633                             "GMX_FORCE_UPDATE_DEFAULT_GPU environment variable, but the following "
634                             "condition(s) were not satisfied:");
635             GMX_LOG(mdlog.warning).asParagraph().appendText(errorMessage.c_str());
636             GMX_LOG(mdlog.warning).asParagraph().appendText("Will use CPU version of update.");
637         }
638         else if (updateTarget == TaskTarget::Gpu)
639         {
640             std::string prefix = gmx::formatString(
641                     "Update task on the GPU was required,\n"
642                     "but the following condition(s) were not satisfied:\n");
643             GMX_THROW(InconsistentInputError((prefix + errorMessage).c_str()));
644         }
645         return false;
646     }
647
648     if (isDomainDecomposition)
649     {
650         return forceGpuUpdateDefault;
651     }
652     else
653     {
654         return (updateTarget == TaskTarget::Gpu || forceGpuUpdateDefault);
655     }
656 }
657
658 } // namespace gmx