Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_cuda / nbnxn_cuda_kernel.cuh
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 "maths.h"
37 /* Note that floating-point constants in CUDA code should be suffixed
38  * with f (e.g. 0.5f), to stop the compiler producing intermediate
39  * code that is in double precision.
40  */
41
42 #if __CUDA_ARCH__ >= 300
43 #define REDUCE_SHUFFLE
44 /* On Kepler pre-loading i-atom types to shmem gives a few %,
45    but on Fermi it does not */
46 #define IATYPE_SHMEM
47 #endif
48
49 /*
50    Kernel launch parameters:
51     - #blocks   = #pair lists, blockId = pair list Id
52     - #threads  = CL_SIZE^2
53     - shmem     = CL_SIZE^2 * sizeof(float)
54
55     Each thread calculates an i force-component taking one pair of i-j atoms.
56  */
57 #ifdef PRUNE_NBL
58 #ifdef CALC_ENERGIES
59 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener_prune)
60 #else
61 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _prune)
62 #endif
63 #else
64 #ifdef CALC_ENERGIES
65 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener)
66 #else
67 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn)
68 #endif
69 #endif
70             (const cu_atomdata_t atdat,
71              const cu_nbparam_t nbparam,
72              const cu_plist_t plist,
73              bool bCalcFshift)
74 {
75     /* convenience variables */
76     const nbnxn_sci_t *pl_sci   = plist.sci;
77 #ifndef PRUNE_NBL
78     const
79 #endif
80     nbnxn_cj4_t *pl_cj4         = plist.cj4;
81     const nbnxn_excl_t *excl    = plist.excl;
82     const int *atom_types       = atdat.atom_types;
83     int ntypes                  = atdat.ntypes;
84     const float4 *xq            = atdat.xq;
85     float3 *f                   = atdat.f;
86     const float3 *shift_vec     = atdat.shift_vec;
87     float rcoulomb_sq           = nbparam.rcoulomb_sq;
88 #ifdef VDW_CUTOFF_CHECK
89     float rvdw_sq               = nbparam.rvdw_sq;
90     float vdw_in_range;
91 #endif
92 #ifdef EL_RF
93     float two_k_rf              = nbparam.two_k_rf;
94 #endif
95 #ifdef EL_EWALD
96     float coulomb_tab_scale     = nbparam.coulomb_tab_scale;
97 #endif
98 #ifdef PRUNE_NBL
99     float rlist_sq              = nbparam.rlist_sq;
100 #endif
101
102 #ifdef CALC_ENERGIES
103     float lj_shift    = nbparam.sh_invrc6;
104 #ifdef EL_EWALD
105     float beta        = nbparam.ewald_beta;
106     float ewald_shift = nbparam.sh_ewald;
107 #else
108     float c_rf        = nbparam.c_rf;
109 #endif
110     float *e_lj       = atdat.e_lj;
111     float *e_el       = atdat.e_el;
112 #endif
113
114     /* thread/block/warp id-s */
115     unsigned int tidxi  = threadIdx.x;
116     unsigned int tidxj  = threadIdx.y;
117     unsigned int tidx   = threadIdx.y * blockDim.x + threadIdx.x;
118     unsigned int bidx   = blockIdx.x;
119     unsigned int widx   = tidx / WARP_SIZE; /* warp index */
120
121     int sci, ci, cj, ci_offset,
122         ai, aj,
123         cij4_start, cij4_end,
124         typei, typej,
125         i, jm, j4, wexcl_idx;
126     float qi, qj_f,
127           r2, inv_r, inv_r2, inv_r6,
128           c6, c12,
129           int_bit,
130 #ifdef CALC_ENERGIES
131           E_lj, E_el, E_lj_p,
132 #endif
133           F_invr;
134     unsigned int wexcl, imask, mask_ji;
135     float4 xqbuf;
136     float3 xi, xj, rv, f_ij, fcj_buf, fshift_buf;
137     float3 fci_buf[NCL_PER_SUPERCL];    /* i force buffer */
138     nbnxn_sci_t nb_sci;
139
140     /* shmem buffer for i x+q pre-loading */
141     extern __shared__  float4 xqib[];
142 #ifdef IATYPE_SHMEM
143     /* shmem buffer for i atom-type pre-loading */
144     int *atib = (int *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
145 #endif
146
147 #ifndef REDUCE_SHUFFLE
148     /* shmem j force buffer */
149 #ifdef IATYPE_SHMEM
150     float *f_buf = (float *)(atib + NCL_PER_SUPERCL * CL_SIZE);
151 #else
152     float *f_buf = (float *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
153 #endif
154 #endif
155
156     nb_sci      = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
157     sci         = nb_sci.sci;           /* super-cluster */
158     cij4_start  = nb_sci.cj4_ind_start; /* first ...*/
159     cij4_end    = nb_sci.cj4_ind_end;   /* and last index of j clusters */
160
161     /* Store the i-atom x and q in shared memory */
162     /* Note: the thread indexing here is inverted with respect to the
163        inner-loop as this results in slightly higher performance */
164     ci = sci * NCL_PER_SUPERCL + tidxi;
165     ai = ci * CL_SIZE + tidxj;
166     xqib[tidxi * CL_SIZE + tidxj] = xq[ai] + shift_vec[nb_sci.shift];
167 #ifdef IATYPE_SHMEM
168     ci = sci * NCL_PER_SUPERCL + tidxj;
169     ai = ci * CL_SIZE + tidxi;
170     atib[tidxj * CL_SIZE + tidxi] = atom_types[ai];
171 #endif
172     __syncthreads();
173
174     for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
175     {
176         fci_buf[ci_offset] = make_float3(0.0f);
177     }
178
179 #ifdef CALC_ENERGIES
180     E_lj = 0.0f;
181     E_el = 0.0f;
182
183 #if defined EL_EWALD || defined EL_RF
184     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
185     {
186         /* we have the diagonal: add the charge self interaction energy term */
187         for (i = 0; i < NCL_PER_SUPERCL; i++)
188         {
189             qi    = xqib[i * CL_SIZE + tidxi].w;
190             E_el += qi*qi;
191         }
192         /* divide the self term equally over the j-threads */
193         E_el /= CL_SIZE;
194 #ifdef EL_RF
195         E_el *= -nbparam.epsfac*0.5f*c_rf;
196 #else
197         E_el *= -nbparam.epsfac*beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
198 #endif
199     }
200 #endif
201 #endif
202
203     /* skip central shifts when summing shift forces */
204     if (nb_sci.shift == CENTRAL)
205     {
206         bCalcFshift = false;
207     }
208
209     fshift_buf = make_float3(0.0f);
210
211     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
212     for (j4 = cij4_start; j4 < cij4_end; j4++)
213     {
214         wexcl_idx   = pl_cj4[j4].imei[widx].excl_ind;
215         imask       = pl_cj4[j4].imei[widx].imask;
216         wexcl       = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
217
218 #ifndef PRUNE_NBL
219         if (imask)
220 #endif
221         {
222             /* Unrolling this loop
223                - with pruning leads to register spilling;
224                - on Kepler is much slower;
225                - doesn't work on CUDA <v4.1
226                Tested with nvcc 3.2 - 5.0.7 */
227 #if !defined PRUNE_NBL && __CUDA_ARCH__ < 300 && CUDA_VERSION >= 4010
228 #pragma unroll 4
229 #endif
230             for (jm = 0; jm < 4; jm++)
231             {
232                 if (imask & (255U << (jm * NCL_PER_SUPERCL)))
233                 {
234                     mask_ji = (1U << (jm * NCL_PER_SUPERCL));
235
236                     cj      = pl_cj4[j4].cj[jm];
237                     aj      = cj * CL_SIZE + tidxj;
238
239                     /* load j atom data */
240                     xqbuf   = xq[aj];
241                     xj      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
242                     qj_f    = nbparam.epsfac * xqbuf.w;
243                     typej   = atom_types[aj];
244
245                     fcj_buf = make_float3(0.0f);
246
247                     /* The PME and RF kernels don't unroll with CUDA <v4.1. */
248 #if !defined PRUNE_NBL && !(CUDA_VERSION < 4010 && (defined EL_EWALD || defined EL_RF))
249 #pragma unroll 8
250 #endif
251                     for(i = 0; i < NCL_PER_SUPERCL; i++)
252                     {
253                         if (imask & mask_ji)
254                         {
255                             ci_offset   = i;    /* i force buffer offset */
256
257                             ci      = sci * NCL_PER_SUPERCL + i; /* i cluster index */
258                             ai      = ci * CL_SIZE + tidxi;      /* i atom index */
259
260                             /* all threads load an atom from i cluster ci into shmem! */
261                             xqbuf   = xqib[i * CL_SIZE + tidxi];
262                             xi      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
263
264                             /* distance between i and j atoms */
265                             rv      = xi - xj;
266                             r2      = norm2(rv);
267
268 #ifdef PRUNE_NBL
269                             /* If _none_ of the atoms pairs are in cutoff range,
270                                the bit corresponding to the current
271                                cluster-pair in imask gets set to 0. */
272                             if (!__any(r2 < rlist_sq))
273                             {
274                                 imask &= ~mask_ji;
275                             }
276 #endif
277
278                             int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
279
280                             /* cutoff & exclusion check */
281 #if defined EL_EWALD || defined EL_RF
282                             if (r2 < rcoulomb_sq *
283                                 (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
284 #else
285                             if (r2 < rcoulomb_sq * int_bit)
286 #endif
287                             {
288                                 /* load the rest of the i-atom parameters */
289                                 qi      = xqbuf.w;
290 #ifdef IATYPE_SHMEM
291                                 typei   = atib[i * CL_SIZE + tidxi];
292 #else
293                                 typei   = atom_types[ai];
294 #endif
295
296                                 /* LJ 6*C6 and 12*C12 */
297                                 c6      = tex1Dfetch(tex_nbfp, 2 * (ntypes * typei + typej));
298                                 c12     = tex1Dfetch(tex_nbfp, 2 * (ntypes * typei + typej) + 1);
299
300                                 /* avoid NaN for excluded pairs at r=0 */
301                                 r2      += (1.0f - int_bit) * NBNXN_AVOID_SING_R2_INC;
302
303                                 inv_r   = rsqrt(r2);
304                                 inv_r2  = inv_r * inv_r;
305                                 inv_r6  = inv_r2 * inv_r2 * inv_r2;
306 #if defined EL_EWALD || defined EL_RF
307                                 /* We could mask inv_r2, but with Ewald
308                                  * masking both inv_r6 and F_invr is faster */
309                                 inv_r6  *= int_bit;
310 #endif
311
312                                 F_invr  = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
313
314 #ifdef CALC_ENERGIES
315                                 E_lj_p  = int_bit * (c12 * (inv_r6 * inv_r6 - lj_shift * lj_shift) * 0.08333333f - c6 * (inv_r6 - lj_shift) * 0.16666667f);
316 #endif
317
318 #ifdef VDW_CUTOFF_CHECK
319                                 /* this enables twin-range cut-offs (rvdw < rcoulomb <= rlist) */
320                                 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
321                                 F_invr  *= vdw_in_range;
322 #ifdef CALC_ENERGIES
323                                 E_lj_p  *= vdw_in_range;
324 #endif
325 #endif
326 #ifdef CALC_ENERGIES
327                                 E_lj    += E_lj_p;
328 #endif
329
330
331 #ifdef EL_CUTOFF
332                                 F_invr  += qi * qj_f * inv_r2 * inv_r;
333 #endif
334 #ifdef EL_RF
335                                 F_invr  += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
336 #endif
337 #ifdef EL_EWALD
338                                 F_invr  += qi * qj_f * (int_bit*inv_r2 - interpolate_coulomb_force_r(r2 * inv_r, coulomb_tab_scale)) * inv_r;
339 #endif
340
341 #ifdef CALC_ENERGIES
342 #ifdef EL_CUTOFF
343                                 E_el    += qi * qj_f * (inv_r - c_rf);
344 #endif
345 #ifdef EL_RF
346                                 E_el    += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
347 #endif
348 #ifdef EL_EWALD
349                                 /* 1.0f - erff is faster than erfcf */
350                                 E_el    += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
351 #endif
352 #endif
353                                 f_ij    = rv * F_invr;
354
355                                 /* accumulate j forces in registers */
356                                 fcj_buf -= f_ij;
357
358                                 /* accumulate i forces in registers */
359                                 fci_buf[ci_offset] += f_ij;
360                             }
361                         }
362
363                         /* shift the mask bit by 1 */
364                         mask_ji += mask_ji;
365                     }
366
367                     /* reduce j forces */
368 #ifdef REDUCE_SHUFFLE
369                     reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj);
370 #else
371                     /* store j forces in shmem */
372                     f_buf[                  tidx] = fcj_buf.x;
373                     f_buf[    FBUF_STRIDE + tidx] = fcj_buf.y;
374                     f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
375
376                     reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
377 #endif
378                 }
379             }
380 #ifdef PRUNE_NBL
381             /* Update the imask with the new one which does not contain the
382                out of range clusters anymore. */
383             pl_cj4[j4].imei[widx].imask = imask;
384 #endif
385         }
386     }
387
388     /* reduce i forces */
389     for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
390     {
391         ai  = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
392 #ifdef REDUCE_SHUFFLE
393         reduce_force_i_warp_shfl(fci_buf[ci_offset], f,
394                                  &fshift_buf, bCalcFshift,
395                                  tidxj, ai);
396 #else
397         f_buf[                  tidx] = fci_buf[ci_offset].x;
398         f_buf[    FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
399         f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
400         __syncthreads();
401         reduce_force_i(f_buf, f,
402                        &fshift_buf, bCalcFshift,
403                        tidxi, tidxj, ai);
404         __syncthreads();
405 #endif
406     }
407
408     /* add up local shift forces into global mem */
409 #ifdef REDUCE_SHUFFLE
410     if (bCalcFshift && (tidxj == 0 || tidxj == 4))
411 #else
412     if (bCalcFshift && tidxj == 0)
413 #endif
414     {
415         atomicAdd(&atdat.fshift[nb_sci.shift].x, fshift_buf.x);
416         atomicAdd(&atdat.fshift[nb_sci.shift].y, fshift_buf.y);
417         atomicAdd(&atdat.fshift[nb_sci.shift].z, fshift_buf.z);
418     }
419
420 #ifdef CALC_ENERGIES
421 #ifdef REDUCE_SHUFFLE
422     /* reduce the energies over warps and store into global memory */
423     reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx);
424 #else
425     /* flush the energies to shmem and reduce them */
426     f_buf[              tidx] = E_lj;
427     f_buf[FBUF_STRIDE + tidx] = E_el;
428     reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);
429 #endif
430 #endif
431 }