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