2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, 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.
48 #include "types/simple.h"
49 #include "types/nbnxn_pairlist.h"
50 #include "types/nb_verlet.h"
51 #include "types/ishift.h"
52 #include "types/force_flags.h"
53 #include "../nbnxn_consts.h"
56 #include "thread_mpi/atomic.h"
59 #include "nbnxn_cuda_types.h"
60 #include "../../gmxlib/cuda_tools/cudautils.cuh"
61 #include "nbnxn_cuda.h"
62 #include "nbnxn_cuda_data_mgmt.h"
64 #if defined TEXOBJ_SUPPORTED && __CUDA_ARCH__ >= 300
68 /*! Texture reference for LJ C6/C12 parameters; bound to cu_nbparam_t.nbfp */
69 texture<float, 1, cudaReadModeElementType> nbfp_texref;
71 /*! Texture reference for LJ-PME parameters; bound to cu_nbparam_t.nbfp_comb */
72 texture<float, 1, cudaReadModeElementType> nbfp_comb_texref;
74 /*! Texture reference for Ewald coulomb force table; bound to cu_nbparam_t.coulomb_tab */
75 texture<float, 1, cudaReadModeElementType> coulomb_tab_texref;
77 /* Convenience defines */
78 #define NCL_PER_SUPERCL (NBNXN_GPU_NCLUSTER_PER_SUPERCLUSTER)
79 #define CL_SIZE (NBNXN_GPU_CLUSTER_SIZE)
81 /***** The kernels come here *****/
82 #include "nbnxn_cuda_kernel_utils.cuh"
84 /* Top-level kernel generation: will generate through multiple inclusion the
85 * following flavors for all kernels:
86 * - force-only output;
87 * - force and energy output;
88 * - force-only with pair list pruning;
89 * - force and energy output with pair list pruning.
92 #include "nbnxn_cuda_kernels.cuh"
93 /** Force & energy **/
95 #include "nbnxn_cuda_kernels.cuh"
98 /*** Pair-list pruning kernels ***/
101 #include "nbnxn_cuda_kernels.cuh"
102 /** Force & energy **/
103 #define CALC_ENERGIES
104 #include "nbnxn_cuda_kernels.cuh"
108 /*! Nonbonded kernel function pointer type */
109 typedef void (*nbnxn_cu_kfunc_ptr_t)(const cu_atomdata_t,
114 /*********************************/
116 /* XXX always/never run the energy/pruning kernels -- only for benchmarking purposes */
117 static bool always_ener = (getenv("GMX_GPU_ALWAYS_ENER") != NULL);
118 static bool never_ener = (getenv("GMX_GPU_NEVER_ENER") != NULL);
119 static bool always_prune = (getenv("GMX_GPU_ALWAYS_PRUNE") != NULL);
122 /* Bit-pattern used for polling-based GPU synchronization. It is used as a float
123 * and corresponds to having the exponent set to the maximum (127 -- single
124 * precision) and the mantissa to 0.
126 static unsigned int poll_wait_pattern = (0x7FU << 23);
128 /*! Returns the number of blocks to be used for the nonbonded GPU kernel. */
129 static inline int calc_nb_kernel_nblock(int nwork_units, cuda_dev_info_t *dinfo)
135 max_grid_x_size = dinfo->prop.maxGridSize[0];
137 /* do we exceed the grid x dimension limit? */
138 if (nwork_units > max_grid_x_size)
140 gmx_fatal(FARGS, "Watch out, the input system is too large to simulate!\n"
141 "The number of nonbonded work units (=number of super-clusters) exceeds the"
142 "maximum grid size in x dimension (%d > %d)!", nwork_units, max_grid_x_size);
149 /* Constant arrays listing all kernel function pointers and enabling selection
150 of a kernel in an elegant manner. */
152 /*! Pointers to the non-bonded kernels organized in 2-dim arrays by:
153 * electrostatics and VDW type.
155 * Note that the row- and column-order of function pointers has to match the
156 * order of corresponding enumerated electrostatics and vdw types, resp.,
157 * defined in nbnxn_cuda_types.h.
160 /*! Force-only kernel function pointers. */
161 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_noprune_ptr[eelCuNR][evdwCuNR] =
163 { nbnxn_kernel_ElecCut_VdwLJ_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 },
164 { nbnxn_kernel_ElecRF_VdwLJ_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 },
165 { nbnxn_kernel_ElecEwQSTab_VdwLJ_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 },
166 { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_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 },
167 { nbnxn_kernel_ElecEw_VdwLJ_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 },
168 { nbnxn_kernel_ElecEwTwinCut_VdwLJ_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 }
171 /*! Force + energy kernel function pointers. */
172 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_noprune_ptr[eelCuNR][evdwCuNR] =
174 { nbnxn_kernel_ElecCut_VdwLJ_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 },
175 { nbnxn_kernel_ElecRF_VdwLJ_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 },
176 { nbnxn_kernel_ElecEwQSTab_VdwLJ_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 },
177 { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_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 },
178 { nbnxn_kernel_ElecEw_VdwLJ_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 },
179 { nbnxn_kernel_ElecEwTwinCut_VdwLJ_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 }
182 /*! Force + pruning kernel function pointers. */
183 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_prune_ptr[eelCuNR][evdwCuNR] =
185 { nbnxn_kernel_ElecCut_VdwLJ_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 },
186 { nbnxn_kernel_ElecRF_VdwLJ_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 },
187 { nbnxn_kernel_ElecEwQSTab_VdwLJ_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 },
188 { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_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 },
189 { nbnxn_kernel_ElecEw_VdwLJ_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 },
190 { nbnxn_kernel_ElecEwTwinCut_VdwLJ_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 }
193 /*! Force + energy + pruning kernel function pointers. */
194 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_prune_ptr[eelCuNR][evdwCuNR] =
196 { nbnxn_kernel_ElecCut_VdwLJ_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 },
197 { nbnxn_kernel_ElecRF_VdwLJ_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 },
198 { nbnxn_kernel_ElecEwQSTab_VdwLJ_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 },
199 { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_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 },
200 { nbnxn_kernel_ElecEw_VdwLJ_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 },
201 { nbnxn_kernel_ElecEwTwinCut_VdwLJ_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 }
204 /*! Return a pointer to the kernel version to be executed at the current step. */
205 static inline nbnxn_cu_kfunc_ptr_t select_nbnxn_kernel(int eeltype,
210 nbnxn_cu_kfunc_ptr_t res;
212 assert(eeltype < eelCuNR);
213 assert(evdwtype < eelCuNR);
219 res = nb_kfunc_ener_prune_ptr[eeltype][evdwtype];
223 res = nb_kfunc_ener_noprune_ptr[eeltype][evdwtype];
230 res = nb_kfunc_noener_prune_ptr[eeltype][evdwtype];
234 res = nb_kfunc_noener_noprune_ptr[eeltype][evdwtype];
241 /*! Calculates the amount of shared memory required by the CUDA kernel in use. */
242 static inline int calc_shmem_required()
246 /* size of shmem (force-buffers/xq/atom type preloading) */
247 /* NOTE: with the default kernel on sm3.0 we need shmem only for pre-loading */
248 /* i-atom x+q in shared memory */
249 shmem = NCL_PER_SUPERCL * CL_SIZE * sizeof(float4);
250 /* cj in shared memory, for both warps separately */
251 shmem += 2 * NBNXN_GPU_JGROUP_SIZE * sizeof(int);
253 /* i-atom types in shared memory */
254 shmem += NCL_PER_SUPERCL * CL_SIZE * sizeof(int);
256 #if __CUDA_ARCH__ < 300
257 /* force reduction buffers in shared memory */
258 shmem += CL_SIZE * CL_SIZE * 3 * sizeof(float);
264 /*! As we execute nonbonded workload in separate streams, before launching
265 the kernel we need to make sure that he following operations have completed:
266 - atomdata allocation and related H2D transfers (every nstlist step);
267 - pair list H2D transfer (every nstlist step);
268 - shift vector H2D transfer (every nstlist step);
269 - force (+shift force and energy) output clearing (every step).
271 These operations are issued in the local stream at the beginning of the step
272 and therefore always complete before the local kernel launch. The non-local
273 kernel is launched after the local on the same device/context, so this is
274 inherently scheduled after the operations in the local stream (including the
276 However, for the sake of having a future-proof implementation, we use the
277 misc_ops_done event to record the point in time when the above operations
278 are finished and synchronize with this event in the non-local stream.
280 void nbnxn_cuda_launch_kernel(nbnxn_cuda_ptr_t cu_nb,
281 const nbnxn_atomdata_t *nbatom,
286 int adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
287 /* CUDA kernel launch-related stuff */
289 dim3 dim_block, dim_grid;
290 nbnxn_cu_kfunc_ptr_t nb_kernel = NULL; /* fn pointer to the nonbonded kernel */
292 cu_atomdata_t *adat = cu_nb->atdat;
293 cu_nbparam_t *nbp = cu_nb->nbparam;
294 cu_plist_t *plist = cu_nb->plist[iloc];
295 cu_timers_t *t = cu_nb->timers;
296 cudaStream_t stream = cu_nb->stream[iloc];
298 bool bCalcEner = flags & GMX_FORCE_VIRIAL;
299 bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
300 bool bDoTime = cu_nb->bDoTime;
302 /* turn energy calculation always on/off (for debugging/testing only) */
303 bCalcEner = (bCalcEner || always_ener) && !never_ener;
305 /* don't launch the kernel if there is no work to do */
306 if (plist->nsci == 0)
311 /* calculate the atom data index range based on locality */
315 adat_len = adat->natoms_local;
319 adat_begin = adat->natoms_local;
320 adat_len = adat->natoms - adat->natoms_local;
323 /* When we get here all misc operations issues in the local stream are done,
324 so we record that in the local stream and wait for it in the nonlocal one. */
325 if (cu_nb->bUseTwoStreams)
327 if (iloc == eintLocal)
329 stat = cudaEventRecord(cu_nb->misc_ops_done, stream);
330 CU_RET_ERR(stat, "cudaEventRecord on misc_ops_done failed");
334 stat = cudaStreamWaitEvent(stream, cu_nb->misc_ops_done, 0);
335 CU_RET_ERR(stat, "cudaStreamWaitEvent on misc_ops_done failed");
339 /* beginning of timed HtoD section */
342 stat = cudaEventRecord(t->start_nb_h2d[iloc], stream);
343 CU_RET_ERR(stat, "cudaEventRecord failed");
347 cu_copy_H2D_async(adat->xq + adat_begin, nbatom->x + adat_begin * 4,
348 adat_len * sizeof(*adat->xq), stream);
352 stat = cudaEventRecord(t->stop_nb_h2d[iloc], stream);
353 CU_RET_ERR(stat, "cudaEventRecord failed");
356 /* beginning of timed nonbonded calculation section */
359 stat = cudaEventRecord(t->start_nb_k[iloc], stream);
360 CU_RET_ERR(stat, "cudaEventRecord failed");
363 /* get the pointer to the kernel flavor we need to use */
364 nb_kernel = select_nbnxn_kernel(nbp->eeltype,
367 plist->bDoPrune || always_prune);
369 /* kernel launch config */
370 nblock = calc_nb_kernel_nblock(plist->nsci, cu_nb->dev_info);
371 dim_block = dim3(CL_SIZE, CL_SIZE, 1);
372 dim_grid = dim3(nblock, 1, 1);
373 shmem = calc_shmem_required();
377 fprintf(debug, "GPU launch configuration:\n\tThread block: %dx%dx%d\n\t"
378 "Grid: %dx%d\n\t#Super-clusters/clusters: %d/%d (%d)\n",
379 dim_block.x, dim_block.y, dim_block.z,
380 dim_grid.x, dim_grid.y, plist->nsci*NCL_PER_SUPERCL,
381 NCL_PER_SUPERCL, plist->na_c);
384 nb_kernel<<< dim_grid, dim_block, shmem, stream>>> (*adat, *nbp, *plist, bCalcFshift);
385 CU_LAUNCH_ERR("k_calc_nb");
389 stat = cudaEventRecord(t->stop_nb_k[iloc], stream);
390 CU_RET_ERR(stat, "cudaEventRecord failed");
394 void nbnxn_cuda_launch_cpyback(nbnxn_cuda_ptr_t cu_nb,
395 const nbnxn_atomdata_t *nbatom,
400 int adat_begin, adat_len, adat_end; /* local/nonlocal offset and length used for xq and f */
403 /* determine interaction locality from atom locality */
408 else if (NONLOCAL_A(aloc))
415 sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
416 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
420 cu_atomdata_t *adat = cu_nb->atdat;
421 cu_timers_t *t = cu_nb->timers;
422 bool bDoTime = cu_nb->bDoTime;
423 cudaStream_t stream = cu_nb->stream[iloc];
425 bool bCalcEner = flags & GMX_FORCE_VIRIAL;
426 bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
428 /* don't launch copy-back if there was no work to do */
429 if (cu_nb->plist[iloc]->nsci == 0)
434 /* calculate the atom data index range based on locality */
438 adat_len = adat->natoms_local;
439 adat_end = cu_nb->atdat->natoms_local;
443 adat_begin = adat->natoms_local;
444 adat_len = adat->natoms - adat->natoms_local;
445 adat_end = cu_nb->atdat->natoms;
448 /* beginning of timed D2H section */
451 stat = cudaEventRecord(t->start_nb_d2h[iloc], stream);
452 CU_RET_ERR(stat, "cudaEventRecord failed");
455 if (!cu_nb->bUseStreamSync)
457 /* For safety reasons set a few (5%) forces to NaN. This way even if the
458 polling "hack" fails with some future NVIDIA driver we'll get a crash. */
459 for (int i = adat_begin; i < 3*adat_end + 2; i += adat_len/20)
462 nbatom->out[0].f[i] = NAN;
465 if (numeric_limits<float>::has_quiet_NaN)
467 nbatom->out[0].f[i] = numeric_limits<float>::quiet_NaN();
472 nbatom->out[0].f[i] = GMX_REAL_MAX;
477 /* Set the last four bytes of the force array to a bit pattern
478 which can't be the result of the force calculation:
479 max exponent (127) and zero mantissa. */
480 *(unsigned int*)&nbatom->out[0].f[adat_end*3 - 1] = poll_wait_pattern;
483 /* With DD the local D2H transfer can only start after the non-local
484 has been launched. */
485 if (iloc == eintLocal && cu_nb->bUseTwoStreams)
487 stat = cudaStreamWaitEvent(stream, cu_nb->nonlocal_done, 0);
488 CU_RET_ERR(stat, "cudaStreamWaitEvent on nonlocal_done failed");
492 cu_copy_D2H_async(nbatom->out[0].f + adat_begin * 3, adat->f + adat_begin,
493 (adat_len)*sizeof(*adat->f), stream);
495 /* After the non-local D2H is launched the nonlocal_done event can be
496 recorded which signals that the local D2H can proceed. This event is not
497 placed after the non-local kernel because we first need the non-local
499 if (iloc == eintNonlocal)
501 stat = cudaEventRecord(cu_nb->nonlocal_done, stream);
502 CU_RET_ERR(stat, "cudaEventRecord on nonlocal_done failed");
505 /* only transfer energies in the local stream */
511 cu_copy_D2H_async(cu_nb->nbst.fshift, adat->fshift,
512 SHIFTS * sizeof(*cu_nb->nbst.fshift), stream);
518 cu_copy_D2H_async(cu_nb->nbst.e_lj, adat->e_lj,
519 sizeof(*cu_nb->nbst.e_lj), stream);
520 cu_copy_D2H_async(cu_nb->nbst.e_el, adat->e_el,
521 sizeof(*cu_nb->nbst.e_el), stream);
527 stat = cudaEventRecord(t->stop_nb_d2h[iloc], stream);
528 CU_RET_ERR(stat, "cudaEventRecord failed");
532 /* Atomic compare-exchange operation on unsigned values. It is used in
533 * polling wait for the GPU.
535 static inline bool atomic_cas(volatile unsigned int *ptr,
542 return tMPI_Atomic_cas((tMPI_Atomic_t *)ptr, oldval, newval);
544 gmx_incons("Atomic operations not available, atomic_cas() should not have been called!");
549 void nbnxn_cuda_wait_gpu(nbnxn_cuda_ptr_t cu_nb,
550 const nbnxn_atomdata_t *nbatom,
552 real *e_lj, real *e_el, rvec *fshift)
554 /* NOTE: only implemented for single-precision at this time */
556 int i, adat_end, iloc = -1;
557 volatile unsigned int *poll_word;
559 /* determine interaction locality from atom locality */
564 else if (NONLOCAL_A(aloc))
571 sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
572 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
576 cu_plist_t *plist = cu_nb->plist[iloc];
577 cu_timers_t *timers = cu_nb->timers;
578 wallclock_gpu_t *timings = cu_nb->timings;
579 nb_staging nbst = cu_nb->nbst;
581 bool bCalcEner = flags & GMX_FORCE_VIRIAL;
582 bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
584 /* turn energy calculation always on/off (for debugging/testing only) */
585 bCalcEner = (bCalcEner || always_ener) && !never_ener;
587 /* don't launch wait/update timers & counters if there was no work to do
589 NOTE: if timing with multiple GPUs (streams) becomes possible, the
590 counters could end up being inconsistent due to not being incremented
591 on some of the nodes! */
592 if (cu_nb->plist[iloc]->nsci == 0)
597 /* calculate the atom data index range based on locality */
600 adat_end = cu_nb->atdat->natoms_local;
604 adat_end = cu_nb->atdat->natoms;
607 if (cu_nb->bUseStreamSync)
609 stat = cudaStreamSynchronize(cu_nb->stream[iloc]);
610 CU_RET_ERR(stat, "cudaStreamSynchronize failed in cu_blockwait_nb");
614 /* Busy-wait until we get the signal pattern set in last byte
615 * of the l/nl float vector. This pattern corresponds to a floating
616 * point number which can't be the result of the force calculation
617 * (maximum, 127 exponent and 0 mantissa).
618 * The polling uses atomic compare-exchange.
620 poll_word = (volatile unsigned int*)&nbatom->out[0].f[adat_end*3 - 1];
621 while (atomic_cas(poll_word, poll_wait_pattern, poll_wait_pattern))
626 /* timing data accumulation */
629 /* only increase counter once (at local F wait) */
633 timings->ktime[plist->bDoPrune ? 1 : 0][bCalcEner ? 1 : 0].c += 1;
637 timings->ktime[plist->bDoPrune ? 1 : 0][bCalcEner ? 1 : 0].t +=
638 cu_event_elapsed(timers->start_nb_k[iloc], timers->stop_nb_k[iloc]);
640 /* X/q H2D and F D2H timings */
641 timings->nb_h2d_t += cu_event_elapsed(timers->start_nb_h2d[iloc],
642 timers->stop_nb_h2d[iloc]);
643 timings->nb_d2h_t += cu_event_elapsed(timers->start_nb_d2h[iloc],
644 timers->stop_nb_d2h[iloc]);
646 /* only count atdat and pair-list H2D at pair-search step */
649 /* atdat transfer timing (add only once, at local F wait) */
653 timings->pl_h2d_t += cu_event_elapsed(timers->start_atdat,
657 timings->pl_h2d_t += cu_event_elapsed(timers->start_pl_h2d[iloc],
658 timers->stop_pl_h2d[iloc]);
662 /* add up energies and shift forces (only once at local F wait) */
673 for (i = 0; i < SHIFTS; i++)
675 fshift[i][0] += nbst.fshift[i].x;
676 fshift[i][1] += nbst.fshift[i].y;
677 fshift[i][2] += nbst.fshift[i].z;
682 /* turn off pruning (doesn't matter if this is pair-search step or not) */
683 plist->bDoPrune = false;
686 /*! Return the reference to the nbfp texture. */
687 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_texref()
692 /*! Return the reference to the nbfp_comb texture. */
693 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_comb_texref()
695 return nbfp_comb_texref;
698 /*! Return the reference to the coulomb_tab. */
699 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_coulomb_tab_texref()
701 return coulomb_tab_texref;
704 /*! Set up the cache configuration for the non-bonded kernels,
706 void nbnxn_cuda_set_cacheconfig(cuda_dev_info_t *devinfo)
710 for (int i = 0; i < eelCuNR; i++)
712 for (int j = 0; j < evdwCuNR; j++)
714 if (devinfo->prop.major >= 3)
716 /* Default kernel on sm 3.x 48/16 kB Shared/L1 */
717 stat = cudaFuncSetCacheConfig(nb_kfunc_ener_prune_ptr[i][j], cudaFuncCachePreferShared);
718 stat = cudaFuncSetCacheConfig(nb_kfunc_ener_noprune_ptr[i][j], cudaFuncCachePreferShared);
719 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_prune_ptr[i][j], cudaFuncCachePreferShared);
720 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_noprune_ptr[i][j], cudaFuncCachePreferShared);
724 /* On Fermi prefer L1 gives 2% higher performance */
725 /* Default kernel on sm_2.x 16/48 kB Shared/L1 */
726 stat = cudaFuncSetCacheConfig(nb_kfunc_ener_prune_ptr[i][j], cudaFuncCachePreferL1);
727 stat = cudaFuncSetCacheConfig(nb_kfunc_ener_noprune_ptr[i][j], cudaFuncCachePreferL1);
728 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_prune_ptr[i][j], cudaFuncCachePreferL1);
729 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_noprune_ptr[i][j], cudaFuncCachePreferL1);
731 CU_RET_ERR(stat, "cudaFuncSetCacheConfig failed");