Removed support for NVIDIA CC 2.x devices (codename Fermi)
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_cuda / nbnxn_cuda.cu
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 /*! \file
36  *  \brief Define CUDA implementation of nbnxn_gpu.h
37  *
38  *  \author Szilard Pall <pall.szilard@gmail.com>
39  */
40 #include "gmxpre.h"
41
42 #include "config.h"
43
44 #include <assert.h>
45 #include <stdlib.h>
46
47 #include "gromacs/mdlib/nbnxn_gpu.h"
48
49 #if defined(_MSVC)
50 #include <limits>
51 #endif
52
53
54 #include "nbnxn_cuda.h"
55
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_common_utils.h"
61 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
62 #include "gromacs/mdlib/nbnxn_pairlist.h"
63 #include "gromacs/timing/gpu_timing.h"
64 #include "gromacs/utility/cstringutil.h"
65 #include "gromacs/utility/gmxassert.h"
66
67 #include "nbnxn_cuda_types.h"
68
69
70 /***** The kernel declarations/definitions come here *****/
71
72 /* Top-level kernel declaration generation: will generate through multiple
73  * inclusion the following flavors for all kernel declarations:
74  * - force-only output;
75  * - force and energy output;
76  * - force-only with pair list pruning;
77  * - force and energy output with pair list pruning.
78  */
79 #define FUNCTION_DECLARATION_ONLY
80 /** Force only **/
81 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
82 /** Force & energy **/
83 #define CALC_ENERGIES
84 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
85 #undef CALC_ENERGIES
86
87 /*** Pair-list pruning kernels ***/
88 /** Force only **/
89 #define PRUNE_NBL
90 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
91 /** Force & energy **/
92 #define CALC_ENERGIES
93 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
94 #undef CALC_ENERGIES
95 #undef PRUNE_NBL
96
97 /* Prune-only kernels */
98 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_pruneonly.cuh"
99 #undef FUNCTION_DECLARATION_ONLY
100
101 /* Now generate the function definitions if we are using a single compilation unit. */
102 #if GMX_CUDA_NB_SINGLE_COMPILATION_UNIT
103 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_F_noprune.cu"
104 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_F_prune.cu"
105 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_VF_noprune.cu"
106 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_VF_prune.cu"
107 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_pruneonly.cu"
108 #endif /* GMX_CUDA_NB_SINGLE_COMPILATION_UNIT */
109
110
111 /*! Nonbonded kernel function pointer type */
112 typedef void (*nbnxn_cu_kfunc_ptr_t)(const cu_atomdata_t,
113                                      const cu_nbparam_t,
114                                      const cu_plist_t,
115                                      bool);
116
117 /*********************************/
118
119 /*! Returns the number of blocks to be used for the nonbonded GPU kernel. */
120 static inline int calc_nb_kernel_nblock(int nwork_units, const gmx_device_info_t *dinfo)
121 {
122     int max_grid_x_size;
123
124     assert(dinfo);
125     /* CUDA does not accept grid dimension of 0 (which can happen e.g. with an
126        empty domain) and that case should be handled before this point. */
127     assert(nwork_units > 0);
128
129     max_grid_x_size = dinfo->prop.maxGridSize[0];
130
131     /* do we exceed the grid x dimension limit? */
132     if (nwork_units > max_grid_x_size)
133     {
134         gmx_fatal(FARGS, "Watch out, the input system is too large to simulate!\n"
135                   "The number of nonbonded work units (=number of super-clusters) exceeds the"
136                   "maximum grid size in x dimension (%d > %d)!", nwork_units, max_grid_x_size);
137     }
138
139     return nwork_units;
140 }
141
142
143 /* Constant arrays listing all kernel function pointers and enabling selection
144    of a kernel in an elegant manner. */
145
146 /*! Pointers to the non-bonded kernels organized in 2-dim arrays by:
147  *  electrostatics and VDW type.
148  *
149  *  Note that the row- and column-order of function pointers has to match the
150  *  order of corresponding enumerated electrostatics and vdw types, resp.,
151  *  defined in nbnxn_cuda_types.h.
152  */
153
154 /*! Force-only kernel function pointers. */
155 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_noprune_ptr[eelCuNR][evdwCuNR] =
156 {
157     { 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            },
158     { 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             },
159     { 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        },
160     { 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 },
161     { 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             },
162     { 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      }
163 };
164
165 /*! Force + energy kernel function pointers. */
166 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_noprune_ptr[eelCuNR][evdwCuNR] =
167 {
168     { 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            },
169     { 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             },
170     { 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        },
171     { 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 },
172     { 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             },
173     { 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      }
174 };
175
176 /*! Force + pruning kernel function pointers. */
177 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_prune_ptr[eelCuNR][evdwCuNR] =
178 {
179     { 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             },
180     { 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              },
181     { 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         },
182     { 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  },
183     { 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              },
184     { 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       }
185 };
186
187 /*! Force + energy + pruning kernel function pointers. */
188 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_prune_ptr[eelCuNR][evdwCuNR] =
189 {
190     { 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            },
191     { 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             },
192     { 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        },
193     { 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 },
194     { 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             },
195     { 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      }
196 };
197
198 /*! Return a pointer to the kernel version to be executed at the current step. */
199 static inline nbnxn_cu_kfunc_ptr_t select_nbnxn_kernel(int                                  eeltype,
200                                                        int                                  evdwtype,
201                                                        bool                                 bDoEne,
202                                                        bool                                 bDoPrune,
203                                                        const gmx_device_info_t gmx_unused  *devInfo)
204 {
205     nbnxn_cu_kfunc_ptr_t res;
206
207     GMX_ASSERT(eeltype < eelCuNR,
208                "The electrostatics type requested is not implemented in the CUDA kernels.");
209     GMX_ASSERT(evdwtype < evdwCuNR,
210                "The VdW type requested is not implemented in the CUDA kernels.");
211
212     /* assert assumptions made by the kernels */
213     GMX_ASSERT(c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit == devInfo->prop.warpSize,
214                "The CUDA kernels require the cluster_size_i*cluster_size_j/nbnxn_gpu_clusterpair_split to match the warp size of the architecture targeted.");
215
216     if (bDoEne)
217     {
218         if (bDoPrune)
219         {
220             res = nb_kfunc_ener_prune_ptr[eeltype][evdwtype];
221         }
222         else
223         {
224             res = nb_kfunc_ener_noprune_ptr[eeltype][evdwtype];
225         }
226     }
227     else
228     {
229         if (bDoPrune)
230         {
231             res = nb_kfunc_noener_prune_ptr[eeltype][evdwtype];
232         }
233         else
234         {
235             res = nb_kfunc_noener_noprune_ptr[eeltype][evdwtype];
236         }
237     }
238
239     return res;
240 }
241
242 /*! \brief Calculates the amount of shared memory required by the nonbonded kernel in use. */
243 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)
244 {
245     int shmem;
246
247     assert(dinfo);
248
249     /* size of shmem (force-buffers/xq/atom type preloading) */
250     /* NOTE: with the default kernel on sm3.0 we need shmem only for pre-loading */
251     /* i-atom x+q in shared memory */
252     shmem  = c_numClPerSupercl * c_clSize * sizeof(float4);
253     /* cj in shared memory, for each warp separately */
254     shmem += num_threads_z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(int);
255
256     if (nbp->vdwtype == evdwCuCUTCOMBGEOM ||
257         nbp->vdwtype == evdwCuCUTCOMBLB)
258     {
259         /* i-atom LJ combination parameters in shared memory */
260         shmem += c_numClPerSupercl * c_clSize * sizeof(float2);
261     }
262     else
263     {
264         /* i-atom types in shared memory */
265         shmem += c_numClPerSupercl * c_clSize * sizeof(int);
266     }
267
268     return shmem;
269 }
270
271 /*! As we execute nonbonded workload in separate streams, before launching
272    the kernel we need to make sure that he following operations have completed:
273    - atomdata allocation and related H2D transfers (every nstlist step);
274    - pair list H2D transfer (every nstlist step);
275    - shift vector H2D transfer (every nstlist step);
276    - force (+shift force and energy) output clearing (every step).
277
278    These operations are issued in the local stream at the beginning of the step
279    and therefore always complete before the local kernel launch. The non-local
280    kernel is launched after the local on the same device/context hence it is
281    inherently scheduled after the operations in the local stream (including the
282    above "misc_ops") on pre-GK110 devices with single hardware queue, but on later
283    devices with multiple hardware queues the dependency needs to be enforced.
284    We use the misc_ops_and_local_H2D_done event to record the point where
285    the local x+q H2D (and all preceding) tasks are complete and synchronize
286    with this event in the non-local stream before launching the non-bonded kernel.
287  */
288 void nbnxn_gpu_launch_kernel(gmx_nbnxn_cuda_t       *nb,
289                              const nbnxn_atomdata_t *nbatom,
290                              int                     flags,
291                              int                     iloc)
292 {
293     cudaError_t          stat;
294     int                  adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
295     /* CUDA kernel launch-related stuff */
296     int                  nblock;
297     dim3                 dim_block, dim_grid;
298     nbnxn_cu_kfunc_ptr_t nb_kernel = nullptr; /* fn pointer to the nonbonded kernel */
299
300     cu_atomdata_t       *adat    = nb->atdat;
301     cu_nbparam_t        *nbp     = nb->nbparam;
302     cu_plist_t          *plist   = nb->plist[iloc];
303     cu_timers_t         *t       = nb->timers;
304     cudaStream_t         stream  = nb->stream[iloc];
305
306     bool                 bCalcEner   = flags & GMX_FORCE_ENERGY;
307     bool                 bCalcFshift = flags & GMX_FORCE_VIRIAL;
308     bool                 bDoTime     = nb->bDoTime;
309
310     /* Don't launch the non-local kernel if there is no work to do.
311        Doing the same for the local kernel is more complicated, since the
312        local part of the force array also depends on the non-local kernel.
313        So to avoid complicating the code and to reduce the risk of bugs,
314        we always call the local kernel, the local x+q copy and later (not in
315        this function) the stream wait, local f copyback and the f buffer
316        clearing. All these operations, except for the local interaction kernel,
317        are needed for the non-local interactions. The skip of the local kernel
318        call is taken care of later in this function. */
319     if (canSkipWork(nb, iloc))
320     {
321         plist->haveFreshList = false;
322
323         return;
324     }
325
326     /* calculate the atom data index range based on locality */
327     if (LOCAL_I(iloc))
328     {
329         adat_begin  = 0;
330         adat_len    = adat->natoms_local;
331     }
332     else
333     {
334         adat_begin  = adat->natoms_local;
335         adat_len    = adat->natoms - adat->natoms_local;
336     }
337
338     /* beginning of timed HtoD section */
339     if (bDoTime)
340     {
341         t->nb_h2d[iloc].openTimingRegion(stream);
342     }
343
344     /* HtoD x, q */
345     cu_copy_H2D_async(adat->xq + adat_begin, nbatom->x + adat_begin * 4,
346                       adat_len * sizeof(*adat->xq), stream);
347
348     if (bDoTime)
349     {
350         t->nb_h2d[iloc].closeTimingRegion(stream);
351     }
352
353     /* When we get here all misc operations issues in the local stream as well as
354        the local xq H2D are done,
355        so we record that in the local stream and wait for it in the nonlocal one. */
356     if (nb->bUseTwoStreams)
357     {
358         if (iloc == eintLocal)
359         {
360             stat = cudaEventRecord(nb->misc_ops_and_local_H2D_done, stream);
361             CU_RET_ERR(stat, "cudaEventRecord on misc_ops_and_local_H2D_done failed");
362         }
363         else
364         {
365             stat = cudaStreamWaitEvent(stream, nb->misc_ops_and_local_H2D_done, 0);
366             CU_RET_ERR(stat, "cudaStreamWaitEvent on misc_ops_and_local_H2D_done failed");
367         }
368     }
369
370     if (nbp->useDynamicPruning && plist->haveFreshList)
371     {
372         /* Prunes for rlistOuter and rlistInner, sets plist->haveFreshList=false
373            (TODO: ATM that's the way the timing accounting can distinguish between
374            separate prune kernel and combined force+prune, maybe we need a better way?).
375          */
376         nbnxn_gpu_launch_kernel_pruneonly(nb, iloc, 1);
377     }
378
379     if (plist->nsci == 0)
380     {
381         /* Don't launch an empty local kernel (not allowed with CUDA) */
382         return;
383     }
384
385     /* beginning of timed nonbonded calculation section */
386     if (bDoTime)
387     {
388         t->nb_k[iloc].openTimingRegion(stream);
389     }
390
391     /* get the pointer to the kernel flavor we need to use */
392     nb_kernel = select_nbnxn_kernel(nbp->eeltype,
393                                     nbp->vdwtype,
394                                     bCalcEner,
395                                     (plist->haveFreshList && !nb->timers->didPrune[iloc]),
396                                     nb->dev_info);
397
398     /* Kernel launch config:
399      * - The thread block dimensions match the size of i-clusters, j-clusters,
400      *   and j-cluster concurrency, in x, y, and z, respectively.
401      * - The 1D block-grid contains as many blocks as super-clusters.
402      */
403     int num_threads_z = 1;
404     if (nb->dev_info->prop.major == 3 && nb->dev_info->prop.minor == 7)
405     {
406         num_threads_z = 2;
407     }
408     nblock    = calc_nb_kernel_nblock(plist->nsci, nb->dev_info);
409
410     KernelLaunchConfig config;
411     config.blockSize[0]     = c_clSize;
412     config.blockSize[1]     = c_clSize;
413     config.blockSize[2]     = num_threads_z;
414     config.gridSize[0]      = nblock;
415     config.sharedMemorySize = calc_shmem_required_nonbonded(num_threads_z, nb->dev_info, nbp);
416     config.stream           = stream;
417
418     if (debug)
419     {
420         fprintf(debug, "Non-bonded GPU launch configuration:\n\tThread block: %zux%zux%zu\n\t"
421                 "\tGrid: %zux%zu\n\t#Super-clusters/clusters: %d/%d (%d)\n"
422                 "\tShMem: %zu\n",
423                 config.blockSize[0], config.blockSize[1], config.blockSize[2],
424                 config.gridSize[0], config.gridSize[1], plist->nsci*c_numClPerSupercl,
425                 c_numClPerSupercl, plist->na_c,
426                 config.sharedMemorySize);
427     }
428
429     auto      *timingEvent = bDoTime ? t->nb_k[iloc].fetchNextEvent() : nullptr;
430     const auto kernelArgs  = prepareGpuKernelArguments(nb_kernel, config, adat, nbp, plist, &bCalcFshift);
431     launchGpuKernel(nb_kernel, config, timingEvent, "k_calc_nb", kernelArgs);
432
433     if (bDoTime)
434     {
435         t->nb_k[iloc].closeTimingRegion(stream);
436     }
437
438     if (GMX_NATIVE_WINDOWS)
439     {
440         /* Windows: force flushing WDDM queue */
441         cudaStreamQuery(stream);
442     }
443 }
444
445 /*! Calculates the amount of shared memory required by the CUDA kernel in use. */
446 static inline int calc_shmem_required_prune(const int num_threads_z)
447 {
448     int shmem;
449
450     /* i-atom x in shared memory */
451     shmem  = c_numClPerSupercl * c_clSize * sizeof(float4);
452     /* cj in shared memory, for each warp separately */
453     shmem += num_threads_z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(int);
454
455     return shmem;
456 }
457
458 void nbnxn_gpu_launch_kernel_pruneonly(gmx_nbnxn_cuda_t       *nb,
459                                        int                     iloc,
460                                        int                     numParts)
461 {
462     cu_atomdata_t       *adat    = nb->atdat;
463     cu_nbparam_t        *nbp     = nb->nbparam;
464     cu_plist_t          *plist   = nb->plist[iloc];
465     cu_timers_t         *t       = nb->timers;
466     cudaStream_t         stream  = nb->stream[iloc];
467
468     bool                 bDoTime = nb->bDoTime;
469
470     if (plist->haveFreshList)
471     {
472         GMX_ASSERT(numParts == 1, "With first pruning we expect 1 part");
473
474         /* Set rollingPruningNumParts to signal that it is not set */
475         plist->rollingPruningNumParts = 0;
476         plist->rollingPruningPart     = 0;
477     }
478     else
479     {
480         if (plist->rollingPruningNumParts == 0)
481         {
482             plist->rollingPruningNumParts = numParts;
483         }
484         else
485         {
486             GMX_ASSERT(numParts == plist->rollingPruningNumParts, "It is not allowed to change numParts in between list generation steps");
487         }
488     }
489
490     /* Use a local variable for part and update in plist, so we can return here
491      * without duplicating the part increment code.
492      */
493     int part = plist->rollingPruningPart;
494
495     plist->rollingPruningPart++;
496     if (plist->rollingPruningPart >= plist->rollingPruningNumParts)
497     {
498         plist->rollingPruningPart = 0;
499     }
500
501     /* Compute the number of list entries to prune in this pass */
502     int numSciInPart = (plist->nsci - part)/numParts;
503
504     /* Don't launch the kernel if there is no work to do (not allowed with CUDA) */
505     if (numSciInPart <= 0)
506     {
507         plist->haveFreshList = false;
508
509         return;
510     }
511
512     GpuRegionTimer *timer = nullptr;
513     if (bDoTime)
514     {
515         timer = &(plist->haveFreshList ? t->prune_k[iloc] : t->rollingPrune_k[iloc]);
516     }
517
518     /* beginning of timed prune calculation section */
519     if (bDoTime)
520     {
521         timer->openTimingRegion(stream);
522     }
523
524     /* Kernel launch config:
525      * - The thread block dimensions match the size of i-clusters, j-clusters,
526      *   and j-cluster concurrency, in x, y, and z, respectively.
527      * - The 1D block-grid contains as many blocks as super-clusters.
528      */
529     int                num_threads_z  = c_cudaPruneKernelJ4Concurrency;
530     int                nblock         = calc_nb_kernel_nblock(numSciInPart, nb->dev_info);
531     KernelLaunchConfig config;
532     config.blockSize[0]     = c_clSize;
533     config.blockSize[1]     = c_clSize;
534     config.blockSize[2]     = num_threads_z;
535     config.gridSize[0]      = nblock;
536     config.sharedMemorySize = calc_shmem_required_prune(num_threads_z);
537     config.stream           = stream;
538
539     if (debug)
540     {
541         fprintf(debug, "Pruning GPU kernel launch configuration:\n\tThread block: %zux%zux%zu\n\t"
542                 "\tGrid: %zux%zu\n\t#Super-clusters/clusters: %d/%d (%d)\n"
543                 "\tShMem: %zu\n",
544                 config.blockSize[0], config.blockSize[1], config.blockSize[2],
545                 config.gridSize[0], config.gridSize[1], numSciInPart*c_numClPerSupercl,
546                 c_numClPerSupercl, plist->na_c,
547                 config.sharedMemorySize);
548     }
549
550     auto          *timingEvent  = bDoTime ? timer->fetchNextEvent() : nullptr;
551     constexpr char kernelName[] = "k_pruneonly";
552     const auto    &kernel       = plist->haveFreshList ? nbnxn_kernel_prune_cuda<true> : nbnxn_kernel_prune_cuda<false>;
553     const auto     kernelArgs   = prepareGpuKernelArguments(kernel, config, adat, nbp, plist, &numParts, &part);
554     launchGpuKernel(kernel, config, timingEvent, kernelName, kernelArgs);
555
556     /* TODO: consider a more elegant way to track which kernel has been called
557        (combined or separate 1st pass prune, rolling prune). */
558     if (plist->haveFreshList)
559     {
560         plist->haveFreshList         = false;
561         /* Mark that pruning has been done */
562         nb->timers->didPrune[iloc] = true;
563     }
564     else
565     {
566         /* Mark that rolling pruning has been done */
567         nb->timers->didRollingPrune[iloc] = true;
568     }
569
570     if (bDoTime)
571     {
572         timer->closeTimingRegion(stream);
573     }
574
575     if (GMX_NATIVE_WINDOWS)
576     {
577         /* Windows: force flushing WDDM queue */
578         cudaStreamQuery(stream);
579     }
580 }
581
582 void nbnxn_gpu_launch_cpyback(gmx_nbnxn_cuda_t       *nb,
583                               const nbnxn_atomdata_t *nbatom,
584                               int                     flags,
585                               int                     aloc)
586 {
587     cudaError_t stat;
588     int         adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
589
590     /* determine interaction locality from atom locality */
591     int              iloc = gpuAtomToInteractionLocality(aloc);
592
593     cu_atomdata_t   *adat    = nb->atdat;
594     cu_timers_t     *t       = nb->timers;
595     bool             bDoTime = nb->bDoTime;
596     cudaStream_t     stream  = nb->stream[iloc];
597
598     bool             bCalcEner   = flags & GMX_FORCE_ENERGY;
599     bool             bCalcFshift = flags & GMX_FORCE_VIRIAL;
600
601     /* don't launch non-local copy-back if there was no non-local work to do */
602     if (canSkipWork(nb, iloc))
603     {
604         return;
605     }
606
607     getGpuAtomRange(adat, aloc, &adat_begin, &adat_len);
608
609     /* beginning of timed D2H section */
610     if (bDoTime)
611     {
612         t->nb_d2h[iloc].openTimingRegion(stream);
613     }
614
615     /* With DD the local D2H transfer can only start after the non-local
616        kernel has finished. */
617     if (iloc == eintLocal && nb->bUseTwoStreams)
618     {
619         stat = cudaStreamWaitEvent(stream, nb->nonlocal_done, 0);
620         CU_RET_ERR(stat, "cudaStreamWaitEvent on nonlocal_done failed");
621     }
622
623     /* DtoH f */
624     cu_copy_D2H_async(nbatom->out[0].f + adat_begin * 3, adat->f + adat_begin,
625                       (adat_len)*sizeof(*adat->f), stream);
626
627     /* After the non-local D2H is launched the nonlocal_done event can be
628        recorded which signals that the local D2H can proceed. This event is not
629        placed after the non-local kernel because we want the non-local data
630        back first. */
631     if (iloc == eintNonlocal)
632     {
633         stat = cudaEventRecord(nb->nonlocal_done, stream);
634         CU_RET_ERR(stat, "cudaEventRecord on nonlocal_done failed");
635     }
636
637     /* only transfer energies in the local stream */
638     if (LOCAL_I(iloc))
639     {
640         /* DtoH fshift */
641         if (bCalcFshift)
642         {
643             cu_copy_D2H_async(nb->nbst.fshift, adat->fshift,
644                               SHIFTS * sizeof(*nb->nbst.fshift), stream);
645         }
646
647         /* DtoH energies */
648         if (bCalcEner)
649         {
650             cu_copy_D2H_async(nb->nbst.e_lj, adat->e_lj,
651                               sizeof(*nb->nbst.e_lj), stream);
652             cu_copy_D2H_async(nb->nbst.e_el, adat->e_el,
653                               sizeof(*nb->nbst.e_el), stream);
654         }
655     }
656
657     if (bDoTime)
658     {
659         t->nb_d2h[iloc].closeTimingRegion(stream);
660     }
661 }
662
663 void nbnxn_cuda_set_cacheconfig()
664 {
665     cudaError_t stat;
666
667     for (int i = 0; i < eelCuNR; i++)
668     {
669         for (int j = 0; j < evdwCuNR; j++)
670         {
671             /* Default kernel 32/32 kB Shared/L1 */
672             cudaFuncSetCacheConfig(nb_kfunc_ener_prune_ptr[i][j], cudaFuncCachePreferEqual);
673             cudaFuncSetCacheConfig(nb_kfunc_ener_noprune_ptr[i][j], cudaFuncCachePreferEqual);
674             cudaFuncSetCacheConfig(nb_kfunc_noener_prune_ptr[i][j], cudaFuncCachePreferEqual);
675             stat = cudaFuncSetCacheConfig(nb_kfunc_noener_noprune_ptr[i][j], cudaFuncCachePreferEqual);
676             CU_RET_ERR(stat, "cudaFuncSetCacheConfig failed");
677         }
678     }
679 }