Enable separate PME GPU rank
[alexxy/gromacs.git] / src / gromacs / taskassignment / resourcedivision.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2015,2016,2017, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  * \brief Defines utility functionality for dividing resources and
37  * checking for consistency and usefulness.
38  *
39  * \author Mark Abraham <mark.j.abraham@gmail.com>
40  * \ingroup module_taskassignment
41  */
42
43 #include "gmxpre.h"
44
45 #include "resourcedivision.h"
46
47 #include "config.h"
48
49 #include <stdlib.h>
50 #include <string.h>
51
52 #include <algorithm>
53
54 #include "gromacs/ewald/pme.h"
55 #include "gromacs/hardware/cpuinfo.h"
56 #include "gromacs/hardware/detecthardware.h"
57 #include "gromacs/hardware/hardwaretopology.h"
58 #include "gromacs/hardware/hw_info.h"
59 #include "gromacs/mdlib/gmx_omp_nthreads.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/topology/mtop_util.h"
64 #include "gromacs/topology/topology.h"
65 #include "gromacs/utility/baseversion.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/gmxassert.h"
68 #include "gromacs/utility/logger.h"
69 #include "gromacs/utility/stringutil.h"
70
71
72 /* DISCLAIMER: All the atom count and thread numbers below are heuristic.
73  * The real switching points will depend on the system simulation,
74  * the algorithms used and the hardware it's running on, as well as if there
75  * are other jobs running on the same machine. We try to take into account
76  * factors that have a large influence, such as recent Intel CPUs being
77  * much better at wide multi-threading. The remaining factors should
78  * (hopefully) have a small influence, such that the performance just before
79  * and after a switch point doesn't change too much.
80  */
81
82 //! Constant used to help minimize preprocessed code
83 static const bool bHasOmpSupport = GMX_OPENMP;
84
85 /*! \brief The minimum number of atoms per thread-MPI thread when GPUs
86  * are present. With fewer atoms than this, the number of thread-MPI
87  * ranks will get lowered.
88  */
89 static const int min_atoms_per_mpi_thread =  90;
90 /*! \brief The minimum number of atoms per GPU with thread-MPI
91  * active. With fewer atoms than this, the number of thread-MPI ranks
92  * will get lowered.
93  */
94 static const int min_atoms_per_gpu        = 900;
95
96 /**@{*/
97 /*! \brief Constants for implementing default divisions of threads */
98
99 /* TODO choose nthreads_omp based on hardware topology
100    when we have a hardware topology detection library */
101 /* First we consider the case of no MPI (1 MPI rank).
102  * In general, when running up to 8 threads, OpenMP should be faster.
103  * Note: on AMD Bulldozer we should avoid running OpenMP over two dies.
104  * On Intel>=Nehalem running OpenMP on a single CPU is always faster,
105  * even on two CPUs it's usually faster (but with many OpenMP threads
106  * it could be faster not to use HT, currently we always use HT).
107  * On Nehalem/Westmere we want to avoid running 16 threads over
108  * two CPUs with HT, so we need a limit<16; thus we use 12.
109  * A reasonable limit for Intel Sandy and Ivy bridge,
110  * not knowing the topology, is 16 threads.
111  * Below we check for Intel and AVX, which for now includes
112  * Sandy/Ivy Bridge, Has/Broadwell. By checking for AVX instead of
113  * model numbers we ensure also future Intel CPUs are covered.
114  */
115 const int nthreads_omp_faster_default   =  8;
116 const int nthreads_omp_faster_Nehalem   = 12;
117 const int nthreads_omp_faster_Intel_AVX = 16;
118 const int nthreads_omp_faster_AMD_Ryzen = 16;
119 /* For CPU only runs the fastest options are usually MPI or OpenMP only.
120  * With one GPU, using MPI only is almost never optimal, so we need to
121  * compare running pure OpenMP with combined MPI+OpenMP. This means higher
122  * OpenMP threads counts can still be ok. Multiplying the numbers above
123  * by a factor of 2 seems to be a good estimate.
124  */
125 const int nthreads_omp_faster_gpu_fac   =  2;
126
127 /* This is the case with MPI (2 or more MPI PP ranks).
128  * By default we will terminate with a fatal error when more than 8
129  * OpenMP thread are (indirectly) requested, since using less threads
130  * nearly always results in better performance.
131  * With thread-mpi and multiple GPUs or one GPU and too many threads
132  * we first try 6 OpenMP threads and then less until the number of MPI ranks
133  * is divisible by the number of GPUs.
134  */
135 #if GMX_OPENMP && GMX_MPI
136 const int nthreads_omp_mpi_ok_max              =  8;
137 const int nthreads_omp_mpi_ok_min_cpu          =  1;
138 #endif
139 const int nthreads_omp_mpi_ok_min_gpu          =  2;
140 const int nthreads_omp_mpi_target_max          =  6;
141
142 /**@}*/
143
144 /*! \brief Returns the maximum OpenMP thread count for which using a single MPI rank
145  * should be faster than using multiple ranks with the same total thread count.
146  */
147 static int nthreads_omp_faster(const gmx::CpuInfo &cpuInfo, gmx_bool bUseGPU)
148 {
149     int nth;
150
151     if (cpuInfo.vendor() == gmx::CpuInfo::Vendor::Intel &&
152         cpuInfo.feature(gmx::CpuInfo::Feature::X86_Avx))
153     {
154         nth = nthreads_omp_faster_Intel_AVX;
155     }
156     else if (gmx::cpuIsX86Nehalem(cpuInfo))
157     {
158         // Intel Nehalem
159         nth = nthreads_omp_faster_Nehalem;
160     }
161     else if (cpuInfo.vendor() == gmx::CpuInfo::Vendor::Amd && cpuInfo.family() >= 23)
162     {
163         // AMD Ryzen
164         nth = nthreads_omp_faster_AMD_Ryzen;
165     }
166     else
167     {
168         nth = nthreads_omp_faster_default;
169     }
170
171     if (bUseGPU)
172     {
173         nth *= nthreads_omp_faster_gpu_fac;
174     }
175
176     nth = std::min(nth, GMX_OPENMP_MAX_THREADS);
177
178     return nth;
179 }
180
181 /*! \brief Returns that maximum OpenMP thread count that passes the efficiency check */
182 gmx_unused static int nthreads_omp_efficient_max(int gmx_unused       nrank,
183                                                  const gmx::CpuInfo  &cpuInfo,
184                                                  gmx_bool             bUseGPU)
185 {
186 #if GMX_OPENMP && GMX_MPI
187     if (nrank > 1)
188     {
189         return nthreads_omp_mpi_ok_max;
190     }
191     else
192 #endif
193     {
194         return nthreads_omp_faster(cpuInfo, bUseGPU);
195     }
196 }
197
198 /*! \brief Return the number of thread-MPI ranks to use.
199  * This is chosen such that we can always obey our own efficiency checks.
200  */
201 gmx_unused static int get_tmpi_omp_thread_division(const gmx_hw_info_t *hwinfo,
202                                                    const gmx_hw_opt_t  &hw_opt,
203                                                    int                  nthreads_tot,
204                                                    int                  ngpu)
205 {
206     int                 nrank;
207     const gmx::CpuInfo &cpuInfo = *hwinfo->cpuInfo;
208
209     GMX_RELEASE_ASSERT(nthreads_tot > 0, "There must be at least one thread per rank");
210
211     /* There are no separate PME nodes here, as we ensured in
212      * check_and_update_hw_opt that nthreads_tmpi>0 with PME nodes
213      * and a conditional ensures we would not have ended up here.
214      * Note that separate PME nodes might be switched on later.
215      */
216     if (ngpu > 0)
217     {
218         nrank = ngpu;
219
220         /* When the user sets nthreads_omp, we can end up oversubscribing CPU cores
221          * if we simply start as many ranks as GPUs. To avoid this, we start as few
222          * tMPI ranks as necessary to avoid oversubscription and instead leave GPUs idle.
223          * If the user does not set the number of OpenMP threads, nthreads_omp==0 and
224          * this code has no effect.
225          */
226         GMX_RELEASE_ASSERT(hw_opt.nthreads_omp >= 0, "nthreads_omp is negative, but previous checks should have prevented this");
227         while (nrank*hw_opt.nthreads_omp > hwinfo->nthreads_hw_avail && nrank > 1)
228         {
229             nrank--;
230         }
231
232         if (nthreads_tot < nrank)
233         {
234             /* #thread < #gpu is very unlikely, but if so: waste gpu(s) */
235             nrank = nthreads_tot;
236         }
237         else if (nthreads_tot > nthreads_omp_faster(cpuInfo, ngpu > 0) ||
238                  (ngpu > 1 && nthreads_tot/ngpu > nthreads_omp_mpi_target_max))
239         {
240             /* The high OpenMP thread count will likely result in sub-optimal
241              * performance. Increase the rank count to reduce the thread count
242              * per rank. This will lead to GPU sharing by MPI ranks/threads.
243              */
244             int nshare;
245
246             /* Increase the rank count as long as have we more than 6 OpenMP
247              * threads per rank or the number of hardware threads is not
248              * divisible by the rank count. Don't go below 2 OpenMP threads.
249              */
250             nshare = 1;
251             do
252             {
253                 nshare++;
254                 nrank = ngpu*nshare;
255             }
256             while (nthreads_tot/nrank > nthreads_omp_mpi_target_max ||
257                    (nthreads_tot/(ngpu*(nshare + 1)) >= nthreads_omp_mpi_ok_min_gpu && nthreads_tot % nrank != 0));
258         }
259     }
260     else if (hw_opt.nthreads_omp > 0)
261     {
262         /* Here we could oversubscribe, when we do, we issue a warning later */
263         nrank = std::max(1, nthreads_tot/hw_opt.nthreads_omp);
264     }
265     else
266     {
267         if (nthreads_tot <= nthreads_omp_faster(cpuInfo, ngpu > 0))
268         {
269             /* Use pure OpenMP parallelization */
270             nrank = 1;
271         }
272         else
273         {
274             /* Don't use OpenMP parallelization */
275             nrank = nthreads_tot;
276         }
277     }
278
279     return nrank;
280 }
281
282 //! Return whether hyper threading is enabled.
283 static bool
284 gmxSmtIsEnabled(const gmx::HardwareTopology &hwTop)
285 {
286     return (hwTop.supportLevel() >= gmx::HardwareTopology::SupportLevel::Basic && hwTop.machine().sockets[0].cores[0].hwThreads.size() > 1);
287 }
288
289 namespace
290 {
291
292 //! Handles checks for algorithms that must use a single rank.
293 class SingleRankChecker
294 {
295     public:
296         //! Constructor
297         SingleRankChecker() : value_(false), reasons_() {}
298         /*! \brief Call this function for each possible condition
299             under which a single rank is required, along with a string
300             describing the constraint when it is applied. */
301         void applyConstraint(bool condition, const char *description)
302         {
303             if (condition)
304             {
305                 value_ = true;
306                 reasons_.push_back(gmx::formatString("%s only supports a single rank.", description));
307             }
308         }
309         //! After applying any conditions, is a single rank required?
310         bool mustUseOneRank() const
311         {
312             return value_;
313         }
314         /*! \brief Return a formatted string to use when writing a
315             message when a single rank is required, (or empty if no
316             constraint exists.) */
317         std::string getMessage() const
318         {
319             return formatAndJoin(reasons_, "\n", gmx::IdentityFormatter());
320         }
321     private:
322         bool                     value_;
323         std::vector<std::string> reasons_;
324 };
325
326 } // namespace
327
328 /* Get the number of MPI ranks to use for thread-MPI based on how many
329  * were requested, which algorithms we're using,
330  * and how many particles there are.
331  * At the point we have already called check_and_update_hw_opt.
332  * Thus all options should be internally consistent and consistent
333  * with the hardware, except that ntmpi could be larger than #GPU.
334  */
335 int get_nthreads_mpi(const gmx_hw_info_t    *hwinfo,
336                      gmx_hw_opt_t           *hw_opt,
337                      const std::vector<int> &gpuIdsToUse,
338                      bool                    nonbondedOnGpu,
339                      bool                    pmeOnGpu,
340                      const t_inputrec       *inputrec,
341                      const gmx_mtop_t       *mtop,
342                      const gmx::MDLogger    &mdlog,
343                      bool                    doMembed)
344 {
345     int                          nthreads_hw, nthreads_tot_max, nrank, ngpu;
346     int                          min_atoms_per_mpi_rank;
347
348     const gmx::CpuInfo          &cpuInfo = *hwinfo->cpuInfo;
349     const gmx::HardwareTopology &hwTop   = *hwinfo->hardwareTopology;
350
351     if (pmeOnGpu)
352     {
353         GMX_RELEASE_ASSERT((EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)) &&
354                            pme_gpu_supports_input(inputrec, nullptr),
355                            "PME can't be on GPUs unless we are using PME");
356
357         // PME on GPUs supports a single PME rank with PP running on the same or few other ranks.
358         // For now, let's treat separate PME GPU rank as opt-in.
359         if (hw_opt->nthreads_tmpi < 1)
360         {
361             return 1;
362         }
363     }
364
365     {
366         /* Check if an algorithm does not support parallel simulation.  */
367         // TODO This might work better if e.g. implemented algorithms
368         // had to define a function that returns such requirements,
369         // and a description string.
370         SingleRankChecker checker;
371         checker.applyConstraint(inputrec->eI == eiLBFGS, "L-BFGS minimization");
372         checker.applyConstraint(inputrec->coulombtype == eelEWALD, "Plain Ewald electrostatics");
373         checker.applyConstraint(doMembed, "Membrane embedding");
374         bool useOrientationRestraints = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
375         checker.applyConstraint(useOrientationRestraints, "Orientation restraints");
376         if (checker.mustUseOneRank())
377         {
378             std::string message = checker.getMessage();
379             if (hw_opt->nthreads_tmpi > 1)
380             {
381                 gmx_fatal(FARGS, "%s However, you asked for more than 1 thread-MPI rank, so mdrun cannot continue. Choose a single rank, or a different algorithm.", message.c_str());
382             }
383             GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted("%s Choosing to use only a single thread-MPI rank.", message.c_str());
384             return 1;
385         }
386     }
387
388     if (hw_opt->nthreads_tmpi > 0)
389     {
390         /* Trivial, return the user's choice right away */
391         return hw_opt->nthreads_tmpi;
392     }
393
394     // Now implement automatic selection of number of thread-MPI ranks
395     nthreads_hw = hwinfo->nthreads_hw_avail;
396
397     if (nthreads_hw <= 0)
398     {
399         /* This should normally not happen, but if it does, we handle it */
400         gmx_fatal(FARGS, "The number of available hardware threads can not be detected, please specify the number of MPI ranks and the number of OpenMP threads (if supported) manually with options -ntmpi and -ntomp, respectively");
401     }
402
403     /* How many total (#tMPI*#OpenMP) threads can we start? */
404     if (hw_opt->nthreads_tot > 0)
405     {
406         nthreads_tot_max = hw_opt->nthreads_tot;
407     }
408     else
409     {
410         nthreads_tot_max = nthreads_hw;
411     }
412
413     /* nonbondedOnGpu might be false e.g. because this simulation uses
414      * the group scheme, or is a rerun with energy groups. */
415     ngpu = (nonbondedOnGpu ? static_cast<int>(gpuIdsToUse.size()) : 0);
416
417     if (inputrec->cutoff_scheme == ecutsGROUP)
418     {
419         /* We checked this before, but it doesn't hurt to do it once more */
420         GMX_RELEASE_ASSERT(hw_opt->nthreads_omp == 1, "The group scheme only supports one OpenMP thread per rank");
421     }
422
423     nrank =
424         get_tmpi_omp_thread_division(hwinfo, *hw_opt, nthreads_tot_max, ngpu);
425
426     if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
427     {
428         /* Dims/steps are divided over the nodes iso splitting the atoms.
429          * With NM we can't have more ranks than #atoms*#dim. With TPI it's
430          * unlikely we have fewer atoms than ranks, and if so, communication
431          * would become a bottleneck, so we set the limit to 1 atom/rank.
432          */
433         min_atoms_per_mpi_rank = 1;
434     }
435     else
436     {
437         if (ngpu >= 1)
438         {
439             min_atoms_per_mpi_rank = min_atoms_per_gpu;
440         }
441         else
442         {
443             min_atoms_per_mpi_rank = min_atoms_per_mpi_thread;
444         }
445     }
446
447     if (mtop->natoms/nrank < min_atoms_per_mpi_rank)
448     {
449         int nrank_new;
450
451         /* the rank number was chosen automatically, but there are too few
452            atoms per rank, so we need to reduce the rank count */
453         nrank_new = std::max(1, mtop->natoms/min_atoms_per_mpi_rank);
454
455         /* Avoid partial use of Hyper-Threading */
456         if (gmxSmtIsEnabled(hwTop) &&
457             nrank_new > nthreads_hw/2 && nrank_new < nthreads_hw)
458         {
459             nrank_new = nthreads_hw/2;
460         }
461
462         /* If the user specified the total thread count, ensure this is
463          * divisible by the number of ranks.
464          * It is quite likely that we have too many total threads compared
465          * to the size of the system, but if the user asked for this many
466          * threads we should respect that.
467          */
468         while (hw_opt->nthreads_tot > 0 &&
469                hw_opt->nthreads_tot % nrank_new != 0)
470         {
471             nrank_new--;
472         }
473
474         /* Avoid large prime numbers in the rank count */
475         if (nrank_new >= 6)
476         {
477             /* Use only 6,8,10 with additional factors of 2 */
478             int fac;
479
480             fac = 2;
481             while (3*fac*2 <= nrank_new)
482             {
483                 fac *= 2;
484             }
485
486             nrank_new = (nrank_new/fac)*fac;
487         }
488         else
489         {
490             /* Avoid 5, since small system won't fit 5 domains along
491              * a dimension. This might lead to waisting some cores, but this
492              * will have a small impact in this regime of very small systems.
493              */
494             if (nrank_new == 5)
495             {
496                 nrank_new = 4;
497             }
498         }
499
500         nrank = nrank_new;
501
502         /* We reduced the number of tMPI ranks, which means we might violate
503          * our own efficiency checks if we simply use all hardware threads.
504          */
505         if (bHasOmpSupport && hw_opt->nthreads_omp <= 0 && hw_opt->nthreads_tot <= 0)
506         {
507             /* The user set neither the total nor the OpenMP thread count,
508              * we should use all hardware threads, unless we will violate
509              * our own efficiency limitation on the thread count.
510              */
511             int  nt_omp_max;
512
513             nt_omp_max = nthreads_omp_efficient_max(nrank, cpuInfo, ngpu >= 1);
514
515             if (nrank*nt_omp_max < hwinfo->nthreads_hw_avail)
516             {
517                 /* Limit the number of OpenMP threads to start */
518                 hw_opt->nthreads_omp = nt_omp_max;
519             }
520         }
521
522         fprintf(stderr, "\n");
523         fprintf(stderr, "NOTE: Parallelization is limited by the small number of atoms,\n");
524         fprintf(stderr, "      only starting %d thread-MPI ranks.\n", nrank);
525         fprintf(stderr, "      You can use the -nt and/or -ntmpi option to optimize the number of threads.\n\n");
526     }
527
528     return nrank;
529 }
530
531
532 void check_resource_division_efficiency(const gmx_hw_info_t *hwinfo,
533                                         int                  numTotalThreads,
534                                         bool                 willUsePhysicalGpu,
535                                         gmx_bool             bNtOmpOptionSet,
536                                         t_commrec           *cr,
537                                         const gmx::MDLogger &mdlog)
538 {
539 #if GMX_OPENMP && GMX_MPI
540     int         nth_omp_min, nth_omp_max;
541     char        buf[1000];
542 #if GMX_THREAD_MPI
543     const char *mpi_option = " (option -ntmpi)";
544 #else
545     const char *mpi_option = "";
546 #endif
547
548     /* This function should be called after thread-MPI (when configured) and
549      * OpenMP have been initialized. Check that here.
550      */
551 #if GMX_THREAD_MPI
552     GMX_RELEASE_ASSERT(nthreads_omp_faster_default >= nthreads_omp_mpi_ok_max, "Inconsistent OpenMP thread count default values");
553 #endif
554     GMX_RELEASE_ASSERT(gmx_omp_nthreads_get(emntDefault) >= 1, "Must have at least one OpenMP thread");
555
556     nth_omp_min = gmx_omp_nthreads_get(emntDefault);
557     nth_omp_max = gmx_omp_nthreads_get(emntDefault);
558
559     bool anyRankIsUsingGpus = willUsePhysicalGpu;
560     /* Thread-MPI seems to have a bug with reduce on 1 node, so use a cond. */
561     if (cr->nnodes + cr->npmenodes > 1)
562     {
563         int count[3], count_max[3];
564
565         count[0] = -nth_omp_min;
566         count[1] =  nth_omp_max;
567         count[2] =  willUsePhysicalGpu;
568
569         MPI_Allreduce(count, count_max, 3, MPI_INT, MPI_MAX, cr->mpi_comm_mysim);
570
571         /* In case of an inhomogeneous run setup we use the maximum counts */
572         nth_omp_min        = -count_max[0];
573         nth_omp_max        =  count_max[1];
574         anyRankIsUsingGpus = count_max[2] > 0;
575     }
576
577     int nthreads_omp_mpi_ok_min;
578
579     if (!anyRankIsUsingGpus)
580     {
581         nthreads_omp_mpi_ok_min = nthreads_omp_mpi_ok_min_cpu;
582     }
583     else
584     {
585         /* With GPUs we set the minimum number of OpenMP threads to 2 to catch
586          * cases where the user specifies #ranks == #cores.
587          */
588         nthreads_omp_mpi_ok_min = nthreads_omp_mpi_ok_min_gpu;
589     }
590
591     if (DOMAINDECOMP(cr) && cr->nnodes > 1)
592     {
593         if (nth_omp_max < nthreads_omp_mpi_ok_min ||
594             nth_omp_max > nthreads_omp_mpi_ok_max)
595         {
596             /* Note that we print target_max here, not ok_max */
597             sprintf(buf, "Your choice of number of MPI ranks and amount of resources results in using %d OpenMP threads per rank, which is most likely inefficient. The optimum is usually between %d and %d threads per rank.",
598                     nth_omp_max,
599                     nthreads_omp_mpi_ok_min,
600                     nthreads_omp_mpi_target_max);
601
602             if (bNtOmpOptionSet)
603             {
604                 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted("NOTE: %s", buf);
605             }
606             else
607             {
608                 /* This fatal error, and the one below, is nasty, but it's
609                  * probably the only way to ensure that all users don't waste
610                  * a lot of resources, since many users don't read logs/stderr.
611                  */
612                 gmx_fatal(FARGS, "%s If you want to run with this setup, specify the -ntomp option. But we suggest to change the number of MPI ranks%s.", buf, mpi_option);
613             }
614         }
615     }
616     else
617     {
618         const gmx::CpuInfo &cpuInfo = *hwinfo->cpuInfo;
619
620         /* No domain decomposition (or only one domain) */
621         if (nth_omp_max > nthreads_omp_faster(cpuInfo, anyRankIsUsingGpus))
622         {
623             /* To arrive here, the user/system set #ranks and/or #OMPthreads */
624             gmx_bool bEnvSet;
625             char     buf2[256];
626
627             bEnvSet = (getenv("OMP_NUM_THREADS") != nullptr);
628
629             if (bNtOmpOptionSet || bEnvSet)
630             {
631                 sprintf(buf2, "You requested %d OpenMP threads", nth_omp_max);
632             }
633             else
634             {
635                 sprintf(buf2, "Your choice of %d MPI rank%s and the use of %d total threads %sleads to the use of %d OpenMP threads",
636                         cr->nnodes + cr->npmenodes,
637                         cr->nnodes + cr->npmenodes == 1 ? "" : "s",
638                         numTotalThreads > 0 ? numTotalThreads : hwinfo->nthreads_hw_avail,
639                         hwinfo->nphysicalnode > 1 ? "on a node " : "",
640                         nth_omp_max);
641             }
642             sprintf(buf, "%s, whereas we expect the optimum to be with more MPI ranks with %d to %d OpenMP threads.",
643                     buf2, nthreads_omp_mpi_ok_min, nthreads_omp_mpi_target_max);
644
645             /* We can not quit with a fatal error when OMP_NUM_THREADS is set
646              * with different values per rank or node, since in that case
647              * the user can not set -ntomp to override the error.
648              */
649             if (bNtOmpOptionSet || (bEnvSet && nth_omp_min != nth_omp_max))
650             {
651                 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted("NOTE: %s", buf);
652             }
653             else
654             {
655                 gmx_fatal(FARGS, "%s If you want to run with this many OpenMP threads, specify the -ntomp option. But we suggest to increase the number of MPI ranks%s.", buf, mpi_option);
656             }
657         }
658     }
659 #else /* GMX_OPENMP && GMX_MPI */
660       /* No OpenMP and/or MPI: it doesn't make much sense to check */
661     GMX_UNUSED_VALUE(bNtOmpOptionSet);
662     GMX_UNUSED_VALUE(numTotalThreads);
663     GMX_UNUSED_VALUE(willUsePhysicalGpu);
664     GMX_UNUSED_VALUE(cr);
665     /* Check if we have more than 1 physical core, if detected,
666      * or more than 1 hardware thread if physical cores were not detected.
667      */
668     if (!GMX_OPENMP && !GMX_MPI && hwinfo->hardwareTopology->numberOfCores() > 1)
669     {
670         GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: GROMACS was compiled without OpenMP and (thread-)MPI support, can only use a single CPU core");
671     }
672 #endif /* GMX_OPENMP && GMX_MPI */
673 }
674
675
676 //! Dump a \c hw_opt to \c fp.
677 static void print_hw_opt(FILE *fp, const gmx_hw_opt_t *hw_opt)
678 {
679     fprintf(fp, "hw_opt: nt %d ntmpi %d ntomp %d ntomp_pme %d gpu_id '%s' gputasks '%s'\n",
680             hw_opt->nthreads_tot,
681             hw_opt->nthreads_tmpi,
682             hw_opt->nthreads_omp,
683             hw_opt->nthreads_omp_pme,
684             hw_opt->gpuIdsAvailable.c_str(),
685             hw_opt->userGpuTaskAssignment.c_str());
686 }
687
688 void check_and_update_hw_opt_1(gmx_hw_opt_t    *hw_opt,
689                                const t_commrec *cr,
690                                int              nPmeRanks)
691 {
692     /* Currently hw_opt only contains default settings or settings supplied
693      * by the user on the command line.
694      */
695     if (hw_opt->nthreads_omp < 0)
696     {
697         gmx_fatal(FARGS, "The number of OpenMP threads supplied on the command line is %d, which is negative and not allowed", hw_opt->nthreads_omp);
698     }
699
700     /* Check for OpenMP settings stored in environment variables, which can
701      * potentially be different on different MPI ranks.
702      */
703     gmx_omp_nthreads_read_env(&hw_opt->nthreads_omp, SIMMASTER(cr));
704
705     /* Check restrictions on the user supplied options before modifying them.
706      * TODO: Put the user values in a const struct and preserve them.
707      */
708 #if !GMX_THREAD_MPI
709     if (hw_opt->nthreads_tot > 0)
710     {
711         gmx_fatal(FARGS, "Setting the total number of threads is only supported with thread-MPI and GROMACS was compiled without thread-MPI");
712     }
713     if (hw_opt->nthreads_tmpi > 0)
714     {
715         gmx_fatal(FARGS, "Setting the number of thread-MPI ranks is only supported with thread-MPI and GROMACS was compiled without thread-MPI");
716     }
717 #endif
718
719     if (bHasOmpSupport)
720     {
721         /* Check restrictions on PME thread related options set by the user */
722
723         if (hw_opt->nthreads_omp_pme > 0 && hw_opt->nthreads_omp <= 0)
724         {
725             gmx_fatal(FARGS, "You need to specify -ntomp in addition to -ntomp_pme");
726         }
727
728         if (hw_opt->nthreads_omp_pme >= 1 &&
729             hw_opt->nthreads_omp_pme != hw_opt->nthreads_omp &&
730             nPmeRanks <= 0)
731         {
732             /* This can result in a fatal error on many MPI ranks,
733              * but since the thread count can differ per rank,
734              * we can't easily avoid this.
735              */
736             gmx_fatal(FARGS, "You need to explicitly specify the number of PME ranks (-npme) when using different number of OpenMP threads for PP and PME ranks");
737         }
738     }
739     else
740     {
741         /* GROMACS was configured without OpenMP support */
742
743         if (hw_opt->nthreads_omp > 1 || hw_opt->nthreads_omp_pme > 1)
744         {
745             gmx_fatal(FARGS, "More than 1 OpenMP thread requested, but GROMACS was compiled without OpenMP support");
746         }
747         hw_opt->nthreads_omp     = 1;
748         hw_opt->nthreads_omp_pme = 1;
749     }
750
751     if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp_pme <= 0)
752     {
753         /* We have the same number of OpenMP threads for PP and PME ranks,
754          * thus we can perform several consistency checks.
755          */
756         if (hw_opt->nthreads_tmpi > 0 &&
757             hw_opt->nthreads_omp > 0 &&
758             hw_opt->nthreads_tot != hw_opt->nthreads_tmpi*hw_opt->nthreads_omp)
759         {
760             gmx_fatal(FARGS, "The total number of threads requested (%d) does not match the thread-MPI ranks (%d) times the OpenMP threads (%d) requested",
761                       hw_opt->nthreads_tot, hw_opt->nthreads_tmpi, hw_opt->nthreads_omp);
762         }
763
764         if (hw_opt->nthreads_tmpi > 0 &&
765             hw_opt->nthreads_tot % hw_opt->nthreads_tmpi != 0)
766         {
767             gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of thread-MPI ranks requested (%d)",
768                       hw_opt->nthreads_tot, hw_opt->nthreads_tmpi);
769         }
770
771         if (hw_opt->nthreads_omp > 0 &&
772             hw_opt->nthreads_tot % hw_opt->nthreads_omp != 0)
773         {
774             gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of OpenMP threads requested (%d)",
775                       hw_opt->nthreads_tot, hw_opt->nthreads_omp);
776         }
777     }
778
779     if (hw_opt->nthreads_tot > 0)
780     {
781         if (hw_opt->nthreads_omp > hw_opt->nthreads_tot)
782         {
783             gmx_fatal(FARGS, "You requested %d OpenMP threads with %d total threads. Choose a total number of threads that is a multiple of the number of OpenMP threads.",
784                       hw_opt->nthreads_omp, hw_opt->nthreads_tot);
785         }
786
787         if (hw_opt->nthreads_tmpi > hw_opt->nthreads_tot)
788         {
789             gmx_fatal(FARGS, "You requested %d thread-MPI ranks with %d total threads. Choose a total number of threads that is a multiple of the number of thread-MPI ranks.",
790                       hw_opt->nthreads_tmpi, hw_opt->nthreads_tot);
791         }
792     }
793
794     if (debug)
795     {
796         print_hw_opt(debug, hw_opt);
797     }
798
799     /* Asserting this simplifies the hardware resource division later
800      * on. */
801     GMX_RELEASE_ASSERT(!(hw_opt->nthreads_omp_pme >= 1 && hw_opt->nthreads_omp <= 0),
802                        "PME thread count should only be set when the normal thread count is also set");
803 }
804
805 void check_and_update_hw_opt_2(gmx_hw_opt_t *hw_opt,
806                                int           cutoff_scheme)
807 {
808     if (cutoff_scheme == ecutsGROUP)
809     {
810         /* We only have OpenMP support for PME only nodes */
811         if (hw_opt->nthreads_omp > 1)
812         {
813             gmx_fatal(FARGS, "OpenMP threads have been requested with cut-off scheme %s, but these are only supported with cut-off scheme %s",
814                       ecutscheme_names[cutoff_scheme],
815                       ecutscheme_names[ecutsVERLET]);
816         }
817         hw_opt->nthreads_omp = 1;
818     }
819 }
820
821 void checkAndUpdateRequestedNumOpenmpThreads(gmx_hw_opt_t        *hw_opt,
822                                              const gmx_hw_info_t &hwinfo,
823                                              const t_commrec     *cr,
824                                              PmeRunMode           pmeRunMode,
825                                              const gmx_mtop_t    &mtop)
826 {
827 #if GMX_THREAD_MPI
828     GMX_RELEASE_ASSERT(hw_opt->nthreads_tmpi >= 1, "Must have at least one thread-MPI rank");
829
830     /* If the user set the total number of threads on the command line
831      * and did not specify the number of OpenMP threads, set the latter here.
832      */
833     if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp <= 0)
834     {
835         hw_opt->nthreads_omp = hw_opt->nthreads_tot/hw_opt->nthreads_tmpi;
836
837         if (!bHasOmpSupport && hw_opt->nthreads_omp > 1)
838         {
839             gmx_fatal(FARGS, "You (indirectly) asked for OpenMP threads by setting -nt > -ntmpi, but GROMACS was compiled without OpenMP support");
840         }
841     }
842 #endif
843
844     /* With both non-bonded and PME on GPU, the work left on the CPU is often
845      * (much) slower with SMT than without SMT. This is mostly the case with
846      * few atoms per core. Thus, if the number of threads is set to auto,
847      * we turn off SMT in that case. Note that PME on GPU implies that also
848      * the non-bonded are computed on the GPU.
849      * We only need to do this when the number of hardware theads is larger
850      * than the number of cores. Note that a queuing system could limit
851      * the number of hardware threads available, but we are not trying to be
852      * too smart here in that case.
853      */
854     /* The number 10000 atoms per core is a reasonable threshold
855      * for Intel x86 and AMD Threadripper. Note that around this regime
856      * the performance is not very sensitive to the SMT choice.
857      */
858     constexpr int c_numAtomsPerCoreSmtThreshold = 10000;
859
860     if (bHasOmpSupport &&
861         hw_opt->nthreads_omp <= 0 &&
862         hwinfo.nthreads_hw_avail > hwinfo.hardwareTopology->numberOfCores() &&
863         pmeRunMode == PmeRunMode::GPU)
864     {
865         /* Note that the queing system might have limited us from using
866          * all detected ncore_tot physical cores. We are currently not
867          * checking for that here.
868          */
869         int numAtomsPerPhysicalCore = mtop.natoms/hwinfo.ncore_tot;
870         if (numAtomsPerPhysicalCore < c_numAtomsPerCoreSmtThreshold)
871         {
872             int numRanksInThisNode = (cr ? cr->nrank_intranode : 1);
873             /* Choose one OpenMP thread per physical core */
874             hw_opt->nthreads_omp = std::max(1, hwinfo.hardwareTopology->numberOfCores()/numRanksInThisNode);
875         }
876     }
877
878     GMX_RELEASE_ASSERT(bHasOmpSupport || hw_opt->nthreads_omp == 1, "Without OpenMP support, only one thread per rank can be used");
879
880     /* We are done with updating nthreads_omp, we can set nthreads_omp_pme */
881     if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
882     {
883         hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
884     }
885
886     if (debug)
887     {
888         print_hw_opt(debug, hw_opt);
889     }
890 }