2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,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 * \brief Define CUDA implementation of nbnxn_gpu.h
38 * \author Szilard Pall <pall.szilard@gmail.com>
47 #include "gromacs/mdlib/nbnxn_gpu.h"
54 #include "nbnxn_cuda.h"
56 #include "gromacs/gpu_utils/cudautils.cuh"
57 #include "gromacs/mdlib/force_flags.h"
58 #include "gromacs/mdlib/nb_verlet.h"
59 #include "gromacs/mdlib/nbnxn_gpu_common.h"
60 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
61 #include "gromacs/mdlib/nbnxn_pairlist.h"
62 #include "gromacs/timing/gpu_timing.h"
63 #include "gromacs/utility/cstringutil.h"
64 #include "gromacs/utility/gmxassert.h"
66 #include "nbnxn_cuda_types.h"
69 * Texture references are created at compile-time and need to be declared
70 * at file scope as global variables (see http://docs.nvidia.com/cuda/cuda-c-programming-guide/index.html#texture-reference-api).
71 * The texture references below are used in two translation units;
72 * we declare them here along the kernels that use them (when compiling legacy Fermi kernels),
73 * and provide getters (see below) used by the data_mgmt module where the
74 * textures are bound/unbound.
75 * (In principle we could do it the other way arond, but that would likely require
76 * device linking and we'd rather avoid technical hurdles.)
78 /*! Texture reference for LJ C6/C12 parameters; bound to cu_nbparam_t.nbfp */
79 texture<float, 1, cudaReadModeElementType> nbfp_texref;
81 /*! Texture reference for LJ-PME parameters; bound to cu_nbparam_t.nbfp_comb */
82 texture<float, 1, cudaReadModeElementType> nbfp_comb_texref;
84 /*! Texture reference for Ewald coulomb force table; bound to cu_nbparam_t.coulomb_tab */
85 texture<float, 1, cudaReadModeElementType> coulomb_tab_texref;
88 /***** The kernel declarations/definitions come here *****/
90 /* Top-level kernel declaration generation: will generate through multiple
91 * inclusion the following flavors for all kernel declarations:
92 * - force-only output;
93 * - force and energy output;
94 * - force-only with pair list pruning;
95 * - force and energy output with pair list pruning.
97 #define FUNCTION_DECLARATION_ONLY
99 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
100 /** Force & energy **/
101 #define CALC_ENERGIES
102 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
105 /*** Pair-list pruning kernels ***/
108 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
109 /** Force & energy **/
110 #define CALC_ENERGIES
111 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
115 /* Prune-only kernels */
116 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_pruneonly.cuh"
117 #undef FUNCTION_DECLARATION_ONLY
119 /* Now generate the function definitions if we are using a single compilation unit. */
120 #if GMX_CUDA_NB_SINGLE_COMPILATION_UNIT
121 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_F_noprune.cu"
122 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_F_prune.cu"
123 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_VF_noprune.cu"
124 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_VF_prune.cu"
125 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_pruneonly.cu"
127 /* Prevent compilation in multiple compilation unit mode for CC 2.x. Although we have
128 * build-time checks to prevent this, the user could manually tweaks nvcc flags
129 * which would lead to buggy kernels getting compiled.
131 #if GMX_PTX_ARCH > 0 && GMX_PTX_ARCH <= 210 && !defined(__clang__)
132 #error Due to an CUDA nvcc compiler bug, the CUDA non-bonded module can not be compiled with multiple compilation units for CC 2.x devices. If you have changed the nvcc flags manually, either use the GMX_CUDA_TARGET_* variables instead or set GMX_CUDA_NB_SINGLE_COMPILATION_UNIT=ON CMake option.
134 #endif /* GMX_CUDA_NB_SINGLE_COMPILATION_UNIT */
137 /*! Nonbonded kernel function pointer type */
138 typedef void (*nbnxn_cu_kfunc_ptr_t)(const cu_atomdata_t,
143 /*********************************/
145 /* XXX switch between chevron and cudaLaunch (supported only in CUDA >=7.0)
146 -- only for benchmarking purposes */
147 static const bool bUseCudaLaunchKernel =
148 (GMX_CUDA_VERSION >= 7000) && (getenv("GMX_DISABLE_CUDALAUNCH") == NULL);
150 /* XXX always/never run the energy/pruning kernels -- only for benchmarking purposes */
151 static bool always_ener = (getenv("GMX_GPU_ALWAYS_ENER") != NULL);
152 static bool never_ener = (getenv("GMX_GPU_NEVER_ENER") != NULL);
153 static bool always_prune = (getenv("GMX_GPU_ALWAYS_PRUNE") != NULL);
156 /*! Returns the number of blocks to be used for the nonbonded GPU kernel. */
157 static inline int calc_nb_kernel_nblock(int nwork_units, const gmx_device_info_t *dinfo)
162 /* CUDA does not accept grid dimension of 0 (which can happen e.g. with an
163 empty domain) and that case should be handled before this point. */
164 assert(nwork_units > 0);
166 max_grid_x_size = dinfo->prop.maxGridSize[0];
168 /* do we exceed the grid x dimension limit? */
169 if (nwork_units > max_grid_x_size)
171 gmx_fatal(FARGS, "Watch out, the input system is too large to simulate!\n"
172 "The number of nonbonded work units (=number of super-clusters) exceeds the"
173 "maximum grid size in x dimension (%d > %d)!", nwork_units, max_grid_x_size);
180 /* Constant arrays listing all kernel function pointers and enabling selection
181 of a kernel in an elegant manner. */
183 /*! Pointers to the non-bonded kernels organized in 2-dim arrays by:
184 * electrostatics and VDW type.
186 * Note that the row- and column-order of function pointers has to match the
187 * order of corresponding enumerated electrostatics and vdw types, resp.,
188 * defined in nbnxn_cuda_types.h.
191 /*! Force-only kernel function pointers. */
192 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_noprune_ptr[eelCuNR][evdwCuNR] =
194 { nbnxn_kernel_ElecCut_VdwLJ_F_cuda, nbnxn_kernel_ElecCut_VdwLJCombGeom_F_cuda, nbnxn_kernel_ElecCut_VdwLJCombLB_F_cuda, nbnxn_kernel_ElecCut_VdwLJFsw_F_cuda, nbnxn_kernel_ElecCut_VdwLJPsw_F_cuda, nbnxn_kernel_ElecCut_VdwLJEwCombGeom_F_cuda, nbnxn_kernel_ElecCut_VdwLJEwCombLB_F_cuda },
195 { nbnxn_kernel_ElecRF_VdwLJ_F_cuda, nbnxn_kernel_ElecRF_VdwLJCombGeom_F_cuda, nbnxn_kernel_ElecRF_VdwLJCombLB_F_cuda, nbnxn_kernel_ElecRF_VdwLJFsw_F_cuda, nbnxn_kernel_ElecRF_VdwLJPsw_F_cuda, nbnxn_kernel_ElecRF_VdwLJEwCombGeom_F_cuda, nbnxn_kernel_ElecRF_VdwLJEwCombLB_F_cuda },
196 { nbnxn_kernel_ElecEwQSTab_VdwLJ_F_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_F_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_F_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJFsw_F_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJPsw_F_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_F_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_F_cuda },
197 { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_F_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_F_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_F_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_F_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_F_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_F_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_F_cuda },
198 { nbnxn_kernel_ElecEw_VdwLJ_F_cuda, nbnxn_kernel_ElecEw_VdwLJCombGeom_F_cuda, nbnxn_kernel_ElecEw_VdwLJCombLB_F_cuda, nbnxn_kernel_ElecEw_VdwLJFsw_F_cuda, nbnxn_kernel_ElecEw_VdwLJPsw_F_cuda, nbnxn_kernel_ElecEw_VdwLJEwCombGeom_F_cuda, nbnxn_kernel_ElecEw_VdwLJEwCombLB_F_cuda },
199 { nbnxn_kernel_ElecEwTwinCut_VdwLJ_F_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_F_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_F_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_F_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_F_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_F_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_F_cuda }
202 /*! Force + energy kernel function pointers. */
203 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_noprune_ptr[eelCuNR][evdwCuNR] =
205 { nbnxn_kernel_ElecCut_VdwLJ_VF_cuda, nbnxn_kernel_ElecCut_VdwLJCombGeom_VF_cuda, nbnxn_kernel_ElecCut_VdwLJCombLB_VF_cuda, nbnxn_kernel_ElecCut_VdwLJFsw_VF_cuda, nbnxn_kernel_ElecCut_VdwLJPsw_VF_cuda, nbnxn_kernel_ElecCut_VdwLJEwCombGeom_VF_cuda, nbnxn_kernel_ElecCut_VdwLJEwCombLB_VF_cuda },
206 { nbnxn_kernel_ElecRF_VdwLJ_VF_cuda, nbnxn_kernel_ElecRF_VdwLJCombGeom_VF_cuda, nbnxn_kernel_ElecRF_VdwLJCombLB_VF_cuda, nbnxn_kernel_ElecRF_VdwLJFsw_VF_cuda, nbnxn_kernel_ElecRF_VdwLJPsw_VF_cuda, nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VF_cuda, nbnxn_kernel_ElecRF_VdwLJEwCombLB_VF_cuda },
207 { nbnxn_kernel_ElecEwQSTab_VdwLJ_VF_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_VF_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_VF_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJFsw_VF_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJPsw_VF_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_VF_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_VF_cuda },
208 { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_VF_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_VF_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_VF_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_VF_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_VF_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_VF_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_VF_cuda },
209 { nbnxn_kernel_ElecEw_VdwLJ_VF_cuda, nbnxn_kernel_ElecEw_VdwLJCombGeom_VF_cuda, nbnxn_kernel_ElecEw_VdwLJCombLB_VF_cuda, nbnxn_kernel_ElecEw_VdwLJFsw_VF_cuda, nbnxn_kernel_ElecEw_VdwLJPsw_VF_cuda, nbnxn_kernel_ElecEw_VdwLJEwCombGeom_VF_cuda, nbnxn_kernel_ElecEw_VdwLJEwCombLB_VF_cuda },
210 { nbnxn_kernel_ElecEwTwinCut_VdwLJ_VF_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_VF_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_VF_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_VF_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_VF_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_VF_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_VF_cuda }
213 /*! Force + pruning kernel function pointers. */
214 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_prune_ptr[eelCuNR][evdwCuNR] =
216 { nbnxn_kernel_ElecCut_VdwLJ_F_prune_cuda, nbnxn_kernel_ElecCut_VdwLJCombGeom_F_prune_cuda, nbnxn_kernel_ElecCut_VdwLJCombLB_F_prune_cuda, nbnxn_kernel_ElecCut_VdwLJFsw_F_prune_cuda, nbnxn_kernel_ElecCut_VdwLJPsw_F_prune_cuda, nbnxn_kernel_ElecCut_VdwLJEwCombGeom_F_prune_cuda, nbnxn_kernel_ElecCut_VdwLJEwCombLB_F_prune_cuda },
217 { nbnxn_kernel_ElecRF_VdwLJ_F_prune_cuda, nbnxn_kernel_ElecRF_VdwLJCombGeom_F_prune_cuda, nbnxn_kernel_ElecRF_VdwLJCombLB_F_prune_cuda, nbnxn_kernel_ElecRF_VdwLJFsw_F_prune_cuda, nbnxn_kernel_ElecRF_VdwLJPsw_F_prune_cuda, nbnxn_kernel_ElecRF_VdwLJEwCombGeom_F_prune_cuda, nbnxn_kernel_ElecRF_VdwLJEwCombLB_F_prune_cuda },
218 { nbnxn_kernel_ElecEwQSTab_VdwLJ_F_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_F_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_F_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJFsw_F_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJPsw_F_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_F_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_F_prune_cuda },
219 { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_F_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_F_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_F_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_F_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_F_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_F_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_F_prune_cuda },
220 { nbnxn_kernel_ElecEw_VdwLJ_F_prune_cuda, nbnxn_kernel_ElecEw_VdwLJCombGeom_F_prune_cuda, nbnxn_kernel_ElecEw_VdwLJCombLB_F_prune_cuda, nbnxn_kernel_ElecEw_VdwLJFsw_F_prune_cuda, nbnxn_kernel_ElecEw_VdwLJPsw_F_prune_cuda, nbnxn_kernel_ElecEw_VdwLJEwCombGeom_F_prune_cuda, nbnxn_kernel_ElecEw_VdwLJEwCombLB_F_prune_cuda },
221 { nbnxn_kernel_ElecEwTwinCut_VdwLJ_F_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_F_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_F_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_F_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_F_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_F_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_F_prune_cuda }
224 /*! Force + energy + pruning kernel function pointers. */
225 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_prune_ptr[eelCuNR][evdwCuNR] =
227 { nbnxn_kernel_ElecCut_VdwLJ_VF_prune_cuda, nbnxn_kernel_ElecCut_VdwLJCombGeom_VF_prune_cuda, nbnxn_kernel_ElecCut_VdwLJCombLB_VF_prune_cuda, nbnxn_kernel_ElecCut_VdwLJFsw_VF_prune_cuda, nbnxn_kernel_ElecCut_VdwLJPsw_VF_prune_cuda, nbnxn_kernel_ElecCut_VdwLJEwCombGeom_VF_prune_cuda, nbnxn_kernel_ElecCut_VdwLJEwCombLB_VF_prune_cuda },
228 { nbnxn_kernel_ElecRF_VdwLJ_VF_prune_cuda, nbnxn_kernel_ElecRF_VdwLJCombGeom_VF_prune_cuda, nbnxn_kernel_ElecRF_VdwLJCombLB_VF_prune_cuda, nbnxn_kernel_ElecRF_VdwLJFsw_VF_prune_cuda, nbnxn_kernel_ElecRF_VdwLJPsw_VF_prune_cuda, nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VF_prune_cuda, nbnxn_kernel_ElecRF_VdwLJEwCombLB_VF_prune_cuda },
229 { nbnxn_kernel_ElecEwQSTab_VdwLJ_VF_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJCombGeom_VF_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJCombLB_VF_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJFsw_VF_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJPsw_VF_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJEwCombGeom_VF_prune_cuda, nbnxn_kernel_ElecEwQSTab_VdwLJEwCombLB_VF_prune_cuda },
230 { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_VF_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombGeom_VF_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJCombLB_VF_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_VF_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_VF_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombGeom_VF_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJEwCombLB_VF_prune_cuda },
231 { nbnxn_kernel_ElecEw_VdwLJ_VF_prune_cuda, nbnxn_kernel_ElecEw_VdwLJCombGeom_VF_prune_cuda, nbnxn_kernel_ElecEw_VdwLJCombLB_VF_prune_cuda, nbnxn_kernel_ElecEw_VdwLJFsw_VF_prune_cuda, nbnxn_kernel_ElecEw_VdwLJPsw_VF_prune_cuda, nbnxn_kernel_ElecEw_VdwLJEwCombGeom_VF_prune_cuda, nbnxn_kernel_ElecEw_VdwLJEwCombLB_VF_prune_cuda },
232 { nbnxn_kernel_ElecEwTwinCut_VdwLJ_VF_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_VF_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_VF_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_VF_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_VF_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_VF_prune_cuda, nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombLB_VF_prune_cuda }
235 /*! Return a pointer to the kernel version to be executed at the current step. */
236 static inline nbnxn_cu_kfunc_ptr_t select_nbnxn_kernel(int eeltype,
240 const gmx_device_info_t gmx_unused *devInfo)
242 nbnxn_cu_kfunc_ptr_t res;
244 GMX_ASSERT(eeltype < eelCuNR,
245 "The electrostatics type requested is not implemented in the CUDA kernels.");
246 GMX_ASSERT(evdwtype < evdwCuNR,
247 "The VdW type requested is not implemented in the CUDA kernels.");
249 /* assert assumptions made by the kernels */
250 GMX_ASSERT(c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit == devInfo->prop.warpSize,
251 "The CUDA kernels require the cluster_size_i*cluster_size_j/nbnxn_gpu_clusterpair_split to match the warp size of the architecture targeted.");
257 res = nb_kfunc_ener_prune_ptr[eeltype][evdwtype];
261 res = nb_kfunc_ener_noprune_ptr[eeltype][evdwtype];
268 res = nb_kfunc_noener_prune_ptr[eeltype][evdwtype];
272 res = nb_kfunc_noener_noprune_ptr[eeltype][evdwtype];
279 /*! \brief Calculates the amount of shared memory required by the nonbonded kernel in use. */
280 static inline int calc_shmem_required_nonbonded(const int num_threads_z, const gmx_device_info_t gmx_unused *dinfo, const cu_nbparam_t *nbp)
286 /* size of shmem (force-buffers/xq/atom type preloading) */
287 /* NOTE: with the default kernel on sm3.0 we need shmem only for pre-loading */
288 /* i-atom x+q in shared memory */
289 shmem = c_numClPerSupercl * c_clSize * sizeof(float4);
290 /* cj in shared memory, for each warp separately */
291 shmem += num_threads_z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(int);
292 if (dinfo->prop.major >= 3)
294 if (nbp->vdwtype == evdwCuCUTCOMBGEOM ||
295 nbp->vdwtype == evdwCuCUTCOMBLB)
297 /* i-atom LJ combination parameters in shared memory */
298 shmem += c_numClPerSupercl * c_clSize * sizeof(float2);
302 /* i-atom types in shared memory */
303 shmem += c_numClPerSupercl * c_clSize * sizeof(int);
306 if (dinfo->prop.major < 3)
308 /* force reduction buffers in shared memory */
309 shmem += c_clSize * c_clSize * 3 * sizeof(float);
314 /*! As we execute nonbonded workload in separate streams, before launching
315 the kernel we need to make sure that he following operations have completed:
316 - atomdata allocation and related H2D transfers (every nstlist step);
317 - pair list H2D transfer (every nstlist step);
318 - shift vector H2D transfer (every nstlist step);
319 - force (+shift force and energy) output clearing (every step).
321 These operations are issued in the local stream at the beginning of the step
322 and therefore always complete before the local kernel launch. The non-local
323 kernel is launched after the local on the same device/context hence it is
324 inherently scheduled after the operations in the local stream (including the
325 above "misc_ops") on pre-GK110 devices with single hardware queue, but on later
326 devices with multiple hardware queues the dependency needs to be enforced.
327 We use the misc_ops_and_local_H2D_done event to record the point where
328 the local x+q H2D (and all preceding) tasks are complete and synchronize
329 with this event in the non-local stream before launching the non-bonded kernel.
331 void nbnxn_gpu_launch_kernel(gmx_nbnxn_cuda_t *nb,
332 const nbnxn_atomdata_t *nbatom,
337 int adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
338 /* CUDA kernel launch-related stuff */
340 dim3 dim_block, dim_grid;
341 nbnxn_cu_kfunc_ptr_t nb_kernel = NULL; /* fn pointer to the nonbonded kernel */
343 cu_atomdata_t *adat = nb->atdat;
344 cu_nbparam_t *nbp = nb->nbparam;
345 cu_plist_t *plist = nb->plist[iloc];
346 cu_timers_t *t = nb->timers;
347 cudaStream_t stream = nb->stream[iloc];
349 bool bCalcEner = flags & GMX_FORCE_ENERGY;
350 bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
351 bool bDoTime = nb->bDoTime;
353 /* turn energy calculation always on/off (for debugging/testing only) */
354 bCalcEner = (bCalcEner || always_ener) && !never_ener;
356 /* Don't launch the non-local kernel if there is no work to do.
357 Doing the same for the local kernel is more complicated, since the
358 local part of the force array also depends on the non-local kernel.
359 So to avoid complicating the code and to reduce the risk of bugs,
360 we always call the local kernel, the local x+q copy and later (not in
361 this function) the stream wait, local f copyback and the f buffer
362 clearing. All these operations, except for the local interaction kernel,
363 are needed for the non-local interactions. The skip of the local kernel
364 call is taken care of later in this function. */
365 if (canSkipWork(nb, iloc))
367 plist->haveFreshList = false;
372 /* calculate the atom data index range based on locality */
376 adat_len = adat->natoms_local;
380 adat_begin = adat->natoms_local;
381 adat_len = adat->natoms - adat->natoms_local;
384 /* beginning of timed HtoD section */
387 stat = cudaEventRecord(t->start_nb_h2d[iloc], stream);
388 CU_RET_ERR(stat, "cudaEventRecord failed");
392 cu_copy_H2D_async(adat->xq + adat_begin, nbatom->x + adat_begin * 4,
393 adat_len * sizeof(*adat->xq), stream);
397 stat = cudaEventRecord(t->stop_nb_h2d[iloc], stream);
398 CU_RET_ERR(stat, "cudaEventRecord failed");
401 /* When we get here all misc operations issues in the local stream as well as
402 the local xq H2D are done,
403 so we record that in the local stream and wait for it in the nonlocal one. */
404 if (nb->bUseTwoStreams)
406 if (iloc == eintLocal)
408 stat = cudaEventRecord(nb->misc_ops_and_local_H2D_done, stream);
409 CU_RET_ERR(stat, "cudaEventRecord on misc_ops_and_local_H2D_done failed");
413 stat = cudaStreamWaitEvent(stream, nb->misc_ops_and_local_H2D_done, 0);
414 CU_RET_ERR(stat, "cudaStreamWaitEvent on misc_ops_and_local_H2D_done failed");
418 if (nbp->useDynamicPruning && plist->haveFreshList)
420 /* Prunes for rlistOuter and rlistInner, sets plist->haveFreshList=false
421 (TODO: ATM that's the way the timing accounting can distinguish between
422 separate prune kernel and combined force+prune, maybe we need a better way?).
424 nbnxn_gpu_launch_kernel_pruneonly(nb, iloc, 1);
427 if (plist->nsci == 0)
429 /* Don't launch an empty local kernel (not allowed with CUDA) */
433 /* beginning of timed nonbonded calculation section */
436 stat = cudaEventRecord(t->start_nb_k[iloc], stream);
437 CU_RET_ERR(stat, "cudaEventRecord failed");
440 /* get the pointer to the kernel flavor we need to use */
441 nb_kernel = select_nbnxn_kernel(nbp->eeltype,
444 (plist->haveFreshList && !nb->timers->didPrune[iloc]) || always_prune,
447 /* Kernel launch config:
448 * - The thread block dimensions match the size of i-clusters, j-clusters,
449 * and j-cluster concurrency, in x, y, and z, respectively.
450 * - The 1D block-grid contains as many blocks as super-clusters.
452 int num_threads_z = 1;
453 if (nb->dev_info->prop.major == 3 && nb->dev_info->prop.minor == 7)
457 nblock = calc_nb_kernel_nblock(plist->nsci, nb->dev_info);
458 dim_block = dim3(c_clSize, c_clSize, num_threads_z);
459 dim_grid = dim3(nblock, 1, 1);
460 shmem = calc_shmem_required_nonbonded(num_threads_z, nb->dev_info, nbp);
464 fprintf(debug, "Non-bonded GPU launch configuration:\n\tThread block: %ux%ux%u\n\t"
465 "\tGrid: %ux%u\n\t#Super-clusters/clusters: %d/%d (%d)\n"
467 dim_block.x, dim_block.y, dim_block.z,
468 dim_grid.x, dim_grid.y, plist->nsci*c_numClPerSupercl,
469 c_numClPerSupercl, plist->na_c,
473 if (bUseCudaLaunchKernel)
475 gmx_unused void* kernel_args[4];
476 kernel_args[0] = adat;
477 kernel_args[1] = nbp;
478 kernel_args[2] = plist;
479 kernel_args[3] = &bCalcFshift;
481 #if GMX_CUDA_VERSION >= 7000
482 cudaLaunchKernel((void *)nb_kernel, dim_grid, dim_block, kernel_args, shmem, stream);
487 nb_kernel<<< dim_grid, dim_block, shmem, stream>>> (*adat, *nbp, *plist, bCalcFshift);
489 CU_LAUNCH_ERR("k_calc_nb");
493 stat = cudaEventRecord(t->stop_nb_k[iloc], stream);
494 CU_RET_ERR(stat, "cudaEventRecord failed");
497 #if (defined(WIN32) || defined( _WIN32 ))
498 /* Windows: force flushing WDDM queue */
499 stat = cudaStreamQuery(stream);
503 /*! Calculates the amount of shared memory required by the CUDA kernel in use. */
504 static inline int calc_shmem_required_prune(const int num_threads_z)
508 /* i-atom x in shared memory */
509 shmem = c_numClPerSupercl * c_clSize * sizeof(float4);
510 /* cj in shared memory, for each warp separately */
511 shmem += num_threads_z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(int);
516 void nbnxn_gpu_launch_kernel_pruneonly(gmx_nbnxn_cuda_t *nb,
522 cu_atomdata_t *adat = nb->atdat;
523 cu_nbparam_t *nbp = nb->nbparam;
524 cu_plist_t *plist = nb->plist[iloc];
525 cu_timers_t *t = nb->timers;
526 cudaStream_t stream = nb->stream[iloc];
528 bool bDoTime = nb->bDoTime;
530 if (plist->haveFreshList)
532 GMX_ASSERT(numParts == 1, "With first pruning we expect 1 part");
534 /* Set rollingPruningNumParts to signal that it is not set */
535 plist->rollingPruningNumParts = 0;
536 plist->rollingPruningPart = 0;
540 if (plist->rollingPruningNumParts == 0)
542 plist->rollingPruningNumParts = numParts;
546 GMX_ASSERT(numParts == plist->rollingPruningNumParts, "It is not allowed to change numParts in between list generation steps");
550 /* Use a local variable for part and update in plist, so we can return here
551 * without duplicating the part increment code.
553 int part = plist->rollingPruningPart;
555 plist->rollingPruningPart++;
556 if (plist->rollingPruningPart >= plist->rollingPruningNumParts)
558 plist->rollingPruningPart = 0;
561 /* Compute the number of list entries to prune in this pass */
562 int numSciInPart = (plist->nsci - part)/numParts;
564 /* Don't launch the kernel if there is no work to do (not allowed with CUDA) */
565 if (numSciInPart <= 0)
567 plist->haveFreshList = false;
572 cudaEvent_t startEvent, stopEvent;
575 startEvent = (plist->haveFreshList ? t->start_prune_k[iloc] : t->start_rollingPrune_k[iloc]);
576 stopEvent = (plist->haveFreshList ? t->stop_prune_k[iloc] : t->stop_rollingPrune_k[iloc]);
579 /* beginning of timed prune calculation section */
582 stat = cudaEventRecord(startEvent, stream);
583 CU_RET_ERR(stat, "cudaEventRecord failed");
586 /* Kernel launch config:
587 * - The thread block dimensions match the size of i-clusters, j-clusters,
588 * and j-cluster concurrency, in x, y, and z, respectively.
589 * - The 1D block-grid contains as many blocks as super-clusters.
591 int num_threads_z = c_cudaPruneKernelJ4Concurrency;
592 int nblock = calc_nb_kernel_nblock(numSciInPart, nb->dev_info);
593 dim3 dim_block = dim3(c_clSize, c_clSize, num_threads_z);
594 dim3 dim_grid = dim3(nblock, 1, 1);
595 int shmem = calc_shmem_required_prune(num_threads_z);
599 fprintf(debug, "Pruning GPU kernel launch configuration:\n\tThread block: %dx%dx%d\n\t"
600 "\tGrid: %dx%d\n\t#Super-clusters/clusters: %d/%d (%d)\n"
602 dim_block.x, dim_block.y, dim_block.z,
603 dim_grid.x, dim_grid.y, numSciInPart*c_numClPerSupercl,
604 c_numClPerSupercl, plist->na_c,
608 if (bUseCudaLaunchKernel)
610 gmx_unused void* kernel_args[5];
611 kernel_args[0] = adat;
612 kernel_args[1] = nbp;
613 kernel_args[2] = plist;
614 kernel_args[3] = &numParts;
615 kernel_args[4] = ∂
617 #if GMX_CUDA_VERSION >= 7000
618 if (plist->haveFreshList)
620 cudaLaunchKernel((void *)nbnxn_kernel_prune_cuda<true>, dim_grid, dim_block, kernel_args, shmem, stream);
624 cudaLaunchKernel((void *)nbnxn_kernel_prune_cuda<false>, dim_grid, dim_block, kernel_args, shmem, stream);
630 if (plist->haveFreshList)
632 nbnxn_kernel_prune_cuda<true><<< dim_grid, dim_block, shmem, stream>>> (*adat, *nbp, *plist, numParts, part);
636 nbnxn_kernel_prune_cuda<false><<< dim_grid, dim_block, shmem, stream>>> (*adat, *nbp, *plist, numParts, part);
639 CU_LAUNCH_ERR("k_pruneonly");
641 /* TODO: consider a more elegant way to track which kernel has been called
642 (combined or separate 1st pass prune, rolling prune). */
643 if (plist->haveFreshList)
645 plist->haveFreshList = false;
646 /* Mark that pruning has been done */
647 nb->timers->didPrune[iloc] = true;
651 /* Mark that rolling pruning has been done */
652 nb->timers->didRollingPrune[iloc] = true;
657 stat = cudaEventRecord(stopEvent, stream);
658 CU_RET_ERR(stat, "cudaEventRecord failed");
661 #if (defined(WIN32) || defined( _WIN32 ))
662 /* Windows: force flushing WDDM queue */
663 stat = cudaStreamQuery(stream);
667 void nbnxn_gpu_launch_cpyback(gmx_nbnxn_cuda_t *nb,
668 const nbnxn_atomdata_t *nbatom,
673 int adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
676 /* determine interaction locality from atom locality */
681 else if (NONLOCAL_A(aloc))
688 sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
689 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
693 cu_atomdata_t *adat = nb->atdat;
694 cu_timers_t *t = nb->timers;
695 bool bDoTime = nb->bDoTime;
696 cudaStream_t stream = nb->stream[iloc];
698 bool bCalcEner = flags & GMX_FORCE_ENERGY;
699 bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
701 /* don't launch non-local copy-back if there was no non-local work to do */
702 if (canSkipWork(nb, iloc))
707 /* calculate the atom data index range based on locality */
711 adat_len = adat->natoms_local;
715 adat_begin = adat->natoms_local;
716 adat_len = adat->natoms - adat->natoms_local;
719 /* beginning of timed D2H section */
722 stat = cudaEventRecord(t->start_nb_d2h[iloc], stream);
723 CU_RET_ERR(stat, "cudaEventRecord failed");
726 /* With DD the local D2H transfer can only start after the non-local
727 kernel has finished. */
728 if (iloc == eintLocal && nb->bUseTwoStreams)
730 stat = cudaStreamWaitEvent(stream, nb->nonlocal_done, 0);
731 CU_RET_ERR(stat, "cudaStreamWaitEvent on nonlocal_done failed");
735 cu_copy_D2H_async(nbatom->out[0].f + adat_begin * 3, adat->f + adat_begin,
736 (adat_len)*sizeof(*adat->f), stream);
738 /* After the non-local D2H is launched the nonlocal_done event can be
739 recorded which signals that the local D2H can proceed. This event is not
740 placed after the non-local kernel because we want the non-local data
742 if (iloc == eintNonlocal)
744 stat = cudaEventRecord(nb->nonlocal_done, stream);
745 CU_RET_ERR(stat, "cudaEventRecord on nonlocal_done failed");
748 /* only transfer energies in the local stream */
754 cu_copy_D2H_async(nb->nbst.fshift, adat->fshift,
755 SHIFTS * sizeof(*nb->nbst.fshift), stream);
761 cu_copy_D2H_async(nb->nbst.e_lj, adat->e_lj,
762 sizeof(*nb->nbst.e_lj), stream);
763 cu_copy_D2H_async(nb->nbst.e_el, adat->e_el,
764 sizeof(*nb->nbst.e_el), stream);
770 stat = cudaEventRecord(t->stop_nb_d2h[iloc], stream);
771 CU_RET_ERR(stat, "cudaEventRecord failed");
775 /*! \brief Count pruning kernel time if either kernel has been triggered
777 * We do the accounting for either of the two pruning kernel flavors:
778 * - 1st pass prune: ran during the current step (prior to the force kernel);
779 * - rolling prune: ran at the end of the previous step (prior to the current step H2D xq);
781 * Note that the resetting of cu_timers_t::didPrune and cu_timers_t::didRollingPrune should happen
782 * after calling this function.
784 * \param[in] timers structs with CUDA timer objects
785 * \param[inout] timings GPU task timing data
786 * \param[in] iloc interaction locality
788 static void countPruneKernelTime(const cu_timers_t *timers,
789 gmx_wallclock_gpu_t *timings,
792 // We might have not done any pruning (e.g. if we skipped with empty domains).
793 if (!timers->didPrune[iloc] && !timers->didRollingPrune[iloc])
798 if (timers->didPrune[iloc])
800 timings->pruneTime.c++;
801 timings->pruneTime.t += cu_event_elapsed(timers->start_prune_k[iloc],
802 timers->stop_prune_k[iloc]);
804 if (timers->didRollingPrune[iloc])
806 timings->dynamicPruneTime.c++;
807 timings->dynamicPruneTime.t += cu_event_elapsed(timers->start_rollingPrune_k[iloc],
808 timers->stop_rollingPrune_k[iloc]);
812 void nbnxn_gpu_wait_for_gpu(gmx_nbnxn_cuda_t *nb,
814 real *e_lj, real *e_el, rvec *fshift)
816 /* NOTE: only implemented for single-precision at this time */
820 /* determine interaction locality from atom locality */
825 else if (NONLOCAL_A(aloc))
832 sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
833 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
837 cu_plist_t *plist = nb->plist[iloc];
838 cu_timers_t *timers = nb->timers;
839 struct gmx_wallclock_gpu_t *timings = nb->timings;
840 nb_staging nbst = nb->nbst;
842 bool bCalcEner = flags & GMX_FORCE_ENERGY;
843 bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
845 /* turn energy calculation always on/off (for debugging/testing only) */
846 bCalcEner = (bCalcEner || always_ener) && !never_ener;
848 /* Launch wait/update timers & counters and do reduction into staging buffers
849 BUT skip it when during the non-local phase there was actually no work to do.
850 This is consistent with nbnxn_gpu_launch_kernel.
852 NOTE: if timing with multiple GPUs (streams) becomes possible, the
853 counters could end up being inconsistent due to not being incremented
854 on some of the nodes! */
855 if (!canSkipWork(nb, iloc))
857 stat = cudaStreamSynchronize(nb->stream[iloc]);
858 CU_RET_ERR(stat, "cudaStreamSynchronize failed in cu_blockwait_nb");
860 /* timing data accumulation */
863 /* only increase counter once (at local F wait) */
867 timings->ktime[plist->haveFreshList ? 1 : 0][bCalcEner ? 1 : 0].c += 1;
871 timings->ktime[plist->haveFreshList ? 1 : 0][bCalcEner ? 1 : 0].t +=
872 cu_event_elapsed(timers->start_nb_k[iloc], timers->stop_nb_k[iloc]);
874 /* X/q H2D and F D2H timings */
875 timings->nb_h2d_t += cu_event_elapsed(timers->start_nb_h2d[iloc],
876 timers->stop_nb_h2d[iloc]);
877 timings->nb_d2h_t += cu_event_elapsed(timers->start_nb_d2h[iloc],
878 timers->stop_nb_d2h[iloc]);
880 /* Count the pruning kernel times for both cases:1st pass (at search step)
881 and rolling pruning (if called at the previous step).
882 We do the accounting here as this is the only sync point where we
883 know (without checking or additional sync-ing) that prune tasks in
884 in the current stream have completed (having just blocking-waited
885 for the force D2H). */
886 countPruneKernelTime(timers, timings, iloc);
888 /* only count atdat and pair-list H2D at pair-search step */
889 if (timers->didPairlistH2D[iloc])
891 /* atdat transfer timing (add only once, at local F wait) */
895 timings->pl_h2d_t += cu_event_elapsed(timers->start_atdat,
899 timings->pl_h2d_t += cu_event_elapsed(timers->start_pl_h2d[iloc],
900 timers->stop_pl_h2d[iloc]);
902 /* Clear the timing flag for the next step */
903 timers->didPairlistH2D[iloc] = false;
907 /* add up energies and shift forces (only once at local F wait) */
918 for (int i = 0; i < SHIFTS; i++)
920 fshift[i][0] += nbst.fshift[i].x;
921 fshift[i][1] += nbst.fshift[i].y;
922 fshift[i][2] += nbst.fshift[i].z;
928 /* Always reset both pruning flags (doesn't hurt doing it even when timing is off). */
929 timers->didPrune[iloc] = timers->didRollingPrune[iloc] = false;
931 /* Turn off initial list pruning (doesn't hurt if this is not pair-search step). */
932 plist->haveFreshList = false;
935 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_texref()
940 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_comb_texref()
942 return nbfp_comb_texref;
945 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_coulomb_tab_texref()
947 return coulomb_tab_texref;
950 void nbnxn_cuda_set_cacheconfig(const gmx_device_info_t *devinfo)
954 for (int i = 0; i < eelCuNR; i++)
956 for (int j = 0; j < evdwCuNR; j++)
958 if (devinfo->prop.major >= 3)
960 /* Default kernel on sm 3.x and later 32/32 kB Shared/L1 */
961 cudaFuncSetCacheConfig(nb_kfunc_ener_prune_ptr[i][j], cudaFuncCachePreferEqual);
962 cudaFuncSetCacheConfig(nb_kfunc_ener_noprune_ptr[i][j], cudaFuncCachePreferEqual);
963 cudaFuncSetCacheConfig(nb_kfunc_noener_prune_ptr[i][j], cudaFuncCachePreferEqual);
964 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_noprune_ptr[i][j], cudaFuncCachePreferEqual);
968 /* On Fermi prefer L1 gives 2% higher performance */
969 /* Default kernel on sm_2.x 16/48 kB Shared/L1 */
970 cudaFuncSetCacheConfig(nb_kfunc_ener_prune_ptr[i][j], cudaFuncCachePreferL1);
971 cudaFuncSetCacheConfig(nb_kfunc_ener_noprune_ptr[i][j], cudaFuncCachePreferL1);
972 cudaFuncSetCacheConfig(nb_kfunc_noener_prune_ptr[i][j], cudaFuncCachePreferL1);
973 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_noprune_ptr[i][j], cudaFuncCachePreferL1);
975 CU_RET_ERR(stat, "cudaFuncSetCacheConfig failed");