Clean up PME rank variables
[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, 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 <cstdlib>
50 #include <cstring>
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/math/functions.h"
60 #include "gromacs/mdlib/gmx_omp_nthreads.h"
61 #include "gromacs/mdrunutility/multisim.h"
62 #include "gromacs/mdtypes/commrec.h"
63 #include "gromacs/mdtypes/inputrec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/topology/mtop_util.h"
66 #include "gromacs/topology/topology.h"
67 #include "gromacs/utility/baseversion.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/gmxassert.h"
70 #include "gromacs/utility/logger.h"
71 #include "gromacs/utility/physicalnodecommunicator.h"
72 #include "gromacs/utility/stringutil.h"
73
74
75 /* DISCLAIMER: All the atom count and thread numbers below are heuristic.
76  * The real switching points will depend on the system simulation,
77  * the algorithms used and the hardware it's running on, as well as if there
78  * are other jobs running on the same machine. We try to take into account
79  * factors that have a large influence, such as recent Intel CPUs being
80  * much better at wide multi-threading. The remaining factors should
81  * (hopefully) have a small influence, such that the performance just before
82  * and after a switch point doesn't change too much.
83  */
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 constexpr 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 constexpr 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 constexpr int nthreads_omp_faster_default   =  8;
116 constexpr int nthreads_omp_faster_Nehalem   = 12;
117 constexpr int nthreads_omp_faster_Intel_AVX = 16;
118 constexpr 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 constexpr 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 constexpr int nthreads_omp_mpi_ok_max              =  8;
136 constexpr int nthreads_omp_mpi_ok_min_cpu          =  1;
137 constexpr int nthreads_omp_mpi_ok_min_gpu          =  2;
138 constexpr int nthreads_omp_mpi_target_max          =  6;
139
140 /**@}*/
141
142 /*! \brief Returns the maximum OpenMP thread count for which using a single MPI rank
143  * should be faster than using multiple ranks with the same total thread count.
144  */
145 static int nthreads_omp_faster(const gmx::CpuInfo &cpuInfo, gmx_bool bUseGPU)
146 {
147     int nth;
148
149     if (cpuInfo.vendor() == gmx::CpuInfo::Vendor::Intel &&
150         cpuInfo.feature(gmx::CpuInfo::Feature::X86_Avx))
151     {
152         nth = nthreads_omp_faster_Intel_AVX;
153     }
154     else if (gmx::cpuIsX86Nehalem(cpuInfo))
155     {
156         // Intel Nehalem
157         nth = nthreads_omp_faster_Nehalem;
158     }
159     else if ((cpuInfo.vendor() == gmx::CpuInfo::Vendor::Amd && cpuInfo.family() >= 23) ||
160              cpuInfo.vendor() == gmx::CpuInfo::Vendor::Hygon)
161     {
162         // AMD Ryzen || Hygon Dhyana
163         nth = nthreads_omp_faster_AMD_Ryzen;
164     }
165     else
166     {
167         nth = nthreads_omp_faster_default;
168     }
169
170     if (bUseGPU)
171     {
172         nth *= nthreads_omp_faster_gpu_fac;
173     }
174
175     nth = std::min(nth, GMX_OPENMP_MAX_THREADS);
176
177     return nth;
178 }
179
180 /*! \brief Returns that maximum OpenMP thread count that passes the efficiency check */
181 gmx_unused static int nthreads_omp_efficient_max(int gmx_unused       nrank,
182                                                  const gmx::CpuInfo  &cpuInfo,
183                                                  gmx_bool             bUseGPU)
184 {
185     if (GMX_OPENMP && GMX_MPI && (nrank > 1))
186     {
187         return nthreads_omp_mpi_ok_max;
188     }
189     else
190     {
191         return nthreads_omp_faster(cpuInfo, bUseGPU);
192     }
193 }
194
195 /*! \brief Return the number of thread-MPI ranks to use.
196  * This is chosen such that we can always obey our own efficiency checks.
197  */
198 gmx_unused static int get_tmpi_omp_thread_division(const gmx_hw_info_t *hwinfo,
199                                                    const gmx_hw_opt_t  &hw_opt,
200                                                    int                  nthreads_tot,
201                                                    int                  ngpu)
202 {
203     int                 nrank;
204     const gmx::CpuInfo &cpuInfo = *hwinfo->cpuInfo;
205
206     GMX_RELEASE_ASSERT(nthreads_tot > 0, "There must be at least one thread per rank");
207
208     /* There are no separate PME nodes here, as we ensured in
209      * check_and_update_hw_opt that nthreads_tmpi>0 with PME nodes
210      * and a conditional ensures we would not have ended up here.
211      * Note that separate PME nodes might be switched on later.
212      */
213     if (ngpu > 0)
214     {
215         if (hw_opt.nthreads_omp > 0)
216         {
217             /* In this case it is unclear if we should use 1 rank per GPU
218              * or more or less, so we require also setting the number of ranks.
219              */
220             gmx_fatal(FARGS, "When using GPUs, setting the number of OpenMP threads without specifying the number "
221                       "of ranks can lead to conflicting demands. Please specify the number of thread-MPI ranks "
222                       "as well (option -ntmpi).");
223         }
224
225         nrank = ngpu;
226
227         /* When the user sets nthreads_omp, we can end up oversubscribing CPU cores
228          * if we simply start as many ranks as GPUs. To avoid this, we start as few
229          * tMPI ranks as necessary to avoid oversubscription and instead leave GPUs idle.
230          * If the user does not set the number of OpenMP threads, nthreads_omp==0 and
231          * this code has no effect.
232          */
233         GMX_RELEASE_ASSERT(hw_opt.nthreads_omp >= 0, "nthreads_omp is negative, but previous checks should "
234                            "have prevented this");
235         while (nrank*hw_opt.nthreads_omp > hwinfo->nthreads_hw_avail && nrank > 1)
236         {
237             nrank--;
238         }
239
240         if (nthreads_tot < nrank)
241         {
242             /* #thread < #gpu is very unlikely, but if so: waste gpu(s) */
243             nrank = nthreads_tot;
244         }
245         else if (nthreads_tot > nthreads_omp_faster(cpuInfo, ngpu > 0) ||
246                  (ngpu > 1 && nthreads_tot/ngpu > nthreads_omp_mpi_target_max))
247         {
248             /* The high OpenMP thread count will likely result in sub-optimal
249              * performance. Increase the rank count to reduce the thread count
250              * per rank. This will lead to GPU sharing by MPI ranks/threads.
251              */
252             int nshare;
253
254             /* Increase the rank count as long as have we more than 6 OpenMP
255              * threads per rank or the number of hardware threads is not
256              * divisible by the rank count. Don't go below 2 OpenMP threads.
257              */
258             nshare = 1;
259             do
260             {
261                 nshare++;
262                 nrank = ngpu*nshare;
263             }
264             while (nthreads_tot/nrank > nthreads_omp_mpi_target_max ||
265                    (nthreads_tot/(ngpu*(nshare + 1)) >= nthreads_omp_mpi_ok_min_gpu && nthreads_tot % nrank != 0));
266         }
267     }
268     else if (hw_opt.nthreads_omp > 0)
269     {
270         /* Here we could oversubscribe, when we do, we issue a warning later */
271         nrank = std::max(1, nthreads_tot/hw_opt.nthreads_omp);
272     }
273     else
274     {
275         if (nthreads_tot <= nthreads_omp_faster(cpuInfo, ngpu > 0))
276         {
277             /* Use pure OpenMP parallelization */
278             nrank = 1;
279         }
280         else
281         {
282             /* Don't use OpenMP parallelization */
283             nrank = nthreads_tot;
284         }
285     }
286
287     return nrank;
288 }
289
290 //! Return whether hyper threading is enabled.
291 static bool
292 gmxSmtIsEnabled(const gmx::HardwareTopology &hwTop)
293 {
294     return (hwTop.supportLevel() >= gmx::HardwareTopology::SupportLevel::Basic && hwTop.machine().sockets[0].cores[0].hwThreads.size() > 1);
295 }
296
297 namespace
298 {
299
300 //! Handles checks for algorithms that must use a single rank.
301 class SingleRankChecker
302 {
303     public:
304         SingleRankChecker() : value_(false) {}
305         /*! \brief Call this function for each possible condition
306             under which a single rank is required, along with a string
307             describing the constraint when it is applied. */
308         void applyConstraint(bool condition, const char *description)
309         {
310             if (condition)
311             {
312                 value_ = true;
313                 reasons_.push_back(gmx::formatString("%s only supports a single rank.", description));
314             }
315         }
316         //! After applying any conditions, is a single rank required?
317         bool mustUseOneRank() const
318         {
319             return value_;
320         }
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     private:
329         bool                     value_;
330         std::vector<std::string> reasons_;
331 };
332
333 } // namespace
334
335 /* Get the number of MPI ranks to use for thread-MPI based on how many
336  * were requested, which algorithms we're using,
337  * and how many particles there are.
338  * At the point we have already called check_and_update_hw_opt.
339  * Thus all options should be internally consistent and consistent
340  * with the hardware, except that ntmpi could be larger than #GPU.
341  */
342 int get_nthreads_mpi(const gmx_hw_info_t    *hwinfo,
343                      gmx_hw_opt_t           *hw_opt,
344                      const std::vector<int> &gpuIdsToUse,
345                      bool                    nonbondedOnGpu,
346                      bool                    pmeOnGpu,
347                      const t_inputrec       *inputrec,
348                      const gmx_mtop_t       *mtop,
349                      const gmx::MDLogger    &mdlog,
350                      bool                    doMembed)
351 {
352     int                          nthreads_hw, nthreads_tot_max, nrank, ngpu;
353     int                          min_atoms_per_mpi_rank;
354
355     const gmx::CpuInfo          &cpuInfo = *hwinfo->cpuInfo;
356     const gmx::HardwareTopology &hwTop   = *hwinfo->hardwareTopology;
357
358     if (pmeOnGpu)
359     {
360         GMX_RELEASE_ASSERT((EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)) &&
361                            pme_gpu_supports_build(nullptr) &&
362                            pme_gpu_supports_hardware(*hwinfo, nullptr) &&
363                            pme_gpu_supports_input(*inputrec, *mtop, nullptr),
364                            "PME can't be on GPUs unless we are using PME");
365
366         // PME on GPUs supports a single PME rank with PP running on the same or few other ranks.
367         // For now, let's treat separate PME GPU rank as opt-in.
368         if (hw_opt->nthreads_tmpi < 1)
369         {
370             return 1;
371         }
372     }
373
374     {
375         /* Check if an algorithm does not support parallel simulation.  */
376         // TODO This might work better if e.g. implemented algorithms
377         // had to define a function that returns such requirements,
378         // and a description string.
379         SingleRankChecker checker;
380         checker.applyConstraint(inputrec->eI == eiLBFGS, "L-BFGS minimization");
381         checker.applyConstraint(inputrec->coulombtype == eelEWALD, "Plain Ewald electrostatics");
382         checker.applyConstraint(doMembed, "Membrane embedding");
383         bool useOrientationRestraints = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
384         checker.applyConstraint(useOrientationRestraints, "Orientation restraints");
385         if (checker.mustUseOneRank())
386         {
387             std::string message = checker.getMessage();
388             if (hw_opt->nthreads_tmpi > 1)
389             {
390                 gmx_fatal(FARGS, "%s However, you asked for more than 1 thread-MPI rank, so mdrun cannot continue. "
391                           "Choose a single rank, or a different algorithm.", message.c_str());
392             }
393             GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted("%s Choosing to use only a single thread-MPI rank.", message.c_str());
394             return 1;
395         }
396     }
397
398     if (hw_opt->nthreads_tmpi > 0)
399     {
400         /* Trivial, return the user's choice right away */
401         return hw_opt->nthreads_tmpi;
402     }
403
404     // Now implement automatic selection of number of thread-MPI ranks
405     nthreads_hw = hwinfo->nthreads_hw_avail;
406
407     if (nthreads_hw <= 0)
408     {
409         /* This should normally not happen, but if it does, we handle it */
410         gmx_fatal(FARGS, "The number of available hardware threads can not be detected, please specify the number of "
411                   "MPI ranks and the number of OpenMP threads (if supported) manually with options "
412                   "-ntmpi and -ntomp, respectively");
413     }
414
415     /* How many total (#tMPI*#OpenMP) threads can we start? */
416     if (hw_opt->nthreads_tot > 0)
417     {
418         nthreads_tot_max = hw_opt->nthreads_tot;
419     }
420     else
421     {
422         nthreads_tot_max = nthreads_hw;
423     }
424
425     /* nonbondedOnGpu might be false e.g. because this simulation uses
426      * the group scheme, or is a rerun with energy groups. */
427     ngpu = (nonbondedOnGpu ? gmx::ssize(gpuIdsToUse) : 0);
428
429     nrank =
430         get_tmpi_omp_thread_division(hwinfo, *hw_opt, nthreads_tot_max, ngpu);
431
432     if (inputrec->eI == eiNM || EI_TPI(inputrec->eI))
433     {
434         /* Dims/steps are divided over the nodes iso splitting the atoms.
435          * With NM we can't have more ranks than #atoms*#dim. With TPI it's
436          * unlikely we have fewer atoms than ranks, and if so, communication
437          * would become a bottleneck, so we set the limit to 1 atom/rank.
438          */
439         min_atoms_per_mpi_rank = 1;
440     }
441     else
442     {
443         if (ngpu >= 1)
444         {
445             min_atoms_per_mpi_rank = min_atoms_per_gpu;
446         }
447         else
448         {
449             min_atoms_per_mpi_rank = min_atoms_per_mpi_thread;
450         }
451     }
452
453     if (mtop->natoms/nrank < min_atoms_per_mpi_rank)
454     {
455         int nrank_new;
456
457         /* the rank number was chosen automatically, but there are too few
458            atoms per rank, so we need to reduce the rank count */
459         nrank_new = std::max(1, mtop->natoms/min_atoms_per_mpi_rank);
460
461         /* Avoid partial use of Hyper-Threading */
462         if (gmxSmtIsEnabled(hwTop) &&
463             nrank_new > nthreads_hw/2 && nrank_new < nthreads_hw)
464         {
465             nrank_new = nthreads_hw/2;
466         }
467
468         /* If the user specified the total thread count, ensure this is
469          * divisible by the number of ranks.
470          * It is quite likely that we have too many total threads compared
471          * to the size of the system, but if the user asked for this many
472          * threads we should respect that.
473          */
474         while (hw_opt->nthreads_tot > 0 &&
475                hw_opt->nthreads_tot % nrank_new != 0)
476         {
477             nrank_new--;
478         }
479
480         /* Avoid large prime numbers in the rank count */
481         if (nrank_new >= 6)
482         {
483             /* Use only 6,8,10 with additional factors of 2 */
484             int fac;
485
486             fac = 2;
487             while (3*fac*2 <= nrank_new)
488             {
489                 fac *= 2;
490             }
491
492             nrank_new = (nrank_new/fac)*fac;
493         }
494         else
495         {
496             /* Avoid 5, since small system won't fit 5 domains along
497              * a dimension. This might lead to waisting some cores, but this
498              * will have a small impact in this regime of very small systems.
499              */
500             if (nrank_new == 5)
501             {
502                 nrank_new = 4;
503             }
504         }
505
506         nrank = nrank_new;
507
508         /* We reduced the number of tMPI ranks, which means we might violate
509          * our own efficiency checks if we simply use all hardware threads.
510          */
511         if (GMX_OPENMP && hw_opt->nthreads_omp <= 0 && hw_opt->nthreads_tot <= 0)
512         {
513             /* The user set neither the total nor the OpenMP thread count,
514              * we should use all hardware threads, unless we will violate
515              * our own efficiency limitation on the thread count.
516              */
517             int  nt_omp_max;
518
519             nt_omp_max = nthreads_omp_efficient_max(nrank, cpuInfo, ngpu >= 1);
520
521             if (nrank*nt_omp_max < hwinfo->nthreads_hw_avail)
522             {
523                 /* Limit the number of OpenMP threads to start */
524                 hw_opt->nthreads_omp = nt_omp_max;
525             }
526         }
527
528         fprintf(stderr, "\n");
529         fprintf(stderr, "NOTE: Parallelization is limited by the small number of atoms,\n");
530         fprintf(stderr, "      only starting %d thread-MPI ranks.\n", nrank);
531         fprintf(stderr, "      You can use the -nt and/or -ntmpi option to optimize the number of threads.\n\n");
532     }
533
534     return nrank;
535 }
536
537
538 void check_resource_division_efficiency(const gmx_hw_info_t *hwinfo,
539                                         bool                 willUsePhysicalGpu,
540                                         gmx_bool             bNtOmpOptionSet,
541                                         t_commrec           *cr,
542                                         const gmx::MDLogger &mdlog)
543 {
544 #if GMX_OPENMP && GMX_MPI
545     GMX_UNUSED_VALUE(hwinfo);
546
547     int         nth_omp_min, nth_omp_max;
548     char        buf[1000];
549     const char *mpi_option = GMX_THREAD_MPI ?  " (option -ntmpi)" : "";
550
551     /* This function should be called after thread-MPI (when configured) and
552      * OpenMP have been initialized. Check that here.
553      */
554     if (GMX_THREAD_MPI)
555     {
556         GMX_RELEASE_ASSERT(nthreads_omp_faster_default >= nthreads_omp_mpi_ok_max,
557                            "Inconsistent OpenMP thread count default values");
558     }
559     GMX_RELEASE_ASSERT(gmx_omp_nthreads_get(emntDefault) >= 1, "Must have at least one OpenMP thread");
560
561     nth_omp_min = gmx_omp_nthreads_get(emntDefault);
562     nth_omp_max = gmx_omp_nthreads_get(emntDefault);
563
564     bool anyRankIsUsingGpus = willUsePhysicalGpu;
565     /* Thread-MPI seems to have a bug with reduce on 1 node, so use a cond. */
566     if (cr->nnodes > 1)
567     {
568         int count[3], count_max[3];
569
570         count[0] = -nth_omp_min;
571         count[1] = nth_omp_max;
572         count[2] = int(willUsePhysicalGpu);
573
574         MPI_Allreduce(count, count_max, 3, MPI_INT, MPI_MAX, cr->mpi_comm_mysim);
575
576         /* In case of an inhomogeneous run setup we use the maximum counts */
577         nth_omp_min        = -count_max[0];
578         nth_omp_max        = count_max[1];
579         anyRankIsUsingGpus = count_max[2] > 0;
580     }
581
582     int nthreads_omp_mpi_ok_min;
583
584     if (!anyRankIsUsingGpus)
585     {
586         nthreads_omp_mpi_ok_min = nthreads_omp_mpi_ok_min_cpu;
587     }
588     else
589     {
590         /* With GPUs we set the minimum number of OpenMP threads to 2 to catch
591          * cases where the user specifies #ranks == #cores.
592          */
593         nthreads_omp_mpi_ok_min = nthreads_omp_mpi_ok_min_gpu;
594     }
595
596     if (DOMAINDECOMP(cr))
597     {
598         if (nth_omp_max < nthreads_omp_mpi_ok_min ||
599             nth_omp_max > nthreads_omp_mpi_ok_max)
600         {
601             /* Note that we print target_max here, not ok_max */
602             sprintf(buf, "Your choice of number of MPI ranks and amount of resources results in using %d OpenMP "
603                     "threads per rank, which is most likely inefficient. The optimum is usually between %d and"
604                     " %d threads per rank.",
605                     nth_omp_max,
606                     nthreads_omp_mpi_ok_min,
607                     nthreads_omp_mpi_target_max);
608
609             if (bNtOmpOptionSet)
610             {
611                 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted("NOTE: %s", buf);
612             }
613             else
614             {
615                 /* This fatal error, and the one below, is nasty, but it's
616                  * probably the only way to ensure that all users don't waste
617                  * a lot of resources, since many users don't read logs/stderr.
618                  */
619                 gmx_fatal(FARGS, "%s If you want to run with this setup, specify the -ntomp option. But we suggest to "
620                           "change the number of MPI ranks%s.",
621                           buf, mpi_option);
622             }
623         }
624     }
625 #else // !GMX_OPENMP || ! GMX_MPI
626     GMX_UNUSED_VALUE(bNtOmpOptionSet);
627     GMX_UNUSED_VALUE(willUsePhysicalGpu);
628     GMX_UNUSED_VALUE(cr);
629     GMX_UNUSED_VALUE(nthreads_omp_mpi_ok_max);
630     GMX_UNUSED_VALUE(nthreads_omp_mpi_ok_min_cpu);
631     /* Check if we have more than 1 physical core, if detected,
632      * or more than 1 hardware thread if physical cores were not detected.
633      */
634     if (!GMX_OPENMP && !GMX_MPI && hwinfo->hardwareTopology->numberOfCores() > 1)
635     {
636         GMX_LOG(mdlog.warning).asParagraph().appendText(
637                 "NOTE: GROMACS was compiled without OpenMP and (thread-)MPI support, can only use a single CPU core");
638     }
639 #endif // end GMX_OPENMP && GMX_MPI
640 }
641
642
643 //! Dump a \c hw_opt to \c fp.
644 static void print_hw_opt(FILE *fp, const gmx_hw_opt_t *hw_opt)
645 {
646     fprintf(fp, "hw_opt: nt %d ntmpi %d ntomp %d ntomp_pme %d gpu_id '%s' gputasks '%s'\n",
647             hw_opt->nthreads_tot,
648             hw_opt->nthreads_tmpi,
649             hw_opt->nthreads_omp,
650             hw_opt->nthreads_omp_pme,
651             hw_opt->gpuIdsAvailable.c_str(),
652             hw_opt->userGpuTaskAssignment.c_str());
653 }
654
655 void checkAndUpdateHardwareOptions(const gmx::MDLogger &mdlog,
656                                    gmx_hw_opt_t        *hw_opt,
657                                    const bool           isSimulationMasterRank,
658                                    const int            nPmeRanks)
659 {
660     /* Currently hw_opt only contains default settings or settings supplied
661      * by the user on the command line.
662      */
663     if (hw_opt->nthreads_omp < 0)
664     {
665         gmx_fatal(FARGS, "The number of OpenMP threads supplied on the command line is %d, which is negative "
666                   "and not allowed", hw_opt->nthreads_omp);
667     }
668
669     /* Check for OpenMP settings stored in environment variables, which can
670      * potentially be different on different MPI ranks.
671      */
672     gmx_omp_nthreads_read_env(mdlog, &hw_opt->nthreads_omp);
673
674     /* Check restrictions on the user supplied options before modifying them.
675      * TODO: Put the user values in a const struct and preserve them.
676      */
677     if (!GMX_THREAD_MPI)
678     {
679
680         if (hw_opt->nthreads_tot > 0)
681         {
682             gmx_fatal(FARGS, "Setting the total number of threads is only supported with thread-MPI and GROMACS was "
683                       "compiled without thread-MPI");
684         }
685         if (hw_opt->nthreads_tmpi > 0)
686         {
687             gmx_fatal(FARGS, "Setting the number of thread-MPI ranks is only supported with thread-MPI and GROMACS was "
688                       "compiled without thread-MPI");
689         }
690     }
691     /* With thread-MPI the master thread sets hw_opt->totNumThreadsIsAuto.
692      * The other threads receive a partially processed hw_opt from the master
693      * thread and should not set hw_opt->totNumThreadsIsAuto again.
694      */
695     if (!GMX_THREAD_MPI || isSimulationMasterRank)
696     {
697         /* Check if mdrun is free to choose the total number of threads */
698         hw_opt->totNumThreadsIsAuto = (hw_opt->nthreads_omp == 0 && hw_opt->nthreads_omp_pme == 0 && hw_opt->nthreads_tot == 0);
699     }
700
701     if (GMX_OPENMP)
702     {
703         /* Check restrictions on PME thread related options set by the user */
704
705         if (hw_opt->nthreads_omp_pme > 0 && hw_opt->nthreads_omp <= 0)
706         {
707             gmx_fatal(FARGS, "You need to specify -ntomp in addition to -ntomp_pme");
708         }
709
710         if (hw_opt->nthreads_omp_pme >= 1 &&
711             hw_opt->nthreads_omp_pme != hw_opt->nthreads_omp &&
712             nPmeRanks <= 0)
713         {
714             /* This can result in a fatal error on many MPI ranks,
715              * but since the thread count can differ per rank,
716              * we can't easily avoid this.
717              */
718             gmx_fatal(FARGS, "You need to explicitly specify the number of PME ranks (-npme) when using "
719                       "different numbers of OpenMP threads for PP and PME ranks");
720         }
721     }
722     else
723     {
724         /* GROMACS was configured without OpenMP support */
725
726         if (hw_opt->nthreads_omp > 1 || hw_opt->nthreads_omp_pme > 1)
727         {
728             gmx_fatal(FARGS, "More than 1 OpenMP thread requested, but GROMACS was compiled without OpenMP support");
729         }
730         hw_opt->nthreads_omp     = 1;
731         hw_opt->nthreads_omp_pme = 1;
732     }
733
734     if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp_pme <= 0)
735     {
736         /* We have the same number of OpenMP threads for PP and PME ranks,
737          * thus we can perform several consistency checks.
738          */
739         if (hw_opt->nthreads_tmpi > 0 &&
740             hw_opt->nthreads_omp > 0 &&
741             hw_opt->nthreads_tot != hw_opt->nthreads_tmpi*hw_opt->nthreads_omp)
742         {
743             gmx_fatal(FARGS, "The total number of threads requested (%d) does not match the thread-MPI ranks (%d) "
744                       "times the OpenMP threads (%d) requested",
745                       hw_opt->nthreads_tot, hw_opt->nthreads_tmpi, hw_opt->nthreads_omp);
746         }
747
748         if (hw_opt->nthreads_tmpi > 0 &&
749             hw_opt->nthreads_tot % hw_opt->nthreads_tmpi != 0)
750         {
751             gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of thread-MPI "
752                       "ranks requested (%d)",
753                       hw_opt->nthreads_tot, hw_opt->nthreads_tmpi);
754         }
755
756         if (hw_opt->nthreads_omp > 0 &&
757             hw_opt->nthreads_tot % hw_opt->nthreads_omp != 0)
758         {
759             gmx_fatal(FARGS, "The total number of threads requested (%d) is not divisible by the number of OpenMP "
760                       "threads requested (%d)",
761                       hw_opt->nthreads_tot, hw_opt->nthreads_omp);
762         }
763     }
764
765     if (hw_opt->nthreads_tot > 0)
766     {
767         if (hw_opt->nthreads_omp > hw_opt->nthreads_tot)
768         {
769             gmx_fatal(FARGS, "You requested %d OpenMP threads with %d total threads. Choose a total number of threads "
770                       "that is a multiple of the number of OpenMP threads.",
771                       hw_opt->nthreads_omp, hw_opt->nthreads_tot);
772         }
773
774         if (hw_opt->nthreads_tmpi > hw_opt->nthreads_tot)
775         {
776             gmx_fatal(FARGS, "You requested %d thread-MPI ranks with %d total threads. Choose a total number of "
777                       "threads that is a multiple of the number of thread-MPI ranks.",
778                       hw_opt->nthreads_tmpi, hw_opt->nthreads_tot);
779         }
780     }
781
782     if (GMX_THREAD_MPI && nPmeRanks > 0 && hw_opt->nthreads_tmpi <= 0)
783     {
784         gmx_fatal(FARGS, "You need to explicitly specify the number of MPI threads (-ntmpi) when using separate PME ranks");
785     }
786
787     if (debug)
788     {
789         print_hw_opt(debug, hw_opt);
790     }
791
792     /* Asserting this simplifies the hardware resource division later
793      * on. */
794     GMX_RELEASE_ASSERT(!(hw_opt->nthreads_omp_pme >= 1 && hw_opt->nthreads_omp <= 0),
795                        "PME thread count should only be set when the normal thread count is also set");
796 }
797
798 void checkAndUpdateRequestedNumOpenmpThreads(gmx_hw_opt_t         *hw_opt,
799                                              const gmx_hw_info_t  &hwinfo,
800                                              const t_commrec      *cr,
801                                              const gmx_multisim_t *ms,
802                                              int                   numRanksOnThisNode,
803                                              PmeRunMode            pmeRunMode,
804                                              const gmx_mtop_t     &mtop)
805 {
806     if (GMX_THREAD_MPI)
807     {
808
809         GMX_RELEASE_ASSERT(hw_opt->nthreads_tmpi >= 1, "Must have at least one thread-MPI rank");
810
811         /* If the user set the total number of threads on the command line
812          * and did not specify the number of OpenMP threads, set the latter here.
813          */
814         if (hw_opt->nthreads_tot > 0 && hw_opt->nthreads_omp <= 0)
815         {
816             hw_opt->nthreads_omp = hw_opt->nthreads_tot / hw_opt->nthreads_tmpi;
817
818             if (!GMX_OPENMP && hw_opt->nthreads_omp > 1)
819             {
820                 gmx_fatal(FARGS, "You (indirectly) asked for OpenMP threads by setting -nt > -ntmpi, but GROMACS was "
821                           "compiled without OpenMP support");
822             }
823         }
824     }
825     /* With both non-bonded and PME on GPU, the work left on the CPU is often
826      * (much) slower with SMT than without SMT. This is mostly the case with
827      * few atoms per core. Thus, if the number of threads is set to auto,
828      * we turn off SMT in that case. Note that PME on GPU implies that also
829      * the non-bonded are computed on the GPU.
830      * We only need to do this when the number of hardware theads is larger
831      * than the number of cores. Note that a queuing system could limit
832      * the number of hardware threads available, but we are not trying to be
833      * too smart here in that case.
834      */
835     /* The thread reduction and synchronization costs go up roughy quadratically
836      * with the threads count, so we apply a threshold quadratic in #cores.
837      * Also more cores per GPU usually means the CPU gets faster than the GPU.
838      * The number 1000 atoms per core^2 is a reasonable threshold
839      * for Intel x86 and AMD Threadripper.
840      */
841     constexpr int c_numAtomsPerCoreSquaredSmtThreshold = 1000;
842
843     /* Prepare conditions for deciding if we should disable SMT.
844      * We currently only limit SMT for simulations using a single rank.
845      * TODO: Consider limiting also for multi-rank simulations.
846      */
847     bool canChooseNumOpenmpThreads      = (GMX_OPENMP && hw_opt->nthreads_omp <= 0);
848     bool haveSmtSupport                 = (hwinfo.hardwareTopology->supportLevel() >= gmx::HardwareTopology::SupportLevel::Basic &&
849                                            hwinfo.hardwareTopology->machine().logicalProcessorCount > hwinfo.hardwareTopology->numberOfCores());
850     bool simRunsSingleRankNBAndPmeOnGpu = (cr->nnodes == 1 && pmeRunMode == PmeRunMode::GPU);
851
852     if (canChooseNumOpenmpThreads && haveSmtSupport &&
853         simRunsSingleRankNBAndPmeOnGpu)
854     {
855         /* Note that the queing system might have limited us from using
856          * all detected ncore_tot physical cores. We are currently not
857          * checking for that here.
858          */
859         int numRanksTot     = cr->nnodes*(isMultiSim(ms) ? ms->nsim : 1);
860         int numAtomsPerRank = mtop.natoms/cr->nnodes;
861         int numCoresPerRank = hwinfo.ncore_tot/numRanksTot;
862         if (numAtomsPerRank < c_numAtomsPerCoreSquaredSmtThreshold*gmx::square(numCoresPerRank))
863         {
864             /* Choose one OpenMP thread per physical core */
865             hw_opt->nthreads_omp = std::max(1, hwinfo.hardwareTopology->numberOfCores()/numRanksOnThisNode);
866         }
867     }
868
869     GMX_RELEASE_ASSERT(GMX_OPENMP || hw_opt->nthreads_omp == 1, "Without OpenMP support, only one thread per rank can be used");
870
871     /* We are done with updating nthreads_omp, we can set nthreads_omp_pme */
872     if (hw_opt->nthreads_omp_pme <= 0 && hw_opt->nthreads_omp > 0)
873     {
874         hw_opt->nthreads_omp_pme = hw_opt->nthreads_omp;
875     }
876
877     if (debug)
878     {
879         print_hw_opt(debug, hw_opt);
880     }
881 }
882
883 namespace gmx
884 {
885
886 void checkHardwareOversubscription(int                             numThreadsOnThisRank,
887                                    int                             rank,
888                                    const HardwareTopology         &hwTop,
889                                    const PhysicalNodeCommunicator &comm,
890                                    const MDLogger                 &mdlog)
891 {
892     if (hwTop.supportLevel() < HardwareTopology::SupportLevel::LogicalProcessorCount)
893     {
894         /* There is nothing we can check */
895         return;
896     }
897
898     int numRanksOnThisNode   = comm.size_;
899     int numThreadsOnThisNode = numThreadsOnThisRank;
900     /* Avoid MPI calls with uninitialized thread-MPI communicators */
901     if (comm.size_ > 1)
902     {
903 #if GMX_MPI
904         /* Count the threads within this physical node */
905         MPI_Allreduce(&numThreadsOnThisRank, &numThreadsOnThisNode, 1, MPI_INT, MPI_SUM, comm.comm_);
906 #endif
907     }
908
909     if (numThreadsOnThisNode > hwTop.machine().logicalProcessorCount)
910     {
911         std::string mesg = "WARNING: ";
912         if (GMX_LIB_MPI)
913         {
914             mesg += formatString("On rank %d: o", rank);
915         }
916         else
917         {
918             mesg += "O";
919         }
920         mesg     += formatString("versubscribing the available %d logical CPU cores", hwTop.machine().logicalProcessorCount);
921         if (GMX_LIB_MPI)
922         {
923             mesg += " per node";
924         }
925         mesg     += formatString(" with %d ", numThreadsOnThisNode);
926         if (numRanksOnThisNode == numThreadsOnThisNode)
927         {
928             if (GMX_THREAD_MPI)
929             {
930                 mesg += "thread-MPI threads.";
931             }
932             else
933             {
934                 mesg += "MPI processes.";
935             }
936         }
937         else
938         {
939             mesg += "threads.";
940         }
941         mesg     += "\n         This will cause considerable performance loss.";
942         /* Note that only the master rank logs to stderr and only ranks
943          * with an open log file write to log.
944          * TODO: When we have a proper parallel logging framework,
945          *       the framework should add the rank and node numbers.
946          */
947         GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted("%s", mesg.c_str());
948     }
949 }
950
951 }  // namespace gmx