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