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