Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / gmxlib / gmx_omp_nthreads.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  * 
4  *                This source code is part of
5  * 
6  *                 G   R   O   M   A   C   S
7  * 
8  *          GROningen MAchine for Chemical Simulations
9  * 
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2010, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  * 
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  * 
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
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 "statutil.h"
50 #include "gmx_omp.h"
51 #include "gmx_omp_nthreads.h"
52 #include "md_logging.h"
53
54 /** Structure with the number of threads for each OpenMP multi-threaded
55  *  algorithmic module in mdrun. */
56 typedef struct
57 {
58     int gnth;               /**< Global num. of threads per PP or PP+PME process/tMPI thread. */
59     int gnth_pme;           /**< Global num. of threads per PME only process/tMPI thread. */
60
61     int nth[emntNR];        /**< Number of threads for each module, indexed with module_nth_t */
62     gmx_bool initialized;   /**< TRUE if the module as been initialized. */
63 } omp_module_nthreads_t;
64
65 /** Names of environment variables to set the per module number of threads.
66  *
67  *  Indexed with the values of module_nth_t.
68  * */
69 static const char *modth_env_var[emntNR] =
70 {
71     "GMX_DEFAULT_NUM_THREADS should never be set",
72     "GMX_DOMDEC_NUM_THREADS", "GMX_PAIRSEARCH_NUM_THREADS",
73     "GMX_NONBONDED_NUM_THREADS", "GMX_BONDED_NUM_THREADS",
74     "GMX_PME_NUM_THREADS", "GMX_UPDATE_NUM_THREADS",
75     "GMX_LINCS_NUM_THREADS", "GMX_SETTLE_NUM_THREADS"
76 };
77
78 /** Names of the modules. */
79 static const char *mod_name[emntNR] =
80 {
81     "default", "domain decomposition", "pair search", "non-bonded",
82     "bonded", "PME", "update", "LINCS", "SETTLE"
83 };
84
85 /** Number of threads for each algorithmic module.
86  *
87  *  File-scope global variable that gets set once in \init_module_nthreads
88  *  and queried via gmx_omp_nthreads_get.
89  *
90  *  All fields are initialized to 0 which should result in errors if
91  *  the init call is omitted.
92  * */
93 static omp_module_nthreads_t modth = { 0, 0, {0, 0, 0, 0, 0, 0, 0, 0}, FALSE};
94
95
96 /** Determine the number of threads for module \mod.
97  *
98  *  \m takes values form the module_nth_t enum and maps these to the
99  *  corresponding value in modth_env_var.
100  *
101  *  Each number of threads per module takes the default value unless
102  *  GMX_*_NUM_THERADS env var is set, case in which its value overrides
103  *  the deafult.
104  *
105  *  The "group" scheme supports OpenMP only in PME and in thise case all but
106  *  the PME nthread values default to 1.
107  */
108 static int pick_module_nthreads(FILE *fplog, int m,
109                                 gmx_bool bSimMaster,
110                                 gmx_bool bFullOmpSupport,
111                                 gmx_bool bSepPME)
112 {
113     char *env;
114     int  nth;
115     char sbuf[STRLEN];
116     gmx_bool bOMP;
117
118 #ifdef GMX_OPENMP
119     bOMP = TRUE;
120 #else
121     bOMP = FALSE;
122 #endif /* GMX_OPENMP */
123
124     /* The default should never be set through a GMX_*_NUM_THREADS env var
125      * as it's always equal with gnth. */
126     if (m == emntDefault)
127     {
128         return modth.nth[emntDefault];
129     }
130
131     /* check the environment variable */
132     if ((env = getenv(modth_env_var[m])) != NULL)
133     {
134         sscanf(env, "%d", &nth);
135
136         if (!bOMP)
137         {
138             gmx_warning("%s=%d is set, but %s is compiled without OpenMP!",
139                         modth_env_var[m], nth, ShortProgram());
140         }
141
142         /* with the verlet codepath, when any GMX_*_NUM_THREADS env var is set,
143          * OMP_NUM_THREADS also has to be set */
144         if (bFullOmpSupport && getenv("OMP_NUM_THREADS") == NULL)
145         {
146             gmx_fatal(FARGS, "%s=%d is set, the default number of threads also "
147                       "needs to be set with OMP_NUM_THREADS!",
148                       modth_env_var[m], nth);
149         }
150
151         /* with the group scheme warn if any env var except PME is set */
152         if (!bFullOmpSupport)
153         {
154             if (m != emntPME)
155             {
156                 gmx_warning("%s=%d is set, but OpenMP multithreading is not "
157                             "supported in %s!",
158                             modth_env_var[m], nth, mod_name[m]);
159                 nth = 1;
160             }
161         }
162
163         /* only babble if we are really overriding with a different value */
164         if ((bSepPME && m == emntPME && nth != modth.gnth_pme) || (nth != modth.gnth))
165         {
166             sprintf(sbuf, "%s=%d set, overriding the default number of %s threads",
167                     modth_env_var[m], nth, mod_name[m]);
168             if (bSimMaster)
169             {
170                 fprintf(stderr, "\n%s\n", sbuf);
171             }
172             if (fplog)
173             {
174                 fprintf(fplog, "%s\n", sbuf);
175             }
176         }
177     }
178     else
179     {
180         /* pick the global PME node nthreads if we are setting the number
181          * of threads in separate PME nodes  */
182         nth = (bSepPME && m == emntPME) ? modth.gnth_pme : modth.gnth;
183     }
184
185     return modth.nth[m] = nth;
186 }
187
188 void gmx_omp_nthreads_read_env(int *nthreads_omp)
189 {
190     char *env;
191
192     assert(nthreads_omp);
193
194     if ((env = getenv("OMP_NUM_THREADS")) != NULL)
195     {
196         int nt_omp;
197
198         sscanf(env,"%d",&nt_omp);
199         if (nt_omp <= 0)
200         {
201             gmx_fatal(FARGS,"OMP_NUM_THREADS is invalid: '%s'",env);
202         }
203
204         if (*nthreads_omp > 0 && nt_omp != *nthreads_omp)
205         {
206             gmx_fatal(FARGS,"OMP_NUM_THREADS (%d) and the number of threads requested on the command line (%d) have different values",nt_omp,*nthreads_omp);
207         }
208
209         /* Setting the number of OpenMP threads.
210          * NOTE: with tMPI this function is only called on the master node,
211          * but with MPI on all nodes which means lots of messages on stderr.
212          */
213         fprintf(stderr,"Getting the number of OpenMP threads from OMP_NUM_THREADS: %d\n",nt_omp);
214         *nthreads_omp = nt_omp;
215     }
216 }
217
218 void gmx_omp_nthreads_init(FILE *fplog, t_commrec *cr,
219                            int nthreads_hw_avail,
220                            int omp_nthreads_req,
221                            int omp_nthreads_pme_req,
222                            gmx_bool bThisNodePMEOnly,
223                            gmx_bool bFullOmpSupport)
224 {
225     int  nth, nth_pmeonly, gmx_maxth, nppn;
226     char *env;
227     gmx_bool bSepPME, bOMP;
228
229 #ifdef GMX_OPENMP
230     bOMP = TRUE;
231 #else
232     bOMP = FALSE;
233 #endif /* GMX_OPENMP */
234
235     /* number of processes per node */
236     nppn = cr->nnodes_intra;
237
238     bSepPME = ( (cr->duty & DUTY_PP) && !(cr->duty & DUTY_PME)) ||
239               (!(cr->duty & DUTY_PP) &&  (cr->duty & DUTY_PME));
240
241 #ifdef GMX_THREAD_MPI
242     /* modth is shared among tMPI threads, so for thread safety do the
243      * detection is done on the master only. It is not thread-safe with
244      * multiple simulations, but that's anyway not supported by tMPI. */
245     if (SIMMASTER(cr))
246 #endif
247     {
248         /* just return if the initialization has already been done */
249         if (modth.initialized)
250         {
251             return;
252         }
253
254         /* With full OpenMP support (verlet scheme) set the number of threads
255          * per process / default:
256          * - 1 if not compiled with OpenMP or
257          * - OMP_NUM_THREADS if the env. var is set, or
258          * - omp_nthreads_req = #of threads requested by the user on the mdrun
259          *   command line, otherwise
260          * - take the max number of available threads and distribute them
261          *   on the processes/tMPI threads.
262          * ~ The GMX_*_NUM_THREADS env var overrides the number of threads of
263          *   the respective module and it has to be used in conjunction with
264          *   OMP_NUM_THREADS.
265          *
266          * With the group scheme OpenMP multithreading is only supported in PME,
267          * for all other modules nthreads is set to 1.
268          * The number of PME threads is equal to:
269          * - 1 if not compiled with OpenMP or
270          * - GMX_PME_NUM_THREADS if defined, otherwise
271          * - OMP_NUM_THREADS if defined, otherwise
272          * - 1
273          */
274         nth = 1;
275         if ((env = getenv("OMP_NUM_THREADS")) != NULL)
276         {
277             if (!bOMP && (strncmp(env, "1", 1) != 0))
278             {
279                 gmx_warning("OMP_NUM_THREADS is set, but %s was compiled without OpenMP support!",
280                             ShortProgram());
281             }
282             else
283             {
284                 nth = gmx_omp_get_max_threads();
285             }
286         }
287         else if (omp_nthreads_req > 0)
288         {
289             nth = omp_nthreads_req;
290         }
291         else if (bFullOmpSupport && bOMP)
292         {
293             /* max available threads per node */
294             nth = nthreads_hw_avail;
295
296             /* divide the threads among the MPI processes/tMPI threads */
297             if (nth >= nppn)
298             {
299                 nth /= nppn;
300             }
301             else
302             {
303                 nth = 1;
304             }
305         }
306
307         /* now we have the global values, set them:
308          * - 1 if not compiled with OpenMP and for the group scheme
309          * - nth for the verlet scheme when compiled with OpenMP
310          */
311         if (bFullOmpSupport && bOMP)
312         {
313             modth.gnth = nth;
314         }
315         else
316         {
317             modth.gnth = 1;
318         }
319
320         if (bSepPME)
321         {
322             if (omp_nthreads_pme_req > 0)
323             {
324                 modth.gnth_pme = omp_nthreads_pme_req;
325             }
326             else
327             {
328                 modth.gnth_pme = nth;
329             }
330         }
331         else
332         {
333             modth.gnth_pme = 0;
334         }
335
336         /* now set the per-module values */
337         modth.nth[emntDefault] = modth.gnth;
338         pick_module_nthreads(fplog, emntDomdec, SIMMASTER(cr), bFullOmpSupport, bSepPME);
339         pick_module_nthreads(fplog, emntPairsearch, SIMMASTER(cr), bFullOmpSupport, bSepPME);
340         pick_module_nthreads(fplog, emntNonbonded, SIMMASTER(cr), bFullOmpSupport, bSepPME);
341         pick_module_nthreads(fplog, emntBonded, SIMMASTER(cr), bFullOmpSupport, bSepPME);
342         pick_module_nthreads(fplog, emntPME, SIMMASTER(cr), bFullOmpSupport, bSepPME);
343         pick_module_nthreads(fplog, emntUpdate, SIMMASTER(cr), bFullOmpSupport, bSepPME);
344         pick_module_nthreads(fplog, emntLINCS, SIMMASTER(cr), bFullOmpSupport, bSepPME);
345         pick_module_nthreads(fplog, emntSETTLE, SIMMASTER(cr), bFullOmpSupport, bSepPME);
346
347         /* set the number of threads globally */
348         if (bOMP)
349         {
350 #ifndef GMX_THREAD_MPI
351             if (bThisNodePMEOnly)
352             {
353                 gmx_omp_set_num_threads(modth.gnth_pme);
354             }
355             else
356 #endif /* GMX_THREAD_MPI */
357             {
358                 if (bFullOmpSupport)
359                 {
360                     gmx_omp_set_num_threads(nth);
361                 }
362                 else
363                 {
364                     gmx_omp_set_num_threads(1);
365                 }
366             }
367         }
368
369         modth.initialized = TRUE;
370     }
371 #ifdef GMX_THREAD_MPI
372     /* Non-master threads have to wait for the detection to be done. */
373     if (PAR(cr))
374     {
375         MPI_Barrier(cr->mpi_comm_mysim);
376     }
377 #endif
378
379     /* inform the user about the settings */
380     if (SIMMASTER(cr) && bOMP)
381     {
382 #ifdef GMX_THREAD_MPI
383         const char *mpi_str="per tMPI thread";
384 #else
385         const char *mpi_str="per MPI process";
386 #endif
387
388         /* for group scheme we print PME threads info only */
389         if (bFullOmpSupport)
390         {
391             fprintf(stderr, "Using %d OpenMP thread%s %s\n",
392                     modth.gnth,modth.gnth > 1 ? "s" : "",
393                     cr->nnodes > 1 ? mpi_str : "");
394         }
395         if (bSepPME && modth.gnth_pme != modth.gnth)
396         {
397             fprintf(stderr, "Using %d OpenMP thread%s %s for PME\n",
398                     modth.gnth_pme,modth.gnth_pme > 1 ? "s" : "",
399                     cr->nnodes > 1 ? mpi_str : "");
400         }
401     }
402
403     /* detect and warn about oversubscription
404      * TODO: enable this for separate PME nodes as well! */
405     if (!bSepPME && cr->nodeid_intra == 0)
406     {
407         char sbuf[STRLEN], sbuf1[STRLEN], sbuf2[STRLEN];
408
409         if (modth.gnth*nppn > nthreads_hw_avail)
410         {
411             sprintf(sbuf, "threads");
412             sbuf1[0] = '\0';
413             sprintf(sbuf2, "O");
414 #ifdef GMX_MPI
415             if (modth.gnth == 1)
416             {
417 #ifdef GMX_THREAD_MPI
418                 sprintf(sbuf, "thread-MPI threads");
419 #else
420                 sprintf(sbuf, "MPI processes");
421                 sprintf(sbuf1, " per node");
422                 sprintf(sbuf2, "On node %d: o", cr->sim_nodeid);
423 #endif
424             }
425 #endif
426             md_print_warn(cr, fplog,
427                           "WARNING: %sversubscribing the available %d logical CPU cores%s with %d %s.\n"
428                           "         This will cause considerable performance loss!",
429                           sbuf2, nthreads_hw_avail, sbuf1, nppn*modth.gnth, sbuf);
430         }
431     }
432 }
433
434 int gmx_omp_nthreads_get(int mod)
435 {
436     if (mod < 0 || mod >= emntNR)
437     {
438         /* invalid module queried */
439         return -1;
440     }
441     else
442     {
443         return modth.nth[mod];
444     }
445 }