Fix CUDA inter-stream synchronization issue
[alexxy/gromacs.git] / src / mdlib / nbnxn_cuda / nbnxn_cuda.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38
39 #include <stdlib.h>
40 #include <assert.h>
41
42 #if defined(_MSVC)
43 #include <limits>
44 #endif
45
46 #include <cuda.h>
47
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"
54
55 #ifdef TMPI_ATOMICS
56 #include "thread_mpi/atomic.h"
57 #endif
58
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"
63
64 #if defined TEXOBJ_SUPPORTED && __CUDA_ARCH__ >= 300
65 #define USE_TEXOBJ
66 #endif
67
68 /*! Texture reference for nonbonded parameters; bound to cu_nbparam_t.nbfp*/
69 texture<float, 1, cudaReadModeElementType> nbfp_texref;
70
71 /*! Texture reference for Ewald coulomb force table; bound to cu_nbparam_t.coulomb_tab */
72 texture<float, 1, cudaReadModeElementType> coulomb_tab_texref;
73
74 /* Convenience defines */
75 #define NCL_PER_SUPERCL         (NBNXN_GPU_NCLUSTER_PER_SUPERCLUSTER)
76 #define CL_SIZE                 (NBNXN_GPU_CLUSTER_SIZE)
77
78 /***** The kernels come here *****/
79 #include "nbnxn_cuda_kernel_utils.cuh"
80
81 /* Top-level kernel generation: will generate through multiple inclusion the
82  * following flavors for all kernels:
83  * - force-only output;
84  * - force and energy output;
85  * - force-only with pair list pruning;
86  * - force and energy output with pair list pruning.
87  */
88 /** Force only **/
89 #include "nbnxn_cuda_kernels.cuh"
90 /** Force & energy **/
91 #define CALC_ENERGIES
92 #include "nbnxn_cuda_kernels.cuh"
93 #undef CALC_ENERGIES
94
95 /*** Pair-list pruning kernels ***/
96 /** Force only **/
97 #define PRUNE_NBL
98 #include "nbnxn_cuda_kernels.cuh"
99 /** Force & energy **/
100 #define CALC_ENERGIES
101 #include "nbnxn_cuda_kernels.cuh"
102 #undef CALC_ENERGIES
103 #undef PRUNE_NBL
104
105 /*! Nonbonded kernel function pointer type */
106 typedef void (*nbnxn_cu_kfunc_ptr_t)(const cu_atomdata_t,
107                                      const cu_nbparam_t,
108                                      const cu_plist_t,
109                                      bool);
110
111 /*********************************/
112
113 /* XXX always/never run the energy/pruning kernels -- only for benchmarking purposes */
114 static bool always_ener  = (getenv("GMX_GPU_ALWAYS_ENER") != NULL);
115 static bool never_ener   = (getenv("GMX_GPU_NEVER_ENER") != NULL);
116 static bool always_prune = (getenv("GMX_GPU_ALWAYS_PRUNE") != NULL);
117
118
119 /* Bit-pattern used for polling-based GPU synchronization. It is used as a float
120  * and corresponds to having the exponent set to the maximum (127 -- single
121  * precision) and the mantissa to 0.
122  */
123 static unsigned int poll_wait_pattern = (0x7FU << 23);
124
125 /*! Returns the number of blocks to be used for the nonbonded GPU kernel. */
126 static inline int calc_nb_kernel_nblock(int nwork_units, cuda_dev_info_t *dinfo)
127 {
128     int max_grid_x_size;
129
130     assert(dinfo);
131
132     max_grid_x_size = dinfo->prop.maxGridSize[0];
133
134     /* do we exceed the grid x dimension limit? */
135     if (nwork_units > max_grid_x_size)
136     {
137         gmx_fatal(FARGS, "Watch out, the input system is too large to simulate!\n"
138                   "The number of nonbonded work units (=number of super-clusters) exceeds the"
139                   "maximum grid size in x dimension (%d > %d)!", nwork_units, max_grid_x_size);
140     }
141
142     return nwork_units;
143 }
144
145
146 /* Constant arrays listing all kernel function pointers and enabling selection
147    of a kernel in an elegant manner. */
148
149 static const int nEnergyKernelTypes = 2; /* 0 - no energy, 1 - energy */
150 static const int nPruneKernelTypes  = 2; /* 0 - no prune, 1 - prune */
151
152 /*! Pointers to the default kernels organized in a 3 dim array by:
153  *  electrostatics type, energy calculation on/off, and pruning on/off.
154  *
155  *  Note that the order of electrostatics (1st dimension) has to match the
156  *  order of corresponding enumerated types defined in nbnxn_cuda_types.h.
157  */
158 static const nbnxn_cu_kfunc_ptr_t
159 nb_default_kfunc_ptr[eelCuNR][nEnergyKernelTypes][nPruneKernelTypes] =
160 {
161     { { k_nbnxn_cutoff,                     k_nbnxn_cutoff_prune },
162       { k_nbnxn_cutoff_ener,                k_nbnxn_cutoff_ener_prune } },
163     { { k_nbnxn_rf,                         k_nbnxn_rf_prune },
164       { k_nbnxn_rf_ener,                    k_nbnxn_rf_ener_prune } },
165     { { k_nbnxn_ewald_tab,                  k_nbnxn_ewald_tab_prune },
166       { k_nbnxn_ewald_tab_ener,             k_nbnxn_ewald_tab_ener_prune } },
167     { { k_nbnxn_ewald_tab_twin,             k_nbnxn_ewald_tab_twin_prune },
168       { k_nbnxn_ewald_tab_twin_ener,        k_nbnxn_ewald_twin_ener_prune } },
169     { { k_nbnxn_ewald,                      k_nbnxn_ewald_prune },
170       { k_nbnxn_ewald_ener,                 k_nbnxn_ewald_ener_prune } },
171     { { k_nbnxn_ewald_twin,                 k_nbnxn_ewald_twin_prune },
172       { k_nbnxn_ewald_twin_ener,            k_nbnxn_ewald_twin_ener_prune } },
173 };
174
175 /*! Pointers to the legacy kernels organized in a 3 dim array by:
176  *  electrostatics type, energy calculation on/off, and pruning on/off.
177  *
178  *  Note that the order of electrostatics (1st dimension) has to match the
179  *  order of corresponding enumerated types defined in nbnxn_cuda_types.h.
180  */
181 static const nbnxn_cu_kfunc_ptr_t
182 nb_legacy_kfunc_ptr[eelCuNR][nEnergyKernelTypes][nPruneKernelTypes] =
183 {
184     { { k_nbnxn_cutoff_legacy,              k_nbnxn_cutoff_prune_legacy },
185       { k_nbnxn_cutoff_ener_legacy,         k_nbnxn_cutoff_ener_prune_legacy } },
186     { { k_nbnxn_rf_legacy,                  k_nbnxn_rf_prune_legacy },
187       { k_nbnxn_rf_ener_legacy,             k_nbnxn_rf_ener_prune_legacy } },
188     { { k_nbnxn_ewald_tab_legacy,           k_nbnxn_ewald_tab_prune_legacy },
189       { k_nbnxn_ewald_tab_ener_legacy,      k_nbnxn_ewald_tab_ener_prune_legacy } },
190     { { k_nbnxn_ewald_tab_twin_legacy,      k_nbnxn_ewald_tab_twin_prune_legacy },
191       { k_nbnxn_ewald_tab_twin_ener_legacy, k_nbnxn_ewald_tab_twin_ener_prune_legacy } },
192 };
193
194 /*! Return a pointer to the kernel version to be executed at the current step. */
195 static inline nbnxn_cu_kfunc_ptr_t select_nbnxn_kernel(int kver, int eeltype,
196                                                        bool bDoEne, bool bDoPrune)
197 {
198     assert(kver < eNbnxnCuKNR);
199     assert(eeltype < eelCuNR);
200
201     if (NBNXN_KVER_LEGACY(kver))
202     {
203         /* no analytical Ewald with legacy kernels */
204         assert(eeltype <= eelCuEWALD_TAB_TWIN);
205
206         return nb_legacy_kfunc_ptr[eeltype][bDoEne][bDoPrune];
207     }
208     else
209     {
210         return nb_default_kfunc_ptr[eeltype][bDoEne][bDoPrune];
211     }
212 }
213
214 /*! Calculates the amount of shared memory required for kernel version in use. */
215 static inline int calc_shmem_required(int kver)
216 {
217     int shmem;
218
219     /* size of shmem (force-buffers/xq/atom type preloading) */
220     if (NBNXN_KVER_LEGACY(kver))
221     {
222         /* i-atom x+q in shared memory */
223         shmem =  NCL_PER_SUPERCL * CL_SIZE * sizeof(float4);
224         /* force reduction buffers in shared memory */
225         shmem += CL_SIZE * CL_SIZE * 3 * sizeof(float);
226     }
227     else
228     {
229         /* NOTE: with the default kernel on sm3.0 we need shmem only for pre-loading */
230         /* i-atom x+q in shared memory */
231         shmem  = NCL_PER_SUPERCL * CL_SIZE * sizeof(float4);
232         /* cj in shared memory, for both warps separately */
233         shmem += 2 * NBNXN_GPU_JGROUP_SIZE * sizeof(int);
234 #ifdef IATYPE_SHMEM
235         /* i-atom types in shared memory */
236         shmem += NCL_PER_SUPERCL * CL_SIZE * sizeof(int);
237 #endif
238 #if __CUDA_ARCH__ < 300
239         /* force reduction buffers in shared memory */
240         shmem += CL_SIZE * CL_SIZE * 3 * sizeof(float);
241 #endif
242     }
243
244     return shmem;
245 }
246
247 /*! As we execute nonbonded workload in separate streams, before launching 
248    the kernel we need to make sure that he following operations have completed:
249    - atomdata allocation and related H2D transfers (every nstlist step);
250    - pair list H2D transfer (every nstlist step);
251    - shift vector H2D transfer (every nstlist step);
252    - force (+shift force and energy) output clearing (every step).
253
254    These operations are issued in the local stream at the beginning of the step
255    and therefore always complete before the local kernel launch. The non-local
256    kernel is launched after the local on the same device/context hence it is
257    inherently scheduled after the operations in the local stream (including the
258    above "misc_ops") on pre-GK110 devices with single hardware queue, but on later
259    devices with multiple hardware queues the dependency needs to be enforced.
260    We use the misc_ops_and_local_H2D_done event to record the point where
261    the local x+q H2D (and all preceding) tasks are complete and synchronize
262    with this event in the non-local stream before launching the non-bonded kernel.
263 */
264 void nbnxn_cuda_launch_kernel(nbnxn_cuda_ptr_t cu_nb,
265                               const nbnxn_atomdata_t *nbatom,
266                               int flags,
267                               int iloc)
268 {
269     cudaError_t stat;
270     int adat_begin, adat_len;  /* local/nonlocal offset and length used for xq and f */
271     /* CUDA kernel launch-related stuff */
272     int  shmem, nblock;
273     dim3 dim_block, dim_grid;
274     nbnxn_cu_kfunc_ptr_t nb_kernel = NULL; /* fn pointer to the nonbonded kernel */
275
276     cu_atomdata_t   *adat   = cu_nb->atdat;
277     cu_nbparam_t    *nbp    = cu_nb->nbparam;
278     cu_plist_t      *plist  = cu_nb->plist[iloc];
279     cu_timers_t     *t      = cu_nb->timers;
280     cudaStream_t    stream  = cu_nb->stream[iloc];
281
282     bool bCalcEner   = flags & GMX_FORCE_VIRIAL;
283     bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
284     bool bDoTime     = cu_nb->bDoTime;
285
286     /* turn energy calculation always on/off (for debugging/testing only) */
287     bCalcEner = (bCalcEner || always_ener) && !never_ener;
288
289     /* Don't launch the non-local kernel if there is no work to do.
290        Doing the same for the local kernel is more complicated, since the
291        local part of the force array also depends on the non-local kernel.
292        So to avoid complicating the code and to reduce the risk of bugs,
293        we always call the local kernel, the local x+q copy and later (not in
294        this function) the stream wait, local f copyback and the f buffer
295        clearing. All these operations, except for the local interaction kernel,
296        are needed for the non-local interactions. The skip of the local kernel
297        call is taken care of later in this function. */
298     if (iloc == eintNonlocal && plist->nsci == 0)
299     {
300         return;
301     }
302
303     /* calculate the atom data index range based on locality */
304     if (LOCAL_I(iloc))
305     {
306         adat_begin  = 0;
307         adat_len    = adat->natoms_local;
308     }
309     else
310     {
311         adat_begin  = adat->natoms_local;
312         adat_len    = adat->natoms - adat->natoms_local;
313     }
314
315     /* beginning of timed HtoD section */
316     if (bDoTime)
317     {
318         stat = cudaEventRecord(t->start_nb_h2d[iloc], stream);
319         CU_RET_ERR(stat, "cudaEventRecord failed");
320     }
321
322     /* HtoD x, q */
323     cu_copy_H2D_async(adat->xq + adat_begin, nbatom->x + adat_begin * 4,
324                       adat_len * sizeof(*adat->xq), stream); 
325
326     /* When we get here all misc operations issues in the local stream as well as
327        the local xq H2D are done,
328        so we record that in the local stream and wait for it in the nonlocal one. */
329     if (cu_nb->bUseTwoStreams)
330     {
331         if (iloc == eintLocal)
332         {
333             stat = cudaEventRecord(cu_nb->misc_ops_and_local_H2D_done, stream);
334             CU_RET_ERR(stat, "cudaEventRecord on misc_ops_and_local_H2D_done failed");
335         }
336         else
337         {
338             stat = cudaStreamWaitEvent(stream, cu_nb->misc_ops_and_local_H2D_done, 0);
339             CU_RET_ERR(stat, "cudaStreamWaitEvent on misc_ops_and_local_H2D_done failed");
340         }
341     }
342
343     if (bDoTime)
344     {
345         stat = cudaEventRecord(t->stop_nb_h2d[iloc], stream);
346         CU_RET_ERR(stat, "cudaEventRecord failed");
347     }
348
349     if (plist->nsci == 0)
350     {
351         /* Don't launch an empty local kernel (not allowed with CUDA) */
352         return;
353     }
354
355     /* beginning of timed nonbonded calculation section */
356     if (bDoTime)
357     {
358         stat = cudaEventRecord(t->start_nb_k[iloc], stream);
359         CU_RET_ERR(stat, "cudaEventRecord failed");
360     }
361
362     /* get the pointer to the kernel flavor we need to use */
363     nb_kernel = select_nbnxn_kernel(cu_nb->kernel_ver, nbp->eeltype, bCalcEner,
364                                     plist->bDoPrune || always_prune);
365
366     /* kernel launch config */
367     nblock    = calc_nb_kernel_nblock(plist->nsci, cu_nb->dev_info);
368     dim_block = dim3(CL_SIZE, CL_SIZE, 1);
369     dim_grid  = dim3(nblock, 1, 1);
370     shmem     = calc_shmem_required(cu_nb->kernel_ver);
371
372     if (debug)
373     {
374         fprintf(debug, "GPU launch configuration:\n\tThread block: %dx%dx%d\n\t"
375                 "Grid: %dx%d\n\t#Super-clusters/clusters: %d/%d (%d)\n",
376                 dim_block.x, dim_block.y, dim_block.z,
377                 dim_grid.x, dim_grid.y, plist->nsci*NCL_PER_SUPERCL,
378                 NCL_PER_SUPERCL, plist->na_c);
379     }
380
381     nb_kernel<<<dim_grid, dim_block, shmem, stream>>>(*adat, *nbp, *plist, bCalcFshift);
382     CU_LAUNCH_ERR("k_calc_nb");
383
384     if (bDoTime)
385     {
386         stat = cudaEventRecord(t->stop_nb_k[iloc], stream);
387         CU_RET_ERR(stat, "cudaEventRecord failed");
388     }
389 }
390
391 void nbnxn_cuda_launch_cpyback(nbnxn_cuda_ptr_t cu_nb,
392                                const nbnxn_atomdata_t *nbatom,
393                                int flags,
394                                int aloc)
395 {
396     cudaError_t stat;
397     int adat_begin, adat_len, adat_end;  /* local/nonlocal offset and length used for xq and f */
398     int iloc = -1;
399
400     /* determine interaction locality from atom locality */
401     if (LOCAL_A(aloc))
402     {
403         iloc = eintLocal;
404     }
405     else if (NONLOCAL_A(aloc))
406     {
407         iloc = eintNonlocal;
408     }
409     else
410     {
411         char stmp[STRLEN];
412         sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
413                 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
414         gmx_incons(stmp);
415     }
416
417     cu_atomdata_t   *adat   = cu_nb->atdat;
418     cu_timers_t     *t      = cu_nb->timers;
419     bool            bDoTime = cu_nb->bDoTime;
420     cudaStream_t    stream  = cu_nb->stream[iloc];
421
422     bool bCalcEner   = flags & GMX_FORCE_VIRIAL;
423     bool bCalcFshift = flags & GMX_FORCE_VIRIAL;
424
425     /* don't launch non-local copy-back if there was no non-local work to do */
426     if (iloc == eintNonlocal && cu_nb->plist[iloc]->nsci == 0)
427     {
428         return;
429     }
430
431     /* calculate the atom data index range based on locality */
432     if (LOCAL_A(aloc))
433     {
434         adat_begin  = 0;
435         adat_len    = adat->natoms_local;
436         adat_end    = cu_nb->atdat->natoms_local;
437     }
438     else
439     {
440         adat_begin  = adat->natoms_local;
441         adat_len    = adat->natoms - adat->natoms_local;
442         adat_end    = cu_nb->atdat->natoms;
443     }
444
445     /* beginning of timed D2H section */
446     if (bDoTime)
447     {
448         stat = cudaEventRecord(t->start_nb_d2h[iloc], stream);
449         CU_RET_ERR(stat, "cudaEventRecord failed");
450     }
451
452     if (!cu_nb->bUseStreamSync)
453     {
454         /* For safety reasons set a few (5%) forces to NaN. This way even if the
455            polling "hack" fails with some future NVIDIA driver we'll get a crash. */
456         for (int i = adat_begin; i < 3*adat_end + 2; i += adat_len/20)
457         {
458 #ifdef NAN
459             nbatom->out[0].f[i] = NAN;
460 #else
461 #  ifdef _MSVC
462             if (numeric_limits<float>::has_quiet_NaN)
463             {
464                 nbatom->out[0].f[i] = numeric_limits<float>::quiet_NaN();
465             }
466             else
467 #  endif
468             {
469                 nbatom->out[0].f[i] = GMX_REAL_MAX;
470             }
471 #endif
472         }
473
474         /* Set the last four bytes of the force array to a bit pattern
475            which can't be the result of the force calculation:
476            max exponent (127) and zero mantissa. */
477         *(unsigned int*)&nbatom->out[0].f[adat_end*3 - 1] = poll_wait_pattern;
478     }
479
480     /* With DD the local D2H transfer can only start after the non-local 
481        kernel has finished. */
482     if (iloc == eintLocal && cu_nb->bUseTwoStreams)
483     {
484         stat = cudaStreamWaitEvent(stream, cu_nb->nonlocal_done, 0);
485         CU_RET_ERR(stat, "cudaStreamWaitEvent on nonlocal_done failed");
486     }
487
488     /* DtoH f */
489     cu_copy_D2H_async(nbatom->out[0].f + adat_begin * 3, adat->f + adat_begin, 
490                       (adat_len)*sizeof(*adat->f), stream);
491
492     /* After the non-local D2H is launched the nonlocal_done event can be
493        recorded which signals that the local D2H can proceed. This event is not
494        placed after the non-local kernel because we want the non-local data
495        back first. */
496     if (iloc == eintNonlocal)
497     {
498         stat = cudaEventRecord(cu_nb->nonlocal_done, stream);
499         CU_RET_ERR(stat, "cudaEventRecord on nonlocal_done failed");
500     }
501
502     /* only transfer energies in the local stream */
503     if (LOCAL_I(iloc))
504     {
505         /* DtoH fshift */
506         if (bCalcFshift)
507         {
508             cu_copy_D2H_async(cu_nb->nbst.fshift, adat->fshift,
509                               SHIFTS * sizeof(*cu_nb->nbst.fshift), stream);
510         }
511
512         /* DtoH energies */
513         if (bCalcEner)
514         {
515             cu_copy_D2H_async(cu_nb->nbst.e_lj, adat->e_lj,
516                               sizeof(*cu_nb->nbst.e_lj), stream);
517             cu_copy_D2H_async(cu_nb->nbst.e_el, adat->e_el,
518                               sizeof(*cu_nb->nbst.e_el), stream);
519         }
520     }
521
522     if (bDoTime)
523     {
524         stat = cudaEventRecord(t->stop_nb_d2h[iloc], stream);
525         CU_RET_ERR(stat, "cudaEventRecord failed");
526     }
527 }
528
529 /* Atomic compare-exchange operation on unsigned values. It is used in
530  * polling wait for the GPU.
531  */
532 static inline bool atomic_cas(volatile unsigned int *ptr,
533                               unsigned int oldval,
534                               unsigned int newval)
535 {
536     assert(ptr);
537
538 #ifdef TMPI_ATOMICS
539     return tMPI_Atomic_cas((tMPI_Atomic_t *)ptr, oldval, newval);
540 #else
541     gmx_incons("Atomic operations not available, atomic_cas() should not have been called!");
542     return true;
543 #endif
544 }
545
546 void nbnxn_cuda_wait_gpu(nbnxn_cuda_ptr_t cu_nb,
547                          const nbnxn_atomdata_t *nbatom,
548                          int flags, int aloc,
549                          real *e_lj, real *e_el, rvec *fshift)
550 {
551     /* NOTE:  only implemented for single-precision at this time */
552     cudaError_t stat;
553     int i, adat_end, iloc = -1;
554     volatile unsigned int *poll_word;
555
556     /* determine interaction locality from atom locality */
557     if (LOCAL_A(aloc))
558     {
559         iloc = eintLocal;
560     }
561     else if (NONLOCAL_A(aloc))
562     {
563         iloc = eintNonlocal;
564     }
565     else
566     {
567         char stmp[STRLEN];
568         sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
569                 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
570         gmx_incons(stmp);
571     }
572
573     cu_plist_t      *plist   = cu_nb->plist[iloc];
574     cu_timers_t     *timers  = cu_nb->timers;
575     wallclock_gpu_t *timings = cu_nb->timings;
576     nb_staging      nbst     = cu_nb->nbst;
577
578     bool    bCalcEner   = flags & GMX_FORCE_VIRIAL;
579     bool    bCalcFshift = flags & GMX_FORCE_VIRIAL;
580
581     /* turn energy calculation always on/off (for debugging/testing only) */
582     bCalcEner = (bCalcEner || always_ener) && !never_ener;
583
584     /* Launch wait/update timers & counters, unless doing the non-local phase
585        when there is not actually work to do. This is consistent with
586        nbnxn_cuda_launch_kernel.
587
588        NOTE: if timing with multiple GPUs (streams) becomes possible, the
589        counters could end up being inconsistent due to not being incremented
590        on some of the nodes! */
591     if (iloc == eintNonlocal && cu_nb->plist[iloc]->nsci == 0)
592     {
593         return;
594     }
595
596     /* calculate the atom data index range based on locality */
597     if (LOCAL_A(aloc))
598     {
599         adat_end = cu_nb->atdat->natoms_local;
600     }
601     else
602     {
603         adat_end = cu_nb->atdat->natoms;
604     }
605
606     if (cu_nb->bUseStreamSync)
607     {
608         stat = cudaStreamSynchronize(cu_nb->stream[iloc]);
609         CU_RET_ERR(stat, "cudaStreamSynchronize failed in cu_blockwait_nb");
610     }
611     else 
612     {
613         /* Busy-wait until we get the signal pattern set in last byte
614          * of the l/nl float vector. This pattern corresponds to a floating
615          * point number which can't be the result of the force calculation
616          * (maximum, 127 exponent and 0 mantissa).
617          * The polling uses atomic compare-exchange.
618          */
619         poll_word = (volatile unsigned int*)&nbatom->out[0].f[adat_end*3 - 1];
620         while (atomic_cas(poll_word, poll_wait_pattern, poll_wait_pattern)) {}
621     }
622
623     /* timing data accumulation */
624     if (cu_nb->bDoTime)
625     {
626         /* only increase counter once (at local F wait) */
627         if (LOCAL_I(iloc))
628         {
629             timings->nb_c++;
630             timings->ktime[plist->bDoPrune ? 1 : 0][bCalcEner ? 1 : 0].c += 1;
631         }
632
633         /* kernel timings */
634         timings->ktime[plist->bDoPrune ? 1 : 0][bCalcEner ? 1 : 0].t +=
635             cu_event_elapsed(timers->start_nb_k[iloc], timers->stop_nb_k[iloc]);
636
637         /* X/q H2D and F D2H timings */
638         timings->nb_h2d_t += cu_event_elapsed(timers->start_nb_h2d[iloc],
639                                                  timers->stop_nb_h2d[iloc]);
640         timings->nb_d2h_t += cu_event_elapsed(timers->start_nb_d2h[iloc],
641                                                  timers->stop_nb_d2h[iloc]);
642
643         /* only count atdat and pair-list H2D at pair-search step */
644         if (plist->bDoPrune)
645         {
646             /* atdat transfer timing (add only once, at local F wait) */
647             if (LOCAL_A(aloc))
648             {
649                 timings->pl_h2d_c++;
650                 timings->pl_h2d_t += cu_event_elapsed(timers->start_atdat,
651                                                          timers->stop_atdat);
652             }
653
654             timings->pl_h2d_t += cu_event_elapsed(timers->start_pl_h2d[iloc],
655                                                      timers->stop_pl_h2d[iloc]);
656         }
657     }
658
659     /* add up energies and shift forces (only once at local F wait) */
660     if (LOCAL_I(iloc))
661     {
662         if (bCalcEner)
663         {
664             *e_lj += *nbst.e_lj;
665             *e_el += *nbst.e_el;
666         }
667
668         if (bCalcFshift)
669         {
670             for (i = 0; i < SHIFTS; i++)
671             {
672                 fshift[i][0] += nbst.fshift[i].x;
673                 fshift[i][1] += nbst.fshift[i].y;
674                 fshift[i][2] += nbst.fshift[i].z;
675             }
676         }
677     }
678
679     /* turn off pruning (doesn't matter if this is pair-search step or not) */
680     plist->bDoPrune = false;
681 }
682
683 /*! Return the reference to the nbfp texture. */
684 const struct texture<float, 1, cudaReadModeElementType>& nbnxn_cuda_get_nbfp_texref()
685 {
686     return nbfp_texref;
687 }
688
689 /*! Return the reference to the coulomb_tab. */
690 const struct texture<float, 1, cudaReadModeElementType>& nbnxn_cuda_get_coulomb_tab_texref()
691 {
692     return coulomb_tab_texref;
693 }
694
695 /*! Set up the cache configuration for the non-bonded kernels,
696  */
697 void nbnxn_cuda_set_cacheconfig(cuda_dev_info_t *devinfo)
698 {
699     cudaError_t stat;
700
701     for (int i = 0; i < eelCuNR; i++)
702     {
703         for (int j = 0; j < nEnergyKernelTypes; j++)
704         {
705             for (int k = 0; k < nPruneKernelTypes; k++)
706             {
707                 /* Legacy kernel 16/48 kB Shared/L1
708                  * No analytical Ewald!
709                  */
710                 if (i != eelCuEWALD_ANA && i != eelCuEWALD_ANA_TWIN)
711                 {
712                     stat = cudaFuncSetCacheConfig(nb_legacy_kfunc_ptr[i][j][k], cudaFuncCachePreferL1);
713                     CU_RET_ERR(stat, "cudaFuncSetCacheConfig failed");
714                 }
715
716                 if (devinfo->prop.major >= 3)
717                 {
718                     /* Default kernel on sm 3.x 48/16 kB Shared/L1 */
719                     stat = cudaFuncSetCacheConfig(nb_default_kfunc_ptr[i][j][k], cudaFuncCachePreferShared);
720                 }
721                 else
722                 {
723                     /* On Fermi prefer L1 gives 2% higher performance */
724                     /* Default kernel on sm_2.x 16/48 kB Shared/L1 */
725                     stat = cudaFuncSetCacheConfig(nb_default_kfunc_ptr[i][j][k], cudaFuncCachePreferL1);
726                 }
727                 CU_RET_ERR(stat, "cudaFuncSetCacheConfig failed");
728             }
729         }
730     }
731 }