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