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