Merge release-4-6 into release-5-0
[alexxy/gromacs.git] / src / gromacs / gmxlib / gmx_omp_nthreads.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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 <stdio.h>
37 #include <stdlib.h>
38 #include <string.h>
39 #include <assert.h>
40
41 #ifdef HAVE_CONFIG_H
42 #include <config.h>
43 #endif
44
45 #include "gmx_fatal.h"
46 #include "typedefs.h"
47 #include "macros.h"
48 #include "network.h"
49 #include "copyrite.h"
50 #include "gmx_omp_nthreads.h"
51 #include "md_logging.h"
52
53 #include "gromacs/utility/gmxomp.h"
54
55 /** Structure with the number of threads for each OpenMP multi-threaded
56  *  algorithmic module in mdrun. */
57 typedef struct
58 {
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. */
61
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;
65
66 /** Names of environment variables to set the per module number of threads.
67  *
68  *  Indexed with the values of module_nth_t.
69  * */
70 static const char *modth_env_var[emntNR] =
71 {
72     "GMX_DEFAULT_NUM_THREADS should never be set",
73     "GMX_DOMDEC_NUM_THREADS", "GMX_PAIRSEARCH_NUM_THREADS",
74     "GMX_NONBONDED_NUM_THREADS", "GMX_BONDED_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"
78 };
79
80 /** Names of the modules. */
81 static const char *mod_name[emntNR] =
82 {
83     "default", "domain decomposition", "pair search", "non-bonded",
84     "bonded", "PME", "update", "LINCS", "SETTLE"
85 };
86
87 /** Number of threads for each algorithmic module.
88  *
89  *  File-scope global variable that gets set once in pick_module_nthreads()
90  *  and queried via gmx_omp_nthreads_get().
91  *
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};
96
97
98 /** Determine the number of threads for module \p mod.
99  *
100  *  \p m takes values form the module_nth_t enum and maps these to the
101  *  corresponding value in modth_env_var.
102  *
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.
106  *
107  *  The "group" scheme supports OpenMP only in PME and in thise case all but
108  *  the PME nthread values default to 1.
109  */
110 static void pick_module_nthreads(FILE *fplog, int m,
111                                  gmx_bool bSimMaster,
112                                  gmx_bool bFullOmpSupport,
113                                  gmx_bool bSepPME)
114 {
115     char    *env;
116     int      nth;
117     char     sbuf[STRLEN];
118     gmx_bool bOMP;
119
120 #ifdef GMX_OPENMP
121     bOMP = TRUE;
122 #else
123     bOMP = FALSE;
124 #endif /* GMX_OPENMP */
125
126     /* The default should never be set through a GMX_*_NUM_THREADS env var
127      * as it's always equal with gnth. */
128     if (m == emntDefault)
129     {
130         return;
131     }
132
133     /* check the environment variable */
134     if ((env = getenv(modth_env_var[m])) != NULL)
135     {
136         sscanf(env, "%d", &nth);
137
138         if (!bOMP)
139         {
140             gmx_warning("%s=%d is set, but %s is compiled without OpenMP!",
141                         modth_env_var[m], nth, ShortProgram());
142         }
143
144         /* with the verlet codepath, when any GMX_*_NUM_THREADS env var is set,
145          * OMP_NUM_THREADS also has to be set */
146         if (bFullOmpSupport && getenv("OMP_NUM_THREADS") == NULL)
147         {
148             gmx_fatal(FARGS, "%s=%d is set, the default number of threads also "
149                       "needs to be set with OMP_NUM_THREADS!",
150                       modth_env_var[m], nth);
151         }
152
153         /* with the group scheme warn if any env var except PME is set */
154         if (!bFullOmpSupport)
155         {
156             if (m != emntPME)
157             {
158                 gmx_warning("%s=%d is set, but OpenMP multithreading is not "
159                             "supported in %s!",
160                             modth_env_var[m], nth, mod_name[m]);
161                 nth = 1;
162             }
163         }
164
165         /* only babble if we are really overriding with a different value */
166         if ((bSepPME && m == emntPME && nth != modth.gnth_pme) || (nth != modth.gnth))
167         {
168             sprintf(sbuf, "%s=%d set, overriding the default number of %s threads",
169                     modth_env_var[m], nth, mod_name[m]);
170             if (bSimMaster)
171             {
172                 fprintf(stderr, "\n%s\n", sbuf);
173             }
174             if (fplog)
175             {
176                 fprintf(fplog, "%s\n", sbuf);
177             }
178         }
179     }
180     else
181     {
182         /* pick the global PME node nthreads if we are setting the number
183          * of threads in separate PME nodes  */
184         nth = (bSepPME && m == emntPME) ? modth.gnth_pme : modth.gnth;
185     }
186
187     gmx_omp_nthreads_set(m, nth);
188 }
189
190 void gmx_omp_nthreads_read_env(int     *nthreads_omp,
191                                gmx_bool bIsSimMaster)
192 {
193     char    *env;
194     gmx_bool bCommandLineSetNthreadsOMP = *nthreads_omp > 0;
195     char     buffer[STRLEN];
196
197     assert(nthreads_omp);
198
199     if ((env = getenv("OMP_NUM_THREADS")) != NULL)
200     {
201         int nt_omp;
202
203         sscanf(env, "%d", &nt_omp);
204         if (nt_omp <= 0)
205         {
206             gmx_fatal(FARGS, "OMP_NUM_THREADS is invalid: '%s'", env);
207         }
208
209         if (bCommandLineSetNthreadsOMP && nt_omp != *nthreads_omp)
210         {
211             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);
212         }
213
214         /* Setting the number of OpenMP threads. */
215         *nthreads_omp = nt_omp;
216
217         /* Output the results */
218         sprintf(buffer,
219                 "The number of OpenMP threads was set by environment variable OMP_NUM_THREADS to %d%s\n",
220                 nt_omp,
221                 bCommandLineSetNthreadsOMP ? " (and the command-line setting agreed with that)" : "");
222         if (bIsSimMaster)
223         {
224             /* This prints once per simulation for multi-simulations,
225              * which might help diagnose issues with inhomogenous
226              * cluster setups. */
227             fputs(buffer, stderr);
228         }
229         if (debug)
230         {
231             /* This prints once per process for real MPI (i.e. once
232              * per debug file), and once per simulation for thread MPI
233              * (because of logic in the calling function). */
234             fputs(buffer, debug);
235         }
236     }
237 }
238
239 void gmx_omp_nthreads_init(FILE *fplog, t_commrec *cr,
240                            int nthreads_hw_avail,
241                            int omp_nthreads_req,
242                            int omp_nthreads_pme_req,
243                            gmx_bool gmx_unused bThisNodePMEOnly,
244                            gmx_bool bFullOmpSupport)
245 {
246     int      nth, nth_pmeonly, gmx_maxth, nppn;
247     char    *env;
248     gmx_bool bSepPME, bOMP;
249
250 #ifdef GMX_OPENMP
251     bOMP = TRUE;
252 #else
253     bOMP = FALSE;
254 #endif /* GMX_OPENMP */
255
256     /* number of MPI processes/threads per physical node */
257     nppn = cr->nrank_intranode;
258
259     bSepPME = ( (cr->duty & DUTY_PP) && !(cr->duty & DUTY_PME)) ||
260         (!(cr->duty & DUTY_PP) &&  (cr->duty & DUTY_PME));
261
262 #ifdef GMX_THREAD_MPI
263     /* modth is shared among tMPI threads, so for thread safety do the
264      * detection is done on the master only. It is not thread-safe with
265      * multiple simulations, but that's anyway not supported by tMPI. */
266     if (SIMMASTER(cr))
267 #endif
268     {
269         /* just return if the initialization has already been done */
270         if (modth.initialized)
271         {
272             return;
273         }
274
275         /* With full OpenMP support (verlet scheme) set the number of threads
276          * per process / default:
277          * - 1 if not compiled with OpenMP or
278          * - OMP_NUM_THREADS if the env. var is set, or
279          * - omp_nthreads_req = #of threads requested by the user on the mdrun
280          *   command line, otherwise
281          * - take the max number of available threads and distribute them
282          *   on the processes/tMPI threads.
283          * ~ The GMX_*_NUM_THREADS env var overrides the number of threads of
284          *   the respective module and it has to be used in conjunction with
285          *   OMP_NUM_THREADS.
286          *
287          * With the group scheme OpenMP multithreading is only supported in PME,
288          * for all other modules nthreads is set to 1.
289          * The number of PME threads is equal to:
290          * - 1 if not compiled with OpenMP or
291          * - GMX_PME_NUM_THREADS if defined, otherwise
292          * - OMP_NUM_THREADS if defined, otherwise
293          * - 1
294          */
295         nth = 1;
296         if ((env = getenv("OMP_NUM_THREADS")) != NULL)
297         {
298             if (!bOMP && (strncmp(env, "1", 1) != 0))
299             {
300                 gmx_warning("OMP_NUM_THREADS is set, but %s was compiled without OpenMP support!",
301                             ShortProgram());
302             }
303             else
304             {
305                 nth = gmx_omp_get_max_threads();
306             }
307         }
308         else if (omp_nthreads_req > 0)
309         {
310             nth = omp_nthreads_req;
311         }
312         else if (bFullOmpSupport && bOMP)
313         {
314             /* max available threads per node */
315             nth = nthreads_hw_avail;
316
317             /* divide the threads among the MPI processes/tMPI threads */
318             if (nth >= nppn)
319             {
320                 nth /= nppn;
321             }
322             else
323             {
324                 nth = 1;
325             }
326         }
327
328         /* now we have the global values, set them:
329          * - 1 if not compiled with OpenMP and for the group scheme
330          * - nth for the verlet scheme when compiled with OpenMP
331          */
332         if (bFullOmpSupport && bOMP)
333         {
334             modth.gnth = nth;
335         }
336         else
337         {
338             modth.gnth = 1;
339         }
340
341         if (bSepPME)
342         {
343             if (omp_nthreads_pme_req > 0)
344             {
345                 modth.gnth_pme = omp_nthreads_pme_req;
346             }
347             else
348             {
349                 modth.gnth_pme = nth;
350             }
351         }
352         else
353         {
354             modth.gnth_pme = 0;
355         }
356
357         /* now set the per-module values */
358         modth.nth[emntDefault] = modth.gnth;
359         pick_module_nthreads(fplog, emntDomdec, SIMMASTER(cr), bFullOmpSupport, bSepPME);
360         pick_module_nthreads(fplog, emntPairsearch, SIMMASTER(cr), bFullOmpSupport, bSepPME);
361         pick_module_nthreads(fplog, emntNonbonded, SIMMASTER(cr), bFullOmpSupport, bSepPME);
362         pick_module_nthreads(fplog, emntBonded, SIMMASTER(cr), bFullOmpSupport, bSepPME);
363         pick_module_nthreads(fplog, emntPME, SIMMASTER(cr), bFullOmpSupport, bSepPME);
364         pick_module_nthreads(fplog, emntUpdate, SIMMASTER(cr), bFullOmpSupport, bSepPME);
365         pick_module_nthreads(fplog, emntVSITE, SIMMASTER(cr), bFullOmpSupport, bSepPME);
366         pick_module_nthreads(fplog, emntLINCS, SIMMASTER(cr), bFullOmpSupport, bSepPME);
367         pick_module_nthreads(fplog, emntSETTLE, SIMMASTER(cr), bFullOmpSupport, bSepPME);
368
369         /* set the number of threads globally */
370         if (bOMP)
371         {
372 #ifndef GMX_THREAD_MPI
373             if (bThisNodePMEOnly)
374             {
375                 gmx_omp_set_num_threads(modth.gnth_pme);
376             }
377             else
378 #endif      /* GMX_THREAD_MPI */
379             {
380                 if (bFullOmpSupport)
381                 {
382                     gmx_omp_set_num_threads(nth);
383                 }
384                 else
385                 {
386                     gmx_omp_set_num_threads(1);
387                 }
388             }
389         }
390
391         modth.initialized = TRUE;
392     }
393 #ifdef GMX_THREAD_MPI
394     /* Non-master threads have to wait for the detection to be done. */
395     if (PAR(cr))
396     {
397         MPI_Barrier(cr->mpi_comm_mysim);
398     }
399 #endif
400
401     /* inform the user about the settings */
402     if (bOMP)
403     {
404 #ifdef GMX_THREAD_MPI
405         const char *mpi_str = "per tMPI thread";
406 #else
407         const char *mpi_str = "per MPI process";
408 #endif
409
410         /* for group scheme we print PME threads info only */
411         if (bFullOmpSupport)
412         {
413             md_print_info(cr, fplog, "Using %d OpenMP thread%s %s\n",
414                           modth.gnth, modth.gnth > 1 ? "s" : "",
415                           cr->nnodes > 1 ? mpi_str : "");
416         }
417         if (bSepPME && modth.gnth_pme != modth.gnth)
418         {
419             md_print_info(cr, fplog, "Using %d OpenMP thread%s %s for PME\n",
420                           modth.gnth_pme, modth.gnth_pme > 1 ? "s" : "",
421                           cr->nnodes > 1 ? mpi_str : "");
422         }
423     }
424
425     /* detect and warn about oversubscription
426      * TODO: enable this for separate PME nodes as well! */
427     if (!bSepPME && cr->rank_pp_intranode == 0)
428     {
429         char sbuf[STRLEN], sbuf1[STRLEN], sbuf2[STRLEN];
430
431         if (modth.gnth*nppn > nthreads_hw_avail)
432         {
433             sprintf(sbuf, "threads");
434             sbuf1[0] = '\0';
435             sprintf(sbuf2, "O");
436 #ifdef GMX_MPI
437             if (modth.gnth == 1)
438             {
439 #ifdef GMX_THREAD_MPI
440                 sprintf(sbuf, "thread-MPI threads");
441 #else
442                 sprintf(sbuf, "MPI processes");
443                 sprintf(sbuf1, " per node");
444                 sprintf(sbuf2, "On node %d: o", cr->sim_nodeid);
445 #endif
446             }
447 #endif
448             md_print_warn(cr, fplog,
449                           "WARNING: %sversubscribing the available %d logical CPU cores%s with %d %s.\n"
450                           "         This will cause considerable performance loss!",
451                           sbuf2, nthreads_hw_avail, sbuf1, nppn*modth.gnth, sbuf);
452         }
453     }
454 }
455
456 int gmx_omp_nthreads_get(int mod)
457 {
458     if (mod < 0 || mod >= emntNR)
459     {
460         /* invalid module queried */
461         return -1;
462     }
463     else
464     {
465         return modth.nth[mod];
466     }
467 }
468
469 void
470 gmx_omp_nthreads_set(int mod, int nthreads)
471 {
472     /* Catch an attempt to set the number of threads on an invalid
473      * OpenMP module. */
474     assert(mod >= 0 && mod < emntNR);
475
476     modth.nth[mod] = nthreads;
477 }