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