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