Split lines with many copyright years
[gromacs.git] / src / gromacs / mdlib / gmx_omp_nthreads.cpp
blobaa3327c8221a9703360160b55b42cdc1731ba4c4
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5 * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "gmxpre.h"
39 #include "gmx_omp_nthreads.h"
41 #include "config.h"
43 #include <cstdio>
44 #include <cstdlib>
45 #include <cstring>
47 #include "gromacs/gmxlib/network.h"
48 #include "gromacs/mdtypes/commrec.h"
49 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/gmxassert.h"
52 #include "gromacs/utility/gmxomp.h"
53 #include "gromacs/utility/logger.h"
54 #include "gromacs/utility/programcontext.h"
56 /** Structure with the number of threads for each OpenMP multi-threaded
57 * algorithmic module in mdrun. */
58 typedef struct
60 int gnth; /**< Global num. of threads per PP or PP+PME process/tMPI thread. */
61 int gnth_pme; /**< Global num. of threads per PME only process/tMPI thread. */
63 int nth[emntNR]; /**< Number of threads for each module, indexed with module_nth_t */
64 gmx_bool initialized; /**< TRUE if the module as been initialized. */
65 } omp_module_nthreads_t;
67 /** Names of environment variables to set the per module number of threads.
69 * Indexed with the values of module_nth_t.
70 * */
71 static const char* modth_env_var[emntNR] = { "GMX_DEFAULT_NUM_THREADS should never be set",
72 "GMX_DOMDEC_NUM_THREADS",
73 "GMX_PAIRSEARCH_NUM_THREADS",
74 "GMX_NONBONDED_NUM_THREADS",
75 "GMX_LISTED_FORCES_NUM_THREADS",
76 "GMX_PME_NUM_THREADS",
77 "GMX_UPDATE_NUM_THREADS",
78 "GMX_VSITE_NUM_THREADS",
79 "GMX_LINCS_NUM_THREADS",
80 "GMX_SETTLE_NUM_THREADS" };
82 /** Names of the modules. */
83 static const char* mod_name[emntNR] = { "default", "domain decomposition",
84 "pair search", "non-bonded",
85 "bonded", "PME",
86 "update", "LINCS",
87 "SETTLE" };
89 /** Number of threads for each algorithmic module.
91 * File-scope global variable that gets set once in pick_module_nthreads()
92 * and queried via gmx_omp_nthreads_get().
94 * All fields are initialized to 0 which should result in errors if
95 * the init call is omitted.
96 * */
97 static omp_module_nthreads_t modth = { 0, 0, { 0, 0, 0, 0, 0, 0, 0, 0, 0 }, FALSE };
100 /** Determine the number of threads for module \p mod.
102 * \p m takes values form the module_nth_t enum and maps these to the
103 * corresponding value in modth_env_var.
105 * Each number of threads per module takes the default value unless
106 * GMX_*_NUM_THERADS env var is set, case in which its value overrides
107 * the default.
109 static void pick_module_nthreads(const gmx::MDLogger& mdlog, int m, gmx_bool bSepPME)
111 char* env;
112 int nth;
114 const bool bOMP = GMX_OPENMP;
116 /* The default should never be set through a GMX_*_NUM_THREADS env var
117 * as it's always equal with gnth. */
118 if (m == emntDefault)
120 return;
123 /* check the environment variable */
124 if ((env = getenv(modth_env_var[m])) != nullptr)
126 sscanf(env, "%d", &nth);
128 if (!bOMP)
130 gmx_warning("%s=%d is set, but %s is compiled without OpenMP!", modth_env_var[m], nth,
131 gmx::getProgramContext().displayName());
134 /* with the verlet codepath, when any GMX_*_NUM_THREADS env var is set,
135 * OMP_NUM_THREADS also has to be set */
136 if (getenv("OMP_NUM_THREADS") == nullptr)
138 gmx_warning(
139 "%s=%d is set, the default number of threads also "
140 "needs to be set with OMP_NUM_THREADS!",
141 modth_env_var[m], nth);
144 /* only babble if we are really overriding with a different value */
145 if ((bSepPME && m == emntPME && nth != modth.gnth_pme) || (nth != modth.gnth))
147 GMX_LOG(mdlog.warning)
148 .asParagraph()
149 .appendTextFormatted("%s=%d set, overriding the default number of %s threads",
150 modth_env_var[m], nth, mod_name[m]);
153 else
155 /* pick the global PME node nthreads if we are setting the number
156 * of threads in separate PME nodes */
157 nth = (bSepPME && m == emntPME) ? modth.gnth_pme : modth.gnth;
160 gmx_omp_nthreads_set(m, nth);
163 void gmx_omp_nthreads_read_env(const gmx::MDLogger& mdlog, int* nthreads_omp)
165 char* env;
166 gmx_bool bCommandLineSetNthreadsOMP = *nthreads_omp > 0;
167 char buffer[STRLEN];
169 GMX_RELEASE_ASSERT(nthreads_omp, "nthreads_omp must be a non-NULL pointer");
171 if ((env = getenv("OMP_NUM_THREADS")) != nullptr)
173 int nt_omp;
175 sscanf(env, "%d", &nt_omp);
176 if (nt_omp <= 0)
178 gmx_fatal(FARGS, "OMP_NUM_THREADS is invalid: '%s'", env);
181 if (bCommandLineSetNthreadsOMP && nt_omp != *nthreads_omp)
183 gmx_fatal(FARGS,
184 "Environment variable OMP_NUM_THREADS (%d) and the number of threads "
185 "requested on the command line (%d) have different values. Either omit one, "
186 "or set them both to the same value.",
187 nt_omp, *nthreads_omp);
190 /* Setting the number of OpenMP threads. */
191 *nthreads_omp = nt_omp;
193 /* Output the results */
194 sprintf(buffer,
195 "\nThe number of OpenMP threads was set by environment variable OMP_NUM_THREADS to "
196 "%d%s\n\n",
197 nt_omp,
198 bCommandLineSetNthreadsOMP ? " (and the command-line setting agreed with that)" : "");
200 /* This prints once per simulation for multi-simulations,
201 * which might help diagnose issues with inhomogenous
202 * cluster setups. */
203 GMX_LOG(mdlog.info).appendTextFormatted("%s", buffer);
204 if (debug)
206 /* This prints once per process for real MPI (i.e. once
207 * per debug file), and once per simulation for thread MPI
208 * (because of logic in the calling function). */
209 fputs(buffer, debug);
214 /*! \brief Helper function for parsing various input about the number
215 of OpenMP threads to use in various modules and deciding what to
216 do about it. */
217 static void manage_number_of_openmp_threads(const gmx::MDLogger& mdlog,
218 const t_commrec* cr,
219 bool bOMP,
220 int nthreads_hw_avail,
221 int omp_nthreads_req,
222 int omp_nthreads_pme_req,
223 gmx_bool gmx_unused bThisNodePMEOnly,
224 int numRanksOnThisNode,
225 gmx_bool bSepPME)
227 int nth;
228 char* env;
230 #if GMX_THREAD_MPI
231 /* modth is shared among tMPI threads, so for thread safety, the
232 * detection is done on the master only. It is not thread-safe
233 * with multiple simulations, but that's anyway not supported by
234 * tMPI. */
235 if (!SIMMASTER(cr))
237 return;
239 #else
240 GMX_UNUSED_VALUE(cr);
241 #endif
243 if (modth.initialized)
245 /* Just return if the initialization has already been
246 done. This could only happen if gmx_omp_nthreads_init() has
247 already been called. */
248 return;
251 /* With full OpenMP support (verlet scheme) set the number of threads
252 * per process / default:
253 * - 1 if not compiled with OpenMP or
254 * - OMP_NUM_THREADS if the env. var is set, or
255 * - omp_nthreads_req = #of threads requested by the user on the mdrun
256 * command line, otherwise
257 * - take the max number of available threads and distribute them
258 * on the processes/tMPI threads.
259 * ~ The GMX_*_NUM_THREADS env var overrides the number of threads of
260 * the respective module and it has to be used in conjunction with
261 * OMP_NUM_THREADS.
263 * The number of PME threads is equal to:
264 * - 1 if not compiled with OpenMP or
265 * - GMX_PME_NUM_THREADS if defined, otherwise
266 * - OMP_NUM_THREADS if defined, otherwise
267 * - 1
269 nth = 1;
270 if ((env = getenv("OMP_NUM_THREADS")) != nullptr)
272 if (!bOMP && (std::strncmp(env, "1", 1) != 0))
274 gmx_warning("OMP_NUM_THREADS is set, but %s was compiled without OpenMP support!",
275 gmx::getProgramContext().displayName());
277 else
279 nth = gmx_omp_get_max_threads();
282 else if (omp_nthreads_req > 0)
284 nth = omp_nthreads_req;
286 else if (bOMP)
288 /* max available threads per node */
289 nth = nthreads_hw_avail;
291 /* divide the threads among the MPI ranks */
292 if (nth >= numRanksOnThisNode)
294 nth /= numRanksOnThisNode;
296 else
298 nth = 1;
302 /* now we have the global values, set them:
303 * - 1 if not compiled with OpenMP
304 * - nth for the verlet scheme when compiled with OpenMP
306 if (bOMP)
308 modth.gnth = nth;
310 else
312 modth.gnth = 1;
315 if (bSepPME)
317 if (omp_nthreads_pme_req > 0)
319 modth.gnth_pme = omp_nthreads_pme_req;
321 else
323 modth.gnth_pme = nth;
326 else
328 modth.gnth_pme = 0;
331 /* now set the per-module values */
332 modth.nth[emntDefault] = modth.gnth;
333 pick_module_nthreads(mdlog, emntDomdec, bSepPME);
334 pick_module_nthreads(mdlog, emntPairsearch, bSepPME);
335 pick_module_nthreads(mdlog, emntNonbonded, bSepPME);
336 pick_module_nthreads(mdlog, emntBonded, bSepPME);
337 pick_module_nthreads(mdlog, emntPME, bSepPME);
338 pick_module_nthreads(mdlog, emntUpdate, bSepPME);
339 pick_module_nthreads(mdlog, emntVSITE, bSepPME);
340 pick_module_nthreads(mdlog, emntLINCS, bSepPME);
341 pick_module_nthreads(mdlog, emntSETTLE, bSepPME);
343 /* set the number of threads globally */
344 if (bOMP)
346 #if !GMX_THREAD_MPI
347 if (bThisNodePMEOnly)
349 gmx_omp_set_num_threads(modth.gnth_pme);
351 else
352 #endif /* GMX_THREAD_MPI */
354 gmx_omp_set_num_threads(nth);
358 modth.initialized = TRUE;
361 /*! \brief Report on the OpenMP settings that will be used */
362 static void reportOpenmpSettings(const gmx::MDLogger& mdlog, const t_commrec* cr, gmx_bool bOMP, gmx_bool bSepPME)
364 #if GMX_THREAD_MPI
365 const char* mpi_str = "per tMPI thread";
366 #else
367 const char* mpi_str = "per MPI process";
368 #endif
369 int nth_min, nth_max, nth_pme_min, nth_pme_max;
371 /* inform the user about the settings */
372 if (!bOMP)
374 return;
377 #if GMX_MPI
378 if (cr->nnodes > 1)
380 /* Get the min and max thread counts over the MPI ranks */
381 int buf_in[4], buf_out[4];
383 buf_in[0] = -modth.gnth;
384 buf_in[1] = modth.gnth;
385 buf_in[2] = -modth.gnth_pme;
386 buf_in[3] = modth.gnth_pme;
388 MPI_Allreduce(buf_in, buf_out, 4, MPI_INT, MPI_MAX, cr->mpi_comm_mysim);
390 nth_min = -buf_out[0];
391 nth_max = buf_out[1];
392 nth_pme_min = -buf_out[2];
393 nth_pme_max = buf_out[3];
395 else
396 #endif
398 nth_min = modth.gnth;
399 nth_max = modth.gnth;
400 nth_pme_min = modth.gnth_pme;
401 nth_pme_max = modth.gnth_pme;
405 if (nth_max == nth_min)
407 GMX_LOG(mdlog.warning)
408 .appendTextFormatted("Using %d OpenMP thread%s %s", nth_min, nth_min > 1 ? "s" : "",
409 cr->nnodes > 1 ? mpi_str : "");
411 else
413 GMX_LOG(mdlog.warning).appendTextFormatted("Using %d - %d OpenMP threads %s", nth_min, nth_max, mpi_str);
416 if (bSepPME && (nth_pme_min != nth_min || nth_pme_max != nth_max))
418 if (nth_pme_max == nth_pme_min)
420 GMX_LOG(mdlog.warning)
421 .appendTextFormatted("Using %d OpenMP thread%s %s for PME", nth_pme_min,
422 nth_pme_min > 1 ? "s" : "", cr->nnodes > 1 ? mpi_str : "");
424 else
426 GMX_LOG(mdlog.warning)
427 .appendTextFormatted("Using %d - %d OpenMP threads %s for PME", nth_pme_min,
428 nth_pme_max, mpi_str);
431 GMX_LOG(mdlog.warning);
434 void gmx_omp_nthreads_init(const gmx::MDLogger& mdlog,
435 t_commrec* cr,
436 int nthreads_hw_avail,
437 int numRanksOnThisNode,
438 int omp_nthreads_req,
439 int omp_nthreads_pme_req,
440 gmx_bool bThisNodePMEOnly)
442 gmx_bool bSepPME;
444 const bool bOMP = GMX_OPENMP;
446 bSepPME = (thisRankHasDuty(cr, DUTY_PP) != thisRankHasDuty(cr, DUTY_PME));
448 manage_number_of_openmp_threads(mdlog, cr, bOMP, nthreads_hw_avail, omp_nthreads_req,
449 omp_nthreads_pme_req, bThisNodePMEOnly, numRanksOnThisNode, bSepPME);
450 #if GMX_THREAD_MPI
451 /* Non-master threads have to wait for the OpenMP management to be
452 * done, so that code elsewhere that uses OpenMP can be certain
453 * the setup is complete. */
454 if (PAR(cr))
456 MPI_Barrier(cr->mpi_comm_mysim);
458 #endif
460 reportOpenmpSettings(mdlog, cr, bOMP, bSepPME);
463 int gmx_omp_nthreads_get(int mod)
465 if (mod < 0 || mod >= emntNR)
467 /* invalid module queried */
468 return -1;
470 else
472 return modth.nth[mod];
476 void gmx_omp_nthreads_set(int mod, int nthreads)
478 /* Catch an attempt to set the number of threads on an invalid
479 * OpenMP module. */
480 GMX_RELEASE_ASSERT(mod >= 0 && mod < emntNR, "Trying to set nthreads on invalid OpenMP module");
482 modth.nth[mod] = nthreads;