Prepared t_mdatoms for using vector
[gromacs.git] / src / gromacs / mdlib / gmx_omp_nthreads.cpp
blob6c3b06b0992affb830af5c6d7c06184739cadbdd
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,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.
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.
36 #include "gmxpre.h"
38 #include "gmx_omp_nthreads.h"
40 #include "config.h"
42 #include <cstdio>
43 #include <cstdlib>
44 #include <cstring>
46 #include "gromacs/gmxlib/network.h"
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/utility/cstringutil.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/gmxassert.h"
51 #include "gromacs/utility/gmxomp.h"
52 #include "gromacs/utility/logger.h"
53 #include "gromacs/utility/programcontext.h"
55 /** Structure with the number of threads for each OpenMP multi-threaded
56 * algorithmic module in mdrun. */
57 typedef struct
59 int gnth; /**< Global num. of threads per PP or PP+PME process/tMPI thread. */
60 int gnth_pme; /**< Global num. of threads per PME only process/tMPI thread. */
62 int nth[emntNR]; /**< Number of threads for each module, indexed with module_nth_t */
63 gmx_bool initialized; /**< TRUE if the module as been initialized. */
64 } omp_module_nthreads_t;
66 /** Names of environment variables to set the per module number of threads.
68 * Indexed with the values of module_nth_t.
69 * */
70 static const char *modth_env_var[emntNR] =
72 "GMX_DEFAULT_NUM_THREADS should never be set",
73 "GMX_DOMDEC_NUM_THREADS", "GMX_PAIRSEARCH_NUM_THREADS",
74 "GMX_NONBONDED_NUM_THREADS", "GMX_LISTED_FORCES_NUM_THREADS",
75 "GMX_PME_NUM_THREADS", "GMX_UPDATE_NUM_THREADS",
76 "GMX_VSITE_NUM_THREADS",
77 "GMX_LINCS_NUM_THREADS", "GMX_SETTLE_NUM_THREADS"
80 /** Names of the modules. */
81 static const char *mod_name[emntNR] =
83 "default", "domain decomposition", "pair search", "non-bonded",
84 "bonded", "PME", "update", "LINCS", "SETTLE"
87 /** Number of threads for each algorithmic module.
89 * File-scope global variable that gets set once in pick_module_nthreads()
90 * and queried via gmx_omp_nthreads_get().
92 * All fields are initialized to 0 which should result in errors if
93 * the init call is omitted.
94 * */
95 static omp_module_nthreads_t modth = { 0, 0, {0, 0, 0, 0, 0, 0, 0, 0, 0}, FALSE};
98 /** Determine the number of threads for module \p mod.
100 * \p m takes values form the module_nth_t enum and maps these to the
101 * corresponding value in modth_env_var.
103 * Each number of threads per module takes the default value unless
104 * GMX_*_NUM_THERADS env var is set, case in which its value overrides
105 * the deafult.
107 * The "group" scheme supports OpenMP only in PME and in thise case all but
108 * the PME nthread values default to 1.
110 static void pick_module_nthreads(const gmx::MDLogger &mdlog, int m,
111 gmx_bool bFullOmpSupport,
112 gmx_bool bSepPME)
114 char *env;
115 int nth;
117 const bool bOMP = GMX_OPENMP;
119 /* The default should never be set through a GMX_*_NUM_THREADS env var
120 * as it's always equal with gnth. */
121 if (m == emntDefault)
123 return;
126 /* check the environment variable */
127 if ((env = getenv(modth_env_var[m])) != nullptr)
129 sscanf(env, "%d", &nth);
131 // cppcheck-suppress knownConditionTrueFalse
132 if (!bOMP)
134 gmx_warning("%s=%d is set, but %s is compiled without OpenMP!",
135 modth_env_var[m], nth,
136 gmx::getProgramContext().displayName());
139 /* with the verlet codepath, when any GMX_*_NUM_THREADS env var is set,
140 * OMP_NUM_THREADS also has to be set */
141 if (bFullOmpSupport && getenv("OMP_NUM_THREADS") == nullptr)
143 gmx_warning("%s=%d is set, the default number of threads also "
144 "needs to be set with OMP_NUM_THREADS!",
145 modth_env_var[m], nth);
148 /* with the group scheme warn if any env var except PME is set */
149 if (!bFullOmpSupport)
151 if (m != emntPME)
153 gmx_warning("%s=%d is set, but OpenMP multithreading is not "
154 "supported in %s!",
155 modth_env_var[m], nth, mod_name[m]);
156 nth = 1;
160 /* only babble if we are really overriding with a different value */
161 if ((bSepPME && m == emntPME && nth != modth.gnth_pme) || (nth != modth.gnth))
163 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
164 "%s=%d set, overriding the default number of %s threads",
165 modth_env_var[m], nth, mod_name[m]);
168 else
170 /* pick the global PME node nthreads if we are setting the number
171 * of threads in separate PME nodes */
172 nth = (bSepPME && m == emntPME) ? modth.gnth_pme : modth.gnth;
175 gmx_omp_nthreads_set(m, nth);
178 void gmx_omp_nthreads_read_env(int *nthreads_omp,
179 gmx_bool bIsSimMaster)
181 char *env;
182 gmx_bool bCommandLineSetNthreadsOMP = *nthreads_omp > 0;
183 char buffer[STRLEN];
185 GMX_RELEASE_ASSERT(nthreads_omp, "nthreads_omp must be a non-NULL pointer");
187 if ((env = getenv("OMP_NUM_THREADS")) != nullptr)
189 int nt_omp;
191 sscanf(env, "%d", &nt_omp);
192 if (nt_omp <= 0)
194 gmx_fatal(FARGS, "OMP_NUM_THREADS is invalid: '%s'", env);
197 if (bCommandLineSetNthreadsOMP && nt_omp != *nthreads_omp)
199 gmx_fatal(FARGS, "Environment variable OMP_NUM_THREADS (%d) and the number of threads requested on the command line (%d) have different values. Either omit one, or set them both to the same value.", nt_omp, *nthreads_omp);
202 /* Setting the number of OpenMP threads. */
203 *nthreads_omp = nt_omp;
205 /* Output the results */
206 sprintf(buffer,
207 "The number of OpenMP threads was set by environment variable OMP_NUM_THREADS to %d%s\n",
208 nt_omp,
209 bCommandLineSetNthreadsOMP ? " (and the command-line setting agreed with that)" : "");
210 if (bIsSimMaster)
212 /* This prints once per simulation for multi-simulations,
213 * which might help diagnose issues with inhomogenous
214 * cluster setups. */
215 fputs(buffer, stderr);
217 if (debug)
219 /* This prints once per process for real MPI (i.e. once
220 * per debug file), and once per simulation for thread MPI
221 * (because of logic in the calling function). */
222 fputs(buffer, debug);
227 /*! \brief Helper function for parsing various input about the number
228 of OpenMP threads to use in various modules and deciding what to
229 do about it. */
230 static void manage_number_of_openmp_threads(const gmx::MDLogger &mdlog,
231 const t_commrec *cr,
232 bool bOMP,
233 int nthreads_hw_avail,
234 int omp_nthreads_req,
235 int omp_nthreads_pme_req,
236 gmx_bool gmx_unused bThisNodePMEOnly,
237 gmx_bool bFullOmpSupport,
238 int nppn,
239 gmx_bool bSepPME)
241 int nth;
242 char *env;
244 #if GMX_THREAD_MPI
245 /* modth is shared among tMPI threads, so for thread safety, the
246 * detection is done on the master only. It is not thread-safe
247 * with multiple simulations, but that's anyway not supported by
248 * tMPI. */
249 if (!SIMMASTER(cr))
251 return;
253 #else
254 GMX_UNUSED_VALUE(cr);
255 #endif
257 if (modth.initialized)
259 /* Just return if the initialization has already been
260 done. This could only happen if gmx_omp_nthreads_init() has
261 already been called. */
262 return;
265 /* With full OpenMP support (verlet scheme) set the number of threads
266 * per process / default:
267 * - 1 if not compiled with OpenMP or
268 * - OMP_NUM_THREADS if the env. var is set, or
269 * - omp_nthreads_req = #of threads requested by the user on the mdrun
270 * command line, otherwise
271 * - take the max number of available threads and distribute them
272 * on the processes/tMPI threads.
273 * ~ The GMX_*_NUM_THREADS env var overrides the number of threads of
274 * the respective module and it has to be used in conjunction with
275 * OMP_NUM_THREADS.
277 * With the group scheme OpenMP multithreading is only supported in PME,
278 * for all other modules nthreads is set to 1.
279 * The number of PME threads is equal to:
280 * - 1 if not compiled with OpenMP or
281 * - GMX_PME_NUM_THREADS if defined, otherwise
282 * - OMP_NUM_THREADS if defined, otherwise
283 * - 1
285 nth = 1;
286 if ((env = getenv("OMP_NUM_THREADS")) != nullptr)
288 if (!bOMP && (std::strncmp(env, "1", 1) != 0))
290 gmx_warning("OMP_NUM_THREADS is set, but %s was compiled without OpenMP support!",
291 gmx::getProgramContext().displayName());
293 else
295 nth = gmx_omp_get_max_threads();
298 else if (omp_nthreads_req > 0)
300 nth = omp_nthreads_req;
302 else if (bFullOmpSupport && bOMP)
304 /* max available threads per node */
305 nth = nthreads_hw_avail;
307 /* divide the threads among the MPI processes/tMPI threads */
308 if (nth >= nppn)
310 nth /= nppn;
312 else
314 nth = 1;
318 /* now we have the global values, set them:
319 * - 1 if not compiled with OpenMP and for the group scheme
320 * - nth for the verlet scheme when compiled with OpenMP
322 if (bFullOmpSupport && bOMP)
324 modth.gnth = nth;
326 else
328 modth.gnth = 1;
331 if (bSepPME)
333 if (omp_nthreads_pme_req > 0)
335 modth.gnth_pme = omp_nthreads_pme_req;
337 else
339 modth.gnth_pme = nth;
342 else
344 modth.gnth_pme = 0;
347 /* now set the per-module values */
348 modth.nth[emntDefault] = modth.gnth;
349 pick_module_nthreads(mdlog, emntDomdec, bFullOmpSupport, bSepPME);
350 pick_module_nthreads(mdlog, emntPairsearch, bFullOmpSupport, bSepPME);
351 pick_module_nthreads(mdlog, emntNonbonded, bFullOmpSupport, bSepPME);
352 pick_module_nthreads(mdlog, emntBonded, bFullOmpSupport, bSepPME);
353 pick_module_nthreads(mdlog, emntPME, bFullOmpSupport, bSepPME);
354 pick_module_nthreads(mdlog, emntUpdate, bFullOmpSupport, bSepPME);
355 pick_module_nthreads(mdlog, emntVSITE, bFullOmpSupport, bSepPME);
356 pick_module_nthreads(mdlog, emntLINCS, bFullOmpSupport, bSepPME);
357 pick_module_nthreads(mdlog, emntSETTLE, bFullOmpSupport, bSepPME);
359 /* set the number of threads globally */
360 if (bOMP)
362 #if !GMX_THREAD_MPI
363 if (bThisNodePMEOnly)
365 gmx_omp_set_num_threads(modth.gnth_pme);
367 else
368 #endif /* GMX_THREAD_MPI */
370 if (bFullOmpSupport)
372 gmx_omp_set_num_threads(nth);
374 else
376 gmx_omp_set_num_threads(1);
381 modth.initialized = TRUE;
384 /*! \brief Report on the OpenMP settings that will be used */
385 static void
386 reportOpenmpSettings(const gmx::MDLogger &mdlog,
387 const t_commrec *cr,
388 gmx_bool bOMP,
389 gmx_bool bFullOmpSupport,
390 gmx_bool bSepPME)
392 #if GMX_THREAD_MPI
393 const char *mpi_str = "per tMPI thread";
394 #else
395 const char *mpi_str = "per MPI process";
396 #endif
397 int nth_min, nth_max, nth_pme_min, nth_pme_max;
399 /* inform the user about the settings */
400 if (!bOMP)
402 return;
405 #if GMX_MPI
406 if (cr->nnodes + cr->npmenodes > 1)
408 /* Get the min and max thread counts over the MPI ranks */
409 int buf_in[4], buf_out[4];
411 buf_in[0] = -modth.gnth;
412 buf_in[1] = modth.gnth;
413 buf_in[2] = -modth.gnth_pme;
414 buf_in[3] = modth.gnth_pme;
416 MPI_Allreduce(buf_in, buf_out, 4, MPI_INT, MPI_MAX, cr->mpi_comm_mysim);
418 nth_min = -buf_out[0];
419 nth_max = buf_out[1];
420 nth_pme_min = -buf_out[2];
421 nth_pme_max = buf_out[3];
423 else
424 #endif
426 nth_min = modth.gnth;
427 nth_max = modth.gnth;
428 nth_pme_min = modth.gnth_pme;
429 nth_pme_max = modth.gnth_pme;
432 /* for group scheme we print PME threads info only */
433 if (bFullOmpSupport)
435 if (nth_max == nth_min)
437 GMX_LOG(mdlog.warning).appendTextFormatted(
438 "Using %d OpenMP thread%s %s",
439 nth_min, nth_min > 1 ? "s" : "",
440 cr->nnodes > 1 ? mpi_str : "");
442 else
444 GMX_LOG(mdlog.warning).appendTextFormatted(
445 "Using %d - %d OpenMP threads %s",
446 nth_min, nth_max, mpi_str);
449 if (bSepPME && (nth_pme_min != nth_min || nth_pme_max != nth_max))
451 if (nth_pme_max == nth_pme_min)
453 GMX_LOG(mdlog.warning).appendTextFormatted(
454 "Using %d OpenMP thread%s %s for PME",
455 nth_pme_min, nth_pme_min > 1 ? "s" : "",
456 cr->nnodes > 1 ? mpi_str : "");
458 else
460 GMX_LOG(mdlog.warning).appendTextFormatted(
461 "Using %d - %d OpenMP threads %s for PME",
462 nth_pme_min, nth_pme_max, mpi_str);
465 GMX_LOG(mdlog.warning);
468 /*! \brief Detect and warn about oversubscription of cores.
470 * \todo This could probably live elsewhere, since it is not specifc
471 * to OpenMP, and only needs modth.gnth.
473 * \todo Enable this for separate PME nodes as well! */
474 static void
475 issueOversubscriptionWarning(const gmx::MDLogger &mdlog,
476 const t_commrec *cr,
477 int nthreads_hw_avail,
478 int nppn,
479 gmx_bool bSepPME)
481 char sbuf[STRLEN], sbuf1[STRLEN], sbuf2[STRLEN];
483 if (bSepPME || 0 != cr->rank_pp_intranode)
485 return;
488 if (modth.gnth*nppn > nthreads_hw_avail)
490 sprintf(sbuf, "threads");
491 sbuf1[0] = '\0';
492 sprintf(sbuf2, "O");
493 #if GMX_MPI
494 if (modth.gnth == 1)
496 #if GMX_THREAD_MPI
497 sprintf(sbuf, "thread-MPI threads");
498 #else
499 sprintf(sbuf, "MPI processes");
500 sprintf(sbuf1, " per rank");
501 sprintf(sbuf2, "On rank %d: o", cr->sim_nodeid);
502 #endif
504 #endif
505 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
506 "WARNING: %sversubscribing the available %d logical CPU cores%s with %d %s.\n"
507 " This will cause considerable performance loss!",
508 sbuf2, nthreads_hw_avail, sbuf1, nppn*modth.gnth, sbuf);
512 void gmx_omp_nthreads_init(const gmx::MDLogger &mdlog, t_commrec *cr,
513 int nthreads_hw_avail,
514 int omp_nthreads_req,
515 int omp_nthreads_pme_req,
516 gmx_bool bThisNodePMEOnly,
517 gmx_bool bFullOmpSupport)
519 int nppn;
520 gmx_bool bSepPME;
522 const bool bOMP = GMX_OPENMP;
524 /* number of MPI processes/threads per physical node */
525 nppn = cr->nrank_intranode;
527 bSepPME = (thisRankHasDuty(cr, DUTY_PP) != thisRankHasDuty(cr, DUTY_PME));
529 manage_number_of_openmp_threads(mdlog, cr, bOMP,
530 nthreads_hw_avail,
531 omp_nthreads_req, omp_nthreads_pme_req,
532 bThisNodePMEOnly, bFullOmpSupport,
533 nppn, bSepPME);
534 #if GMX_THREAD_MPI
535 /* Non-master threads have to wait for the OpenMP management to be
536 * done, so that code elsewhere that uses OpenMP can be certain
537 * the setup is complete. */
538 if (PAR(cr))
540 MPI_Barrier(cr->mpi_comm_mysim);
542 #endif
544 reportOpenmpSettings(mdlog, cr, bOMP, bFullOmpSupport, bSepPME);
545 issueOversubscriptionWarning(mdlog, cr, nthreads_hw_avail, nppn, bSepPME);
548 int gmx_omp_nthreads_get(int mod)
550 if (mod < 0 || mod >= emntNR)
552 /* invalid module queried */
553 return -1;
555 else
557 return modth.nth[mod];
561 void
562 gmx_omp_nthreads_set(int mod, int nthreads)
564 /* Catch an attempt to set the number of threads on an invalid
565 * OpenMP module. */
566 GMX_RELEASE_ASSERT(mod >= 0 && mod < emntNR, "Trying to set nthreads on invalid OpenMP module");
568 modth.nth[mod] = nthreads;