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