2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2012, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * 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"
65 /*! Texture reference for nonbonded parameters; bound to cu_nbparam_t.nbfp*/
66 texture<float, 1, cudaReadModeElementType> tex_nbfp;
68 /*! Texture reference for Ewald coulomb force table; bound to cu_nbparam_t.coulomb_tab */
69 texture<float, 1, cudaReadModeElementType> tex_coulomb_tab;
71 /* Convenience defines */
72 #define NCL_PER_SUPERCL (NBNXN_GPU_NCLUSTER_PER_SUPERCLUSTER)
73 #define CL_SIZE (NBNXN_GPU_CLUSTER_SIZE)
75 /***** The kernels come here *****/
76 #include "nbnxn_cuda_kernel_utils.cuh"
78 /* Generate all combinations of kernels through multiple inclusion:
79 F, F + E, F + prune, F + E + prune. */
81 #include "nbnxn_cuda_kernels.cuh"
82 /** Force & energy **/
84 #include "nbnxn_cuda_kernels.cuh"
87 /*** Pair-list pruning kernels ***/
90 #include "nbnxn_cuda_kernels.cuh"
91 /** Force & energy **/
93 #include "nbnxn_cuda_kernels.cuh"
97 /*! Nonbonded kernel function pointer type */
98 typedef void (*nbnxn_cu_kfunc_ptr_t)(const cu_atomdata_t,
103 /*********************************/
105 /* XXX always/never run the energy/pruning kernels -- only for benchmarking purposes */
106 static bool always_ener = (getenv("GMX_GPU_ALWAYS_ENER") != NULL);
107 static bool never_ener = (getenv("GMX_GPU_NEVER_ENER") != NULL);
108 static bool always_prune = (getenv("GMX_GPU_ALWAYS_PRUNE") != NULL);
111 /* Bit-pattern used for polling-based GPU synchronization. It is used as a float
112 * and corresponds to having the exponent set to the maximum (127 -- single
113 * precision) and the mantissa to 0.
115 static unsigned int poll_wait_pattern = (0x7FU << 23);
117 /*! Returns the number of blocks to be used for the nonbonded GPU kernel. */
118 static inline int calc_nb_kernel_nblock(int nwork_units, cuda_dev_info_t *dinfo)
124 max_grid_x_size = dinfo->prop.maxGridSize[0];
126 /* do we exceed the grid x dimension limit? */
127 if (nwork_units > max_grid_x_size)
129 gmx_fatal(FARGS, "Watch out system too large to simulate!\n"
130 "The number of nonbonded work units (=number of super-clusters) exceeds the"
131 "maximum grid size in x dimension (%d > %d)!", nwork_units, max_grid_x_size);
138 /* Constant arrays listing all kernel function pointers and enabling selection
139 of a kernel in an elegant manner. */
141 static const int nEnergyKernelTypes = 2; /* 0 - no energy, 1 - energy */
142 static const int nPruneKernelTypes = 2; /* 0 - no prune, 1 - prune */
144 /* Default kernels */
145 static const nbnxn_cu_kfunc_ptr_t
146 nb_default_kfunc_ptr[eelCuNR][nEnergyKernelTypes][nPruneKernelTypes] =
148 { { k_nbnxn_ewald, k_nbnxn_ewald_prune },
149 { k_nbnxn_ewald_ener, k_nbnxn_ewald_ener_prune } },
150 { { k_nbnxn_ewald_twin, k_nbnxn_ewald_twin_prune },
151 { k_nbnxn_ewald_twin_ener, k_nbnxn_ewald_twin_ener_prune } },
152 { { k_nbnxn_rf, k_nbnxn_rf_prune },
153 { k_nbnxn_rf_ener, k_nbnxn_rf_ener_prune } },
154 { { k_nbnxn_cutoff, k_nbnxn_cutoff_prune },
155 { k_nbnxn_cutoff_ener, k_nbnxn_cutoff_ener_prune } },
159 static const nbnxn_cu_kfunc_ptr_t
160 nb_legacy_kfunc_ptr[eelCuNR][nEnergyKernelTypes][nPruneKernelTypes] =
162 { { k_nbnxn_ewald_legacy, k_nbnxn_ewald_prune_legacy },
163 { k_nbnxn_ewald_ener_legacy, k_nbnxn_ewald_ener_prune_legacy } },
164 { { k_nbnxn_ewald_twin_legacy, k_nbnxn_ewald_twin_prune_legacy },
165 { k_nbnxn_ewald_twin_ener_legacy, k_nbnxn_ewald_twin_ener_prune_legacy } },
166 { { k_nbnxn_rf_legacy, k_nbnxn_rf_prune_legacy },
167 { k_nbnxn_rf_ener_legacy, k_nbnxn_rf_ener_prune_legacy } },
168 { { k_nbnxn_cutoff_legacy, k_nbnxn_cutoff_prune_legacy },
169 { k_nbnxn_cutoff_ener_legacy, k_nbnxn_cutoff_ener_prune_legacy } },
172 /*! Return a pointer to the kernel version to be executed at the current step. */
173 static inline nbnxn_cu_kfunc_ptr_t select_nbnxn_kernel(int kver, int eeltype,
174 bool bDoEne, bool bDoPrune)
176 assert(kver < eNbnxnCuKNR);
177 assert(eeltype < eelCuNR);
179 if (NBNXN_KVER_LEGACY(kver))
181 return nb_legacy_kfunc_ptr[eeltype][bDoEne][bDoPrune];
185 return nb_default_kfunc_ptr[eeltype][bDoEne][bDoPrune];
189 /*! Calculates the amount of shared memory required for kernel version in use. */
190 static inline int calc_shmem_required(int kver)
194 /* size of shmem (force-buffers/xq/atom type preloading) */
195 if (NBNXN_KVER_LEGACY(kver))
197 /* i-atom x+q in shared memory */
198 shmem = NCL_PER_SUPERCL * CL_SIZE * sizeof(float4);
199 /* force reduction buffers in shared memory */
200 shmem += CL_SIZE * CL_SIZE * 3 * sizeof(float);
204 /* NOTE: with the default kernel on sm3.0 we need shmem only for pre-loading */
205 /* i-atom x+q in shared memory */
206 shmem = NCL_PER_SUPERCL * CL_SIZE * sizeof(float4);
207 /* cj in shared memory, for both warps separately */
208 shmem += 2 * NBNXN_GPU_JGROUP_SIZE * sizeof(int);
210 /* i-atom types in shared memory */
211 shmem += NCL_PER_SUPERCL * CL_SIZE * sizeof(int);
213 #if __CUDA_ARCH__ < 300
214 /* force reduction buffers in shared memory */
215 shmem += CL_SIZE * CL_SIZE * 3 * sizeof(float);
222 /*! As we execute nonbonded workload in separate streams, before launching
223 the kernel we need to make sure that he following operations have completed:
224 - atomdata allocation and related H2D transfers (every nstlist step);
225 - pair list H2D transfer (every nstlist step);
226 - shift vector H2D transfer (every nstlist step);
227 - force (+shift force and energy) output clearing (every step).
229 These operations are issued in the local stream at the beginning of the step
230 and therefore always complete before the local kernel launch. The non-local
231 kernel is launched after the local on the same device/context, so this is
232 inherently scheduled after the operations in the local stream (including the
234 However, for the sake of having a future-proof implementation, we use the
235 misc_ops_done event to record the point in time when the above operations
236 are finished and synchronize with this event in the non-local stream.
238 void nbnxn_cuda_launch_kernel(nbnxn_cuda_ptr_t cu_nb,
239 const nbnxn_atomdata_t *nbatom,
244 int adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
245 /* CUDA kernel launch-related stuff */
247 dim3 dim_block, dim_grid;
248 nbnxn_cu_kfunc_ptr_t nb_kernel = NULL; /* fn pointer to the nonbonded kernel */
250 cu_atomdata_t *adat = cu_nb->atdat;
251 cu_nbparam_t *nbp = cu_nb->nbparam;
252 cu_plist_t *plist = cu_nb->plist[iloc];
253 cu_timers_t *t = cu_nb->timers;
254 cudaStream_t stream = cu_nb->stream[iloc];
256 bool bCalcEner = flags & GMX_FORCE_VIRIAL;
257 bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
258 bool bDoTime = cu_nb->bDoTime;
260 /* turn energy calculation always on/off (for debugging/testing only) */
261 bCalcEner = (bCalcEner || always_ener) && !never_ener;
263 /* don't launch the kernel if there is no work to do */
264 if (plist->nsci == 0)
269 /* calculate the atom data index range based on locality */
273 adat_len = adat->natoms_local;
277 adat_begin = adat->natoms_local;
278 adat_len = adat->natoms - adat->natoms_local;
281 /* When we get here all misc operations issues in the local stream are done,
282 so we record that in the local stream and wait for it in the nonlocal one. */
283 if (cu_nb->bUseTwoStreams)
285 if (iloc == eintLocal)
287 stat = cudaEventRecord(cu_nb->misc_ops_done, stream);
288 CU_RET_ERR(stat, "cudaEventRecord on misc_ops_done failed");
292 stat = cudaStreamWaitEvent(stream, cu_nb->misc_ops_done, 0);
293 CU_RET_ERR(stat, "cudaStreamWaitEvent on misc_ops_done failed");
297 /* beginning of timed HtoD section */
300 stat = cudaEventRecord(t->start_nb_h2d[iloc], stream);
301 CU_RET_ERR(stat, "cudaEventRecord failed");
305 cu_copy_H2D_async(adat->xq + adat_begin, nbatom->x + adat_begin * 4,
306 adat_len * sizeof(*adat->xq), stream);
310 stat = cudaEventRecord(t->stop_nb_h2d[iloc], stream);
311 CU_RET_ERR(stat, "cudaEventRecord failed");
314 /* beginning of timed nonbonded calculation section */
317 stat = cudaEventRecord(t->start_nb_k[iloc], stream);
318 CU_RET_ERR(stat, "cudaEventRecord failed");
321 /* get the pointer to the kernel flavor we need to use */
322 nb_kernel = select_nbnxn_kernel(cu_nb->kernel_ver, nbp->eeltype, bCalcEner,
323 plist->bDoPrune || always_prune);
325 /* kernel launch config */
326 nblock = calc_nb_kernel_nblock(plist->nsci, cu_nb->dev_info);
327 dim_block = dim3(CL_SIZE, CL_SIZE, 1);
328 dim_grid = dim3(nblock, 1, 1);
329 shmem = calc_shmem_required(cu_nb->kernel_ver);
333 fprintf(debug, "GPU launch configuration:\n\tThread block: %dx%dx%d\n\t"
334 "Grid: %dx%d\n\t#Super-clusters/clusters: %d/%d (%d)\n",
335 dim_block.x, dim_block.y, dim_block.z,
336 dim_grid.x, dim_grid.y, plist->nsci*NCL_PER_SUPERCL,
337 NCL_PER_SUPERCL, plist->na_c);
340 nb_kernel<<<dim_grid, dim_block, shmem, stream>>>(*adat, *nbp, *plist, bCalcFshift);
341 CU_LAUNCH_ERR("k_calc_nb");
345 stat = cudaEventRecord(t->stop_nb_k[iloc], stream);
346 CU_RET_ERR(stat, "cudaEventRecord failed");
350 void nbnxn_cuda_launch_cpyback(nbnxn_cuda_ptr_t cu_nb,
351 const nbnxn_atomdata_t *nbatom,
356 int adat_begin, adat_len, adat_end; /* local/nonlocal offset and length used for xq and f */
359 /* determine interaction locality from atom locality */
364 else if (NONLOCAL_A(aloc))
371 sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
372 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
376 cu_atomdata_t *adat = cu_nb->atdat;
377 cu_timers_t *t = cu_nb->timers;
378 bool bDoTime = cu_nb->bDoTime;
379 cudaStream_t stream = cu_nb->stream[iloc];
381 bool bCalcEner = flags & GMX_FORCE_VIRIAL;
382 bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
384 /* don't launch copy-back if there was no work to do */
385 if (cu_nb->plist[iloc]->nsci == 0)
390 /* calculate the atom data index range based on locality */
394 adat_len = adat->natoms_local;
395 adat_end = cu_nb->atdat->natoms_local;
399 adat_begin = adat->natoms_local;
400 adat_len = adat->natoms - adat->natoms_local;
401 adat_end = cu_nb->atdat->natoms;
404 /* beginning of timed D2H section */
407 stat = cudaEventRecord(t->start_nb_d2h[iloc], stream);
408 CU_RET_ERR(stat, "cudaEventRecord failed");
411 if (!cu_nb->bUseStreamSync)
413 /* For safety reasons set a few (5%) forces to NaN. This way even if the
414 polling "hack" fails with some future NVIDIA driver we'll get a crash. */
415 for (int i = adat_begin; i < 3*adat_end + 2; i += adat_len/20)
418 nbatom->out[0].f[i] = NAN;
421 if (numeric_limits<float>::has_quiet_NaN)
423 nbatom->out[0].f[i] = numeric_limits<float>::quiet_NaN();
428 nbatom->out[0].f[i] = GMX_REAL_MAX;
433 /* Set the last four bytes of the force array to a bit pattern
434 which can't be the result of the force calculation:
435 max exponent (127) and zero mantissa. */
436 *(unsigned int*)&nbatom->out[0].f[adat_end*3 - 1] = poll_wait_pattern;
439 /* With DD the local D2H transfer can only start after the non-local
440 has been launched. */
441 if (iloc == eintLocal && cu_nb->bUseTwoStreams)
443 stat = cudaStreamWaitEvent(stream, cu_nb->nonlocal_done, 0);
444 CU_RET_ERR(stat, "cudaStreamWaitEvent on nonlocal_done failed");
448 cu_copy_D2H_async(nbatom->out[0].f + adat_begin * 3, adat->f + adat_begin,
449 (adat_len)*sizeof(*adat->f), stream);
451 /* After the non-local D2H is launched the nonlocal_done event can be
452 recorded which signals that the local D2H can proceed. This event is not
453 placed after the non-local kernel because we first need the non-local
455 if (iloc == eintNonlocal)
457 stat = cudaEventRecord(cu_nb->nonlocal_done, stream);
458 CU_RET_ERR(stat, "cudaEventRecord on nonlocal_done failed");
461 /* only transfer energies in the local stream */
467 cu_copy_D2H_async(cu_nb->nbst.fshift, adat->fshift,
468 SHIFTS * sizeof(*cu_nb->nbst.fshift), stream);
474 cu_copy_D2H_async(cu_nb->nbst.e_lj, adat->e_lj,
475 sizeof(*cu_nb->nbst.e_lj), stream);
476 cu_copy_D2H_async(cu_nb->nbst.e_el, adat->e_el,
477 sizeof(*cu_nb->nbst.e_el), stream);
483 stat = cudaEventRecord(t->stop_nb_d2h[iloc], stream);
484 CU_RET_ERR(stat, "cudaEventRecord failed");
488 /* Atomic compare-exchange operation on unsigned values. It is used in
489 * polling wait for the GPU.
491 static inline bool atomic_cas(volatile unsigned int *ptr,
498 return tMPI_Atomic_cas((tMPI_Atomic_t *)ptr, oldval, newval);
500 gmx_incons("Atomic operations not available, atomic_cas() should not have been called!");
505 void nbnxn_cuda_wait_gpu(nbnxn_cuda_ptr_t cu_nb,
506 const nbnxn_atomdata_t *nbatom,
508 real *e_lj, real *e_el, rvec *fshift)
510 /* NOTE: only implemented for single-precision at this time */
512 int i, adat_end, iloc = -1;
513 volatile unsigned int *poll_word;
515 /* determine interaction locality from atom locality */
520 else if (NONLOCAL_A(aloc))
527 sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
528 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
532 cu_plist_t *plist = cu_nb->plist[iloc];
533 cu_timers_t *timers = cu_nb->timers;
534 wallclock_gpu_t *timings = cu_nb->timings;
535 nb_staging nbst = cu_nb->nbst;
537 bool bCalcEner = flags & GMX_FORCE_VIRIAL;
538 bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
540 /* turn energy calculation always on/off (for debugging/testing only) */
541 bCalcEner = (bCalcEner || always_ener) && !never_ener;
543 /* don't launch wait/update timers & counters if there was no work to do
545 NOTE: if timing with multiple GPUs (streams) becomes possible, the
546 counters could end up being inconsistent due to not being incremented
547 on some of the nodes! */
548 if (cu_nb->plist[iloc]->nsci == 0)
553 /* calculate the atom data index range based on locality */
556 adat_end = cu_nb->atdat->natoms_local;
560 adat_end = cu_nb->atdat->natoms;
563 if (cu_nb->bUseStreamSync)
565 stat = cudaStreamSynchronize(cu_nb->stream[iloc]);
566 CU_RET_ERR(stat, "cudaStreamSynchronize failed in cu_blockwait_nb");
570 /* Busy-wait until we get the signal pattern set in last byte
571 * of the l/nl float vector. This pattern corresponds to a floating
572 * point number which can't be the result of the force calculation
573 * (maximum, 127 exponent and 0 mantissa).
574 * The polling uses atomic compare-exchange.
576 poll_word = (volatile unsigned int*)&nbatom->out[0].f[adat_end*3 - 1];
577 while (atomic_cas(poll_word, poll_wait_pattern, poll_wait_pattern)) {}
580 /* timing data accumulation */
583 /* only increase counter once (at local F wait) */
587 timings->ktime[plist->bDoPrune ? 1 : 0][bCalcEner ? 1 : 0].c += 1;
591 timings->ktime[plist->bDoPrune ? 1 : 0][bCalcEner ? 1 : 0].t +=
592 cu_event_elapsed(timers->start_nb_k[iloc], timers->stop_nb_k[iloc]);
594 /* X/q H2D and F D2H timings */
595 timings->nb_h2d_t += cu_event_elapsed(timers->start_nb_h2d[iloc],
596 timers->stop_nb_h2d[iloc]);
597 timings->nb_d2h_t += cu_event_elapsed(timers->start_nb_d2h[iloc],
598 timers->stop_nb_d2h[iloc]);
600 /* only count atdat and pair-list H2D at pair-search step */
603 /* atdat transfer timing (add only once, at local F wait) */
607 timings->pl_h2d_t += cu_event_elapsed(timers->start_atdat,
611 timings->pl_h2d_t += cu_event_elapsed(timers->start_pl_h2d[iloc],
612 timers->stop_pl_h2d[iloc]);
616 /* add up energies and shift forces (only once at local F wait) */
627 for (i = 0; i < SHIFTS; i++)
629 fshift[i][0] += nbst.fshift[i].x;
630 fshift[i][1] += nbst.fshift[i].y;
631 fshift[i][2] += nbst.fshift[i].z;
636 /* turn off pruning (doesn't matter if this is pair-search step or not) */
637 plist->bDoPrune = false;
640 /*! Return the reference to the nbfp texture. */
641 const struct texture<float, 1, cudaReadModeElementType>& nbnxn_cuda_get_nbfp_texref()
646 /*! Return the reference to the coulomb_tab. */
647 const struct texture<float, 1, cudaReadModeElementType>& nbnxn_cuda_get_coulomb_tab_texref()
649 return tex_coulomb_tab;
652 /*! Set up the cache configuration for the non-bonded kernels,
654 void nbnxn_cuda_set_cacheconfig(cuda_dev_info_t *devinfo)
658 for (int i = 0; i < eelCuNR; i++)
659 for (int j = 0; j < nEnergyKernelTypes; j++)
660 for (int k = 0; k < nPruneKernelTypes; k++)
662 /* Legacy kernel 16/48 kB Shared/L1 */
663 stat = cudaFuncSetCacheConfig(nb_legacy_kfunc_ptr[i][j][k], cudaFuncCachePreferL1);
664 CU_RET_ERR(stat, "cudaFuncSetCacheConfig failed");
666 if (devinfo->prop.major >= 3)
668 /* Default kernel on sm 3.x 48/16 kB Shared/L1 */
669 stat = cudaFuncSetCacheConfig(nb_default_kfunc_ptr[i][j][k], cudaFuncCachePreferShared);
673 /* On Fermi prefer L1 gives 2% higher performance */
674 /* Default kernel on sm_2.x 16/48 kB Shared/L1 */
675 stat = cudaFuncSetCacheConfig(nb_default_kfunc_ptr[i][j][k], cudaFuncCachePreferL1);
677 CU_RET_ERR(stat, "cudaFuncSetCacheConfig failed");