Merge release-5-0 into master
[gromacs.git] / src / gromacs / gmxlib / gmx_omp_nthreads.c
blob76824dee392284ae76893bd3c5a40f5b5389553e
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, 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 <stdio.h>
39 #include <stdlib.h>
40 #include <string.h>
41 #include <assert.h>
43 #include "config.h"
45 #include "gromacs/legacyheaders/typedefs.h"
46 #include "gromacs/legacyheaders/types/commrec.h"
47 #include "gromacs/legacyheaders/macros.h"
48 #include "gromacs/legacyheaders/network.h"
49 #include "gromacs/legacyheaders/copyrite.h"
50 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
51 #include "gromacs/legacyheaders/md_logging.h"
53 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/gmxomp.h"
57 /** Structure with the number of threads for each OpenMP multi-threaded
58 * algorithmic module in mdrun. */
59 typedef struct
61 int gnth; /**< Global num. of threads per PP or PP+PME process/tMPI thread. */
62 int gnth_pme; /**< Global num. of threads per PME only process/tMPI thread. */
64 int nth[emntNR]; /**< Number of threads for each module, indexed with module_nth_t */
65 gmx_bool initialized; /**< TRUE if the module as been initialized. */
66 } omp_module_nthreads_t;
68 /** Names of environment variables to set the per module number of threads.
70 * Indexed with the values of module_nth_t.
71 * */
72 static const char *modth_env_var[emntNR] =
74 "GMX_DEFAULT_NUM_THREADS should never be set",
75 "GMX_DOMDEC_NUM_THREADS", "GMX_PAIRSEARCH_NUM_THREADS",
76 "GMX_NONBONDED_NUM_THREADS", "GMX_BONDED_NUM_THREADS",
77 "GMX_PME_NUM_THREADS", "GMX_UPDATE_NUM_THREADS",
78 "GMX_VSITE_NUM_THREADS",
79 "GMX_LINCS_NUM_THREADS", "GMX_SETTLE_NUM_THREADS"
82 /** Names of the modules. */
83 static const char *mod_name[emntNR] =
85 "default", "domain decomposition", "pair search", "non-bonded",
86 "bonded", "PME", "update", "LINCS", "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 deafult.
109 * The "group" scheme supports OpenMP only in PME and in thise case all but
110 * the PME nthread values default to 1.
112 static void pick_module_nthreads(FILE *fplog, int m,
113 gmx_bool bSimMaster,
114 gmx_bool bFullOmpSupport,
115 gmx_bool bSepPME)
117 char *env;
118 int nth;
119 char sbuf[STRLEN];
120 gmx_bool bOMP;
122 #ifdef GMX_OPENMP
123 bOMP = TRUE;
124 #else
125 bOMP = FALSE;
126 #endif /* GMX_OPENMP */
128 /* The default should never be set through a GMX_*_NUM_THREADS env var
129 * as it's always equal with gnth. */
130 if (m == emntDefault)
132 return;
135 /* check the environment variable */
136 if ((env = getenv(modth_env_var[m])) != NULL)
138 sscanf(env, "%d", &nth);
140 if (!bOMP)
142 gmx_warning("%s=%d is set, but %s is compiled without OpenMP!",
143 modth_env_var[m], nth, ShortProgram());
146 /* with the verlet codepath, when any GMX_*_NUM_THREADS env var is set,
147 * OMP_NUM_THREADS also has to be set */
148 if (bFullOmpSupport && getenv("OMP_NUM_THREADS") == NULL)
150 gmx_warning("%s=%d is set, the default number of threads also "
151 "needs to be set with OMP_NUM_THREADS!",
152 modth_env_var[m], nth);
155 /* with the group scheme warn if any env var except PME is set */
156 if (!bFullOmpSupport)
158 if (m != emntPME)
160 gmx_warning("%s=%d is set, but OpenMP multithreading is not "
161 "supported in %s!",
162 modth_env_var[m], nth, mod_name[m]);
163 nth = 1;
167 /* only babble if we are really overriding with a different value */
168 if ((bSepPME && m == emntPME && nth != modth.gnth_pme) || (nth != modth.gnth))
170 sprintf(sbuf, "%s=%d set, overriding the default number of %s threads",
171 modth_env_var[m], nth, mod_name[m]);
172 if (bSimMaster)
174 fprintf(stderr, "\n%s\n", sbuf);
176 if (fplog)
178 fprintf(fplog, "%s\n", sbuf);
182 else
184 /* pick the global PME node nthreads if we are setting the number
185 * of threads in separate PME nodes */
186 nth = (bSepPME && m == emntPME) ? modth.gnth_pme : modth.gnth;
189 gmx_omp_nthreads_set(m, nth);
192 void gmx_omp_nthreads_read_env(int *nthreads_omp,
193 gmx_bool bIsSimMaster)
195 char *env;
196 gmx_bool bCommandLineSetNthreadsOMP = *nthreads_omp > 0;
197 char buffer[STRLEN];
199 assert(nthreads_omp);
201 if ((env = getenv("OMP_NUM_THREADS")) != NULL)
203 int nt_omp;
205 sscanf(env, "%d", &nt_omp);
206 if (nt_omp <= 0)
208 gmx_fatal(FARGS, "OMP_NUM_THREADS is invalid: '%s'", env);
211 if (bCommandLineSetNthreadsOMP && nt_omp != *nthreads_omp)
213 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);
216 /* Setting the number of OpenMP threads. */
217 *nthreads_omp = nt_omp;
219 /* Output the results */
220 sprintf(buffer,
221 "The number of OpenMP threads was set by environment variable OMP_NUM_THREADS to %d%s\n",
222 nt_omp,
223 bCommandLineSetNthreadsOMP ? " (and the command-line setting agreed with that)" : "");
224 if (bIsSimMaster)
226 /* This prints once per simulation for multi-simulations,
227 * which might help diagnose issues with inhomogenous
228 * cluster setups. */
229 fputs(buffer, stderr);
231 if (debug)
233 /* This prints once per process for real MPI (i.e. once
234 * per debug file), and once per simulation for thread MPI
235 * (because of logic in the calling function). */
236 fputs(buffer, debug);
241 void gmx_omp_nthreads_init(FILE *fplog, t_commrec *cr,
242 int nthreads_hw_avail,
243 int omp_nthreads_req,
244 int omp_nthreads_pme_req,
245 gmx_bool gmx_unused bThisNodePMEOnly,
246 gmx_bool bFullOmpSupport)
248 int nth, nth_pmeonly, gmx_maxth, nppn;
249 char *env;
250 gmx_bool bSepPME, bOMP;
252 #ifdef GMX_OPENMP
253 bOMP = TRUE;
254 #else
255 bOMP = FALSE;
256 #endif /* GMX_OPENMP */
258 /* number of MPI processes/threads per physical node */
259 nppn = cr->nrank_intranode;
261 bSepPME = ( (cr->duty & DUTY_PP) && !(cr->duty & DUTY_PME)) ||
262 (!(cr->duty & DUTY_PP) && (cr->duty & DUTY_PME));
264 #ifdef GMX_THREAD_MPI
265 /* modth is shared among tMPI threads, so for thread safety do the
266 * detection is done on the master only. It is not thread-safe with
267 * multiple simulations, but that's anyway not supported by tMPI. */
268 if (SIMMASTER(cr))
269 #endif
271 /* just return if the initialization has already been done */
272 if (modth.initialized)
274 return;
277 /* With full OpenMP support (verlet scheme) set the number of threads
278 * per process / default:
279 * - 1 if not compiled with OpenMP or
280 * - OMP_NUM_THREADS if the env. var is set, or
281 * - omp_nthreads_req = #of threads requested by the user on the mdrun
282 * command line, otherwise
283 * - take the max number of available threads and distribute them
284 * on the processes/tMPI threads.
285 * ~ The GMX_*_NUM_THREADS env var overrides the number of threads of
286 * the respective module and it has to be used in conjunction with
287 * OMP_NUM_THREADS.
289 * With the group scheme OpenMP multithreading is only supported in PME,
290 * for all other modules nthreads is set to 1.
291 * The number of PME threads is equal to:
292 * - 1 if not compiled with OpenMP or
293 * - GMX_PME_NUM_THREADS if defined, otherwise
294 * - OMP_NUM_THREADS if defined, otherwise
295 * - 1
297 nth = 1;
298 if ((env = getenv("OMP_NUM_THREADS")) != NULL)
300 if (!bOMP && (strncmp(env, "1", 1) != 0))
302 gmx_warning("OMP_NUM_THREADS is set, but %s was compiled without OpenMP support!",
303 ShortProgram());
305 else
307 nth = gmx_omp_get_max_threads();
310 else if (omp_nthreads_req > 0)
312 nth = omp_nthreads_req;
314 else if (bFullOmpSupport && bOMP)
316 /* max available threads per node */
317 nth = nthreads_hw_avail;
319 /* divide the threads among the MPI processes/tMPI threads */
320 if (nth >= nppn)
322 nth /= nppn;
324 else
326 nth = 1;
330 /* now we have the global values, set them:
331 * - 1 if not compiled with OpenMP and for the group scheme
332 * - nth for the verlet scheme when compiled with OpenMP
334 if (bFullOmpSupport && bOMP)
336 modth.gnth = nth;
338 else
340 modth.gnth = 1;
343 if (bSepPME)
345 if (omp_nthreads_pme_req > 0)
347 modth.gnth_pme = omp_nthreads_pme_req;
349 else
351 modth.gnth_pme = nth;
354 else
356 modth.gnth_pme = 0;
359 /* now set the per-module values */
360 modth.nth[emntDefault] = modth.gnth;
361 pick_module_nthreads(fplog, emntDomdec, SIMMASTER(cr), bFullOmpSupport, bSepPME);
362 pick_module_nthreads(fplog, emntPairsearch, SIMMASTER(cr), bFullOmpSupport, bSepPME);
363 pick_module_nthreads(fplog, emntNonbonded, SIMMASTER(cr), bFullOmpSupport, bSepPME);
364 pick_module_nthreads(fplog, emntBonded, SIMMASTER(cr), bFullOmpSupport, bSepPME);
365 pick_module_nthreads(fplog, emntPME, SIMMASTER(cr), bFullOmpSupport, bSepPME);
366 pick_module_nthreads(fplog, emntUpdate, SIMMASTER(cr), bFullOmpSupport, bSepPME);
367 pick_module_nthreads(fplog, emntVSITE, SIMMASTER(cr), bFullOmpSupport, bSepPME);
368 pick_module_nthreads(fplog, emntLINCS, SIMMASTER(cr), bFullOmpSupport, bSepPME);
369 pick_module_nthreads(fplog, emntSETTLE, SIMMASTER(cr), bFullOmpSupport, bSepPME);
371 /* set the number of threads globally */
372 if (bOMP)
374 #ifndef GMX_THREAD_MPI
375 if (bThisNodePMEOnly)
377 gmx_omp_set_num_threads(modth.gnth_pme);
379 else
380 #endif /* GMX_THREAD_MPI */
382 if (bFullOmpSupport)
384 gmx_omp_set_num_threads(nth);
386 else
388 gmx_omp_set_num_threads(1);
393 modth.initialized = TRUE;
395 #ifdef GMX_THREAD_MPI
396 /* Non-master threads have to wait for the detection to be done. */
397 if (PAR(cr))
399 MPI_Barrier(cr->mpi_comm_mysim);
401 #endif
403 /* inform the user about the settings */
404 if (bOMP)
406 #ifdef GMX_THREAD_MPI
407 const char *mpi_str = "per tMPI thread";
408 #else
409 const char *mpi_str = "per MPI process";
410 #endif
412 /* for group scheme we print PME threads info only */
413 if (bFullOmpSupport)
415 md_print_info(cr, fplog, "Using %d OpenMP thread%s %s\n",
416 modth.gnth, modth.gnth > 1 ? "s" : "",
417 cr->nnodes > 1 ? mpi_str : "");
419 if (bSepPME && modth.gnth_pme != modth.gnth)
421 md_print_info(cr, fplog, "Using %d OpenMP thread%s %s for PME\n",
422 modth.gnth_pme, modth.gnth_pme > 1 ? "s" : "",
423 cr->nnodes > 1 ? mpi_str : "");
427 /* detect and warn about oversubscription
428 * TODO: enable this for separate PME nodes as well! */
429 if (!bSepPME && cr->rank_pp_intranode == 0)
431 char sbuf[STRLEN], sbuf1[STRLEN], sbuf2[STRLEN];
433 if (modth.gnth*nppn > nthreads_hw_avail)
435 sprintf(sbuf, "threads");
436 sbuf1[0] = '\0';
437 sprintf(sbuf2, "O");
438 #ifdef GMX_MPI
439 if (modth.gnth == 1)
441 #ifdef GMX_THREAD_MPI
442 sprintf(sbuf, "thread-MPI threads");
443 #else
444 sprintf(sbuf, "MPI processes");
445 sprintf(sbuf1, " per rank");
446 sprintf(sbuf2, "On rank %d: o", cr->sim_nodeid);
447 #endif
449 #endif
450 md_print_warn(cr, fplog,
451 "WARNING: %sversubscribing the available %d logical CPU cores%s with %d %s.\n"
452 " This will cause considerable performance loss!",
453 sbuf2, nthreads_hw_avail, sbuf1, nppn*modth.gnth, sbuf);
458 int gmx_omp_nthreads_get(int mod)
460 if (mod < 0 || mod >= emntNR)
462 /* invalid module queried */
463 return -1;
465 else
467 return modth.nth[mod];
471 void
472 gmx_omp_nthreads_set(int mod, int nthreads)
474 /* Catch an attempt to set the number of threads on an invalid
475 * OpenMP module. */
476 assert(mod >= 0 && mod < emntNR);
478 modth.nth[mod] = nthreads;