940d7b5b3471c2c06aa48d8545f07ba1519ac3c4
[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, 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 "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"
62
63 #include "nbnxn_cuda_types.h"
64
65 /*
66  * Texture references are created at compile-time and need to be declared
67  * at file scope as global variables (see http://docs.nvidia.com/cuda/cuda-c-programming-guide/index.html#texture-reference-api).
68  * The texture references below are used in two translation units;
69  * we declare them here along the kernels that use them (when compiling legacy Fermi kernels),
70  * and provide getters (see below) used by the data_mgmt module where the
71  * textures are bound/unbound.
72  * (In principle we could do it the other way arond, but that would likely require
73  * device linking and we'd rather avoid technical hurdles.)
74  */
75 /*! Texture reference for LJ C6/C12 parameters; bound to cu_nbparam_t.nbfp */
76 texture<float, 1, cudaReadModeElementType> nbfp_texref;
77
78 /*! Texture reference for LJ-PME parameters; bound to cu_nbparam_t.nbfp_comb */
79 texture<float, 1, cudaReadModeElementType> nbfp_comb_texref;
80
81 /*! Texture reference for Ewald coulomb force table; bound to cu_nbparam_t.coulomb_tab */
82 texture<float, 1, cudaReadModeElementType> coulomb_tab_texref;
83
84
85 /***** The kernel declarations/definitions come here *****/
86
87 /* Top-level kernel declaration generation: will generate through multiple
88  * inclusion the following flavors for all kernel declarations:
89  * - force-only output;
90  * - force and energy output;
91  * - force-only with pair list pruning;
92  * - force and energy output with pair list pruning.
93  */
94 #define FUNCTION_DECLARATION_ONLY
95 /** Force only **/
96 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
97 /** Force & energy **/
98 #define CALC_ENERGIES
99 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
100 #undef CALC_ENERGIES
101
102 /*** Pair-list pruning kernels ***/
103 /** Force only **/
104 #define PRUNE_NBL
105 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
106 /** Force & energy **/
107 #define CALC_ENERGIES
108 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernels.cuh"
109 #undef CALC_ENERGIES
110 #undef PRUNE_NBL
111
112 /* Prune-only kernels */
113 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_pruneonly.cuh"
114 #undef FUNCTION_DECLARATION_ONLY
115
116 /* Now generate the function definitions if we are using a single compilation unit. */
117 #if GMX_CUDA_NB_SINGLE_COMPILATION_UNIT
118 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_F_noprune.cu"
119 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_F_prune.cu"
120 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_VF_noprune.cu"
121 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_VF_prune.cu"
122 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_pruneonly.cu"
123 #else
124 /* Prevent compilation in multiple compilation unit mode for CC 2.x. Although we have
125  * build-time checks to prevent this, the user could manually tweaks nvcc flags
126  * which would lead to buggy kernels getting compiled.
127  */
128 #if GMX_PTX_ARCH > 0 && GMX_PTX_ARCH <= 210
129 #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.
130 #endif
131 #endif /* GMX_CUDA_NB_SINGLE_COMPILATION_UNIT */
132
133
134 /*! Nonbonded kernel function pointer type */
135 typedef void (*nbnxn_cu_kfunc_ptr_t)(const cu_atomdata_t,
136                                      const cu_nbparam_t,
137                                      const cu_plist_t,
138                                      bool);
139
140 /*********************************/
141
142 /* XXX switch between chevron and cudaLaunch (supported only in CUDA >=7.0)
143    -- only for benchmarking purposes */
144 static const bool bUseCudaLaunchKernel =
145     (GMX_CUDA_VERSION >= 7000) && (getenv("GMX_DISABLE_CUDALAUNCH") == NULL);
146
147 /* XXX always/never run the energy/pruning kernels -- only for benchmarking purposes */
148 static bool always_ener  = (getenv("GMX_GPU_ALWAYS_ENER") != NULL);
149 static bool never_ener   = (getenv("GMX_GPU_NEVER_ENER") != NULL);
150 static bool always_prune = (getenv("GMX_GPU_ALWAYS_PRUNE") != NULL);
151
152
153 /*! Returns the number of blocks to be used for the nonbonded GPU kernel. */
154 static inline int calc_nb_kernel_nblock(int nwork_units, const gmx_device_info_t *dinfo)
155 {
156     int max_grid_x_size;
157
158     assert(dinfo);
159     /* CUDA does not accept grid dimension of 0 (which can happen e.g. with an
160        empty domain) and that case should be handled before this point. */
161     assert(nwork_units > 0);
162
163     max_grid_x_size = dinfo->prop.maxGridSize[0];
164
165     /* do we exceed the grid x dimension limit? */
166     if (nwork_units > max_grid_x_size)
167     {
168         gmx_fatal(FARGS, "Watch out, the input system is too large to simulate!\n"
169                   "The number of nonbonded work units (=number of super-clusters) exceeds the"
170                   "maximum grid size in x dimension (%d > %d)!", nwork_units, max_grid_x_size);
171     }
172
173     return nwork_units;
174 }
175
176
177 /* Constant arrays listing all kernel function pointers and enabling selection
178    of a kernel in an elegant manner. */
179
180 /*! Pointers to the non-bonded kernels organized in 2-dim arrays by:
181  *  electrostatics and VDW type.
182  *
183  *  Note that the row- and column-order of function pointers has to match the
184  *  order of corresponding enumerated electrostatics and vdw types, resp.,
185  *  defined in nbnxn_cuda_types.h.
186  */
187
188 /*! Force-only kernel function pointers. */
189 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_noprune_ptr[eelCuNR][evdwCuNR] =
190 {
191     { 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            },
192     { 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             },
193     { 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        },
194     { 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 },
195     { 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             },
196     { 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      }
197 };
198
199 /*! Force + energy kernel function pointers. */
200 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_noprune_ptr[eelCuNR][evdwCuNR] =
201 {
202     { 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            },
203     { 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             },
204     { 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        },
205     { 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 },
206     { 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             },
207     { 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      }
208 };
209
210 /*! Force + pruning kernel function pointers. */
211 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_prune_ptr[eelCuNR][evdwCuNR] =
212 {
213     { 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             },
214     { 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              },
215     { 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         },
216     { 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  },
217     { 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              },
218     { 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       }
219 };
220
221 /*! Force + energy + pruning kernel function pointers. */
222 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_prune_ptr[eelCuNR][evdwCuNR] =
223 {
224     { 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            },
225     { 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             },
226     { 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        },
227     { 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 },
228     { 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             },
229     { 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      }
230 };
231
232 /*! Return a pointer to the kernel version to be executed at the current step. */
233 static inline nbnxn_cu_kfunc_ptr_t select_nbnxn_kernel(int                                  eeltype,
234                                                        int                                  evdwtype,
235                                                        bool                                 bDoEne,
236                                                        bool                                 bDoPrune,
237                                                        const gmx_device_info_t gmx_unused  *devInfo)
238 {
239     nbnxn_cu_kfunc_ptr_t res;
240
241     GMX_ASSERT(eeltype < eelCuNR,
242                "The electrostatics type requested is not implemented in the CUDA kernels.");
243     GMX_ASSERT(evdwtype < evdwCuNR,
244                "The VdW type requested is not implemented in the CUDA kernels.");
245
246     /* assert assumptions made by the kernels */
247     GMX_ASSERT(c_nbnxnGpuClusterSize*c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit == devInfo->prop.warpSize,
248                "The CUDA kernels require the cluster_size_i*cluster_size_j/nbnxn_gpu_clusterpair_split to match the warp size of the architecture targeted.");
249
250     if (bDoEne)
251     {
252         if (bDoPrune)
253         {
254             res = nb_kfunc_ener_prune_ptr[eeltype][evdwtype];
255         }
256         else
257         {
258             res = nb_kfunc_ener_noprune_ptr[eeltype][evdwtype];
259         }
260     }
261     else
262     {
263         if (bDoPrune)
264         {
265             res = nb_kfunc_noener_prune_ptr[eeltype][evdwtype];
266         }
267         else
268         {
269             res = nb_kfunc_noener_noprune_ptr[eeltype][evdwtype];
270         }
271     }
272
273     return res;
274 }
275
276 /*! \brief Calculates the amount of shared memory required by the nonbonded kernel in use. */
277 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)
278 {
279     int shmem;
280
281     assert(dinfo);
282
283     /* size of shmem (force-buffers/xq/atom type preloading) */
284     /* NOTE: with the default kernel on sm3.0 we need shmem only for pre-loading */
285     /* i-atom x+q in shared memory */
286     shmem  = c_numClPerSupercl * c_clSize * sizeof(float4);
287     /* cj in shared memory, for each warp separately */
288     shmem += num_threads_z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(int);
289     if (dinfo->prop.major >= 3)
290     {
291         if (nbp->vdwtype == evdwCuCUTCOMBGEOM ||
292             nbp->vdwtype == evdwCuCUTCOMBLB)
293         {
294             /* i-atom LJ combination parameters in shared memory */
295             shmem += c_numClPerSupercl * c_clSize * sizeof(float2);
296         }
297         else
298         {
299             /* i-atom types in shared memory */
300             shmem += c_numClPerSupercl * c_clSize * sizeof(int);
301         }
302     }
303     if (dinfo->prop.major < 3)
304     {
305         /* force reduction buffers in shared memory */
306         shmem += c_clSize * c_clSize * 3 * sizeof(float);
307     }
308     return shmem;
309 }
310
311 /*! As we execute nonbonded workload in separate streams, before launching
312    the kernel we need to make sure that he following operations have completed:
313    - atomdata allocation and related H2D transfers (every nstlist step);
314    - pair list H2D transfer (every nstlist step);
315    - shift vector H2D transfer (every nstlist step);
316    - force (+shift force and energy) output clearing (every step).
317
318    These operations are issued in the local stream at the beginning of the step
319    and therefore always complete before the local kernel launch. The non-local
320    kernel is launched after the local on the same device/context hence it is
321    inherently scheduled after the operations in the local stream (including the
322    above "misc_ops") on pre-GK110 devices with single hardware queue, but on later
323    devices with multiple hardware queues the dependency needs to be enforced.
324    We use the misc_ops_and_local_H2D_done event to record the point where
325    the local x+q H2D (and all preceding) tasks are complete and synchronize
326    with this event in the non-local stream before launching the non-bonded kernel.
327  */
328 void nbnxn_gpu_launch_kernel(gmx_nbnxn_cuda_t       *nb,
329                              const nbnxn_atomdata_t *nbatom,
330                              int                     flags,
331                              int                     iloc)
332 {
333     cudaError_t          stat;
334     int                  adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
335     /* CUDA kernel launch-related stuff */
336     int                  shmem, nblock;
337     dim3                 dim_block, dim_grid;
338     nbnxn_cu_kfunc_ptr_t nb_kernel = NULL; /* fn pointer to the nonbonded kernel */
339
340     cu_atomdata_t       *adat    = nb->atdat;
341     cu_nbparam_t        *nbp     = nb->nbparam;
342     cu_plist_t          *plist   = nb->plist[iloc];
343     cu_timers_t         *t       = nb->timers;
344     cudaStream_t         stream  = nb->stream[iloc];
345
346     bool                 bCalcEner   = flags & GMX_FORCE_ENERGY;
347     bool                 bCalcFshift = flags & GMX_FORCE_VIRIAL;
348     bool                 bDoTime     = nb->bDoTime;
349
350     /* turn energy calculation always on/off (for debugging/testing only) */
351     bCalcEner = (bCalcEner || always_ener) && !never_ener;
352
353     /* Don't launch the non-local kernel if there is no work to do.
354        Doing the same for the local kernel is more complicated, since the
355        local part of the force array also depends on the non-local kernel.
356        So to avoid complicating the code and to reduce the risk of bugs,
357        we always call the local kernel, the local x+q copy and later (not in
358        this function) the stream wait, local f copyback and the f buffer
359        clearing. All these operations, except for the local interaction kernel,
360        are needed for the non-local interactions. The skip of the local kernel
361        call is taken care of later in this function. */
362     if (iloc == eintNonlocal && plist->nsci == 0)
363     {
364         plist->haveFreshList = false;
365
366         return;
367     }
368
369     /* calculate the atom data index range based on locality */
370     if (LOCAL_I(iloc))
371     {
372         adat_begin  = 0;
373         adat_len    = adat->natoms_local;
374     }
375     else
376     {
377         adat_begin  = adat->natoms_local;
378         adat_len    = adat->natoms - adat->natoms_local;
379     }
380
381     /* beginning of timed HtoD section */
382     if (bDoTime)
383     {
384         stat = cudaEventRecord(t->start_nb_h2d[iloc], stream);
385         CU_RET_ERR(stat, "cudaEventRecord failed");
386     }
387
388     /* HtoD x, q */
389     cu_copy_H2D_async(adat->xq + adat_begin, nbatom->x + adat_begin * 4,
390                       adat_len * sizeof(*adat->xq), stream);
391
392     /* When we get here all misc operations issues in the local stream as well as
393        the local xq H2D are done,
394        so we record that in the local stream and wait for it in the nonlocal one. */
395     if (nb->bUseTwoStreams)
396     {
397         if (iloc == eintLocal)
398         {
399             stat = cudaEventRecord(nb->misc_ops_and_local_H2D_done, stream);
400             CU_RET_ERR(stat, "cudaEventRecord on misc_ops_and_local_H2D_done failed");
401         }
402         else
403         {
404             stat = cudaStreamWaitEvent(stream, nb->misc_ops_and_local_H2D_done, 0);
405             CU_RET_ERR(stat, "cudaStreamWaitEvent on misc_ops_and_local_H2D_done failed");
406         }
407     }
408
409     if (bDoTime)
410     {
411         stat = cudaEventRecord(t->stop_nb_h2d[iloc], stream);
412         CU_RET_ERR(stat, "cudaEventRecord failed");
413     }
414
415     if (nbp->useDynamicPruning && plist->haveFreshList)
416     {
417         /* Prunes for rlistOuter and rlistInner, sets plist->haveFreshList=false
418            (TODO: ATM that's the way the timing accounting can distinguish between
419            separate prune kernel and combined force+prune, maybe we need a better way?).
420          */
421         nbnxn_gpu_launch_kernel_pruneonly(nb, iloc, 1);
422     }
423
424     if (plist->nsci == 0)
425     {
426         /* Don't launch an empty local kernel (not allowed with CUDA) */
427         return;
428     }
429
430     /* beginning of timed nonbonded calculation section */
431     if (bDoTime)
432     {
433         stat = cudaEventRecord(t->start_nb_k[iloc], stream);
434         CU_RET_ERR(stat, "cudaEventRecord failed");
435     }
436
437     /* get the pointer to the kernel flavor we need to use */
438     nb_kernel = select_nbnxn_kernel(nbp->eeltype,
439                                     nbp->vdwtype,
440                                     bCalcEner,
441                                     (plist->haveFreshList && !nb->timers->didPrune[iloc]) || always_prune,
442                                     nb->dev_info);
443
444     /* Kernel launch config:
445      * - The thread block dimensions match the size of i-clusters, j-clusters,
446      *   and j-cluster concurrency, in x, y, and z, respectively.
447      * - The 1D block-grid contains as many blocks as super-clusters.
448      */
449     int num_threads_z = 1;
450     if (nb->dev_info->prop.major == 3 && nb->dev_info->prop.minor == 7)
451     {
452         num_threads_z = 2;
453     }
454     nblock    = calc_nb_kernel_nblock(plist->nsci, nb->dev_info);
455     dim_block = dim3(c_clSize, c_clSize, num_threads_z);
456     dim_grid  = dim3(nblock, 1, 1);
457     shmem     = calc_shmem_required_nonbonded(num_threads_z, nb->dev_info, nbp);
458
459     if (debug)
460     {
461         fprintf(debug, "Non-bonded GPU launch configuration:\n\tThread block: %ux%ux%u\n\t"
462                 "\tGrid: %ux%u\n\t#Super-clusters/clusters: %d/%d (%d)\n"
463                 "\tShMem: %d\n",
464                 dim_block.x, dim_block.y, dim_block.z,
465                 dim_grid.x, dim_grid.y, plist->nsci*c_numClPerSupercl,
466                 c_numClPerSupercl, plist->na_c,
467                 shmem);
468     }
469
470     if (bUseCudaLaunchKernel)
471     {
472         gmx_unused void* kernel_args[4];
473         kernel_args[0] = adat;
474         kernel_args[1] = nbp;
475         kernel_args[2] = plist;
476         kernel_args[3] = &bCalcFshift;
477
478 #if GMX_CUDA_VERSION >= 7000
479         cudaLaunchKernel((void *)nb_kernel, dim_grid, dim_block, kernel_args, shmem, stream);
480 #endif
481     }
482     else
483     {
484         nb_kernel<<< dim_grid, dim_block, shmem, stream>>> (*adat, *nbp, *plist, bCalcFshift);
485     }
486     CU_LAUNCH_ERR("k_calc_nb");
487
488     if (bDoTime)
489     {
490         stat = cudaEventRecord(t->stop_nb_k[iloc], stream);
491         CU_RET_ERR(stat, "cudaEventRecord failed");
492     }
493
494 #if (defined(WIN32) || defined( _WIN32 ))
495     /* Windows: force flushing WDDM queue */
496     stat = cudaStreamQuery(stream);
497 #endif
498 }
499
500 /*! Calculates the amount of shared memory required by the CUDA kernel in use. */
501 static inline int calc_shmem_required_prune(const int num_threads_z)
502 {
503     int shmem;
504
505     /* i-atom x in shared memory */
506     shmem  = c_numClPerSupercl * c_clSize * sizeof(float4);
507     /* cj in shared memory, for each warp separately */
508     shmem += num_threads_z * c_nbnxnGpuClusterpairSplit * c_nbnxnGpuJgroupSize * sizeof(int);
509
510     return shmem;
511 }
512
513 void nbnxn_gpu_launch_kernel_pruneonly(gmx_nbnxn_cuda_t       *nb,
514                                        int                     iloc,
515                                        int                     numParts)
516 {
517     cudaError_t          stat;
518
519     cu_atomdata_t       *adat    = nb->atdat;
520     cu_nbparam_t        *nbp     = nb->nbparam;
521     cu_plist_t          *plist   = nb->plist[iloc];
522     cu_timers_t         *t       = nb->timers;
523     cudaStream_t         stream  = nb->stream[iloc];
524
525     bool                 bDoTime = nb->bDoTime;
526
527     if (plist->haveFreshList)
528     {
529         GMX_ASSERT(numParts == 1, "With first pruning we expect 1 part");
530
531         /* Set rollingPruningNumParts to signal that it is not set */
532         plist->rollingPruningNumParts = 0;
533         plist->rollingPruningPart     = 0;
534     }
535     else
536     {
537         if (plist->rollingPruningNumParts == 0)
538         {
539             plist->rollingPruningNumParts = numParts;
540         }
541         else
542         {
543             GMX_ASSERT(numParts == plist->rollingPruningNumParts, "It is not allowed to change numParts in between list generation steps");
544         }
545     }
546
547     /* Use a local variable for part and update in plist, so we can return here
548      * without duplicating the part increment code.
549      */
550     int part = plist->rollingPruningPart;
551
552     plist->rollingPruningPart++;
553     if (plist->rollingPruningPart >= plist->rollingPruningNumParts)
554     {
555         plist->rollingPruningPart = 0;
556     }
557
558     /* Compute the number of list entries to prune in this pass */
559     int numSciInPart = (plist->nsci - part)/numParts;
560
561     /* Don't launch the kernel if there is no work to do (not allowed with CUDA) */
562     if (numSciInPart <= 0)
563     {
564         plist->haveFreshList = false;
565
566         return;
567     }
568
569     cudaEvent_t startEvent, stopEvent;
570     if (bDoTime)
571     {
572         startEvent = (plist->haveFreshList ? t->start_prune_k[iloc] : t->start_rollingPrune_k[iloc]);
573         stopEvent  = (plist->haveFreshList ? t->stop_prune_k[iloc]  : t->stop_rollingPrune_k[iloc]);
574     }
575
576     /* beginning of timed prune calculation section */
577     if (bDoTime)
578     {
579         stat = cudaEventRecord(startEvent, stream);
580         CU_RET_ERR(stat, "cudaEventRecord failed");
581     }
582
583     /* Kernel launch config:
584      * - The thread block dimensions match the size of i-clusters, j-clusters,
585      *   and j-cluster concurrency, in x, y, and z, respectively.
586      * - The 1D block-grid contains as many blocks as super-clusters.
587      */
588     int  num_threads_z  = c_cudaPruneKernelJ4Concurrency;
589     int  nblock         = calc_nb_kernel_nblock(numSciInPart, nb->dev_info);
590     dim3 dim_block      = dim3(c_clSize, c_clSize, num_threads_z);
591     dim3 dim_grid       = dim3(nblock, 1, 1);
592     int  shmem          = calc_shmem_required_prune(num_threads_z);
593
594     if (debug)
595     {
596         fprintf(debug, "Pruning GPU kernel launch configuration:\n\tThread block: %dx%dx%d\n\t"
597                 "\tGrid: %dx%d\n\t#Super-clusters/clusters: %d/%d (%d)\n"
598                 "\tShMem: %d\n",
599                 dim_block.x, dim_block.y, dim_block.z,
600                 dim_grid.x, dim_grid.y, numSciInPart*c_numClPerSupercl,
601                 c_numClPerSupercl, plist->na_c,
602                 shmem);
603     }
604
605     if (bUseCudaLaunchKernel)
606     {
607         gmx_unused void* kernel_args[5];
608         kernel_args[0] = adat;
609         kernel_args[1] = nbp;
610         kernel_args[2] = plist;
611         kernel_args[3] = &numParts;
612         kernel_args[4] = &part;
613
614 #if GMX_CUDA_VERSION >= 7000
615         if (plist->haveFreshList)
616         {
617             cudaLaunchKernel((void *)nbnxn_kernel_prune_cuda<true>, dim_grid, dim_block, kernel_args, shmem, stream);
618         }
619         else
620         {
621             cudaLaunchKernel((void *)nbnxn_kernel_prune_cuda<false>, dim_grid, dim_block, kernel_args, shmem, stream);
622         }
623 #endif
624     }
625     else
626     {
627         if (plist->haveFreshList)
628         {
629             nbnxn_kernel_prune_cuda<true><<< dim_grid, dim_block, shmem, stream>>> (*adat, *nbp, *plist, numParts, part);
630         }
631         else
632         {
633             nbnxn_kernel_prune_cuda<false><<< dim_grid, dim_block, shmem, stream>>> (*adat, *nbp, *plist, numParts, part);
634         }
635     }
636     CU_LAUNCH_ERR("k_pruneonly");
637
638     /* TODO: consider a more elegant way to track which kernel has been called
639        (combined or separate 1st pass prune, rolling prune). */
640     if (plist->haveFreshList)
641     {
642         plist->haveFreshList         = false;
643         /* Mark that pruning has been done */
644         nb->timers->didPrune[iloc] = true;
645     }
646     else
647     {
648         /* Mark that rolling pruning has been done */
649         nb->timers->didRollingPrune[iloc] = true;
650     }
651
652     if (bDoTime)
653     {
654         stat = cudaEventRecord(stopEvent, stream);
655         CU_RET_ERR(stat, "cudaEventRecord failed");
656     }
657
658 #if (defined(WIN32) || defined( _WIN32 ))
659     /* Windows: force flushing WDDM queue */
660     stat = cudaStreamQuery(stream);
661 #endif
662 }
663
664 void nbnxn_gpu_launch_cpyback(gmx_nbnxn_cuda_t       *nb,
665                               const nbnxn_atomdata_t *nbatom,
666                               int                     flags,
667                               int                     aloc)
668 {
669     cudaError_t stat;
670     int         adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
671     int         iloc = -1;
672
673     /* determine interaction locality from atom locality */
674     if (LOCAL_A(aloc))
675     {
676         iloc = eintLocal;
677     }
678     else if (NONLOCAL_A(aloc))
679     {
680         iloc = eintNonlocal;
681     }
682     else
683     {
684         char stmp[STRLEN];
685         sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
686                 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
687         gmx_incons(stmp);
688     }
689
690     cu_atomdata_t   *adat    = nb->atdat;
691     cu_timers_t     *t       = nb->timers;
692     bool             bDoTime = nb->bDoTime;
693     cudaStream_t     stream  = nb->stream[iloc];
694
695     bool             bCalcEner   = flags & GMX_FORCE_ENERGY;
696     bool             bCalcFshift = flags & GMX_FORCE_VIRIAL;
697
698     /* don't launch non-local copy-back if there was no non-local work to do */
699     if (iloc == eintNonlocal && nb->plist[iloc]->nsci == 0)
700     {
701         return;
702     }
703
704     /* calculate the atom data index range based on locality */
705     if (LOCAL_A(aloc))
706     {
707         adat_begin  = 0;
708         adat_len    = adat->natoms_local;
709     }
710     else
711     {
712         adat_begin  = adat->natoms_local;
713         adat_len    = adat->natoms - adat->natoms_local;
714     }
715
716     /* beginning of timed D2H section */
717     if (bDoTime)
718     {
719         stat = cudaEventRecord(t->start_nb_d2h[iloc], stream);
720         CU_RET_ERR(stat, "cudaEventRecord failed");
721     }
722
723     /* With DD the local D2H transfer can only start after the non-local
724        kernel has finished. */
725     if (iloc == eintLocal && nb->bUseTwoStreams)
726     {
727         stat = cudaStreamWaitEvent(stream, nb->nonlocal_done, 0);
728         CU_RET_ERR(stat, "cudaStreamWaitEvent on nonlocal_done failed");
729     }
730
731     /* DtoH f */
732     cu_copy_D2H_async(nbatom->out[0].f + adat_begin * 3, adat->f + adat_begin,
733                       (adat_len)*sizeof(*adat->f), stream);
734
735     /* After the non-local D2H is launched the nonlocal_done event can be
736        recorded which signals that the local D2H can proceed. This event is not
737        placed after the non-local kernel because we want the non-local data
738        back first. */
739     if (iloc == eintNonlocal)
740     {
741         stat = cudaEventRecord(nb->nonlocal_done, stream);
742         CU_RET_ERR(stat, "cudaEventRecord on nonlocal_done failed");
743     }
744
745     /* only transfer energies in the local stream */
746     if (LOCAL_I(iloc))
747     {
748         /* DtoH fshift */
749         if (bCalcFshift)
750         {
751             cu_copy_D2H_async(nb->nbst.fshift, adat->fshift,
752                               SHIFTS * sizeof(*nb->nbst.fshift), stream);
753         }
754
755         /* DtoH energies */
756         if (bCalcEner)
757         {
758             cu_copy_D2H_async(nb->nbst.e_lj, adat->e_lj,
759                               sizeof(*nb->nbst.e_lj), stream);
760             cu_copy_D2H_async(nb->nbst.e_el, adat->e_el,
761                               sizeof(*nb->nbst.e_el), stream);
762         }
763     }
764
765     if (bDoTime)
766     {
767         stat = cudaEventRecord(t->stop_nb_d2h[iloc], stream);
768         CU_RET_ERR(stat, "cudaEventRecord failed");
769     }
770 }
771
772 /*! \brief Count pruning kernel time if either kernel has been triggered
773  *
774  *  We do the accounting for either of the two pruning kernel flavors:
775  *   - 1st pass prune: ran during the current step (prior to the force kernel);
776  *   - rolling prune:  ran at the end of the previous step (prior to the current step H2D xq);
777  *
778  * Note that the resetting of cu_timers_t::didPrune and cu_timers_t::didRollingPrune should happen
779  * after calling this function.
780  *
781  * \param[in] timers      structs with CUDA timer objects
782  * \param[inout] timings  GPU task timing data
783  * \param[in] iloc        interaction locality
784  */
785 static void countPruneKernelTime(const cu_timers_t   *timers,
786                                  gmx_wallclock_gpu_t *timings,
787                                  const int            iloc)
788 {
789     // We might have not done any pruning (e.g. if we skipped with empty domains).
790     if (!timers->didPrune[iloc] && !timers->didRollingPrune[iloc])
791     {
792         return;
793     }
794
795     if (timers->didPrune[iloc])
796     {
797         timings->pruneTime.c++;
798         timings->pruneTime.t += cu_event_elapsed(timers->start_prune_k[iloc],
799                                                  timers->stop_prune_k[iloc]);
800     }
801     if (timers->didRollingPrune[iloc])
802     {
803         timings->dynamicPruneTime.c++;
804         timings->dynamicPruneTime.t += cu_event_elapsed(timers->start_rollingPrune_k[iloc],
805                                                         timers->stop_rollingPrune_k[iloc]);
806     }
807 }
808
809 void nbnxn_gpu_wait_for_gpu(gmx_nbnxn_cuda_t *nb,
810                             int flags, int aloc,
811                             real *e_lj, real *e_el, rvec *fshift)
812 {
813     /* NOTE:  only implemented for single-precision at this time */
814     cudaError_t stat;
815     int         iloc = -1;
816
817     /* determine interaction locality from atom locality */
818     if (LOCAL_A(aloc))
819     {
820         iloc = eintLocal;
821     }
822     else if (NONLOCAL_A(aloc))
823     {
824         iloc = eintNonlocal;
825     }
826     else
827     {
828         char stmp[STRLEN];
829         sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
830                 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
831         gmx_incons(stmp);
832     }
833
834     cu_plist_t                 *plist    = nb->plist[iloc];
835     cu_timers_t                *timers   = nb->timers;
836     struct gmx_wallclock_gpu_t *timings  = nb->timings;
837     nb_staging                  nbst     = nb->nbst;
838
839     bool                        bCalcEner   = flags & GMX_FORCE_ENERGY;
840     bool                        bCalcFshift = flags & GMX_FORCE_VIRIAL;
841
842     /* turn energy calculation always on/off (for debugging/testing only) */
843     bCalcEner = (bCalcEner || always_ener) && !never_ener;
844
845     /* Launch wait/update timers & counters and do reduction into staging buffers
846        BUT skip it when during the non-local phase there was actually no work to do.
847        This is consistent with nbnxn_gpu_launch_kernel.
848
849        NOTE: if timing with multiple GPUs (streams) becomes possible, the
850        counters could end up being inconsistent due to not being incremented
851        on some of the nodes! */
852     if (!(iloc == eintNonlocal && nb->plist[iloc]->nsci == 0))
853     {
854         stat = cudaStreamSynchronize(nb->stream[iloc]);
855         CU_RET_ERR(stat, "cudaStreamSynchronize failed in cu_blockwait_nb");
856
857         /* timing data accumulation */
858         if (nb->bDoTime)
859         {
860             /* only increase counter once (at local F wait) */
861             if (LOCAL_I(iloc))
862             {
863                 timings->nb_c++;
864                 timings->ktime[plist->haveFreshList ? 1 : 0][bCalcEner ? 1 : 0].c += 1;
865             }
866
867             /* kernel timings */
868             timings->ktime[plist->haveFreshList ? 1 : 0][bCalcEner ? 1 : 0].t +=
869                 cu_event_elapsed(timers->start_nb_k[iloc], timers->stop_nb_k[iloc]);
870
871             /* X/q H2D and F D2H timings */
872             timings->nb_h2d_t += cu_event_elapsed(timers->start_nb_h2d[iloc],
873                                                   timers->stop_nb_h2d[iloc]);
874             timings->nb_d2h_t += cu_event_elapsed(timers->start_nb_d2h[iloc],
875                                                   timers->stop_nb_d2h[iloc]);
876
877             /* Count the pruning kernel times for both cases:1st pass (at search step)
878                and rolling pruning (if called at the previous step).
879                We do the accounting here as this is the only sync point where we
880                know (without checking or additional sync-ing) that prune tasks in
881                in the current stream have completed (having just blocking-waited
882                for the force D2H). */
883             countPruneKernelTime(timers, timings, iloc);
884
885             /* only count atdat and pair-list H2D at pair-search step */
886             if (timers->didPairlistH2D[iloc])
887             {
888                 /* atdat transfer timing (add only once, at local F wait) */
889                 if (LOCAL_A(aloc))
890                 {
891                     timings->pl_h2d_c++;
892                     timings->pl_h2d_t += cu_event_elapsed(timers->start_atdat,
893                                                           timers->stop_atdat);
894                 }
895
896                 timings->pl_h2d_t += cu_event_elapsed(timers->start_pl_h2d[iloc],
897                                                       timers->stop_pl_h2d[iloc]);
898
899                 /* Clear the timing flag for the next step */
900                 timers->didPairlistH2D[iloc] = false;
901             }
902         }
903
904         /* add up energies and shift forces (only once at local F wait) */
905         if (LOCAL_I(iloc))
906         {
907             if (bCalcEner)
908             {
909                 *e_lj += *nbst.e_lj;
910                 *e_el += *nbst.e_el;
911             }
912
913             if (bCalcFshift)
914             {
915                 for (int i = 0; i < SHIFTS; i++)
916                 {
917                     fshift[i][0] += nbst.fshift[i].x;
918                     fshift[i][1] += nbst.fshift[i].y;
919                     fshift[i][2] += nbst.fshift[i].z;
920                 }
921             }
922         }
923     }
924
925     /* Always reset both pruning flags (doesn't hurt doing it even when timing is off). */
926     timers->didPrune[iloc] = timers->didRollingPrune[iloc] = false;
927
928     /* Turn off initial list pruning (doesn't hurt if this is not pair-search step). */
929     plist->haveFreshList = false;
930 }
931
932 /*! Return the reference to the nbfp texture. */
933 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_texref()
934 {
935     assert(!c_disableCudaTextures);
936     return nbfp_texref;
937 }
938
939 /*! Return the reference to the nbfp_comb texture. */
940 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_comb_texref()
941 {
942     assert(!c_disableCudaTextures);
943     return nbfp_comb_texref;
944 }
945
946 /*! Return the reference to the coulomb_tab. */
947 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_coulomb_tab_texref()
948 {
949     assert(!c_disableCudaTextures);
950     return coulomb_tab_texref;
951 }
952
953 /*! Set up the cache configuration for the non-bonded kernels,
954  */
955 void nbnxn_cuda_set_cacheconfig(const gmx_device_info_t *devinfo)
956 {
957     cudaError_t stat;
958
959     for (int i = 0; i < eelCuNR; i++)
960     {
961         for (int j = 0; j < evdwCuNR; j++)
962         {
963             if (devinfo->prop.major >= 3)
964             {
965                 /* Default kernel on sm 3.x and later 32/32 kB Shared/L1 */
966                 cudaFuncSetCacheConfig(nb_kfunc_ener_prune_ptr[i][j], cudaFuncCachePreferEqual);
967                 cudaFuncSetCacheConfig(nb_kfunc_ener_noprune_ptr[i][j], cudaFuncCachePreferEqual);
968                 cudaFuncSetCacheConfig(nb_kfunc_noener_prune_ptr[i][j], cudaFuncCachePreferEqual);
969                 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_noprune_ptr[i][j], cudaFuncCachePreferEqual);
970             }
971             else
972             {
973                 /* On Fermi prefer L1 gives 2% higher performance */
974                 /* Default kernel on sm_2.x 16/48 kB Shared/L1 */
975                 cudaFuncSetCacheConfig(nb_kfunc_ener_prune_ptr[i][j], cudaFuncCachePreferL1);
976                 cudaFuncSetCacheConfig(nb_kfunc_ener_noprune_ptr[i][j], cudaFuncCachePreferL1);
977                 cudaFuncSetCacheConfig(nb_kfunc_noener_prune_ptr[i][j], cudaFuncCachePreferL1);
978                 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_noprune_ptr[i][j], cudaFuncCachePreferL1);
979             }
980             CU_RET_ERR(stat, "cudaFuncSetCacheConfig failed");
981         }
982     }
983 }