2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016, 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 "gromacs/gpu_utils/cudautils.cuh"
55 #include "gromacs/mdlib/force_flags.h"
56 #include "gromacs/mdlib/nb_verlet.h"
57 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
58 #include "gromacs/mdlib/nbnxn_pairlist.h"
59 #include "gromacs/timing/gpu_timing.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/gmxassert.h"
63 #include "nbnxn_cuda_types.h"
65 /*! Texture reference for LJ C6/C12 parameters; bound to cu_nbparam_t.nbfp */
66 texture<float, 1, cudaReadModeElementType> nbfp_texref;
68 /*! Texture reference for LJ-PME parameters; bound to cu_nbparam_t.nbfp_comb */
69 texture<float, 1, cudaReadModeElementType> nbfp_comb_texref;
71 /*! Texture reference for Ewald coulomb force table; bound to cu_nbparam_t.coulomb_tab */
72 texture<float, 1, cudaReadModeElementType> coulomb_tab_texref;
75 /***** The kernel declarations/definitions come here *****/
77 /* Top-level kernel declaration generation: will generate through multiple
78 * inclusion the following flavors for all kernel declarations:
79 * - force-only output;
80 * - force and energy output;
81 * - force-only with pair list pruning;
82 * - force and energy output with pair list pruning.
84 #define FUNCTION_DECLARATION_ONLY
86 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
87 /** Force & energy **/
89 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
92 /*** Pair-list pruning kernels ***/
95 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
96 /** Force & energy **/
98 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
101 #undef FUNCTION_DECLARATION_ONLY
103 /* Now generate the function definitions if we are using a single compilation unit. */
104 #if GMX_CUDA_NB_SINGLE_COMPILATION_UNIT
105 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_F_noprune.cu"
106 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_F_prune.cu"
107 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_VF_noprune.cu"
108 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_VF_prune.cu"
110 /* Prevent compilation in multiple compilation unit mode for CC 2.x. Although we have
111 * build-time checks to prevent this, the user could manually tweaks nvcc flags
112 * which would lead to buggy kernels getting compiled.
114 #if GMX_PTX_ARCH > 0 && GMX_PTX_ARCH <= 210
115 #error Due to an CUDA 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.
117 #endif /* GMX_CUDA_NB_SINGLE_COMPILATION_UNIT */
121 /*! Nonbonded kernel function pointer type */
122 typedef void (*nbnxn_cu_kfunc_ptr_t)(const cu_atomdata_t,
127 /*********************************/
129 /* XXX switch between chevron and cudaLaunch (supported only in CUDA >=7.0)
130 -- only for benchmarking purposes */
131 static const bool bUseCudaLaunchKernel =
132 (GMX_CUDA_VERSION >= 7000) && (getenv("GMX_DISABLE_CUDALAUNCH") == NULL);
134 /* XXX always/never run the energy/pruning kernels -- only for benchmarking purposes */
135 static bool always_ener = (getenv("GMX_GPU_ALWAYS_ENER") != NULL);
136 static bool never_ener = (getenv("GMX_GPU_NEVER_ENER") != NULL);
137 static bool always_prune = (getenv("GMX_GPU_ALWAYS_PRUNE") != NULL);
140 /*! Returns the number of blocks to be used for the nonbonded GPU kernel. */
141 static inline int calc_nb_kernel_nblock(int nwork_units, gmx_device_info_t *dinfo)
146 /* CUDA does not accept grid dimension of 0 (which can happen e.g. with an
147 empty domain) and that case should be handled before this point. */
148 assert(nwork_units > 0);
150 max_grid_x_size = dinfo->prop.maxGridSize[0];
152 /* do we exceed the grid x dimension limit? */
153 if (nwork_units > max_grid_x_size)
155 gmx_fatal(FARGS, "Watch out, the input system is too large to simulate!\n"
156 "The number of nonbonded work units (=number of super-clusters) exceeds the"
157 "maximum grid size in x dimension (%d > %d)!", nwork_units, max_grid_x_size);
164 /* Constant arrays listing all kernel function pointers and enabling selection
165 of a kernel in an elegant manner. */
167 /*! Pointers to the non-bonded kernels organized in 2-dim arrays by:
168 * electrostatics and VDW type.
170 * Note that the row- and column-order of function pointers has to match the
171 * order of corresponding enumerated electrostatics and vdw types, resp.,
172 * defined in nbnxn_cuda_types.h.
175 /*! Force-only kernel function pointers. */
176 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_noprune_ptr[eelCuNR][evdwCuNR] =
178 { 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 },
179 { 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 },
180 { 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 },
181 { 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 },
182 { 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 },
183 { 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 }
186 /*! Force + energy kernel function pointers. */
187 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_noprune_ptr[eelCuNR][evdwCuNR] =
189 { 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 },
190 { 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 },
191 { 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 },
192 { 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 },
193 { 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 },
194 { 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 }
197 /*! Force + pruning kernel function pointers. */
198 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_prune_ptr[eelCuNR][evdwCuNR] =
200 { 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 },
201 { 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 },
202 { 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 },
203 { 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 },
204 { 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 },
205 { 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 }
208 /*! Force + energy + pruning kernel function pointers. */
209 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_prune_ptr[eelCuNR][evdwCuNR] =
211 { 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 },
212 { 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 },
213 { 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 },
214 { 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 },
215 { 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 },
216 { 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 }
219 /*! Return a pointer to the kernel version to be executed at the current step. */
220 static inline nbnxn_cu_kfunc_ptr_t select_nbnxn_kernel(int eeltype,
224 struct gmx_device_info_t gmx_unused *devInfo)
226 nbnxn_cu_kfunc_ptr_t res;
228 GMX_ASSERT(eeltype < eelCuNR,
229 "The electrostatics type requested is not implemented in the CUDA kernels.");
230 GMX_ASSERT(evdwtype < evdwCuNR,
231 "The VdW type requested is not implemented in the CUDA kernels.");
233 /* assert assumptions made by the kernels */
234 GMX_ASSERT(c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit == devInfo->prop.warpSize,
235 "The CUDA kernels require the cluster_size_i*cluster_size_j/nbnxn_gpu_clusterpair_split to match the warp size of the architecture targeted.");
241 res = nb_kfunc_ener_prune_ptr[eeltype][evdwtype];
245 res = nb_kfunc_ener_noprune_ptr[eeltype][evdwtype];
252 res = nb_kfunc_noener_prune_ptr[eeltype][evdwtype];
256 res = nb_kfunc_noener_noprune_ptr[eeltype][evdwtype];
263 /*! Calculates the amount of shared memory required by the CUDA kernel in use. */
264 static inline int calc_shmem_required(const int num_threads_z, gmx_device_info_t gmx_unused *dinfo, const cu_nbparam_t *nbp)
270 /* size of shmem (force-buffers/xq/atom type preloading) */
271 /* NOTE: with the default kernel on sm3.0 we need shmem only for pre-loading */
272 /* i-atom x+q in shared memory */
273 shmem = c_numClPerSupercl * c_clSize * sizeof(float4);
274 /* cj in shared memory, for each warp separately */
275 shmem += num_threads_z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(int);
276 if (dinfo->prop.major >= 3)
278 if (nbp->vdwtype == evdwCuCUTCOMBGEOM ||
279 nbp->vdwtype == evdwCuCUTCOMBLB)
281 /* i-atom LJ combination parameters in shared memory */
282 shmem += c_numClPerSupercl * c_clSize * sizeof(float2);
286 /* i-atom types in shared memory */
287 shmem += c_numClPerSupercl * c_clSize * sizeof(int);
290 if (dinfo->prop.major < 3)
292 /* force reduction buffers in shared memory */
293 shmem += c_clSize * c_clSize * 3 * sizeof(float);
298 /*! As we execute nonbonded workload in separate streams, before launching
299 the kernel we need to make sure that he following operations have completed:
300 - atomdata allocation and related H2D transfers (every nstlist step);
301 - pair list H2D transfer (every nstlist step);
302 - shift vector H2D transfer (every nstlist step);
303 - force (+shift force and energy) output clearing (every step).
305 These operations are issued in the local stream at the beginning of the step
306 and therefore always complete before the local kernel launch. The non-local
307 kernel is launched after the local on the same device/context hence it is
308 inherently scheduled after the operations in the local stream (including the
309 above "misc_ops") on pre-GK110 devices with single hardware queue, but on later
310 devices with multiple hardware queues the dependency needs to be enforced.
311 We use the misc_ops_and_local_H2D_done event to record the point where
312 the local x+q H2D (and all preceding) tasks are complete and synchronize
313 with this event in the non-local stream before launching the non-bonded kernel.
315 void nbnxn_gpu_launch_kernel(gmx_nbnxn_cuda_t *nb,
316 const nbnxn_atomdata_t *nbatom,
321 int adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
322 /* CUDA kernel launch-related stuff */
324 dim3 dim_block, dim_grid;
325 nbnxn_cu_kfunc_ptr_t nb_kernel = NULL; /* fn pointer to the nonbonded kernel */
327 cu_atomdata_t *adat = nb->atdat;
328 cu_nbparam_t *nbp = nb->nbparam;
329 cu_plist_t *plist = nb->plist[iloc];
330 cu_timers_t *t = nb->timers;
331 cudaStream_t stream = nb->stream[iloc];
333 bool bCalcEner = flags & GMX_FORCE_ENERGY;
334 bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
335 bool bDoTime = nb->bDoTime;
337 /* turn energy calculation always on/off (for debugging/testing only) */
338 bCalcEner = (bCalcEner || always_ener) && !never_ener;
340 /* Don't launch the non-local kernel if there is no work to do.
341 Doing the same for the local kernel is more complicated, since the
342 local part of the force array also depends on the non-local kernel.
343 So to avoid complicating the code and to reduce the risk of bugs,
344 we always call the local kernel, the local x+q copy and later (not in
345 this function) the stream wait, local f copyback and the f buffer
346 clearing. All these operations, except for the local interaction kernel,
347 are needed for the non-local interactions. The skip of the local kernel
348 call is taken care of later in this function. */
349 if (iloc == eintNonlocal && plist->nsci == 0)
354 /* calculate the atom data index range based on locality */
358 adat_len = adat->natoms_local;
362 adat_begin = adat->natoms_local;
363 adat_len = adat->natoms - adat->natoms_local;
366 /* beginning of timed HtoD section */
369 stat = cudaEventRecord(t->start_nb_h2d[iloc], stream);
370 CU_RET_ERR(stat, "cudaEventRecord failed");
374 cu_copy_H2D_async(adat->xq + adat_begin, nbatom->x + adat_begin * 4,
375 adat_len * sizeof(*adat->xq), stream);
377 /* When we get here all misc operations issues in the local stream as well as
378 the local xq H2D are done,
379 so we record that in the local stream and wait for it in the nonlocal one. */
380 if (nb->bUseTwoStreams)
382 if (iloc == eintLocal)
384 stat = cudaEventRecord(nb->misc_ops_and_local_H2D_done, stream);
385 CU_RET_ERR(stat, "cudaEventRecord on misc_ops_and_local_H2D_done failed");
389 stat = cudaStreamWaitEvent(stream, nb->misc_ops_and_local_H2D_done, 0);
390 CU_RET_ERR(stat, "cudaStreamWaitEvent on misc_ops_and_local_H2D_done failed");
396 stat = cudaEventRecord(t->stop_nb_h2d[iloc], stream);
397 CU_RET_ERR(stat, "cudaEventRecord failed");
400 if (plist->nsci == 0)
402 /* Don't launch an empty local kernel (not allowed with CUDA) */
406 /* beginning of timed nonbonded calculation section */
409 stat = cudaEventRecord(t->start_nb_k[iloc], stream);
410 CU_RET_ERR(stat, "cudaEventRecord failed");
413 /* get the pointer to the kernel flavor we need to use */
414 nb_kernel = select_nbnxn_kernel(nbp->eeltype,
417 plist->bDoPrune || always_prune,
420 /* Kernel launch config:
421 * - The thread block dimensions match the size of i-clusters, j-clusters,
422 * and j-cluster concurrency, in x, y, and z, respectively.
423 * - The 1D block-grid contains as many blocks as super-clusters.
425 int num_threads_z = 1;
426 if (nb->dev_info->prop.major == 3 && nb->dev_info->prop.minor == 7)
430 nblock = calc_nb_kernel_nblock(plist->nsci, nb->dev_info);
431 dim_block = dim3(c_clSize, c_clSize, num_threads_z);
432 dim_grid = dim3(nblock, 1, 1);
433 shmem = calc_shmem_required(num_threads_z, nb->dev_info, nbp);
437 fprintf(debug, "GPU launch configuration:\n\tThread block: %dx%dx%d\n\t"
438 "\tGrid: %dx%d\n\t#Super-clusters/clusters: %d/%d (%d)\n"
440 dim_block.x, dim_block.y, dim_block.z,
441 dim_grid.x, dim_grid.y, plist->nsci*c_numClPerSupercl,
442 c_numClPerSupercl, plist->na_c,
446 if (bUseCudaLaunchKernel)
448 gmx_unused void* kernel_args[4];
449 kernel_args[0] = adat;
450 kernel_args[1] = nbp;
451 kernel_args[2] = plist;
452 kernel_args[3] = &bCalcFshift;
454 #if GMX_CUDA_VERSION >= 7000
455 cudaLaunchKernel((void *)nb_kernel, dim_grid, dim_block, kernel_args, shmem, stream);
460 nb_kernel<<< dim_grid, dim_block, shmem, stream>>> (*adat, *nbp, *plist, bCalcFshift);
462 CU_LAUNCH_ERR("k_calc_nb");
466 stat = cudaEventRecord(t->stop_nb_k[iloc], stream);
467 CU_RET_ERR(stat, "cudaEventRecord failed");
470 #if (defined(WIN32) || defined( _WIN32 ))
471 /* Windows: force flushing WDDM queue */
472 stat = cudaStreamQuery(stream);
476 void nbnxn_gpu_launch_cpyback(gmx_nbnxn_cuda_t *nb,
477 const nbnxn_atomdata_t *nbatom,
482 int adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
485 /* determine interaction locality from atom locality */
490 else if (NONLOCAL_A(aloc))
497 sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
498 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
502 cu_atomdata_t *adat = nb->atdat;
503 cu_timers_t *t = nb->timers;
504 bool bDoTime = nb->bDoTime;
505 cudaStream_t stream = nb->stream[iloc];
507 bool bCalcEner = flags & GMX_FORCE_ENERGY;
508 bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
510 /* don't launch non-local copy-back if there was no non-local work to do */
511 if (iloc == eintNonlocal && nb->plist[iloc]->nsci == 0)
516 /* calculate the atom data index range based on locality */
520 adat_len = adat->natoms_local;
524 adat_begin = adat->natoms_local;
525 adat_len = adat->natoms - adat->natoms_local;
528 /* beginning of timed D2H section */
531 stat = cudaEventRecord(t->start_nb_d2h[iloc], stream);
532 CU_RET_ERR(stat, "cudaEventRecord failed");
535 /* With DD the local D2H transfer can only start after the non-local
536 kernel has finished. */
537 if (iloc == eintLocal && nb->bUseTwoStreams)
539 stat = cudaStreamWaitEvent(stream, nb->nonlocal_done, 0);
540 CU_RET_ERR(stat, "cudaStreamWaitEvent on nonlocal_done failed");
544 cu_copy_D2H_async(nbatom->out[0].f + adat_begin * 3, adat->f + adat_begin,
545 (adat_len)*sizeof(*adat->f), stream);
547 /* After the non-local D2H is launched the nonlocal_done event can be
548 recorded which signals that the local D2H can proceed. This event is not
549 placed after the non-local kernel because we want the non-local data
551 if (iloc == eintNonlocal)
553 stat = cudaEventRecord(nb->nonlocal_done, stream);
554 CU_RET_ERR(stat, "cudaEventRecord on nonlocal_done failed");
557 /* only transfer energies in the local stream */
563 cu_copy_D2H_async(nb->nbst.fshift, adat->fshift,
564 SHIFTS * sizeof(*nb->nbst.fshift), stream);
570 cu_copy_D2H_async(nb->nbst.e_lj, adat->e_lj,
571 sizeof(*nb->nbst.e_lj), stream);
572 cu_copy_D2H_async(nb->nbst.e_el, adat->e_el,
573 sizeof(*nb->nbst.e_el), stream);
579 stat = cudaEventRecord(t->stop_nb_d2h[iloc], stream);
580 CU_RET_ERR(stat, "cudaEventRecord failed");
584 void nbnxn_gpu_wait_for_gpu(gmx_nbnxn_cuda_t *nb,
586 real *e_lj, real *e_el, rvec *fshift)
588 /* NOTE: only implemented for single-precision at this time */
592 /* determine interaction locality from atom locality */
597 else if (NONLOCAL_A(aloc))
604 sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
605 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
609 cu_plist_t *plist = nb->plist[iloc];
610 cu_timers_t *timers = nb->timers;
611 struct gmx_wallclock_gpu_t *timings = nb->timings;
612 nb_staging nbst = nb->nbst;
614 bool bCalcEner = flags & GMX_FORCE_ENERGY;
615 bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
617 /* turn energy calculation always on/off (for debugging/testing only) */
618 bCalcEner = (bCalcEner || always_ener) && !never_ener;
620 /* Launch wait/update timers & counters, unless doing the non-local phase
621 when there is not actually work to do. This is consistent with
622 nbnxn_cuda_launch_kernel.
624 NOTE: if timing with multiple GPUs (streams) becomes possible, the
625 counters could end up being inconsistent due to not being incremented
626 on some of the nodes! */
627 if (iloc == eintNonlocal && nb->plist[iloc]->nsci == 0)
632 stat = cudaStreamSynchronize(nb->stream[iloc]);
633 CU_RET_ERR(stat, "cudaStreamSynchronize failed in cu_blockwait_nb");
635 /* timing data accumulation */
638 /* only increase counter once (at local F wait) */
642 timings->ktime[plist->bDoPrune ? 1 : 0][bCalcEner ? 1 : 0].c += 1;
646 timings->ktime[plist->bDoPrune ? 1 : 0][bCalcEner ? 1 : 0].t +=
647 cu_event_elapsed(timers->start_nb_k[iloc], timers->stop_nb_k[iloc]);
649 /* X/q H2D and F D2H timings */
650 timings->nb_h2d_t += cu_event_elapsed(timers->start_nb_h2d[iloc],
651 timers->stop_nb_h2d[iloc]);
652 timings->nb_d2h_t += cu_event_elapsed(timers->start_nb_d2h[iloc],
653 timers->stop_nb_d2h[iloc]);
655 /* only count atdat and pair-list H2D at pair-search step */
658 /* atdat transfer timing (add only once, at local F wait) */
662 timings->pl_h2d_t += cu_event_elapsed(timers->start_atdat,
666 timings->pl_h2d_t += cu_event_elapsed(timers->start_pl_h2d[iloc],
667 timers->stop_pl_h2d[iloc]);
671 /* add up energies and shift forces (only once at local F wait) */
682 for (int i = 0; i < SHIFTS; i++)
684 fshift[i][0] += nbst.fshift[i].x;
685 fshift[i][1] += nbst.fshift[i].y;
686 fshift[i][2] += nbst.fshift[i].z;
691 /* turn off pruning (doesn't matter if this is pair-search step or not) */
692 plist->bDoPrune = false;
695 /*! Return the reference to the nbfp texture. */
696 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_texref()
701 /*! Return the reference to the nbfp_comb texture. */
702 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_comb_texref()
704 return nbfp_comb_texref;
707 /*! Return the reference to the coulomb_tab. */
708 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_coulomb_tab_texref()
710 return coulomb_tab_texref;
713 /*! Set up the cache configuration for the non-bonded kernels,
715 void nbnxn_cuda_set_cacheconfig(gmx_device_info_t *devinfo)
719 for (int i = 0; i < eelCuNR; i++)
721 for (int j = 0; j < evdwCuNR; j++)
723 if (devinfo->prop.major >= 3)
725 /* Default kernel on sm 3.x and later 32/32 kB Shared/L1 */
726 cudaFuncSetCacheConfig(nb_kfunc_ener_prune_ptr[i][j], cudaFuncCachePreferEqual);
727 cudaFuncSetCacheConfig(nb_kfunc_ener_noprune_ptr[i][j], cudaFuncCachePreferEqual);
728 cudaFuncSetCacheConfig(nb_kfunc_noener_prune_ptr[i][j], cudaFuncCachePreferEqual);
729 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_noprune_ptr[i][j], cudaFuncCachePreferEqual);
733 /* On Fermi prefer L1 gives 2% higher performance */
734 /* Default kernel on sm_2.x 16/48 kB Shared/L1 */
735 cudaFuncSetCacheConfig(nb_kfunc_ener_prune_ptr[i][j], cudaFuncCachePreferL1);
736 cudaFuncSetCacheConfig(nb_kfunc_ener_noprune_ptr[i][j], cudaFuncCachePreferL1);
737 cudaFuncSetCacheConfig(nb_kfunc_noener_prune_ptr[i][j], cudaFuncCachePreferL1);
738 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_noprune_ptr[i][j], cudaFuncCachePreferL1);
740 CU_RET_ERR(stat, "cudaFuncSetCacheConfig failed");