6ae2e122ff49fc04a425496b68c8bcee476f3623
[alexxy/gromacs.git] / src / gromacs / mdlib / gmx_omp_nthreads.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, 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 "gmxpre.h"
37
38 #include "gmx_omp_nthreads.h"
39
40 #include "config.h"
41
42 #include <cstdio>
43 #include <cstdlib>
44 #include <cstring>
45
46 #include "gromacs/gmxlib/network.h"
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/utility/cstringutil.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/gmxassert.h"
51 #include "gromacs/utility/gmxomp.h"
52 #include "gromacs/utility/logger.h"
53 #include "gromacs/utility/programcontext.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] = { "GMX_DEFAULT_NUM_THREADS should never be set",
71                                              "GMX_DOMDEC_NUM_THREADS",
72                                              "GMX_PAIRSEARCH_NUM_THREADS",
73                                              "GMX_NONBONDED_NUM_THREADS",
74                                              "GMX_LISTED_FORCES_NUM_THREADS",
75                                              "GMX_PME_NUM_THREADS",
76                                              "GMX_UPDATE_NUM_THREADS",
77                                              "GMX_VSITE_NUM_THREADS",
78                                              "GMX_LINCS_NUM_THREADS",
79                                              "GMX_SETTLE_NUM_THREADS" };
80
81 /** Names of the modules. */
82 static const char* mod_name[emntNR] = { "default",     "domain decomposition",
83                                         "pair search", "non-bonded",
84                                         "bonded",      "PME",
85                                         "update",      "LINCS",
86                                         "SETTLE" };
87
88 /** Number of threads for each algorithmic module.
89  *
90  *  File-scope global variable that gets set once in pick_module_nthreads()
91  *  and queried via gmx_omp_nthreads_get().
92  *
93  *  All fields are initialized to 0 which should result in errors if
94  *  the init call is omitted.
95  * */
96 static omp_module_nthreads_t modth = { 0, 0, { 0, 0, 0, 0, 0, 0, 0, 0, 0 }, FALSE };
97
98
99 /** Determine the number of threads for module \p mod.
100  *
101  *  \p m takes values form the module_nth_t enum and maps these to the
102  *  corresponding value in modth_env_var.
103  *
104  *  Each number of threads per module takes the default value unless
105  *  GMX_*_NUM_THERADS env var is set, case in which its value overrides
106  *  the default.
107  */
108 static void pick_module_nthreads(const gmx::MDLogger& mdlog, int m, gmx_bool bSepPME)
109 {
110     char* env;
111     int   nth;
112
113     const bool bOMP = GMX_OPENMP;
114
115     /* The default should never be set through a GMX_*_NUM_THREADS env var
116      * as it's always equal with gnth. */
117     if (m == emntDefault)
118     {
119         return;
120     }
121
122     /* check the environment variable */
123     if ((env = getenv(modth_env_var[m])) != nullptr)
124     {
125         sscanf(env, "%d", &nth);
126
127         if (!bOMP)
128         {
129             gmx_warning("%s=%d is set, but %s is compiled without OpenMP!", modth_env_var[m], nth,
130                         gmx::getProgramContext().displayName());
131         }
132
133         /* with the verlet codepath, when any GMX_*_NUM_THREADS env var is set,
134          * OMP_NUM_THREADS also has to be set */
135         if (getenv("OMP_NUM_THREADS") == nullptr)
136         {
137             gmx_warning(
138                     "%s=%d is set, the default number of threads also "
139                     "needs to be set with OMP_NUM_THREADS!",
140                     modth_env_var[m], nth);
141         }
142
143         /* only babble if we are really overriding with a different value */
144         if ((bSepPME && m == emntPME && nth != modth.gnth_pme) || (nth != modth.gnth))
145         {
146             GMX_LOG(mdlog.warning)
147                     .asParagraph()
148                     .appendTextFormatted("%s=%d set, overriding the default number of %s threads",
149                                          modth_env_var[m], nth, mod_name[m]);
150         }
151     }
152     else
153     {
154         /* pick the global PME node nthreads if we are setting the number
155          * of threads in separate PME nodes  */
156         nth = (bSepPME && m == emntPME) ? modth.gnth_pme : modth.gnth;
157     }
158
159     gmx_omp_nthreads_set(m, nth);
160 }
161
162 void gmx_omp_nthreads_read_env(const gmx::MDLogger& mdlog, int* nthreads_omp)
163 {
164     char*    env;
165     gmx_bool bCommandLineSetNthreadsOMP = *nthreads_omp > 0;
166     char     buffer[STRLEN];
167
168     GMX_RELEASE_ASSERT(nthreads_omp, "nthreads_omp must be a non-NULL pointer");
169
170     if ((env = getenv("OMP_NUM_THREADS")) != nullptr)
171     {
172         int nt_omp;
173
174         sscanf(env, "%d", &nt_omp);
175         if (nt_omp <= 0)
176         {
177             gmx_fatal(FARGS, "OMP_NUM_THREADS is invalid: '%s'", env);
178         }
179
180         if (bCommandLineSetNthreadsOMP && nt_omp != *nthreads_omp)
181         {
182             gmx_fatal(FARGS,
183                       "Environment variable OMP_NUM_THREADS (%d) and the number of threads "
184                       "requested on the command line (%d) have different values. Either omit one, "
185                       "or set them both to the same value.",
186                       nt_omp, *nthreads_omp);
187         }
188
189         /* Setting the number of OpenMP threads. */
190         *nthreads_omp = nt_omp;
191
192         /* Output the results */
193         sprintf(buffer,
194                 "\nThe number of OpenMP threads was set by environment variable OMP_NUM_THREADS to "
195                 "%d%s\n\n",
196                 nt_omp,
197                 bCommandLineSetNthreadsOMP ? " (and the command-line setting agreed with that)" : "");
198
199         /* This prints once per simulation for multi-simulations,
200          * which might help diagnose issues with inhomogenous
201          * cluster setups. */
202         GMX_LOG(mdlog.info).appendTextFormatted("%s", buffer);
203         if (debug)
204         {
205             /* This prints once per process for real MPI (i.e. once
206              * per debug file), and once per simulation for thread MPI
207              * (because of logic in the calling function). */
208             fputs(buffer, debug);
209         }
210     }
211 }
212
213 /*! \brief Helper function for parsing various input about the number
214     of OpenMP threads to use in various modules and deciding what to
215     do about it. */
216 static void manage_number_of_openmp_threads(const gmx::MDLogger& mdlog,
217                                             const t_commrec*     cr,
218                                             bool                 bOMP,
219                                             int                  nthreads_hw_avail,
220                                             int                  omp_nthreads_req,
221                                             int                  omp_nthreads_pme_req,
222                                             gmx_bool gmx_unused bThisNodePMEOnly,
223                                             int                 numRanksOnThisNode,
224                                             gmx_bool            bSepPME)
225 {
226     int   nth;
227     char* env;
228
229 #if GMX_THREAD_MPI
230     /* modth is shared among tMPI threads, so for thread safety, the
231      * detection is done on the master only. It is not thread-safe
232      * with multiple simulations, but that's anyway not supported by
233      * tMPI. */
234     if (!SIMMASTER(cr))
235     {
236         return;
237     }
238 #else
239     GMX_UNUSED_VALUE(cr);
240 #endif
241
242     if (modth.initialized)
243     {
244         /* Just return if the initialization has already been
245            done. This could only happen if gmx_omp_nthreads_init() has
246            already been called. */
247         return;
248     }
249
250     /* With full OpenMP support (verlet scheme) set the number of threads
251      * per process / default:
252      * - 1 if not compiled with OpenMP or
253      * - OMP_NUM_THREADS if the env. var is set, or
254      * - omp_nthreads_req = #of threads requested by the user on the mdrun
255      *   command line, otherwise
256      * - take the max number of available threads and distribute them
257      *   on the processes/tMPI threads.
258      * ~ The GMX_*_NUM_THREADS env var overrides the number of threads of
259      *   the respective module and it has to be used in conjunction with
260      *   OMP_NUM_THREADS.
261      *
262      * The number of PME threads is equal to:
263      * - 1 if not compiled with OpenMP or
264      * - GMX_PME_NUM_THREADS if defined, otherwise
265      * - OMP_NUM_THREADS if defined, otherwise
266      * - 1
267      */
268     nth = 1;
269     if ((env = getenv("OMP_NUM_THREADS")) != nullptr)
270     {
271         if (!bOMP && (std::strncmp(env, "1", 1) != 0))
272         {
273             gmx_warning("OMP_NUM_THREADS is set, but %s was compiled without OpenMP support!",
274                         gmx::getProgramContext().displayName());
275         }
276         else
277         {
278             nth = gmx_omp_get_max_threads();
279         }
280     }
281     else if (omp_nthreads_req > 0)
282     {
283         nth = omp_nthreads_req;
284     }
285     else if (bOMP)
286     {
287         /* max available threads per node */
288         nth = nthreads_hw_avail;
289
290         /* divide the threads among the MPI ranks */
291         if (nth >= numRanksOnThisNode)
292         {
293             nth /= numRanksOnThisNode;
294         }
295         else
296         {
297             nth = 1;
298         }
299     }
300
301     /* now we have the global values, set them:
302      * - 1 if not compiled with OpenMP
303      * - nth for the verlet scheme when compiled with OpenMP
304      */
305     if (bOMP)
306     {
307         modth.gnth = nth;
308     }
309     else
310     {
311         modth.gnth = 1;
312     }
313
314     if (bSepPME)
315     {
316         if (omp_nthreads_pme_req > 0)
317         {
318             modth.gnth_pme = omp_nthreads_pme_req;
319         }
320         else
321         {
322             modth.gnth_pme = nth;
323         }
324     }
325     else
326     {
327         modth.gnth_pme = 0;
328     }
329
330     /* now set the per-module values */
331     modth.nth[emntDefault] = modth.gnth;
332     pick_module_nthreads(mdlog, emntDomdec, bSepPME);
333     pick_module_nthreads(mdlog, emntPairsearch, bSepPME);
334     pick_module_nthreads(mdlog, emntNonbonded, bSepPME);
335     pick_module_nthreads(mdlog, emntBonded, bSepPME);
336     pick_module_nthreads(mdlog, emntPME, bSepPME);
337     pick_module_nthreads(mdlog, emntUpdate, bSepPME);
338     pick_module_nthreads(mdlog, emntVSITE, bSepPME);
339     pick_module_nthreads(mdlog, emntLINCS, bSepPME);
340     pick_module_nthreads(mdlog, emntSETTLE, bSepPME);
341
342     /* set the number of threads globally */
343     if (bOMP)
344     {
345 #if !GMX_THREAD_MPI
346         if (bThisNodePMEOnly)
347         {
348             gmx_omp_set_num_threads(modth.gnth_pme);
349         }
350         else
351 #endif /* GMX_THREAD_MPI */
352         {
353             gmx_omp_set_num_threads(nth);
354         }
355     }
356
357     modth.initialized = TRUE;
358 }
359
360 /*! \brief Report on the OpenMP settings that will be used */
361 static void reportOpenmpSettings(const gmx::MDLogger& mdlog, const t_commrec* cr, gmx_bool bOMP, gmx_bool bSepPME)
362 {
363 #if GMX_THREAD_MPI
364     const char* mpi_str = "per tMPI thread";
365 #else
366     const char* mpi_str = "per MPI process";
367 #endif
368     int nth_min, nth_max, nth_pme_min, nth_pme_max;
369
370     /* inform the user about the settings */
371     if (!bOMP)
372     {
373         return;
374     }
375
376 #if GMX_MPI
377     if (cr->nnodes > 1)
378     {
379         /* Get the min and max thread counts over the MPI ranks */
380         int buf_in[4], buf_out[4];
381
382         buf_in[0] = -modth.gnth;
383         buf_in[1] = modth.gnth;
384         buf_in[2] = -modth.gnth_pme;
385         buf_in[3] = modth.gnth_pme;
386
387         MPI_Allreduce(buf_in, buf_out, 4, MPI_INT, MPI_MAX, cr->mpi_comm_mysim);
388
389         nth_min     = -buf_out[0];
390         nth_max     = buf_out[1];
391         nth_pme_min = -buf_out[2];
392         nth_pme_max = buf_out[3];
393     }
394     else
395 #endif
396     {
397         nth_min     = modth.gnth;
398         nth_max     = modth.gnth;
399         nth_pme_min = modth.gnth_pme;
400         nth_pme_max = modth.gnth_pme;
401     }
402
403
404     if (nth_max == nth_min)
405     {
406         GMX_LOG(mdlog.warning)
407                 .appendTextFormatted("Using %d OpenMP thread%s %s", nth_min, nth_min > 1 ? "s" : "",
408                                      cr->nnodes > 1 ? mpi_str : "");
409     }
410     else
411     {
412         GMX_LOG(mdlog.warning).appendTextFormatted("Using %d - %d OpenMP threads %s", nth_min, nth_max, mpi_str);
413     }
414
415     if (bSepPME && (nth_pme_min != nth_min || nth_pme_max != nth_max))
416     {
417         if (nth_pme_max == nth_pme_min)
418         {
419             GMX_LOG(mdlog.warning)
420                     .appendTextFormatted("Using %d OpenMP thread%s %s for PME", nth_pme_min,
421                                          nth_pme_min > 1 ? "s" : "", cr->nnodes > 1 ? mpi_str : "");
422         }
423         else
424         {
425             GMX_LOG(mdlog.warning)
426                     .appendTextFormatted("Using %d - %d OpenMP threads %s for PME", nth_pme_min,
427                                          nth_pme_max, mpi_str);
428         }
429     }
430     GMX_LOG(mdlog.warning);
431 }
432
433 void gmx_omp_nthreads_init(const gmx::MDLogger& mdlog,
434                            t_commrec*           cr,
435                            int                  nthreads_hw_avail,
436                            int                  numRanksOnThisNode,
437                            int                  omp_nthreads_req,
438                            int                  omp_nthreads_pme_req,
439                            gmx_bool             bThisNodePMEOnly)
440 {
441     gmx_bool bSepPME;
442
443     const bool bOMP = GMX_OPENMP;
444
445     bSepPME = (thisRankHasDuty(cr, DUTY_PP) != thisRankHasDuty(cr, DUTY_PME));
446
447     manage_number_of_openmp_threads(mdlog, cr, bOMP, nthreads_hw_avail, omp_nthreads_req,
448                                     omp_nthreads_pme_req, bThisNodePMEOnly, numRanksOnThisNode, bSepPME);
449 #if GMX_THREAD_MPI
450     /* Non-master threads have to wait for the OpenMP management to be
451      * done, so that code elsewhere that uses OpenMP can be certain
452      * the setup is complete. */
453     if (PAR(cr))
454     {
455         MPI_Barrier(cr->mpi_comm_mysim);
456     }
457 #endif
458
459     reportOpenmpSettings(mdlog, cr, bOMP, bSepPME);
460 }
461
462 int gmx_omp_nthreads_get(int mod)
463 {
464     if (mod < 0 || mod >= emntNR)
465     {
466         /* invalid module queried */
467         return -1;
468     }
469     else
470     {
471         return modth.nth[mod];
472     }
473 }
474
475 void gmx_omp_nthreads_set(int mod, int nthreads)
476 {
477     /* Catch an attempt to set the number of threads on an invalid
478      * OpenMP module. */
479     GMX_RELEASE_ASSERT(mod >= 0 && mod < emntNR, "Trying to set nthreads on invalid OpenMP module");
480
481     modth.nth[mod] = nthreads;
482 }