Apply re-formatting to C++ in src/ tree.
[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 by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 #include "gmxpre.h"
38
39 #include "gmx_omp_nthreads.h"
40
41 #include "config.h"
42
43 #include <cstdio>
44 #include <cstdlib>
45 #include <cstring>
46
47 #include "gromacs/gmxlib/network.h"
48 #include "gromacs/mdtypes/commrec.h"
49 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/gmxassert.h"
52 #include "gromacs/utility/gmxomp.h"
53 #include "gromacs/utility/logger.h"
54 #include "gromacs/utility/programcontext.h"
55
56 /** Structure with the number of threads for each OpenMP multi-threaded
57  *  algorithmic module in mdrun. */
58 typedef struct
59 {
60     int gnth;     /**< Global num. of threads per PP or PP+PME process/tMPI thread. */
61     int gnth_pme; /**< Global num. of threads per PME only process/tMPI thread. */
62
63     int      nth[emntNR]; /**< Number of threads for each module, indexed with module_nth_t */
64     gmx_bool initialized; /**< TRUE if the module as been initialized. */
65 } omp_module_nthreads_t;
66
67 /** Names of environment variables to set the per module number of threads.
68  *
69  *  Indexed with the values of module_nth_t.
70  * */
71 static const char* modth_env_var[emntNR] = { "GMX_DEFAULT_NUM_THREADS should never be set",
72                                              "GMX_DOMDEC_NUM_THREADS",
73                                              "GMX_PAIRSEARCH_NUM_THREADS",
74                                              "GMX_NONBONDED_NUM_THREADS",
75                                              "GMX_LISTED_FORCES_NUM_THREADS",
76                                              "GMX_PME_NUM_THREADS",
77                                              "GMX_UPDATE_NUM_THREADS",
78                                              "GMX_VSITE_NUM_THREADS",
79                                              "GMX_LINCS_NUM_THREADS",
80                                              "GMX_SETTLE_NUM_THREADS" };
81
82 /** Names of the modules. */
83 static const char* mod_name[emntNR] = { "default",     "domain decomposition",
84                                         "pair search", "non-bonded",
85                                         "bonded",      "PME",
86                                         "update",      "LINCS",
87                                         "SETTLE" };
88
89 /** Number of threads for each algorithmic module.
90  *
91  *  File-scope global variable that gets set once in pick_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 \p mod.
101  *
102  *  \p 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 default.
108  */
109 static void pick_module_nthreads(const gmx::MDLogger& mdlog, int m, gmx_bool bSepPME)
110 {
111     char* env;
112     int   nth;
113
114     const bool bOMP = GMX_OPENMP;
115
116     /* The default should never be set through a GMX_*_NUM_THREADS env var
117      * as it's always equal with gnth. */
118     if (m == emntDefault)
119     {
120         return;
121     }
122
123     /* check the environment variable */
124     if ((env = getenv(modth_env_var[m])) != nullptr)
125     {
126         sscanf(env, "%d", &nth);
127
128         if (!bOMP)
129         {
130             gmx_warning("%s=%d is set, but %s is compiled without OpenMP!",
131                         modth_env_var[m],
132                         nth,
133                         gmx::getProgramContext().displayName());
134         }
135
136         /* with the verlet codepath, when any GMX_*_NUM_THREADS env var is set,
137          * OMP_NUM_THREADS also has to be set */
138         if (getenv("OMP_NUM_THREADS") == nullptr)
139         {
140             gmx_warning(
141                     "%s=%d is set, the default number of threads also "
142                     "needs to be set with OMP_NUM_THREADS!",
143                     modth_env_var[m],
144                     nth);
145         }
146
147         /* only babble if we are really overriding with a different value */
148         if ((bSepPME && m == emntPME && nth != modth.gnth_pme) || (nth != modth.gnth))
149         {
150             GMX_LOG(mdlog.warning)
151                     .asParagraph()
152                     .appendTextFormatted("%s=%d set, overriding the default number of %s threads",
153                                          modth_env_var[m],
154                                          nth,
155                                          mod_name[m]);
156         }
157     }
158     else
159     {
160         /* pick the global PME node nthreads if we are setting the number
161          * of threads in separate PME nodes  */
162         nth = (bSepPME && m == emntPME) ? modth.gnth_pme : modth.gnth;
163     }
164
165     gmx_omp_nthreads_set(m, nth);
166 }
167
168 void gmx_omp_nthreads_read_env(const gmx::MDLogger& mdlog, int* nthreads_omp)
169 {
170     char*    env;
171     gmx_bool bCommandLineSetNthreadsOMP = *nthreads_omp > 0;
172     char     buffer[STRLEN];
173
174     GMX_RELEASE_ASSERT(nthreads_omp, "nthreads_omp must be a non-NULL pointer");
175
176     if ((env = getenv("OMP_NUM_THREADS")) != nullptr)
177     {
178         int nt_omp;
179
180         sscanf(env, "%d", &nt_omp);
181         if (nt_omp <= 0)
182         {
183             gmx_fatal(FARGS, "OMP_NUM_THREADS is invalid: '%s'", env);
184         }
185
186         if (bCommandLineSetNthreadsOMP && nt_omp != *nthreads_omp)
187         {
188             gmx_fatal(FARGS,
189                       "Environment variable OMP_NUM_THREADS (%d) and the number of threads "
190                       "requested on the command line (%d) have different values. Either omit one, "
191                       "or set them both to the same value.",
192                       nt_omp,
193                       *nthreads_omp);
194         }
195
196         /* Setting the number of OpenMP threads. */
197         *nthreads_omp = nt_omp;
198
199         /* Output the results */
200         sprintf(buffer,
201                 "\nThe number of OpenMP threads was set by environment variable OMP_NUM_THREADS to "
202                 "%d%s\n\n",
203                 nt_omp,
204                 bCommandLineSetNthreadsOMP ? " (and the command-line setting agreed with that)" : "");
205
206         /* This prints once per simulation for multi-simulations,
207          * which might help diagnose issues with inhomogenous
208          * cluster setups. */
209         GMX_LOG(mdlog.info).appendTextFormatted("%s", buffer);
210         if (debug)
211         {
212             /* This prints once per process for real MPI (i.e. once
213              * per debug file), and once per simulation for thread MPI
214              * (because of logic in the calling function). */
215             fputs(buffer, debug);
216         }
217     }
218 }
219
220 /*! \brief Helper function for parsing various input about the number
221     of OpenMP threads to use in various modules and deciding what to
222     do about it. */
223 static void manage_number_of_openmp_threads(const gmx::MDLogger& mdlog,
224                                             const t_commrec*     cr,
225                                             bool                 bOMP,
226                                             int                  nthreads_hw_avail,
227                                             int                  omp_nthreads_req,
228                                             int                  omp_nthreads_pme_req,
229                                             gmx_bool gmx_unused bThisNodePMEOnly,
230                                             int                 numRanksOnThisNode,
231                                             gmx_bool            bSepPME)
232 {
233     int   nth;
234     char* env;
235
236 #if GMX_THREAD_MPI
237     /* modth is shared among tMPI threads, so for thread safety, the
238      * detection is done on the master only. It is not thread-safe
239      * with multiple simulations, but that's anyway not supported by
240      * tMPI. */
241     if (!SIMMASTER(cr))
242     {
243         return;
244     }
245 #else
246     GMX_UNUSED_VALUE(cr);
247 #endif
248
249     if (modth.initialized)
250     {
251         /* Just return if the initialization has already been
252            done. This could only happen if gmx_omp_nthreads_init() has
253            already been called. */
254         return;
255     }
256
257     /* With full OpenMP support (verlet scheme) set the number of threads
258      * per process / default:
259      * - 1 if not compiled with OpenMP or
260      * - OMP_NUM_THREADS if the env. var is set, or
261      * - omp_nthreads_req = #of threads requested by the user on the mdrun
262      *   command line, otherwise
263      * - take the max number of available threads and distribute them
264      *   on the processes/tMPI threads.
265      * ~ The GMX_*_NUM_THREADS env var overrides the number of threads of
266      *   the respective module and it has to be used in conjunction with
267      *   OMP_NUM_THREADS.
268      *
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")) != nullptr)
277     {
278         if (!bOMP && (std::strncmp(env, "1", 1) != 0))
279         {
280             gmx_warning("OMP_NUM_THREADS is set, but %s was compiled without OpenMP support!",
281                         gmx::getProgramContext().displayName());
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 (bOMP)
293     {
294         /* max available threads per node */
295         nth = nthreads_hw_avail;
296
297         /* divide the threads among the MPI ranks */
298         if (nth >= numRanksOnThisNode)
299         {
300             nth /= numRanksOnThisNode;
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
310      * - nth for the verlet scheme when compiled with OpenMP
311      */
312     if (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(mdlog, emntDomdec, bSepPME);
340     pick_module_nthreads(mdlog, emntPairsearch, bSepPME);
341     pick_module_nthreads(mdlog, emntNonbonded, bSepPME);
342     pick_module_nthreads(mdlog, emntBonded, bSepPME);
343     pick_module_nthreads(mdlog, emntPME, bSepPME);
344     pick_module_nthreads(mdlog, emntUpdate, bSepPME);
345     pick_module_nthreads(mdlog, emntVSITE, bSepPME);
346     pick_module_nthreads(mdlog, emntLINCS, bSepPME);
347     pick_module_nthreads(mdlog, emntSETTLE, bSepPME);
348
349     /* set the number of threads globally */
350     if (bOMP)
351     {
352 #if !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             gmx_omp_set_num_threads(nth);
361         }
362     }
363
364     modth.initialized = TRUE;
365 }
366
367 /*! \brief Report on the OpenMP settings that will be used */
368 static void reportOpenmpSettings(const gmx::MDLogger& mdlog, const t_commrec* cr, gmx_bool bOMP, gmx_bool bSepPME)
369 {
370 #if GMX_THREAD_MPI
371     const char* mpi_str = "per tMPI thread";
372 #else
373     const char* mpi_str = "per MPI process";
374 #endif
375     int nth_min, nth_max, nth_pme_min, nth_pme_max;
376
377     /* inform the user about the settings */
378     if (!bOMP)
379     {
380         return;
381     }
382
383 #if GMX_MPI
384     if (cr->nnodes > 1)
385     {
386         /* Get the min and max thread counts over the MPI ranks */
387         int buf_in[4], buf_out[4];
388
389         buf_in[0] = -modth.gnth;
390         buf_in[1] = modth.gnth;
391         buf_in[2] = -modth.gnth_pme;
392         buf_in[3] = modth.gnth_pme;
393
394         MPI_Allreduce(buf_in, buf_out, 4, MPI_INT, MPI_MAX, cr->mpi_comm_mysim);
395
396         nth_min     = -buf_out[0];
397         nth_max     = buf_out[1];
398         nth_pme_min = -buf_out[2];
399         nth_pme_max = buf_out[3];
400     }
401     else
402 #endif
403     {
404         nth_min     = modth.gnth;
405         nth_max     = modth.gnth;
406         nth_pme_min = modth.gnth_pme;
407         nth_pme_max = modth.gnth_pme;
408     }
409
410
411     if (nth_max == nth_min)
412     {
413         GMX_LOG(mdlog.warning)
414                 .appendTextFormatted("Using %d OpenMP thread%s %s",
415                                      nth_min,
416                                      nth_min > 1 ? "s" : "",
417                                      cr->nnodes > 1 ? mpi_str : "");
418     }
419     else
420     {
421         GMX_LOG(mdlog.warning).appendTextFormatted("Using %d - %d OpenMP threads %s", nth_min, nth_max, mpi_str);
422     }
423
424     if (bSepPME && (nth_pme_min != nth_min || nth_pme_max != nth_max))
425     {
426         if (nth_pme_max == nth_pme_min)
427         {
428             GMX_LOG(mdlog.warning)
429                     .appendTextFormatted("Using %d OpenMP thread%s %s for PME",
430                                          nth_pme_min,
431                                          nth_pme_min > 1 ? "s" : "",
432                                          cr->nnodes > 1 ? mpi_str : "");
433         }
434         else
435         {
436             GMX_LOG(mdlog.warning)
437                     .appendTextFormatted(
438                             "Using %d - %d OpenMP threads %s for PME", nth_pme_min, nth_pme_max, mpi_str);
439         }
440     }
441     GMX_LOG(mdlog.warning);
442 }
443
444 void gmx_omp_nthreads_init(const gmx::MDLogger& mdlog,
445                            t_commrec*           cr,
446                            int                  nthreads_hw_avail,
447                            int                  numRanksOnThisNode,
448                            int                  omp_nthreads_req,
449                            int                  omp_nthreads_pme_req,
450                            gmx_bool             bThisNodePMEOnly)
451 {
452     gmx_bool bSepPME;
453
454     const bool bOMP = GMX_OPENMP;
455
456     bSepPME = (thisRankHasDuty(cr, DUTY_PP) != thisRankHasDuty(cr, DUTY_PME));
457
458     manage_number_of_openmp_threads(mdlog,
459                                     cr,
460                                     bOMP,
461                                     nthreads_hw_avail,
462                                     omp_nthreads_req,
463                                     omp_nthreads_pme_req,
464                                     bThisNodePMEOnly,
465                                     numRanksOnThisNode,
466                                     bSepPME);
467 #if GMX_THREAD_MPI
468     /* Non-master threads have to wait for the OpenMP management to be
469      * done, so that code elsewhere that uses OpenMP can be certain
470      * the setup is complete. */
471     if (PAR(cr))
472     {
473         MPI_Barrier(cr->mpi_comm_mysim);
474     }
475 #endif
476
477     reportOpenmpSettings(mdlog, cr, bOMP, bSepPME);
478 }
479
480 int gmx_omp_nthreads_get(int mod)
481 {
482     if (mod < 0 || mod >= emntNR)
483     {
484         /* invalid module queried */
485         return -1;
486     }
487     else
488     {
489         return modth.nth[mod];
490     }
491 }
492
493 void gmx_omp_nthreads_set(int mod, int nthreads)
494 {
495     /* Catch an attempt to set the number of threads on an invalid
496      * OpenMP module. */
497     GMX_RELEASE_ASSERT(mod >= 0 && mod < emntNR, "Trying to set nthreads on invalid OpenMP module");
498
499     modth.nth[mod] = nthreads;
500 }