e7dd9907b5256e44d1291f4b6926f7e4ff80ffac
[alexxy/gromacs.git] / src / gromacs / mdlib / gmx_omp_nthreads.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, 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 #ifndef GMX_OMP_NTHREADS_H
37 #define GMX_OMP_NTHREADS_H
38
39 #include <stdio.h>
40
41 #include "gromacs/utility/basedefinitions.h"
42
43 struct t_commrec;
44
45 namespace gmx
46 {
47 class MDLogger;
48 }
49
50 /** Enum values corresponding to multithreaded algorithmic modules. */
51 typedef enum module_nth
52 {
53     /* Default is meant to be used in OMP regions outside the named
54      * algorithmic modules listed below. */
55     emntDefault, emntDomdec, emntPairsearch, emntNonbonded,
56     emntBonded, emntPME,  emntUpdate, emntVSITE, emntLINCS, emntSETTLE,
57     emntNR
58 } module_nth_t;
59
60 /*! \brief
61  * Initializes the per-module thread count.
62  *
63  * It is compatible with tMPI, thread-safety is ensured (for the features
64  * available with tMPI).
65  * This function should caled only once during the initialization of mdrun. */
66 void gmx_omp_nthreads_init(const gmx::MDLogger &fplog, t_commrec *cr,
67                            int nthreads_hw_avail,
68                            int numRanksOnThisNode,
69                            int omp_nthreads_req,
70                            int omp_nthreads_pme_req,
71                            gmx_bool bCurrNodePMEOnly,
72                            gmx_bool bFullOmpSupport);
73
74 /*! \brief
75  * Returns the number of threads to be used in the given module \p mod. */
76 int gmx_omp_nthreads_get(int mod);
77
78 /*! \brief
79  * Returns the number of threads to be used in the given module \p mod for simple rvec operations.
80  *
81  * When the, potentially, parallel task only consists of a loop of clear_rvec
82  * or rvec_inc for nrvec elements, the OpenMP overhead might be higher than
83  * the reduction in computional cost due to parallelization. This routine
84  * returns 1 when the overhead is expected to be higher than the gain.
85  */
86 static inline int gmx_omp_nthreads_get_simple_rvec_task(int mod, int nrvec)
87 {
88     /* There can be a relatively large overhead to an OpenMP parallel for loop.
89      * This overhead increases, slowly, with the numbe of threads used.
90      * The computational gain goes as 1/#threads. The two effects combined
91      * lead to a cross-over point for a (non-)parallel loop at loop count
92      * that is not strongly dependent on the thread count.
93      * Note that a (non-)parallel loop can have benefit later in the code
94      * due to generating more cache hits, depending on how the next lask
95      * that accesses the same data is (not) parallelized over threads.
96      *
97      * A value of 2000 is the switch-over point for Haswell without
98      * hyper-threading. With hyper-threading it is about a factor 1.5 higher.
99      */
100     const int nrvec_omp = 2000;
101
102     if (nrvec < nrvec_omp)
103     {
104         return 1;
105     }
106     else
107     {
108         return gmx_omp_nthreads_get(mod);
109     }
110 }
111
112 /*! \brief Sets the number of threads to be used in module.
113  *
114  * Intended for use in testing. */
115 void gmx_omp_nthreads_set(int mod, int nthreads);
116
117 /*! \brief
118  * Read the OMP_NUM_THREADS env. var. and check against the value set on the
119  * command line. */
120 void gmx_omp_nthreads_read_env(const gmx::MDLogger &mdlog,
121                                int                 *nthreads_omp);
122
123 #endif