7477ac7863a2de3f1fd48b72a9574319cf0638fc
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_cuda / nbnxn_cuda.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
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 /*! Pointers to the non-bonded kernels organized in 2-dim arrays by:
147  *  electrostatics and VDW type.
148  *
149  *  Note that the row- and column-order of function pointers has to match the
150  *  order of corresponding enumerated electrostatics and vdw types, resp.,
151  *  defined in nbnxn_cuda_types.h.
152  */
153
154 /*! Force-only kernel function pointers. */
155 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_noprune_ptr[eelCuNR][evdwCuNR] =
156 {
157     { nbnxn_kernel_ElecCut_VdwLJ_F_cuda,                nbnxn_kernel_ElecCut_VdwLJFsw_F_cuda,                    nbnxn_kernel_ElecCut_VdwLJPsw_F_cuda            },
158     { nbnxn_kernel_ElecRF_VdwLJ_F_cuda,                 nbnxn_kernel_ElecRF_VdwLJFsw_F_cuda,                     nbnxn_kernel_ElecRF_VdwLJPsw_F_cuda             },
159     { nbnxn_kernel_ElecEwQSTab_VdwLJ_F_cuda,            nbnxn_kernel_ElecEwQSTab_VdwLJFsw_F_cuda,                nbnxn_kernel_ElecEwQSTab_VdwLJPsw_F_cuda        },
160     { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_F_cuda,     nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_F_cuda,         nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_F_cuda },
161     { nbnxn_kernel_ElecEw_VdwLJ_F_cuda,                 nbnxn_kernel_ElecEw_VdwLJFsw_F_cuda,                     nbnxn_kernel_ElecEw_VdwLJPsw_F_cuda             },
162     { nbnxn_kernel_ElecEwTwinCut_VdwLJ_F_cuda,          nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_F_cuda,              nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_F_cuda      }
163 };
164
165 /*! Force + energy kernel function pointers. */
166 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_noprune_ptr[eelCuNR][evdwCuNR] =
167 {
168     { nbnxn_kernel_ElecCut_VdwLJ_VF_cuda,                 nbnxn_kernel_ElecCut_VdwLJFsw_VF_cuda,                 nbnxn_kernel_ElecCut_VdwLJPsw_VF_cuda            },
169     { nbnxn_kernel_ElecRF_VdwLJ_VF_cuda,                  nbnxn_kernel_ElecRF_VdwLJFsw_VF_cuda,                  nbnxn_kernel_ElecRF_VdwLJPsw_VF_cuda             },
170     { nbnxn_kernel_ElecEwQSTab_VdwLJ_VF_cuda,             nbnxn_kernel_ElecEwQSTab_VdwLJFsw_VF_cuda,             nbnxn_kernel_ElecEwQSTab_VdwLJPsw_VF_cuda        },
171     { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_VF_cuda,      nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_VF_cuda,      nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_VF_cuda },
172     { nbnxn_kernel_ElecEw_VdwLJ_VF_cuda,                  nbnxn_kernel_ElecEw_VdwLJFsw_VF_cuda,                  nbnxn_kernel_ElecEw_VdwLJPsw_VF_cuda             },
173     { nbnxn_kernel_ElecEwTwinCut_VdwLJ_VF_cuda,           nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_VF_cuda,           nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_VF_cuda      }
174 };
175
176 /*! Force + pruning kernel function pointers. */
177 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_noener_prune_ptr[eelCuNR][evdwCuNR] =
178 {
179     { nbnxn_kernel_ElecCut_VdwLJ_F_prune_cuda,            nbnxn_kernel_ElecCut_VdwLJFsw_F_prune_cuda,            nbnxn_kernel_ElecCut_VdwLJPsw_F_prune_cuda            },
180     { nbnxn_kernel_ElecRF_VdwLJ_F_prune_cuda,             nbnxn_kernel_ElecRF_VdwLJFsw_F_prune_cuda,             nbnxn_kernel_ElecRF_VdwLJPsw_F_prune_cuda             },
181     { nbnxn_kernel_ElecEwQSTab_VdwLJ_F_prune_cuda,        nbnxn_kernel_ElecEwQSTab_VdwLJFsw_F_prune_cuda,        nbnxn_kernel_ElecEwQSTab_VdwLJPsw_F_prune_cuda        },
182     { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_F_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_F_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_F_prune_cuda },
183     { nbnxn_kernel_ElecEw_VdwLJ_F_prune_cuda,             nbnxn_kernel_ElecEw_VdwLJFsw_F_prune_cuda,             nbnxn_kernel_ElecEw_VdwLJPsw_F_prune_cuda             },
184     { nbnxn_kernel_ElecEwTwinCut_VdwLJ_F_prune_cuda,      nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_F_prune_cuda,      nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_F_prune_cuda      }
185 };
186
187 /*! Force + energy + pruning kernel function pointers. */
188 static const nbnxn_cu_kfunc_ptr_t nb_kfunc_ener_prune_ptr[eelCuNR][evdwCuNR] =
189 {
190     { nbnxn_kernel_ElecCut_VdwLJ_VF_prune_cuda,            nbnxn_kernel_ElecCut_VdwLJFsw_VF_prune_cuda,           nbnxn_kernel_ElecCut_VdwLJPsw_VF_prune_cuda             },
191     { nbnxn_kernel_ElecRF_VdwLJ_VF_prune_cuda,             nbnxn_kernel_ElecRF_VdwLJFsw_VF_prune_cuda,            nbnxn_kernel_ElecRF_VdwLJPsw_VF_prune_cuda              },
192     { nbnxn_kernel_ElecEwQSTab_VdwLJ_VF_prune_cuda,        nbnxn_kernel_ElecEwQSTab_VdwLJFsw_VF_prune_cuda,        nbnxn_kernel_ElecEwQSTab_VdwLJPsw_VF_prune_cuda        },
193     { nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJ_VF_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJFsw_VF_prune_cuda, nbnxn_kernel_ElecEwQSTabTwinCut_VdwLJPsw_VF_prune_cuda },
194     { nbnxn_kernel_ElecEw_VdwLJ_VF_prune_cuda,             nbnxn_kernel_ElecEw_VdwLJFsw_VF_prune_cuda,             nbnxn_kernel_ElecEw_VdwLJPsw_VF_prune_cuda             },
195     { nbnxn_kernel_ElecEwTwinCut_VdwLJ_VF_prune_cuda,      nbnxn_kernel_ElecEwTwinCut_VdwLJFsw_VF_prune_cuda,      nbnxn_kernel_ElecEwTwinCut_VdwLJPsw_VF_prune_cuda      }
196 };
197
198 /*! Return a pointer to the kernel version to be executed at the current step. */
199 static inline nbnxn_cu_kfunc_ptr_t select_nbnxn_kernel(int  eeltype,
200                                                        int  evdwtype,
201                                                        bool bDoEne,
202                                                        bool bDoPrune)
203 {
204     nbnxn_cu_kfunc_ptr_t res;
205
206     assert(eeltype < eelCuNR);
207     assert(evdwtype < eelCuNR);
208
209     if (bDoEne)
210     {
211         if (bDoPrune)
212         {
213             res = nb_kfunc_ener_prune_ptr[eeltype][evdwtype];
214         }
215         else
216         {
217             res = nb_kfunc_ener_noprune_ptr[eeltype][evdwtype];
218         }
219     }
220     else
221     {
222         if (bDoPrune)
223         {
224             res = nb_kfunc_noener_prune_ptr[eeltype][evdwtype];
225         }
226         else
227         {
228             res = nb_kfunc_noener_noprune_ptr[eeltype][evdwtype];
229         }
230     }
231
232     return res;
233 }
234
235 /*! Calculates the amount of shared memory required by the CUDA kernel in use. */
236 static inline int calc_shmem_required()
237 {
238     int shmem;
239
240     /* size of shmem (force-buffers/xq/atom type preloading) */
241     /* NOTE: with the default kernel on sm3.0 we need shmem only for pre-loading */
242     /* i-atom x+q in shared memory */
243     shmem  = NCL_PER_SUPERCL * CL_SIZE * sizeof(float4);
244     /* cj in shared memory, for both warps separately */
245     shmem += 2 * NBNXN_GPU_JGROUP_SIZE * sizeof(int);
246 #ifdef IATYPE_SHMEM
247     /* i-atom types in shared memory */
248     shmem += NCL_PER_SUPERCL * CL_SIZE * sizeof(int);
249 #endif
250 #if __CUDA_ARCH__ < 300
251     /* force reduction buffers in shared memory */
252     shmem += CL_SIZE * CL_SIZE * 3 * sizeof(float);
253 #endif
254
255     return shmem;
256 }
257
258 /*! As we execute nonbonded workload in separate streams, before launching
259    the kernel we need to make sure that he following operations have completed:
260    - atomdata allocation and related H2D transfers (every nstlist step);
261    - pair list H2D transfer (every nstlist step);
262    - shift vector H2D transfer (every nstlist step);
263    - force (+shift force and energy) output clearing (every step).
264
265    These operations are issued in the local stream at the beginning of the step
266    and therefore always complete before the local kernel launch. The non-local
267    kernel is launched after the local on the same device/context, so this is
268    inherently scheduled after the operations in the local stream (including the
269    above "misc_ops").
270    However, for the sake of having a future-proof implementation, we use the
271    misc_ops_done event to record the point in time when the above  operations
272    are finished and synchronize with this event in the non-local stream.
273  */
274 void nbnxn_cuda_launch_kernel(nbnxn_cuda_ptr_t        cu_nb,
275                               const nbnxn_atomdata_t *nbatom,
276                               int                     flags,
277                               int                     iloc)
278 {
279     cudaError_t          stat;
280     int                  adat_begin, adat_len; /* local/nonlocal offset and length used for xq and f */
281     /* CUDA kernel launch-related stuff */
282     int                  shmem, nblock;
283     dim3                 dim_block, dim_grid;
284     nbnxn_cu_kfunc_ptr_t nb_kernel = NULL; /* fn pointer to the nonbonded kernel */
285
286     cu_atomdata_t       *adat    = cu_nb->atdat;
287     cu_nbparam_t        *nbp     = cu_nb->nbparam;
288     cu_plist_t          *plist   = cu_nb->plist[iloc];
289     cu_timers_t         *t       = cu_nb->timers;
290     cudaStream_t         stream  = cu_nb->stream[iloc];
291
292     bool                 bCalcEner   = flags & GMX_FORCE_VIRIAL;
293     bool                 bCalcFshift = flags & GMX_FORCE_VIRIAL;
294     bool                 bDoTime     = cu_nb->bDoTime;
295
296     /* turn energy calculation always on/off (for debugging/testing only) */
297     bCalcEner = (bCalcEner || always_ener) && !never_ener;
298
299     /* don't launch the kernel if there is no work to do */
300     if (plist->nsci == 0)
301     {
302         return;
303     }
304
305     /* calculate the atom data index range based on locality */
306     if (LOCAL_I(iloc))
307     {
308         adat_begin  = 0;
309         adat_len    = adat->natoms_local;
310     }
311     else
312     {
313         adat_begin  = adat->natoms_local;
314         adat_len    = adat->natoms - adat->natoms_local;
315     }
316
317     /* When we get here all misc operations issues in the local stream are done,
318        so we record that in the local stream and wait for it in the nonlocal one. */
319     if (cu_nb->bUseTwoStreams)
320     {
321         if (iloc == eintLocal)
322         {
323             stat = cudaEventRecord(cu_nb->misc_ops_done, stream);
324             CU_RET_ERR(stat, "cudaEventRecord on misc_ops_done failed");
325         }
326         else
327         {
328             stat = cudaStreamWaitEvent(stream, cu_nb->misc_ops_done, 0);
329             CU_RET_ERR(stat, "cudaStreamWaitEvent on misc_ops_done failed");
330         }
331     }
332
333     /* beginning of timed HtoD section */
334     if (bDoTime)
335     {
336         stat = cudaEventRecord(t->start_nb_h2d[iloc], stream);
337         CU_RET_ERR(stat, "cudaEventRecord failed");
338     }
339
340     /* HtoD x, q */
341     cu_copy_H2D_async(adat->xq + adat_begin, nbatom->x + adat_begin * 4,
342                       adat_len * sizeof(*adat->xq), stream);
343
344     if (bDoTime)
345     {
346         stat = cudaEventRecord(t->stop_nb_h2d[iloc], stream);
347         CU_RET_ERR(stat, "cudaEventRecord failed");
348     }
349
350     /* beginning of timed nonbonded calculation section */
351     if (bDoTime)
352     {
353         stat = cudaEventRecord(t->start_nb_k[iloc], stream);
354         CU_RET_ERR(stat, "cudaEventRecord failed");
355     }
356
357     /* get the pointer to the kernel flavor we need to use */
358     nb_kernel = select_nbnxn_kernel(nbp->eeltype,
359                                     nbp->vdwtype,
360                                     bCalcEner,
361                                     plist->bDoPrune || always_prune);
362
363     /* kernel launch config */
364     nblock    = calc_nb_kernel_nblock(plist->nsci, cu_nb->dev_info);
365     dim_block = dim3(CL_SIZE, CL_SIZE, 1);
366     dim_grid  = dim3(nblock, 1, 1);
367     shmem     = calc_shmem_required();
368
369     if (debug)
370     {
371         fprintf(debug, "GPU launch configuration:\n\tThread block: %dx%dx%d\n\t"
372                 "Grid: %dx%d\n\t#Super-clusters/clusters: %d/%d (%d)\n",
373                 dim_block.x, dim_block.y, dim_block.z,
374                 dim_grid.x, dim_grid.y, plist->nsci*NCL_PER_SUPERCL,
375                 NCL_PER_SUPERCL, plist->na_c);
376     }
377
378     nb_kernel<<< dim_grid, dim_block, shmem, stream>>> (*adat, *nbp, *plist, bCalcFshift);
379     CU_LAUNCH_ERR("k_calc_nb");
380
381     if (bDoTime)
382     {
383         stat = cudaEventRecord(t->stop_nb_k[iloc], stream);
384         CU_RET_ERR(stat, "cudaEventRecord failed");
385     }
386 }
387
388 void nbnxn_cuda_launch_cpyback(nbnxn_cuda_ptr_t        cu_nb,
389                                const nbnxn_atomdata_t *nbatom,
390                                int                     flags,
391                                int                     aloc)
392 {
393     cudaError_t stat;
394     int         adat_begin, adat_len, adat_end; /* local/nonlocal offset and length used for xq and f */
395     int         iloc = -1;
396
397     /* determine interaction locality from atom locality */
398     if (LOCAL_A(aloc))
399     {
400         iloc = eintLocal;
401     }
402     else if (NONLOCAL_A(aloc))
403     {
404         iloc = eintNonlocal;
405     }
406     else
407     {
408         char stmp[STRLEN];
409         sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
410                 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
411         gmx_incons(stmp);
412     }
413
414     cu_atomdata_t   *adat    = cu_nb->atdat;
415     cu_timers_t     *t       = cu_nb->timers;
416     bool             bDoTime = cu_nb->bDoTime;
417     cudaStream_t     stream  = cu_nb->stream[iloc];
418
419     bool             bCalcEner   = flags & GMX_FORCE_VIRIAL;
420     bool             bCalcFshift = flags & GMX_FORCE_VIRIAL;
421
422     /* don't launch copy-back if there was no work to do */
423     if (cu_nb->plist[iloc]->nsci == 0)
424     {
425         return;
426     }
427
428     /* calculate the atom data index range based on locality */
429     if (LOCAL_A(aloc))
430     {
431         adat_begin  = 0;
432         adat_len    = adat->natoms_local;
433         adat_end    = cu_nb->atdat->natoms_local;
434     }
435     else
436     {
437         adat_begin  = adat->natoms_local;
438         adat_len    = adat->natoms - adat->natoms_local;
439         adat_end    = cu_nb->atdat->natoms;
440     }
441
442     /* beginning of timed D2H section */
443     if (bDoTime)
444     {
445         stat = cudaEventRecord(t->start_nb_d2h[iloc], stream);
446         CU_RET_ERR(stat, "cudaEventRecord failed");
447     }
448
449     if (!cu_nb->bUseStreamSync)
450     {
451         /* For safety reasons set a few (5%) forces to NaN. This way even if the
452            polling "hack" fails with some future NVIDIA driver we'll get a crash. */
453         for (int i = adat_begin; i < 3*adat_end + 2; i += adat_len/20)
454         {
455 #ifdef NAN
456             nbatom->out[0].f[i] = NAN;
457 #else
458 #  ifdef _MSVC
459             if (numeric_limits<float>::has_quiet_NaN)
460             {
461                 nbatom->out[0].f[i] = numeric_limits<float>::quiet_NaN();
462             }
463             else
464 #  endif
465             {
466                 nbatom->out[0].f[i] = GMX_REAL_MAX;
467             }
468 #endif
469         }
470
471         /* Set the last four bytes of the force array to a bit pattern
472            which can't be the result of the force calculation:
473            max exponent (127) and zero mantissa. */
474         *(unsigned int*)&nbatom->out[0].f[adat_end*3 - 1] = poll_wait_pattern;
475     }
476
477     /* With DD the local D2H transfer can only start after the non-local
478        has been launched. */
479     if (iloc == eintLocal && cu_nb->bUseTwoStreams)
480     {
481         stat = cudaStreamWaitEvent(stream, cu_nb->nonlocal_done, 0);
482         CU_RET_ERR(stat, "cudaStreamWaitEvent on nonlocal_done failed");
483     }
484
485     /* DtoH f */
486     cu_copy_D2H_async(nbatom->out[0].f + adat_begin * 3, adat->f + adat_begin,
487                       (adat_len)*sizeof(*adat->f), stream);
488
489     /* After the non-local D2H is launched the nonlocal_done event can be
490        recorded which signals that the local D2H can proceed. This event is not
491        placed after the non-local kernel because we first need the non-local
492        data back first. */
493     if (iloc == eintNonlocal)
494     {
495         stat = cudaEventRecord(cu_nb->nonlocal_done, stream);
496         CU_RET_ERR(stat, "cudaEventRecord on nonlocal_done failed");
497     }
498
499     /* only transfer energies in the local stream */
500     if (LOCAL_I(iloc))
501     {
502         /* DtoH fshift */
503         if (bCalcFshift)
504         {
505             cu_copy_D2H_async(cu_nb->nbst.fshift, adat->fshift,
506                               SHIFTS * sizeof(*cu_nb->nbst.fshift), stream);
507         }
508
509         /* DtoH energies */
510         if (bCalcEner)
511         {
512             cu_copy_D2H_async(cu_nb->nbst.e_lj, adat->e_lj,
513                               sizeof(*cu_nb->nbst.e_lj), stream);
514             cu_copy_D2H_async(cu_nb->nbst.e_el, adat->e_el,
515                               sizeof(*cu_nb->nbst.e_el), stream);
516         }
517     }
518
519     if (bDoTime)
520     {
521         stat = cudaEventRecord(t->stop_nb_d2h[iloc], stream);
522         CU_RET_ERR(stat, "cudaEventRecord failed");
523     }
524 }
525
526 /* Atomic compare-exchange operation on unsigned values. It is used in
527  * polling wait for the GPU.
528  */
529 static inline bool atomic_cas(volatile unsigned int *ptr,
530                               unsigned int           oldval,
531                               unsigned int           newval)
532 {
533     assert(ptr);
534
535 #ifdef TMPI_ATOMICS
536     return tMPI_Atomic_cas((tMPI_Atomic_t *)ptr, oldval, newval);
537 #else
538     gmx_incons("Atomic operations not available, atomic_cas() should not have been called!");
539     return true;
540 #endif
541 }
542
543 void nbnxn_cuda_wait_gpu(nbnxn_cuda_ptr_t cu_nb,
544                          const nbnxn_atomdata_t *nbatom,
545                          int flags, int aloc,
546                          real *e_lj, real *e_el, rvec *fshift)
547 {
548     /* NOTE:  only implemented for single-precision at this time */
549     cudaError_t            stat;
550     int                    i, adat_end, iloc = -1;
551     volatile unsigned int *poll_word;
552
553     /* determine interaction locality from atom locality */
554     if (LOCAL_A(aloc))
555     {
556         iloc = eintLocal;
557     }
558     else if (NONLOCAL_A(aloc))
559     {
560         iloc = eintNonlocal;
561     }
562     else
563     {
564         char stmp[STRLEN];
565         sprintf(stmp, "Invalid atom locality passed (%d); valid here is only "
566                 "local (%d) or nonlocal (%d)", aloc, eatLocal, eatNonlocal);
567         gmx_incons(stmp);
568     }
569
570     cu_plist_t      *plist    = cu_nb->plist[iloc];
571     cu_timers_t     *timers   = cu_nb->timers;
572     wallclock_gpu_t *timings  = cu_nb->timings;
573     nb_staging       nbst     = cu_nb->nbst;
574
575     bool             bCalcEner   = flags & GMX_FORCE_VIRIAL;
576     bool             bCalcFshift = flags & GMX_FORCE_VIRIAL;
577
578     /* turn energy calculation always on/off (for debugging/testing only) */
579     bCalcEner = (bCalcEner || always_ener) && !never_ener;
580
581     /* don't launch wait/update timers & counters if there was no work to do
582
583        NOTE: if timing with multiple GPUs (streams) becomes possible, the
584        counters could end up being inconsistent due to not being incremented
585        on some of the nodes! */
586     if (cu_nb->plist[iloc]->nsci == 0)
587     {
588         return;
589     }
590
591     /* calculate the atom data index range based on locality */
592     if (LOCAL_A(aloc))
593     {
594         adat_end = cu_nb->atdat->natoms_local;
595     }
596     else
597     {
598         adat_end = cu_nb->atdat->natoms;
599     }
600
601     if (cu_nb->bUseStreamSync)
602     {
603         stat = cudaStreamSynchronize(cu_nb->stream[iloc]);
604         CU_RET_ERR(stat, "cudaStreamSynchronize failed in cu_blockwait_nb");
605     }
606     else
607     {
608         /* Busy-wait until we get the signal pattern set in last byte
609          * of the l/nl float vector. This pattern corresponds to a floating
610          * point number which can't be the result of the force calculation
611          * (maximum, 127 exponent and 0 mantissa).
612          * The polling uses atomic compare-exchange.
613          */
614         poll_word = (volatile unsigned int*)&nbatom->out[0].f[adat_end*3 - 1];
615         while (atomic_cas(poll_word, poll_wait_pattern, poll_wait_pattern))
616         {
617         }
618     }
619
620     /* timing data accumulation */
621     if (cu_nb->bDoTime)
622     {
623         /* only increase counter once (at local F wait) */
624         if (LOCAL_I(iloc))
625         {
626             timings->nb_c++;
627             timings->ktime[plist->bDoPrune ? 1 : 0][bCalcEner ? 1 : 0].c += 1;
628         }
629
630         /* kernel timings */
631         timings->ktime[plist->bDoPrune ? 1 : 0][bCalcEner ? 1 : 0].t +=
632             cu_event_elapsed(timers->start_nb_k[iloc], timers->stop_nb_k[iloc]);
633
634         /* X/q H2D and F D2H timings */
635         timings->nb_h2d_t += cu_event_elapsed(timers->start_nb_h2d[iloc],
636                                               timers->stop_nb_h2d[iloc]);
637         timings->nb_d2h_t += cu_event_elapsed(timers->start_nb_d2h[iloc],
638                                               timers->stop_nb_d2h[iloc]);
639
640         /* only count atdat and pair-list H2D at pair-search step */
641         if (plist->bDoPrune)
642         {
643             /* atdat transfer timing (add only once, at local F wait) */
644             if (LOCAL_A(aloc))
645             {
646                 timings->pl_h2d_c++;
647                 timings->pl_h2d_t += cu_event_elapsed(timers->start_atdat,
648                                                       timers->stop_atdat);
649             }
650
651             timings->pl_h2d_t += cu_event_elapsed(timers->start_pl_h2d[iloc],
652                                                   timers->stop_pl_h2d[iloc]);
653         }
654     }
655
656     /* add up energies and shift forces (only once at local F wait) */
657     if (LOCAL_I(iloc))
658     {
659         if (bCalcEner)
660         {
661             *e_lj += *nbst.e_lj;
662             *e_el += *nbst.e_el;
663         }
664
665         if (bCalcFshift)
666         {
667             for (i = 0; i < SHIFTS; i++)
668             {
669                 fshift[i][0] += nbst.fshift[i].x;
670                 fshift[i][1] += nbst.fshift[i].y;
671                 fshift[i][2] += nbst.fshift[i].z;
672             }
673         }
674     }
675
676     /* turn off pruning (doesn't matter if this is pair-search step or not) */
677     plist->bDoPrune = false;
678 }
679
680 /*! Return the reference to the nbfp texture. */
681 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_texref()
682 {
683     return nbfp_texref;
684 }
685
686 /*! Return the reference to the coulomb_tab. */
687 const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_coulomb_tab_texref()
688 {
689     return coulomb_tab_texref;
690 }
691
692 /*! Set up the cache configuration for the non-bonded kernels,
693  */
694 void nbnxn_cuda_set_cacheconfig(cuda_dev_info_t *devinfo)
695 {
696     cudaError_t stat;
697
698     for (int i = 0; i < eelCuNR; i++)
699     {
700         for (int j = 0; j < evdwCuNR; j++)
701         {
702             if (devinfo->prop.major >= 3)
703             {
704                 /* Default kernel on sm 3.x 48/16 kB Shared/L1 */
705                 stat = cudaFuncSetCacheConfig(nb_kfunc_ener_prune_ptr[i][j], cudaFuncCachePreferShared);
706                 stat = cudaFuncSetCacheConfig(nb_kfunc_ener_noprune_ptr[i][j], cudaFuncCachePreferShared);
707                 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_prune_ptr[i][j], cudaFuncCachePreferShared);
708                 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_noprune_ptr[i][j], cudaFuncCachePreferShared);
709             }
710             else
711             {
712                 /* On Fermi prefer L1 gives 2% higher performance */
713                 /* Default kernel on sm_2.x 16/48 kB Shared/L1 */
714                 stat = cudaFuncSetCacheConfig(nb_kfunc_ener_prune_ptr[i][j], cudaFuncCachePreferL1);
715                 stat = cudaFuncSetCacheConfig(nb_kfunc_ener_noprune_ptr[i][j], cudaFuncCachePreferL1);
716                 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_prune_ptr[i][j], cudaFuncCachePreferL1);
717                 stat = cudaFuncSetCacheConfig(nb_kfunc_noener_noprune_ptr[i][j], cudaFuncCachePreferL1);
718             }
719             CU_RET_ERR(stat, "cudaFuncSetCacheConfig failed");
720         }
721     }
722 }