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     /* shmem buffer for cj, for both warps separately */
143     int *cjs     = (int *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
144 #ifdef IATYPE_SHMEM
145     /* shmem buffer for i atom-type pre-loading */
146     int *atib = (int *)(cjs + 2 * NBNXN_GPU_JGROUP_SIZE);
147 #endif
148
149 #ifndef REDUCE_SHUFFLE
150     /* shmem j force buffer */
151 #ifdef IATYPE_SHMEM
152     float *f_buf = (float *)(atib + NCL_PER_SUPERCL * CL_SIZE);
153 #else
154     float *f_buf = (float *)(cjs + 2 * NBNXN_GPU_JGROUP_SIZE);
155 #endif
156 #endif
157
158     nb_sci      = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
159     sci         = nb_sci.sci;           /* super-cluster */
160     cij4_start  = nb_sci.cj4_ind_start; /* first ...*/
161     cij4_end    = nb_sci.cj4_ind_end;   /* and last index of j clusters */
162
163     /* Store the i-atom x and q in shared memory */
164     /* Note: the thread indexing here is inverted with respect to the
165        inner-loop as this results in slightly higher performance */
166     ci = sci * NCL_PER_SUPERCL + tidxi;
167     ai = ci * CL_SIZE + tidxj;
168     xqib[tidxi * CL_SIZE + tidxj] = xq[ai] + shift_vec[nb_sci.shift];
169 #ifdef IATYPE_SHMEM
170     ci = sci * NCL_PER_SUPERCL + tidxj;
171     ai = ci * CL_SIZE + tidxi;
172     atib[tidxj * CL_SIZE + tidxi] = atom_types[ai];
173 #endif
174     __syncthreads();
175
176     for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
177     {
178         fci_buf[ci_offset] = make_float3(0.0f);
179     }
180
181 #ifdef CALC_ENERGIES
182     E_lj = 0.0f;
183     E_el = 0.0f;
184
185 #if defined EL_EWALD || defined EL_RF
186     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
187     {
188         /* we have the diagonal: add the charge self interaction energy term */
189         for (i = 0; i < NCL_PER_SUPERCL; i++)
190         {
191             qi    = xqib[i * CL_SIZE + tidxi].w;
192             E_el += qi*qi;
193         }
194         /* divide the self term equally over the j-threads */
195         E_el /= CL_SIZE;
196 #ifdef EL_RF
197         E_el *= -nbparam.epsfac*0.5f*c_rf;
198 #else
199         E_el *= -nbparam.epsfac*beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
200 #endif
201     }
202 #endif
203 #endif
204
205     /* skip central shifts when summing shift forces */
206     if (nb_sci.shift == CENTRAL)
207     {
208         bCalcFshift = false;
209     }
210
211     fshift_buf = make_float3(0.0f);
212
213     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
214     for (j4 = cij4_start; j4 < cij4_end; j4++)
215     {
216         wexcl_idx   = pl_cj4[j4].imei[widx].excl_ind;
217         imask       = pl_cj4[j4].imei[widx].imask;
218         wexcl       = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
219
220 #ifndef PRUNE_NBL
221         if (imask)
222 #endif
223         {
224             /* Pre-load cj into shared memory on both warps separately */
225             if ((tidxj == 0 || tidxj == 4) && tidxi < NBNXN_GPU_JGROUP_SIZE)
226             {
227                 cjs[tidxi + tidxj * NBNXN_GPU_JGROUP_SIZE / 4] = pl_cj4[j4].cj[tidxi];
228             }
229
230             /* Unrolling this loop
231                - with pruning leads to register spilling;
232                - on Kepler is much slower;
233                - doesn't work on CUDA <v4.1
234                Tested with nvcc 3.2 - 5.0.7 */
235 #if !defined PRUNE_NBL && __CUDA_ARCH__ < 300 && CUDA_VERSION >= 4010
236 #pragma unroll 4
237 #endif
238             for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
239             {
240                 if (imask & (supercl_interaction_mask << (jm * NCL_PER_SUPERCL)))
241                 {
242                     mask_ji = (1U << (jm * NCL_PER_SUPERCL));
243
244                     cj      = cjs[jm + (tidxj & 4) * NBNXN_GPU_JGROUP_SIZE / 4];
245                     aj      = cj * CL_SIZE + tidxj;
246
247                     /* load j atom data */
248                     xqbuf   = xq[aj];
249                     xj      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
250                     qj_f    = nbparam.epsfac * xqbuf.w;
251                     typej   = atom_types[aj];
252
253                     fcj_buf = make_float3(0.0f);
254
255                     /* The PME and RF kernels don't unroll with CUDA <v4.1. */
256 #if !defined PRUNE_NBL && !(CUDA_VERSION < 4010 && (defined EL_EWALD || defined EL_RF))
257 #pragma unroll 8
258 #endif
259                     for(i = 0; i < NCL_PER_SUPERCL; i++)
260                     {
261                         if (imask & mask_ji)
262                         {
263                             ci_offset   = i;    /* i force buffer offset */
264
265                             ci      = sci * NCL_PER_SUPERCL + i; /* i cluster index */
266                             ai      = ci * CL_SIZE + tidxi;      /* i atom index */
267
268                             /* all threads load an atom from i cluster ci into shmem! */
269                             xqbuf   = xqib[i * CL_SIZE + tidxi];
270                             xi      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
271
272                             /* distance between i and j atoms */
273                             rv      = xi - xj;
274                             r2      = norm2(rv);
275
276 #ifdef PRUNE_NBL
277                             /* If _none_ of the atoms pairs are in cutoff range,
278                                the bit corresponding to the current
279                                cluster-pair in imask gets set to 0. */
280                             if (!__any(r2 < rlist_sq))
281                             {
282                                 imask &= ~mask_ji;
283                             }
284 #endif
285
286                             int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
287
288                             /* cutoff & exclusion check */
289 #if defined EL_EWALD || defined EL_RF
290                             if (r2 < rcoulomb_sq *
291                                 (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
292 #else
293                             if (r2 < rcoulomb_sq * int_bit)
294 #endif
295                             {
296                                 /* load the rest of the i-atom parameters */
297                                 qi      = xqbuf.w;
298 #ifdef IATYPE_SHMEM
299                                 typei   = atib[i * CL_SIZE + tidxi];
300 #else
301                                 typei   = atom_types[ai];
302 #endif
303
304                                 /* LJ 6*C6 and 12*C12 */
305                                 c6      = tex1Dfetch(tex_nbfp, 2 * (ntypes * typei + typej));
306                                 c12     = tex1Dfetch(tex_nbfp, 2 * (ntypes * typei + typej) + 1);
307
308                                 /* avoid NaN for excluded pairs at r=0 */
309                                 r2      += (1.0f - int_bit) * NBNXN_AVOID_SING_R2_INC;
310
311                                 inv_r   = rsqrt(r2);
312                                 inv_r2  = inv_r * inv_r;
313                                 inv_r6  = inv_r2 * inv_r2 * inv_r2;
314 #if defined EL_EWALD || defined EL_RF
315                                 /* We could mask inv_r2, but with Ewald
316                                  * masking both inv_r6 and F_invr is faster */
317                                 inv_r6  *= int_bit;
318 #endif
319
320                                 F_invr  = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
321
322 #ifdef CALC_ENERGIES
323                                 E_lj_p  = int_bit * (c12 * (inv_r6 * inv_r6 - lj_shift * lj_shift) * 0.08333333f - c6 * (inv_r6 - lj_shift) * 0.16666667f);
324 #endif
325
326 #ifdef VDW_CUTOFF_CHECK
327                                 /* this enables twin-range cut-offs (rvdw < rcoulomb <= rlist) */
328                                 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
329                                 F_invr  *= vdw_in_range;
330 #ifdef CALC_ENERGIES
331                                 E_lj_p  *= vdw_in_range;
332 #endif
333 #endif
334 #ifdef CALC_ENERGIES
335                                 E_lj    += E_lj_p;
336 #endif
337
338
339 #ifdef EL_CUTOFF
340                                 F_invr  += qi * qj_f * inv_r2 * inv_r;
341 #endif
342 #ifdef EL_RF
343                                 F_invr  += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
344 #endif
345 #ifdef EL_EWALD
346                                 F_invr  += qi * qj_f * (int_bit*inv_r2 - interpolate_coulomb_force_r(r2 * inv_r, coulomb_tab_scale)) * inv_r;
347 #endif
348
349 #ifdef CALC_ENERGIES
350 #ifdef EL_CUTOFF
351                                 E_el    += qi * qj_f * (inv_r - c_rf);
352 #endif
353 #ifdef EL_RF
354                                 E_el    += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
355 #endif
356 #ifdef EL_EWALD
357                                 /* 1.0f - erff is faster than erfcf */
358                                 E_el    += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
359 #endif
360 #endif
361                                 f_ij    = rv * F_invr;
362
363                                 /* accumulate j forces in registers */
364                                 fcj_buf -= f_ij;
365
366                                 /* accumulate i forces in registers */
367                                 fci_buf[ci_offset] += f_ij;
368                             }
369                         }
370
371                         /* shift the mask bit by 1 */
372                         mask_ji += mask_ji;
373                     }
374
375                     /* reduce j forces */
376 #ifdef REDUCE_SHUFFLE
377                     reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj);
378 #else
379                     /* store j forces in shmem */
380                     f_buf[                  tidx] = fcj_buf.x;
381                     f_buf[    FBUF_STRIDE + tidx] = fcj_buf.y;
382                     f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
383
384                     reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
385 #endif
386                 }
387             }
388 #ifdef PRUNE_NBL
389             /* Update the imask with the new one which does not contain the
390                out of range clusters anymore. */
391             pl_cj4[j4].imei[widx].imask = imask;
392 #endif
393         }
394     }
395
396     /* reduce i forces */
397     for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
398     {
399         ai  = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
400 #ifdef REDUCE_SHUFFLE
401         reduce_force_i_warp_shfl(fci_buf[ci_offset], f,
402                                  &fshift_buf, bCalcFshift,
403                                  tidxj, ai);
404 #else
405         f_buf[                  tidx] = fci_buf[ci_offset].x;
406         f_buf[    FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
407         f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
408         __syncthreads();
409         reduce_force_i(f_buf, f,
410                        &fshift_buf, bCalcFshift,
411                        tidxi, tidxj, ai);
412         __syncthreads();
413 #endif
414     }
415
416     /* add up local shift forces into global mem */
417 #ifdef REDUCE_SHUFFLE
418     if (bCalcFshift && (tidxj == 0 || tidxj == 4))
419 #else
420     if (bCalcFshift && tidxj == 0)
421 #endif
422     {
423         atomicAdd(&atdat.fshift[nb_sci.shift].x, fshift_buf.x);
424         atomicAdd(&atdat.fshift[nb_sci.shift].y, fshift_buf.y);
425         atomicAdd(&atdat.fshift[nb_sci.shift].z, fshift_buf.z);
426     }
427
428 #ifdef CALC_ENERGIES
429 #ifdef REDUCE_SHUFFLE
430     /* reduce the energies over warps and store into global memory */
431     reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx);
432 #else
433     /* flush the energies to shmem and reduce them */
434     f_buf[              tidx] = E_lj;
435     f_buf[FBUF_STRIDE + tidx] = E_el;
436     reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);
437 #endif
438 #endif
439 }