Merge "Extended genbox to insert molecules at given positions"
[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_VSITE_NUM_THREADS",
76     "GMX_LINCS_NUM_THREADS", "GMX_SETTLE_NUM_THREADS"
77 };
78
79 /** Names of the modules. */
80 static const char *mod_name[emntNR] =
81 {
82     "default", "domain decomposition", "pair search", "non-bonded",
83     "bonded", "PME", "update", "LINCS", "SETTLE"
84 };
85
86 /** Number of threads for each algorithmic module.
87  *
88  *  File-scope global variable that gets set once in \init_module_nthreads
89  *  and queried via gmx_omp_nthreads_get.
90  *
91  *  All fields are initialized to 0 which should result in errors if
92  *  the init call is omitted.
93  * */
94 static omp_module_nthreads_t modth = { 0, 0, {0, 0, 0, 0, 0, 0, 0, 0, 0}, FALSE};
95
96
97 /** Determine the number of threads for module \mod.
98  *
99  *  \m takes values form the module_nth_t enum and maps these to the
100  *  corresponding value in modth_env_var.
101  *
102  *  Each number of threads per module takes the default value unless
103  *  GMX_*_NUM_THERADS env var is set, case in which its value overrides
104  *  the deafult.
105  *
106  *  The "group" scheme supports OpenMP only in PME and in thise case all but
107  *  the PME nthread values default to 1.
108  */
109 static int pick_module_nthreads(FILE *fplog, int m,
110                                 gmx_bool bSimMaster,
111                                 gmx_bool bFullOmpSupport,
112                                 gmx_bool bSepPME)
113 {
114     char *env;
115     int  nth;
116     char sbuf[STRLEN];
117     gmx_bool bOMP;
118
119 #ifdef GMX_OPENMP
120     bOMP = TRUE;
121 #else
122     bOMP = FALSE;
123 #endif /* GMX_OPENMP */
124
125     /* The default should never be set through a GMX_*_NUM_THREADS env var
126      * as it's always equal with gnth. */
127     if (m == emntDefault)
128     {
129         return modth.nth[emntDefault];
130     }
131
132     /* check the environment variable */
133     if ((env = getenv(modth_env_var[m])) != NULL)
134     {
135         sscanf(env, "%d", &nth);
136
137         if (!bOMP)
138         {
139             gmx_warning("%s=%d is set, but %s is compiled without OpenMP!",
140                         modth_env_var[m], nth, ShortProgram());
141         }
142
143         /* with the verlet codepath, when any GMX_*_NUM_THREADS env var is set,
144          * OMP_NUM_THREADS also has to be set */
145         if (bFullOmpSupport && getenv("OMP_NUM_THREADS") == NULL)
146         {
147             gmx_fatal(FARGS, "%s=%d is set, the default number of threads also "
148                       "needs to be set with OMP_NUM_THREADS!",
149                       modth_env_var[m], nth);
150         }
151
152         /* with the group scheme warn if any env var except PME is set */
153         if (!bFullOmpSupport)
154         {
155             if (m != emntPME)
156             {
157                 gmx_warning("%s=%d is set, but OpenMP multithreading is not "
158                             "supported in %s!",
159                             modth_env_var[m], nth, mod_name[m]);
160                 nth = 1;
161             }
162         }
163
164         /* only babble if we are really overriding with a different value */
165         if ((bSepPME && m == emntPME && nth != modth.gnth_pme) || (nth != modth.gnth))
166         {
167             sprintf(sbuf, "%s=%d set, overriding the default number of %s threads",
168                     modth_env_var[m], nth, mod_name[m]);
169             if (bSimMaster)
170             {
171                 fprintf(stderr, "\n%s\n", sbuf);
172             }
173             if (fplog)
174             {
175                 fprintf(fplog, "%s\n", sbuf);
176             }
177         }
178     }
179     else
180     {
181         /* pick the global PME node nthreads if we are setting the number
182          * of threads in separate PME nodes  */
183         nth = (bSepPME && m == emntPME) ? modth.gnth_pme : modth.gnth;
184     }
185
186     return modth.nth[m] = nth;
187 }
188
189 void gmx_omp_nthreads_read_env(int *nthreads_omp)
190 {
191     char *env;
192
193     assert(nthreads_omp);
194
195     if ((env = getenv("OMP_NUM_THREADS")) != NULL)
196     {
197         int nt_omp;
198
199         sscanf(env,"%d",&nt_omp);
200         if (nt_omp <= 0)
201         {
202             gmx_fatal(FARGS,"OMP_NUM_THREADS is invalid: '%s'",env);
203         }
204
205         if (*nthreads_omp > 0 && nt_omp != *nthreads_omp)
206         {
207             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);
208         }
209
210         /* Setting the number of OpenMP threads.
211          * NOTE: with tMPI this function is only called on the master node,
212          * but with MPI on all nodes which means lots of messages on stderr.
213          */
214         fprintf(stderr,"Getting the number of OpenMP threads from OMP_NUM_THREADS: %d\n",nt_omp);
215         *nthreads_omp = nt_omp;
216     }
217 }
218
219 void gmx_omp_nthreads_init(FILE *fplog, t_commrec *cr,
220                            int nthreads_hw_avail,
221                            int omp_nthreads_req,
222                            int omp_nthreads_pme_req,
223                            gmx_bool bThisNodePMEOnly,
224                            gmx_bool bFullOmpSupport)
225 {
226     int  nth, nth_pmeonly, gmx_maxth, nppn;
227     char *env;
228     gmx_bool bSepPME, bOMP;
229
230 #ifdef GMX_OPENMP
231     bOMP = TRUE;
232 #else
233     bOMP = FALSE;
234 #endif /* GMX_OPENMP */
235
236     /* number of MPI processes/threads per physical node */
237     nppn = cr->nrank_intranode;
238
239     bSepPME = ( (cr->duty & DUTY_PP) && !(cr->duty & DUTY_PME)) ||
240               (!(cr->duty & DUTY_PP) &&  (cr->duty & DUTY_PME));
241
242 #ifdef GMX_THREAD_MPI
243     /* modth is shared among tMPI threads, so for thread safety do the
244      * detection is done on the master only. It is not thread-safe with
245      * multiple simulations, but that's anyway not supported by tMPI. */
246     if (SIMMASTER(cr))
247 #endif
248     {
249         /* just return if the initialization has already been done */
250         if (modth.initialized)
251         {
252             return;
253         }
254
255         /* With full OpenMP support (verlet scheme) set the number of threads
256          * per process / default:
257          * - 1 if not compiled with OpenMP or
258          * - OMP_NUM_THREADS if the env. var is set, or
259          * - omp_nthreads_req = #of threads requested by the user on the mdrun
260          *   command line, otherwise
261          * - take the max number of available threads and distribute them
262          *   on the processes/tMPI threads.
263          * ~ The GMX_*_NUM_THREADS env var overrides the number of threads of
264          *   the respective module and it has to be used in conjunction with
265          *   OMP_NUM_THREADS.
266          *
267          * With the group scheme OpenMP multithreading is only supported in PME,
268          * for all other modules nthreads is set to 1.
269          * The number of PME threads is equal to:
270          * - 1 if not compiled with OpenMP or
271          * - GMX_PME_NUM_THREADS if defined, otherwise
272          * - OMP_NUM_THREADS if defined, otherwise
273          * - 1
274          */
275         nth = 1;
276         if ((env = getenv("OMP_NUM_THREADS")) != NULL)
277         {
278             if (!bOMP && (strncmp(env, "1", 1) != 0))
279             {
280                 gmx_warning("OMP_NUM_THREADS is set, but %s was compiled without OpenMP support!",
281                             ShortProgram());
282             }
283             else
284             {
285                 nth = gmx_omp_get_max_threads();
286             }
287         }
288         else if (omp_nthreads_req > 0)
289         {
290             nth = omp_nthreads_req;
291         }
292         else if (bFullOmpSupport && bOMP)
293         {
294             /* max available threads per node */
295             nth = nthreads_hw_avail;
296
297             /* divide the threads among the MPI processes/tMPI threads */
298             if (nth >= nppn)
299             {
300                 nth /= nppn;
301             }
302             else
303             {
304                 nth = 1;
305             }
306         }
307
308         /* now we have the global values, set them:
309          * - 1 if not compiled with OpenMP and for the group scheme
310          * - nth for the verlet scheme when compiled with OpenMP
311          */
312         if (bFullOmpSupport && bOMP)
313         {
314             modth.gnth = nth;
315         }
316         else
317         {
318             modth.gnth = 1;
319         }
320
321         if (bSepPME)
322         {
323             if (omp_nthreads_pme_req > 0)
324             {
325                 modth.gnth_pme = omp_nthreads_pme_req;
326             }
327             else
328             {
329                 modth.gnth_pme = nth;
330             }
331         }
332         else
333         {
334             modth.gnth_pme = 0;
335         }
336
337         /* now set the per-module values */
338         modth.nth[emntDefault] = modth.gnth;
339         pick_module_nthreads(fplog, emntDomdec, SIMMASTER(cr), bFullOmpSupport, bSepPME);
340         pick_module_nthreads(fplog, emntPairsearch, SIMMASTER(cr), bFullOmpSupport, bSepPME);
341         pick_module_nthreads(fplog, emntNonbonded, SIMMASTER(cr), bFullOmpSupport, bSepPME);
342         pick_module_nthreads(fplog, emntBonded, SIMMASTER(cr), bFullOmpSupport, bSepPME);
343         pick_module_nthreads(fplog, emntPME, SIMMASTER(cr), bFullOmpSupport, bSepPME);
344         pick_module_nthreads(fplog, emntUpdate, SIMMASTER(cr), bFullOmpSupport, bSepPME);
345         pick_module_nthreads(fplog, emntVSITE, SIMMASTER(cr), bFullOmpSupport, bSepPME);
346         pick_module_nthreads(fplog, emntLINCS, SIMMASTER(cr), bFullOmpSupport, bSepPME);
347         pick_module_nthreads(fplog, emntSETTLE, SIMMASTER(cr), bFullOmpSupport, bSepPME);
348
349         /* set the number of threads globally */
350         if (bOMP)
351         {
352 #ifndef GMX_THREAD_MPI
353             if (bThisNodePMEOnly)
354             {
355                 gmx_omp_set_num_threads(modth.gnth_pme);
356             }
357             else
358 #endif /* GMX_THREAD_MPI */
359             {
360                 if (bFullOmpSupport)
361                 {
362                     gmx_omp_set_num_threads(nth);
363                 }
364                 else
365                 {
366                     gmx_omp_set_num_threads(1);
367                 }
368             }
369         }
370
371         modth.initialized = TRUE;
372     }
373 #ifdef GMX_THREAD_MPI
374     /* Non-master threads have to wait for the detection to be done. */
375     if (PAR(cr))
376     {
377         MPI_Barrier(cr->mpi_comm_mysim);
378     }
379 #endif
380
381     /* inform the user about the settings */
382     if (SIMMASTER(cr) && bOMP)
383     {
384 #ifdef GMX_THREAD_MPI
385         const char *mpi_str="per tMPI thread";
386 #else
387         const char *mpi_str="per MPI process";
388 #endif
389
390         /* for group scheme we print PME threads info only */
391         if (bFullOmpSupport)
392         {
393             fprintf(stderr, "Using %d OpenMP thread%s %s\n",
394                     modth.gnth,modth.gnth > 1 ? "s" : "",
395                     cr->nnodes > 1 ? mpi_str : "");
396         }
397         if (bSepPME && modth.gnth_pme != modth.gnth)
398         {
399             fprintf(stderr, "Using %d OpenMP thread%s %s for PME\n",
400                     modth.gnth_pme,modth.gnth_pme > 1 ? "s" : "",
401                     cr->nnodes > 1 ? mpi_str : "");
402         }
403     }
404
405     /* detect and warn about oversubscription
406      * TODO: enable this for separate PME nodes as well! */
407     if (!bSepPME && cr->rank_pp_intranode == 0)
408     {
409         char sbuf[STRLEN], sbuf1[STRLEN], sbuf2[STRLEN];
410
411         if (modth.gnth*nppn > nthreads_hw_avail)
412         {
413             sprintf(sbuf, "threads");
414             sbuf1[0] = '\0';
415             sprintf(sbuf2, "O");
416 #ifdef GMX_MPI
417             if (modth.gnth == 1)
418             {
419 #ifdef GMX_THREAD_MPI
420                 sprintf(sbuf, "thread-MPI threads");
421 #else
422                 sprintf(sbuf, "MPI processes");
423                 sprintf(sbuf1, " per node");
424                 sprintf(sbuf2, "On node %d: o", cr->sim_nodeid);
425 #endif
426             }
427 #endif
428             md_print_warn(cr, fplog,
429                           "WARNING: %sversubscribing the available %d logical CPU cores%s with %d %s.\n"
430                           "         This will cause considerable performance loss!",
431                           sbuf2, nthreads_hw_avail, sbuf1, nppn*modth.gnth, sbuf);
432         }
433     }
434 }
435
436 int gmx_omp_nthreads_get(int mod)
437 {
438     if (mod < 0 || mod >= emntNR)
439     {
440         /* invalid module queried */
441         return -1;
442     }
443     else
444     {
445         return modth.nth[mod];
446     }
447 }