Introduce type-safe GPU emulation variables
[gromacs.git] / src / gromacs / mdlib / nb_verlet.h
blob9d744b7f5924f523805e232cf407a1b499fccf9e
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2017, 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.
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.
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.
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.
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.
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.
36 // FIXME: remove the "__" prefix in front of the group def when we move the
37 // nonbonded code into separate dir.
39 /*! \libinternal \defgroup __module_nb_verlet Short-range non-bonded interaction module
40 * \ingroup group_mdrun
42 * \brief Computes forces and energies for short-range pair-interactions
43 * based on the Verlet algorithm. The algorithm uses pair-lists generated
44 * at fixed intervals as well as various flavors of pair interaction kernels
45 * implemented for a wide range of CPU and GPU architectures.
47 * The module includes support for flavors of Coulomb and Lennard-Jones interaction
48 * treatment implemented for a large range of SIMD instruction sets for CPU
49 * architectures as well as in CUDA and OpenCL for GPU architectures.
50 * Additionally there is a reference CPU non-SIMD and a reference CPU
51 * for GPU pair-list setup interaction kernel.
53 * The implementation of the kernels is based on the cluster non-bonded algorithm
54 * which in the code is referred to as the NxN algorithms ("nbnxn_" prefix);
55 * for details of the algorithm see DOI:10.1016/j.cpc.2013.06.003.
57 * Algorithmically, the non-bonded computation has two different modes:
58 * A "classical" mode: generate a list every nstlist steps containing at least
59 * all atom pairs up to a distance of rlistOuter and compute pair interactions
60 * for all pairs that are within the interaction cut-off.
61 * A "dynamic pruning" mode: generate an "outer-list" up to cut-off rlistOuter
62 * every nstlist steps and prune the outer-list using a cut-off of rlistInner
63 * every nstlistPrune steps to obtain a, smaller, "inner-list". This
64 * results in fewer interaction computations and allows for a larger nstlist.
65 * On a GPU, this dynamic pruning is performed in a rolling fashion, pruning
66 * only a sub-part of the list each (second) step. This way it can often
67 * overlap with integration and constraints on the CPU.
68 * Currently a simple heuristic determines which mode will be used.
70 * TODO: add a summary list and brief descriptions of the different submodules:
71 * search, CPU kernels, GPU glue code + kernels.
73 * \author Berk Hess <hess@kth.se>
74 * \author Szilárd Páll <pall.szilard@gmail.com>
75 * \author Mark Abraham <mark.j.abraham@gmail.com>
76 * \author Anca Hamuraru <anca@streamcomputing.eu>
77 * \author Teemu Virolainen <teemu@streamcomputing.eu>
78 * \author Dimitrios Karkoulis <dimitris.karkoulis@gmail.com>
80 * TODO: add more authors!
83 /*! \libinternal \file
85 * \brief This file contains the public interface of the non-bonded Verlet module
86 * that implements the NxN cluster non-bonded algorithm to efficiently compute
87 * pair forces.
90 * \author Berk Hess <hess@kth.se>
91 * \author Szilárd Páll <pall.szilard@gmail.com>
93 * \inlibraryapi
94 * \ingroup __module_nb_verlet
98 #ifndef NB_VERLET_H
99 #define NB_VERLET_H
101 #include <memory>
103 #include "gromacs/mdlib/nbnxn_gpu_types.h"
104 #include "gromacs/mdlib/nbnxn_pairlist.h"
106 //! Help pass GPU-emulation parameters with type safety.
107 enum class EmulateGpuNonbonded : bool
109 //! Do not emulate GPUs.
111 //! Do emulate GPUs.
115 #ifdef __cplusplus
116 extern "C" {
117 #endif
120 /*! \brief Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation */
121 typedef enum
123 nbnxnkNotSet = 0,
124 nbnxnk4x4_PlainC,
125 nbnxnk4xN_SIMD_4xN,
126 nbnxnk4xN_SIMD_2xNN,
127 nbnxnk8x8x8_GPU,
128 nbnxnk8x8x8_PlainC,
129 nbnxnkNR
130 } nbnxn_kernel_type;
132 /*! \brief Return a string identifying the kernel type.
134 * \param [in] kernel_type nonbonded kernel types, takes values from the nbnxn_kernel_type enum
135 * \returns a string identifying the kernel corresponding to the type passed as argument
137 const char *lookup_nbnxn_kernel_name(int kernel_type);
139 /*! \brief Ewald exclusion types */
140 enum {
141 ewaldexclTable, ewaldexclAnalytical
144 /*! \brief Atom locality indicator: local, non-local, all.
146 * Used for calls to:
147 * gridding, pair-search, force calculation, x/f buffer operations
148 * */
149 enum {
150 eatLocal = 0, eatNonlocal = 1, eatAll
153 /*! \brief Tests for local atom range */
154 #define LOCAL_A(x) ((x) == eatLocal)
155 /*! \brief Tests for non-local atom range */
156 #define NONLOCAL_A(x) ((x) == eatNonlocal)
157 /*! \brief Tests for either local or non-local atom range */
158 #define LOCAL_OR_NONLOCAL_A(x) (LOCAL_A(x) || NONLOCAL_A(x))
160 /*! \brief Interaction locality indicator
162 * Used in pair-list search/calculations in the following manner:
163 * - local interactions require local atom data and affect local output only;
164 * - non-local interactions require both local and non-local atom data and
165 * affect both local- and non-local output.
167 enum {
168 eintLocal = 0, eintNonlocal = 1
171 /*! \brief Tests for local interaction indicator */
172 #define LOCAL_I(x) ((x) == eintLocal)
173 /*! \brief Tests for non-local interaction indicator */
174 #define NONLOCAL_I(x) ((x) == eintNonlocal)
176 /*! \brief Flag to tell the nonbonded kernels whether to clear the force output buffers */
177 enum {
178 enbvClearFNo, enbvClearFYes
181 /*! \libinternal
182 * \brief Non-bonded interaction group data structure. */
183 typedef struct nonbonded_verlet_group_t {
184 nbnxn_pairlist_set_t nbl_lists; /**< pair list(s) */
185 nbnxn_atomdata_t *nbat; /**< atom data */
186 int kernel_type; /**< non-bonded kernel - see enum above */
187 int ewald_excl; /**< Ewald exclusion - see enum above */
188 } nonbonded_verlet_group_t;
190 /*! \libinternal
191 * \brief Top-level non-bonded data structure for the Verlet-type cut-off scheme. */
192 typedef struct nonbonded_verlet_t {
193 std::unique_ptr<NbnxnListParameters> listParams; /**< Parameters for the search and list pruning setup */
194 nbnxn_search_t nbs; /**< n vs n atom pair searching data */
195 int ngrp; /**< number of interaction groups */
196 nonbonded_verlet_group_t grp[2]; /**< local and non-local interaction group */
198 gmx_bool bUseGPU; /**< TRUE when non-bonded interactions are computed on a physical GPU */
199 EmulateGpuNonbonded emulateGpu; /**< true when non-bonded interactions are computed on the CPU using GPU-style pair lists */
200 gmx_nbnxn_gpu_t *gpu_nbv; /**< pointer to GPU nb verlet data */
201 int min_ci_balanced; /**< pair list balancing parameter
202 used for the 8x8x8 GPU kernels */
203 } nonbonded_verlet_t;
205 /*! \brief Getter for bUseGPU */
206 gmx_bool
207 usingGpu(nonbonded_verlet_t *nbv);
209 #ifdef __cplusplus
211 #endif
213 #endif /* NB_VERLET_H */