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