2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016, 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.
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.
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.
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.
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.
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.
38 #include "resource-division.h"
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/topology.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/gmxassert.h"
62 /* DISCLAIMER: All the atom count and thread numbers below are heuristic.
63 * The real switching points will depend on the system simulation,
64 * the algorithms used and the hardware it's running on, as well as if there
65 * are other jobs running on the same machine. We try to take into account
66 * factors that have a large influence, such as recent Intel CPUs being
67 * much better at wide multi-threading. The remaining factors should
68 * (hopefully) have a small influence, such that the performance just before
69 * and after a switch point doesn't change too much.
72 static const bool bHasOmpSupport
= GMX_OPENMP
;
75 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
76 * the number of threads will get lowered.
78 static const int min_atoms_per_mpi_thread
= 90;
79 static const int min_atoms_per_gpu
= 900;
80 #endif /* GMX_THREAD_MPI */
82 /* TODO choose nthreads_omp based on hardware topology
83 when we have a hardware topology detection library */
84 /* First we consider the case of no MPI (1 MPI rank).
85 * In general, when running up to 8 threads, OpenMP should be faster.
86 * Note: on AMD Bulldozer we should avoid running OpenMP over two dies.
87 * On Intel>=Nehalem running OpenMP on a single CPU is always faster,
88 * even on two CPUs it's usually faster (but with many OpenMP threads
89 * it could be faster not to use HT, currently we always use HT).
90 * On Nehalem/Westmere we want to avoid running 16 threads over
91 * two CPUs with HT, so we need a limit<16; thus we use 12.
92 * A reasonable limit for Intel Sandy and Ivy bridge,
93 * not knowing the topology, is 16 threads.
94 * Below we check for Intel and AVX, which for now includes
95 * Sandy/Ivy Bridge, Has/Broadwell. By checking for AVX instead of
96 * model numbers we ensure also future Intel CPUs are covered.
98 const int nthreads_omp_faster_default
= 8;
99 const int nthreads_omp_faster_Nehalem
= 12;
100 const int nthreads_omp_faster_Intel_AVX
= 16;
101 /* For CPU only runs the fastest options are usually MPI or OpenMP only.
102 * With one GPU, using MPI only is almost never optimal, so we need to
103 * compare running pure OpenMP with combined MPI+OpenMP. This means higher
104 * OpenMP threads counts can still be ok. Multiplying the numbers above
105 * by a factor of 2 seems to be a good estimate.
107 const int nthreads_omp_faster_gpu_fac
= 2;
109 /* This is the case with MPI (2 or more MPI PP ranks).
110 * By default we will terminate with a fatal error when more than 8
111 * OpenMP thread are (indirectly) requested, since using less threads
112 * nearly always results in better performance.
113 * With thread-mpi and multiple GPUs or one GPU and too many threads
114 * we first try 6 OpenMP threads and then less until the number of MPI ranks
115 * is divisible by the number of GPUs.
117 #if GMX_OPENMP && GMX_MPI
118 const int nthreads_omp_mpi_ok_max
= 8;
119 const int nthreads_omp_mpi_ok_min_cpu
= 1;
121 const int nthreads_omp_mpi_ok_min_gpu
= 2;
122 const int nthreads_omp_mpi_target_max
= 6;
125 /* Returns the maximum OpenMP thread count for which using a single MPI rank
126 * should be faster than using multiple ranks with the same total thread count.
128 static int nthreads_omp_faster(const gmx::CpuInfo
&cpuInfo
, gmx_bool bUseGPU
)
132 if (cpuInfo
.vendor() == gmx::CpuInfo::Vendor::Intel
&&
133 cpuInfo
.feature(gmx::CpuInfo::Feature::X86_Avx
))
135 nth
= nthreads_omp_faster_Intel_AVX
;
137 else if (gmx::cpuIsX86Nehalem(cpuInfo
))
140 nth
= nthreads_omp_faster_Nehalem
;
144 nth
= nthreads_omp_faster_default
;
149 nth
*= nthreads_omp_faster_gpu_fac
;
152 nth
= std::min(nth
, GMX_OPENMP_MAX_THREADS
);
157 /* Returns that maximum OpenMP thread count that passes the efficiency check */
158 static int nthreads_omp_efficient_max(int gmx_unused nrank
,
159 const gmx::CpuInfo
&cpuInfo
,
162 #if GMX_OPENMP && GMX_MPI
165 return nthreads_omp_mpi_ok_max
;
170 return nthreads_omp_faster(cpuInfo
, bUseGPU
);
174 /* Return the number of thread-MPI ranks to use.
175 * This is chosen such that we can always obey our own efficiency checks.
177 static int get_tmpi_omp_thread_division(const gmx_hw_info_t
*hwinfo
,
178 const gmx_hw_opt_t
*hw_opt
,
183 const gmx::CpuInfo
&cpuInfo
= *hwinfo
->cpuInfo
;
185 GMX_RELEASE_ASSERT(nthreads_tot
> 0, "There must be at least one thread per rank");
187 /* There are no separate PME nodes here, as we ensured in
188 * check_and_update_hw_opt that nthreads_tmpi>0 with PME nodes
189 * and a conditional ensures we would not have ended up here.
190 * Note that separate PME nodes might be switched on later.
196 /* When the user sets nthreads_omp, we can end up oversubscribing CPU cores
197 * if we simply start as many ranks as GPUs. To avoid this, we start as few
198 * tMPI ranks as necessary to avoid oversubscription and instead leave GPUs idle.
199 * If the user does not set the number of OpenMP threads, nthreads_omp==0 and
200 * this code has no effect.
202 GMX_RELEASE_ASSERT(hw_opt
->nthreads_omp
>= 0, "nthreads_omp is negative, but previous checks should have prevented this");
203 while (nrank
*hw_opt
->nthreads_omp
> hwinfo
->nthreads_hw_avail
&& nrank
> 1)
208 if (nthreads_tot
< nrank
)
210 /* #thread < #gpu is very unlikely, but if so: waste gpu(s) */
211 nrank
= nthreads_tot
;
213 else if (gmx_gpu_sharing_supported() &&
214 (nthreads_tot
> nthreads_omp_faster(cpuInfo
, ngpu
> 0) ||
215 (ngpu
> 1 && nthreads_tot
/ngpu
> nthreads_omp_mpi_target_max
)))
217 /* The high OpenMP thread count will likely result in sub-optimal
218 * performance. Increase the rank count to reduce the thread count
219 * per rank. This will lead to GPU sharing by MPI ranks/threads.
223 /* Increase the rank count as long as have we more than 6 OpenMP
224 * threads per rank or the number of hardware threads is not
225 * divisible by the rank count. Don't go below 2 OpenMP threads.
233 while (nthreads_tot
/nrank
> nthreads_omp_mpi_target_max
||
234 (nthreads_tot
/(ngpu
*(nshare
+ 1)) >= nthreads_omp_mpi_ok_min_gpu
&& nthreads_tot
% nrank
!= 0));
237 else if (hw_opt
->nthreads_omp
> 0)
239 /* Here we could oversubscribe, when we do, we issue a warning later */
240 nrank
= std::max(1, nthreads_tot
/hw_opt
->nthreads_omp
);
244 if (nthreads_tot
<= nthreads_omp_faster(cpuInfo
, ngpu
> 0))
246 /* Use pure OpenMP parallelization */
251 /* Don't use OpenMP parallelization */
252 nrank
= nthreads_tot
;
260 static int getMaxGpuUsable(FILE *fplog
, const t_commrec
*cr
, const gmx_hw_info_t
*hwinfo
,
261 int cutoff_scheme
, gmx_bool bUseGpu
)
263 /* This code relies on the fact that GPU are not detected when GPU
264 * acceleration was disabled at run time by the user.
266 if (cutoff_scheme
== ecutsVERLET
&&
268 hwinfo
->gpu_info
.n_dev_compatible
> 0)
270 if (gmx_multiple_gpu_per_node_supported())
272 return hwinfo
->gpu_info
.n_dev_compatible
;
276 if (hwinfo
->gpu_info
.n_dev_compatible
> 1)
278 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");
293 gmxSmtIsEnabled(const gmx::HardwareTopology
&hwTop
)
295 return (hwTop
.supportLevel() >= gmx::HardwareTopology::SupportLevel::Basic
&& hwTop
.machine().sockets
[0].cores
[0].hwThreads
.size() > 1);
298 /* Get the number of MPI ranks to use for thread-MPI based on how many
299 * were requested, which algorithms we're using,
300 * and how many particles there are.
301 * At the point we have already called check_and_update_hw_opt.
302 * Thus all options should be internally consistent and consistent
303 * with the hardware, except that ntmpi could be larger than #GPU.
305 int get_nthreads_mpi(const gmx_hw_info_t
*hwinfo
,
306 gmx_hw_opt_t
*hw_opt
,
307 const t_inputrec
*inputrec
,
308 const gmx_mtop_t
*mtop
,
313 int nthreads_hw
, nthreads_tot_max
, nrank
, ngpu
;
314 int min_atoms_per_mpi_rank
;
316 const gmx::CpuInfo
&cpuInfo
= *hwinfo
->cpuInfo
;
317 const gmx::HardwareTopology
&hwTop
= *hwinfo
->hardwareTopology
;
319 /* Check if an algorithm does not support parallel simulation. */
320 if (inputrec
->eI
== eiLBFGS
||
321 inputrec
->coulombtype
== eelEWALD
)
323 md_print_warn(cr
, fplog
, "The integration or electrostatics algorithm doesn't support parallel runs. Using a single thread-MPI rank.\n");
324 if (hw_opt
->nthreads_tmpi
> 1)
326 gmx_fatal(FARGS
, "You asked for more than 1 thread-MPI rank, but an algorithm doesn't support that");
332 if (hw_opt
->nthreads_tmpi
> 0)
334 /* Trivial, return right away */
335 return hw_opt
->nthreads_tmpi
;
338 nthreads_hw
= hwinfo
->nthreads_hw_avail
;
340 if (nthreads_hw
<= 0)
342 /* This should normally not happen, but if it does, we handle it */
343 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");
346 /* How many total (#tMPI*#OpenMP) threads can we start? */
347 if (hw_opt
->nthreads_tot
> 0)
349 nthreads_tot_max
= hw_opt
->nthreads_tot
;
353 nthreads_tot_max
= nthreads_hw
;
356 ngpu
= getMaxGpuUsable(fplog
, cr
, hwinfo
, inputrec
->cutoff_scheme
, bUseGpu
);
358 if (inputrec
->cutoff_scheme
== ecutsGROUP
)
360 /* We checked this before, but it doesn't hurt to do it once more */
361 GMX_RELEASE_ASSERT(hw_opt
->nthreads_omp
== 1, "The group scheme only supports one OpenMP thread per rank");
365 get_tmpi_omp_thread_division(hwinfo
, hw_opt
, nthreads_tot_max
, ngpu
);
367 if (inputrec
->eI
== eiNM
|| EI_TPI(inputrec
->eI
))
369 /* Dims/steps are divided over the nodes iso splitting the atoms.
370 * With NM we can't have more ranks than #atoms*#dim. With TPI it's
371 * unlikely we have fewer atoms than ranks, and if so, communication
372 * would become a bottleneck, so we set the limit to 1 atom/rank.
374 min_atoms_per_mpi_rank
= 1;
380 min_atoms_per_mpi_rank
= min_atoms_per_gpu
;
384 min_atoms_per_mpi_rank
= min_atoms_per_mpi_thread
;
388 if (mtop
->natoms
/nrank
< min_atoms_per_mpi_rank
)
392 /* the rank number was chosen automatically, but there are too few
393 atoms per rank, so we need to reduce the rank count */
394 nrank_new
= std::max(1, mtop
->natoms
/min_atoms_per_mpi_rank
);
396 /* Avoid partial use of Hyper-Threading */
397 if (gmxSmtIsEnabled(hwTop
) &&
398 nrank_new
> nthreads_hw
/2 && nrank_new
< nthreads_hw
)
400 nrank_new
= nthreads_hw
/2;
403 /* If the user specified the total thread count, ensure this is
404 * divisible by the number of ranks.
405 * It is quite likely that we have too many total threads compared
406 * to the size of the system, but if the user asked for this many
407 * threads we should respect that.
409 while (hw_opt
->nthreads_tot
> 0 &&
410 hw_opt
->nthreads_tot
% nrank_new
!= 0)
415 /* Avoid large prime numbers in the rank count */
418 /* Use only 6,8,10 with additional factors of 2 */
422 while (3*fac
*2 <= nrank_new
)
427 nrank_new
= (nrank_new
/fac
)*fac
;
431 /* Avoid 5, since small system won't fit 5 domains along
432 * a dimension. This might lead to waisting some cores, but this
433 * will have a small impact in this regime of very small systems.
443 /* We reduced the number of tMPI ranks, which means we might violate
444 * our own efficiency checks if we simply use all hardware threads.
446 if (bHasOmpSupport
&& hw_opt
->nthreads_omp
<= 0 && hw_opt
->nthreads_tot
<= 0)
448 /* The user set neither the total nor the OpenMP thread count,
449 * we should use all hardware threads, unless we will violate
450 * our own efficiency limitation on the thread count.
454 nt_omp_max
= nthreads_omp_efficient_max(nrank
, cpuInfo
, ngpu
>= 1);
456 if (nrank
*nt_omp_max
< hwinfo
->nthreads_hw_avail
)
458 /* Limit the number of OpenMP threads to start */
459 hw_opt
->nthreads_omp
= nt_omp_max
;
463 fprintf(stderr
, "\n");
464 fprintf(stderr
, "NOTE: Parallelization is limited by the small number of atoms,\n");
465 fprintf(stderr
, " only starting %d thread-MPI ranks.\n", nrank
);
466 fprintf(stderr
, " You can use the -nt and/or -ntmpi option to optimize the number of threads.\n\n");
471 #endif /* GMX_THREAD_MPI */
474 void check_resource_division_efficiency(const gmx_hw_info_t
*hwinfo
,
475 const gmx_hw_opt_t
*hw_opt
,
476 gmx_bool bNtOmpOptionSet
,
480 #if GMX_OPENMP && GMX_MPI
481 int nth_omp_min
, nth_omp_max
, ngpu
;
484 const char *mpi_option
= " (option -ntmpi)";
486 const char *mpi_option
= "";
489 /* This function should be called after thread-MPI (when configured) and
490 * OpenMP have been initialized. Check that here.
493 GMX_RELEASE_ASSERT(nthreads_omp_faster_default
>= nthreads_omp_mpi_ok_max
, "Inconsistent OpenMP thread count default values");
494 GMX_RELEASE_ASSERT(hw_opt
->nthreads_tmpi
>= 1, "Must have at least one thread-MPI rank");
496 GMX_RELEASE_ASSERT(gmx_omp_nthreads_get(emntDefault
) >= 1, "Must have at least one OpenMP thread");
498 nth_omp_min
= gmx_omp_nthreads_get(emntDefault
);
499 nth_omp_max
= gmx_omp_nthreads_get(emntDefault
);
500 ngpu
= hw_opt
->gpu_opt
.n_dev_use
;
502 /* Thread-MPI seems to have a bug with reduce on 1 node, so use a cond. */
503 if (cr
->nnodes
+ cr
->npmenodes
> 1)
505 int count
[3], count_max
[3];
507 count
[0] = -nth_omp_min
;
508 count
[1] = nth_omp_max
;
511 MPI_Allreduce(count
, count_max
, 3, MPI_INT
, MPI_MAX
, cr
->mpi_comm_mysim
);
513 /* In case of an inhomogeneous run setup we use the maximum counts */
514 nth_omp_min
= -count_max
[0];
515 nth_omp_max
= count_max
[1];
519 int nthreads_omp_mpi_ok_min
;
523 nthreads_omp_mpi_ok_min
= nthreads_omp_mpi_ok_min_cpu
;
527 /* With GPUs we set the minimum number of OpenMP threads to 2 to catch
528 * cases where the user specifies #ranks == #cores.
530 nthreads_omp_mpi_ok_min
= nthreads_omp_mpi_ok_min_gpu
;
533 if (DOMAINDECOMP(cr
) && cr
->nnodes
> 1)
535 if (nth_omp_max
< nthreads_omp_mpi_ok_min
||
536 (!(ngpu
> 0 && !gmx_gpu_sharing_supported()) &&
537 nth_omp_max
> nthreads_omp_mpi_ok_max
))
539 /* Note that we print target_max here, not ok_max */
540 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.",
542 nthreads_omp_mpi_ok_min
,
543 nthreads_omp_mpi_target_max
);
547 md_print_warn(cr
, fplog
, "NOTE: %s\n", buf
);
551 /* This fatal error, and the one below, is nasty, but it's
552 * probably the only way to ensure that all users don't waste
553 * a lot of resources, since many users don't read logs/stderr.
555 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
);
561 const gmx::CpuInfo
&cpuInfo
= *hwinfo
->cpuInfo
;
563 /* No domain decomposition (or only one domain) */
564 if (!(ngpu
> 0 && !gmx_gpu_sharing_supported()) &&
565 nth_omp_max
> nthreads_omp_faster(cpuInfo
, ngpu
> 0))
567 /* To arrive here, the user/system set #ranks and/or #OMPthreads */
571 bEnvSet
= (getenv("OMP_NUM_THREADS") != NULL
);
573 if (bNtOmpOptionSet
|| bEnvSet
)
575 sprintf(buf2
, "You requested %d OpenMP threads", nth_omp_max
);
579 sprintf(buf2
, "Your choice of %d MPI rank%s and the use of %d total threads %sleads to the use of %d OpenMP threads",
580 cr
->nnodes
+ cr
->npmenodes
,
581 cr
->nnodes
+ cr
->npmenodes
== 1 ? "" : "s",
582 hw_opt
->nthreads_tot
> 0 ? hw_opt
->nthreads_tot
: hwinfo
->nthreads_hw_avail
,
583 hwinfo
->nphysicalnode
> 1 ? "on a node " : "",
586 sprintf(buf
, "%s, whereas we expect the optimum to be with more MPI ranks with %d to %d OpenMP threads.",
587 buf2
, nthreads_omp_mpi_ok_min
, nthreads_omp_mpi_target_max
);
589 /* We can not quit with a fatal error when OMP_NUM_THREADS is set
590 * with different values per rank or node, since in that case
591 * the user can not set -ntomp to override the error.
593 if (bNtOmpOptionSet
|| (bEnvSet
&& nth_omp_min
!= nth_omp_max
))
595 md_print_warn(cr
, fplog
, "NOTE: %s\n", buf
);
599 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
);
603 #else /* GMX_OPENMP && GMX_MPI */
604 /* No OpenMP and/or MPI: it doesn't make much sense to check */
605 GMX_UNUSED_VALUE(hw_opt
);
606 GMX_UNUSED_VALUE(bNtOmpOptionSet
);
607 /* Check if we have more than 1 physical core, if detected,
608 * or more than 1 hardware thread if physical cores were not detected.
610 #if !GMX_OPENMP && !GMX_MPI
611 if (hwinfo
->hardwareTopology
->numberOfCores() > 1)
613 md_print_warn(cr
, fplog
, "NOTE: GROMACS was compiled without OpenMP and (thread-)MPI support, can only use a single CPU core\n");
616 GMX_UNUSED_VALUE(hwinfo
);
617 GMX_UNUSED_VALUE(cr
);
618 GMX_UNUSED_VALUE(fplog
);
621 #endif /* GMX_OPENMP && GMX_MPI */
625 static void print_hw_opt(FILE *fp
, const gmx_hw_opt_t
*hw_opt
)
627 fprintf(fp
, "hw_opt: nt %d ntmpi %d ntomp %d ntomp_pme %d gpu_id '%s'\n",
628 hw_opt
->nthreads_tot
,
629 hw_opt
->nthreads_tmpi
,
630 hw_opt
->nthreads_omp
,
631 hw_opt
->nthreads_omp_pme
,
632 hw_opt
->gpu_opt
.gpu_id
!= NULL
? hw_opt
->gpu_opt
.gpu_id
: "");
635 /* Checks we can do when we don't (yet) know the cut-off scheme */
636 void check_and_update_hw_opt_1(gmx_hw_opt_t
*hw_opt
,
640 /* Currently hw_opt only contains default settings or settings supplied
641 * by the user on the command line.
643 if (hw_opt
->nthreads_omp
< 0)
645 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
);
648 /* Check for OpenMP settings stored in environment variables, which can
649 * potentially be different on different MPI ranks.
651 gmx_omp_nthreads_read_env(&hw_opt
->nthreads_omp
, SIMMASTER(cr
));
653 /* Check restrictions on the user supplied options before modifying them.
654 * TODO: Put the user values in a const struct and preserve them.
657 if (hw_opt
->nthreads_tot
> 0)
659 gmx_fatal(FARGS
, "Setting the total number of threads is only supported with thread-MPI and GROMACS was compiled without thread-MPI");
661 if (hw_opt
->nthreads_tmpi
> 0)
663 gmx_fatal(FARGS
, "Setting the number of thread-MPI ranks is only supported with thread-MPI and GROMACS was compiled without thread-MPI");
669 /* Check restrictions on PME thread related options set by the user */
671 if (hw_opt
->nthreads_omp_pme
> 0 && hw_opt
->nthreads_omp
<= 0)
673 gmx_fatal(FARGS
, "You need to specify -ntomp in addition to -ntomp_pme");
676 if (hw_opt
->nthreads_omp_pme
>= 1 &&
677 hw_opt
->nthreads_omp_pme
!= hw_opt
->nthreads_omp
&&
680 /* This can result in a fatal error on many MPI ranks,
681 * but since the thread count can differ per rank,
682 * we can't easily avoid this.
684 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");
689 /* GROMACS was configured without OpenMP support */
691 if (hw_opt
->nthreads_omp
> 1 || hw_opt
->nthreads_omp_pme
> 1)
693 gmx_fatal(FARGS
, "More than 1 OpenMP thread requested, but GROMACS was compiled without OpenMP support");
695 hw_opt
->nthreads_omp
= 1;
696 hw_opt
->nthreads_omp_pme
= 1;
699 if (hw_opt
->nthreads_tot
> 0 && hw_opt
->nthreads_omp_pme
<= 0)
701 /* We have the same number of OpenMP threads for PP and PME ranks,
702 * thus we can perform several consistency checks.
704 if (hw_opt
->nthreads_tmpi
> 0 &&
705 hw_opt
->nthreads_omp
> 0 &&
706 hw_opt
->nthreads_tot
!= hw_opt
->nthreads_tmpi
*hw_opt
->nthreads_omp
)
708 gmx_fatal(FARGS
, "The total number of threads requested (%d) does not match the thread-MPI ranks (%d) times the OpenMP threads (%d) requested",
709 hw_opt
->nthreads_tot
, hw_opt
->nthreads_tmpi
, hw_opt
->nthreads_omp
);
712 if (hw_opt
->nthreads_tmpi
> 0 &&
713 hw_opt
->nthreads_tot
% hw_opt
->nthreads_tmpi
!= 0)
715 gmx_fatal(FARGS
, "The total number of threads requested (%d) is not divisible by the number of thread-MPI ranks requested (%d)",
716 hw_opt
->nthreads_tot
, hw_opt
->nthreads_tmpi
);
719 if (hw_opt
->nthreads_omp
> 0 &&
720 hw_opt
->nthreads_tot
% hw_opt
->nthreads_omp
!= 0)
722 gmx_fatal(FARGS
, "The total number of threads requested (%d) is not divisible by the number of OpenMP threads requested (%d)",
723 hw_opt
->nthreads_tot
, hw_opt
->nthreads_omp
);
727 if (hw_opt
->nthreads_tot
== 1)
729 hw_opt
->nthreads_tmpi
= 1;
731 if (hw_opt
->nthreads_omp
> 1)
733 gmx_fatal(FARGS
, "You requested %d OpenMP threads with %d total threads",
734 hw_opt
->nthreads_tmpi
, hw_opt
->nthreads_tot
);
736 hw_opt
->nthreads_omp
= 1;
737 hw_opt
->nthreads_omp_pme
= 1;
740 /* Parse GPU IDs, if provided.
741 * We check consistency with the tMPI thread count later.
743 gmx_parse_gpu_ids(&hw_opt
->gpu_opt
);
746 if (hw_opt
->gpu_opt
.n_dev_use
> 0 && hw_opt
->nthreads_tmpi
== 0)
748 /* Set the number of MPI threads equal to the number of GPUs */
749 hw_opt
->nthreads_tmpi
= hw_opt
->gpu_opt
.n_dev_use
;
751 if (hw_opt
->nthreads_tot
> 0 &&
752 hw_opt
->nthreads_tmpi
> hw_opt
->nthreads_tot
)
754 /* We have more GPUs than total threads requested.
755 * We choose to (later) generate a mismatch error,
756 * instead of launching more threads than requested.
758 hw_opt
->nthreads_tmpi
= hw_opt
->nthreads_tot
;
765 print_hw_opt(debug
, hw_opt
);
768 /* Asserting this simplifies the hardware resource division later
770 GMX_RELEASE_ASSERT(!(hw_opt
->nthreads_omp_pme
>= 1 && hw_opt
->nthreads_omp
<= 0),
771 "PME thread count should only be set when the normal thread count is also set");
774 /* Checks we can do when we know the cut-off scheme */
775 void check_and_update_hw_opt_2(gmx_hw_opt_t
*hw_opt
,
778 if (cutoff_scheme
== ecutsGROUP
)
780 /* We only have OpenMP support for PME only nodes */
781 if (hw_opt
->nthreads_omp
> 1)
783 gmx_fatal(FARGS
, "OpenMP threads have been requested with cut-off scheme %s, but these are only supported with cut-off scheme %s",
784 ecutscheme_names
[cutoff_scheme
],
785 ecutscheme_names
[ecutsVERLET
]);
787 hw_opt
->nthreads_omp
= 1;
791 /* Checks we can do when we know the thread-MPI rank count */
792 void check_and_update_hw_opt_3(gmx_hw_opt_t gmx_unused
*hw_opt
)
795 GMX_RELEASE_ASSERT(hw_opt
->nthreads_tmpi
>= 1, "Must have at least one thread-MPI rank");
797 /* If the user set the total number of threads on the command line
798 * and did not specify the number of OpenMP threads, set the latter here.
800 if (hw_opt
->nthreads_tot
> 0 && hw_opt
->nthreads_omp
<= 0)
802 hw_opt
->nthreads_omp
= hw_opt
->nthreads_tot
/hw_opt
->nthreads_tmpi
;
804 if (!bHasOmpSupport
&& hw_opt
->nthreads_omp
> 1)
806 gmx_fatal(FARGS
, "You (indirectly) asked for OpenMP threads by setting -nt > -ntmpi, but GROMACS was compiled without OpenMP support");
811 GMX_RELEASE_ASSERT(bHasOmpSupport
|| hw_opt
->nthreads_omp
== 1, "Without OpenMP support, only one thread per rank can be used");
813 /* We are done with updating nthreads_omp, we can set nthreads_omp_pme */
814 if (hw_opt
->nthreads_omp_pme
<= 0 && hw_opt
->nthreads_omp
> 0)
816 hw_opt
->nthreads_omp_pme
= hw_opt
->nthreads_omp
;
821 print_hw_opt(debug
, hw_opt
);